xref: /petsc/src/snes/impls/ntrdc/ntrdc.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
141ba4c6cSHeeho Park #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I   "petscsnes.h"   I*/
241ba4c6cSHeeho Park 
341ba4c6cSHeeho Park typedef struct {
441ba4c6cSHeeho Park   SNES snes;
541ba4c6cSHeeho Park   /*  Information on the regular SNES convergence test; which may have been user provided
641ba4c6cSHeeho Park       Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho
741ba4c6cSHeeho Park       Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private
841ba4c6cSHeeho Park  */
941ba4c6cSHeeho Park 
104d4d2bdcSBarry Smith   KSPConvergenceTestFn *convtest;
114d4d2bdcSBarry Smith   PetscCtxDestroyFn    *convdestroy;
1241ba4c6cSHeeho Park   void                 *convctx;
1341ba4c6cSHeeho Park } SNES_TRDC_KSPConverged_Ctx;
1441ba4c6cSHeeho Park 
SNESNewtonTRSetTolerances_TRDC(SNES snes,PetscReal delta_min,PetscReal delta_max,PetscReal delta_0)153201ab8dSStefano Zampini static PetscErrorCode SNESNewtonTRSetTolerances_TRDC(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0)
163201ab8dSStefano Zampini {
173201ab8dSStefano Zampini   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
183201ab8dSStefano Zampini 
193201ab8dSStefano Zampini   PetscFunctionBegin;
203201ab8dSStefano Zampini   if (delta_min == PETSC_DETERMINE) delta_min = 1.e-12;
213201ab8dSStefano Zampini   if (delta_max == PETSC_DETERMINE) delta_max = 0.5;
223201ab8dSStefano Zampini   if (delta_0 == PETSC_DETERMINE) delta_0 = 0.1;
233201ab8dSStefano Zampini   if (delta_min != PETSC_CURRENT) tr->deltatol = delta_min;
243201ab8dSStefano Zampini   if (delta_max != PETSC_CURRENT) tr->deltaM = delta_max;
253201ab8dSStefano Zampini   if (delta_0 != PETSC_CURRENT) tr->delta0 = delta_0;
263201ab8dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
273201ab8dSStefano Zampini }
283201ab8dSStefano Zampini 
SNESTRDC_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason * reason,void * cctx)29d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx)
30d71ae5a4SJacob Faibussowitsch {
3141ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx  = (SNES_TRDC_KSPConverged_Ctx *)cctx;
3241ba4c6cSHeeho Park   SNES                        snes = ctx->snes;
3341ba4c6cSHeeho Park   SNES_NEWTONTRDC            *neP  = (SNES_NEWTONTRDC *)snes->data;
3441ba4c6cSHeeho Park   Vec                         x;
3541ba4c6cSHeeho Park   PetscReal                   nrm;
3641ba4c6cSHeeho Park 
3741ba4c6cSHeeho Park   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
3948a46eb9SPierre Jolivet   if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm));
4041ba4c6cSHeeho Park   /* Determine norm of solution */
419566063dSJacob Faibussowitsch   PetscCall(KSPBuildSolution(ksp, NULL, &x));
429566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &nrm));
4341ba4c6cSHeeho Park   if (nrm >= neP->delta) {
449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
4541ba4c6cSHeeho Park     *reason = KSP_CONVERGED_STEP_LENGTH;
4641ba4c6cSHeeho Park   }
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4841ba4c6cSHeeho Park }
4941ba4c6cSHeeho Park 
SNESTRDC_KSPConverged_Destroy(PetscCtxRt cctx)50*2a8381b2SBarry Smith static PetscErrorCode SNESTRDC_KSPConverged_Destroy(PetscCtxRt cctx)
51d71ae5a4SJacob Faibussowitsch {
52*2a8381b2SBarry Smith   SNES_TRDC_KSPConverged_Ctx *ctx = *(SNES_TRDC_KSPConverged_Ctx **)cctx;
5341ba4c6cSHeeho Park 
5441ba4c6cSHeeho Park   PetscFunctionBegin;
554d4d2bdcSBarry Smith   PetscCall((*ctx->convdestroy)(&ctx->convctx));
569566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5841ba4c6cSHeeho Park }
5941ba4c6cSHeeho Park 
6041ba4c6cSHeeho Park /*
614d4d2bdcSBarry Smith    SNESTRDC_Converged_Private -test convergence JUST for the trust region tolerance.
6241ba4c6cSHeeho Park */
SNESTRDC_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason * reason,void * dummy)63d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
64d71ae5a4SJacob Faibussowitsch {
6541ba4c6cSHeeho Park   SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
6641ba4c6cSHeeho Park 
6741ba4c6cSHeeho Park   PetscFunctionBegin;
6841ba4c6cSHeeho Park   *reason = SNES_CONVERGED_ITERATING;
693201ab8dSStefano Zampini   if (neP->delta < xnorm * neP->deltatol) {
703201ab8dSStefano Zampini     PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)neP->deltatol));
7141ba4c6cSHeeho Park     *reason = SNES_DIVERGED_TR_DELTA;
7241ba4c6cSHeeho Park   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
739566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
7441ba4c6cSHeeho Park     *reason = SNES_DIVERGED_FUNCTION_COUNT;
7541ba4c6cSHeeho Park   }
763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7741ba4c6cSHeeho Park }
7841ba4c6cSHeeho Park 
7941ba4c6cSHeeho Park /*@
80f6dfbefdSBarry Smith   SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region.
8141ba4c6cSHeeho Park 
8220f4b53cSBarry Smith   Logically Collective
8320f4b53cSBarry Smith 
84f6dfbefdSBarry Smith   Input Parameter:
8541ba4c6cSHeeho Park . snes - the nonlinear solver object
8641ba4c6cSHeeho Park 
87f6dfbefdSBarry Smith   Output Parameter:
88420bcc1bSBarry Smith . rho_flag - `PETSC_FALSE` or `PETSC_TRUE`
8941ba4c6cSHeeho Park 
9041ba4c6cSHeeho Park   Level: developer
9141ba4c6cSHeeho Park 
92420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPreCheck()`,
93f6dfbefdSBarry Smith           `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
9441ba4c6cSHeeho Park @*/
SNESNewtonTRDCGetRhoFlag(SNES snes,PetscBool * rho_flag)95d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag)
96d71ae5a4SJacob Faibussowitsch {
9741ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
9841ba4c6cSHeeho Park 
9941ba4c6cSHeeho Park   PetscFunctionBegin;
10041ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1014f572ea9SToby Isaac   PetscAssertPointer(rho_flag, 2);
10241ba4c6cSHeeho Park   *rho_flag = tr->rho_satisfied;
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10441ba4c6cSHeeho Park }
10541ba4c6cSHeeho Park 
10641ba4c6cSHeeho Park /*@C
10741ba4c6cSHeeho Park   SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined.
10841ba4c6cSHeeho Park   Allows the user a chance to change or override the trust region decision.
10941ba4c6cSHeeho Park 
110c3339decSBarry Smith   Logically Collective
11141ba4c6cSHeeho Park 
11241ba4c6cSHeeho Park   Input Parameters:
11341ba4c6cSHeeho Park + snes - the nonlinear solver object
11420f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()`
11520f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
11641ba4c6cSHeeho Park 
11741ba4c6cSHeeho Park   Level: intermediate
11841ba4c6cSHeeho Park 
119f6dfbefdSBarry Smith   Note:
120f6dfbefdSBarry Smith   This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver.
12141ba4c6cSHeeho Park 
122420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`,
123f6dfbefdSBarry Smith           `SNESNewtonTRDCGetRhoFlag()`
12441ba4c6cSHeeho Park @*/
SNESNewtonTRDCSetPreCheck(SNES snes,PetscErrorCode (* func)(SNES,Vec,Vec,PetscBool *,void *),PetscCtx ctx)125*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), PetscCtx ctx)
126d71ae5a4SJacob Faibussowitsch {
12741ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
12841ba4c6cSHeeho Park 
12941ba4c6cSHeeho Park   PetscFunctionBegin;
13041ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
13141ba4c6cSHeeho Park   if (func) tr->precheck = func;
13241ba4c6cSHeeho Park   if (ctx) tr->precheckctx = ctx;
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13441ba4c6cSHeeho Park }
13541ba4c6cSHeeho Park 
13641ba4c6cSHeeho Park /*@C
137420bcc1bSBarry Smith   SNESNewtonTRDCGetPreCheck - Gets the pre-check function optionally set with `SNESNewtonTRDCSetPreCheck()`
13841ba4c6cSHeeho Park 
13920f4b53cSBarry Smith   Not Collective
14041ba4c6cSHeeho Park 
14141ba4c6cSHeeho Park   Input Parameter:
14241ba4c6cSHeeho Park . snes - the nonlinear solver context
14341ba4c6cSHeeho Park 
14441ba4c6cSHeeho Park   Output Parameters:
14520f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()`
14620f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
14741ba4c6cSHeeho Park 
14841ba4c6cSHeeho Park   Level: intermediate
14941ba4c6cSHeeho Park 
150420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()`
15141ba4c6cSHeeho Park @*/
SNESNewtonTRDCGetPreCheck(SNES snes,PetscErrorCode (** func)(SNES,Vec,Vec,PetscBool *,void *),PetscCtxRt ctx)152*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), PetscCtxRt ctx)
153d71ae5a4SJacob Faibussowitsch {
15441ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
15541ba4c6cSHeeho Park 
15641ba4c6cSHeeho Park   PetscFunctionBegin;
15741ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
15841ba4c6cSHeeho Park   if (func) *func = tr->precheck;
159*2a8381b2SBarry Smith   if (ctx) *(void **)ctx = tr->precheckctx;
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16141ba4c6cSHeeho Park }
16241ba4c6cSHeeho Park 
16341ba4c6cSHeeho Park /*@C
16441ba4c6cSHeeho Park   SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
16541ba4c6cSHeeho Park   function evaluation. Allows the user a chance to change or override the decision of the line search routine
16641ba4c6cSHeeho Park 
167c3339decSBarry Smith   Logically Collective
16841ba4c6cSHeeho Park 
16941ba4c6cSHeeho Park   Input Parameters:
17041ba4c6cSHeeho Park + snes - the nonlinear solver object
17120f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()`
17220f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
17341ba4c6cSHeeho Park 
17441ba4c6cSHeeho Park   Level: intermediate
17541ba4c6cSHeeho Park 
176f6dfbefdSBarry Smith   Note:
177f6dfbefdSBarry Smith   This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in
178f6dfbefdSBarry Smith   `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation.
17941ba4c6cSHeeho Park 
180420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`
18141ba4c6cSHeeho Park @*/
SNESNewtonTRDCSetPostCheck(SNES snes,PetscErrorCode (* func)(SNES,Vec,Vec,Vec,PetscBool *,PetscBool *,void *),PetscCtx ctx)182*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtx ctx)
183d71ae5a4SJacob Faibussowitsch {
18441ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
18541ba4c6cSHeeho Park 
18641ba4c6cSHeeho Park   PetscFunctionBegin;
18741ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
18841ba4c6cSHeeho Park   if (func) tr->postcheck = func;
18941ba4c6cSHeeho Park   if (ctx) tr->postcheckctx = ctx;
1903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19141ba4c6cSHeeho Park }
19241ba4c6cSHeeho Park 
19341ba4c6cSHeeho Park /*@C
194420bcc1bSBarry Smith   SNESNewtonTRDCGetPostCheck - Gets the post-check function optionally set with `SNESNewtonTRDCSetPostCheck()`
19541ba4c6cSHeeho Park 
19620f4b53cSBarry Smith   Not Collective
19741ba4c6cSHeeho Park 
19841ba4c6cSHeeho Park   Input Parameter:
19941ba4c6cSHeeho Park . snes - the nonlinear solver context
20041ba4c6cSHeeho Park 
20141ba4c6cSHeeho Park   Output Parameters:
20220f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()`
20320f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
20441ba4c6cSHeeho Park 
20541ba4c6cSHeeho Park   Level: intermediate
20641ba4c6cSHeeho Park 
207420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`
20841ba4c6cSHeeho Park @*/
SNESNewtonTRDCGetPostCheck(SNES snes,PetscErrorCode (** func)(SNES,Vec,Vec,Vec,PetscBool *,PetscBool *,void *),PetscCtxRt ctx)209*2a8381b2SBarry Smith PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtxRt ctx)
210d71ae5a4SJacob Faibussowitsch {
21141ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
21241ba4c6cSHeeho Park 
21341ba4c6cSHeeho Park   PetscFunctionBegin;
21441ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
21541ba4c6cSHeeho Park   if (func) *func = tr->postcheck;
216*2a8381b2SBarry Smith   if (ctx) *(void **)ctx = tr->postcheckctx;
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21841ba4c6cSHeeho Park }
21941ba4c6cSHeeho Park 
220ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-internal-linkage
22141ba4c6cSHeeho Park /*@C
222f6dfbefdSBarry Smith    SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC`
22341ba4c6cSHeeho Park 
224c3339decSBarry Smith    Logically Collective
22541ba4c6cSHeeho Park 
22641ba4c6cSHeeho Park    Input Parameters:
22741ba4c6cSHeeho Park +  snes - the solver
22841ba4c6cSHeeho Park .  X - The last solution
22941ba4c6cSHeeho Park -  Y - The step direction
23041ba4c6cSHeeho Park 
2312fe279fdSBarry Smith    Output Parameter:
23220f4b53cSBarry Smith .  changed_Y - Indicator that the step direction `Y` has been changed.
23341ba4c6cSHeeho Park 
23441ba4c6cSHeeho Park    Level: developer
23541ba4c6cSHeeho Park 
236420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()`
23741ba4c6cSHeeho Park @*/
SNESNewtonTRDCPreCheck(SNES snes,Vec X,Vec Y,PetscBool * changed_Y)238d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y)
239d71ae5a4SJacob Faibussowitsch {
24041ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
24141ba4c6cSHeeho Park 
24241ba4c6cSHeeho Park   PetscFunctionBegin;
24341ba4c6cSHeeho Park   *changed_Y = PETSC_FALSE;
24441ba4c6cSHeeho Park   if (tr->precheck) {
2459566063dSJacob Faibussowitsch     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
24641ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_Y, 4);
24741ba4c6cSHeeho Park   }
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24941ba4c6cSHeeho Park }
25041ba4c6cSHeeho Park 
251ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-internal-linkage
25241ba4c6cSHeeho Park /*@C
253f6dfbefdSBarry Smith    SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step
25441ba4c6cSHeeho Park 
255c3339decSBarry Smith    Logically Collective
25641ba4c6cSHeeho Park 
25741ba4c6cSHeeho Park    Input Parameters:
258f1a722f8SMatthew G. Knepley +  snes - the solver
259f1a722f8SMatthew G. Knepley .  X - The last solution
26041ba4c6cSHeeho Park .  Y - The full step direction
26141ba4c6cSHeeho Park -  W - The updated solution, W = X - Y
26241ba4c6cSHeeho Park 
26341ba4c6cSHeeho Park    Output Parameters:
26441ba4c6cSHeeho Park +  changed_Y - indicator if step has been changed
26520f4b53cSBarry Smith -  changed_W - Indicator if the new candidate solution `W` has been changed.
26641ba4c6cSHeeho Park 
2672fe279fdSBarry Smith    Level: developer
2682fe279fdSBarry Smith 
269f6dfbefdSBarry Smith    Note:
27020f4b53cSBarry Smith      If `Y` is changed then `W` is recomputed as `X` - `Y`
27141ba4c6cSHeeho Park 
272420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck()
27341ba4c6cSHeeho Park @*/
SNESNewtonTRDCPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool * changed_Y,PetscBool * changed_W)274d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
275d71ae5a4SJacob Faibussowitsch {
27641ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
27741ba4c6cSHeeho Park 
27841ba4c6cSHeeho Park   PetscFunctionBegin;
27941ba4c6cSHeeho Park   *changed_Y = PETSC_FALSE;
28041ba4c6cSHeeho Park   *changed_W = PETSC_FALSE;
28141ba4c6cSHeeho Park   if (tr->postcheck) {
2829566063dSJacob Faibussowitsch     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
28341ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_Y, 5);
28441ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_W, 6);
28541ba4c6cSHeeho Park   }
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28741ba4c6cSHeeho Park }
28841ba4c6cSHeeho Park 
28941ba4c6cSHeeho Park /*
29041ba4c6cSHeeho Park    SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
29141ba4c6cSHeeho Park    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
29241ba4c6cSHeeho Park    nonlinear equations
29341ba4c6cSHeeho Park 
29441ba4c6cSHeeho Park */
SNESSolve_NEWTONTRDC(SNES snes)295d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes)
296d71ae5a4SJacob Faibussowitsch {
29741ba4c6cSHeeho Park   SNES_NEWTONTRDC            *neP = (SNES_NEWTONTRDC *)snes->data;
29841ba4c6cSHeeho Park   Vec                         X, F, Y, G, W, GradF, YNtmp;
29941ba4c6cSHeeho Park   Vec                         YCtmp;
30041ba4c6cSHeeho Park   Mat                         jac;
30141ba4c6cSHeeho Park   PetscInt                    maxits, i, j, lits, inner_count, bs;
30241ba4c6cSHeeho Park   PetscReal                   rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */
30341ba4c6cSHeeho Park   PetscReal                   inorms[99];                                                         /* need to make it dynamic eventually, fixed max block size of 99 for now */
30441ba4c6cSHeeho Park   PetscReal                   deltaM, ynnorm, f0, mp, gTy, g, yTHy;                               /* rho calculation */
30541ba4c6cSHeeho Park   PetscReal                   auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg;       /* Cauchy Point */
30641ba4c6cSHeeho Park   KSP                         ksp;
30741ba4c6cSHeeho Park   SNESConvergedReason         reason   = SNES_CONVERGED_ITERATING;
30841ba4c6cSHeeho Park   PetscBool                   breakout = PETSC_FALSE;
30941ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx;
3104d4d2bdcSBarry Smith   KSPConvergenceTestFn       *convtest;
3114d4d2bdcSBarry Smith   PetscCtxDestroyFn          *convdestroy;
31241ba4c6cSHeeho Park   void                       *convctx;
31341ba4c6cSHeeho Park 
31441ba4c6cSHeeho Park   PetscFunctionBegin;
31541ba4c6cSHeeho Park   maxits = snes->max_its;  /* maximum number of iterations */
31641ba4c6cSHeeho Park   X      = snes->vec_sol;  /* solution vector */
31741ba4c6cSHeeho Park   F      = snes->vec_func; /* residual vector */
31841ba4c6cSHeeho Park   Y      = snes->work[0];  /* update vector */
31941ba4c6cSHeeho Park   G      = snes->work[1];  /* updated residual */
32041ba4c6cSHeeho Park   W      = snes->work[2];  /* temporary vector */
32141ba4c6cSHeeho Park   GradF  = snes->work[3];  /* grad f = J^T F */
32241ba4c6cSHeeho Park   YNtmp  = snes->work[4];  /* Newton solution */
32341ba4c6cSHeeho Park   YCtmp  = snes->work[5];  /* Cauchy solution */
32441ba4c6cSHeeho Park 
3250b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
32641ba4c6cSHeeho Park 
3279566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(YNtmp, &bs));
32841ba4c6cSHeeho Park 
3299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
33041ba4c6cSHeeho Park   snes->iter = 0;
3319566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
33241ba4c6cSHeeho Park 
33341ba4c6cSHeeho Park   /* Set the linear stopping criteria to use the More' trick. From tr.c */
3349566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
3359566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
33641ba4c6cSHeeho Park   if (convtest != SNESTRDC_KSPConverged_Private) {
3379566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
33841ba4c6cSHeeho Park     ctx->snes = snes;
3399566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
3409566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy));
3419566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n"));
34241ba4c6cSHeeho Park   }
34341ba4c6cSHeeho Park 
34441ba4c6cSHeeho Park   if (!snes->vec_func_init_set) {
3459566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
34641ba4c6cSHeeho Park   } else snes->vec_func_init_set = PETSC_FALSE;
34741ba4c6cSHeeho Park 
3489566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
34976c63389SBarry Smith   SNESCheckFunctionDomainError(snes, fnorm);
3509566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
35141ba4c6cSHeeho Park 
3529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
35341ba4c6cSHeeho Park   snes->norm = fnorm;
3549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
35541ba4c6cSHeeho Park   delta      = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */
35641ba4c6cSHeeho Park   deltaM     = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */
35741ba4c6cSHeeho Park   neP->delta = delta;
3589566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3599566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
36041ba4c6cSHeeho Park 
36141ba4c6cSHeeho Park   neP->rho_satisfied = PETSC_FALSE;
36241ba4c6cSHeeho Park 
36341ba4c6cSHeeho Park   /* test convergence */
364dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
3653ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
36641ba4c6cSHeeho Park 
36741ba4c6cSHeeho Park   for (i = 0; i < maxits; i++) {
36841ba4c6cSHeeho Park     PetscBool changed_y;
36941ba4c6cSHeeho Park     PetscBool changed_w;
37041ba4c6cSHeeho Park 
37141ba4c6cSHeeho Park     /* dogleg method */
3729566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
37376c63389SBarry Smith     SNESCheckJacobianDomainError(snes);
3749566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian));
3759566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */
37641ba4c6cSHeeho Park     SNESCheckKSPSolve(snes);                  /* this is necessary but old tr.c did not have it*/
3779566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
3789566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
37941ba4c6cSHeeho Park 
38041ba4c6cSHeeho Park     /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable)
38141ba4c6cSHeeho Park        for inner iteration and Cauchy direction calculation
38241ba4c6cSHeeho Park     */
38341ba4c6cSHeeho Park     if (bs > 1 && neP->auto_scale_multiphase) {
3849566063dSJacob Faibussowitsch       PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms));
38541ba4c6cSHeeho Park       for (j = 0; j < bs; j++) {
38641ba4c6cSHeeho Park         if (neP->auto_scale_max > 1.0) {
387ad540459SPierre Jolivet           if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max;
38841ba4c6cSHeeho Park         }
3899566063dSJacob Faibussowitsch         PetscCall(VecStrideSet(W, j, inorms[j]));
3909566063dSJacob Faibussowitsch         PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j]));
3919566063dSJacob Faibussowitsch         PetscCall(VecStrideScale(X, j, 1.0 / inorms[j]));
39241ba4c6cSHeeho Park       }
3939566063dSJacob Faibussowitsch       PetscCall(VecNorm(X, NORM_2, &xnorm));
39441ba4c6cSHeeho Park       if (i == 0) {
39541ba4c6cSHeeho Park         delta = neP->delta0 * xnorm;
39641ba4c6cSHeeho Park       } else {
39741ba4c6cSHeeho Park         delta = neP->delta * xnorm;
39841ba4c6cSHeeho Park       }
39941ba4c6cSHeeho Park       deltaM = neP->deltaM * xnorm;
400f3fa974cSJacob Faibussowitsch       PetscCall(MatDiagonalScale(jac, NULL, W));
40141ba4c6cSHeeho Park     }
40241ba4c6cSHeeho Park 
40341ba4c6cSHeeho Park     /* calculating GradF of minimization function */
4049566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */
4059566063dSJacob Faibussowitsch     PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */
40641ba4c6cSHeeho Park 
40741ba4c6cSHeeho Park     inner_count        = 0;
40841ba4c6cSHeeho Park     neP->rho_satisfied = PETSC_FALSE;
40941ba4c6cSHeeho Park     while (1) {
41041ba4c6cSHeeho Park       if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */
4119566063dSJacob Faibussowitsch         PetscCall(VecCopy(YNtmp, Y));
41241ba4c6cSHeeho Park       } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */
4139566063dSJacob Faibussowitsch         PetscCall(MatMult(jac, GradF, W));
4149566063dSJacob Faibussowitsch         PetscCall(VecDotRealPart(W, W, &gTBg));     /* completes GradF^T J^T J GradF */
4159566063dSJacob Faibussowitsch         PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */
41641ba4c6cSHeeho Park         if (gTBg <= 0.0) {
41741ba4c6cSHeeho Park           auk = PETSC_MAX_REAL;
41841ba4c6cSHeeho Park         } else {
41941ba4c6cSHeeho Park           auk = PetscSqr(gfnorm) / gTBg;
42041ba4c6cSHeeho Park         }
42141ba4c6cSHeeho Park         auk = PetscMin(delta / gfnorm, auk);
4229566063dSJacob Faibussowitsch         PetscCall(VecCopy(GradF, YCtmp));           /* this could be improved */
4239566063dSJacob Faibussowitsch         PetscCall(VecScale(YCtmp, auk));            /* YCtmp, Cauchy solution*/
4249566063dSJacob Faibussowitsch         PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */
42541ba4c6cSHeeho Park         if (ycnorm >= delta) {                      /* see if the Cauchy solution meets the criteria */
4269566063dSJacob Faibussowitsch           PetscCall(VecCopy(YCtmp, Y));
4279566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm));
42841ba4c6cSHeeho Park         } else {                                  /* take ratio, tau, of Cauchy and Newton direction and step */
4299566063dSJacob Faibussowitsch           PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */
4309566063dSJacob Faibussowitsch           PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */
43141ba4c6cSHeeho Park           c0 = PetscSqr(c0);
4329566063dSJacob Faibussowitsch           PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1));
43341ba4c6cSHeeho Park           c1 = 2.0 * c1;
4349566063dSJacob Faibussowitsch           PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */
43541ba4c6cSHeeho Park           c2      = PetscSqr(c2) - PetscSqr(delta);
43641ba4c6cSHeeho Park           tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */
43741ba4c6cSHeeho Park           tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0);
43841ba4c6cSHeeho Park           tau     = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */
4399566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm));
4409566063dSJacob Faibussowitsch           PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp));
4419566063dSJacob Faibussowitsch           PetscCall(VecAXPY(W, -tau, YCtmp));
4429566063dSJacob Faibussowitsch           PetscCall(VecCopy(W, Y)); /* this could be improved */
44341ba4c6cSHeeho Park         }
44441ba4c6cSHeeho Park       } else {
44541ba4c6cSHeeho Park         /* if Cauchy is disabled, only use Newton direction */
44641ba4c6cSHeeho Park         auk = delta / ynnorm;
4479566063dSJacob Faibussowitsch         PetscCall(VecScale(YNtmp, auk));
4489566063dSJacob Faibussowitsch         PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/
44941ba4c6cSHeeho Park       }
45041ba4c6cSHeeho Park 
4519566063dSJacob Faibussowitsch       PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm  */
45241ba4c6cSHeeho Park       f0 = 0.5 * PetscSqr(fnorm);            /* minimizing function f(X) */
4539566063dSJacob Faibussowitsch       PetscCall(MatMult(jac, Y, W));
4549566063dSJacob Faibussowitsch       PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */
4559566063dSJacob Faibussowitsch       PetscCall(VecDotRealPart(GradF, Y, &gTy));
45641ba4c6cSHeeho Park       mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/
45741ba4c6cSHeeho Park 
45841ba4c6cSHeeho Park       /* scale back solution update */
45941ba4c6cSHeeho Park       if (bs > 1 && neP->auto_scale_multiphase) {
46041ba4c6cSHeeho Park         for (j = 0; j < bs; j++) {
4619566063dSJacob Faibussowitsch           PetscCall(VecStrideScale(Y, j, inorms[j]));
46241ba4c6cSHeeho Park           if (inner_count == 0) {
46341ba4c6cSHeeho Park             /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */
46441ba4c6cSHeeho Park             /* need to scale back X to match Y and provide proper update to the external code */
4659566063dSJacob Faibussowitsch             PetscCall(VecStrideScale(X, j, inorms[j]));
46641ba4c6cSHeeho Park           }
46741ba4c6cSHeeho Park         }
4689566063dSJacob Faibussowitsch         if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */
4699566063dSJacob Faibussowitsch         PetscCall(VecNorm(Y, NORM_2, &temp_ynorm));
47041ba4c6cSHeeho Park       } else {
47141ba4c6cSHeeho Park         temp_xnorm = xnorm;
47241ba4c6cSHeeho Park         temp_ynorm = ynorm;
47341ba4c6cSHeeho Park       }
47441ba4c6cSHeeho Park       inner_count++;
47541ba4c6cSHeeho Park 
47641ba4c6cSHeeho Park       /* Evaluate the solution to meet the improvement ratio criteria */
4779566063dSJacob Faibussowitsch       PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y));
4789566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(W, -1.0, Y, X));
4799566063dSJacob Faibussowitsch       PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w));
4809566063dSJacob Faibussowitsch       if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X));
4819566063dSJacob Faibussowitsch       PetscCall(VecCopy(Y, snes->vec_sol_update));
4829566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, W, G)); /*  F(X-Y) = G */
4839566063dSJacob Faibussowitsch       PetscCall(VecNorm(G, NORM_2, &gnorm));      /* gnorm <- || g || */
48476c63389SBarry Smith       SNESCheckFunctionDomainError(snes, gnorm);
48541ba4c6cSHeeho Park       g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */
48641ba4c6cSHeeho Park       if (f0 == mp) rho = 0.0;
48741ba4c6cSHeeho Park       else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */
48841ba4c6cSHeeho Park 
48941ba4c6cSHeeho Park       if (rho < neP->eta2) {
49041ba4c6cSHeeho Park         delta *= neP->t1; /* shrink the region */
49141ba4c6cSHeeho Park       } else if (rho > neP->eta3) {
49241ba4c6cSHeeho Park         delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */
49341ba4c6cSHeeho Park       }
49441ba4c6cSHeeho Park 
49541ba4c6cSHeeho Park       neP->delta = delta;
49641ba4c6cSHeeho Park       if (rho >= neP->eta1) {
49741ba4c6cSHeeho Park         /* unscale delta and xnorm before going to the next outer iteration */
49841ba4c6cSHeeho Park         if (bs > 1 && neP->auto_scale_multiphase) {
49941ba4c6cSHeeho Park           neP->delta = delta / xnorm;
50041ba4c6cSHeeho Park           xnorm      = temp_xnorm;
50141ba4c6cSHeeho Park           ynorm      = temp_ynorm;
50241ba4c6cSHeeho Park         }
50341ba4c6cSHeeho Park         neP->rho_satisfied = PETSC_TRUE;
50441ba4c6cSHeeho Park         break; /* the improvement ratio is satisfactory */
50541ba4c6cSHeeho Park       }
5069566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
50741ba4c6cSHeeho Park 
50841ba4c6cSHeeho Park       /* check to see if progress is hopeless */
50941ba4c6cSHeeho Park       neP->itflag = PETSC_FALSE;
51041ba4c6cSHeeho Park       /* both delta, ynorm, and xnorm are either scaled or unscaled */
5119566063dSJacob Faibussowitsch       PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP));
51241ba4c6cSHeeho Park       /* if multiphase state changes, break out inner iteration */
51341ba4c6cSHeeho Park       if (reason == SNES_BREAKOUT_INNER_ITER) {
51441ba4c6cSHeeho Park         if (bs > 1 && neP->auto_scale_multiphase) {
51541ba4c6cSHeeho Park           /* unscale delta and xnorm before going to the next outer iteration */
51641ba4c6cSHeeho Park           neP->delta = delta / xnorm;
51741ba4c6cSHeeho Park           xnorm      = temp_xnorm;
51841ba4c6cSHeeho Park           ynorm      = temp_ynorm;
51941ba4c6cSHeeho Park         }
52041ba4c6cSHeeho Park         reason = SNES_CONVERGED_ITERATING;
52141ba4c6cSHeeho Park         break;
52241ba4c6cSHeeho Park       }
52341ba4c6cSHeeho Park       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
52441ba4c6cSHeeho Park       if (reason) {
52541ba4c6cSHeeho Park         if (reason < 0) {
52641ba4c6cSHeeho Park           /* We're not progressing, so return with the current iterate */
5279566063dSJacob Faibussowitsch           PetscCall(SNESMonitor(snes, i + 1, fnorm));
52841ba4c6cSHeeho Park           breakout = PETSC_TRUE;
52941ba4c6cSHeeho Park           break;
53041ba4c6cSHeeho Park         } else if (reason > 0) {
53141ba4c6cSHeeho Park           /* We're converged, so return with the current iterate and update solution */
5329566063dSJacob Faibussowitsch           PetscCall(SNESMonitor(snes, i + 1, fnorm));
53341ba4c6cSHeeho Park           breakout = PETSC_FALSE;
53441ba4c6cSHeeho Park           break;
53541ba4c6cSHeeho Park         }
53641ba4c6cSHeeho Park       }
53741ba4c6cSHeeho Park       snes->numFailures++;
53841ba4c6cSHeeho Park     }
53941ba4c6cSHeeho Park     if (!breakout) {
54041ba4c6cSHeeho Park       /* Update function and solution vectors */
54141ba4c6cSHeeho Park       fnorm = gnorm;
5429566063dSJacob Faibussowitsch       PetscCall(VecCopy(G, F));
5439566063dSJacob Faibussowitsch       PetscCall(VecCopy(W, X));
54441ba4c6cSHeeho Park       /* Monitor convergence */
5459566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
54641ba4c6cSHeeho Park       snes->iter  = i + 1;
54741ba4c6cSHeeho Park       snes->norm  = fnorm;
54841ba4c6cSHeeho Park       snes->xnorm = xnorm;
54941ba4c6cSHeeho Park       snes->ynorm = ynorm;
5509566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5519566063dSJacob Faibussowitsch       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
5529566063dSJacob Faibussowitsch       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
55341ba4c6cSHeeho Park       /* Test for convergence, xnorm = || X || */
55441ba4c6cSHeeho Park       neP->itflag = PETSC_TRUE;
5559566063dSJacob Faibussowitsch       if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
556dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
55741ba4c6cSHeeho Park       if (reason) break;
55841ba4c6cSHeeho Park     } else break;
55941ba4c6cSHeeho Park   }
56041ba4c6cSHeeho Park 
5619566063dSJacob Faibussowitsch   /* PetscCall(PetscFree(inorms)); */
56241ba4c6cSHeeho Park   if (i == maxits) {
5639566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
56441ba4c6cSHeeho Park     if (!reason) reason = SNES_DIVERGED_MAX_IT;
56541ba4c6cSHeeho Park   }
5669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
56741ba4c6cSHeeho Park   snes->reason = reason;
5689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
56941ba4c6cSHeeho Park   if (convtest != SNESTRDC_KSPConverged_Private) {
5709566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
5719566063dSJacob Faibussowitsch     PetscCall(PetscFree(ctx));
5729566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
57341ba4c6cSHeeho Park   }
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57541ba4c6cSHeeho Park }
57641ba4c6cSHeeho Park 
SNESSetUp_NEWTONTRDC(SNES snes)577d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes)
578d71ae5a4SJacob Faibussowitsch {
57941ba4c6cSHeeho Park   PetscFunctionBegin;
5809566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 6));
5819566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
5823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58341ba4c6cSHeeho Park }
58441ba4c6cSHeeho Park 
SNESDestroy_NEWTONTRDC(SNES snes)585d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes)
586d71ae5a4SJacob Faibussowitsch {
58741ba4c6cSHeeho Park   PetscFunctionBegin;
5883201ab8dSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL));
5899566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59141ba4c6cSHeeho Park }
59241ba4c6cSHeeho Park 
SNESSetFromOptions_NEWTONTRDC(SNES snes,PetscOptionItems PetscOptionsObject)593ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems PetscOptionsObject)
594d71ae5a4SJacob Faibussowitsch {
59541ba4c6cSHeeho Park   SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data;
59641ba4c6cSHeeho Park 
59741ba4c6cSHeeho Park   PetscFunctionBegin;
598d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
5993201ab8dSStefano Zampini   PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESNewtonTRSetTolerances", ctx->deltatol, &ctx->deltatol, NULL));
6009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
6019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
6029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
6039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
6049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
6059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
6069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
6079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL));
6089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL));
6099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_trdc_auto_scale_multiphase", "auto_scale_multiphase", "Auto scaling for proper cauchy direction", ctx->auto_scale_multiphase, &ctx->auto_scale_multiphase, NULL));
610d0609cedSBarry Smith   PetscOptionsHeadEnd();
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61241ba4c6cSHeeho Park }
61341ba4c6cSHeeho Park 
SNESView_NEWTONTRDC(SNES snes,PetscViewer viewer)614d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer)
615d71ae5a4SJacob Faibussowitsch {
61641ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
6179f196a02SMartin Diehl   PetscBool        isascii;
61841ba4c6cSHeeho Park 
61941ba4c6cSHeeho Park   PetscFunctionBegin;
6209f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
6219f196a02SMartin Diehl   if (isascii) {
6223201ab8dSStefano Zampini     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region tolerance %g\n", (double)tr->deltatol));
6239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
6249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM));
62541ba4c6cSHeeho Park   }
6263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62741ba4c6cSHeeho Park }
628f6dfbefdSBarry Smith 
62941ba4c6cSHeeho Park /*MC
63041ba4c6cSHeeho Park       SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
63141ba4c6cSHeeho Park 
632f6dfbefdSBarry Smith    Options Database Keys:
63341ba4c6cSHeeho Park +   -snes_trdc_tol <tol>                                     - trust region tolerance
63441ba4c6cSHeeho Park .   -snes_trdc_eta1 <eta1>                                   - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
63541ba4c6cSHeeho Park .   -snes_trdc_eta2 <eta2>                                   - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
63641ba4c6cSHeeho Park .   -snes_trdc_eta3 <eta3>                                   - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
63741ba4c6cSHeeho Park .   -snes_trdc_t1 <t1>                                       - trust region parameter, shrinking factor of trust region (default: 0.25)
63841ba4c6cSHeeho Park .   -snes_trdc_t2 <t2>                                       - trust region parameter, expanding factor of trust region (default: 2.0)
6391d27aa22SBarry Smith .   -snes_trdc_deltaM <deltaM>                               - trust region parameter, max size of trust region, $deltaM*norm2(x)$ (default: 0.5)
6401d27aa22SBarry Smith .   -snes_trdc_delta0 <delta0>                               - trust region parameter, initial size of trust region, $delta0*norm2(x)$ (default: 0.1)
64141ba4c6cSHeeho Park .   -snes_trdc_auto_scale_max <auto_scale_max>               - used with auto_scale_multiphase, caps the maximum auto-scaling factor
64241ba4c6cSHeeho Park .   -snes_trdc_use_cauchy <use_cauchy>                       - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm
64341ba4c6cSHeeho Park -   -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region
64441ba4c6cSHeeho Park 
6453201ab8dSStefano Zampini    Level: advanced
64620f4b53cSBarry Smith 
6473201ab8dSStefano Zampini    Notes:
6483201ab8dSStefano Zampini    `SNESNEWTONTRDC` only works for root-finding problems and does not support objective functions.
6493201ab8dSStefano Zampini    The main difference with respect to `SNESNEWTONTR` is that `SNESNEWTONTRDC` scales the trust region by the norm of the current linearization point.
6503201ab8dSStefano Zampini    Future version may extend the `SNESNEWTONTR` code and deprecate `SNESNEWTONTRDC`.
65141ba4c6cSHeeho Park 
6523201ab8dSStefano Zampini    For details, see {cite}`park2021linear`
6533201ab8dSStefano Zampini 
6543201ab8dSStefano Zampini .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNewtonTRSetTolerances()`,
655f6dfbefdSBarry Smith           `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`,
656f6dfbefdSBarry Smith           `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()`
65741ba4c6cSHeeho Park M*/
SNESCreate_NEWTONTRDC(SNES snes)658d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes)
659d71ae5a4SJacob Faibussowitsch {
66041ba4c6cSHeeho Park   SNES_NEWTONTRDC *neP;
66141ba4c6cSHeeho Park 
66241ba4c6cSHeeho Park   PetscFunctionBegin;
66341ba4c6cSHeeho Park   snes->ops->setup          = SNESSetUp_NEWTONTRDC;
66441ba4c6cSHeeho Park   snes->ops->solve          = SNESSolve_NEWTONTRDC;
66541ba4c6cSHeeho Park   snes->ops->destroy        = SNESDestroy_NEWTONTRDC;
66641ba4c6cSHeeho Park   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC;
66741ba4c6cSHeeho Park   snes->ops->view           = SNESView_NEWTONTRDC;
66841ba4c6cSHeeho Park 
66941ba4c6cSHeeho Park   snes->usesksp = PETSC_TRUE;
67041ba4c6cSHeeho Park   snes->usesnpc = PETSC_FALSE;
67141ba4c6cSHeeho Park 
67241ba4c6cSHeeho Park   snes->alwayscomputesfinalresidual = PETSC_TRUE;
67341ba4c6cSHeeho Park 
67477e5a1f9SBarry Smith   PetscCall(SNESParametersInitialize(snes));
67577e5a1f9SBarry Smith 
6764dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
67741ba4c6cSHeeho Park   snes->data                 = (void *)neP;
67841ba4c6cSHeeho Park   neP->eta1                  = 0.001;
67941ba4c6cSHeeho Park   neP->eta2                  = 0.25;
68041ba4c6cSHeeho Park   neP->eta3                  = 0.75;
68141ba4c6cSHeeho Park   neP->t1                    = 0.25;
68241ba4c6cSHeeho Park   neP->t2                    = 2.0;
68341ba4c6cSHeeho Park   neP->sigma                 = 0.0001;
68441ba4c6cSHeeho Park   neP->itflag                = PETSC_FALSE;
68541ba4c6cSHeeho Park   neP->rnorm0                = 0.0;
68641ba4c6cSHeeho Park   neP->ttol                  = 0.0;
68741ba4c6cSHeeho Park   neP->use_cauchy            = PETSC_TRUE;
68841ba4c6cSHeeho Park   neP->auto_scale_multiphase = PETSC_FALSE;
68941ba4c6cSHeeho Park   neP->auto_scale_max        = -1.0;
69041ba4c6cSHeeho Park   neP->rho_satisfied         = PETSC_FALSE;
6913201ab8dSStefano Zampini   neP->delta                 = 0.0;
6923201ab8dSStefano Zampini   neP->deltaM                = 0.5;
6933201ab8dSStefano Zampini   neP->delta0                = 0.1;
6943201ab8dSStefano Zampini   neP->deltatol              = 1.e-12;
69541ba4c6cSHeeho Park 
69641ba4c6cSHeeho Park   /* for multiphase (multivariable) scaling */
69741ba4c6cSHeeho Park   /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13
69841ba4c6cSHeeho Park      on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now.
6999566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(snes->work[0],&neP->bs));
7009566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(neP->bs,&neP->inorms));
70141ba4c6cSHeeho Park   */
7023201ab8dSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TRDC));
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70441ba4c6cSHeeho Park }
705