1 2 #include <../src/snes/impls/richardson/snesrichardsonimpl.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESReset_NRichardson" 6 PetscErrorCode SNESReset_NRichardson(SNES snes) 7 { 8 PetscErrorCode ierr; 9 10 PetscFunctionBegin; 11 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 12 PetscFunctionReturn(0); 13 } 14 15 /* 16 SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson(). 17 18 Input Parameter: 19 . snes - the SNES context 20 21 Application Interface Routine: SNESDestroy() 22 */ 23 #undef __FUNCT__ 24 #define __FUNCT__ "SNESDestroy_NRichardson" 25 PetscErrorCode SNESDestroy_NRichardson(SNES snes) 26 { 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr); 31 ierr = PetscFree(snes->data);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 /* 36 SNESSetUp_NRichardson - Sets up the internal data structures for the later use 37 of the SNESNRICHARDSON nonlinear solver. 38 39 Input Parameters: 40 + snes - the SNES context 41 - x - the solution vector 42 43 Application Interface Routine: SNESSetUp() 44 */ 45 #undef __FUNCT__ 46 #define __FUNCT__ "SNESSetUp_NRichardson" 47 PetscErrorCode SNESSetUp_NRichardson(SNES snes) 48 { 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr); 53 PetscFunctionReturn(0); 54 } 55 56 /* 57 SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method. 58 59 Input Parameter: 60 . snes - the SNES context 61 62 Application Interface Routine: SNESSetFromOptions() 63 */ 64 #undef __FUNCT__ 65 #define __FUNCT__ "SNESSetFromOptions_NRichardson" 66 static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes) 67 { 68 PetscErrorCode ierr; 69 PetscFunctionBegin; 70 ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr); 71 ierr = PetscOptionsTail();CHKERRQ(ierr); 72 PetscFunctionReturn(0); 73 } 74 75 /* 76 SNESView_NRichardson - Prints info from the SNESRichardson data structure. 77 78 Input Parameters: 79 + SNES - the SNES context 80 - viewer - visualization context 81 82 Application Interface Routine: SNESView() 83 */ 84 #undef __FUNCT__ 85 #define __FUNCT__ "SNESView_NRichardson" 86 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 87 { 88 PetscBool iascii; 89 PetscErrorCode ierr; 90 91 PetscFunctionBegin; 92 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 93 if (iascii) { 94 ierr = PetscViewerASCIIPrintf(viewer," richardson variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr); 95 } 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "SNESLineSearchNo_NRichardson" 101 PetscErrorCode SNESLineSearchNo_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 102 { 103 PetscErrorCode ierr; 104 105 PetscFunctionBegin; 106 ierr = VecAXPY(X, snes->damping, Y);CHKERRQ(ierr); 107 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 108 ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr); 109 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm"); 110 PetscFunctionReturn(0); 111 } 112 113 #undef __FUNCT__ 114 #define __FUNCT__ "SNESLineSearchNoNorms_NRichardson" 115 PetscErrorCode SNESLineSearchNoNorms_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 116 { 117 PetscErrorCode ierr; 118 119 PetscFunctionBegin; 120 ierr = VecAXPY(X, snes->damping, Y);CHKERRQ(ierr); 121 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "SNESLineSearchQuadratic_NRichardson" 127 PetscErrorCode SNESLineSearchQuadratic_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 128 { 129 PetscInt i; 130 PetscReal alphas[3] = {0.0, 0.5, 1.0}; 131 PetscReal norms[3]; 132 PetscReal alpha,a,b; 133 PetscErrorCode ierr; 134 135 PetscFunctionBegin; 136 norms[0] = fnorm; 137 for(i=1; i < 3; ++i) { 138 ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 139 ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 140 ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr); 141 } 142 for(i = 0; i < 3; ++i) { 143 norms[i] = PetscSqr(norms[i]); 144 } 145 /* Fit a quadratic: 146 If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 147 a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1) 148 b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1) 149 c = y_0 150 x_min = -b/2a 151 152 If we let x_{0,1,2} = 0, 0.5, 1.0 153 a = 2 y_2 - 4 y_1 + 2 y_0 154 b = -y_2 + 4 y_1 - 3 y_0 155 c = y_0 156 */ 157 a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1])); 158 b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]); 159 /* Check for positive a (concave up) */ 160 if (a >= 0.0) { 161 alpha = -b/(2.0*a); 162 alpha = PetscMin(alpha, alphas[2]); 163 alpha = PetscMax(alpha, alphas[0]); 164 } else { 165 alpha = 1.0; 166 } 167 if (snes->ls_monitor) { 168 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 169 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr); 170 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 171 } 172 ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr); 173 if (alpha != 1.0) { 174 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 175 ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr); 176 } else { 177 *gnorm = PetscSqrtReal(norms[2]); 178 } 179 if (alpha == 0.0) *flag = PETSC_FALSE; 180 else *flag = PETSC_TRUE; 181 PetscFunctionReturn(0); 182 } 183 184 /* 185 SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 186 187 Input Parameters: 188 . snes - the SNES context 189 190 Output Parameter: 191 . outits - number of iterations until termination 192 193 Application Interface Routine: SNESSolve() 194 */ 195 #undef __FUNCT__ 196 #define __FUNCT__ "SNESSolve_NRichardson" 197 PetscErrorCode SNESSolve_NRichardson(SNES snes) 198 { 199 Vec X, Y, F, W; 200 PetscReal fnorm; 201 PetscInt maxits, i; 202 PetscErrorCode ierr; 203 SNESConvergedReason reason; 204 205 PetscFunctionBegin; 206 snes->reason = SNES_CONVERGED_ITERATING; 207 208 maxits = snes->max_its; /* maximum number of iterations */ 209 X = snes->vec_sol; /* X^n */ 210 Y = snes->vec_sol_update; /* \tilde X */ 211 F = snes->vec_func; /* residual vector */ 212 W = snes->work[0]; /* work vector */ 213 214 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 215 snes->iter = 0; 216 snes->norm = 0.; 217 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 218 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 219 if (snes->domainerror) { 220 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 221 PetscFunctionReturn(0); 222 } 223 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 224 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 225 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 226 snes->norm = fnorm; 227 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 228 SNESLogConvHistory(snes,fnorm,0); 229 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 230 231 /* set parameter for default relative tolerance convergence test */ 232 snes->ttol = fnorm*snes->rtol; 233 /* test convergence */ 234 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 235 if (snes->reason) PetscFunctionReturn(0); 236 237 for(i = 0; i < maxits; i++) { 238 PetscBool lsSuccess = PETSC_TRUE; 239 PetscReal dummyNorm; 240 241 /* Call general purpose update function */ 242 if (snes->ops->update) { 243 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 244 } 245 if (!snes->pc) { 246 ierr = VecCopy(F,Y);CHKERRQ(ierr); 247 ierr = VecScale(Y,-1.0);CHKERRQ(ierr); 248 } else { 249 ierr = VecCopy(X,Y);CHKERRQ(ierr); 250 ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr); 251 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 252 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 253 snes->reason = SNES_DIVERGED_INNER; 254 PetscFunctionReturn(0); 255 } 256 ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 257 } 258 259 ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, Y, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr); 260 if (!lsSuccess) { 261 if (++snes->numFailures >= snes->maxFailures) { 262 snes->reason = SNES_DIVERGED_LINE_SEARCH; 263 break; 264 } 265 } 266 if (snes->nfuncs >= snes->max_funcs) { 267 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 268 break; 269 } 270 if (snes->domainerror) { 271 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 272 PetscFunctionReturn(0); 273 } 274 /* Monitor convergence */ 275 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 276 snes->iter = i+1; 277 snes->norm = fnorm; 278 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 279 SNESLogConvHistory(snes,snes->norm,0); 280 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 281 /* Test for convergence */ 282 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 283 if (snes->reason) break; 284 } 285 if (i == maxits) { 286 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 287 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 288 } 289 PetscFunctionReturn(0); 290 } 291 292 293 EXTERN_C_BEGIN 294 #undef __FUNCT__ 295 #define __FUNCT__ "SNESLineSearchSetType_NRichardson" 296 PetscErrorCode SNESLineSearchSetType_NRichardson(SNES snes, SNESLineSearchType type) 297 { 298 PetscErrorCode ierr; 299 PetscFunctionBegin; 300 301 switch (type) { 302 case SNES_LS_BASIC: 303 ierr = SNESLineSearchSet(snes,SNESLineSearchNo_NRichardson,PETSC_NULL);CHKERRQ(ierr); 304 break; 305 case SNES_LS_BASIC_NONORMS: 306 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_NRichardson,PETSC_NULL);CHKERRQ(ierr); 307 break; 308 case SNES_LS_QUADRATIC: 309 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NRichardson,PETSC_NULL);CHKERRQ(ierr); 310 break; 311 default: 312 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 313 break; 314 } 315 snes->ls_type = type; 316 PetscFunctionReturn(0); 317 } 318 EXTERN_C_END 319 320 /*MC 321 SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 322 323 Level: beginner 324 325 Options Database: 326 + -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms) 327 - -snes_ls <basic,basicnormnorms,quadratic> 328 329 Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda (F(x^n) - b) 330 where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. 331 If an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner solver is called an initial solution x^n and the 332 nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 333 334 This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 335 336 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN 337 M*/ 338 EXTERN_C_BEGIN 339 #undef __FUNCT__ 340 #define __FUNCT__ "SNESCreate_NRichardson" 341 PetscErrorCode SNESCreate_NRichardson(SNES snes) 342 { 343 PetscErrorCode ierr; 344 PetscFunctionBegin; 345 snes->ops->destroy = SNESDestroy_NRichardson; 346 snes->ops->setup = SNESSetUp_NRichardson; 347 snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 348 snes->ops->view = SNESView_NRichardson; 349 snes->ops->solve = SNESSolve_NRichardson; 350 snes->ops->reset = SNESReset_NRichardson; 351 352 snes->usesksp = PETSC_FALSE; 353 snes->usespc = PETSC_TRUE; 354 355 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NRichardson",SNESLineSearchSetType_NRichardson);CHKERRQ(ierr); 356 ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 357 358 PetscFunctionReturn(0); 359 } 360 EXTERN_C_END 361