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 92f6dfbefdSBarry Smith . func - [optional] function evaluation routine, see `SNESNewtonTRPreCheck()` for the calling sequence 93c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 94c9368356SGlenn Hammond 95c9368356SGlenn Hammond Level: intermediate 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 120c9368356SGlenn Hammond Not collective 121c9368356SGlenn Hammond 122c9368356SGlenn Hammond Input Parameter: 123c9368356SGlenn Hammond . snes - the nonlinear solver context 124c9368356SGlenn Hammond 125c9368356SGlenn Hammond Output Parameters: 126f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPreCheck()` 127c9368356SGlenn Hammond - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 128c9368356SGlenn Hammond 129c9368356SGlenn Hammond Level: intermediate 130c9368356SGlenn Hammond 1314a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 132c9368356SGlenn Hammond @*/ 133d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 134d71ae5a4SJacob Faibussowitsch { 135c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1364a221d59SStefano Zampini PetscBool flg; 137c9368356SGlenn Hammond 138c9368356SGlenn Hammond PetscFunctionBegin; 139c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1404a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1414a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 142c9368356SGlenn Hammond if (func) *func = tr->precheck; 143c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145c9368356SGlenn Hammond } 146c9368356SGlenn Hammond 147c9368356SGlenn Hammond /*@C 1487cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 1494a221d59SStefano Zampini function evaluation. Allows the user a chance to change or override the internal decision of the solver 150f6dfbefdSBarry Smith 151c3339decSBarry Smith Logically Collective 1527cb011f5SBarry Smith 1537cb011f5SBarry Smith Input Parameters: 1547cb011f5SBarry Smith + snes - the nonlinear solver object 155f6dfbefdSBarry Smith . func - [optional] function evaluation routine, see `SNESNewtonTRPostCheck()` for the calling sequence 1567cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1577cb011f5SBarry Smith 1587cb011f5SBarry Smith Level: intermediate 1597cb011f5SBarry Smith 160f6dfbefdSBarry Smith Note: 1614a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver while the function set in 162f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 1637cb011f5SBarry Smith 1644a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 1657cb011f5SBarry Smith @*/ 166d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 167d71ae5a4SJacob Faibussowitsch { 1687cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1694a221d59SStefano Zampini PetscBool flg; 1707cb011f5SBarry Smith 1717cb011f5SBarry Smith PetscFunctionBegin; 1727cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1734a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1744a221d59SStefano Zampini if (flg) { 1757cb011f5SBarry Smith if (func) tr->postcheck = func; 1767cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 1774a221d59SStefano Zampini } 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1797cb011f5SBarry Smith } 1807cb011f5SBarry Smith 1817cb011f5SBarry Smith /*@C 1827cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 1837cb011f5SBarry Smith 1847cb011f5SBarry Smith Not collective 1857cb011f5SBarry Smith 1867cb011f5SBarry Smith Input Parameter: 1877cb011f5SBarry Smith . snes - the nonlinear solver context 1887cb011f5SBarry Smith 1897cb011f5SBarry Smith Output Parameters: 190f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPostCheck()` 1917cb011f5SBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 1927cb011f5SBarry Smith 1937cb011f5SBarry Smith Level: intermediate 1947cb011f5SBarry Smith 1954a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 1967cb011f5SBarry Smith @*/ 197d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 198d71ae5a4SJacob Faibussowitsch { 1997cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2004a221d59SStefano Zampini PetscBool flg; 2017cb011f5SBarry Smith 2027cb011f5SBarry Smith PetscFunctionBegin; 2037cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2044a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2054a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 2067cb011f5SBarry Smith if (func) *func = tr->postcheck; 2077cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2097cb011f5SBarry Smith } 2107cb011f5SBarry Smith 2117cb011f5SBarry Smith /*@C 2124a221d59SStefano Zampini SNESNewtonTRPreCheck - Runs the precheck routine 213c9368356SGlenn Hammond 214c3339decSBarry Smith Logically Collective 215c9368356SGlenn Hammond 216c9368356SGlenn Hammond Input Parameters: 217c9368356SGlenn Hammond + snes - the solver 218c9368356SGlenn Hammond . X - The last solution 219c9368356SGlenn Hammond - Y - The step direction 220c9368356SGlenn Hammond 221c9368356SGlenn Hammond Output Parameters: 222c9368356SGlenn Hammond . changed_Y - Indicator that the step direction Y has been changed. 223c9368356SGlenn Hammond 2244a221d59SStefano Zampini Level: intermediate 225c9368356SGlenn Hammond 2264a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 227c9368356SGlenn Hammond @*/ 2284a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 229d71ae5a4SJacob Faibussowitsch { 230c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2314a221d59SStefano Zampini PetscBool flg; 232c9368356SGlenn Hammond 233c9368356SGlenn Hammond PetscFunctionBegin; 2344a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2354a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2364a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 237c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 238c9368356SGlenn Hammond if (tr->precheck) { 2399566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 240c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 241c9368356SGlenn Hammond } 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 243c9368356SGlenn Hammond } 244c9368356SGlenn Hammond 245c9368356SGlenn Hammond /*@C 2464a221d59SStefano Zampini SNESNewtonTRPostCheck - Runs the postcheck routine 2477cb011f5SBarry Smith 248c3339decSBarry Smith Logically Collective 2497cb011f5SBarry Smith 2507cb011f5SBarry Smith Input Parameters: 2516b867d5aSJose E. Roman + snes - the solver 2526b867d5aSJose E. Roman . X - The last solution 2537cb011f5SBarry Smith . Y - The full step direction 2543312a946SBarry Smith - W - The updated solution, W = X - Y 2557cb011f5SBarry Smith 2567cb011f5SBarry Smith Output Parameters: 2573312a946SBarry Smith + changed_Y - indicator if step has been changed 2583312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 2597cb011f5SBarry Smith 260f6dfbefdSBarry Smith Note: 2613312a946SBarry Smith If Y is changed then W is recomputed as X - Y 2627cb011f5SBarry Smith 2634a221d59SStefano Zampini Level: intermediate 2647cb011f5SBarry Smith 2654a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck() 2667cb011f5SBarry Smith @*/ 2674a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 268d71ae5a4SJacob Faibussowitsch { 2697cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2704a221d59SStefano Zampini PetscBool flg; 2717cb011f5SBarry Smith 2727cb011f5SBarry Smith PetscFunctionBegin; 2734a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2744a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2754a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 276c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 2777cb011f5SBarry Smith *changed_W = PETSC_FALSE; 2787cb011f5SBarry Smith if (tr->postcheck) { 2799566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 280c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 2817cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 2827cb011f5SBarry Smith } 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2847cb011f5SBarry Smith } 28585385478SLisandro Dalcin 2864a221d59SStefano Zampini static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 2874a221d59SStefano Zampini { 2884a221d59SStefano Zampini PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 2894a221d59SStefano Zampini PetscReal x1 = temp / a; 2904a221d59SStefano Zampini PetscReal x2 = c / temp; 2914a221d59SStefano Zampini *xm = PetscMin(x1, x2); 2924a221d59SStefano Zampini *xp = PetscMax(x1, x2); 2934a221d59SStefano Zampini } 2944a221d59SStefano Zampini 295*7aa289d8SStefano Zampini /* Computes the quadratic model difference */ 296*7aa289d8SStefano Zampini static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy, PetscReal *gTy, PetscReal *deltaqm) 297*7aa289d8SStefano Zampini { 298*7aa289d8SStefano Zampini PetscFunctionBegin; 299*7aa289d8SStefano Zampini PetscCall(MatMult(snes->jacobian, Y, W)); 300*7aa289d8SStefano Zampini if (has_objective) PetscCall(VecDotRealPart(Y, W, yTHy)); 301*7aa289d8SStefano Zampini else PetscCall(VecDotRealPart(W, W, yTHy)); /* Gauss-Newton approximation J^t * J */ 302*7aa289d8SStefano Zampini PetscCall(VecDotRealPart(GradF, Y, gTy)); 303*7aa289d8SStefano Zampini *deltaqm = -(-(*gTy) + 0.5 * (*yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */ 304*7aa289d8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 305*7aa289d8SStefano Zampini } 306*7aa289d8SStefano Zampini 307df60cc22SBarry Smith /* 3084a221d59SStefano Zampini SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 3094a221d59SStefano Zampini (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 3104a221d59SStefano Zampini nonlinear equations 3114800dd8cSBarry Smith 3124800dd8cSBarry Smith */ 313d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 314d71ae5a4SJacob Faibussowitsch { 31504d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 3164a221d59SStefano Zampini Vec X, F, Y, G, W, GradF, YU; 3174a221d59SStefano Zampini PetscInt maxits, lits; 3184a221d59SStefano Zampini PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm; 3194a221d59SStefano Zampini PetscReal deltaM, fk, fkp1, deltaqm, gTy, yTHy; 3204a221d59SStefano Zampini PetscReal auk, gfnorm, ycnorm, gTBg; 321df60cc22SBarry Smith KSP ksp; 3224a221d59SStefano Zampini PetscBool already_done = PETSC_FALSE; 323*7aa289d8SStefano Zampini PetscBool clear_converged_test, rho_satisfied, has_objective; 3244a221d59SStefano Zampini PetscVoidFunction ksp_has_radius; 325df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 3265e28bcb6Sprj- void *convctx; 3274a221d59SStefano Zampini PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 3284a221d59SStefano Zampini PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 3294800dd8cSBarry Smith 3303a40ed3dSBarry Smith PetscFunctionBegin; 3314a221d59SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 332*7aa289d8SStefano Zampini has_objective = objective ? PETSC_TRUE : PETSC_FALSE; 333c579b300SPatrick Farrell 334fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 335fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 33639e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 3374a221d59SStefano Zampini Y = snes->vec_sol_update; /* update vector */ 3384a221d59SStefano Zampini G = snes->work[0]; /* updated residual */ 3394a221d59SStefano Zampini W = snes->work[1]; /* temporary vector */ 340*7aa289d8SStefano Zampini GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 3414a221d59SStefano Zampini YU = snes->work[3]; /* work vector for dogleg method */ 3424a221d59SStefano Zampini 3434a221d59SStefano 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); 3444800dd8cSBarry Smith 3459566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 3464c49b128SBarry Smith snes->iter = 0; 3479566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3484800dd8cSBarry Smith 3494a221d59SStefano Zampini /* Set the linear stopping criteria to use the More' trick if needed */ 3504a221d59SStefano Zampini clear_converged_test = PETSC_FALSE; 3519566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3529566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 3534a221d59SStefano Zampini PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &ksp_has_radius)); 3544a221d59SStefano Zampini if (convtest != SNESTR_KSPConverged_Private && !ksp_has_radius) { 3554a221d59SStefano Zampini clear_converged_test = PETSC_TRUE; 3569566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 357df8705c3SBarry Smith ctx->snes = snes; 3589566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 3599566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 3609566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 361df8705c3SBarry Smith } 362df8705c3SBarry Smith 363e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3649566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 3651aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 366e4ed7901SPeter Brune 3679566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 368422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 3699566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 3704a221d59SStefano Zampini 3719566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 372fbe28522SBarry Smith snes->norm = fnorm; 3739566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3744a221d59SStefano Zampini delta = neP->delta0; 3754a221d59SStefano Zampini deltaM = neP->deltaM; 3764800dd8cSBarry Smith neP->delta = delta; 3779566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3789566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 379b37302c6SLois Curfman McInnes 38085385478SLisandro Dalcin /* test convergence */ 3814a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 382dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 3833ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 3843f149594SLisandro Dalcin 385*7aa289d8SStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 3864a221d59SStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 38799a96b7cSMatthew Knepley 3884a221d59SStefano Zampini while (snes->iter < maxits) { 389c9368356SGlenn Hammond PetscBool changed_y; 3907cb011f5SBarry Smith PetscBool changed_w; 3916b5873e3SBarry Smith 392*7aa289d8SStefano Zampini /* calculating GradF of minimization function only once */ 3934a221d59SStefano Zampini if (!already_done) { 3944a221d59SStefano Zampini PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 3954a221d59SStefano Zampini SNESCheckJacobianDomainerror(snes); 396*7aa289d8SStefano Zampini if (has_objective) gfnorm = fnorm; 397*7aa289d8SStefano Zampini else { 398*7aa289d8SStefano Zampini PetscCall(MatMultTranspose(snes->jacobian, F, GradF)); /* grad f = J^T F */ 399*7aa289d8SStefano Zampini PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 400fbe28522SBarry Smith } 401*7aa289d8SStefano Zampini } 402*7aa289d8SStefano Zampini already_done = PETSC_TRUE; 403*7aa289d8SStefano Zampini 404*7aa289d8SStefano Zampini /* solve trust-region subproblem */ 4054a221d59SStefano Zampini PetscCall(KSPCGSetRadius(snes->ksp, delta)); 4064a221d59SStefano Zampini PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 4074a221d59SStefano Zampini PetscCall(KSPSolve(snes->ksp, F, Y)); 4084a221d59SStefano Zampini SNESCheckKSPSolve(snes); 4094a221d59SStefano Zampini PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 4104800dd8cSBarry Smith 4114a221d59SStefano Zampini /* decide what to do when the update is outside of trust region */ 4124a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 4134a221d59SStefano Zampini if (ynorm > delta) { 4144a221d59SStefano Zampini switch (neP->fallback) { 4154a221d59SStefano Zampini case SNES_TR_FALLBACK_NEWTON: 4164a221d59SStefano Zampini auk = delta / ynorm; 4174a221d59SStefano Zampini PetscCall(VecScale(Y, auk)); 4184a221d59SStefano Zampini break; 4194a221d59SStefano Zampini case SNES_TR_FALLBACK_CAUCHY: 4204a221d59SStefano Zampini case SNES_TR_FALLBACK_DOGLEG: 4214a221d59SStefano Zampini PetscCall(MatMult(snes->jacobian, GradF, W)); 422*7aa289d8SStefano Zampini if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 4234a221d59SStefano Zampini else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 4244a221d59SStefano Zampini /* Eqs 4.7 and 4.8 in Nocedal and Wright */ 4254a221d59SStefano Zampini auk = delta / gfnorm; 4264a221d59SStefano Zampini if (gTBg > 0.0) auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1); 4274a221d59SStefano Zampini ycnorm = auk * gfnorm; 4284a221d59SStefano Zampini if (neP->fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 4294a221d59SStefano Zampini /* Cauchy solution */ 4304a221d59SStefano Zampini PetscCall(VecAXPBY(Y, auk, 0.0, GradF)); 4314a221d59SStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 4324a221d59SStefano Zampini } else { /* take linear combination of Cauchy and Newton direction and step */ 4334a221d59SStefano Zampini PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 4344a221d59SStefano Zampini PetscBool noroots; 435284fb49fSHeeho Park 4364a221d59SStefano Zampini auk = gfnorm * gfnorm / gTBg; 4374a221d59SStefano Zampini PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 4384a221d59SStefano Zampini PetscCall(VecAXPY(Y, -1.0, YU)); 4394a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); 4404a221d59SStefano Zampini PetscCall(VecDotRealPart(YU, Y, &c1)); 4414a221d59SStefano Zampini c0 = PetscSqr(c0); 4424a221d59SStefano Zampini c2 = PetscSqr(ycnorm) - PetscSqr(delta); 4434a221d59SStefano Zampini PetscQuadraticRoots(c0, c1, c2, &tneg, &tpos); 4444a221d59SStefano Zampini 4454a221d59SStefano Zampini noroots = PetscIsInfOrNanReal(tneg); 4464a221d59SStefano Zampini if (noroots) { /* No roots, select Cauchy point */ 4474a221d59SStefano Zampini auk = delta / gfnorm; 4484a221d59SStefano Zampini auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1); 4494a221d59SStefano Zampini PetscCall(VecAXPBY(Y, auk, 0.0, GradF)); 4504a221d59SStefano Zampini } else { /* Here roots corresponds to tau-1 in Nocedal and Wright */ 4514a221d59SStefano Zampini tpos += 1.0; 4524a221d59SStefano Zampini tneg += 1.0; 4534a221d59SStefano Zampini tau = PetscClipInterval(tpos, 0.0, 2.0); /* clip to tau [0,2] */ 4544a221d59SStefano Zampini if (tau < 1.0) { 4554a221d59SStefano Zampini PetscCall(VecAXPBY(Y, tau, 0.0, YU)); 4564a221d59SStefano Zampini } else { 4574a221d59SStefano Zampini PetscCall(VecAXPBY(Y, 1.0, tau - 1, YU)); 4584a221d59SStefano Zampini } 4594a221d59SStefano Zampini } 4604a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); /* this norm will be cached and reused later */ 4614a221d59SStefano 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)); 4624a221d59SStefano Zampini } 4634a221d59SStefano Zampini break; 4644a221d59SStefano Zampini default: 4654a221d59SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 466454a90a3SBarry Smith break; 46752392280SLois Curfman McInnes } 4684800dd8cSBarry Smith } 4694a221d59SStefano Zampini 4704a221d59SStefano Zampini /* Evaluate the solution to meet the improvement ratio criteria */ 4714a221d59SStefano Zampini 4724a221d59SStefano Zampini /* compute the final ynorm */ 4734a221d59SStefano Zampini PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 4744a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 4754a221d59SStefano Zampini 476*7aa289d8SStefano Zampini /* compute the quadratic model difference */ 477*7aa289d8SStefano Zampini PetscCall(SNESNewtonTRQuadraticDelta(snes, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm)); 4784a221d59SStefano Zampini 4794a221d59SStefano Zampini /* update */ 4804a221d59SStefano Zampini PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 4814a221d59SStefano Zampini PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 4824a221d59SStefano Zampini if (changed_y) { 4834a221d59SStefano Zampini /* Need to recompute the quadratic model difference */ 484*7aa289d8SStefano Zampini PetscCall(SNESNewtonTRQuadraticDelta(snes, has_objective, Y, GradF, YU, &yTHy, &gTy, &deltaqm)); 4854a221d59SStefano Zampini /* User changed Y but not W */ 4864a221d59SStefano Zampini if (!changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 4874a221d59SStefano Zampini } 4884a221d59SStefano Zampini 4894a221d59SStefano Zampini /* Compute new objective function */ 4904a221d59SStefano Zampini PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 4914a221d59SStefano Zampini PetscCall(VecNorm(G, NORM_2, &gnorm)); 492*7aa289d8SStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, W, &fkp1)); 4934a221d59SStefano Zampini else fkp1 = 0.5 * PetscSqr(gnorm); 4944a221d59SStefano Zampini SNESCheckFunctionNorm(snes, fkp1); 4954a221d59SStefano Zampini 4964a221d59SStefano Zampini if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 4974a221d59SStefano Zampini else rho = -1.0; /* no reduction in quadratic model, step must be rejected */ 4984a221d59SStefano 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)); 4994a221d59SStefano Zampini 5004a221d59SStefano Zampini if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 5014a221d59SStefano Zampini else if (rho > neP->eta3) delta *= neP->t2; /* expand the region */ 5024a221d59SStefano Zampini delta = PetscMin(delta, deltaM); /* but not greater than deltaM */ 5034a221d59SStefano Zampini 5044a221d59SStefano Zampini neP->delta = delta; 5054a221d59SStefano Zampini if (rho >= neP->eta1) { 5064a221d59SStefano Zampini rho_satisfied = PETSC_TRUE; 5074a221d59SStefano Zampini } else { 5084a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 5094a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 5104a221d59SStefano Zampini /* check to see if progress is hopeless */ 5114a221d59SStefano Zampini PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 5124a221d59SStefano Zampini if (!snes->reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 5134a221d59SStefano Zampini if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_INNER; 5144a221d59SStefano Zampini snes->numFailures++; 5154a221d59SStefano Zampini /* We're not progressing, so return with the current iterate */ 5164a221d59SStefano Zampini if (snes->reason) break; 5174a221d59SStefano Zampini } 5184a221d59SStefano Zampini if (rho_satisfied) { 5194a221d59SStefano Zampini /* Update function values */ 5204a221d59SStefano Zampini already_done = PETSC_FALSE; 5214800dd8cSBarry Smith fnorm = gnorm; 5224a221d59SStefano Zampini fk = fkp1; 5234a221d59SStefano Zampini 5244a221d59SStefano Zampini /* New residual and linearization point */ 5259566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 5269566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 5274a221d59SStefano Zampini 52885385478SLisandro Dalcin /* Monitor convergence */ 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 5304a221d59SStefano Zampini snes->iter++; 531fbe28522SBarry Smith snes->norm = fnorm; 532c1e67a49SFande Kong snes->xnorm = xnorm; 533c1e67a49SFande Kong snes->ynorm = ynorm; 5349566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5359566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 5369566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 5374a221d59SStefano Zampini 53885385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 5394a221d59SStefano Zampini PetscCall(VecNorm(X, NORM_2, &xnorm)); 5404a221d59SStefano Zampini PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 5414a221d59SStefano Zampini if (snes->reason) break; 5424a221d59SStefano Zampini } 54338442cffSBarry Smith } 544284fb49fSHeeho Park 5454a221d59SStefano Zampini if (snes->iter == maxits) { 5469566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 5474a221d59SStefano Zampini if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 54852392280SLois Curfman McInnes } 5494a221d59SStefano Zampini if (clear_converged_test) { 5509566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 5519566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 5529566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 5535e28bcb6Sprj- } 5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5554800dd8cSBarry Smith } 556284fb49fSHeeho Park 557d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 558d71ae5a4SJacob Faibussowitsch { 5593a40ed3dSBarry Smith PetscFunctionBegin; 5609566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 4)); 5619566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5634800dd8cSBarry Smith } 5646b8b9a38SLisandro Dalcin 565d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 566d71ae5a4SJacob Faibussowitsch { 5673a40ed3dSBarry Smith PetscFunctionBegin; 5689566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 5693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5704800dd8cSBarry Smith } 5714800dd8cSBarry Smith 572d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 573d71ae5a4SJacob Faibussowitsch { 57404d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 5754800dd8cSBarry Smith 5763a40ed3dSBarry Smith PetscFunctionBegin; 577d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 5784a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 5794a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 5804a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 5814a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 5824a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 5834a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 5844a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 5859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 5864a221d59SStefano 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)); 587d0609cedSBarry Smith PetscOptionsHeadEnd(); 5883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5894800dd8cSBarry Smith } 5904800dd8cSBarry Smith 591d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 592d71ae5a4SJacob Faibussowitsch { 59304d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 594ace3abfcSBarry Smith PetscBool iascii; 595a935fc98SLois Curfman McInnes 5963a40ed3dSBarry Smith PetscFunctionBegin; 5979566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 59832077d6dSBarry Smith if (iascii) { 5994a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)snes->deltatol)); 6004a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 6014a221d59SStefano 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)); 6024a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 60319bcc07fSBarry Smith } 6043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 605a935fc98SLois Curfman McInnes } 606f6dfbefdSBarry Smith 607ebe3b25bSBarry Smith /*MC 6084a221d59SStefano Zampini SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 609f6dfbefdSBarry Smith 610f6dfbefdSBarry Smith Options Database Keys: 6114a221d59SStefano Zampini + -snes_tr_tol <tol> - trust region tolerance 6124a221d59SStefano Zampini . -snes_tr_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 6134a221d59SStefano Zampini . -snes_tr_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 6144a221d59SStefano Zampini . -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 6154a221d59SStefano Zampini . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 6164a221d59SStefano Zampini . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 6174a221d59SStefano Zampini . -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL) 6184a221d59SStefano Zampini . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 6194a221d59SStefano 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. 620ebe3b25bSBarry Smith 621f6dfbefdSBarry Smith Reference: 6224a221d59SStefano Zampini . * - "Numerical Optimization" by Nocedal and Wright, chapter 4. 623ebe3b25bSBarry Smith 624ee3001cbSBarry Smith Level: intermediate 625ee3001cbSBarry Smith 6264a221d59SStefano Zampini .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 6274a221d59SStefano Zampini `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 6284a221d59SStefano Zampini `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()` 629ebe3b25bSBarry Smith M*/ 630d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 631d71ae5a4SJacob Faibussowitsch { 63204d7464bSBarry Smith SNES_NEWTONTR *neP; 6334800dd8cSBarry Smith 6343a40ed3dSBarry Smith PetscFunctionBegin; 63504d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 63604d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 63704d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 63804d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 63904d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 640fbe28522SBarry Smith 64142f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 642efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 64342f4f86dSBarry Smith 6444fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 6454fc747eaSLawrence Mitchell 6464dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 647fbe28522SBarry Smith snes->data = (void *)neP; 648fbe28522SBarry Smith neP->delta = 0.0; 649fbe28522SBarry Smith neP->delta0 = 0.2; 6504a221d59SStefano Zampini neP->eta1 = 0.001; 6514a221d59SStefano Zampini neP->eta2 = 0.25; 6524a221d59SStefano Zampini neP->eta3 = 0.75; 6534a221d59SStefano Zampini neP->t1 = 0.25; 6544a221d59SStefano Zampini neP->t2 = 2.0; 6554a221d59SStefano Zampini neP->deltaM = 1.e10; 6564a221d59SStefano Zampini neP->fallback = SNES_TR_FALLBACK_NEWTON; 6573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6584800dd8cSBarry Smith } 659