1 #include <../src/snes/impls/richardson/snesrichardsonimpl.h> 2 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESReset_NRichardson" 6 PetscErrorCode SNESReset_NRichardson(SNES snes) 7 { 8 PetscFunctionBegin; 9 PetscFunctionReturn(0); 10 } 11 12 /* 13 SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson(). 14 15 Input Parameter: 16 . snes - the SNES context 17 18 Application Interface Routine: SNESDestroy() 19 */ 20 #undef __FUNCT__ 21 #define __FUNCT__ "SNESDestroy_NRichardson" 22 PetscErrorCode SNESDestroy_NRichardson(SNES snes) 23 { 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr); 28 ierr = PetscFree(snes->data);CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 /* 33 SNESSetUp_NRichardson - Sets up the internal data structures for the later use 34 of the SNESNRICHARDSON nonlinear solver. 35 36 Input Parameters: 37 + snes - the SNES context 38 - x - the solution vector 39 40 Application Interface Routine: SNESSetUp() 41 */ 42 #undef __FUNCT__ 43 #define __FUNCT__ "SNESSetUp_NRichardson" 44 PetscErrorCode SNESSetUp_NRichardson(SNES snes) 45 { 46 PetscFunctionBegin; 47 if (snes->pcside == PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");} 48 if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 49 PetscFunctionReturn(0); 50 } 51 52 /* 53 SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method. 54 55 Input Parameter: 56 . snes - the SNES context 57 58 Application Interface Routine: SNESSetFromOptions() 59 */ 60 #undef __FUNCT__ 61 #define __FUNCT__ "SNESSetFromOptions_NRichardson" 62 static PetscErrorCode SNESSetFromOptions_NRichardson(PetscOptions *PetscOptionsObject,SNES snes) 63 { 64 PetscErrorCode ierr; 65 SNESLineSearch linesearch; 66 67 PetscFunctionBegin; 68 ierr = PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");CHKERRQ(ierr); 69 ierr = PetscOptionsTail();CHKERRQ(ierr); 70 if (!snes->linesearch) { 71 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 72 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 73 } 74 PetscFunctionReturn(0); 75 } 76 77 /* 78 SNESView_NRichardson - Prints info from the SNESRichardson data structure. 79 80 Input Parameters: 81 + SNES - the SNES context 82 - viewer - visualization context 83 84 Application Interface Routine: SNESView() 85 */ 86 #undef __FUNCT__ 87 #define __FUNCT__ "SNESView_NRichardson" 88 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 89 { 90 PetscBool iascii; 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 95 if (iascii) { 96 } 97 PetscFunctionReturn(0); 98 } 99 100 /* 101 SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 102 103 Input Parameters: 104 . snes - the SNES context 105 106 Output Parameter: 107 . outits - number of iterations until termination 108 109 Application Interface Routine: SNESSolve() 110 */ 111 #undef __FUNCT__ 112 #define __FUNCT__ "SNESSolve_NRichardson" 113 PetscErrorCode SNESSolve_NRichardson(SNES snes) 114 { 115 Vec X, Y, F; 116 PetscReal xnorm, fnorm, ynorm; 117 PetscInt maxits, i; 118 PetscErrorCode ierr; 119 SNESLineSearchReason lsresult; 120 SNESConvergedReason reason; 121 122 PetscFunctionBegin; 123 124 if (snes->xl || snes->xu || snes->ops->computevariablebounds) { 125 SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 126 } 127 128 snes->reason = SNES_CONVERGED_ITERATING; 129 130 maxits = snes->max_its; /* maximum number of iterations */ 131 X = snes->vec_sol; /* X^n */ 132 Y = snes->vec_sol_update; /* \tilde X */ 133 F = snes->vec_func; /* residual vector */ 134 135 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 136 snes->iter = 0; 137 snes->norm = 0.; 138 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 139 140 if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 141 ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 142 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 143 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 144 snes->reason = SNES_DIVERGED_INNER; 145 PetscFunctionReturn(0); 146 } 147 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 148 } else { 149 if (!snes->vec_func_init_set) { 150 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 151 } else snes->vec_func_init_set = PETSC_FALSE; 152 153 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 154 SNESCheckFunctionNorm(snes,fnorm); 155 } 156 if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 157 ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr); 158 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 159 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 160 snes->reason = SNES_DIVERGED_INNER; 161 PetscFunctionReturn(0); 162 } 163 } else { 164 ierr = VecCopy(F,Y);CHKERRQ(ierr); 165 } 166 167 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 168 snes->norm = fnorm; 169 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 170 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 171 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 172 173 /* test convergence */ 174 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 175 if (snes->reason) PetscFunctionReturn(0); 176 177 /* Call general purpose update function */ 178 if (snes->ops->update) { 179 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 180 } 181 182 /* set parameter for default relative tolerance convergence test */ 183 snes->ttol = fnorm*snes->rtol; 184 /* test convergence */ 185 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 186 if (snes->reason) PetscFunctionReturn(0); 187 188 for (i = 1; i < maxits+1; i++) { 189 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 190 ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr); 191 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 192 if (lsresult) { 193 if (++snes->numFailures >= snes->maxFailures) { 194 snes->reason = SNES_DIVERGED_LINE_SEARCH; 195 break; 196 } 197 } 198 if (snes->nfuncs >= snes->max_funcs) { 199 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 200 break; 201 } 202 203 /* Monitor convergence */ 204 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 205 snes->iter = i; 206 snes->norm = fnorm; 207 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 208 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 209 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 210 /* Test for convergence */ 211 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 212 if (snes->reason) break; 213 214 /* Call general purpose update function */ 215 if (snes->ops->update) { 216 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 217 } 218 219 if (snes->pc) { 220 if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 221 ierr = SNESApplyNPC(snes,X,NULL,Y);CHKERRQ(ierr); 222 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 223 ierr = VecCopy(Y,F);CHKERRQ(ierr); 224 } else { 225 ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr); 226 } 227 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 228 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 229 snes->reason = SNES_DIVERGED_INNER; 230 PetscFunctionReturn(0); 231 } 232 } else { 233 ierr = VecCopy(F,Y);CHKERRQ(ierr); 234 } 235 } 236 if (i == maxits+1) { 237 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 238 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 239 } 240 PetscFunctionReturn(0); 241 } 242 243 /*MC 244 SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 245 246 Level: beginner 247 248 Options Database: 249 + -snes_linesearch_type <l2,cp,basic> Line search type. 250 - -snes_linesearch_damping<1.0> Damping for the line search. 251 252 Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 253 (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 254 an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner 255 solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 256 where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 257 258 The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 259 linesearch, one may have to scale the update with -snes_linesearch_damping 260 261 This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 262 263 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG 264 M*/ 265 #undef __FUNCT__ 266 #define __FUNCT__ "SNESCreate_NRichardson" 267 PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) 268 { 269 PetscErrorCode ierr; 270 SNES_NRichardson *neP; 271 272 PetscFunctionBegin; 273 snes->ops->destroy = SNESDestroy_NRichardson; 274 snes->ops->setup = SNESSetUp_NRichardson; 275 snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 276 snes->ops->view = SNESView_NRichardson; 277 snes->ops->solve = SNESSolve_NRichardson; 278 snes->ops->reset = SNESReset_NRichardson; 279 280 snes->usesksp = PETSC_FALSE; 281 snes->usespc = PETSC_TRUE; 282 283 snes->pcside = PC_LEFT; 284 285 ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 286 snes->data = (void*) neP; 287 288 if (!snes->tolerancesset) { 289 snes->max_funcs = 30000; 290 snes->max_its = 10000; 291 snes->stol = 1e-20; 292 } 293 PetscFunctionReturn(0); 294 } 295