1c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 24800dd8cSBarry Smith 3971273eeSBarry Smith typedef struct { 4971273eeSBarry Smith SNES snes; 5df8705c3SBarry Smith /* Information on the regular SNES convergence test; which may have been user provided */ 6df8705c3SBarry Smith PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 7df8705c3SBarry Smith PetscErrorCode (*convdestroy)(void *); 8df8705c3SBarry Smith void *convctx; 9971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx; 10971273eeSBarry Smith 119371c9d4SSatish Balay static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) { 12971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 13971273eeSBarry Smith SNES snes = ctx->snes; 1404d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 15df60cc22SBarry Smith Vec x; 16064f8208SBarry Smith PetscReal nrm; 17df60cc22SBarry Smith 183a40ed3dSBarry Smith PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 20*48a46eb9SPierre Jolivet if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 21a935fc98SLois Curfman McInnes /* Determine norm of solution */ 229566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 239566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &nrm)); 24064f8208SBarry Smith if (nrm >= neP->delta) { 259566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 26329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 27df60cc22SBarry Smith } 283a40ed3dSBarry Smith PetscFunctionReturn(0); 29df60cc22SBarry Smith } 3082bf6240SBarry Smith 319371c9d4SSatish Balay static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) { 32971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 33971273eeSBarry Smith 34971273eeSBarry Smith PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 369566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 37971273eeSBarry Smith PetscFunctionReturn(0); 38971273eeSBarry Smith } 39971273eeSBarry Smith 4085385478SLisandro Dalcin /* ---------------------------------------------------------------- */ 4185385478SLisandro Dalcin /* 42470880abSPatrick Sanan SNESTR_Converged_Private -test convergence JUST for 4385385478SLisandro Dalcin the trust region tolerance. 4485385478SLisandro Dalcin 4585385478SLisandro Dalcin */ 469371c9d4SSatish Balay static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) { 4704d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 4885385478SLisandro Dalcin 4985385478SLisandro Dalcin PetscFunctionBegin; 5085385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 5185385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 529566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to trust region param %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol)); 531c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 54e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 5563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 5685385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 5785385478SLisandro Dalcin } 5885385478SLisandro Dalcin PetscFunctionReturn(0); 5985385478SLisandro Dalcin } 6085385478SLisandro Dalcin 617cb011f5SBarry Smith /*@C 62c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 63c9368356SGlenn Hammond Allows the user a chance to change or override the decision of the line search routine. 64c9368356SGlenn Hammond 65c9368356SGlenn Hammond Logically Collective on snes 66c9368356SGlenn Hammond 67c9368356SGlenn Hammond Input Parameters: 68c9368356SGlenn Hammond + snes - the nonlinear solver object 69c9368356SGlenn Hammond . func - [optional] function evaluation routine, see SNESNewtonTRPreCheck() for the calling sequence 70c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 71c9368356SGlenn Hammond 72c9368356SGlenn Hammond Level: intermediate 73c9368356SGlenn Hammond 74c9368356SGlenn Hammond Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver. 75c9368356SGlenn Hammond 76db781477SPatrick Sanan .seealso: `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 77c9368356SGlenn Hammond @*/ 789371c9d4SSatish Balay PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) { 79c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 80c9368356SGlenn Hammond 81c9368356SGlenn Hammond PetscFunctionBegin; 82c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 83c9368356SGlenn Hammond if (func) tr->precheck = func; 84c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 85c9368356SGlenn Hammond PetscFunctionReturn(0); 86c9368356SGlenn Hammond } 87c9368356SGlenn Hammond 88c9368356SGlenn Hammond /*@C 89c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 90c9368356SGlenn Hammond 91c9368356SGlenn Hammond Not collective 92c9368356SGlenn Hammond 93c9368356SGlenn Hammond Input Parameter: 94c9368356SGlenn Hammond . snes - the nonlinear solver context 95c9368356SGlenn Hammond 96c9368356SGlenn Hammond Output Parameters: 97c9368356SGlenn Hammond + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck() 98c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 99c9368356SGlenn Hammond 100c9368356SGlenn Hammond Level: intermediate 101c9368356SGlenn Hammond 102db781477SPatrick Sanan .seealso: `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 103c9368356SGlenn Hammond @*/ 1049371c9d4SSatish Balay PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) { 105c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 106c9368356SGlenn Hammond 107c9368356SGlenn Hammond PetscFunctionBegin; 108c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 109c9368356SGlenn Hammond if (func) *func = tr->precheck; 110c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 111c9368356SGlenn Hammond PetscFunctionReturn(0); 112c9368356SGlenn Hammond } 113c9368356SGlenn Hammond 114c9368356SGlenn Hammond /*@C 1157cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 1167cb011f5SBarry Smith function evaluation. Allows the user a chance to change or override the decision of the line search routine 1177cb011f5SBarry Smith 1187cb011f5SBarry Smith Logically Collective on snes 1197cb011f5SBarry Smith 1207cb011f5SBarry Smith Input Parameters: 1217cb011f5SBarry Smith + snes - the nonlinear solver object 1227cb011f5SBarry Smith . func - [optional] function evaluation routine, see SNESNewtonTRPostCheck() for the calling sequence 1237cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1247cb011f5SBarry Smith 1257cb011f5SBarry Smith Level: intermediate 1267cb011f5SBarry Smith 127c9368356SGlenn Hammond Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in 1287cb011f5SBarry Smith SNESLineSearchSetPostCheck() is called AFTER the function evaluation. 1297cb011f5SBarry Smith 130db781477SPatrick Sanan .seealso: `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()` 1317cb011f5SBarry Smith @*/ 1329371c9d4SSatish Balay PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) { 1337cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1347cb011f5SBarry Smith 1357cb011f5SBarry Smith PetscFunctionBegin; 1367cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1377cb011f5SBarry Smith if (func) tr->postcheck = func; 1387cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 1397cb011f5SBarry Smith PetscFunctionReturn(0); 1407cb011f5SBarry Smith } 1417cb011f5SBarry Smith 1427cb011f5SBarry Smith /*@C 1437cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 1447cb011f5SBarry Smith 1457cb011f5SBarry Smith Not collective 1467cb011f5SBarry Smith 1477cb011f5SBarry Smith Input Parameter: 1487cb011f5SBarry Smith . snes - the nonlinear solver context 1497cb011f5SBarry Smith 1507cb011f5SBarry Smith Output Parameters: 1517cb011f5SBarry Smith + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck() 1527cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1537cb011f5SBarry Smith 1547cb011f5SBarry Smith Level: intermediate 1557cb011f5SBarry Smith 156db781477SPatrick Sanan .seealso: `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 1577cb011f5SBarry Smith @*/ 1589371c9d4SSatish Balay PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) { 1597cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1607cb011f5SBarry Smith 1617cb011f5SBarry Smith PetscFunctionBegin; 1627cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1637cb011f5SBarry Smith if (func) *func = tr->postcheck; 1647cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 1657cb011f5SBarry Smith PetscFunctionReturn(0); 1667cb011f5SBarry Smith } 1677cb011f5SBarry Smith 1687cb011f5SBarry Smith /*@C 169c9368356SGlenn Hammond SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR 170c9368356SGlenn Hammond 171c9368356SGlenn Hammond Logically Collective on snes 172c9368356SGlenn Hammond 173c9368356SGlenn Hammond Input Parameters: 174c9368356SGlenn Hammond + snes - the solver 175c9368356SGlenn Hammond . X - The last solution 176c9368356SGlenn Hammond - Y - The step direction 177c9368356SGlenn Hammond 178c9368356SGlenn Hammond Output Parameters: 179c9368356SGlenn Hammond . changed_Y - Indicator that the step direction Y has been changed. 180c9368356SGlenn Hammond 181c9368356SGlenn Hammond Level: developer 182c9368356SGlenn Hammond 183db781477SPatrick Sanan .seealso: `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 184c9368356SGlenn Hammond @*/ 1859371c9d4SSatish Balay static PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) { 186c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 187c9368356SGlenn Hammond 188c9368356SGlenn Hammond PetscFunctionBegin; 189c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 190c9368356SGlenn Hammond if (tr->precheck) { 1919566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 192c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 193c9368356SGlenn Hammond } 194c9368356SGlenn Hammond PetscFunctionReturn(0); 195c9368356SGlenn Hammond } 196c9368356SGlenn Hammond 197c9368356SGlenn Hammond /*@C 1987cb011f5SBarry Smith SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation 1997cb011f5SBarry Smith 2007cb011f5SBarry Smith Logically Collective on snes 2017cb011f5SBarry Smith 2027cb011f5SBarry Smith Input Parameters: 2036b867d5aSJose E. Roman + snes - the solver 2046b867d5aSJose E. Roman . X - The last solution 2057cb011f5SBarry Smith . Y - The full step direction 2063312a946SBarry Smith - W - The updated solution, W = X - Y 2077cb011f5SBarry Smith 2087cb011f5SBarry Smith Output Parameters: 2093312a946SBarry Smith + changed_Y - indicator if step has been changed 2103312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 2117cb011f5SBarry Smith 2127cb011f5SBarry Smith Notes: 2133312a946SBarry Smith If Y is changed then W is recomputed as X - Y 2147cb011f5SBarry Smith 2157cb011f5SBarry Smith Level: developer 2167cb011f5SBarry Smith 217db781477SPatrick Sanan .seealso: `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 2187cb011f5SBarry Smith @*/ 2199371c9d4SSatish Balay static PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) { 2207cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2217cb011f5SBarry Smith 2227cb011f5SBarry Smith PetscFunctionBegin; 223c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 2247cb011f5SBarry Smith *changed_W = PETSC_FALSE; 2257cb011f5SBarry Smith if (tr->postcheck) { 2269566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 227c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 2287cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 2297cb011f5SBarry Smith } 2307cb011f5SBarry Smith PetscFunctionReturn(0); 2317cb011f5SBarry Smith } 23285385478SLisandro Dalcin 233df60cc22SBarry Smith /* 23404d7464bSBarry Smith SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 235ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 2364800dd8cSBarry Smith 2374800dd8cSBarry Smith */ 2389371c9d4SSatish Balay static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) { 23904d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 240c9368356SGlenn Hammond Vec X, F, Y, G, Ytmp, W; 241a7cc72afSBarry Smith PetscInt maxits, i, lits; 24285385478SLisandro Dalcin PetscReal rho, fnorm, gnorm, gpnorm, xnorm = 0, delta, nrm, ynorm, norm1; 24379f36460SBarry Smith PetscScalar cnorm; 244df60cc22SBarry Smith KSP ksp; 2453f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 246df8705c3SBarry Smith PetscBool breakout = PETSC_FALSE; 247df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 2485e28bcb6Sprj- PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 2495e28bcb6Sprj- void *convctx; 2504800dd8cSBarry Smith 2513a40ed3dSBarry Smith PetscFunctionBegin; 2520b121fc5SBarry Smith 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); 253c579b300SPatrick Farrell 254fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 255fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 25639e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 257fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 258fbe28522SBarry Smith G = snes->work[1]; 2596b5873e3SBarry Smith Ytmp = snes->work[2]; 260c9368356SGlenn Hammond W = snes->work[3]; 2614800dd8cSBarry Smith 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 2634c49b128SBarry Smith snes->iter = 0; 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2654800dd8cSBarry Smith 266df8705c3SBarry Smith /* Set the linear stopping criteria to use the More' trick. */ 2679566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2689566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 269df8705c3SBarry Smith if (convtest != SNESTR_KSPConverged_Private) { 2709566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 271df8705c3SBarry Smith ctx->snes = snes; 2729566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 2739566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 2749566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 275df8705c3SBarry Smith } 276df8705c3SBarry Smith 277e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2789566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 2791aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 280e4ed7901SPeter Brune 2819566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 282422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 2839566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 2849566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 285fbe28522SBarry Smith snes->norm = fnorm; 2869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2878b84755bSBarry Smith delta = xnorm ? neP->delta0 * xnorm : neP->delta0; 2884800dd8cSBarry Smith neP->delta = delta; 2899566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 2909566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 291b37302c6SLois Curfman McInnes 29285385478SLisandro Dalcin /* test convergence */ 293dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 29485385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 2953f149594SLisandro Dalcin 2964800dd8cSBarry Smith for (i = 0; i < maxits; i++) { 29799a96b7cSMatthew Knepley /* Call general purpose update function */ 298dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 29999a96b7cSMatthew Knepley 30085385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 3019566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 30207b62357SFande Kong SNESCheckJacobianDomainerror(snes); 3039566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 3049566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Ytmp)); 3059566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 306329f5518SBarry Smith snes->linear_its += lits; 3071aa26658SKarl Rupp 3089566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 3099566063dSJacob Faibussowitsch PetscCall(VecNorm(Ytmp, NORM_2, &nrm)); 310064f8208SBarry Smith norm1 = nrm; 311284fb49fSHeeho Park 3126b5873e3SBarry Smith while (1) { 313c9368356SGlenn Hammond PetscBool changed_y; 3147cb011f5SBarry Smith PetscBool changed_w; 3159566063dSJacob Faibussowitsch PetscCall(VecCopy(Ytmp, Y)); 316064f8208SBarry Smith nrm = norm1; 3176b5873e3SBarry Smith 3182deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 319064f8208SBarry Smith if (nrm >= delta) { 320064f8208SBarry Smith nrm = delta / nrm; 321064f8208SBarry Smith gpnorm = (1.0 - nrm) * fnorm; 322064f8208SBarry Smith cnorm = nrm; 3239566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Scaling direction by %g\n", (double)nrm)); 3249566063dSJacob Faibussowitsch PetscCall(VecScale(Y, cnorm)); 325064f8208SBarry Smith nrm = gpnorm; 326fbe28522SBarry Smith ynorm = delta; 327fbe28522SBarry Smith } else { 328fbe28522SBarry Smith gpnorm = 0.0; 3299566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Direction is in Trust Region\n")); 330064f8208SBarry Smith ynorm = nrm; 331fbe28522SBarry Smith } 332c9368356SGlenn Hammond /* PreCheck() allows for updates to Y prior to W <- X - Y */ 3339566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 3349566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -1.0, Y, X)); /* W <- X - Y */ 3359566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 3369566063dSJacob Faibussowitsch if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X)); 3379566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 3389566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */ 3399566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */ 3407551ef16SBarry Smith SNESCheckFunctionNorm(snes, gnorm); 3414800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 342fbe28522SBarry Smith else rho = (fnorm * fnorm - gnorm * gnorm) / (fnorm * fnorm - gpnorm * gpnorm); 3434800dd8cSBarry Smith 3444800dd8cSBarry Smith /* Update size of trust region */ 3454800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 3464800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 3474800dd8cSBarry Smith else delta *= neP->delta3; 3489566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%g, gnorm=%g, ynorm=%g\n", (double)fnorm, (double)gnorm, (double)ynorm)); 3499566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "gpred=%g, rho=%g, delta=%g\n", (double)gpnorm, (double)rho, (double)delta)); 3501aa26658SKarl Rupp 3514800dd8cSBarry Smith neP->delta = delta; 3524800dd8cSBarry Smith if (rho > neP->sigma) break; 3539566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 354284fb49fSHeeho Park 3556b5873e3SBarry Smith /* check to see if progress is hopeless */ 356ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 3579566063dSJacob Faibussowitsch PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP)); 358dbbe0bcdSBarry Smith if (!reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 3597cb011f5SBarry Smith if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 360184914b5SBarry Smith if (reason) { 36152392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 3629566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 363a7cc72afSBarry Smith breakout = PETSC_TRUE; 364454a90a3SBarry Smith break; 36552392280SLois Curfman McInnes } 36650ffb88aSMatthew Knepley snes->numFailures++; 3674800dd8cSBarry Smith } 3681acabf8cSLois Curfman McInnes if (!breakout) { 36985385478SLisandro Dalcin /* Update function and solution vectors */ 3704800dd8cSBarry Smith fnorm = gnorm; 3719566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 3729566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 37385385478SLisandro Dalcin /* Monitor convergence */ 3749566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 375c509a162SBarry Smith snes->iter = i + 1; 376fbe28522SBarry Smith snes->norm = fnorm; 377c1e67a49SFande Kong snes->xnorm = xnorm; 378c1e67a49SFande Kong snes->ynorm = ynorm; 3799566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3809566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 3819566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 38285385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 383ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 3849566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 385dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 3863f149594SLisandro Dalcin if (reason) break; 3871aa26658SKarl Rupp } else break; 38838442cffSBarry Smith } 389284fb49fSHeeho Park 39052392280SLois Curfman McInnes if (i == maxits) { 3919566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 39285385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 39352392280SLois Curfman McInnes } 3949566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 395184914b5SBarry Smith snes->reason = reason; 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3975e28bcb6Sprj- if (convtest != SNESTR_KSPConverged_Private) { 3989566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 3999566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 4009566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 4015e28bcb6Sprj- } 4023a40ed3dSBarry Smith PetscFunctionReturn(0); 4034800dd8cSBarry Smith } 404284fb49fSHeeho Park 4054800dd8cSBarry Smith /*------------------------------------------------------------*/ 4069371c9d4SSatish Balay static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) { 4073a40ed3dSBarry Smith PetscFunctionBegin; 4089566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 4)); 4099566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 4103a40ed3dSBarry Smith PetscFunctionReturn(0); 4114800dd8cSBarry Smith } 4126b8b9a38SLisandro Dalcin 4139371c9d4SSatish Balay PetscErrorCode SNESReset_NEWTONTR(SNES snes) { 4146b8b9a38SLisandro Dalcin PetscFunctionBegin; 4156b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 4166b8b9a38SLisandro Dalcin } 4176b8b9a38SLisandro Dalcin 4189371c9d4SSatish Balay static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) { 4193a40ed3dSBarry Smith PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONTR(snes)); 4219566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4223a40ed3dSBarry Smith PetscFunctionReturn(0); 4234800dd8cSBarry Smith } 4244800dd8cSBarry Smith /*------------------------------------------------------------*/ 4254800dd8cSBarry Smith 4269371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) { 42704d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 4284800dd8cSBarry Smith 4293a40ed3dSBarry Smith PetscFunctionBegin; 430d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 4319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trtol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 4329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_mu", "mu", "None", ctx->mu, &ctx->mu, NULL)); 4339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_eta", "eta", "None", ctx->eta, &ctx->eta, NULL)); 4349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_sigma", "sigma", "None", ctx->sigma, &ctx->sigma, NULL)); 4359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 4369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta1", "delta1", "None", ctx->delta1, &ctx->delta1, NULL)); 4379566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta2", "delta2", "None", ctx->delta2, &ctx->delta2, NULL)); 4389566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta3", "delta3", "None", ctx->delta3, &ctx->delta3, NULL)); 439d0609cedSBarry Smith PetscOptionsHeadEnd(); 4403a40ed3dSBarry Smith PetscFunctionReturn(0); 4414800dd8cSBarry Smith } 4424800dd8cSBarry Smith 4439371c9d4SSatish Balay static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) { 44404d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 445ace3abfcSBarry Smith PetscBool iascii; 446a935fc98SLois Curfman McInnes 4473a40ed3dSBarry Smith PetscFunctionBegin; 4489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 44932077d6dSBarry Smith if (iascii) { 45063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol)); 4519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " mu=%g, eta=%g, sigma=%g\n", (double)tr->mu, (double)tr->eta, (double)tr->sigma)); 4529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, delta1=%g, delta2=%g, delta3=%g\n", (double)tr->delta0, (double)tr->delta1, (double)tr->delta2, (double)tr->delta3)); 45319bcc07fSBarry Smith } 4543a40ed3dSBarry Smith PetscFunctionReturn(0); 455a935fc98SLois Curfman McInnes } 45652392280SLois Curfman McInnes /* ------------------------------------------------------------ */ 457ebe3b25bSBarry Smith /*MC 45804d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 459ebe3b25bSBarry Smith 460ebe3b25bSBarry Smith Options Database: 4618e434772SBarry Smith + -snes_trtol <tol> - trust region tolerance 4628e434772SBarry Smith . -snes_tr_mu <mu> - trust region parameter 4638e434772SBarry Smith . -snes_tr_eta <eta> - trust region parameter 4648e434772SBarry Smith . -snes_tr_sigma <sigma> - trust region parameter 4658b84755bSBarry Smith . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 4668e434772SBarry Smith . -snes_tr_delta1 <delta1> - trust region parameter 4678e434772SBarry Smith . -snes_tr_delta2 <delta2> - trust region parameter 4688e434772SBarry Smith - -snes_tr_delta3 <delta3> - trust region parameter 469ebe3b25bSBarry Smith 470ebe3b25bSBarry Smith The basic algorithm is taken from "The Minpack Project", by More', 471ebe3b25bSBarry Smith Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 472ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 473ebe3b25bSBarry Smith 474ee3001cbSBarry Smith Level: intermediate 475ee3001cbSBarry Smith 476db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()` 477ebe3b25bSBarry Smith 478ebe3b25bSBarry Smith M*/ 4799371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) { 48004d7464bSBarry Smith SNES_NEWTONTR *neP; 4814800dd8cSBarry Smith 4823a40ed3dSBarry Smith PetscFunctionBegin; 48304d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 48404d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 48504d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 48604d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 48704d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 48804d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 489fbe28522SBarry Smith 49042f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 491efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 49242f4f86dSBarry Smith 4934fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 4944fc747eaSLawrence Mitchell 4959566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes, &neP)); 496fbe28522SBarry Smith snes->data = (void *)neP; 497fbe28522SBarry Smith neP->mu = 0.25; 498fbe28522SBarry Smith neP->eta = 0.75; 499fbe28522SBarry Smith neP->delta = 0.0; 500fbe28522SBarry Smith neP->delta0 = 0.2; 501fbe28522SBarry Smith neP->delta1 = 0.3; 502fbe28522SBarry Smith neP->delta2 = 0.75; 503fbe28522SBarry Smith neP->delta3 = 2.0; 504fbe28522SBarry Smith neP->sigma = 0.0001; 505ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 50675567043SBarry Smith neP->rnorm0 = 0.0; 50775567043SBarry Smith neP->ttol = 0.0; 5083a40ed3dSBarry Smith PetscFunctionReturn(0); 5094800dd8cSBarry Smith } 510