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,2);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__ "SNESLineSearchQuadratic_NRichardson" 101 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) 102 { 103 PetscInt i; 104 PetscReal alphas[3] = {0.0, 0.5, 1.0}; 105 PetscReal norms[3]; 106 PetscReal alpha,a,b; 107 PetscErrorCode ierr; 108 109 PetscFunctionBegin; 110 norms[0] = fnorm; 111 for(i=1; i < 3; ++i) { 112 ierr = VecWAXPY(W, -alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 113 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 114 ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr); 115 } 116 for(i = 0; i < 3; ++i) { 117 norms[i] = PetscSqr(norms[i]); 118 } 119 /* Fit a quadratic: 120 If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 121 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) 122 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) 123 c = y_0 124 x_min = -b/2a 125 126 If we let x_{0,1,2} = 0, 0.5, 1.0 127 a = 2 y_2 - 4 y_1 + 2 y_0 128 b = -y_2 + 4 y_1 - 3 y_0 129 c = y_0 130 */ 131 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])); 132 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]); 133 /* Check for positive a (concave up) */ 134 if (a >= 0.0) { 135 alpha = -b/(2.0*a); 136 alpha = PetscMin(alpha, alphas[2]); 137 alpha = PetscMax(alpha, alphas[0]); 138 } else { 139 alpha = 1.0; 140 } 141 if (snes->ls_monitor) { 142 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 143 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", 144 PetscSqrtReal(norms[0]),PetscSqrtReal(norms[1]),PetscSqrtReal(norms[2]),alpha);CHKERRQ(ierr); 145 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 146 } 147 ierr = VecCopy(X, W);CHKERRQ(ierr); 148 ierr = VecAXPY(W, -alpha, Y);CHKERRQ(ierr); 149 if (alpha != 1.0) { 150 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 151 ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 152 } else { 153 *gnorm = PetscSqrtReal(norms[2]); 154 } 155 if (alpha == 0.0) *flag = PETSC_FALSE; 156 else *flag = PETSC_TRUE; 157 PetscFunctionReturn(0); 158 } 159 160 /* 161 SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 162 163 Input Parameters: 164 . snes - the SNES context 165 166 Output Parameter: 167 . outits - number of iterations until termination 168 169 Application Interface Routine: SNESSolve() 170 */ 171 #undef __FUNCT__ 172 #define __FUNCT__ "SNESSolve_NRichardson" 173 PetscErrorCode SNESSolve_NRichardson(SNES snes) 174 { 175 Vec X, Y, F, W, G; 176 PetscReal fnorm; 177 PetscInt maxits, i; 178 PetscErrorCode ierr; 179 SNESConvergedReason reason; 180 181 PetscFunctionBegin; 182 snes->reason = SNES_CONVERGED_ITERATING; 183 184 maxits = snes->max_its; /* maximum number of iterations */ 185 X = snes->vec_sol; /* X^n */ 186 Y = snes->vec_sol_update; /* \tilde X */ 187 F = snes->vec_func; /* residual vector */ 188 W = snes->work[0]; /* work vector */ 189 G = snes->work[1]; /* line search function vector */ 190 191 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 192 snes->iter = 0; 193 snes->norm = 0.; 194 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 195 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 196 if (snes->domainerror) { 197 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 198 PetscFunctionReturn(0); 199 } 200 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 201 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 202 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 203 snes->norm = fnorm; 204 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 205 SNESLogConvHistory(snes,fnorm,0); 206 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 207 208 /* set parameter for default relative tolerance convergence test */ 209 snes->ttol = fnorm*snes->rtol; 210 /* test convergence */ 211 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 212 if (snes->reason) PetscFunctionReturn(0); 213 214 for(i = 0; i < maxits; i++) { 215 PetscBool lsSuccess = PETSC_TRUE; 216 PetscReal dummyNorm; 217 218 /* Call general purpose update function */ 219 if (snes->ops->update) { 220 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 221 } 222 else if (snes->pc) { 223 ierr = VecCopy(X,Y);CHKERRQ(ierr); 224 ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr); 225 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 226 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 227 snes->reason = SNES_DIVERGED_INNER; 228 PetscFunctionReturn(0); 229 } 230 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 231 } else { 232 ierr = VecCopy(F,Y);CHKERRQ(ierr); 233 } 234 ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, Y, fnorm, 0.0, G, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr); 235 if (!lsSuccess) { 236 if (++snes->numFailures >= snes->maxFailures) { 237 snes->reason = SNES_DIVERGED_LINE_SEARCH; 238 break; 239 } 240 } 241 if (snes->nfuncs >= snes->max_funcs) { 242 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 243 break; 244 } 245 if (snes->domainerror) { 246 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 247 PetscFunctionReturn(0); 248 } 249 ierr = VecCopy(G, F);CHKERRQ(ierr); 250 ierr = VecCopy(W, X);CHKERRQ(ierr); 251 /* Monitor convergence */ 252 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 253 snes->iter = i+1; 254 snes->norm = fnorm; 255 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 256 SNESLogConvHistory(snes,snes->norm,0); 257 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 258 /* Test for convergence */ 259 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 260 if (snes->reason) break; 261 } 262 if (i == maxits) { 263 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 264 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 265 } 266 PetscFunctionReturn(0); 267 } 268 269 270 EXTERN_C_BEGIN 271 #undef __FUNCT__ 272 #define __FUNCT__ "SNESLineSearchSetType_NRichardson" 273 PetscErrorCode SNESLineSearchSetType_NRichardson(SNES snes, SNESLineSearchType type) 274 { 275 PetscErrorCode ierr; 276 PetscFunctionBegin; 277 278 switch (type) { 279 case SNES_LS_BASIC: 280 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 281 break; 282 case SNES_LS_BASIC_NONORMS: 283 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 284 break; 285 case SNES_LS_QUADRATIC: 286 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NRichardson,PETSC_NULL);CHKERRQ(ierr); 287 break; 288 case SNES_LS_SECANT: 289 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 290 break; 291 default: 292 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 293 break; 294 } 295 snes->ls_type = type; 296 PetscFunctionReturn(0); 297 } 298 EXTERN_C_END 299 300 /*MC 301 SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 302 303 Level: beginner 304 305 Options Database: 306 + -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms) 307 - -snes_ls <basic,basicnormnorms,quadratic> 308 309 Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 310 (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 311 an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 312 solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 313 where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 314 315 This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 316 317 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN 318 M*/ 319 EXTERN_C_BEGIN 320 #undef __FUNCT__ 321 #define __FUNCT__ "SNESCreate_NRichardson" 322 PetscErrorCode SNESCreate_NRichardson(SNES snes) 323 { 324 PetscErrorCode ierr; 325 PetscFunctionBegin; 326 snes->ops->destroy = SNESDestroy_NRichardson; 327 snes->ops->setup = SNESSetUp_NRichardson; 328 snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 329 snes->ops->view = SNESView_NRichardson; 330 snes->ops->solve = SNESSolve_NRichardson; 331 snes->ops->reset = SNESReset_NRichardson; 332 333 snes->usesksp = PETSC_FALSE; 334 snes->usespc = PETSC_TRUE; 335 336 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NRichardson",SNESLineSearchSetType_NRichardson);CHKERRQ(ierr); 337 ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr); 338 339 PetscFunctionReturn(0); 340 } 341 EXTERN_C_END 342