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}; 1124fb275aSStefano Zampini const char *const SNESNewtonTRQNTypes[] = {"NONE", "SAME", "DIFFERENT", "SNESNewtonTRQNType", "SNES_TR_QN_", NULL}; 1224fb275aSStefano Zampini 133201ab8dSStefano Zampini static PetscErrorCode SNESNewtonTRSetTolerances_TR(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0) 143201ab8dSStefano Zampini { 153201ab8dSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 163201ab8dSStefano Zampini 173201ab8dSStefano Zampini PetscFunctionBegin; 183201ab8dSStefano Zampini if (delta_min == PETSC_DETERMINE) delta_min = tr->default_deltam; 193201ab8dSStefano Zampini if (delta_max == PETSC_DETERMINE) delta_max = tr->default_deltaM; 203201ab8dSStefano Zampini if (delta_0 == PETSC_DETERMINE) delta_0 = tr->default_delta0; 213201ab8dSStefano Zampini if (delta_min != PETSC_CURRENT) tr->deltam = delta_min; 223201ab8dSStefano Zampini if (delta_max != PETSC_CURRENT) tr->deltaM = delta_max; 233201ab8dSStefano Zampini if (delta_0 != PETSC_CURRENT) tr->delta0 = delta_0; 243201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 253201ab8dSStefano Zampini } 263201ab8dSStefano Zampini 273201ab8dSStefano Zampini static PetscErrorCode SNESNewtonTRGetTolerances_TR(SNES snes, PetscReal *delta_min, PetscReal *delta_max, PetscReal *delta_0) 283201ab8dSStefano Zampini { 293201ab8dSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 303201ab8dSStefano Zampini 313201ab8dSStefano Zampini PetscFunctionBegin; 323201ab8dSStefano Zampini if (delta_min) *delta_min = tr->deltam; 333201ab8dSStefano Zampini if (delta_max) *delta_max = tr->deltaM; 343201ab8dSStefano Zampini if (delta_0) *delta_0 = tr->delta0; 353201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 363201ab8dSStefano Zampini } 373201ab8dSStefano Zampini 3824fb275aSStefano Zampini static PetscErrorCode SNESComputeJacobian_MATLMVM(SNES snes, Vec X, Mat J, Mat B, void *dummy) 3924fb275aSStefano Zampini { 4024fb275aSStefano Zampini PetscFunctionBegin; 4124fb275aSStefano Zampini // PetscCall(MatLMVMSymBroydenSetDelta(B, _some_delta)); 4224fb275aSStefano Zampini PetscCall(MatLMVMUpdate(B, X, snes->vec_func)); 4324fb275aSStefano Zampini PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 4424fb275aSStefano Zampini PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 4524fb275aSStefano Zampini if (J != B) { 4624fb275aSStefano Zampini // PetscCall(MatLMVMSymBroydenSetDelta(J, _some_delta)); 4724fb275aSStefano Zampini PetscCall(MatLMVMUpdate(J, X, snes->vec_func)); 4824fb275aSStefano Zampini PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4924fb275aSStefano Zampini PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 5024fb275aSStefano Zampini } 5124fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 5224fb275aSStefano Zampini } 534a221d59SStefano Zampini 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 55d71ae5a4SJacob Faibussowitsch { 56971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 57971273eeSBarry Smith SNES snes = ctx->snes; 5804d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 59df60cc22SBarry Smith Vec x; 60064f8208SBarry Smith PetscReal nrm; 61df60cc22SBarry Smith 623a40ed3dSBarry Smith PetscFunctionBegin; 63a935fc98SLois Curfman McInnes /* Determine norm of solution */ 649566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 6524fb275aSStefano Zampini PetscCall(VecNorm(x, neP->norm, &nrm)); 66064f8208SBarry Smith if (nrm >= neP->delta) { 676d51442aSBarry Smith PetscCall(PetscInfo(snes, "Ending linear iteration early due to exiting trust region, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 68329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 696d51442aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 70df60cc22SBarry Smith } 716d51442aSBarry Smith PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 726d51442aSBarry Smith if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74df60cc22SBarry Smith } 7582bf6240SBarry Smith 76d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 77d71ae5a4SJacob Faibussowitsch { 78971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 79971273eeSBarry Smith 80971273eeSBarry Smith PetscFunctionBegin; 819566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 829566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84971273eeSBarry Smith } 85971273eeSBarry Smith 86d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 87d71ae5a4SJacob Faibussowitsch { 8804d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 8985385478SLisandro Dalcin 9085385478SLisandro Dalcin PetscFunctionBegin; 9185385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 923201ab8dSStefano Zampini if (neP->delta < neP->deltam) { 933201ab8dSStefano Zampini PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)neP->deltam)); 941c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 95e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 9663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 9785385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 9885385478SLisandro Dalcin } 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10085385478SLisandro Dalcin } 10185385478SLisandro Dalcin 1024a221d59SStefano Zampini /*@ 10324fb275aSStefano Zampini SNESNewtonTRSetNormType - Specify the type of norm to use for the computation of the trust region. 10424fb275aSStefano Zampini 10524fb275aSStefano Zampini Input Parameters: 10624fb275aSStefano Zampini + snes - the nonlinear solver object 10724fb275aSStefano Zampini - norm - the norm type 10824fb275aSStefano Zampini 10924fb275aSStefano Zampini Level: intermediate 11024fb275aSStefano Zampini 11124fb275aSStefano Zampini .seealso: `SNESNEWTONTR`, `NormType` 11224fb275aSStefano Zampini @*/ 11324fb275aSStefano Zampini PetscErrorCode SNESNewtonTRSetNormType(SNES snes, NormType norm) 11424fb275aSStefano Zampini { 11524fb275aSStefano Zampini PetscBool flg; 11624fb275aSStefano Zampini 11724fb275aSStefano Zampini PetscFunctionBegin; 11824fb275aSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 11924fb275aSStefano Zampini PetscValidLogicalCollectiveEnum(snes, norm, 2); 12024fb275aSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 12124fb275aSStefano Zampini if (flg) { 12224fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 12324fb275aSStefano Zampini 12424fb275aSStefano Zampini tr->norm = norm; 12524fb275aSStefano Zampini } 12624fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 12724fb275aSStefano Zampini } 12824fb275aSStefano Zampini 12924fb275aSStefano Zampini /*@ 13024fb275aSStefano Zampini SNESNewtonTRSetQNType - Specify to use a quasi-Newton model. 13124fb275aSStefano Zampini 13224fb275aSStefano Zampini Input Parameters: 13324fb275aSStefano Zampini + snes - the nonlinear solver object 13424fb275aSStefano Zampini - use - the type of approximations to be used 13524fb275aSStefano Zampini 13624fb275aSStefano Zampini Level: intermediate 13724fb275aSStefano Zampini 13824fb275aSStefano Zampini Notes: 139*ac84dfd5SBarry Smith Options for the approximations can be set with the `snes_tr_qn_` and `snes_tr_qn_pre_` prefixes. 14024fb275aSStefano Zampini 14124fb275aSStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRQNType`, `MATLMVM` 14224fb275aSStefano Zampini @*/ 14324fb275aSStefano Zampini PetscErrorCode SNESNewtonTRSetQNType(SNES snes, SNESNewtonTRQNType use) 14424fb275aSStefano Zampini { 14524fb275aSStefano Zampini PetscBool flg; 14624fb275aSStefano Zampini 14724fb275aSStefano Zampini PetscFunctionBegin; 14824fb275aSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 14924fb275aSStefano Zampini PetscValidLogicalCollectiveEnum(snes, use, 2); 15024fb275aSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 15124fb275aSStefano Zampini if (flg) { 15224fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 15324fb275aSStefano Zampini 15424fb275aSStefano Zampini tr->qn = use; 15524fb275aSStefano Zampini } 15624fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 15724fb275aSStefano Zampini } 15824fb275aSStefano Zampini 15924fb275aSStefano Zampini /*@ 160420bcc1bSBarry Smith SNESNewtonTRSetFallbackType - Set the type of fallback to use if the solution of the trust region subproblem is outside the radius 1614a221d59SStefano Zampini 1624a221d59SStefano Zampini Input Parameters: 1634a221d59SStefano Zampini + snes - the nonlinear solver object 1644a221d59SStefano Zampini - ftype - the fallback type, see `SNESNewtonTRFallbackType` 1654a221d59SStefano Zampini 1664a221d59SStefano Zampini Level: intermediate 1674a221d59SStefano Zampini 168420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`, 1694a221d59SStefano Zampini `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 1704a221d59SStefano Zampini @*/ 1714a221d59SStefano Zampini PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype) 1724a221d59SStefano Zampini { 1734a221d59SStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1744a221d59SStefano Zampini PetscBool flg; 1754a221d59SStefano Zampini 1764a221d59SStefano Zampini PetscFunctionBegin; 1774a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1784a221d59SStefano Zampini PetscValidLogicalCollectiveEnum(snes, ftype, 2); 1794a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1804a221d59SStefano Zampini if (flg) tr->fallback = ftype; 1814a221d59SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 1824a221d59SStefano Zampini } 1834a221d59SStefano Zampini 1847cb011f5SBarry Smith /*@C 185c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 1864a221d59SStefano Zampini Allows the user a chance to change or override the trust region decision. 187f6dfbefdSBarry Smith 188c3339decSBarry Smith Logically Collective 189c9368356SGlenn Hammond 190c9368356SGlenn Hammond Input Parameters: 191c9368356SGlenn Hammond + snes - the nonlinear solver object 19220f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 19320f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 194c9368356SGlenn Hammond 195acba1f63SStefano Zampini Level: intermediate 196c9368356SGlenn Hammond 197f6dfbefdSBarry Smith Note: 1984a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver. 199c9368356SGlenn Hammond 200420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 201c9368356SGlenn Hammond @*/ 202d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 203d71ae5a4SJacob Faibussowitsch { 204c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2054a221d59SStefano Zampini PetscBool flg; 206c9368356SGlenn Hammond 207c9368356SGlenn Hammond PetscFunctionBegin; 208c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2094a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2104a221d59SStefano Zampini if (flg) { 211c9368356SGlenn Hammond if (func) tr->precheck = func; 212c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 2134a221d59SStefano Zampini } 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215c9368356SGlenn Hammond } 216c9368356SGlenn Hammond 217c9368356SGlenn Hammond /*@C 218c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 219c9368356SGlenn Hammond 22020f4b53cSBarry Smith Not Collective 221c9368356SGlenn Hammond 222c9368356SGlenn Hammond Input Parameter: 223c9368356SGlenn Hammond . snes - the nonlinear solver context 224c9368356SGlenn Hammond 225c9368356SGlenn Hammond Output Parameters: 22620f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 22720f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 228c9368356SGlenn Hammond 229acba1f63SStefano Zampini Level: intermediate 230c9368356SGlenn Hammond 231420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 232c9368356SGlenn Hammond @*/ 233d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 234d71ae5a4SJacob Faibussowitsch { 235c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2364a221d59SStefano Zampini PetscBool flg; 237c9368356SGlenn Hammond 238c9368356SGlenn Hammond PetscFunctionBegin; 239c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2404a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2414a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 242c9368356SGlenn Hammond if (func) *func = tr->precheck; 243c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245c9368356SGlenn Hammond } 246c9368356SGlenn Hammond 247c9368356SGlenn Hammond /*@C 2487cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 2494a221d59SStefano Zampini function evaluation. Allows the user a chance to change or override the internal decision of the solver 250f6dfbefdSBarry Smith 251c3339decSBarry Smith Logically Collective 2527cb011f5SBarry Smith 2537cb011f5SBarry Smith Input Parameters: 2547cb011f5SBarry Smith + snes - the nonlinear solver object 25520f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 25620f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 2577cb011f5SBarry Smith 258acba1f63SStefano Zampini Level: intermediate 2597cb011f5SBarry Smith 260f6dfbefdSBarry Smith Note: 2614a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver while the function set in 262f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 2637cb011f5SBarry Smith 264420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 2657cb011f5SBarry Smith @*/ 266d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 267d71ae5a4SJacob Faibussowitsch { 2687cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2694a221d59SStefano Zampini PetscBool flg; 2707cb011f5SBarry Smith 2717cb011f5SBarry Smith PetscFunctionBegin; 2727cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2734a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2744a221d59SStefano Zampini if (flg) { 2757cb011f5SBarry Smith if (func) tr->postcheck = func; 2767cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 2774a221d59SStefano Zampini } 2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2797cb011f5SBarry Smith } 2807cb011f5SBarry Smith 2817cb011f5SBarry Smith /*@C 2827cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 2837cb011f5SBarry Smith 28420f4b53cSBarry Smith Not Collective 2857cb011f5SBarry Smith 2867cb011f5SBarry Smith Input Parameter: 2877cb011f5SBarry Smith . snes - the nonlinear solver context 2887cb011f5SBarry Smith 2897cb011f5SBarry Smith Output Parameters: 29020f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 29120f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 2927cb011f5SBarry Smith 2937cb011f5SBarry Smith Level: intermediate 2947cb011f5SBarry Smith 295420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 2967cb011f5SBarry Smith @*/ 297d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 298d71ae5a4SJacob Faibussowitsch { 2997cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3004a221d59SStefano Zampini PetscBool flg; 3017cb011f5SBarry Smith 3027cb011f5SBarry Smith PetscFunctionBegin; 3037cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3044a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3054a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 3067cb011f5SBarry Smith if (func) *func = tr->postcheck; 3077cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3097cb011f5SBarry Smith } 3107cb011f5SBarry Smith 3117cb011f5SBarry Smith /*@C 3124a221d59SStefano Zampini SNESNewtonTRPreCheck - Runs the precheck routine 313c9368356SGlenn Hammond 314c3339decSBarry Smith Logically Collective 315c9368356SGlenn Hammond 316c9368356SGlenn Hammond Input Parameters: 317c9368356SGlenn Hammond + snes - the solver 318c9368356SGlenn Hammond . X - The last solution 319c9368356SGlenn Hammond - Y - The step direction 320c9368356SGlenn Hammond 3212fe279fdSBarry Smith Output Parameter: 3222fe279fdSBarry Smith . changed_Y - Indicator that the step direction `Y` has been changed. 323c9368356SGlenn Hammond 3244a221d59SStefano Zampini Level: intermediate 325c9368356SGlenn Hammond 326420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 327c9368356SGlenn Hammond @*/ 3284a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 329d71ae5a4SJacob Faibussowitsch { 330c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3314a221d59SStefano Zampini PetscBool flg; 332c9368356SGlenn Hammond 333c9368356SGlenn Hammond PetscFunctionBegin; 3344a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3354a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3364a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 337c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 338c9368356SGlenn Hammond if (tr->precheck) { 3399566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 340c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 341c9368356SGlenn Hammond } 3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 343c9368356SGlenn Hammond } 344c9368356SGlenn Hammond 345c9368356SGlenn Hammond /*@C 3464a221d59SStefano Zampini SNESNewtonTRPostCheck - Runs the postcheck routine 3477cb011f5SBarry Smith 348c3339decSBarry Smith Logically Collective 3497cb011f5SBarry Smith 3507cb011f5SBarry Smith Input Parameters: 3516b867d5aSJose E. Roman + snes - the solver 3526b867d5aSJose E. Roman . X - The last solution 3537cb011f5SBarry Smith . Y - The full step direction 3543312a946SBarry Smith - W - The updated solution, W = X - Y 3557cb011f5SBarry Smith 3567cb011f5SBarry Smith Output Parameters: 3573312a946SBarry Smith + changed_Y - indicator if step has been changed 3583312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 3597cb011f5SBarry Smith 360f6dfbefdSBarry Smith Note: 3613312a946SBarry Smith If Y is changed then W is recomputed as X - Y 3627cb011f5SBarry Smith 3634a221d59SStefano Zampini Level: intermediate 3647cb011f5SBarry Smith 365420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()` 3667cb011f5SBarry Smith @*/ 3674a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 368d71ae5a4SJacob Faibussowitsch { 3697cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3704a221d59SStefano Zampini PetscBool flg; 3717cb011f5SBarry Smith 3727cb011f5SBarry Smith PetscFunctionBegin; 3734a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3744a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3754a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 376c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 3777cb011f5SBarry Smith *changed_W = PETSC_FALSE; 3787cb011f5SBarry Smith if (tr->postcheck) { 3799566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 380c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 3817cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 3827cb011f5SBarry Smith } 3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3847cb011f5SBarry Smith } 38585385478SLisandro Dalcin 38624fb275aSStefano Zampini /* stable implementation of roots of a*x^2 + b*x + c = 0 */ 3874a221d59SStefano Zampini static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 3884a221d59SStefano Zampini { 3894a221d59SStefano Zampini PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 3904a221d59SStefano Zampini PetscReal x1 = temp / a; 3914a221d59SStefano Zampini PetscReal x2 = c / temp; 3924a221d59SStefano Zampini *xm = PetscMin(x1, x2); 3934a221d59SStefano Zampini *xp = PetscMax(x1, x2); 3944a221d59SStefano Zampini } 3954a221d59SStefano Zampini 3967aa289d8SStefano Zampini /* Computes the quadratic model difference */ 39724fb275aSStefano Zampini static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, Mat J, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy_, PetscReal *gTy_, PetscReal *deltaqm) 3987aa289d8SStefano Zampini { 39924fb275aSStefano Zampini PetscReal yTHy, gTy; 40024fb275aSStefano Zampini 4017aa289d8SStefano Zampini PetscFunctionBegin; 40224fb275aSStefano Zampini PetscCall(MatMult(J, Y, W)); 40324fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 40424fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */ 40524fb275aSStefano Zampini PetscCall(VecDotRealPart(GradF, Y, &gTy)); 40624fb275aSStefano Zampini *deltaqm = -(-(gTy) + 0.5 * (yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */ 40724fb275aSStefano Zampini if (yTHy_) *yTHy_ = yTHy; 40824fb275aSStefano Zampini if (gTy_) *gTy_ = gTy; 40924fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 41024fb275aSStefano Zampini } 41124fb275aSStefano Zampini 41224fb275aSStefano Zampini /* Computes the new objective given X = Xk, Y = direction 41324fb275aSStefano Zampini W work vector, on output W = X - Y 41424fb275aSStefano Zampini G work vector, on output G = SNESFunction(W) */ 41524fb275aSStefano Zampini static PetscErrorCode SNESNewtonTRObjective(SNES snes, PetscBool has_objective, Vec X, Vec Y, Vec W, Vec G, PetscReal *gnorm, PetscReal *fkp1) 41624fb275aSStefano Zampini { 41724fb275aSStefano Zampini PetscBool changed_y, changed_w; 41824fb275aSStefano Zampini 41924fb275aSStefano Zampini PetscFunctionBegin; 42024fb275aSStefano Zampini /* TODO: we can add a linesearch here */ 42124fb275aSStefano Zampini PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 42224fb275aSStefano Zampini PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 42324fb275aSStefano Zampini PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 42424fb275aSStefano Zampini if (changed_y && !changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 42524fb275aSStefano Zampini 42624fb275aSStefano Zampini PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 42724fb275aSStefano Zampini PetscCall(VecNorm(G, NORM_2, gnorm)); 42824fb275aSStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, W, fkp1)); 42924fb275aSStefano Zampini else *fkp1 = 0.5 * PetscSqr(*gnorm); 43024fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 43124fb275aSStefano Zampini } 43224fb275aSStefano Zampini 43324fb275aSStefano Zampini static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes) 43424fb275aSStefano Zampini { 43524fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 43624fb275aSStefano Zampini 43724fb275aSStefano Zampini PetscFunctionBegin; 43824fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB)); 43924fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre)); 44024fb275aSStefano Zampini if (tr->qn) { 44124fb275aSStefano Zampini PetscInt n, N; 44224fb275aSStefano Zampini const char *optionsprefix; 44324fb275aSStefano Zampini Mat B; 44424fb275aSStefano Zampini 44524fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 44624fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 44724fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_")); 44824fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 44924fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS)); 45024fb275aSStefano Zampini PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 45124fb275aSStefano Zampini PetscCall(VecGetSize(snes->vec_sol, &N)); 45224fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 45324fb275aSStefano Zampini PetscCall(MatSetUp(B)); 45424fb275aSStefano Zampini PetscCall(MatSetFromOptions(B)); 45524fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 45624fb275aSStefano Zampini tr->qnB = B; 45724fb275aSStefano Zampini if (tr->qn == SNES_TR_QN_DIFFERENT) { 45824fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 45924fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 46024fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_")); 46124fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 46224fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS)); 46324fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 46424fb275aSStefano Zampini PetscCall(MatSetUp(B)); 46524fb275aSStefano Zampini PetscCall(MatSetFromOptions(B)); 46624fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 46724fb275aSStefano Zampini tr->qnB_pre = B; 46824fb275aSStefano Zampini } else { 46924fb275aSStefano Zampini PetscCall(PetscObjectReference((PetscObject)tr->qnB)); 47024fb275aSStefano Zampini tr->qnB_pre = tr->qnB; 47124fb275aSStefano Zampini } 47224fb275aSStefano Zampini } 4737aa289d8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 4747aa289d8SStefano Zampini } 4757aa289d8SStefano Zampini 476df60cc22SBarry Smith /* 4774a221d59SStefano Zampini SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 4784a221d59SStefano Zampini (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 4794a221d59SStefano Zampini nonlinear equations 4804800dd8cSBarry Smith 4814800dd8cSBarry Smith */ 482d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 483d71ae5a4SJacob Faibussowitsch { 48404d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 48524fb275aSStefano Zampini Vec X, F, Y, G, W, GradF, YU, Yc; 4864a221d59SStefano Zampini PetscInt maxits, lits; 4874b0a5c37SStefano Zampini PetscReal rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm; 4888b630c82SStefano Zampini PetscReal fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0; 48924fb275aSStefano Zampini PetscReal auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0; 490a0254a93SStefano Zampini PC pc; 49124fb275aSStefano Zampini Mat J, Jp; 49237d4d9b4SStefano Zampini PetscBool already_done = PETSC_FALSE, on_boundary, use_cauchy; 4937aa289d8SStefano Zampini PetscBool clear_converged_test, rho_satisfied, has_objective; 494df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 4955e28bcb6Sprj- void *convctx; 4966b72add0SBarry Smith SNESObjectiveFn *objective; 4974a221d59SStefano Zampini PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 4984800dd8cSBarry Smith 4993a40ed3dSBarry Smith PetscFunctionBegin; 5004a221d59SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 5017aa289d8SStefano Zampini has_objective = objective ? PETSC_TRUE : PETSC_FALSE; 502c579b300SPatrick Farrell 503fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 504fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 50539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 5064a221d59SStefano Zampini Y = snes->vec_sol_update; /* update vector */ 5074a221d59SStefano Zampini G = snes->work[0]; /* updated residual */ 5084a221d59SStefano Zampini W = snes->work[1]; /* temporary vector */ 5097aa289d8SStefano Zampini GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 5104a221d59SStefano Zampini YU = snes->work[3]; /* work vector for dogleg method */ 51124fb275aSStefano Zampini Yc = snes->work[4]; /* Cauchy point */ 5124a221d59SStefano Zampini 5134a221d59SStefano 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); 5144800dd8cSBarry Smith 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 5164c49b128SBarry Smith snes->iter = 0; 5179566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5184800dd8cSBarry Smith 51924fb275aSStefano Zampini /* setup QN matrices if needed */ 52024fb275aSStefano Zampini PetscCall(SNESSetUpQN_NEWTONTR(snes)); 52124fb275aSStefano Zampini 5224a221d59SStefano Zampini /* Set the linear stopping criteria to use the More' trick if needed */ 5234a221d59SStefano Zampini clear_converged_test = PETSC_FALSE; 524a0254a93SStefano Zampini PetscCall(SNESGetKSP(snes, &snes->ksp)); 525a0254a93SStefano Zampini PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy)); 526fcc61681SStefano Zampini if (convtest != SNESTR_KSPConverged_Private) { 5274a221d59SStefano Zampini clear_converged_test = PETSC_TRUE; 5289566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 529df8705c3SBarry Smith ctx->snes = snes; 530a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 531a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 5329566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 533df8705c3SBarry Smith } 534df8705c3SBarry Smith 535e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 5369566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 5371aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 538e4ed7901SPeter Brune 5399566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 540422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 5419566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 5424a221d59SStefano Zampini 5439566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 544fbe28522SBarry Smith snes->norm = fnorm; 5459566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5464a221d59SStefano Zampini delta = neP->delta0; 5474800dd8cSBarry Smith neP->delta = delta; 5489566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 549b37302c6SLois Curfman McInnes 55085385478SLisandro Dalcin /* test convergence */ 5514a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 5522d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 5532d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 5543ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 5553f149594SLisandro Dalcin 5567aa289d8SStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 5574a221d59SStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 55899a96b7cSMatthew Knepley 559a0254a93SStefano Zampini /* hook state vector to BFGS preconditioner */ 560a0254a93SStefano Zampini PetscCall(KSPGetPC(snes->ksp, &pc)); 561a0254a93SStefano Zampini PetscCall(PCLMVMSetUpdateVec(pc, X)); 562a0254a93SStefano Zampini 56324fb275aSStefano Zampini if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE)); 5646b5873e3SBarry Smith 56524fb275aSStefano Zampini while (snes->iter < maxits) { 56612d0050eSStefano Zampini /* calculating Jacobian and GradF of minimization function only once */ 5674a221d59SStefano Zampini if (!already_done) { 56812d0050eSStefano Zampini /* Call general purpose update function */ 56912d0050eSStefano Zampini PetscTryTypeMethod(snes, update, snes->iter); 57012d0050eSStefano Zampini 5714b0a5c37SStefano Zampini /* apply the nonlinear preconditioner */ 5724b0a5c37SStefano Zampini if (snes->npc && snes->npcside == PC_RIGHT) { 5734b0a5c37SStefano Zampini SNESConvergedReason reason; 5744b0a5c37SStefano Zampini 5754b0a5c37SStefano Zampini PetscCall(SNESSetInitialFunction(snes->npc, F)); 5764b0a5c37SStefano Zampini PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 5774b0a5c37SStefano Zampini PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 5784b0a5c37SStefano Zampini PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 5794b0a5c37SStefano Zampini PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 5804b0a5c37SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 5814b0a5c37SStefano Zampini snes->reason = SNES_DIVERGED_INNER; 5824b0a5c37SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 5834b0a5c37SStefano Zampini } 5844b0a5c37SStefano Zampini // XXX 5854b0a5c37SStefano Zampini PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 58612d0050eSStefano Zampini } 58712d0050eSStefano Zampini 58812d0050eSStefano Zampini /* Jacobian */ 58924fb275aSStefano Zampini J = NULL; 59024fb275aSStefano Zampini Jp = NULL; 59124fb275aSStefano Zampini if (!neP->qnB) { 5924a221d59SStefano Zampini PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 59324fb275aSStefano Zampini J = snes->jacobian; 59424fb275aSStefano Zampini Jp = snes->jacobian_pre; 59524fb275aSStefano Zampini } else { /* QN model */ 59624fb275aSStefano Zampini PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL)); 59724fb275aSStefano Zampini J = neP->qnB; 59824fb275aSStefano Zampini Jp = neP->qnB_pre; 59924fb275aSStefano Zampini } 6004a221d59SStefano Zampini SNESCheckJacobianDomainerror(snes); 60112d0050eSStefano Zampini 60224fb275aSStefano Zampini /* objective function */ 60324fb275aSStefano Zampini PetscCall(VecNorm(F, NORM_2, &fnorm)); 60424fb275aSStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 60524fb275aSStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 60624fb275aSStefano Zampini 60712d0050eSStefano Zampini /* GradF */ 6087aa289d8SStefano Zampini if (has_objective) gfnorm = fnorm; 6097aa289d8SStefano Zampini else { 61024fb275aSStefano Zampini PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */ 6117aa289d8SStefano Zampini PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 612fbe28522SBarry Smith } 61324fb275aSStefano Zampini PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k)); 6147aa289d8SStefano Zampini } 6157aa289d8SStefano Zampini already_done = PETSC_TRUE; 6167aa289d8SStefano Zampini 6174b0a5c37SStefano Zampini /* solve trust-region subproblem */ 6184b0a5c37SStefano Zampini 61924fb275aSStefano Zampini /* first compute Cauchy Point */ 62024fb275aSStefano Zampini PetscCall(MatMult(J, GradF, W)); 62124fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 62224fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 62324fb275aSStefano Zampini /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */ 62424fb275aSStefano Zampini auk = delta / gfnorm_k; 62524fb275aSStefano Zampini if (gTBg < 0.0) tauk = 1.0; 62624fb275aSStefano Zampini else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1); 62724fb275aSStefano Zampini auk *= tauk; 62824fb275aSStefano Zampini ycnorm = auk * gfnorm; 62924fb275aSStefano Zampini PetscCall(VecAXPBY(Yc, auk, 0.0, GradF)); 63024fb275aSStefano Zampini 63124fb275aSStefano Zampini on_boundary = PETSC_FALSE; 63237d4d9b4SStefano Zampini use_cauchy = (PetscBool)(tauk == 1.0 && has_objective); 63337d4d9b4SStefano Zampini if (!use_cauchy) { 63424fb275aSStefano Zampini KSPConvergedReason reason; 63524fb275aSStefano Zampini 6364b0a5c37SStefano Zampini /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods") 63724fb275aSStefano Zampini beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */ 63824fb275aSStefano Zampini objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta); 639fb01098fSStefano Zampini PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin)); 6404b0a5c37SStefano Zampini 64124fb275aSStefano Zampini /* specify radius if looking for Newton step and trust region norm is the l2 norm */ 64224fb275aSStefano Zampini PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0)); 64324fb275aSStefano Zampini PetscCall(KSPSetOperators(snes->ksp, J, Jp)); 6444a221d59SStefano Zampini PetscCall(KSPSolve(snes->ksp, F, Y)); 6454a221d59SStefano Zampini SNESCheckKSPSolve(snes); 6464a221d59SStefano Zampini PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 64724fb275aSStefano Zampini PetscCall(KSPGetConvergedReason(snes->ksp, &reason)); 64824fb275aSStefano Zampini on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH); 6494b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 65024fb275aSStefano Zampini if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */ 65124fb275aSStefano Zampini PetscReal emax, emin; 65224fb275aSStefano Zampini PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin)); 65324fb275aSStefano Zampini if (emax > 0.0) beta_k = emax + 1; 65424fb275aSStefano Zampini } 65524fb275aSStefano Zampini } else { /* Cauchy point is on the boundary, accept it */ 65624fb275aSStefano Zampini on_boundary = PETSC_TRUE; 65724fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 65824fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg)); 65924fb275aSStefano Zampini } 66024fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm)); 6614800dd8cSBarry Smith 6624a221d59SStefano Zampini /* decide what to do when the update is outside of trust region */ 66337d4d9b4SStefano Zampini if (!use_cauchy && (ynorm > delta || ynorm == 0.0)) { 6645ec2728bSStefano Zampini SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY; 6655ec2728bSStefano Zampini 66624fb275aSStefano Zampini PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented"); 6675ec2728bSStefano Zampini switch (fallback) { 6684a221d59SStefano Zampini case SNES_TR_FALLBACK_NEWTON: 6694a221d59SStefano Zampini auk = delta / ynorm; 6704a221d59SStefano Zampini PetscCall(VecScale(Y, auk)); 6714b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm)); 6724a221d59SStefano Zampini break; 6734a221d59SStefano Zampini case SNES_TR_FALLBACK_CAUCHY: 6744a221d59SStefano Zampini case SNES_TR_FALLBACK_DOGLEG: 6755ec2728bSStefano Zampini if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 67624fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 6774a221d59SStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 6784a221d59SStefano Zampini } else { /* take linear combination of Cauchy and Newton direction and step */ 67924fb275aSStefano Zampini auk = gfnorm * gfnorm / gTBg; 68024fb275aSStefano Zampini if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */ 68124fb275aSStefano Zampini PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF)); 68224fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm)); 68324fb275aSStefano Zampini } else { /* second leg */ 6844a221d59SStefano Zampini PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 6854a221d59SStefano Zampini PetscBool noroots; 686284fb49fSHeeho Park 68724fb275aSStefano Zampini /* Find solutions of (Eq. 4.16 in Nocedal and Wright) 68824fb275aSStefano Zampini ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0, 68924fb275aSStefano Zampini where p_U the Cauchy direction, p_B the Newton direction */ 6904a221d59SStefano Zampini PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 6914a221d59SStefano Zampini PetscCall(VecAXPY(Y, -1.0, YU)); 6924a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); 6934a221d59SStefano Zampini PetscCall(VecDotRealPart(YU, Y, &c1)); 6944a221d59SStefano Zampini c0 = PetscSqr(c0); 6954a221d59SStefano Zampini c2 = PetscSqr(ycnorm) - PetscSqr(delta); 69624fb275aSStefano Zampini PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos); 6974a221d59SStefano Zampini 69824fb275aSStefano Zampini /* In principle the DL strategy as a unique solution in [0,1] 69924fb275aSStefano Zampini here we check that for some reason we numerically failed 70024fb275aSStefano Zampini to compute it. In that case, we use the Cauchy point */ 7014a221d59SStefano Zampini noroots = PetscIsInfOrNanReal(tneg); 70224fb275aSStefano Zampini if (!noroots) { 70324fb275aSStefano Zampini if (tpos > 1) { 70424fb275aSStefano Zampini if (tneg >= 0 && tneg <= 1) { 70524fb275aSStefano Zampini tau = tneg; 70624fb275aSStefano Zampini } else noroots = PETSC_TRUE; 70724fb275aSStefano Zampini } else if (tpos >= 0) { 70824fb275aSStefano Zampini tau = tpos; 70924fb275aSStefano Zampini } else noroots = PETSC_TRUE; 71024fb275aSStefano Zampini } 7114a221d59SStefano Zampini if (noroots) { /* No roots, select Cauchy point */ 71224fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 7134a221d59SStefano Zampini } else { 71424fb275aSStefano Zampini PetscCall(VecAXPBY(Y, 1.0, tau, YU)); 7154a221d59SStefano Zampini } 71624fb275aSStefano Zampini PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)gTBg)); 7174a221d59SStefano Zampini } 7184a221d59SStefano Zampini } 7194a221d59SStefano Zampini break; 7204a221d59SStefano Zampini default: 7214a221d59SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 722454a90a3SBarry Smith break; 72352392280SLois Curfman McInnes } 7244800dd8cSBarry Smith } 7254a221d59SStefano Zampini 7267aa289d8SStefano Zampini /* compute the quadratic model difference */ 72724fb275aSStefano Zampini PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm)); 7284a221d59SStefano Zampini 7294a221d59SStefano Zampini /* Compute new objective function */ 73024fb275aSStefano Zampini PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1)); 73124fb275aSStefano Zampini if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1; 73224fb275aSStefano Zampini else { 7334a221d59SStefano Zampini if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 73424fb275aSStefano Zampini else rho = neP->eta1; /* no reduction in quadratic model, step must be rejected */ 73524fb275aSStefano Zampini } 7364a221d59SStefano Zampini 73724fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm)); 73824fb275aSStefano Zampini PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g, ynormk=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy, (double)ynorm)); 73924fb275aSStefano Zampini 74024fb275aSStefano Zampini /* update the size of the trust region */ 7414a221d59SStefano Zampini if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 74224fb275aSStefano Zampini else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */ 7438b630c82SStefano Zampini delta = PetscMin(delta, neP->deltaM); /* but not greater than deltaM */ 7444a221d59SStefano Zampini 74524fb275aSStefano Zampini /* log 2-norm of update for moniroting routines */ 74624fb275aSStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 74724fb275aSStefano Zampini 74824fb275aSStefano Zampini /* decide on new step */ 7494a221d59SStefano Zampini neP->delta = delta; 75024fb275aSStefano Zampini if (rho > neP->eta1) { 7514a221d59SStefano Zampini rho_satisfied = PETSC_TRUE; 7524a221d59SStefano Zampini } else { 7534a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 7544a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 7554a221d59SStefano Zampini /* check to see if progress is hopeless */ 7564a221d59SStefano Zampini PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 7572d157150SStefano Zampini if (!snes->reason) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 7584b0a5c37SStefano Zampini if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_TR_DELTA; 7594a221d59SStefano Zampini snes->numFailures++; 7604a221d59SStefano Zampini /* We're not progressing, so return with the current iterate */ 7614a221d59SStefano Zampini if (snes->reason) break; 7624a221d59SStefano Zampini } 7634a221d59SStefano Zampini if (rho_satisfied) { 7644a221d59SStefano Zampini /* Update function values */ 7654a221d59SStefano Zampini already_done = PETSC_FALSE; 7664800dd8cSBarry Smith fnorm = gnorm; 7674a221d59SStefano Zampini fk = fkp1; 7684a221d59SStefano Zampini 7694a221d59SStefano Zampini /* New residual and linearization point */ 7709566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 7719566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 7724a221d59SStefano Zampini 77385385478SLisandro Dalcin /* Monitor convergence */ 7749566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 7754a221d59SStefano Zampini snes->iter++; 776fbe28522SBarry Smith snes->norm = fnorm; 777c1e67a49SFande Kong snes->xnorm = xnorm; 778c1e67a49SFande Kong snes->ynorm = ynorm; 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7809566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 7814a221d59SStefano Zampini 78285385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 7834a221d59SStefano Zampini PetscCall(VecNorm(X, NORM_2, &xnorm)); 7842d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 7852d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 7864a221d59SStefano Zampini if (snes->reason) break; 7874a221d59SStefano Zampini } 78838442cffSBarry Smith } 789284fb49fSHeeho Park 7904a221d59SStefano Zampini if (clear_converged_test) { 791a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 7929566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 793a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy)); 7945e28bcb6Sprj- } 7953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7964800dd8cSBarry Smith } 797284fb49fSHeeho Park 798d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 799d71ae5a4SJacob Faibussowitsch { 8003a40ed3dSBarry Smith PetscFunctionBegin; 80124fb275aSStefano Zampini PetscCall(SNESSetWorkVecs(snes, 5)); 8029566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 8033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8044800dd8cSBarry Smith } 8056b8b9a38SLisandro Dalcin 80624fb275aSStefano Zampini static PetscErrorCode SNESReset_NEWTONTR(SNES snes) 80724fb275aSStefano Zampini { 80824fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 80924fb275aSStefano Zampini 81024fb275aSStefano Zampini PetscFunctionBegin; 81124fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB)); 81224fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre)); 81324fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 81424fb275aSStefano Zampini } 81524fb275aSStefano Zampini 816d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 817d71ae5a4SJacob Faibussowitsch { 8183a40ed3dSBarry Smith PetscFunctionBegin; 81924fb275aSStefano Zampini PetscCall(SNESReset_NEWTONTR(snes)); 8203201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL)); 8213201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", NULL)); 8229566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8244800dd8cSBarry Smith } 8254800dd8cSBarry Smith 826ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems PetscOptionsObject) 827d71ae5a4SJacob Faibussowitsch { 82804d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 82924fb275aSStefano Zampini SNESNewtonTRQNType qn; 83024fb275aSStefano Zampini SNESNewtonTRFallbackType fallback; 83124fb275aSStefano Zampini NormType norm; 83224fb275aSStefano Zampini PetscBool flg; 8334800dd8cSBarry Smith 8343a40ed3dSBarry Smith PetscFunctionBegin; 835d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 8363201ab8dSStefano Zampini PetscCall(PetscOptionsDeprecated("-snes_tr_deltaM", "-snes_tr_deltamax", "3.22", NULL)); 8373201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "SNESNewtonTRSetUpdateParameters", ctx->eta1, &ctx->eta1, NULL)); 8383201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "SNESNewtonTRSetUpdateParameters", ctx->eta2, &ctx->eta2, NULL)); 8393201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "SNESNewtonTRSetUpdateParameters", ctx->eta3, &ctx->eta3, NULL)); 8403201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "SNESNewtonTRSetUpdateParameters", ctx->t1, &ctx->t1, NULL)); 8413201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "SNESNewtonTRSetUpdateParameters", ctx->t2, &ctx->t2, NULL)); 8423201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_delta0", "Initial trust region size", "SNESNewtonTRSetTolerances", ctx->delta0, &ctx->delta0, NULL)); 8433201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltamin", "Minimum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltam, &ctx->deltam, NULL)); 8443201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltamax", "Maximum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltaM, &ctx->deltaM, NULL)); 8454b0a5c37SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL)); 84624fb275aSStefano Zampini 84724fb275aSStefano Zampini fallback = ctx->fallback; 84824fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)fallback, (PetscEnum *)&fallback, &flg)); 84924fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback)); 85024fb275aSStefano Zampini 85124fb275aSStefano Zampini qn = ctx->qn; 85224fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg)); 85324fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn)); 85424fb275aSStefano Zampini 85524fb275aSStefano Zampini norm = ctx->norm; 85624fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg)); 85724fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm)); 85824fb275aSStefano Zampini 859d0609cedSBarry Smith PetscOptionsHeadEnd(); 8603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8614800dd8cSBarry Smith } 8624800dd8cSBarry Smith 863d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 864d71ae5a4SJacob Faibussowitsch { 86504d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 866ace3abfcSBarry Smith PetscBool iascii; 867a935fc98SLois Curfman McInnes 8683a40ed3dSBarry Smith PetscFunctionBegin; 8699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 87032077d6dSBarry Smith if (iascii) { 8713201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region parameters:\n")); 8724a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 8733201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " t1=%g, t2=%g\n", (double)tr->t1, (double)tr->t2)); 8743201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " delta_min=%g, delta_0=%g, delta_max=%g\n", (double)tr->deltam, (double)tr->delta0, (double)tr->deltaM)); 8754b0a5c37SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " kmdc=%g\n", (double)tr->kmdc)); 8764a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 87724fb275aSStefano Zampini if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, " qn=%s\n", SNESNewtonTRQNTypes[tr->qn])); 87824fb275aSStefano Zampini if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, " norm=%s\n", NormTypes[tr->norm])); 87919bcc07fSBarry Smith } 8803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 881a935fc98SLois Curfman McInnes } 882f6dfbefdSBarry Smith 8833201ab8dSStefano Zampini /*@ 8843201ab8dSStefano Zampini SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 8853201ab8dSStefano Zampini 8863201ab8dSStefano Zampini Logically Collective 8873201ab8dSStefano Zampini 8883201ab8dSStefano Zampini Input Parameters: 8893201ab8dSStefano Zampini + snes - the `SNES` context 8903201ab8dSStefano Zampini - tol - tolerance 8913201ab8dSStefano Zampini 8923201ab8dSStefano Zampini Level: deprecated 8933201ab8dSStefano Zampini 8943201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetTolerances()` 8953201ab8dSStefano Zampini @*/ 8963201ab8dSStefano Zampini PetscErrorCode SNESSetTrustRegionTolerance(SNES snes, PetscReal tol) 8973201ab8dSStefano Zampini { 8983201ab8dSStefano Zampini return SNESNewtonTRSetTolerances(snes, tol, PETSC_CURRENT, PETSC_CURRENT); 8993201ab8dSStefano Zampini } 9003201ab8dSStefano Zampini 9013201ab8dSStefano Zampini /*@ 9023201ab8dSStefano Zampini SNESNewtonTRSetTolerances - Sets the trust region parameter tolerances. 9033201ab8dSStefano Zampini 9043201ab8dSStefano Zampini Logically Collective 9053201ab8dSStefano Zampini 9063201ab8dSStefano Zampini Input Parameters: 9073201ab8dSStefano Zampini + snes - the `SNES` context 9083201ab8dSStefano Zampini . delta_min - minimum allowed trust region size 9093201ab8dSStefano Zampini . delta_max - maximum allowed trust region size 9103201ab8dSStefano Zampini - delta_0 - initial trust region size 9113201ab8dSStefano Zampini 9123201ab8dSStefano Zampini Options Database Key: 9133201ab8dSStefano Zampini + -snes_tr_deltamin <tol> - Set minimum size 9143201ab8dSStefano Zampini . -snes_tr_deltamax <tol> - Set maximum size 9153201ab8dSStefano Zampini - -snes_tr_delta0 <tol> - Set initial size 9163201ab8dSStefano Zampini 9173201ab8dSStefano Zampini Note: 9183201ab8dSStefano Zampini Use `PETSC_DETERMINE` to use the default value for the given `SNES`. 9193201ab8dSStefano Zampini Use `PETSC_CURRENT` to retain a value. 9203201ab8dSStefano Zampini 9213201ab8dSStefano Zampini Fortran Note: 9223201ab8dSStefano Zampini Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL` 9233201ab8dSStefano Zampini 9243201ab8dSStefano Zampini Level: intermediate 9253201ab8dSStefano Zampini 9263201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRGetTolerances()` 9273201ab8dSStefano Zampini @*/ 9283201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRSetTolerances(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0) 9293201ab8dSStefano Zampini { 9303201ab8dSStefano Zampini PetscFunctionBegin; 9313201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 9323201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, delta_min, 2); 9333201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, delta_max, 3); 9343201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, delta_0, 4); 9353201ab8dSStefano Zampini PetscTryMethod(snes, "SNESNewtonTRSetTolerances_C", (SNES, PetscReal, PetscReal, PetscReal), (snes, delta_min, delta_max, delta_0)); 9363201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9373201ab8dSStefano Zampini } 9383201ab8dSStefano Zampini 9393201ab8dSStefano Zampini /*@ 9403201ab8dSStefano Zampini SNESNewtonTRGetTolerances - Gets the trust region parameter tolerances. 9413201ab8dSStefano Zampini 9423201ab8dSStefano Zampini Not Collective 9433201ab8dSStefano Zampini 9443201ab8dSStefano Zampini Input Parameter: 9453201ab8dSStefano Zampini . snes - the `SNES` context 9463201ab8dSStefano Zampini 9473201ab8dSStefano Zampini Output Parameters: 9483201ab8dSStefano Zampini + delta_min - minimum allowed trust region size or `NULL` 9493201ab8dSStefano Zampini . delta_max - maximum allowed trust region size or `NULL` 9503201ab8dSStefano Zampini - delta_0 - initial trust region size or `NULL` 9513201ab8dSStefano Zampini 9523201ab8dSStefano Zampini Level: intermediate 9533201ab8dSStefano Zampini 9543201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetTolerances()` 9553201ab8dSStefano Zampini @*/ 9563201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRGetTolerances(SNES snes, PetscReal *delta_min, PetscReal *delta_max, PetscReal *delta_0) 9573201ab8dSStefano Zampini { 9583201ab8dSStefano Zampini PetscFunctionBegin; 9593201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 9603201ab8dSStefano Zampini if (delta_min) PetscAssertPointer(delta_min, 2); 9613201ab8dSStefano Zampini if (delta_max) PetscAssertPointer(delta_max, 3); 9623201ab8dSStefano Zampini if (delta_0) PetscAssertPointer(delta_0, 4); 9633201ab8dSStefano Zampini PetscUseMethod(snes, "SNESNewtonTRGetTolerances_C", (SNES, PetscReal *, PetscReal *, PetscReal *), (snes, delta_min, delta_max, delta_0)); 9643201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 9653201ab8dSStefano Zampini } 9663201ab8dSStefano Zampini 9673201ab8dSStefano Zampini /*@ 9683201ab8dSStefano Zampini SNESNewtonTRSetUpdateParameters - Sets the trust region update parameters. 9693201ab8dSStefano Zampini 9703201ab8dSStefano Zampini Logically Collective 9713201ab8dSStefano Zampini 9723201ab8dSStefano Zampini Input Parameters: 9733201ab8dSStefano Zampini + snes - the `SNES` context 9743201ab8dSStefano Zampini . eta1 - acceptance tolerance 9753201ab8dSStefano Zampini . eta2 - shrinking tolerance 9763201ab8dSStefano Zampini . eta3 - enlarging tolerance 9773201ab8dSStefano Zampini . t1 - shrink factor 9783201ab8dSStefano Zampini - t2 - enlarge factor 9793201ab8dSStefano Zampini 9803201ab8dSStefano Zampini Options Database Key: 981*ac84dfd5SBarry Smith + -snes_tr_eta1 <tol> - Set `eta1` 982*ac84dfd5SBarry Smith . -snes_tr_eta2 <tol> - Set `eta2` 983*ac84dfd5SBarry Smith . -snes_tr_eta3 <tol> - Set `eta3` 984*ac84dfd5SBarry Smith . -snes_tr_t1 <tol> - Set `t1` 985*ac84dfd5SBarry Smith - -snes_tr_t2 <tol> - Set `t2` 9863201ab8dSStefano Zampini 9873201ab8dSStefano Zampini Notes: 9883201ab8dSStefano Zampini Given the ratio $\rho = \frac{f(x_k) - f(x_k+s_k)}{m(0) - m(s_k)}$, with $x_k$ the current iterate, 9893201ab8dSStefano Zampini $s_k$ the computed step, $f$ the objective function, and $m$ the quadratic model, the trust region 9903201ab8dSStefano Zampini radius is modified as follows 9913201ab8dSStefano Zampini $$ 9923201ab8dSStefano Zampini \delta = 9933201ab8dSStefano Zampini \begin{cases} 9943201ab8dSStefano Zampini \delta * t_1 ,& \rho < \eta_2 \\ 9953201ab8dSStefano Zampini \delta * t_2 ,& \rho > \eta_3 \\ 9963201ab8dSStefano Zampini \end{cases} 9973201ab8dSStefano Zampini $$ 9983201ab8dSStefano Zampini The step is accepted if $\rho > \eta_1$. 9993201ab8dSStefano Zampini Use `PETSC_DETERMINE` to use the default value for the given `SNES`. 10003201ab8dSStefano Zampini Use `PETSC_CURRENT` to retain a value. 10013201ab8dSStefano Zampini 10023201ab8dSStefano Zampini Fortran Note: 10033201ab8dSStefano Zampini Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL` 10043201ab8dSStefano Zampini 10053201ab8dSStefano Zampini Level: intermediate 10063201ab8dSStefano Zampini 10073201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetObjective()`, `SNESNewtonTRGetUpdateParameters()` 10083201ab8dSStefano Zampini @*/ 10093201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRSetUpdateParameters(SNES snes, PetscReal eta1, PetscReal eta2, PetscReal eta3, PetscReal t1, PetscReal t2) 10103201ab8dSStefano Zampini { 10113201ab8dSStefano Zampini PetscBool flg; 10123201ab8dSStefano Zampini 10133201ab8dSStefano Zampini PetscFunctionBegin; 10143201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 10153201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, eta1, 2); 10163201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, eta2, 3); 10173201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, eta3, 4); 10183201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, t1, 5); 10193201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, t2, 6); 10203201ab8dSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 10213201ab8dSStefano Zampini if (flg) { 10223201ab8dSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 10233201ab8dSStefano Zampini 10243201ab8dSStefano Zampini if (eta1 == PETSC_DETERMINE) eta1 = tr->default_eta1; 10253201ab8dSStefano Zampini if (eta2 == PETSC_DETERMINE) eta2 = tr->default_eta2; 10263201ab8dSStefano Zampini if (eta3 == PETSC_DETERMINE) eta3 = tr->default_eta3; 10273201ab8dSStefano Zampini if (t1 == PETSC_DETERMINE) t1 = tr->default_t1; 10283201ab8dSStefano Zampini if (t2 == PETSC_DETERMINE) t2 = tr->default_t2; 10293201ab8dSStefano Zampini if (eta1 != PETSC_CURRENT) tr->eta1 = eta1; 10303201ab8dSStefano Zampini if (eta2 != PETSC_CURRENT) tr->eta2 = eta2; 10313201ab8dSStefano Zampini if (eta3 != PETSC_CURRENT) tr->eta3 = eta3; 10323201ab8dSStefano Zampini if (t1 != PETSC_CURRENT) tr->t1 = t1; 10333201ab8dSStefano Zampini if (t2 != PETSC_CURRENT) tr->t2 = t2; 10343201ab8dSStefano Zampini } 10353201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10363201ab8dSStefano Zampini } 10373201ab8dSStefano Zampini 10383201ab8dSStefano Zampini /*@ 10393201ab8dSStefano Zampini SNESNewtonTRGetUpdateParameters - Gets the trust region update parameters. 10403201ab8dSStefano Zampini 10413201ab8dSStefano Zampini Not Collective 10423201ab8dSStefano Zampini 10433201ab8dSStefano Zampini Input Parameter: 10443201ab8dSStefano Zampini . snes - the `SNES` context 10453201ab8dSStefano Zampini 10463201ab8dSStefano Zampini Output Parameters: 10473201ab8dSStefano Zampini + eta1 - acceptance tolerance 10483201ab8dSStefano Zampini . eta2 - shrinking tolerance 10493201ab8dSStefano Zampini . eta3 - enlarging tolerance 10503201ab8dSStefano Zampini . t1 - shrink factor 10513201ab8dSStefano Zampini - t2 - enlarge factor 10523201ab8dSStefano Zampini 10533201ab8dSStefano Zampini Level: intermediate 10543201ab8dSStefano Zampini 10553201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetUpdateParameters()` 10563201ab8dSStefano Zampini @*/ 10573201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRGetUpdateParameters(SNES snes, PetscReal *eta1, PetscReal *eta2, PetscReal *eta3, PetscReal *t1, PetscReal *t2) 10583201ab8dSStefano Zampini { 10593201ab8dSStefano Zampini SNES_NEWTONTR *tr; 10603201ab8dSStefano Zampini PetscBool flg; 10613201ab8dSStefano Zampini 10623201ab8dSStefano Zampini PetscFunctionBegin; 10633201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 10643201ab8dSStefano Zampini if (eta1) PetscAssertPointer(eta1, 2); 10653201ab8dSStefano Zampini if (eta2) PetscAssertPointer(eta2, 3); 10663201ab8dSStefano Zampini if (eta3) PetscAssertPointer(eta3, 4); 10673201ab8dSStefano Zampini if (t1) PetscAssertPointer(t1, 5); 10683201ab8dSStefano Zampini if (t2) PetscAssertPointer(t2, 6); 10693201ab8dSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 10703201ab8dSStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 10713201ab8dSStefano Zampini tr = (SNES_NEWTONTR *)snes->data; 10723201ab8dSStefano Zampini if (eta1) *eta1 = tr->eta1; 10733201ab8dSStefano Zampini if (eta2) *eta2 = tr->eta2; 10743201ab8dSStefano Zampini if (eta3) *eta3 = tr->eta3; 10753201ab8dSStefano Zampini if (t1) *t1 = tr->t1; 10763201ab8dSStefano Zampini if (t2) *t2 = tr->t2; 10773201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10783201ab8dSStefano Zampini } 10793201ab8dSStefano Zampini 1080ebe3b25bSBarry Smith /*MC 10813201ab8dSStefano Zampini SNESNEWTONTR - Newton based nonlinear solver that uses a trust-region strategy 1082f6dfbefdSBarry Smith 1083f6dfbefdSBarry Smith Options Database Keys: 10843201ab8dSStefano Zampini + -snes_tr_deltamin <deltamin> - trust region parameter, minimum size of trust region 10853201ab8dSStefano Zampini . -snes_tr_deltamax <deltamax> - trust region parameter, max size of trust region (default: 1e10) 10863201ab8dSStefano Zampini . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 1087*ac84dfd5SBarry Smith . -snes_tr_eta1 <eta1> - trust region parameter $ eta1 \le eta2$, $rho > eta1 $ breaks out of the inner iteration (default: 0.001) 1088*ac84dfd5SBarry Smith . -snes_tr_eta2 <eta2> - trust region parameter, $ rho \le eta2$ shrinks the trust region (default: 0.25) 1089*ac84dfd5SBarry Smith . -snes_tr_eta3 <eta3> - trust region parameter $eta3 > eta2$, $rho \ge eta3$ expands the trust region (default: 0.75) 10904a221d59SStefano Zampini . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 10914a221d59SStefano Zampini . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 10923201ab8dSStefano Zampini . -snes_tr_norm_type <1,2,infinity> - Type of norm for trust region bounds (default: "2") 10934a221d59SStefano 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. 1094acba1f63SStefano Zampini 1095acba1f63SStefano Zampini Level: beginner 1096b3113221SBarry Smith 10973201ab8dSStefano Zampini Notes: 10983201ab8dSStefano Zampini The code is largely based on the book {cite}`nocedal2006numerical` and supports minimizing objective functions using a quadratic model. 10993201ab8dSStefano Zampini Quasi-Newton models are also supported. 11003201ab8dSStefano Zampini 11013201ab8dSStefano Zampini Default step computation uses the Newton direction, but a dogleg type update is also supported. 11023201ab8dSStefano Zampini The 1- and infinity-norms are also supported when computing the trust region bounds. 11033201ab8dSStefano Zampini 11043201ab8dSStefano Zampini .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetObjective()`, 11053201ab8dSStefano Zampini `SNESNewtonTRSetTolerances()`, `SNESNewtonTRSetUpdateParameters()` 11063201ab8dSStefano Zampini `SNESNewtonTRSetNormType()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()` 11073201ab8dSStefano Zampini `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRSetPreCheck()`, 1108ebe3b25bSBarry Smith M*/ 1109d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 1110d71ae5a4SJacob Faibussowitsch { 111104d7464bSBarry Smith SNES_NEWTONTR *neP; 11124800dd8cSBarry Smith 11133a40ed3dSBarry Smith PetscFunctionBegin; 111404d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 111504d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 111624fb275aSStefano Zampini snes->ops->reset = SNESReset_NEWTONTR; 111704d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 111804d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 111904d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 1120fbe28522SBarry Smith 112177e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 112237d4d9b4SStefano Zampini PetscObjectParameterSetDefault(snes, stol, 0.0); 112337d4d9b4SStefano Zampini 112442f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 11254b0a5c37SStefano Zampini snes->npcside = PC_RIGHT; 11264b0a5c37SStefano Zampini snes->usesnpc = PETSC_TRUE; 112742f4f86dSBarry Smith 11284fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 11294fc747eaSLawrence Mitchell 11304dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 1131fbe28522SBarry Smith snes->data = (void *)neP; 11323201ab8dSStefano Zampini 11333201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, eta1, 0.001); 11343201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, eta2, 0.25); 11353201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, eta3, 0.75); 11363201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, t1, 0.25); 11373201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, t2, 2.0); 11383201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, deltam, PetscDefined(USE_REAL_SINGLE) ? 1.e-6 : 1.e-12); 11393201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, delta0, 0.2); 11403201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, deltaM, 1.e10); 11413201ab8dSStefano Zampini 114224fb275aSStefano Zampini neP->norm = NORM_2; 11434a221d59SStefano Zampini neP->fallback = SNES_TR_FALLBACK_NEWTON; 11444b0a5c37SStefano Zampini neP->kmdc = 0.0; /* by default do not use sufficient decrease */ 11453201ab8dSStefano Zampini 11463201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TR)); 11473201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", SNESNewtonTRGetTolerances_TR)); 11483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11494800dd8cSBarry Smith } 1150