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(SNES snes) 63 { 64 PetscErrorCode ierr; 65 SNESLineSearch linesearch; 66 67 PetscFunctionBegin; 68 ierr = PetscOptionsHead("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 PetscBool lsSuccess; 120 SNESConvergedReason reason; 121 122 PetscFunctionBegin; 123 snes->reason = SNES_CONVERGED_ITERATING; 124 125 maxits = snes->max_its; /* maximum number of iterations */ 126 X = snes->vec_sol; /* X^n */ 127 Y = snes->vec_sol_update; /* \tilde X */ 128 F = snes->vec_func; /* residual vector */ 129 130 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 131 snes->iter = 0; 132 snes->norm = 0.; 133 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 134 135 if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 136 ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 137 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 138 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 139 snes->reason = SNES_DIVERGED_INNER; 140 PetscFunctionReturn(0); 141 } 142 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 143 } else { 144 if (!snes->vec_func_init_set) { 145 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 146 if (snes->domainerror) { 147 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 148 PetscFunctionReturn(0); 149 } 150 } else snes->vec_func_init_set = PETSC_FALSE; 151 if (!snes->norm_init_set) { 152 /* convergence test */ 153 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 154 if (PetscIsInfOrNanReal(fnorm)) { 155 snes->reason = SNES_DIVERGED_FNORM_NAN; 156 PetscFunctionReturn(0); 157 } 158 } else { 159 fnorm = snes->norm_init; 160 snes->norm_init_set = PETSC_FALSE; 161 } 162 } 163 if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 164 ierr = SNESApplyPC(snes,X,F,&fnorm,Y);CHKERRQ(ierr); 165 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 166 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 167 snes->reason = SNES_DIVERGED_INNER; 168 PetscFunctionReturn(0); 169 } 170 } else { 171 ierr = VecCopy(F,Y);CHKERRQ(ierr); 172 } 173 174 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 175 snes->norm = fnorm; 176 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 177 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 178 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 179 180 /* set parameter for default relative tolerance convergence test */ 181 snes->ttol = fnorm*snes->rtol; 182 /* test convergence */ 183 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 184 if (snes->reason) PetscFunctionReturn(0); 185 186 /* Call general purpose update function */ 187 if (snes->ops->update) { 188 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 189 } 190 191 /* set parameter for default relative tolerance convergence test */ 192 snes->ttol = fnorm*snes->rtol; 193 /* test convergence */ 194 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 195 if (snes->reason) PetscFunctionReturn(0); 196 197 for (i = 1; i < maxits+1; i++) { 198 lsSuccess = PETSC_TRUE; 199 200 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 201 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 202 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr); 203 if (!lsSuccess) { 204 if (++snes->numFailures >= snes->maxFailures) { 205 snes->reason = SNES_DIVERGED_LINE_SEARCH; 206 break; 207 } 208 } 209 if (snes->nfuncs >= snes->max_funcs) { 210 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 211 break; 212 } 213 if (snes->domainerror) { 214 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 215 PetscFunctionReturn(0); 216 } 217 218 /* Monitor convergence */ 219 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 220 snes->iter = i; 221 snes->norm = fnorm; 222 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 223 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 224 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 225 /* Test for convergence */ 226 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 227 if (snes->reason) break; 228 229 /* Call general purpose update function */ 230 if (snes->ops->update) { 231 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 232 } 233 234 if (snes->pc) { 235 if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 236 ierr = SNESApplyPC(snes,X,NULL,NULL,Y);CHKERRQ(ierr); 237 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 238 ierr = VecCopy(Y,F);CHKERRQ(ierr); 239 } else { 240 ierr = SNESApplyPC(snes,X,F,&fnorm,Y);CHKERRQ(ierr); 241 } 242 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 243 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 244 snes->reason = SNES_DIVERGED_INNER; 245 PetscFunctionReturn(0); 246 } 247 } else { 248 ierr = VecCopy(F,Y);CHKERRQ(ierr); 249 } 250 } 251 if (i == maxits+1) { 252 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 253 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 254 } 255 PetscFunctionReturn(0); 256 } 257 258 /*MC 259 SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 260 261 Level: beginner 262 263 Options Database: 264 + -snes_linesearch_type <l2,cp,basic> Line search type. 265 - -snes_linesearch_damping<1.0> Damping for the line search. 266 267 Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 268 (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 269 an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 270 solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 271 where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 272 273 The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 274 linesearch, one may have to scale the update with -snes_linesearch_damping 275 276 This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 277 278 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG 279 M*/ 280 #undef __FUNCT__ 281 #define __FUNCT__ "SNESCreate_NRichardson" 282 PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) 283 { 284 PetscErrorCode ierr; 285 SNES_NRichardson *neP; 286 287 PetscFunctionBegin; 288 snes->ops->destroy = SNESDestroy_NRichardson; 289 snes->ops->setup = SNESSetUp_NRichardson; 290 snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 291 snes->ops->view = SNESView_NRichardson; 292 snes->ops->solve = SNESSolve_NRichardson; 293 snes->ops->reset = SNESReset_NRichardson; 294 295 snes->usesksp = PETSC_FALSE; 296 snes->usespc = PETSC_TRUE; 297 298 snes->pcside = PC_LEFT; 299 300 ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr); 301 snes->data = (void*) neP; 302 303 if (!snes->tolerancesset) { 304 snes->max_funcs = 30000; 305 snes->max_its = 10000; 306 snes->stol = 1e-20; 307 } 308 PetscFunctionReturn(0); 309 } 310