1c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/
24800dd8cSBarry Smith
3971273eeSBarry Smith typedef struct {
4971273eeSBarry Smith SNES snes;
54d4d2bdcSBarry Smith KSPConvergenceTestFn *convtest;
64d4d2bdcSBarry Smith PetscCtxDestroyFn *convdestroy;
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
SNESNewtonTRSetTolerances_TR(SNES snes,PetscReal delta_min,PetscReal delta_max,PetscReal delta_0)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
SNESNewtonTRGetTolerances_TR(SNES snes,PetscReal * delta_min,PetscReal * delta_max,PetscReal * delta_0)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
SNESComputeJacobian_MATLMVM(SNES snes,Vec X,Mat J,Mat B,void * dummy)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
SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason * reason,void * cctx)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
SNESTR_KSPConverged_Destroy(PetscCtxRt cctx)76*2a8381b2SBarry Smith static PetscErrorCode SNESTR_KSPConverged_Destroy(PetscCtxRt cctx)
77d71ae5a4SJacob Faibussowitsch {
78*2a8381b2SBarry Smith SNES_TR_KSPConverged_Ctx *ctx = *(SNES_TR_KSPConverged_Ctx **)cctx;
79971273eeSBarry Smith
80971273eeSBarry Smith PetscFunctionBegin;
814d4d2bdcSBarry Smith PetscCall((*ctx->convdestroy)(&ctx->convctx));
829566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx));
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84971273eeSBarry Smith }
85971273eeSBarry Smith
SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason * reason,void * dummy)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 @*/
SNESNewtonTRSetNormType(SNES snes,NormType norm)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:
139ac84dfd5SBarry 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 @*/
SNESNewtonTRSetQNType(SNES snes,SNESNewtonTRQNType use)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 @*/
SNESNewtonTRSetFallbackType(SNES snes,SNESNewtonTRFallbackType ftype)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 @*/
SNESNewtonTRSetPreCheck(SNES snes,PetscErrorCode (* func)(SNES,Vec,Vec,PetscBool *,void *),PetscCtx ctx)202*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), PetscCtx 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 @*/
SNESNewtonTRGetPreCheck(SNES snes,PetscErrorCode (** func)(SNES,Vec,Vec,PetscBool *,void *),PetscCtxRt ctx)233*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), PetscCtxRt 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;
243*2a8381b2SBarry Smith if (ctx) *(void **)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 @*/
SNESNewtonTRSetPostCheck(SNES snes,PetscErrorCode (* func)(SNES,Vec,Vec,Vec,PetscBool *,PetscBool *,void *),PetscCtx ctx)266*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtx 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 @*/
SNESNewtonTRGetPostCheck(SNES snes,PetscErrorCode (** func)(SNES,Vec,Vec,Vec,PetscBool *,PetscBool *,void *),PetscCtxRt ctx)297*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtxRt 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;
307*2a8381b2SBarry Smith if (ctx) *(void **)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 @*/
SNESNewtonTRPreCheck(SNES snes,Vec X,Vec Y,PetscBool * changed_Y)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 @*/
SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool * changed_Y,PetscBool * changed_W)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 */
PetscQuadraticRoots(PetscReal a,PetscReal b,PetscReal c,PetscReal * xm,PetscReal * xp)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 */
SNESNewtonTRQuadraticDelta(SNES snes,Mat J,PetscBool has_objective,Vec Y,Vec GradF,Vec W,PetscReal * yTHy_,PetscReal * gTy_,PetscReal * deltaqm)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) */
SNESNewtonTRObjective(SNES snes,PetscBool has_objective,Vec X,Vec Y,Vec W,Vec G,PetscReal * gnorm,PetscReal * fkp1)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));
42876c63389SBarry Smith SNESCheckFunctionDomainError(snes, *gnorm);
42976c63389SBarry Smith if (has_objective) {
43076c63389SBarry Smith PetscCall(SNESComputeObjective(snes, W, fkp1));
43176c63389SBarry Smith SNESCheckObjectiveDomainError(snes, *fkp1);
43276c63389SBarry Smith } else *fkp1 = 0.5 * PetscSqr(*gnorm);
43324fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
43424fb275aSStefano Zampini }
43524fb275aSStefano Zampini
SNESSetUpQN_NEWTONTR(SNES snes)43624fb275aSStefano Zampini static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes)
43724fb275aSStefano Zampini {
43824fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
43924fb275aSStefano Zampini
44024fb275aSStefano Zampini PetscFunctionBegin;
44124fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB));
44224fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre));
44324fb275aSStefano Zampini if (tr->qn) {
44424fb275aSStefano Zampini PetscInt n, N;
44524fb275aSStefano Zampini const char *optionsprefix;
44624fb275aSStefano Zampini Mat B;
44724fb275aSStefano Zampini
44824fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B));
44924fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
45024fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_"));
45124fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix));
45224fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS));
45324fb275aSStefano Zampini PetscCall(VecGetLocalSize(snes->vec_sol, &n));
45424fb275aSStefano Zampini PetscCall(VecGetSize(snes->vec_sol, &N));
45524fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N));
45624fb275aSStefano Zampini PetscCall(MatSetUp(B));
45724fb275aSStefano Zampini PetscCall(MatSetFromOptions(B));
45824fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func));
45924fb275aSStefano Zampini tr->qnB = B;
46024fb275aSStefano Zampini if (tr->qn == SNES_TR_QN_DIFFERENT) {
46124fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B));
46224fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
46324fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_"));
46424fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix));
46524fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS));
46624fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N));
46724fb275aSStefano Zampini PetscCall(MatSetUp(B));
46824fb275aSStefano Zampini PetscCall(MatSetFromOptions(B));
46924fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func));
47024fb275aSStefano Zampini tr->qnB_pre = B;
47124fb275aSStefano Zampini } else {
47224fb275aSStefano Zampini PetscCall(PetscObjectReference((PetscObject)tr->qnB));
47324fb275aSStefano Zampini tr->qnB_pre = tr->qnB;
47424fb275aSStefano Zampini }
47524fb275aSStefano Zampini }
4767aa289d8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
4777aa289d8SStefano Zampini }
4787aa289d8SStefano Zampini
479df60cc22SBarry Smith /*
4804a221d59SStefano Zampini SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
4814a221d59SStefano Zampini (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
4824a221d59SStefano Zampini nonlinear equations
4834800dd8cSBarry Smith
4844800dd8cSBarry Smith */
SNESSolve_NEWTONTR(SNES snes)485d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
486d71ae5a4SJacob Faibussowitsch {
48704d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;
48824fb275aSStefano Zampini Vec X, F, Y, G, W, GradF, YU, Yc;
4894a221d59SStefano Zampini PetscInt maxits, lits;
4904b0a5c37SStefano Zampini PetscReal rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm;
4918b630c82SStefano Zampini PetscReal fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0;
49224fb275aSStefano Zampini PetscReal auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0;
493a0254a93SStefano Zampini PC pc;
49424fb275aSStefano Zampini Mat J, Jp;
49537d4d9b4SStefano Zampini PetscBool already_done = PETSC_FALSE, on_boundary, use_cauchy;
4967aa289d8SStefano Zampini PetscBool clear_converged_test, rho_satisfied, has_objective;
497df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx;
4985e28bcb6Sprj- void *convctx;
4996b72add0SBarry Smith SNESObjectiveFn *objective;
5004d4d2bdcSBarry Smith KSPConvergenceTestFn *convtest;
5014d4d2bdcSBarry Smith PetscCtxDestroyFn *convdestroy;
5024800dd8cSBarry Smith
5033a40ed3dSBarry Smith PetscFunctionBegin;
5044a221d59SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL));
5057aa289d8SStefano Zampini has_objective = objective ? PETSC_TRUE : PETSC_FALSE;
506c579b300SPatrick Farrell
507fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */
508fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */
50939e2f89bSBarry Smith F = snes->vec_func; /* residual vector */
5104a221d59SStefano Zampini Y = snes->vec_sol_update; /* update vector */
5114a221d59SStefano Zampini G = snes->work[0]; /* updated residual */
5124a221d59SStefano Zampini W = snes->work[1]; /* temporary vector */
5137aa289d8SStefano Zampini GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */
5144a221d59SStefano Zampini YU = snes->work[3]; /* work vector for dogleg method */
51524fb275aSStefano Zampini Yc = snes->work[4]; /* Cauchy point */
5164a221d59SStefano Zampini
5174a221d59SStefano 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);
5184800dd8cSBarry Smith
5199566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
5204c49b128SBarry Smith snes->iter = 0;
5219566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5224800dd8cSBarry Smith
52324fb275aSStefano Zampini /* setup QN matrices if needed */
52424fb275aSStefano Zampini PetscCall(SNESSetUpQN_NEWTONTR(snes));
52524fb275aSStefano Zampini
5264a221d59SStefano Zampini /* Set the linear stopping criteria to use the More' trick if needed */
5274a221d59SStefano Zampini clear_converged_test = PETSC_FALSE;
528a0254a93SStefano Zampini PetscCall(SNESGetKSP(snes, &snes->ksp));
529a0254a93SStefano Zampini PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy));
530fcc61681SStefano Zampini if (convtest != SNESTR_KSPConverged_Private) {
5314a221d59SStefano Zampini clear_converged_test = PETSC_TRUE;
5329566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx));
533df8705c3SBarry Smith ctx->snes = snes;
534a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
535a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy));
5369566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n"));
537df8705c3SBarry Smith }
538df8705c3SBarry Smith
539e4ed7901SPeter Brune if (!snes->vec_func_init_set) {
5409566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
5411aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE;
542e4ed7901SPeter Brune
5439566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
54476c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
5459566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
5464a221d59SStefano Zampini
5479566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
548fbe28522SBarry Smith snes->norm = fnorm;
5499566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5504a221d59SStefano Zampini delta = neP->delta0;
5514800dd8cSBarry Smith neP->delta = delta;
5529566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
553b37302c6SLois Curfman McInnes
55485385478SLisandro Dalcin /* test convergence */
5554a221d59SStefano Zampini rho_satisfied = PETSC_FALSE;
5562d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
5572d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm));
5583ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
5593f149594SLisandro Dalcin
56076c63389SBarry Smith if (has_objective) {
56176c63389SBarry Smith PetscCall(SNESComputeObjective(snes, X, &fk));
56276c63389SBarry Smith SNESCheckObjectiveDomainError(snes, fk);
56376c63389SBarry Smith } else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */
56499a96b7cSMatthew Knepley
565a0254a93SStefano Zampini /* hook state vector to BFGS preconditioner */
566a0254a93SStefano Zampini PetscCall(KSPGetPC(snes->ksp, &pc));
567a0254a93SStefano Zampini PetscCall(PCLMVMSetUpdateVec(pc, X));
568a0254a93SStefano Zampini
56924fb275aSStefano Zampini if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE));
5706b5873e3SBarry Smith
57124fb275aSStefano Zampini while (snes->iter < maxits) {
57212d0050eSStefano Zampini /* calculating Jacobian and GradF of minimization function only once */
5734a221d59SStefano Zampini if (!already_done) {
57412d0050eSStefano Zampini /* Call general purpose update function */
57512d0050eSStefano Zampini PetscTryTypeMethod(snes, update, snes->iter);
57612d0050eSStefano Zampini
5774b0a5c37SStefano Zampini /* apply the nonlinear preconditioner */
5784b0a5c37SStefano Zampini if (snes->npc && snes->npcside == PC_RIGHT) {
5794b0a5c37SStefano Zampini SNESConvergedReason reason;
5804b0a5c37SStefano Zampini
5814b0a5c37SStefano Zampini PetscCall(SNESSetInitialFunction(snes->npc, F));
5824b0a5c37SStefano Zampini PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
5834b0a5c37SStefano Zampini PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
5844b0a5c37SStefano Zampini PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
5854b0a5c37SStefano Zampini PetscCall(SNESGetConvergedReason(snes->npc, &reason));
5864b0a5c37SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
5874b0a5c37SStefano Zampini snes->reason = SNES_DIVERGED_INNER;
5884b0a5c37SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
5894b0a5c37SStefano Zampini }
5904b0a5c37SStefano Zampini // XXX
5914b0a5c37SStefano Zampini PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
59212d0050eSStefano Zampini }
59312d0050eSStefano Zampini
59412d0050eSStefano Zampini /* Jacobian */
59524fb275aSStefano Zampini J = NULL;
59624fb275aSStefano Zampini Jp = NULL;
59724fb275aSStefano Zampini if (!neP->qnB) {
5984a221d59SStefano Zampini PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
59976c63389SBarry Smith SNESCheckJacobianDomainError(snes);
60024fb275aSStefano Zampini J = snes->jacobian;
60124fb275aSStefano Zampini Jp = snes->jacobian_pre;
60224fb275aSStefano Zampini } else { /* QN model */
60324fb275aSStefano Zampini PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL));
60424fb275aSStefano Zampini J = neP->qnB;
60524fb275aSStefano Zampini Jp = neP->qnB_pre;
60624fb275aSStefano Zampini }
60776c63389SBarry Smith SNESCheckJacobianDomainError(snes);
60812d0050eSStefano Zampini
60924fb275aSStefano Zampini /* objective function */
61024fb275aSStefano Zampini PetscCall(VecNorm(F, NORM_2, &fnorm));
61176c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
61276c63389SBarry Smith if (has_objective) {
61376c63389SBarry Smith PetscCall(SNESComputeObjective(snes, X, &fk));
61476c63389SBarry Smith SNESCheckObjectiveDomainError(snes, fk);
61576c63389SBarry Smith } else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */
61624fb275aSStefano Zampini
61712d0050eSStefano Zampini /* GradF */
6187aa289d8SStefano Zampini if (has_objective) gfnorm = fnorm;
6197aa289d8SStefano Zampini else {
62024fb275aSStefano Zampini PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */
6217aa289d8SStefano Zampini PetscCall(VecNorm(GradF, NORM_2, &gfnorm));
622fbe28522SBarry Smith }
62324fb275aSStefano Zampini PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k));
6247aa289d8SStefano Zampini }
6257aa289d8SStefano Zampini already_done = PETSC_TRUE;
6267aa289d8SStefano Zampini
6274b0a5c37SStefano Zampini /* solve trust-region subproblem */
6284b0a5c37SStefano Zampini
62924fb275aSStefano Zampini /* first compute Cauchy Point */
63024fb275aSStefano Zampini PetscCall(MatMult(J, GradF, W));
63124fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg));
63224fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */
63324fb275aSStefano Zampini /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */
63424fb275aSStefano Zampini auk = delta / gfnorm_k;
63524fb275aSStefano Zampini if (gTBg < 0.0) tauk = 1.0;
63624fb275aSStefano Zampini else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1);
63724fb275aSStefano Zampini auk *= tauk;
63824fb275aSStefano Zampini ycnorm = auk * gfnorm;
63924fb275aSStefano Zampini PetscCall(VecAXPBY(Yc, auk, 0.0, GradF));
64024fb275aSStefano Zampini
64124fb275aSStefano Zampini on_boundary = PETSC_FALSE;
64237d4d9b4SStefano Zampini use_cauchy = (PetscBool)(tauk == 1.0 && has_objective);
64337d4d9b4SStefano Zampini if (!use_cauchy) {
64424fb275aSStefano Zampini KSPConvergedReason reason;
64524fb275aSStefano Zampini
6464b0a5c37SStefano Zampini /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods")
64724fb275aSStefano Zampini beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */
64824fb275aSStefano Zampini objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta);
649fb01098fSStefano Zampini PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin));
6504b0a5c37SStefano Zampini
65124fb275aSStefano Zampini /* specify radius if looking for Newton step and trust region norm is the l2 norm */
65224fb275aSStefano Zampini PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0));
65324fb275aSStefano Zampini PetscCall(KSPSetOperators(snes->ksp, J, Jp));
6544a221d59SStefano Zampini PetscCall(KSPSolve(snes->ksp, F, Y));
6554a221d59SStefano Zampini SNESCheckKSPSolve(snes);
6564a221d59SStefano Zampini PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
65724fb275aSStefano Zampini PetscCall(KSPGetConvergedReason(snes->ksp, &reason));
65824fb275aSStefano Zampini on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH);
6594b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
66024fb275aSStefano Zampini if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */
66124fb275aSStefano Zampini PetscReal emax, emin;
66224fb275aSStefano Zampini PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin));
66324fb275aSStefano Zampini if (emax > 0.0) beta_k = emax + 1;
66424fb275aSStefano Zampini }
66524fb275aSStefano Zampini } else { /* Cauchy point is on the boundary, accept it */
66624fb275aSStefano Zampini on_boundary = PETSC_TRUE;
66724fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y));
66824fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg));
66924fb275aSStefano Zampini }
67024fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm));
6714800dd8cSBarry Smith
6724a221d59SStefano Zampini /* decide what to do when the update is outside of trust region */
67337d4d9b4SStefano Zampini if (!use_cauchy && (ynorm > delta || ynorm == 0.0)) {
6745ec2728bSStefano Zampini SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY;
6755ec2728bSStefano Zampini
67624fb275aSStefano Zampini PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented");
6775ec2728bSStefano Zampini switch (fallback) {
6784a221d59SStefano Zampini case SNES_TR_FALLBACK_NEWTON:
6794a221d59SStefano Zampini auk = delta / ynorm;
6804a221d59SStefano Zampini PetscCall(VecScale(Y, auk));
6814b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm));
6824a221d59SStefano Zampini break;
6834a221d59SStefano Zampini case SNES_TR_FALLBACK_CAUCHY:
6844a221d59SStefano Zampini case SNES_TR_FALLBACK_DOGLEG:
6855ec2728bSStefano Zampini if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) {
68624fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y));
6874a221d59SStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg));
6884a221d59SStefano Zampini } else { /* take linear combination of Cauchy and Newton direction and step */
68924fb275aSStefano Zampini auk = gfnorm * gfnorm / gTBg;
69024fb275aSStefano Zampini if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */
69124fb275aSStefano Zampini PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF));
69224fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm));
69324fb275aSStefano Zampini } else { /* second leg */
6944a221d59SStefano Zampini PetscReal c0, c1, c2, tau = 0.0, tpos, tneg;
6954a221d59SStefano Zampini PetscBool noroots;
696284fb49fSHeeho Park
69724fb275aSStefano Zampini /* Find solutions of (Eq. 4.16 in Nocedal and Wright)
69824fb275aSStefano Zampini ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0,
69924fb275aSStefano Zampini where p_U the Cauchy direction, p_B the Newton direction */
7004a221d59SStefano Zampini PetscCall(VecAXPBY(YU, auk, 0.0, GradF));
7014a221d59SStefano Zampini PetscCall(VecAXPY(Y, -1.0, YU));
7024a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0));
7034a221d59SStefano Zampini PetscCall(VecDotRealPart(YU, Y, &c1));
7044a221d59SStefano Zampini c0 = PetscSqr(c0);
7054a221d59SStefano Zampini c2 = PetscSqr(ycnorm) - PetscSqr(delta);
70624fb275aSStefano Zampini PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos);
7074a221d59SStefano Zampini
70824fb275aSStefano Zampini /* In principle the DL strategy as a unique solution in [0,1]
70924fb275aSStefano Zampini here we check that for some reason we numerically failed
71024fb275aSStefano Zampini to compute it. In that case, we use the Cauchy point */
7114a221d59SStefano Zampini noroots = PetscIsInfOrNanReal(tneg);
71224fb275aSStefano Zampini if (!noroots) {
71324fb275aSStefano Zampini if (tpos > 1) {
71424fb275aSStefano Zampini if (tneg >= 0 && tneg <= 1) {
71524fb275aSStefano Zampini tau = tneg;
71624fb275aSStefano Zampini } else noroots = PETSC_TRUE;
71724fb275aSStefano Zampini } else if (tpos >= 0) {
71824fb275aSStefano Zampini tau = tpos;
71924fb275aSStefano Zampini } else noroots = PETSC_TRUE;
72024fb275aSStefano Zampini }
7214a221d59SStefano Zampini if (noroots) { /* No roots, select Cauchy point */
72224fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y));
7234a221d59SStefano Zampini } else {
72424fb275aSStefano Zampini PetscCall(VecAXPBY(Y, 1.0, tau, YU));
7254a221d59SStefano Zampini }
72624fb275aSStefano 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));
7274a221d59SStefano Zampini }
7284a221d59SStefano Zampini }
7294a221d59SStefano Zampini break;
7304a221d59SStefano Zampini default:
7314a221d59SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode");
732454a90a3SBarry Smith break;
73352392280SLois Curfman McInnes }
7344800dd8cSBarry Smith }
7354a221d59SStefano Zampini
7367aa289d8SStefano Zampini /* compute the quadratic model difference */
73724fb275aSStefano Zampini PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm));
7384a221d59SStefano Zampini
7394a221d59SStefano Zampini /* Compute new objective function */
74024fb275aSStefano Zampini PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1));
74124fb275aSStefano Zampini if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1;
74224fb275aSStefano Zampini else {
7434a221d59SStefano Zampini if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */
74424fb275aSStefano Zampini else rho = neP->eta1; /* no reduction in quadratic model, step must be rejected */
74524fb275aSStefano Zampini }
7464a221d59SStefano Zampini
74724fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm));
74824fb275aSStefano 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));
74924fb275aSStefano Zampini
75024fb275aSStefano Zampini /* update the size of the trust region */
7514a221d59SStefano Zampini if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */
75224fb275aSStefano Zampini else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */
7538b630c82SStefano Zampini delta = PetscMin(delta, neP->deltaM); /* but not greater than deltaM */
7544a221d59SStefano Zampini
75524fb275aSStefano Zampini /* log 2-norm of update for moniroting routines */
75624fb275aSStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm));
75724fb275aSStefano Zampini
75824fb275aSStefano Zampini /* decide on new step */
7594a221d59SStefano Zampini neP->delta = delta;
76024fb275aSStefano Zampini if (rho > neP->eta1) {
7614a221d59SStefano Zampini rho_satisfied = PETSC_TRUE;
7624a221d59SStefano Zampini } else {
7634a221d59SStefano Zampini rho_satisfied = PETSC_FALSE;
7644a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
7654a221d59SStefano Zampini /* check to see if progress is hopeless */
7664a221d59SStefano Zampini PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP));
7672d157150SStefano Zampini if (!snes->reason) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
7684b0a5c37SStefano Zampini if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_TR_DELTA;
7694a221d59SStefano Zampini snes->numFailures++;
7704a221d59SStefano Zampini /* We're not progressing, so return with the current iterate */
7714a221d59SStefano Zampini if (snes->reason) break;
7724a221d59SStefano Zampini }
7734a221d59SStefano Zampini if (rho_satisfied) {
7744a221d59SStefano Zampini /* Update function values */
7754a221d59SStefano Zampini already_done = PETSC_FALSE;
7764800dd8cSBarry Smith fnorm = gnorm;
7774a221d59SStefano Zampini fk = fkp1;
7784a221d59SStefano Zampini
7794a221d59SStefano Zampini /* New residual and linearization point */
7809566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F));
7819566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X));
7824a221d59SStefano Zampini
78385385478SLisandro Dalcin /* Monitor convergence */
7849566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
7854a221d59SStefano Zampini snes->iter++;
786fbe28522SBarry Smith snes->norm = fnorm;
787c1e67a49SFande Kong snes->xnorm = xnorm;
788c1e67a49SFande Kong snes->ynorm = ynorm;
7899566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7909566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
7914a221d59SStefano Zampini
79285385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */
7934a221d59SStefano Zampini PetscCall(VecNorm(X, NORM_2, &xnorm));
7942d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
7952d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
7964a221d59SStefano Zampini if (snes->reason) break;
7974a221d59SStefano Zampini }
79838442cffSBarry Smith }
799284fb49fSHeeho Park
8004a221d59SStefano Zampini if (clear_converged_test) {
801a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
8029566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx));
803a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy));
8045e28bcb6Sprj- }
8053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8064800dd8cSBarry Smith }
807284fb49fSHeeho Park
SNESSetUp_NEWTONTR(SNES snes)808d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
809d71ae5a4SJacob Faibussowitsch {
8103a40ed3dSBarry Smith PetscFunctionBegin;
81124fb275aSStefano Zampini PetscCall(SNESSetWorkVecs(snes, 5));
8129566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes));
8133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8144800dd8cSBarry Smith }
8156b8b9a38SLisandro Dalcin
SNESReset_NEWTONTR(SNES snes)81624fb275aSStefano Zampini static PetscErrorCode SNESReset_NEWTONTR(SNES snes)
81724fb275aSStefano Zampini {
81824fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
81924fb275aSStefano Zampini
82024fb275aSStefano Zampini PetscFunctionBegin;
82124fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB));
82224fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre));
82324fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
82424fb275aSStefano Zampini }
82524fb275aSStefano Zampini
SNESDestroy_NEWTONTR(SNES snes)826d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
827d71ae5a4SJacob Faibussowitsch {
8283a40ed3dSBarry Smith PetscFunctionBegin;
82924fb275aSStefano Zampini PetscCall(SNESReset_NEWTONTR(snes));
8303201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL));
8313201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", NULL));
8329566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data));
8333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8344800dd8cSBarry Smith }
8354800dd8cSBarry Smith
SNESSetFromOptions_NEWTONTR(SNES snes,PetscOptionItems PetscOptionsObject)836ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems PetscOptionsObject)
837d71ae5a4SJacob Faibussowitsch {
83804d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data;
83924fb275aSStefano Zampini SNESNewtonTRQNType qn;
84024fb275aSStefano Zampini SNESNewtonTRFallbackType fallback;
84124fb275aSStefano Zampini NormType norm;
84224fb275aSStefano Zampini PetscBool flg;
8434800dd8cSBarry Smith
8443a40ed3dSBarry Smith PetscFunctionBegin;
845d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
8463201ab8dSStefano Zampini PetscCall(PetscOptionsDeprecated("-snes_tr_deltaM", "-snes_tr_deltamax", "3.22", NULL));
8473201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "SNESNewtonTRSetUpdateParameters", ctx->eta1, &ctx->eta1, NULL));
8483201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "SNESNewtonTRSetUpdateParameters", ctx->eta2, &ctx->eta2, NULL));
8493201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "SNESNewtonTRSetUpdateParameters", ctx->eta3, &ctx->eta3, NULL));
8503201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "SNESNewtonTRSetUpdateParameters", ctx->t1, &ctx->t1, NULL));
8513201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "SNESNewtonTRSetUpdateParameters", ctx->t2, &ctx->t2, NULL));
8523201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_delta0", "Initial trust region size", "SNESNewtonTRSetTolerances", ctx->delta0, &ctx->delta0, NULL));
8533201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltamin", "Minimum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltam, &ctx->deltam, NULL));
8543201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltamax", "Maximum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltaM, &ctx->deltaM, NULL));
8554b0a5c37SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL));
85624fb275aSStefano Zampini
85724fb275aSStefano Zampini fallback = ctx->fallback;
85824fb275aSStefano 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));
85924fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback));
86024fb275aSStefano Zampini
86124fb275aSStefano Zampini qn = ctx->qn;
86224fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg));
86324fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn));
86424fb275aSStefano Zampini
86524fb275aSStefano Zampini norm = ctx->norm;
86624fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg));
86724fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm));
86824fb275aSStefano Zampini
869d0609cedSBarry Smith PetscOptionsHeadEnd();
8703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8714800dd8cSBarry Smith }
8724800dd8cSBarry Smith
SNESView_NEWTONTR(SNES snes,PetscViewer viewer)873d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer)
874d71ae5a4SJacob Faibussowitsch {
87504d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
8769f196a02SMartin Diehl PetscBool isascii;
877a935fc98SLois Curfman McInnes
8783a40ed3dSBarry Smith PetscFunctionBegin;
8799f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
8809f196a02SMartin Diehl if (isascii) {
8813201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region parameters:\n"));
8824a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
8833201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " t1=%g, t2=%g\n", (double)tr->t1, (double)tr->t2));
8843201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " delta_min=%g, delta_0=%g, delta_max=%g\n", (double)tr->deltam, (double)tr->delta0, (double)tr->deltaM));
8854b0a5c37SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " kmdc=%g\n", (double)tr->kmdc));
8864a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback]));
88724fb275aSStefano Zampini if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, " qn=%s\n", SNESNewtonTRQNTypes[tr->qn]));
88824fb275aSStefano Zampini if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, " norm=%s\n", NormTypes[tr->norm]));
88919bcc07fSBarry Smith }
8903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
891a935fc98SLois Curfman McInnes }
892f6dfbefdSBarry Smith
8933201ab8dSStefano Zampini /*@
8943201ab8dSStefano Zampini SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
8953201ab8dSStefano Zampini
8963201ab8dSStefano Zampini Logically Collective
8973201ab8dSStefano Zampini
8983201ab8dSStefano Zampini Input Parameters:
8993201ab8dSStefano Zampini + snes - the `SNES` context
9003201ab8dSStefano Zampini - tol - tolerance
9013201ab8dSStefano Zampini
9023201ab8dSStefano Zampini Level: deprecated
9033201ab8dSStefano Zampini
9043201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetTolerances()`
9053201ab8dSStefano Zampini @*/
SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)9063201ab8dSStefano Zampini PetscErrorCode SNESSetTrustRegionTolerance(SNES snes, PetscReal tol)
9073201ab8dSStefano Zampini {
9083201ab8dSStefano Zampini return SNESNewtonTRSetTolerances(snes, tol, PETSC_CURRENT, PETSC_CURRENT);
9093201ab8dSStefano Zampini }
9103201ab8dSStefano Zampini
9113201ab8dSStefano Zampini /*@
9123201ab8dSStefano Zampini SNESNewtonTRSetTolerances - Sets the trust region parameter tolerances.
9133201ab8dSStefano Zampini
9143201ab8dSStefano Zampini Logically Collective
9153201ab8dSStefano Zampini
9163201ab8dSStefano Zampini Input Parameters:
9173201ab8dSStefano Zampini + snes - the `SNES` context
9183201ab8dSStefano Zampini . delta_min - minimum allowed trust region size
9193201ab8dSStefano Zampini . delta_max - maximum allowed trust region size
9203201ab8dSStefano Zampini - delta_0 - initial trust region size
9213201ab8dSStefano Zampini
9223201ab8dSStefano Zampini Options Database Key:
9233201ab8dSStefano Zampini + -snes_tr_deltamin <tol> - Set minimum size
9243201ab8dSStefano Zampini . -snes_tr_deltamax <tol> - Set maximum size
9253201ab8dSStefano Zampini - -snes_tr_delta0 <tol> - Set initial size
9263201ab8dSStefano Zampini
9273201ab8dSStefano Zampini Note:
9283201ab8dSStefano Zampini Use `PETSC_DETERMINE` to use the default value for the given `SNES`.
9293201ab8dSStefano Zampini Use `PETSC_CURRENT` to retain a value.
9303201ab8dSStefano Zampini
9313201ab8dSStefano Zampini Fortran Note:
9323201ab8dSStefano Zampini Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL`
9333201ab8dSStefano Zampini
9343201ab8dSStefano Zampini Level: intermediate
9353201ab8dSStefano Zampini
9363201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRGetTolerances()`
9373201ab8dSStefano Zampini @*/
SNESNewtonTRSetTolerances(SNES snes,PetscReal delta_min,PetscReal delta_max,PetscReal delta_0)9383201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRSetTolerances(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0)
9393201ab8dSStefano Zampini {
9403201ab8dSStefano Zampini PetscFunctionBegin;
9413201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
9423201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, delta_min, 2);
9433201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, delta_max, 3);
9443201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, delta_0, 4);
9453201ab8dSStefano Zampini PetscTryMethod(snes, "SNESNewtonTRSetTolerances_C", (SNES, PetscReal, PetscReal, PetscReal), (snes, delta_min, delta_max, delta_0));
9463201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
9473201ab8dSStefano Zampini }
9483201ab8dSStefano Zampini
9493201ab8dSStefano Zampini /*@
9503201ab8dSStefano Zampini SNESNewtonTRGetTolerances - Gets the trust region parameter tolerances.
9513201ab8dSStefano Zampini
9523201ab8dSStefano Zampini Not Collective
9533201ab8dSStefano Zampini
9543201ab8dSStefano Zampini Input Parameter:
9553201ab8dSStefano Zampini . snes - the `SNES` context
9563201ab8dSStefano Zampini
9573201ab8dSStefano Zampini Output Parameters:
9583201ab8dSStefano Zampini + delta_min - minimum allowed trust region size or `NULL`
9593201ab8dSStefano Zampini . delta_max - maximum allowed trust region size or `NULL`
9603201ab8dSStefano Zampini - delta_0 - initial trust region size or `NULL`
9613201ab8dSStefano Zampini
9623201ab8dSStefano Zampini Level: intermediate
9633201ab8dSStefano Zampini
9643201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetTolerances()`
9653201ab8dSStefano Zampini @*/
SNESNewtonTRGetTolerances(SNES snes,PetscReal * delta_min,PetscReal * delta_max,PetscReal * delta_0)9663201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRGetTolerances(SNES snes, PetscReal *delta_min, PetscReal *delta_max, PetscReal *delta_0)
9673201ab8dSStefano Zampini {
9683201ab8dSStefano Zampini PetscFunctionBegin;
9693201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
9703201ab8dSStefano Zampini if (delta_min) PetscAssertPointer(delta_min, 2);
9713201ab8dSStefano Zampini if (delta_max) PetscAssertPointer(delta_max, 3);
9723201ab8dSStefano Zampini if (delta_0) PetscAssertPointer(delta_0, 4);
9733201ab8dSStefano Zampini PetscUseMethod(snes, "SNESNewtonTRGetTolerances_C", (SNES, PetscReal *, PetscReal *, PetscReal *), (snes, delta_min, delta_max, delta_0));
9743201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
9753201ab8dSStefano Zampini }
9763201ab8dSStefano Zampini
9773201ab8dSStefano Zampini /*@
9783201ab8dSStefano Zampini SNESNewtonTRSetUpdateParameters - Sets the trust region update parameters.
9793201ab8dSStefano Zampini
9803201ab8dSStefano Zampini Logically Collective
9813201ab8dSStefano Zampini
9823201ab8dSStefano Zampini Input Parameters:
9833201ab8dSStefano Zampini + snes - the `SNES` context
9843201ab8dSStefano Zampini . eta1 - acceptance tolerance
9853201ab8dSStefano Zampini . eta2 - shrinking tolerance
9863201ab8dSStefano Zampini . eta3 - enlarging tolerance
9873201ab8dSStefano Zampini . t1 - shrink factor
9883201ab8dSStefano Zampini - t2 - enlarge factor
9893201ab8dSStefano Zampini
9903201ab8dSStefano Zampini Options Database Key:
991ac84dfd5SBarry Smith + -snes_tr_eta1 <tol> - Set `eta1`
992ac84dfd5SBarry Smith . -snes_tr_eta2 <tol> - Set `eta2`
993ac84dfd5SBarry Smith . -snes_tr_eta3 <tol> - Set `eta3`
994ac84dfd5SBarry Smith . -snes_tr_t1 <tol> - Set `t1`
995ac84dfd5SBarry Smith - -snes_tr_t2 <tol> - Set `t2`
9963201ab8dSStefano Zampini
9973201ab8dSStefano Zampini Notes:
9983201ab8dSStefano 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,
9993201ab8dSStefano Zampini $s_k$ the computed step, $f$ the objective function, and $m$ the quadratic model, the trust region
10003201ab8dSStefano Zampini radius is modified as follows
10016b1535e8SBarry Smith
10026b1535e8SBarry Smith $
10033201ab8dSStefano Zampini \delta =
10043201ab8dSStefano Zampini \begin{cases}
10053201ab8dSStefano Zampini \delta * t_1 ,& \rho < \eta_2 \\
10063201ab8dSStefano Zampini \delta * t_2 ,& \rho > \eta_3 \\
10073201ab8dSStefano Zampini \end{cases}
10086b1535e8SBarry Smith $
10096b1535e8SBarry Smith
10103201ab8dSStefano Zampini The step is accepted if $\rho > \eta_1$.
10113201ab8dSStefano Zampini Use `PETSC_DETERMINE` to use the default value for the given `SNES`.
10123201ab8dSStefano Zampini Use `PETSC_CURRENT` to retain a value.
10133201ab8dSStefano Zampini
10143201ab8dSStefano Zampini Fortran Note:
10153201ab8dSStefano Zampini Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL`
10163201ab8dSStefano Zampini
10173201ab8dSStefano Zampini Level: intermediate
10183201ab8dSStefano Zampini
10193201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetObjective()`, `SNESNewtonTRGetUpdateParameters()`
10203201ab8dSStefano Zampini @*/
SNESNewtonTRSetUpdateParameters(SNES snes,PetscReal eta1,PetscReal eta2,PetscReal eta3,PetscReal t1,PetscReal t2)10213201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRSetUpdateParameters(SNES snes, PetscReal eta1, PetscReal eta2, PetscReal eta3, PetscReal t1, PetscReal t2)
10223201ab8dSStefano Zampini {
10233201ab8dSStefano Zampini PetscBool flg;
10243201ab8dSStefano Zampini
10253201ab8dSStefano Zampini PetscFunctionBegin;
10263201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
10273201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, eta1, 2);
10283201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, eta2, 3);
10293201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, eta3, 4);
10303201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, t1, 5);
10313201ab8dSStefano Zampini PetscValidLogicalCollectiveReal(snes, t2, 6);
10323201ab8dSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
10333201ab8dSStefano Zampini if (flg) {
10343201ab8dSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
10353201ab8dSStefano Zampini
10363201ab8dSStefano Zampini if (eta1 == PETSC_DETERMINE) eta1 = tr->default_eta1;
10373201ab8dSStefano Zampini if (eta2 == PETSC_DETERMINE) eta2 = tr->default_eta2;
10383201ab8dSStefano Zampini if (eta3 == PETSC_DETERMINE) eta3 = tr->default_eta3;
10393201ab8dSStefano Zampini if (t1 == PETSC_DETERMINE) t1 = tr->default_t1;
10403201ab8dSStefano Zampini if (t2 == PETSC_DETERMINE) t2 = tr->default_t2;
10413201ab8dSStefano Zampini if (eta1 != PETSC_CURRENT) tr->eta1 = eta1;
10423201ab8dSStefano Zampini if (eta2 != PETSC_CURRENT) tr->eta2 = eta2;
10433201ab8dSStefano Zampini if (eta3 != PETSC_CURRENT) tr->eta3 = eta3;
10443201ab8dSStefano Zampini if (t1 != PETSC_CURRENT) tr->t1 = t1;
10453201ab8dSStefano Zampini if (t2 != PETSC_CURRENT) tr->t2 = t2;
10463201ab8dSStefano Zampini }
10473201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
10483201ab8dSStefano Zampini }
10493201ab8dSStefano Zampini
10503201ab8dSStefano Zampini /*@
10513201ab8dSStefano Zampini SNESNewtonTRGetUpdateParameters - Gets the trust region update parameters.
10523201ab8dSStefano Zampini
10533201ab8dSStefano Zampini Not Collective
10543201ab8dSStefano Zampini
10553201ab8dSStefano Zampini Input Parameter:
10563201ab8dSStefano Zampini . snes - the `SNES` context
10573201ab8dSStefano Zampini
10583201ab8dSStefano Zampini Output Parameters:
10593201ab8dSStefano Zampini + eta1 - acceptance tolerance
10603201ab8dSStefano Zampini . eta2 - shrinking tolerance
10613201ab8dSStefano Zampini . eta3 - enlarging tolerance
10623201ab8dSStefano Zampini . t1 - shrink factor
10633201ab8dSStefano Zampini - t2 - enlarge factor
10643201ab8dSStefano Zampini
10653201ab8dSStefano Zampini Level: intermediate
10663201ab8dSStefano Zampini
10673201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetUpdateParameters()`
10683201ab8dSStefano Zampini @*/
SNESNewtonTRGetUpdateParameters(SNES snes,PetscReal * eta1,PetscReal * eta2,PetscReal * eta3,PetscReal * t1,PetscReal * t2)10693201ab8dSStefano Zampini PetscErrorCode SNESNewtonTRGetUpdateParameters(SNES snes, PetscReal *eta1, PetscReal *eta2, PetscReal *eta3, PetscReal *t1, PetscReal *t2)
10703201ab8dSStefano Zampini {
10713201ab8dSStefano Zampini SNES_NEWTONTR *tr;
10723201ab8dSStefano Zampini PetscBool flg;
10733201ab8dSStefano Zampini
10743201ab8dSStefano Zampini PetscFunctionBegin;
10753201ab8dSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
10763201ab8dSStefano Zampini if (eta1) PetscAssertPointer(eta1, 2);
10773201ab8dSStefano Zampini if (eta2) PetscAssertPointer(eta2, 3);
10783201ab8dSStefano Zampini if (eta3) PetscAssertPointer(eta3, 4);
10793201ab8dSStefano Zampini if (t1) PetscAssertPointer(t1, 5);
10803201ab8dSStefano Zampini if (t2) PetscAssertPointer(t2, 6);
10813201ab8dSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
10823201ab8dSStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
10833201ab8dSStefano Zampini tr = (SNES_NEWTONTR *)snes->data;
10843201ab8dSStefano Zampini if (eta1) *eta1 = tr->eta1;
10853201ab8dSStefano Zampini if (eta2) *eta2 = tr->eta2;
10863201ab8dSStefano Zampini if (eta3) *eta3 = tr->eta3;
10873201ab8dSStefano Zampini if (t1) *t1 = tr->t1;
10883201ab8dSStefano Zampini if (t2) *t2 = tr->t2;
10893201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
10903201ab8dSStefano Zampini }
10913201ab8dSStefano Zampini
1092ebe3b25bSBarry Smith /*MC
10933201ab8dSStefano Zampini SNESNEWTONTR - Newton based nonlinear solver that uses a trust-region strategy
1094f6dfbefdSBarry Smith
1095f6dfbefdSBarry Smith Options Database Keys:
10963201ab8dSStefano Zampini + -snes_tr_deltamin <deltamin> - trust region parameter, minimum size of trust region
10973201ab8dSStefano Zampini . -snes_tr_deltamax <deltamax> - trust region parameter, max size of trust region (default: 1e10)
10983201ab8dSStefano Zampini . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2)
1099ac84dfd5SBarry Smith . -snes_tr_eta1 <eta1> - trust region parameter $eta1 \le eta2$, $rho > eta1$ breaks out of the inner iteration (default: 0.001)
1100ac84dfd5SBarry Smith . -snes_tr_eta2 <eta2> - trust region parameter, $rho \le eta2$ shrinks the trust region (default: 0.25)
1101ac84dfd5SBarry Smith . -snes_tr_eta3 <eta3> - trust region parameter $eta3 > eta2$, $rho \ge eta3$ expands the trust region (default: 0.75)
11024a221d59SStefano Zampini . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
11034a221d59SStefano Zampini . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
11046b1535e8SBarry Smith . -snes_tr_norm_type <1,2,infinity> - Type of norm for trust region bounds (default: 2)
11054a221d59SStefano 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.
1106acba1f63SStefano Zampini
1107acba1f63SStefano Zampini Level: beginner
1108b3113221SBarry Smith
11093201ab8dSStefano Zampini Notes:
11103201ab8dSStefano Zampini The code is largely based on the book {cite}`nocedal2006numerical` and supports minimizing objective functions using a quadratic model.
11113201ab8dSStefano Zampini Quasi-Newton models are also supported.
11123201ab8dSStefano Zampini
11133201ab8dSStefano Zampini Default step computation uses the Newton direction, but a dogleg type update is also supported.
11143201ab8dSStefano Zampini The 1- and infinity-norms are also supported when computing the trust region bounds.
11153201ab8dSStefano Zampini
11163201ab8dSStefano Zampini .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetObjective()`,
11173201ab8dSStefano Zampini `SNESNewtonTRSetTolerances()`, `SNESNewtonTRSetUpdateParameters()`
11183201ab8dSStefano Zampini `SNESNewtonTRSetNormType()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()`
11193201ab8dSStefano Zampini `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRSetPreCheck()`,
1120ebe3b25bSBarry Smith M*/
SNESCreate_NEWTONTR(SNES snes)1121d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
1122d71ae5a4SJacob Faibussowitsch {
112304d7464bSBarry Smith SNES_NEWTONTR *neP;
11244800dd8cSBarry Smith
11253a40ed3dSBarry Smith PetscFunctionBegin;
112604d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR;
112704d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR;
112824fb275aSStefano Zampini snes->ops->reset = SNESReset_NEWTONTR;
112904d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR;
113004d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
113104d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR;
1132fbe28522SBarry Smith
113377e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes));
113437d4d9b4SStefano Zampini PetscObjectParameterSetDefault(snes, stol, 0.0);
113537d4d9b4SStefano Zampini
113642f4f86dSBarry Smith snes->usesksp = PETSC_TRUE;
11374b0a5c37SStefano Zampini snes->npcside = PC_RIGHT;
11384b0a5c37SStefano Zampini snes->usesnpc = PETSC_TRUE;
113942f4f86dSBarry Smith
11404fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE;
11414fc747eaSLawrence Mitchell
11424dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP));
1143fbe28522SBarry Smith snes->data = (void *)neP;
11443201ab8dSStefano Zampini
11453201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, eta1, 0.001);
11463201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, eta2, 0.25);
11473201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, eta3, 0.75);
11483201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, t1, 0.25);
11493201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, t2, 2.0);
11503201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, deltam, PetscDefined(USE_REAL_SINGLE) ? 1.e-6 : 1.e-12);
11513201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, delta0, 0.2);
11523201ab8dSStefano Zampini PetscObjectParameterSetDefault(neP, deltaM, 1.e10);
11533201ab8dSStefano Zampini
115424fb275aSStefano Zampini neP->norm = NORM_2;
11554a221d59SStefano Zampini neP->fallback = SNES_TR_FALLBACK_NEWTON;
11564b0a5c37SStefano Zampini neP->kmdc = 0.0; /* by default do not use sufficient decrease */
11573201ab8dSStefano Zampini
11583201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TR));
11593201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", SNESNewtonTRGetTolerances_TR));
11603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11614800dd8cSBarry Smith }
1162