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 98 PetscFunctionBegin; 99 if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 100 101 maxits = snes->max_its; /* maximum number of iterations */ 102 X = snes->vec_sol; /* solution vector */ 103 F = snes->vec_func; /* residual vector */ 104 Y = snes->work[0]; /* work vectors */ 105 G = snes->work[1]; 106 Ytmp = snes->work[2]; 107 108 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 109 snes->iter = 0; 110 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 111 112 if (!snes->vec_func_init_set) { 113 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 114 } else snes->vec_func_init_set = PETSC_FALSE; 115 116 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 117 SNESCheckFunctionNorm(snes,fnorm); 118 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 119 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 120 snes->norm = fnorm; 121 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 122 delta = xnorm ? neP->delta0*xnorm : neP->delta0; 123 neP->delta = delta; 124 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 125 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 126 127 /* test convergence */ 128 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 129 if (snes->reason) PetscFunctionReturn(0); 130 131 /* Set the stopping criteria to use the More' trick. */ 132 ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr); 133 if (!conv) { 134 SNES_TR_KSPConverged_Ctx *ctx; 135 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 136 ierr = PetscNew(&ctx);CHKERRQ(ierr); 137 ctx->snes = snes; 138 ierr = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr); 139 ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr); 140 ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 141 } 142 143 for (i=0; i<maxits; i++) { 144 145 /* Call general purpose update function */ 146 if (snes->ops->update) { 147 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 148 } 149 150 /* Solve J Y = F, where J is Jacobian matrix */ 151 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 152 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 153 ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 154 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 155 156 snes->linear_its += lits; 157 158 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 159 ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 160 norm1 = nrm; 161 while (1) { 162 ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 163 nrm = norm1; 164 165 /* Scale Y if need be and predict new value of F norm */ 166 if (nrm >= delta) { 167 nrm = delta/nrm; 168 gpnorm = (1.0 - nrm)*fnorm; 169 cnorm = nrm; 170 ierr = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 171 ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 172 nrm = gpnorm; 173 ynorm = delta; 174 } else { 175 gpnorm = 0.0; 176 ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 177 ynorm = nrm; 178 } 179 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 180 ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr); 181 ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 182 ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 183 if (fnorm == gpnorm) rho = 0.0; 184 else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 185 186 /* Update size of trust region */ 187 if (rho < neP->mu) delta *= neP->delta1; 188 else if (rho < neP->eta) delta *= neP->delta2; 189 else delta *= neP->delta3; 190 ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 191 ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 192 193 neP->delta = delta; 194 if (rho > neP->sigma) break; 195 ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 196 /* check to see if progress is hopeless */ 197 neP->itflag = PETSC_FALSE; 198 ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 199 if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); } 200 if (reason) { 201 /* We're not progressing, so return with the current iterate */ 202 ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 203 breakout = PETSC_TRUE; 204 break; 205 } 206 snes->numFailures++; 207 } 208 if (!breakout) { 209 /* Update function and solution vectors */ 210 fnorm = gnorm; 211 ierr = VecCopy(G,F);CHKERRQ(ierr); 212 ierr = VecCopy(Y,X);CHKERRQ(ierr); 213 /* Monitor convergence */ 214 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 215 snes->iter = i+1; 216 snes->norm = fnorm; 217 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 218 ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 219 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 220 /* Test for convergence, xnorm = || X || */ 221 neP->itflag = PETSC_TRUE; 222 if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 223 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 224 if (reason) break; 225 } else break; 226 } 227 if (i == maxits) { 228 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 229 if (!reason) reason = SNES_DIVERGED_MAX_IT; 230 } 231 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 232 snes->reason = reason; 233 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 234 PetscFunctionReturn(0); 235 } 236 /*------------------------------------------------------------*/ 237 #undef __FUNCT__ 238 #define __FUNCT__ "SNESSetUp_NEWTONTR" 239 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 240 { 241 PetscErrorCode ierr; 242 243 PetscFunctionBegin; 244 ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 245 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "SNESReset_NEWTONTR" 251 PetscErrorCode SNESReset_NEWTONTR(SNES snes) 252 { 253 254 PetscFunctionBegin; 255 PetscFunctionReturn(0); 256 } 257 258 #undef __FUNCT__ 259 #define __FUNCT__ "SNESDestroy_NEWTONTR" 260 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 261 { 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 266 ierr = PetscFree(snes->data);CHKERRQ(ierr); 267 PetscFunctionReturn(0); 268 } 269 /*------------------------------------------------------------*/ 270 271 #undef __FUNCT__ 272 #define __FUNCT__ "SNESSetFromOptions_NEWTONTR" 273 static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes) 274 { 275 SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 276 PetscErrorCode ierr; 277 278 PetscFunctionBegin; 279 ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 280 ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 281 ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 282 ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 283 ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 284 ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 285 ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 286 ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 287 ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 288 ierr = PetscOptionsTail();CHKERRQ(ierr); 289 PetscFunctionReturn(0); 290 } 291 292 #undef __FUNCT__ 293 #define __FUNCT__ "SNESView_NEWTONTR" 294 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 295 { 296 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 297 PetscErrorCode ierr; 298 PetscBool iascii; 299 300 PetscFunctionBegin; 301 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 302 if (iascii) { 303 ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 304 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); 305 } 306 PetscFunctionReturn(0); 307 } 308 /* ------------------------------------------------------------ */ 309 /*MC 310 SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 311 312 Options Database: 313 + -snes_trtol <tol> Trust region tolerance 314 . -snes_tr_mu <mu> 315 . -snes_tr_eta <eta> 316 . -snes_tr_sigma <sigma> 317 . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 318 . -snes_tr_delta1 <delta1> 319 . -snes_tr_delta2 <delta2> 320 - -snes_tr_delta3 <delta3> 321 322 The basic algorithm is taken from "The Minpack Project", by More', 323 Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 324 of Mathematical Software", Wayne Cowell, editor. 325 326 Level: intermediate 327 328 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 329 330 M*/ 331 #undef __FUNCT__ 332 #define __FUNCT__ "SNESCreate_NEWTONTR" 333 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 334 { 335 SNES_NEWTONTR *neP; 336 PetscErrorCode ierr; 337 338 PetscFunctionBegin; 339 snes->ops->setup = SNESSetUp_NEWTONTR; 340 snes->ops->solve = SNESSolve_NEWTONTR; 341 snes->ops->destroy = SNESDestroy_NEWTONTR; 342 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 343 snes->ops->view = SNESView_NEWTONTR; 344 snes->ops->reset = SNESReset_NEWTONTR; 345 346 snes->usesksp = PETSC_TRUE; 347 snes->usespc = PETSC_FALSE; 348 349 ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 350 snes->data = (void*)neP; 351 neP->mu = 0.25; 352 neP->eta = 0.75; 353 neP->delta = 0.0; 354 neP->delta0 = 0.2; 355 neP->delta1 = 0.3; 356 neP->delta2 = 0.75; 357 neP->delta3 = 2.0; 358 neP->sigma = 0.0001; 359 neP->itflag = PETSC_FALSE; 360 neP->rnorm0 = 0.0; 361 neP->ttol = 0.0; 362 PetscFunctionReturn(0); 363 } 364 365