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 11d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 12d71ae5a4SJacob Faibussowitsch { 13971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 14971273eeSBarry Smith SNES snes = ctx->snes; 1504d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 16df60cc22SBarry Smith Vec x; 17064f8208SBarry Smith PetscReal nrm; 18df60cc22SBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 2148a46eb9SPierre Jolivet if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 22a935fc98SLois Curfman McInnes /* Determine norm of solution */ 239566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 249566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &nrm)); 25064f8208SBarry Smith if (nrm >= neP->delta) { 269566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 27329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 28df60cc22SBarry Smith } 293a40ed3dSBarry Smith PetscFunctionReturn(0); 30df60cc22SBarry Smith } 3182bf6240SBarry Smith 32d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 33d71ae5a4SJacob Faibussowitsch { 34971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 35971273eeSBarry Smith 36971273eeSBarry Smith PetscFunctionBegin; 379566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 389566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 39971273eeSBarry Smith PetscFunctionReturn(0); 40971273eeSBarry Smith } 41971273eeSBarry Smith 4285385478SLisandro Dalcin /* 43470880abSPatrick Sanan SNESTR_Converged_Private -test convergence JUST for 4485385478SLisandro Dalcin the trust region tolerance. 4585385478SLisandro Dalcin */ 46d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 47d71ae5a4SJacob Faibussowitsch { 4804d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 4985385478SLisandro Dalcin 5085385478SLisandro Dalcin PetscFunctionBegin; 5185385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 5285385478SLisandro Dalcin if (neP->delta < xnorm * snes->deltatol) { 539566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to trust region param %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol)); 541c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 55e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 5663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 5785385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 5885385478SLisandro Dalcin } 5985385478SLisandro Dalcin PetscFunctionReturn(0); 6085385478SLisandro Dalcin } 6185385478SLisandro Dalcin 627cb011f5SBarry Smith /*@C 63c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 64c9368356SGlenn Hammond Allows the user a chance to change or override the decision of the line search routine. 65c9368356SGlenn Hammond 66f6dfbefdSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 67f6dfbefdSBarry Smith 68*c3339decSBarry Smith Logically Collective 69c9368356SGlenn Hammond 70c9368356SGlenn Hammond Input Parameters: 71c9368356SGlenn Hammond + snes - the nonlinear solver object 72f6dfbefdSBarry Smith . func - [optional] function evaluation routine, see `SNESNewtonTRPreCheck()` for the calling sequence 73c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 74c9368356SGlenn Hammond 75c9368356SGlenn Hammond Level: intermediate 76c9368356SGlenn Hammond 77f6dfbefdSBarry Smith Note: 78f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTR` solver. 79c9368356SGlenn Hammond 80f6dfbefdSBarry Smith .seealso: `SNESNEWTONDCTRDC`, `SNESNEWTONDCTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 81c9368356SGlenn Hammond @*/ 82d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 83d71ae5a4SJacob Faibussowitsch { 84c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 85c9368356SGlenn Hammond 86c9368356SGlenn Hammond PetscFunctionBegin; 87c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 88c9368356SGlenn Hammond if (func) tr->precheck = func; 89c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 90c9368356SGlenn Hammond PetscFunctionReturn(0); 91c9368356SGlenn Hammond } 92c9368356SGlenn Hammond 93c9368356SGlenn Hammond /*@C 94c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 95c9368356SGlenn Hammond 96f6dfbefdSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 97f6dfbefdSBarry Smith 98c9368356SGlenn Hammond Not collective 99c9368356SGlenn Hammond 100c9368356SGlenn Hammond Input Parameter: 101c9368356SGlenn Hammond . snes - the nonlinear solver context 102c9368356SGlenn Hammond 103c9368356SGlenn Hammond Output Parameters: 104f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPreCheck()` 105c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 106c9368356SGlenn Hammond 107c9368356SGlenn Hammond Level: intermediate 108c9368356SGlenn Hammond 109f6dfbefdSBarry Smith .seealso: `SNESNEWTONDCTRDC`, `SNESNEWTONDCTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 110c9368356SGlenn Hammond @*/ 111d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 112d71ae5a4SJacob Faibussowitsch { 113c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 114c9368356SGlenn Hammond 115c9368356SGlenn Hammond PetscFunctionBegin; 116c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 117c9368356SGlenn Hammond if (func) *func = tr->precheck; 118c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 119c9368356SGlenn Hammond PetscFunctionReturn(0); 120c9368356SGlenn Hammond } 121c9368356SGlenn Hammond 122c9368356SGlenn Hammond /*@C 1237cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 1247cb011f5SBarry Smith function evaluation. Allows the user a chance to change or override the decision of the line search routine 1257cb011f5SBarry Smith 126f6dfbefdSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 127f6dfbefdSBarry Smith 128*c3339decSBarry Smith Logically Collective 1297cb011f5SBarry Smith 1307cb011f5SBarry Smith Input Parameters: 1317cb011f5SBarry Smith + snes - the nonlinear solver object 132f6dfbefdSBarry Smith . func - [optional] function evaluation routine, see `SNESNewtonTRPostCheck()` for the calling sequence 1337cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1347cb011f5SBarry Smith 1357cb011f5SBarry Smith Level: intermediate 1367cb011f5SBarry Smith 137f6dfbefdSBarry Smith Note: 138f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTR` solver while the function set in 139f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 1407cb011f5SBarry Smith 141f6dfbefdSBarry Smith .seealso: `SNESNEWTONDCTRDC`, `SNESNEWTONDCTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()` 1427cb011f5SBarry Smith @*/ 143d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 144d71ae5a4SJacob Faibussowitsch { 1457cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1467cb011f5SBarry Smith 1477cb011f5SBarry Smith PetscFunctionBegin; 1487cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1497cb011f5SBarry Smith if (func) tr->postcheck = func; 1507cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 1517cb011f5SBarry Smith PetscFunctionReturn(0); 1527cb011f5SBarry Smith } 1537cb011f5SBarry Smith 1547cb011f5SBarry Smith /*@C 1557cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 1567cb011f5SBarry Smith 157f6dfbefdSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 158f6dfbefdSBarry Smith 1597cb011f5SBarry Smith Not collective 1607cb011f5SBarry Smith 1617cb011f5SBarry Smith Input Parameter: 1627cb011f5SBarry Smith . snes - the nonlinear solver context 1637cb011f5SBarry Smith 1647cb011f5SBarry Smith Output Parameters: 165f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPostCheck()` 1667cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1677cb011f5SBarry Smith 1687cb011f5SBarry Smith Level: intermediate 1697cb011f5SBarry Smith 170f6dfbefdSBarry Smith .seealso: `SNESNEWTONDCTRDC`, `SNESNEWTONDCTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 1717cb011f5SBarry Smith @*/ 172d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 173d71ae5a4SJacob Faibussowitsch { 1747cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1757cb011f5SBarry Smith 1767cb011f5SBarry Smith PetscFunctionBegin; 1777cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1787cb011f5SBarry Smith if (func) *func = tr->postcheck; 1797cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 1807cb011f5SBarry Smith PetscFunctionReturn(0); 1817cb011f5SBarry Smith } 1827cb011f5SBarry Smith 1837cb011f5SBarry Smith /*@C 184f6dfbefdSBarry Smith SNESNewtonTRPreCheck - Called before the step has been determined in `SNESNEWTONTR` 185f6dfbefdSBarry Smith 186f6dfbefdSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 187c9368356SGlenn Hammond 188*c3339decSBarry Smith Logically Collective 189c9368356SGlenn Hammond 190c9368356SGlenn Hammond Input Parameters: 191c9368356SGlenn Hammond + snes - the solver 192c9368356SGlenn Hammond . X - The last solution 193c9368356SGlenn Hammond - Y - The step direction 194c9368356SGlenn Hammond 195c9368356SGlenn Hammond Output Parameters: 196c9368356SGlenn Hammond . changed_Y - Indicator that the step direction Y has been changed. 197c9368356SGlenn Hammond 198c9368356SGlenn Hammond Level: developer 199c9368356SGlenn Hammond 200f6dfbefdSBarry Smith .seealso: `SNESNEWTONDCTRDC`, `SNESNEWTONDCTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 201c9368356SGlenn Hammond @*/ 202d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 203d71ae5a4SJacob Faibussowitsch { 204c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 205c9368356SGlenn Hammond 206c9368356SGlenn Hammond PetscFunctionBegin; 207c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 208c9368356SGlenn Hammond if (tr->precheck) { 2099566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 210c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 211c9368356SGlenn Hammond } 212c9368356SGlenn Hammond PetscFunctionReturn(0); 213c9368356SGlenn Hammond } 214c9368356SGlenn Hammond 215c9368356SGlenn Hammond /*@C 216f6dfbefdSBarry Smith SNESNewtonTRPostCheck - Called after the step has been determined in `SNESNEWTONTR` but before the function evaluation 217f6dfbefdSBarry Smith 218f6dfbefdSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 2197cb011f5SBarry Smith 220*c3339decSBarry Smith Logically Collective 2217cb011f5SBarry Smith 2227cb011f5SBarry Smith Input Parameters: 2236b867d5aSJose E. Roman + snes - the solver 2246b867d5aSJose E. Roman . X - The last solution 2257cb011f5SBarry Smith . Y - The full step direction 2263312a946SBarry Smith - W - The updated solution, W = X - Y 2277cb011f5SBarry Smith 2287cb011f5SBarry Smith Output Parameters: 2293312a946SBarry Smith + changed_Y - indicator if step has been changed 2303312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 2317cb011f5SBarry Smith 232f6dfbefdSBarry Smith Note: 2333312a946SBarry Smith If Y is changed then W is recomputed as X - Y 2347cb011f5SBarry Smith 2357cb011f5SBarry Smith Level: developer 2367cb011f5SBarry Smith 237f6dfbefdSBarry Smith .seealso: `SNESNEWTONDCTRDC`, `SNESNEWTONDCTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 2387cb011f5SBarry Smith @*/ 239d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 240d71ae5a4SJacob Faibussowitsch { 2417cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2427cb011f5SBarry Smith 2437cb011f5SBarry Smith PetscFunctionBegin; 244c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 2457cb011f5SBarry Smith *changed_W = PETSC_FALSE; 2467cb011f5SBarry Smith if (tr->postcheck) { 2479566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 248c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 2497cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 2507cb011f5SBarry Smith } 2517cb011f5SBarry Smith PetscFunctionReturn(0); 2527cb011f5SBarry Smith } 25385385478SLisandro Dalcin 254df60cc22SBarry Smith /* 25504d7464bSBarry Smith SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 256ddbbdb52SLois Curfman McInnes region approach for solving systems of nonlinear equations. 2574800dd8cSBarry Smith 2584800dd8cSBarry Smith */ 259d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 260d71ae5a4SJacob Faibussowitsch { 26104d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 262c9368356SGlenn Hammond Vec X, F, Y, G, Ytmp, W; 263a7cc72afSBarry Smith PetscInt maxits, i, lits; 26485385478SLisandro Dalcin PetscReal rho, fnorm, gnorm, gpnorm, xnorm = 0, delta, nrm, ynorm, norm1; 26579f36460SBarry Smith PetscScalar cnorm; 266df60cc22SBarry Smith KSP ksp; 2673f149594SLisandro Dalcin SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 268df8705c3SBarry Smith PetscBool breakout = PETSC_FALSE; 269df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 2705e28bcb6Sprj- PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 2715e28bcb6Sprj- void *convctx; 2724800dd8cSBarry Smith 2733a40ed3dSBarry Smith PetscFunctionBegin; 2740b121fc5SBarry 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); 275c579b300SPatrick Farrell 276fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 277fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 27839e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 279fbe28522SBarry Smith Y = snes->work[0]; /* work vectors */ 280fbe28522SBarry Smith G = snes->work[1]; 2816b5873e3SBarry Smith Ytmp = snes->work[2]; 282c9368356SGlenn Hammond W = snes->work[3]; 2834800dd8cSBarry Smith 2849566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 2854c49b128SBarry Smith snes->iter = 0; 2869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2874800dd8cSBarry Smith 288df8705c3SBarry Smith /* Set the linear stopping criteria to use the More' trick. */ 2899566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2909566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 291df8705c3SBarry Smith if (convtest != SNESTR_KSPConverged_Private) { 2929566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 293df8705c3SBarry Smith ctx->snes = snes; 2949566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 2959566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 2969566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 297df8705c3SBarry Smith } 298df8705c3SBarry Smith 299e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3009566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 3011aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 302e4ed7901SPeter Brune 3039566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 304422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 3059566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 3069566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 307fbe28522SBarry Smith snes->norm = fnorm; 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3098b84755bSBarry Smith delta = xnorm ? neP->delta0 * xnorm : neP->delta0; 3104800dd8cSBarry Smith neP->delta = delta; 3119566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3129566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 313b37302c6SLois Curfman McInnes 31485385478SLisandro Dalcin /* test convergence */ 315dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 31685385478SLisandro Dalcin if (snes->reason) PetscFunctionReturn(0); 3173f149594SLisandro Dalcin 3184800dd8cSBarry Smith for (i = 0; i < maxits; i++) { 31999a96b7cSMatthew Knepley /* Call general purpose update function */ 320dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 32199a96b7cSMatthew Knepley 32285385478SLisandro Dalcin /* Solve J Y = F, where J is Jacobian matrix */ 3239566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 32407b62357SFande Kong SNESCheckJacobianDomainerror(snes); 3259566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 3269566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Ytmp)); 3279566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 328329f5518SBarry Smith snes->linear_its += lits; 3291aa26658SKarl Rupp 3309566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 3319566063dSJacob Faibussowitsch PetscCall(VecNorm(Ytmp, NORM_2, &nrm)); 332064f8208SBarry Smith norm1 = nrm; 333284fb49fSHeeho Park 3346b5873e3SBarry Smith while (1) { 335c9368356SGlenn Hammond PetscBool changed_y; 3367cb011f5SBarry Smith PetscBool changed_w; 3379566063dSJacob Faibussowitsch PetscCall(VecCopy(Ytmp, Y)); 338064f8208SBarry Smith nrm = norm1; 3396b5873e3SBarry Smith 3402deea200SLois Curfman McInnes /* Scale Y if need be and predict new value of F norm */ 341064f8208SBarry Smith if (nrm >= delta) { 342064f8208SBarry Smith nrm = delta / nrm; 343064f8208SBarry Smith gpnorm = (1.0 - nrm) * fnorm; 344064f8208SBarry Smith cnorm = nrm; 3459566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Scaling direction by %g\n", (double)nrm)); 3469566063dSJacob Faibussowitsch PetscCall(VecScale(Y, cnorm)); 347064f8208SBarry Smith nrm = gpnorm; 348fbe28522SBarry Smith ynorm = delta; 349fbe28522SBarry Smith } else { 350fbe28522SBarry Smith gpnorm = 0.0; 3519566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Direction is in Trust Region\n")); 352064f8208SBarry Smith ynorm = nrm; 353fbe28522SBarry Smith } 354c9368356SGlenn Hammond /* PreCheck() allows for updates to Y prior to W <- X - Y */ 3559566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 3569566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -1.0, Y, X)); /* W <- X - Y */ 3579566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 3589566063dSJacob Faibussowitsch if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X)); 3599566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 3609566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */ 3619566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */ 3627551ef16SBarry Smith SNESCheckFunctionNorm(snes, gnorm); 3634800dd8cSBarry Smith if (fnorm == gpnorm) rho = 0.0; 364fbe28522SBarry Smith else rho = (fnorm * fnorm - gnorm * gnorm) / (fnorm * fnorm - gpnorm * gpnorm); 3654800dd8cSBarry Smith 3664800dd8cSBarry Smith /* Update size of trust region */ 3674800dd8cSBarry Smith if (rho < neP->mu) delta *= neP->delta1; 3684800dd8cSBarry Smith else if (rho < neP->eta) delta *= neP->delta2; 3694800dd8cSBarry Smith else delta *= neP->delta3; 3709566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%g, gnorm=%g, ynorm=%g\n", (double)fnorm, (double)gnorm, (double)ynorm)); 3719566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "gpred=%g, rho=%g, delta=%g\n", (double)gpnorm, (double)rho, (double)delta)); 3721aa26658SKarl Rupp 3734800dd8cSBarry Smith neP->delta = delta; 3744800dd8cSBarry Smith if (rho > neP->sigma) break; 3759566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 376284fb49fSHeeho Park 3776b5873e3SBarry Smith /* check to see if progress is hopeless */ 378ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 3799566063dSJacob Faibussowitsch PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP)); 380dbbe0bcdSBarry Smith if (!reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 3817cb011f5SBarry Smith if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 382184914b5SBarry Smith if (reason) { 38352392280SLois Curfman McInnes /* We're not progressing, so return with the current iterate */ 3849566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 385a7cc72afSBarry Smith breakout = PETSC_TRUE; 386454a90a3SBarry Smith break; 38752392280SLois Curfman McInnes } 38850ffb88aSMatthew Knepley snes->numFailures++; 3894800dd8cSBarry Smith } 3901acabf8cSLois Curfman McInnes if (!breakout) { 39185385478SLisandro Dalcin /* Update function and solution vectors */ 3924800dd8cSBarry Smith fnorm = gnorm; 3939566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 3949566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 39585385478SLisandro Dalcin /* Monitor convergence */ 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 397c509a162SBarry Smith snes->iter = i + 1; 398fbe28522SBarry Smith snes->norm = fnorm; 399c1e67a49SFande Kong snes->xnorm = xnorm; 400c1e67a49SFande Kong snes->ynorm = ynorm; 4019566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 4029566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 4039566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 40485385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 405ef87ac0dSKris Buschelman neP->itflag = PETSC_TRUE; 4069566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 407dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 4083f149594SLisandro Dalcin if (reason) break; 4091aa26658SKarl Rupp } else break; 41038442cffSBarry Smith } 411284fb49fSHeeho Park 41252392280SLois Curfman McInnes if (i == maxits) { 4139566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 41485385478SLisandro Dalcin if (!reason) reason = SNES_DIVERGED_MAX_IT; 41552392280SLois Curfman McInnes } 4169566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 417184914b5SBarry Smith snes->reason = reason; 4189566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 4195e28bcb6Sprj- if (convtest != SNESTR_KSPConverged_Private) { 4209566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 4219566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 4229566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 4235e28bcb6Sprj- } 4243a40ed3dSBarry Smith PetscFunctionReturn(0); 4254800dd8cSBarry Smith } 426284fb49fSHeeho Park 427d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 428d71ae5a4SJacob Faibussowitsch { 4293a40ed3dSBarry Smith PetscFunctionBegin; 4309566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 4)); 4319566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 4323a40ed3dSBarry Smith PetscFunctionReturn(0); 4334800dd8cSBarry Smith } 4346b8b9a38SLisandro Dalcin 435d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONTR(SNES snes) 436d71ae5a4SJacob Faibussowitsch { 4376b8b9a38SLisandro Dalcin PetscFunctionBegin; 4386b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 4396b8b9a38SLisandro Dalcin } 4406b8b9a38SLisandro Dalcin 441d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 442d71ae5a4SJacob Faibussowitsch { 4433a40ed3dSBarry Smith PetscFunctionBegin; 4449566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONTR(snes)); 4459566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4463a40ed3dSBarry Smith PetscFunctionReturn(0); 4474800dd8cSBarry Smith } 4484800dd8cSBarry Smith 449d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 450d71ae5a4SJacob Faibussowitsch { 45104d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 4524800dd8cSBarry Smith 4533a40ed3dSBarry Smith PetscFunctionBegin; 454d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 4559566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trtol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 4569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_mu", "mu", "None", ctx->mu, &ctx->mu, NULL)); 4579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_eta", "eta", "None", ctx->eta, &ctx->eta, NULL)); 4589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_sigma", "sigma", "None", ctx->sigma, &ctx->sigma, NULL)); 4599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 4609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta1", "delta1", "None", ctx->delta1, &ctx->delta1, NULL)); 4619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta2", "delta2", "None", ctx->delta2, &ctx->delta2, NULL)); 4629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta3", "delta3", "None", ctx->delta3, &ctx->delta3, NULL)); 463d0609cedSBarry Smith PetscOptionsHeadEnd(); 4643a40ed3dSBarry Smith PetscFunctionReturn(0); 4654800dd8cSBarry Smith } 4664800dd8cSBarry Smith 467d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 468d71ae5a4SJacob Faibussowitsch { 46904d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 470ace3abfcSBarry Smith PetscBool iascii; 471a935fc98SLois Curfman McInnes 4723a40ed3dSBarry Smith PetscFunctionBegin; 4739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 47432077d6dSBarry Smith if (iascii) { 47563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol)); 4769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " mu=%g, eta=%g, sigma=%g\n", (double)tr->mu, (double)tr->eta, (double)tr->sigma)); 4779566063dSJacob 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)); 47819bcc07fSBarry Smith } 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 480a935fc98SLois Curfman McInnes } 481f6dfbefdSBarry Smith 482ebe3b25bSBarry Smith /*MC 48304d7464bSBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 484ebe3b25bSBarry Smith 485f6dfbefdSBarry Smith Deprecated use `SNESNEWTONTRDC` 486f6dfbefdSBarry Smith 487f6dfbefdSBarry Smith Options Database Keys: 4888e434772SBarry Smith + -snes_trtol <tol> - trust region tolerance 4898e434772SBarry Smith . -snes_tr_mu <mu> - trust region parameter 4908e434772SBarry Smith . -snes_tr_eta <eta> - trust region parameter 4918e434772SBarry Smith . -snes_tr_sigma <sigma> - trust region parameter 4928b84755bSBarry Smith . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 4938e434772SBarry Smith . -snes_tr_delta1 <delta1> - trust region parameter 4948e434772SBarry Smith . -snes_tr_delta2 <delta2> - trust region parameter 4958e434772SBarry Smith - -snes_tr_delta3 <delta3> - trust region parameter 496ebe3b25bSBarry Smith 497f6dfbefdSBarry Smith Reference: 498f6dfbefdSBarry Smith . - * "The Minpack Project", by More', Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 499ebe3b25bSBarry Smith of Mathematical Software", Wayne Cowell, editor. 500ebe3b25bSBarry Smith 501ee3001cbSBarry Smith Level: intermediate 502ee3001cbSBarry Smith 503f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()` 504ebe3b25bSBarry Smith M*/ 505d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 506d71ae5a4SJacob Faibussowitsch { 50704d7464bSBarry Smith SNES_NEWTONTR *neP; 5084800dd8cSBarry Smith 5093a40ed3dSBarry Smith PetscFunctionBegin; 51004d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 51104d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 51204d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 51304d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 51404d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 51504d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONTR; 516fbe28522SBarry Smith 51742f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 518efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 51942f4f86dSBarry Smith 5204fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5214fc747eaSLawrence Mitchell 5224dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 523fbe28522SBarry Smith snes->data = (void *)neP; 524fbe28522SBarry Smith neP->mu = 0.25; 525fbe28522SBarry Smith neP->eta = 0.75; 526fbe28522SBarry Smith neP->delta = 0.0; 527fbe28522SBarry Smith neP->delta0 = 0.2; 528fbe28522SBarry Smith neP->delta1 = 0.3; 529fbe28522SBarry Smith neP->delta2 = 0.75; 530fbe28522SBarry Smith neP->delta3 = 2.0; 531fbe28522SBarry Smith neP->sigma = 0.0001; 532ef87ac0dSKris Buschelman neP->itflag = PETSC_FALSE; 53375567043SBarry Smith neP->rnorm0 = 0.0; 53475567043SBarry Smith neP->ttol = 0.0; 5353a40ed3dSBarry Smith PetscFunctionReturn(0); 5364800dd8cSBarry Smith } 537