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