1 2 #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I "petscsnes.h" I*/ 3 4 typedef struct { 5 SNES snes; 6 /* Information on the regular SNES convergence test; which may have been user provided 7 Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho 8 Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private 9 */ 10 11 PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 12 PetscErrorCode (*convdestroy)(void*); 13 void *convctx; 14 } SNES_TRDC_KSPConverged_Ctx; 15 16 static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 17 { 18 SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx*)cctx; 19 SNES snes = ctx->snes; 20 SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC*)snes->data; 21 Vec x; 22 PetscReal nrm; 23 24 PetscFunctionBegin; 25 PetscCall((*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx)); 26 if (*reason) { 27 PetscCall(PetscInfo(snes,"Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n",n,(double)rnorm)); 28 } 29 /* Determine norm of solution */ 30 PetscCall(KSPBuildSolution(ksp,NULL,&x)); 31 PetscCall(VecNorm(x,NORM_2,&nrm)); 32 if (nrm >= neP->delta) { 33 PetscCall(PetscInfo(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm)); 34 *reason = KSP_CONVERGED_STEP_LENGTH; 35 } 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) 40 { 41 SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx*)cctx; 42 43 PetscFunctionBegin; 44 PetscCall((*ctx->convdestroy)(ctx->convctx)); 45 PetscCall(PetscFree(ctx)); 46 47 PetscFunctionReturn(0); 48 } 49 50 /* ---------------------------------------------------------------- */ 51 /* 52 SNESTRDC_Converged_Private -test convergence JUST for 53 the trust region tolerance. 54 55 */ 56 static PetscErrorCode SNESTRDC_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 57 { 58 SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC*)snes->data; 59 60 PetscFunctionBegin; 61 *reason = SNES_CONVERGED_ITERATING; 62 if (neP->delta < xnorm * snes->deltatol) { 63 PetscCall(PetscInfo(snes,"Diverged due to too small a trust region %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol)); 64 *reason = SNES_DIVERGED_TR_DELTA; 65 } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 66 PetscCall(PetscInfo(snes,"Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n",snes->max_funcs)); 67 *reason = SNES_DIVERGED_FUNCTION_COUNT; 68 } 69 PetscFunctionReturn(0); 70 } 71 72 /*@ 73 SNESNewtonTRDCGetRhoFlag - Get whether the solution update is within the trust-region. 74 75 Input Parameters: 76 . snes - the nonlinear solver object 77 78 Output Parameters: 79 . rho_flag: PETSC_TRUE if the solution update is in the trust-region; otherwise, PETSC_FALSE 80 81 Level: developer 82 83 @*/ 84 PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes,PetscBool *rho_flag) 85 { 86 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 90 PetscValidBoolPointer(rho_flag,2); 91 *rho_flag = tr->rho_satisfied; 92 PetscFunctionReturn(0); 93 } 94 95 /*@C 96 SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined. 97 Allows the user a chance to change or override the trust region decision. 98 99 Logically Collective on snes 100 101 Input Parameters: 102 + snes - the nonlinear solver object 103 . func - [optional] function evaluation routine, see SNESNewtonTRDCPreCheck() for the calling sequence 104 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 105 106 Level: intermediate 107 108 Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver. 109 110 .seealso: `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 111 @*/ 112 PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx) 113 { 114 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 115 116 PetscFunctionBegin; 117 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 118 if (func) tr->precheck = func; 119 if (ctx) tr->precheckctx = ctx; 120 PetscFunctionReturn(0); 121 } 122 123 /*@C 124 SNESNewtonTRDCGetPreCheck - Gets the pre-check function 125 126 Not collective 127 128 Input Parameter: 129 . snes - the nonlinear solver context 130 131 Output Parameters: 132 + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPreCheck() 133 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 134 135 Level: intermediate 136 137 .seealso: `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()` 138 @*/ 139 PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx) 140 { 141 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 142 143 PetscFunctionBegin; 144 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 145 if (func) *func = tr->precheck; 146 if (ctx) *ctx = tr->precheckctx; 147 PetscFunctionReturn(0); 148 } 149 150 /*@C 151 SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 152 function evaluation. Allows the user a chance to change or override the decision of the line search routine 153 154 Logically Collective on snes 155 156 Input Parameters: 157 + snes - the nonlinear solver object 158 . func - [optional] function evaluation routine, see SNESNewtonTRDCPostCheck() for the calling sequence 159 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 160 161 Level: intermediate 162 163 Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver while the function set in 164 SNESLineSearchSetPostCheck() is called AFTER the function evaluation. 165 166 .seealso: `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 167 @*/ 168 PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 169 { 170 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 171 172 PetscFunctionBegin; 173 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 174 if (func) tr->postcheck = func; 175 if (ctx) tr->postcheckctx = ctx; 176 PetscFunctionReturn(0); 177 } 178 179 /*@C 180 SNESNewtonTRDCGetPostCheck - Gets the post-check function 181 182 Not collective 183 184 Input Parameter: 185 . snes - the nonlinear solver context 186 187 Output Parameters: 188 + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck() 189 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 190 191 Level: intermediate 192 193 .seealso: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()` 194 @*/ 195 PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes,PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 196 { 197 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 198 199 PetscFunctionBegin; 200 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 201 if (func) *func = tr->postcheck; 202 if (ctx) *ctx = tr->postcheckctx; 203 PetscFunctionReturn(0); 204 } 205 206 /*@C 207 SNESNewtonTRDCPreCheck - Called before the step has been determined in SNESNEWTONTRDC 208 209 Logically Collective on snes 210 211 Input Parameters: 212 + snes - the solver 213 . X - The last solution 214 - Y - The step direction 215 216 Output Parameters: 217 . changed_Y - Indicator that the step direction Y has been changed. 218 219 Level: developer 220 221 .seealso: `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 222 @*/ 223 static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y) 224 { 225 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 226 227 PetscFunctionBegin; 228 *changed_Y = PETSC_FALSE; 229 if (tr->precheck) { 230 PetscCall((*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx)); 231 PetscValidLogicalCollectiveBool(snes,*changed_Y,4); 232 } 233 PetscFunctionReturn(0); 234 } 235 236 /*@C 237 SNESNewtonTRDCPostCheck - Called after the step has been determined in SNESNEWTONTRDC but before the function evaluation 238 239 Logically Collective on snes 240 241 Input Parameters: 242 + snes - the solver 243 . X - The last solution 244 . Y - The full step direction 245 - W - The updated solution, W = X - Y 246 247 Output Parameters: 248 + changed_Y - indicator if step has been changed 249 - changed_W - Indicator if the new candidate solution W has been changed. 250 251 Notes: 252 If Y is changed then W is recomputed as X - Y 253 254 Level: developer 255 256 .seealso: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 257 @*/ 258 static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 259 { 260 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 261 262 PetscFunctionBegin; 263 *changed_Y = PETSC_FALSE; 264 *changed_W = PETSC_FALSE; 265 if (tr->postcheck) { 266 PetscCall((*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx)); 267 PetscValidLogicalCollectiveBool(snes,*changed_Y,5); 268 PetscValidLogicalCollectiveBool(snes,*changed_W,6); 269 } 270 PetscFunctionReturn(0); 271 } 272 273 /* 274 SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 275 (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 276 nonlinear equations 277 278 */ 279 static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) 280 { 281 SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC*)snes->data; 282 Vec X,F,Y,G,W,GradF,YNtmp; 283 Vec YCtmp; 284 Mat jac; 285 PetscInt maxits,i,j,lits,inner_count,bs; 286 PetscReal rho,fnorm,gnorm,xnorm=0,delta,ynorm,temp_xnorm,temp_ynorm; /* TRDC inner iteration */ 287 PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */ 288 PetscReal deltaM,ynnorm,f0,mp,gTy,g,yTHy; /* rho calculation */ 289 PetscReal auk,gfnorm,ycnorm,c0,c1,c2,tau,tau_pos,tau_neg,gTBg; /* Cauchy Point */ 290 KSP ksp; 291 SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 292 PetscBool breakout = PETSC_FALSE; 293 SNES_TRDC_KSPConverged_Ctx *ctx; 294 PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),(*convdestroy)(void*); 295 void *convctx; 296 297 PetscFunctionBegin; 298 maxits = snes->max_its; /* maximum number of iterations */ 299 X = snes->vec_sol; /* solution vector */ 300 F = snes->vec_func; /* residual vector */ 301 Y = snes->work[0]; /* update vector */ 302 G = snes->work[1]; /* updated residual */ 303 W = snes->work[2]; /* temporary vector */ 304 GradF = snes->work[3]; /* grad f = J^T F */ 305 YNtmp = snes->work[4]; /* Newton solution */ 306 YCtmp = snes->work[5]; /* Cauchy solution */ 307 308 PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 309 310 PetscCall(VecGetBlockSize(YNtmp,&bs)); 311 312 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 313 snes->iter = 0; 314 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 315 316 /* Set the linear stopping criteria to use the More' trick. From tr.c */ 317 PetscCall(SNESGetKSP(snes,&ksp)); 318 PetscCall(KSPGetConvergenceTest(ksp,&convtest,&convctx,&convdestroy)); 319 if (convtest != SNESTRDC_KSPConverged_Private) { 320 PetscCall(PetscNew(&ctx)); 321 ctx->snes = snes; 322 PetscCall(KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy)); 323 PetscCall(KSPSetConvergenceTest(ksp,SNESTRDC_KSPConverged_Private,ctx,SNESTRDC_KSPConverged_Destroy)); 324 PetscCall(PetscInfo(snes,"Using Krylov convergence test SNESTRDC_KSPConverged_Private\n")); 325 } 326 327 if (!snes->vec_func_init_set) { 328 PetscCall(SNESComputeFunction(snes,X,F)); /* F(X) */ 329 } else snes->vec_func_init_set = PETSC_FALSE; 330 331 PetscCall(VecNorm(F,NORM_2,&fnorm)); /* fnorm <- || F || */ 332 SNESCheckFunctionNorm(snes,fnorm); 333 PetscCall(VecNorm(X,NORM_2,&xnorm)); /* xnorm <- || X || */ 334 335 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 336 snes->norm = fnorm; 337 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 338 delta = xnorm ? neP->delta0*xnorm : neP->delta0; /* initial trust region size scaled by xnorm */ 339 deltaM = xnorm ? neP->deltaM*xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */ 340 neP->delta = delta; 341 PetscCall(SNESLogConvergenceHistory(snes,fnorm,0)); 342 PetscCall(SNESMonitor(snes,0,fnorm)); 343 344 neP->rho_satisfied = PETSC_FALSE; 345 346 /* test convergence */ 347 PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 348 if (snes->reason) PetscFunctionReturn(0); 349 350 for (i=0; i<maxits; i++) { 351 PetscBool changed_y; 352 PetscBool changed_w; 353 354 /* dogleg method */ 355 PetscCall(SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre)); 356 SNESCheckJacobianDomainerror(snes); 357 PetscCall(KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian)); 358 PetscCall(KSPSolve(snes->ksp,F,YNtmp)); /* Quasi Newton Solution */ 359 SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/ 360 PetscCall(KSPGetIterationNumber(snes->ksp,&lits)); 361 PetscCall(SNESGetJacobian(snes,&jac,NULL,NULL,NULL)); 362 363 /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable) 364 for inner iteration and Cauchy direction calculation 365 */ 366 if (bs > 1 && neP->auto_scale_multiphase) { 367 PetscCall(VecStrideNormAll(YNtmp,NORM_INFINITY,inorms)); 368 for (j=0; j<bs; j++) { 369 if (neP->auto_scale_max > 1.0) { 370 if (inorms[j] < 1.0/neP->auto_scale_max) { 371 inorms[j] = 1.0/neP->auto_scale_max; 372 } 373 } 374 PetscCall(VecStrideSet(W,j,inorms[j])); 375 PetscCall(VecStrideScale(YNtmp,j,1.0/inorms[j])); 376 PetscCall(VecStrideScale(X,j,1.0/inorms[j])); 377 } 378 PetscCall(VecNorm(X,NORM_2,&xnorm)); 379 if (i == 0) { 380 delta = neP->delta0*xnorm; 381 } else { 382 delta = neP->delta*xnorm; 383 } 384 deltaM = neP->deltaM*xnorm; 385 PetscCall(MatDiagonalScale(jac,PETSC_NULL,W)); 386 } 387 388 /* calculating GradF of minimization function */ 389 PetscCall(MatMultTranspose(jac,F,GradF)); /* grad f = J^T F */ 390 PetscCall(VecNorm(YNtmp,NORM_2,&ynnorm)); /* ynnorm <- || Y_newton || */ 391 392 inner_count = 0; 393 neP->rho_satisfied = PETSC_FALSE; 394 while (1) { 395 if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */ 396 PetscCall(VecCopy(YNtmp,Y)); 397 } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */ 398 PetscCall(MatMult(jac,GradF,W)); 399 PetscCall(VecDotRealPart(W,W,&gTBg)); /* completes GradF^T J^T J GradF */ 400 PetscCall(VecNorm(GradF,NORM_2,&gfnorm)); /* grad f norm <- || grad f || */ 401 if (gTBg <= 0.0) { 402 auk = PETSC_MAX_REAL; 403 } else { 404 auk = PetscSqr(gfnorm)/gTBg; 405 } 406 auk = PetscMin(delta/gfnorm,auk); 407 PetscCall(VecCopy(GradF,YCtmp)); /* this could be improved */ 408 PetscCall(VecScale(YCtmp,auk)); /* YCtmp, Cauchy solution*/ 409 PetscCall(VecNorm(YCtmp,NORM_2,&ycnorm)); /* ycnorm <- || Y_cauchy || */ 410 if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */ 411 PetscCall(VecCopy(YCtmp,Y)); 412 PetscCall(PetscInfo(snes,"DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n",(double)delta,(double)ynnorm,(double)ycnorm)); 413 } else { /* take ratio, tau, of Cauchy and Newton direction and step */ 414 PetscCall(VecAXPY(YNtmp,-1.0,YCtmp)); /* YCtmp = A, YNtmp = B */ 415 PetscCall(VecNorm(YNtmp,NORM_2,&c0)); /* this could be improved */ 416 c0 = PetscSqr(c0); 417 PetscCall(VecDotRealPart(YCtmp,YNtmp,&c1)); 418 c1 = 2.0*c1; 419 PetscCall(VecNorm(YCtmp,NORM_2,&c2)); /* this could be improved */ 420 c2 = PetscSqr(c2) - PetscSqr(delta); 421 tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4.*c0*c2))/(2.*c0); /* quadratic formula */ 422 tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4.*c0*c2))/(2.*c0); 423 tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */ 424 PetscCall(PetscInfo(snes,"DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n",(double)tau,(double)ynnorm,(double)ycnorm)); 425 PetscCall(VecWAXPY(W,tau,YNtmp,YCtmp)); 426 PetscCall(VecAXPY(W,-tau,YCtmp)); 427 PetscCall(VecCopy(W, Y)); /* this could be improved */ 428 } 429 } else { 430 /* if Cauchy is disabled, only use Newton direction */ 431 auk = delta/ynnorm; 432 PetscCall(VecScale(YNtmp,auk)); 433 PetscCall(VecCopy(YNtmp,Y)); /* this could be improved (many VecCopy, VecNorm)*/ 434 } 435 436 PetscCall(VecNorm(Y,NORM_2,&ynorm)); /* compute the final ynorm */ 437 f0 = 0.5*PetscSqr(fnorm); /* minimizing function f(X) */ 438 PetscCall(MatMult(jac,Y,W)); 439 PetscCall(VecDotRealPart(W,W,&yTHy)); /* completes GradY^T J^T J GradY */ 440 PetscCall(VecDotRealPart(GradF,Y,&gTy)); 441 mp = f0 - gTy + 0.5*yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/ 442 443 /* scale back solution update */ 444 if (bs > 1 && neP->auto_scale_multiphase) { 445 for (j=0; j<bs; j++) { 446 PetscCall(VecStrideScale(Y,j,inorms[j])); 447 if (inner_count == 0) { 448 /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */ 449 /* need to scale back X to match Y and provide proper update to the external code */ 450 PetscCall(VecStrideScale(X,j,inorms[j])); 451 } 452 } 453 if (inner_count == 0) PetscCall(VecNorm(X,NORM_2,&temp_xnorm)); /* only in the first iteration */ 454 PetscCall(VecNorm(Y,NORM_2,&temp_ynorm)); 455 } else { 456 temp_xnorm = xnorm; 457 temp_ynorm = ynorm; 458 } 459 inner_count++; 460 461 /* Evaluate the solution to meet the improvement ratio criteria */ 462 PetscCall(SNESNewtonTRDCPreCheck(snes,X,Y,&changed_y)); 463 PetscCall(VecWAXPY(W,-1.0,Y,X)); 464 PetscCall(SNESNewtonTRDCPostCheck(snes,X,Y,W,&changed_y,&changed_w)); 465 if (changed_y) PetscCall(VecWAXPY(W,-1.0,Y,X)); 466 PetscCall(VecCopy(Y,snes->vec_sol_update)); 467 PetscCall(SNESComputeFunction(snes,W,G)); /* F(X-Y) = G */ 468 PetscCall(VecNorm(G,NORM_2,&gnorm)); /* gnorm <- || g || */ 469 SNESCheckFunctionNorm(snes,gnorm); 470 g = 0.5*PetscSqr(gnorm); /* minimizing function g(W) */ 471 if (f0 == mp) rho = 0.0; 472 else rho = (f0 - g)/(f0 - mp); /* actual improvement over predicted improvement */ 473 474 if (rho < neP->eta2) { 475 delta *= neP->t1; /* shrink the region */ 476 } else if (rho > neP->eta3) { 477 delta = PetscMin(neP->t2*delta,deltaM); /* expand the region, but not greater than deltaM */ 478 } 479 480 neP->delta = delta; 481 if (rho >= neP->eta1) { 482 /* unscale delta and xnorm before going to the next outer iteration */ 483 if (bs > 1 && neP->auto_scale_multiphase) { 484 neP->delta = delta/xnorm; 485 xnorm = temp_xnorm; 486 ynorm = temp_ynorm; 487 } 488 neP->rho_satisfied = PETSC_TRUE; 489 break; /* the improvement ratio is satisfactory */ 490 } 491 PetscCall(PetscInfo(snes,"Trying again in smaller region\n")); 492 493 /* check to see if progress is hopeless */ 494 neP->itflag = PETSC_FALSE; 495 /* both delta, ynorm, and xnorm are either scaled or unscaled */ 496 PetscCall(SNESTRDC_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP)); 497 if (!reason) { 498 /* temp_xnorm, temp_ynorm is always unscaled */ 499 /* also the inner iteration already calculated the Jacobian and solved the matrix */ 500 /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */ 501 PetscCall((*snes->ops->converged)(snes,snes->iter+1,temp_xnorm,temp_ynorm,fnorm,&reason,snes->cnvP)); 502 } 503 /* if multiphase state changes, break out inner iteration */ 504 if (reason == SNES_BREAKOUT_INNER_ITER) { 505 if (bs > 1 && neP->auto_scale_multiphase) { 506 /* unscale delta and xnorm before going to the next outer iteration */ 507 neP->delta = delta/xnorm; 508 xnorm = temp_xnorm; 509 ynorm = temp_ynorm; 510 } 511 reason = SNES_CONVERGED_ITERATING; 512 break; 513 } 514 if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 515 if (reason) { 516 if (reason < 0) { 517 /* We're not progressing, so return with the current iterate */ 518 PetscCall(SNESMonitor(snes,i+1,fnorm)); 519 breakout = PETSC_TRUE; 520 break; 521 } else if (reason > 0) { 522 /* We're converged, so return with the current iterate and update solution */ 523 PetscCall(SNESMonitor(snes,i+1,fnorm)); 524 breakout = PETSC_FALSE; 525 break; 526 } 527 } 528 snes->numFailures++; 529 } 530 if (!breakout) { 531 /* Update function and solution vectors */ 532 fnorm = gnorm; 533 PetscCall(VecCopy(G,F)); 534 PetscCall(VecCopy(W,X)); 535 /* Monitor convergence */ 536 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 537 snes->iter = i+1; 538 snes->norm = fnorm; 539 snes->xnorm = xnorm; 540 snes->ynorm = ynorm; 541 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 542 PetscCall(SNESLogConvergenceHistory(snes,snes->norm,lits)); 543 PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 544 /* Test for convergence, xnorm = || X || */ 545 neP->itflag = PETSC_TRUE; 546 if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X,NORM_2,&xnorm)); 547 PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP)); 548 if (reason) break; 549 } else break; 550 } 551 552 /* PetscCall(PetscFree(inorms)); */ 553 if (i == maxits) { 554 PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",maxits)); 555 if (!reason) reason = SNES_DIVERGED_MAX_IT; 556 } 557 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 558 snes->reason = reason; 559 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 560 if (convtest != SNESTRDC_KSPConverged_Private) { 561 PetscCall(KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy)); 562 PetscCall(PetscFree(ctx)); 563 PetscCall(KSPSetConvergenceTest(ksp,convtest,convctx,convdestroy)); 564 } 565 PetscFunctionReturn(0); 566 } 567 568 /*------------------------------------------------------------*/ 569 static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) 570 { 571 PetscFunctionBegin; 572 PetscCall(SNESSetWorkVecs(snes,6)); 573 PetscCall(SNESSetUpMatrices(snes)); 574 PetscFunctionReturn(0); 575 } 576 577 PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) 578 { 579 PetscFunctionBegin; 580 PetscFunctionReturn(0); 581 } 582 583 static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) 584 { 585 PetscFunctionBegin; 586 PetscCall(SNESReset_NEWTONTRDC(snes)); 587 PetscCall(PetscFree(snes->data)); 588 PetscFunctionReturn(0); 589 } 590 /*------------------------------------------------------------*/ 591 592 static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(PetscOptionItems *PetscOptionsObject,SNES snes) 593 { 594 SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC*)snes->data; 595 596 PetscFunctionBegin; 597 PetscOptionsHeadBegin(PetscOptionsObject,"SNES trust region options for nonlinear equations"); 598 PetscCall(PetscOptionsReal("-snes_trdc_tol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL)); 599 PetscCall(PetscOptionsReal("-snes_trdc_eta1","eta1","None",ctx->eta1,&ctx->eta1,NULL)); 600 PetscCall(PetscOptionsReal("-snes_trdc_eta2","eta2","None",ctx->eta2,&ctx->eta2,NULL)); 601 PetscCall(PetscOptionsReal("-snes_trdc_eta3","eta3","None",ctx->eta3,&ctx->eta3,NULL)); 602 PetscCall(PetscOptionsReal("-snes_trdc_t1","t1","None",ctx->t1,&ctx->t1,NULL)); 603 PetscCall(PetscOptionsReal("-snes_trdc_t2","t2","None",ctx->t2,&ctx->t2,NULL)); 604 PetscCall(PetscOptionsReal("-snes_trdc_deltaM","deltaM","None",ctx->deltaM,&ctx->deltaM,NULL)); 605 PetscCall(PetscOptionsReal("-snes_trdc_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL)); 606 PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max","auto_scale_max","None",ctx->auto_scale_max,&ctx->auto_scale_max,NULL)); 607 PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy","use_cauchy","use Cauchy step and direction",ctx->use_cauchy,&ctx->use_cauchy,NULL)); 608 PetscCall(PetscOptionsBool("-snes_trdc_auto_scale_multiphase","auto_scale_multiphase","Auto scaling for proper cauchy direction",ctx->auto_scale_multiphase,&ctx->auto_scale_multiphase,NULL)); 609 PetscOptionsHeadEnd(); 610 PetscFunctionReturn(0); 611 } 612 613 static PetscErrorCode SNESView_NEWTONTRDC(SNES snes,PetscViewer viewer) 614 { 615 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC*)snes->data; 616 PetscBool iascii; 617 618 PetscFunctionBegin; 619 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 620 if (iascii) { 621 PetscCall(PetscViewerASCIIPrintf(viewer," Trust region tolerance %g (-snes_trtol)\n",(double)snes->deltatol)); 622 PetscCall(PetscViewerASCIIPrintf(viewer," eta1=%g, eta2=%g, eta3=%g\n",(double)tr->eta1,(double)tr->eta2,(double)tr->eta3)); 623 PetscCall(PetscViewerASCIIPrintf(viewer," delta0=%g, t1=%g, t2=%g, deltaM=%g\n",(double)tr->delta0,(double)tr->t1,(double)tr->t2,(double)tr->deltaM)); 624 } 625 PetscFunctionReturn(0); 626 } 627 /* ------------------------------------------------------------ */ 628 /*MC 629 SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 630 631 Options Database: 632 + -snes_trdc_tol <tol> - trust region tolerance 633 . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 634 . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 635 . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 636 . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 637 . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 638 . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5) 639 . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1) 640 . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor 641 . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm 642 - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region 643 644 Notes: 645 The algorithm is taken from "Linear and Nonlinear Solvers for Simulating Multiphase Flow 646 within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond, 647 Albert J. Valocchi, Tara LaForce. 648 649 Level: intermediate 650 651 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, `SNESNEWTONTRDC` 652 653 M*/ 654 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) 655 { 656 SNES_NEWTONTRDC *neP; 657 658 PetscFunctionBegin; 659 snes->ops->setup = SNESSetUp_NEWTONTRDC; 660 snes->ops->solve = SNESSolve_NEWTONTRDC; 661 snes->ops->destroy = SNESDestroy_NEWTONTRDC; 662 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC; 663 snes->ops->view = SNESView_NEWTONTRDC; 664 snes->ops->reset = SNESReset_NEWTONTRDC; 665 666 snes->usesksp = PETSC_TRUE; 667 snes->usesnpc = PETSC_FALSE; 668 669 snes->alwayscomputesfinalresidual = PETSC_TRUE; 670 671 PetscCall(PetscNewLog(snes,&neP)); 672 snes->data = (void*)neP; 673 neP->delta = 0.0; 674 neP->delta0 = 0.1; 675 neP->eta1 = 0.001; 676 neP->eta2 = 0.25; 677 neP->eta3 = 0.75; 678 neP->t1 = 0.25; 679 neP->t2 = 2.0; 680 neP->deltaM = 0.5; 681 neP->sigma = 0.0001; 682 neP->itflag = PETSC_FALSE; 683 neP->rnorm0 = 0.0; 684 neP->ttol = 0.0; 685 neP->use_cauchy = PETSC_TRUE; 686 neP->auto_scale_multiphase = PETSC_FALSE; 687 neP->auto_scale_max = -1.0; 688 neP->rho_satisfied = PETSC_FALSE; 689 snes->deltatol = 1.e-12; 690 691 /* for multiphase (multivariable) scaling */ 692 /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13 693 on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now. 694 PetscCall(VecGetBlockSize(snes->work[0],&neP->bs)); 695 PetscCall(PetscCalloc1(neP->bs,&neP->inorms)); 696 */ 697 698 PetscFunctionReturn(0); 699 } 700