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 PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 6df8705c3SBarry Smith PetscErrorCode (*convdestroy)(void *); 7df8705c3SBarry Smith void *convctx; 8971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx; 9971273eeSBarry Smith 104a221d59SStefano Zampini const char *const SNESNewtonTRFallbackTypes[] = {"NEWTON", "CAUCHY", "DOGLEG", "SNESNewtonTRFallbackType", "SNES_TR_FALLBACK_", NULL}; 114a221d59SStefano Zampini 12d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 13d71ae5a4SJacob Faibussowitsch { 14971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 15971273eeSBarry Smith SNES snes = ctx->snes; 1604d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 17df60cc22SBarry Smith Vec x; 18064f8208SBarry Smith PetscReal nrm; 19df60cc22SBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 2248a46eb9SPierre Jolivet if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 23a935fc98SLois Curfman McInnes /* Determine norm of solution */ 249566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 259566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &nrm)); 26064f8208SBarry Smith if (nrm >= neP->delta) { 279566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 28329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 29df60cc22SBarry Smith } 303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31df60cc22SBarry Smith } 3282bf6240SBarry Smith 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 34d71ae5a4SJacob Faibussowitsch { 35971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 36971273eeSBarry Smith 37971273eeSBarry Smith PetscFunctionBegin; 389566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 399566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41971273eeSBarry Smith } 42971273eeSBarry Smith 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 44d71ae5a4SJacob Faibussowitsch { 4504d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 4685385478SLisandro Dalcin 4785385478SLisandro Dalcin PetscFunctionBegin; 4885385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 494a221d59SStefano Zampini if (neP->delta < snes->deltatol) { 504a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)snes->deltatol)); 511c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 52e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 5363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 5485385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 5585385478SLisandro Dalcin } 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5785385478SLisandro Dalcin } 5885385478SLisandro Dalcin 594a221d59SStefano Zampini /*@ 604a221d59SStefano Zampini SNESNewtonTRSetFallbackType - Set the type of fallback if the solution of the trust region subproblem is outside the radius 614a221d59SStefano Zampini 624a221d59SStefano Zampini Input Parameters: 634a221d59SStefano Zampini + snes - the nonlinear solver object 644a221d59SStefano Zampini - ftype - the fallback type, see `SNESNewtonTRFallbackType` 654a221d59SStefano Zampini 664a221d59SStefano Zampini Level: intermediate 674a221d59SStefano Zampini 684a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`, 694a221d59SStefano Zampini `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 704a221d59SStefano Zampini @*/ 714a221d59SStefano Zampini PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype) 724a221d59SStefano Zampini { 734a221d59SStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 744a221d59SStefano Zampini PetscBool flg; 754a221d59SStefano Zampini 764a221d59SStefano Zampini PetscFunctionBegin; 774a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 784a221d59SStefano Zampini PetscValidLogicalCollectiveEnum(snes, ftype, 2); 794a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 804a221d59SStefano Zampini if (flg) tr->fallback = ftype; 814a221d59SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 824a221d59SStefano Zampini } 834a221d59SStefano Zampini 847cb011f5SBarry Smith /*@C 85c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 864a221d59SStefano Zampini Allows the user a chance to change or override the trust region decision. 87f6dfbefdSBarry Smith 88c3339decSBarry Smith Logically Collective 89c9368356SGlenn Hammond 90c9368356SGlenn Hammond Input Parameters: 91c9368356SGlenn Hammond + snes - the nonlinear solver object 92*20f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 93*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 94c9368356SGlenn Hammond 95*20f4b53cSBarry Smith Level: deprecated (since 3.19) 96c9368356SGlenn Hammond 97f6dfbefdSBarry Smith Note: 984a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver. 99c9368356SGlenn Hammond 1004a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 101c9368356SGlenn Hammond @*/ 102d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 103d71ae5a4SJacob Faibussowitsch { 104c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1054a221d59SStefano Zampini PetscBool flg; 106c9368356SGlenn Hammond 107c9368356SGlenn Hammond PetscFunctionBegin; 108c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1094a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1104a221d59SStefano Zampini if (flg) { 111c9368356SGlenn Hammond if (func) tr->precheck = func; 112c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 1134a221d59SStefano Zampini } 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115c9368356SGlenn Hammond } 116c9368356SGlenn Hammond 117c9368356SGlenn Hammond /*@C 118c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 119c9368356SGlenn Hammond 120*20f4b53cSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 121*20f4b53cSBarry Smith 122*20f4b53cSBarry Smith Not Collective 123c9368356SGlenn Hammond 124c9368356SGlenn Hammond Input Parameter: 125c9368356SGlenn Hammond . snes - the nonlinear solver context 126c9368356SGlenn Hammond 127c9368356SGlenn Hammond Output Parameters: 128*20f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 129*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 130c9368356SGlenn Hammond 131*20f4b53cSBarry Smith Level: deprecated (since 3.19) 132c9368356SGlenn Hammond 1334a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 134c9368356SGlenn Hammond @*/ 135d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 136d71ae5a4SJacob Faibussowitsch { 137c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1384a221d59SStefano Zampini PetscBool flg; 139c9368356SGlenn Hammond 140c9368356SGlenn Hammond PetscFunctionBegin; 141c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1424a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1434a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 144c9368356SGlenn Hammond if (func) *func = tr->precheck; 145c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147c9368356SGlenn Hammond } 148c9368356SGlenn Hammond 149c9368356SGlenn Hammond /*@C 1507cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 1514a221d59SStefano Zampini function evaluation. Allows the user a chance to change or override the internal decision of the solver 152f6dfbefdSBarry Smith 153c3339decSBarry Smith Logically Collective 1547cb011f5SBarry Smith 1557cb011f5SBarry Smith Input Parameters: 1567cb011f5SBarry Smith + snes - the nonlinear solver object 157*20f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 158*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 1597cb011f5SBarry Smith 160*20f4b53cSBarry Smith Level: deprecated (since 3.19) 1617cb011f5SBarry Smith 162f6dfbefdSBarry Smith Note: 1634a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver while the function set in 164f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 1657cb011f5SBarry Smith 1664a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 1677cb011f5SBarry Smith @*/ 168d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 169d71ae5a4SJacob Faibussowitsch { 1707cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1714a221d59SStefano Zampini PetscBool flg; 1727cb011f5SBarry Smith 1737cb011f5SBarry Smith PetscFunctionBegin; 1747cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1754a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1764a221d59SStefano Zampini if (flg) { 1777cb011f5SBarry Smith if (func) tr->postcheck = func; 1787cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 1794a221d59SStefano Zampini } 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1817cb011f5SBarry Smith } 1827cb011f5SBarry Smith 1837cb011f5SBarry Smith /*@C 1847cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 1857cb011f5SBarry Smith 186*20f4b53cSBarry Smith Not Collective 1877cb011f5SBarry Smith 1887cb011f5SBarry Smith Input Parameter: 1897cb011f5SBarry Smith . snes - the nonlinear solver context 1907cb011f5SBarry Smith 1917cb011f5SBarry Smith Output Parameters: 192*20f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 193*20f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 1947cb011f5SBarry Smith 1957cb011f5SBarry Smith Level: intermediate 1967cb011f5SBarry Smith 1974a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 1987cb011f5SBarry Smith @*/ 199d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 200d71ae5a4SJacob Faibussowitsch { 2017cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2024a221d59SStefano Zampini PetscBool flg; 2037cb011f5SBarry Smith 2047cb011f5SBarry Smith PetscFunctionBegin; 2057cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2064a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2074a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 2087cb011f5SBarry Smith if (func) *func = tr->postcheck; 2097cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2117cb011f5SBarry Smith } 2127cb011f5SBarry Smith 2137cb011f5SBarry Smith /*@C 2144a221d59SStefano Zampini SNESNewtonTRPreCheck - Runs the precheck routine 215c9368356SGlenn Hammond 216c3339decSBarry Smith Logically Collective 217c9368356SGlenn Hammond 218c9368356SGlenn Hammond Input Parameters: 219c9368356SGlenn Hammond + snes - the solver 220c9368356SGlenn Hammond . X - The last solution 221c9368356SGlenn Hammond - Y - The step direction 222c9368356SGlenn Hammond 223c9368356SGlenn Hammond Output Parameters: 224c9368356SGlenn Hammond . changed_Y - Indicator that the step direction Y has been changed. 225c9368356SGlenn Hammond 2264a221d59SStefano Zampini Level: intermediate 227c9368356SGlenn Hammond 2284a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 229c9368356SGlenn Hammond @*/ 2304a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 231d71ae5a4SJacob Faibussowitsch { 232c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2334a221d59SStefano Zampini PetscBool flg; 234c9368356SGlenn Hammond 235c9368356SGlenn Hammond PetscFunctionBegin; 2364a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2374a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2384a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 239c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 240c9368356SGlenn Hammond if (tr->precheck) { 2419566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 242c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 243c9368356SGlenn Hammond } 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245c9368356SGlenn Hammond } 246c9368356SGlenn Hammond 247c9368356SGlenn Hammond /*@C 2484a221d59SStefano Zampini SNESNewtonTRPostCheck - Runs the postcheck routine 2497cb011f5SBarry Smith 250c3339decSBarry Smith Logically Collective 2517cb011f5SBarry Smith 2527cb011f5SBarry Smith Input Parameters: 2536b867d5aSJose E. Roman + snes - the solver 2546b867d5aSJose E. Roman . X - The last solution 2557cb011f5SBarry Smith . Y - The full step direction 2563312a946SBarry Smith - W - The updated solution, W = X - Y 2577cb011f5SBarry Smith 2587cb011f5SBarry Smith Output Parameters: 2593312a946SBarry Smith + changed_Y - indicator if step has been changed 2603312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 2617cb011f5SBarry Smith 262f6dfbefdSBarry Smith Note: 2633312a946SBarry Smith If Y is changed then W is recomputed as X - Y 2647cb011f5SBarry Smith 2654a221d59SStefano Zampini Level: intermediate 2667cb011f5SBarry Smith 2674a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck() 2687cb011f5SBarry Smith @*/ 2694a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 270d71ae5a4SJacob Faibussowitsch { 2717cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2724a221d59SStefano Zampini PetscBool flg; 2737cb011f5SBarry Smith 2747cb011f5SBarry Smith PetscFunctionBegin; 2754a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2764a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2774a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 278c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 2797cb011f5SBarry Smith *changed_W = PETSC_FALSE; 2807cb011f5SBarry Smith if (tr->postcheck) { 2819566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 282c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 2837cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 2847cb011f5SBarry Smith } 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2867cb011f5SBarry Smith } 28785385478SLisandro Dalcin 2884a221d59SStefano Zampini static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 2894a221d59SStefano Zampini { 2904a221d59SStefano Zampini PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 2914a221d59SStefano Zampini PetscReal x1 = temp / a; 2924a221d59SStefano Zampini PetscReal x2 = c / temp; 2934a221d59SStefano Zampini *xm = PetscMin(x1, x2); 2944a221d59SStefano Zampini *xp = PetscMax(x1, x2); 2954a221d59SStefano Zampini } 2964a221d59SStefano Zampini 297df60cc22SBarry Smith /* 2984a221d59SStefano Zampini SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 2994a221d59SStefano Zampini (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 3004a221d59SStefano Zampini nonlinear equations 3014800dd8cSBarry Smith 3024800dd8cSBarry Smith */ 303d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 304d71ae5a4SJacob Faibussowitsch { 30504d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 3064a221d59SStefano Zampini Vec X, F, Y, G, W, GradF, YU; 3074a221d59SStefano Zampini PetscInt maxits, lits; 3084a221d59SStefano Zampini PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm; 3094a221d59SStefano Zampini PetscReal deltaM, fk, fkp1, deltaqm, gTy, yTHy; 3104a221d59SStefano Zampini PetscReal auk, gfnorm, ycnorm, gTBg; 311df60cc22SBarry Smith KSP ksp; 3124a221d59SStefano Zampini PetscBool already_done = PETSC_FALSE; 3134a221d59SStefano Zampini PetscBool clear_converged_test, rho_satisfied; 3144a221d59SStefano Zampini PetscVoidFunction ksp_has_radius; 315df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 3165e28bcb6Sprj- void *convctx; 3174a221d59SStefano Zampini PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 3184a221d59SStefano Zampini PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 3194800dd8cSBarry Smith 3203a40ed3dSBarry Smith PetscFunctionBegin; 3214a221d59SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 322c579b300SPatrick Farrell 323fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 324fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 32539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 3264a221d59SStefano Zampini Y = snes->vec_sol_update; /* update vector */ 3274a221d59SStefano Zampini G = snes->work[0]; /* updated residual */ 3284a221d59SStefano Zampini W = snes->work[1]; /* temporary vector */ 3294a221d59SStefano Zampini GradF = !objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 3304a221d59SStefano Zampini YU = snes->work[3]; /* work vector for dogleg method */ 3314a221d59SStefano Zampini 3324a221d59SStefano Zampini 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); 3334800dd8cSBarry Smith 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 3354c49b128SBarry Smith snes->iter = 0; 3369566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3374800dd8cSBarry Smith 3384a221d59SStefano Zampini /* Set the linear stopping criteria to use the More' trick if needed */ 3394a221d59SStefano Zampini clear_converged_test = PETSC_FALSE; 3409566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3419566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 3424a221d59SStefano Zampini PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &ksp_has_radius)); 3434a221d59SStefano Zampini if (convtest != SNESTR_KSPConverged_Private && !ksp_has_radius) { 3444a221d59SStefano Zampini clear_converged_test = PETSC_TRUE; 3459566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 346df8705c3SBarry Smith ctx->snes = snes; 3479566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 3489566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 3499566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 350df8705c3SBarry Smith } 351df8705c3SBarry Smith 352e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3539566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 3541aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 355e4ed7901SPeter Brune 3569566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 357422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 3589566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 3594a221d59SStefano Zampini 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 361fbe28522SBarry Smith snes->norm = fnorm; 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3634a221d59SStefano Zampini delta = neP->delta0; 3644a221d59SStefano Zampini deltaM = neP->deltaM; 3654800dd8cSBarry Smith neP->delta = delta; 3669566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3679566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 368b37302c6SLois Curfman McInnes 36985385478SLisandro Dalcin /* test convergence */ 3704a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 371dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 3723ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 3733f149594SLisandro Dalcin 3744a221d59SStefano Zampini if (objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 3754a221d59SStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 37699a96b7cSMatthew Knepley 3774a221d59SStefano Zampini while (snes->iter < maxits) { 378c9368356SGlenn Hammond PetscBool changed_y; 3797cb011f5SBarry Smith PetscBool changed_w; 3806b5873e3SBarry Smith 3814a221d59SStefano Zampini /* solve trust-region subproblem */ 3824a221d59SStefano Zampini if (!already_done) { 3834a221d59SStefano Zampini PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 3844a221d59SStefano Zampini SNESCheckJacobianDomainerror(snes); 385fbe28522SBarry Smith } 3864a221d59SStefano Zampini PetscCall(KSPCGSetRadius(snes->ksp, delta)); 3874a221d59SStefano Zampini PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 3884a221d59SStefano Zampini PetscCall(KSPSolve(snes->ksp, F, Y)); 3894a221d59SStefano Zampini SNESCheckKSPSolve(snes); 3904a221d59SStefano Zampini PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 3914800dd8cSBarry Smith 3924a221d59SStefano Zampini /* calculating GradF of minimization function only once */ 3934a221d59SStefano Zampini if (!already_done) { 3944a221d59SStefano Zampini if (objective) gfnorm = fnorm; 3954a221d59SStefano Zampini else { 3964a221d59SStefano Zampini PetscCall(MatMultTranspose(snes->jacobian, F, GradF)); /* grad f = J^T F */ 3974a221d59SStefano Zampini PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 3984a221d59SStefano Zampini } 3994a221d59SStefano Zampini already_done = PETSC_TRUE; 4004a221d59SStefano Zampini } 4011aa26658SKarl Rupp 4024a221d59SStefano Zampini /* decide what to do when the update is outside of trust region */ 4034a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 4044a221d59SStefano Zampini if (ynorm > delta) { 4054a221d59SStefano Zampini switch (neP->fallback) { 4064a221d59SStefano Zampini case SNES_TR_FALLBACK_NEWTON: 4074a221d59SStefano Zampini auk = delta / ynorm; 4084a221d59SStefano Zampini PetscCall(VecScale(Y, auk)); 4094a221d59SStefano Zampini break; 4104a221d59SStefano Zampini case SNES_TR_FALLBACK_CAUCHY: 4114a221d59SStefano Zampini case SNES_TR_FALLBACK_DOGLEG: 4124a221d59SStefano Zampini PetscCall(MatMult(snes->jacobian, GradF, W)); 4134a221d59SStefano Zampini if (objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 4144a221d59SStefano Zampini else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 4154a221d59SStefano Zampini /* Eqs 4.7 and 4.8 in Nocedal and Wright */ 4164a221d59SStefano Zampini auk = delta / gfnorm; 4174a221d59SStefano Zampini if (gTBg > 0.0) auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1); 4184a221d59SStefano Zampini ycnorm = auk * gfnorm; 4194a221d59SStefano Zampini if (neP->fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 4204a221d59SStefano Zampini /* Cauchy solution */ 4214a221d59SStefano Zampini PetscCall(VecAXPBY(Y, auk, 0.0, GradF)); 4224a221d59SStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 4234a221d59SStefano Zampini } else { /* take linear combination of Cauchy and Newton direction and step */ 4244a221d59SStefano Zampini PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 4254a221d59SStefano Zampini PetscBool noroots; 426284fb49fSHeeho Park 4274a221d59SStefano Zampini auk = gfnorm * gfnorm / gTBg; 4284a221d59SStefano Zampini PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 4294a221d59SStefano Zampini PetscCall(VecAXPY(Y, -1.0, YU)); 4304a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); 4314a221d59SStefano Zampini PetscCall(VecDotRealPart(YU, Y, &c1)); 4324a221d59SStefano Zampini c0 = PetscSqr(c0); 4334a221d59SStefano Zampini c2 = PetscSqr(ycnorm) - PetscSqr(delta); 4344a221d59SStefano Zampini PetscQuadraticRoots(c0, c1, c2, &tneg, &tpos); 4354a221d59SStefano Zampini 4364a221d59SStefano Zampini noroots = PetscIsInfOrNanReal(tneg); 4374a221d59SStefano Zampini if (noroots) { /* No roots, select Cauchy point */ 4384a221d59SStefano Zampini auk = delta / gfnorm; 4394a221d59SStefano Zampini auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1); 4404a221d59SStefano Zampini PetscCall(VecAXPBY(Y, auk, 0.0, GradF)); 4414a221d59SStefano Zampini } else { /* Here roots corresponds to tau-1 in Nocedal and Wright */ 4424a221d59SStefano Zampini tpos += 1.0; 4434a221d59SStefano Zampini tneg += 1.0; 4444a221d59SStefano Zampini tau = PetscClipInterval(tpos, 0.0, 2.0); /* clip to tau [0,2] */ 4454a221d59SStefano Zampini if (tau < 1.0) { 4464a221d59SStefano Zampini PetscCall(VecAXPBY(Y, tau, 0.0, YU)); 4474a221d59SStefano Zampini } else { 4484a221d59SStefano Zampini PetscCall(VecAXPBY(Y, 1.0, tau - 1, YU)); 4494a221d59SStefano Zampini } 4504a221d59SStefano Zampini } 4514a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); /* this norm will be cached and reused later */ 4524a221d59SStefano Zampini PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, ydlnorm %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)c0, (double)gTBg)); 4534a221d59SStefano Zampini } 4544a221d59SStefano Zampini break; 4554a221d59SStefano Zampini default: 4564a221d59SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 457454a90a3SBarry Smith break; 45852392280SLois Curfman McInnes } 4594800dd8cSBarry Smith } 4604a221d59SStefano Zampini 4614a221d59SStefano Zampini /* Evaluate the solution to meet the improvement ratio criteria */ 4624a221d59SStefano Zampini 4634a221d59SStefano Zampini /* compute the final ynorm */ 4644a221d59SStefano Zampini PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 4654a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 4664a221d59SStefano Zampini 4674a221d59SStefano Zampini /* the quadratic model difference */ 4684a221d59SStefano Zampini PetscCall(MatMult(snes->jacobian, Y, W)); 4694a221d59SStefano Zampini if (objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 4704a221d59SStefano Zampini else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */ 4714a221d59SStefano Zampini PetscCall(VecDotRealPart(GradF, Y, &gTy)); 4724a221d59SStefano Zampini deltaqm = -(-gTy + 0.5 * yTHy); /* difference in quadratic model, -gTy because SNES solves it this way */ 4734a221d59SStefano Zampini 4744a221d59SStefano Zampini /* update */ 4754a221d59SStefano Zampini PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 4764a221d59SStefano Zampini PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 4774a221d59SStefano Zampini if (changed_y) { 4784a221d59SStefano Zampini /* Need to recompute the quadratic model difference */ 4794a221d59SStefano Zampini PetscCall(MatMult(snes->jacobian, Y, W)); 4804a221d59SStefano Zampini if (objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 4814a221d59SStefano Zampini else PetscCall(VecDotRealPart(W, W, &yTHy)); 4824a221d59SStefano Zampini PetscCall(VecDotRealPart(GradF, Y, &gTy)); 4834a221d59SStefano Zampini deltaqm = -(-gTy + 0.5 * yTHy); 4844a221d59SStefano Zampini /* User changed Y but not W */ 4854a221d59SStefano Zampini if (!changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 4864a221d59SStefano Zampini } 4874a221d59SStefano Zampini 4884a221d59SStefano Zampini /* Compute new objective function */ 4894a221d59SStefano Zampini PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 4904a221d59SStefano Zampini PetscCall(VecNorm(G, NORM_2, &gnorm)); 4914a221d59SStefano Zampini if (objective) PetscCall(SNESComputeObjective(snes, W, &fkp1)); 4924a221d59SStefano Zampini else fkp1 = 0.5 * PetscSqr(gnorm); 4934a221d59SStefano Zampini SNESCheckFunctionNorm(snes, fkp1); 4944a221d59SStefano Zampini 4954a221d59SStefano Zampini if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 4964a221d59SStefano Zampini else rho = -1.0; /* no reduction in quadratic model, step must be rejected */ 4974a221d59SStefano Zampini PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy)); 4984a221d59SStefano Zampini 4994a221d59SStefano Zampini if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 5004a221d59SStefano Zampini else if (rho > neP->eta3) delta *= neP->t2; /* expand the region */ 5014a221d59SStefano Zampini delta = PetscMin(delta, deltaM); /* but not greater than deltaM */ 5024a221d59SStefano Zampini 5034a221d59SStefano Zampini neP->delta = delta; 5044a221d59SStefano Zampini if (rho >= neP->eta1) { 5054a221d59SStefano Zampini rho_satisfied = PETSC_TRUE; 5064a221d59SStefano Zampini } else { 5074a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 5084a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 5094a221d59SStefano Zampini /* check to see if progress is hopeless */ 5104a221d59SStefano Zampini PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 5114a221d59SStefano Zampini if (!snes->reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 5124a221d59SStefano Zampini if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_INNER; 5134a221d59SStefano Zampini snes->numFailures++; 5144a221d59SStefano Zampini /* We're not progressing, so return with the current iterate */ 5154a221d59SStefano Zampini if (snes->reason) break; 5164a221d59SStefano Zampini } 5174a221d59SStefano Zampini if (rho_satisfied) { 5184a221d59SStefano Zampini /* Update function values */ 5194a221d59SStefano Zampini already_done = PETSC_FALSE; 5204800dd8cSBarry Smith fnorm = gnorm; 5214a221d59SStefano Zampini fk = fkp1; 5224a221d59SStefano Zampini 5234a221d59SStefano Zampini /* New residual and linearization point */ 5249566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 5259566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 5264a221d59SStefano Zampini 52785385478SLisandro Dalcin /* Monitor convergence */ 5289566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 5294a221d59SStefano Zampini snes->iter++; 530fbe28522SBarry Smith snes->norm = fnorm; 531c1e67a49SFande Kong snes->xnorm = xnorm; 532c1e67a49SFande Kong snes->ynorm = ynorm; 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5349566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 5359566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 5364a221d59SStefano Zampini 53785385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 5384a221d59SStefano Zampini PetscCall(VecNorm(X, NORM_2, &xnorm)); 5394a221d59SStefano Zampini PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 5404a221d59SStefano Zampini if (snes->reason) break; 5414a221d59SStefano Zampini } 54238442cffSBarry Smith } 543284fb49fSHeeho Park 5444a221d59SStefano Zampini if (snes->iter == maxits) { 5459566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 5464a221d59SStefano Zampini if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 54752392280SLois Curfman McInnes } 5484a221d59SStefano Zampini if (clear_converged_test) { 5499566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 5509566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 5519566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 5525e28bcb6Sprj- } 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5544800dd8cSBarry Smith } 555284fb49fSHeeho Park 556d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 557d71ae5a4SJacob Faibussowitsch { 5583a40ed3dSBarry Smith PetscFunctionBegin; 5599566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 4)); 5609566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 5613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5624800dd8cSBarry Smith } 5636b8b9a38SLisandro Dalcin 564d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 565d71ae5a4SJacob Faibussowitsch { 5663a40ed3dSBarry Smith PetscFunctionBegin; 5679566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5694800dd8cSBarry Smith } 5704800dd8cSBarry Smith 571d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 572d71ae5a4SJacob Faibussowitsch { 57304d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 5744800dd8cSBarry Smith 5753a40ed3dSBarry Smith PetscFunctionBegin; 576d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 5774a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 5784a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 5794a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 5804a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 5814a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 5824a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 5834a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 5849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 5854a221d59SStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)ctx->fallback, (PetscEnum *)&ctx->fallback, NULL)); 586d0609cedSBarry Smith PetscOptionsHeadEnd(); 5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5884800dd8cSBarry Smith } 5894800dd8cSBarry Smith 590d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 591d71ae5a4SJacob Faibussowitsch { 59204d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 593ace3abfcSBarry Smith PetscBool iascii; 594a935fc98SLois Curfman McInnes 5953a40ed3dSBarry Smith PetscFunctionBegin; 5969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 59732077d6dSBarry Smith if (iascii) { 5984a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)snes->deltatol)); 5994a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 6004a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM)); 6014a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 60219bcc07fSBarry Smith } 6033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 604a935fc98SLois Curfman McInnes } 605f6dfbefdSBarry Smith 606ebe3b25bSBarry Smith /*MC 6074a221d59SStefano Zampini SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 608f6dfbefdSBarry Smith 609f6dfbefdSBarry Smith Options Database Keys: 6104a221d59SStefano Zampini + -snes_tr_tol <tol> - trust region tolerance 6114a221d59SStefano Zampini . -snes_tr_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 6124a221d59SStefano Zampini . -snes_tr_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 6134a221d59SStefano Zampini . -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 6144a221d59SStefano Zampini . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 6154a221d59SStefano Zampini . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 6164a221d59SStefano Zampini . -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL) 6174a221d59SStefano Zampini . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 6184a221d59SStefano Zampini - -snes_tr_fallback_type <newton,cauchy,dogleg> - Solution strategy to test reduction when step is outside of trust region. Can use scaled Newton direction, Cauchy point (Steepest Descent direction) or dogleg method. 619ebe3b25bSBarry Smith 620f6dfbefdSBarry Smith Reference: 6214a221d59SStefano Zampini . * - "Numerical Optimization" by Nocedal and Wright, chapter 4. 622ebe3b25bSBarry Smith 623*20f4b53cSBarry Smith Level: deprecated (since 3.19) 624ee3001cbSBarry Smith 6254a221d59SStefano Zampini .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 6264a221d59SStefano Zampini `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 6274a221d59SStefano Zampini `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()` 628ebe3b25bSBarry Smith M*/ 629d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 630d71ae5a4SJacob Faibussowitsch { 63104d7464bSBarry Smith SNES_NEWTONTR *neP; 6324800dd8cSBarry Smith 6333a40ed3dSBarry Smith PetscFunctionBegin; 63404d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 63504d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 63604d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 63704d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 63804d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 639fbe28522SBarry Smith 64042f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 641efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 64242f4f86dSBarry Smith 6434fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 6444fc747eaSLawrence Mitchell 6454dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 646fbe28522SBarry Smith snes->data = (void *)neP; 647fbe28522SBarry Smith neP->delta = 0.0; 648fbe28522SBarry Smith neP->delta0 = 0.2; 6494a221d59SStefano Zampini neP->eta1 = 0.001; 6504a221d59SStefano Zampini neP->eta2 = 0.25; 6514a221d59SStefano Zampini neP->eta3 = 0.75; 6524a221d59SStefano Zampini neP->t1 = 0.25; 6534a221d59SStefano Zampini neP->t2 = 2.0; 6544a221d59SStefano Zampini neP->deltaM = 1.e10; 6554a221d59SStefano Zampini neP->fallback = SNES_TR_FALLBACK_NEWTON; 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6574800dd8cSBarry Smith } 658