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