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