1 2 #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 3 4 typedef struct { 5 void *ctx; 6 SNES snes; 7 } SNES_TR_KSPConverged_Ctx; 8 9 /* 10 This convergence test determines if the two norm of the 11 solution lies outside the trust region, if so it halts. 12 */ 13 #undef __FUNCT__ 14 #define __FUNCT__ "SNES_TR_KSPConverged_Private" 15 PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 16 { 17 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 18 SNES snes = ctx->snes; 19 SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 20 Vec x; 21 PetscReal nrm; 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 ierr = KSPConvergedDefault(ksp,n,rnorm,reason,ctx->ctx);CHKERRQ(ierr); 26 if (*reason) { 27 ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 28 } 29 /* Determine norm of solution */ 30 ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 31 ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 32 if (nrm >= neP->delta) { 33 ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 34 *reason = KSP_CONVERGED_STEP_LENGTH; 35 } 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "SNES_TR_KSPConverged_Destroy" 41 PetscErrorCode SNES_TR_KSPConverged_Destroy(void *cctx) 42 { 43 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 44 PetscErrorCode ierr; 45 46 PetscFunctionBegin; 47 ierr = KSPConvergedDefaultDestroy(ctx->ctx);CHKERRQ(ierr); 48 ierr = PetscFree(ctx);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 /* ---------------------------------------------------------------- */ 53 #undef __FUNCT__ 54 #define __FUNCT__ "SNES_TR_Converged_Private" 55 /* 56 SNES_TR_Converged_Private -test convergence JUST for 57 the trust region tolerance. 58 59 */ 60 static PetscErrorCode SNES_TR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 61 { 62 SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 *reason = SNES_CONVERGED_ITERATING; 67 if (neP->delta < xnorm * snes->deltatol) { 68 ierr = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 69 *reason = SNES_CONVERGED_TR_DELTA; 70 } else if (snes->nfuncs >= snes->max_funcs) { 71 ierr = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 72 *reason = SNES_DIVERGED_FUNCTION_COUNT; 73 } 74 PetscFunctionReturn(0); 75 } 76 77 78 /* 79 SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 80 region approach for solving systems of nonlinear equations. 81 82 83 */ 84 #undef __FUNCT__ 85 #define __FUNCT__ "SNESSolve_NEWTONTR" 86 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 87 { 88 SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 89 Vec X,F,Y,G,Ytmp; 90 PetscErrorCode ierr; 91 PetscInt maxits,i,lits; 92 PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 93 PetscScalar cnorm; 94 KSP ksp; 95 SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 96 PetscBool conv = PETSC_FALSE,breakout = PETSC_FALSE; 97 PetscBool domainerror; 98 99 PetscFunctionBegin; 100 maxits = snes->max_its; /* maximum number of iterations */ 101 X = snes->vec_sol; /* solution vector */ 102 F = snes->vec_func; /* residual vector */ 103 Y = snes->work[0]; /* work vectors */ 104 G = snes->work[1]; 105 Ytmp = snes->work[2]; 106 107 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 108 snes->iter = 0; 109 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 110 111 if (!snes->vec_func_init_set) { 112 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 113 ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 114 if (domainerror) { 115 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 116 PetscFunctionReturn(0); 117 } 118 } else snes->vec_func_init_set = PETSC_FALSE; 119 120 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 121 if (PetscIsInfOrNanReal(fnorm)) { 122 snes->reason = SNES_DIVERGED_FNORM_NAN; 123 PetscFunctionReturn(0); 124 } 125 126 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 127 snes->norm = fnorm; 128 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 129 delta = neP->delta0*fnorm; 130 neP->delta = delta; 131 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 132 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 133 134 /* test convergence */ 135 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 136 if (snes->reason) PetscFunctionReturn(0); 137 138 /* Set the stopping criteria to use the More' trick. */ 139 ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr); 140 if (!conv) { 141 SNES_TR_KSPConverged_Ctx *ctx; 142 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 143 ierr = PetscNew(&ctx);CHKERRQ(ierr); 144 ctx->snes = snes; 145 ierr = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr); 146 ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr); 147 ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 148 } 149 150 for (i=0; i<maxits; i++) { 151 152 /* Call general purpose update function */ 153 if (snes->ops->update) { 154 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 155 } 156 157 /* Solve J Y = F, where J is Jacobian matrix */ 158 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 159 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 160 ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 161 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 162 163 snes->linear_its += lits; 164 165 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 166 ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 167 norm1 = nrm; 168 while (1) { 169 ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 170 nrm = norm1; 171 172 /* Scale Y if need be and predict new value of F norm */ 173 if (nrm >= delta) { 174 nrm = delta/nrm; 175 gpnorm = (1.0 - nrm)*fnorm; 176 cnorm = nrm; 177 ierr = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 178 ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 179 nrm = gpnorm; 180 ynorm = delta; 181 } else { 182 gpnorm = 0.0; 183 ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 184 ynorm = nrm; 185 } 186 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 187 ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr); 188 ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 189 ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 190 if (fnorm == gpnorm) rho = 0.0; 191 else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 192 193 /* Update size of trust region */ 194 if (rho < neP->mu) delta *= neP->delta1; 195 else if (rho < neP->eta) delta *= neP->delta2; 196 else delta *= neP->delta3; 197 ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 198 ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 199 200 neP->delta = delta; 201 if (rho > neP->sigma) break; 202 ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 203 /* check to see if progress is hopeless */ 204 neP->itflag = PETSC_FALSE; 205 ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 206 if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 207 if (reason) { 208 /* We're not progressing, so return with the current iterate */ 209 ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 210 breakout = PETSC_TRUE; 211 break; 212 } 213 snes->numFailures++; 214 } 215 if (!breakout) { 216 /* Update function and solution vectors */ 217 fnorm = gnorm; 218 ierr = VecCopy(G,F);CHKERRQ(ierr); 219 ierr = VecCopy(Y,X);CHKERRQ(ierr); 220 /* Monitor convergence */ 221 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 222 snes->iter = i+1; 223 snes->norm = fnorm; 224 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 225 ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 226 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 227 /* Test for convergence, xnorm = || X || */ 228 neP->itflag = PETSC_TRUE; 229 if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 230 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 231 if (reason) break; 232 } else break; 233 } 234 if (i == maxits) { 235 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 236 if (!reason) reason = SNES_DIVERGED_MAX_IT; 237 } 238 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 239 snes->reason = reason; 240 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 241 PetscFunctionReturn(0); 242 } 243 /*------------------------------------------------------------*/ 244 #undef __FUNCT__ 245 #define __FUNCT__ "SNESSetUp_NEWTONTR" 246 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 247 { 248 PetscErrorCode ierr; 249 250 PetscFunctionBegin; 251 ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 252 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 253 PetscFunctionReturn(0); 254 } 255 256 #undef __FUNCT__ 257 #define __FUNCT__ "SNESReset_NEWTONTR" 258 PetscErrorCode SNESReset_NEWTONTR(SNES snes) 259 { 260 261 PetscFunctionBegin; 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "SNESDestroy_NEWTONTR" 267 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 268 { 269 PetscErrorCode ierr; 270 271 PetscFunctionBegin; 272 ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 273 ierr = PetscFree(snes->data);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 /*------------------------------------------------------------*/ 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "SNESSetFromOptions_NEWTONTR" 280 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes) 281 { 282 SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 287 ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 288 ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 289 ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 290 ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 291 ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 292 ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 293 ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 294 ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 295 ierr = PetscOptionsTail();CHKERRQ(ierr); 296 PetscFunctionReturn(0); 297 } 298 299 #undef __FUNCT__ 300 #define __FUNCT__ "SNESView_NEWTONTR" 301 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 302 { 303 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 304 PetscErrorCode ierr; 305 PetscBool iascii; 306 307 PetscFunctionBegin; 308 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 309 if (iascii) { 310 ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 311 ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);CHKERRQ(ierr); 312 } 313 PetscFunctionReturn(0); 314 } 315 /* ------------------------------------------------------------ */ 316 /*MC 317 SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 318 319 Options Database: 320 + -snes_trtol <tol> Trust region tolerance 321 . -snes_tr_mu <mu> 322 . -snes_tr_eta <eta> 323 . -snes_tr_sigma <sigma> 324 . -snes_tr_delta0 <delta0> 325 . -snes_tr_delta1 <delta1> 326 . -snes_tr_delta2 <delta2> 327 - -snes_tr_delta3 <delta3> 328 329 The basic algorithm is taken from "The Minpack Project", by More', 330 Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 331 of Mathematical Software", Wayne Cowell, editor. 332 333 This is intended as a model implementation, since it does not 334 necessarily have many of the bells and whistles of other 335 implementations. 336 337 Level: intermediate 338 339 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 340 341 M*/ 342 #undef __FUNCT__ 343 #define __FUNCT__ "SNESCreate_NEWTONTR" 344 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 345 { 346 SNES_NEWTONTR *neP; 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 snes->ops->setup = SNESSetUp_NEWTONTR; 351 snes->ops->solve = SNESSolve_NEWTONTR; 352 snes->ops->destroy = SNESDestroy_NEWTONTR; 353 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 354 snes->ops->view = SNESView_NEWTONTR; 355 snes->ops->reset = SNESReset_NEWTONTR; 356 357 snes->usesksp = PETSC_TRUE; 358 snes->usespc = PETSC_FALSE; 359 360 ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 361 snes->data = (void*)neP; 362 neP->mu = 0.25; 363 neP->eta = 0.75; 364 neP->delta = 0.0; 365 neP->delta0 = 0.2; 366 neP->delta1 = 0.3; 367 neP->delta2 = 0.75; 368 neP->delta3 = 2.0; 369 neP->sigma = 0.0001; 370 neP->itflag = PETSC_FALSE; 371 neP->rnorm0 = 0.0; 372 neP->ttol = 0.0; 373 PetscFunctionReturn(0); 374 } 375 376