xref: /petsc/src/snes/impls/ntrdc/ntrdc.c (revision f3fa974cd8ac3307c23c9d7e81b13b5f402f9e6e)
141ba4c6cSHeeho Park 
241ba4c6cSHeeho Park #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I   "petscsnes.h"   I*/
341ba4c6cSHeeho Park 
441ba4c6cSHeeho Park typedef struct {
541ba4c6cSHeeho Park   SNES snes;
641ba4c6cSHeeho Park   /*  Information on the regular SNES convergence test; which may have been user provided
741ba4c6cSHeeho Park       Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho
841ba4c6cSHeeho Park       Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private
941ba4c6cSHeeho Park  */
1041ba4c6cSHeeho Park 
1141ba4c6cSHeeho Park   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
1241ba4c6cSHeeho Park   PetscErrorCode (*convdestroy)(void *);
1341ba4c6cSHeeho Park   void *convctx;
1441ba4c6cSHeeho Park } SNES_TRDC_KSPConverged_Ctx;
1541ba4c6cSHeeho Park 
16d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx)
17d71ae5a4SJacob Faibussowitsch {
1841ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx  = (SNES_TRDC_KSPConverged_Ctx *)cctx;
1941ba4c6cSHeeho Park   SNES                        snes = ctx->snes;
2041ba4c6cSHeeho Park   SNES_NEWTONTRDC            *neP  = (SNES_NEWTONTRDC *)snes->data;
2141ba4c6cSHeeho Park   Vec                         x;
2241ba4c6cSHeeho Park   PetscReal                   nrm;
2341ba4c6cSHeeho Park 
2441ba4c6cSHeeho Park   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
2648a46eb9SPierre Jolivet   if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm));
2741ba4c6cSHeeho Park   /* Determine norm of solution */
289566063dSJacob Faibussowitsch   PetscCall(KSPBuildSolution(ksp, NULL, &x));
299566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &nrm));
3041ba4c6cSHeeho Park   if (nrm >= neP->delta) {
319566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
3241ba4c6cSHeeho Park     *reason = KSP_CONVERGED_STEP_LENGTH;
3341ba4c6cSHeeho Park   }
3441ba4c6cSHeeho Park   PetscFunctionReturn(0);
3541ba4c6cSHeeho Park }
3641ba4c6cSHeeho Park 
37d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx)
38d71ae5a4SJacob Faibussowitsch {
3941ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx;
4041ba4c6cSHeeho Park 
4141ba4c6cSHeeho Park   PetscFunctionBegin;
429566063dSJacob Faibussowitsch   PetscCall((*ctx->convdestroy)(ctx->convctx));
439566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
4441ba4c6cSHeeho Park 
4541ba4c6cSHeeho Park   PetscFunctionReturn(0);
4641ba4c6cSHeeho Park }
4741ba4c6cSHeeho Park 
4841ba4c6cSHeeho Park /*
4941ba4c6cSHeeho Park    SNESTRDC_Converged_Private -test convergence JUST for
5041ba4c6cSHeeho Park    the trust region tolerance.
5141ba4c6cSHeeho Park 
5241ba4c6cSHeeho Park */
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
54d71ae5a4SJacob Faibussowitsch {
5541ba4c6cSHeeho Park   SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
5641ba4c6cSHeeho Park 
5741ba4c6cSHeeho Park   PetscFunctionBegin;
5841ba4c6cSHeeho Park   *reason = SNES_CONVERGED_ITERATING;
5941ba4c6cSHeeho Park   if (neP->delta < xnorm * snes->deltatol) {
609566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol));
6141ba4c6cSHeeho Park     *reason = SNES_DIVERGED_TR_DELTA;
6241ba4c6cSHeeho Park   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
639566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
6441ba4c6cSHeeho Park     *reason = SNES_DIVERGED_FUNCTION_COUNT;
6541ba4c6cSHeeho Park   }
6641ba4c6cSHeeho Park   PetscFunctionReturn(0);
6741ba4c6cSHeeho Park }
6841ba4c6cSHeeho Park 
6941ba4c6cSHeeho Park /*@
70f6dfbefdSBarry Smith   SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region.
7141ba4c6cSHeeho Park 
72f6dfbefdSBarry Smith   Input Parameter:
7341ba4c6cSHeeho Park . snes - the nonlinear solver object
7441ba4c6cSHeeho Park 
75f6dfbefdSBarry Smith   Output Parameter:
76f6dfbefdSBarry Smith . rho_flag: `PETSC_TRUE` if the solution update is in the trust-region; otherwise, `PETSC_FALSE`
7741ba4c6cSHeeho Park 
7841ba4c6cSHeeho Park   Level: developer
7941ba4c6cSHeeho Park 
80f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, , `SNESNewtonTRDCSetPreCheck()`,
81f6dfbefdSBarry Smith           `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
8241ba4c6cSHeeho Park @*/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag)
84d71ae5a4SJacob Faibussowitsch {
8541ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
8641ba4c6cSHeeho Park 
8741ba4c6cSHeeho Park   PetscFunctionBegin;
8841ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
8941ba4c6cSHeeho Park   PetscValidBoolPointer(rho_flag, 2);
9041ba4c6cSHeeho Park   *rho_flag = tr->rho_satisfied;
9141ba4c6cSHeeho Park   PetscFunctionReturn(0);
9241ba4c6cSHeeho Park }
9341ba4c6cSHeeho Park 
9441ba4c6cSHeeho Park /*@C
9541ba4c6cSHeeho Park    SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined.
9641ba4c6cSHeeho Park        Allows the user a chance to change or override the trust region decision.
9741ba4c6cSHeeho Park 
98c3339decSBarry Smith    Logically Collective
9941ba4c6cSHeeho Park 
10041ba4c6cSHeeho Park    Input Parameters:
10141ba4c6cSHeeho Park +  snes - the nonlinear solver object
102f6dfbefdSBarry Smith .  func - [optional] function evaluation routine, see `SNESNewtonTRDCPreCheck()`  for the calling sequence
10341ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
10441ba4c6cSHeeho Park 
10541ba4c6cSHeeho Park    Level: intermediate
10641ba4c6cSHeeho Park 
107f6dfbefdSBarry Smith    Note:
108f6dfbefdSBarry Smith    This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver.
10941ba4c6cSHeeho Park 
110f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`,
111f6dfbefdSBarry Smith           `SNESNewtonTRDCGetRhoFlag()`
11241ba4c6cSHeeho Park @*/
113d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx)
114d71ae5a4SJacob Faibussowitsch {
11541ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
11641ba4c6cSHeeho Park 
11741ba4c6cSHeeho Park   PetscFunctionBegin;
11841ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
11941ba4c6cSHeeho Park   if (func) tr->precheck = func;
12041ba4c6cSHeeho Park   if (ctx) tr->precheckctx = ctx;
12141ba4c6cSHeeho Park   PetscFunctionReturn(0);
12241ba4c6cSHeeho Park }
12341ba4c6cSHeeho Park 
12441ba4c6cSHeeho Park /*@C
12541ba4c6cSHeeho Park    SNESNewtonTRDCGetPreCheck - Gets the pre-check function
12641ba4c6cSHeeho Park 
12741ba4c6cSHeeho Park    Not collective
12841ba4c6cSHeeho Park 
12941ba4c6cSHeeho Park    Input Parameter:
13041ba4c6cSHeeho Park .  snes - the nonlinear solver context
13141ba4c6cSHeeho Park 
13241ba4c6cSHeeho Park    Output Parameters:
133f6dfbefdSBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRDCPreCheck()`
13441ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
13541ba4c6cSHeeho Park 
13641ba4c6cSHeeho Park    Level: intermediate
13741ba4c6cSHeeho Park 
138f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()`
13941ba4c6cSHeeho Park @*/
140d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx)
141d71ae5a4SJacob Faibussowitsch {
14241ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
14341ba4c6cSHeeho Park 
14441ba4c6cSHeeho Park   PetscFunctionBegin;
14541ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
14641ba4c6cSHeeho Park   if (func) *func = tr->precheck;
14741ba4c6cSHeeho Park   if (ctx) *ctx = tr->precheckctx;
14841ba4c6cSHeeho Park   PetscFunctionReturn(0);
14941ba4c6cSHeeho Park }
15041ba4c6cSHeeho Park 
15141ba4c6cSHeeho Park /*@C
15241ba4c6cSHeeho Park    SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
15341ba4c6cSHeeho Park        function evaluation. Allows the user a chance to change or override the decision of the line search routine
15441ba4c6cSHeeho Park 
155c3339decSBarry Smith    Logically Collective
15641ba4c6cSHeeho Park 
15741ba4c6cSHeeho Park    Input Parameters:
15841ba4c6cSHeeho Park +  snes - the nonlinear solver object
159f6dfbefdSBarry Smith .  func - [optional] function evaluation routine, see `SNESNewtonTRDCPostCheck()`  for the calling sequence
16041ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
16141ba4c6cSHeeho Park 
16241ba4c6cSHeeho Park    Level: intermediate
16341ba4c6cSHeeho Park 
164f6dfbefdSBarry Smith    Note:
165f6dfbefdSBarry Smith    This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in
166f6dfbefdSBarry Smith    `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation.
16741ba4c6cSHeeho Park 
168f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`
16941ba4c6cSHeeho Park @*/
170d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx)
171d71ae5a4SJacob Faibussowitsch {
17241ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
17341ba4c6cSHeeho Park 
17441ba4c6cSHeeho Park   PetscFunctionBegin;
17541ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
17641ba4c6cSHeeho Park   if (func) tr->postcheck = func;
17741ba4c6cSHeeho Park   if (ctx) tr->postcheckctx = ctx;
17841ba4c6cSHeeho Park   PetscFunctionReturn(0);
17941ba4c6cSHeeho Park }
18041ba4c6cSHeeho Park 
18141ba4c6cSHeeho Park /*@C
18241ba4c6cSHeeho Park    SNESNewtonTRDCGetPostCheck - Gets the post-check function
18341ba4c6cSHeeho Park 
18441ba4c6cSHeeho Park    Not collective
18541ba4c6cSHeeho Park 
18641ba4c6cSHeeho Park    Input Parameter:
18741ba4c6cSHeeho Park .  snes - the nonlinear solver context
18841ba4c6cSHeeho Park 
18941ba4c6cSHeeho Park    Output Parameters:
19041ba4c6cSHeeho Park +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck()
19141ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
19241ba4c6cSHeeho Park 
19341ba4c6cSHeeho Park    Level: intermediate
19441ba4c6cSHeeho Park 
195f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`
19641ba4c6cSHeeho Park @*/
197d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
198d71ae5a4SJacob Faibussowitsch {
19941ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
20041ba4c6cSHeeho Park 
20141ba4c6cSHeeho Park   PetscFunctionBegin;
20241ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
20341ba4c6cSHeeho Park   if (func) *func = tr->postcheck;
20441ba4c6cSHeeho Park   if (ctx) *ctx = tr->postcheckctx;
20541ba4c6cSHeeho Park   PetscFunctionReturn(0);
20641ba4c6cSHeeho Park }
20741ba4c6cSHeeho Park 
20841ba4c6cSHeeho Park /*@C
209f6dfbefdSBarry Smith    SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC`
21041ba4c6cSHeeho Park 
211c3339decSBarry Smith    Logically Collective
21241ba4c6cSHeeho Park 
21341ba4c6cSHeeho Park    Input Parameters:
21441ba4c6cSHeeho Park +  snes - the solver
21541ba4c6cSHeeho Park .  X - The last solution
21641ba4c6cSHeeho Park -  Y - The step direction
21741ba4c6cSHeeho Park 
21841ba4c6cSHeeho Park    Output Parameters:
21941ba4c6cSHeeho Park .  changed_Y - Indicator that the step direction Y has been changed.
22041ba4c6cSHeeho Park 
22141ba4c6cSHeeho Park    Level: developer
22241ba4c6cSHeeho Park 
223f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()`
22441ba4c6cSHeeho Park @*/
225d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y)
226d71ae5a4SJacob Faibussowitsch {
22741ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
22841ba4c6cSHeeho Park 
22941ba4c6cSHeeho Park   PetscFunctionBegin;
23041ba4c6cSHeeho Park   *changed_Y = PETSC_FALSE;
23141ba4c6cSHeeho Park   if (tr->precheck) {
2329566063dSJacob Faibussowitsch     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
23341ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_Y, 4);
23441ba4c6cSHeeho Park   }
23541ba4c6cSHeeho Park   PetscFunctionReturn(0);
23641ba4c6cSHeeho Park }
23741ba4c6cSHeeho Park 
23841ba4c6cSHeeho Park /*@C
239f6dfbefdSBarry Smith    SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step
24041ba4c6cSHeeho Park 
241c3339decSBarry Smith    Logically Collective
24241ba4c6cSHeeho Park 
24341ba4c6cSHeeho Park    Input Parameters:
244f1a722f8SMatthew G. Knepley +  snes - the solver
245f1a722f8SMatthew G. Knepley .  X - The last solution
24641ba4c6cSHeeho Park .  Y - The full step direction
24741ba4c6cSHeeho Park -  W - The updated solution, W = X - Y
24841ba4c6cSHeeho Park 
24941ba4c6cSHeeho Park    Output Parameters:
25041ba4c6cSHeeho Park +  changed_Y - indicator if step has been changed
25141ba4c6cSHeeho Park -  changed_W - Indicator if the new candidate solution W has been changed.
25241ba4c6cSHeeho Park 
253f6dfbefdSBarry Smith    Note:
25441ba4c6cSHeeho Park      If Y is changed then W is recomputed as X - Y
25541ba4c6cSHeeho Park 
25641ba4c6cSHeeho Park    Level: developer
25741ba4c6cSHeeho Park 
258f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck()
25941ba4c6cSHeeho Park @*/
260d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
261d71ae5a4SJacob Faibussowitsch {
26241ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
26341ba4c6cSHeeho Park 
26441ba4c6cSHeeho Park   PetscFunctionBegin;
26541ba4c6cSHeeho Park   *changed_Y = PETSC_FALSE;
26641ba4c6cSHeeho Park   *changed_W = PETSC_FALSE;
26741ba4c6cSHeeho Park   if (tr->postcheck) {
2689566063dSJacob Faibussowitsch     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
26941ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_Y, 5);
27041ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_W, 6);
27141ba4c6cSHeeho Park   }
27241ba4c6cSHeeho Park   PetscFunctionReturn(0);
27341ba4c6cSHeeho Park }
27441ba4c6cSHeeho Park 
27541ba4c6cSHeeho Park /*
27641ba4c6cSHeeho Park    SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
27741ba4c6cSHeeho Park    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
27841ba4c6cSHeeho Park    nonlinear equations
27941ba4c6cSHeeho Park 
28041ba4c6cSHeeho Park */
281d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes)
282d71ae5a4SJacob Faibussowitsch {
28341ba4c6cSHeeho Park   SNES_NEWTONTRDC            *neP = (SNES_NEWTONTRDC *)snes->data;
28441ba4c6cSHeeho Park   Vec                         X, F, Y, G, W, GradF, YNtmp;
28541ba4c6cSHeeho Park   Vec                         YCtmp;
28641ba4c6cSHeeho Park   Mat                         jac;
28741ba4c6cSHeeho Park   PetscInt                    maxits, i, j, lits, inner_count, bs;
28841ba4c6cSHeeho Park   PetscReal                   rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */
28941ba4c6cSHeeho Park   PetscReal                   inorms[99];                                                         /* need to make it dynamic eventually, fixed max block size of 99 for now */
29041ba4c6cSHeeho Park   PetscReal                   deltaM, ynnorm, f0, mp, gTy, g, yTHy;                               /* rho calculation */
29141ba4c6cSHeeho Park   PetscReal                   auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg;       /* Cauchy Point */
29241ba4c6cSHeeho Park   KSP                         ksp;
29341ba4c6cSHeeho Park   SNESConvergedReason         reason   = SNES_CONVERGED_ITERATING;
29441ba4c6cSHeeho Park   PetscBool                   breakout = PETSC_FALSE;
29541ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx;
29641ba4c6cSHeeho Park   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
29741ba4c6cSHeeho Park   void *convctx;
29841ba4c6cSHeeho Park 
29941ba4c6cSHeeho Park   PetscFunctionBegin;
30041ba4c6cSHeeho Park   maxits = snes->max_its;  /* maximum number of iterations */
30141ba4c6cSHeeho Park   X      = snes->vec_sol;  /* solution vector */
30241ba4c6cSHeeho Park   F      = snes->vec_func; /* residual vector */
30341ba4c6cSHeeho Park   Y      = snes->work[0];  /* update vector */
30441ba4c6cSHeeho Park   G      = snes->work[1];  /* updated residual */
30541ba4c6cSHeeho Park   W      = snes->work[2];  /* temporary vector */
30641ba4c6cSHeeho Park   GradF  = snes->work[3];  /* grad f = J^T F */
30741ba4c6cSHeeho Park   YNtmp  = snes->work[4];  /* Newton solution */
30841ba4c6cSHeeho Park   YCtmp  = snes->work[5];  /* Cauchy solution */
30941ba4c6cSHeeho Park 
3100b121fc5SBarry 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);
31141ba4c6cSHeeho Park 
3129566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(YNtmp, &bs));
31341ba4c6cSHeeho Park 
3149566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
31541ba4c6cSHeeho Park   snes->iter = 0;
3169566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
31741ba4c6cSHeeho Park 
31841ba4c6cSHeeho Park   /* Set the linear stopping criteria to use the More' trick. From tr.c */
3199566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
3209566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
32141ba4c6cSHeeho Park   if (convtest != SNESTRDC_KSPConverged_Private) {
3229566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
32341ba4c6cSHeeho Park     ctx->snes = snes;
3249566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
3259566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy));
3269566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n"));
32741ba4c6cSHeeho Park   }
32841ba4c6cSHeeho Park 
32941ba4c6cSHeeho Park   if (!snes->vec_func_init_set) {
3309566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
33141ba4c6cSHeeho Park   } else snes->vec_func_init_set = PETSC_FALSE;
33241ba4c6cSHeeho Park 
3339566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
33441ba4c6cSHeeho Park   SNESCheckFunctionNorm(snes, fnorm);
3359566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
33641ba4c6cSHeeho Park 
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
33841ba4c6cSHeeho Park   snes->norm = fnorm;
3399566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
34041ba4c6cSHeeho Park   delta      = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */
34141ba4c6cSHeeho Park   deltaM     = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */
34241ba4c6cSHeeho Park   neP->delta = delta;
3439566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3449566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
34541ba4c6cSHeeho Park 
34641ba4c6cSHeeho Park   neP->rho_satisfied = PETSC_FALSE;
34741ba4c6cSHeeho Park 
34841ba4c6cSHeeho Park   /* test convergence */
349dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
35041ba4c6cSHeeho Park   if (snes->reason) PetscFunctionReturn(0);
35141ba4c6cSHeeho Park 
35241ba4c6cSHeeho Park   for (i = 0; i < maxits; i++) {
35341ba4c6cSHeeho Park     PetscBool changed_y;
35441ba4c6cSHeeho Park     PetscBool changed_w;
35541ba4c6cSHeeho Park 
35641ba4c6cSHeeho Park     /* dogleg method */
3579566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
35841ba4c6cSHeeho Park     SNESCheckJacobianDomainerror(snes);
3599566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian));
3609566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */
36141ba4c6cSHeeho Park     SNESCheckKSPSolve(snes);                  /* this is necessary but old tr.c did not have it*/
3629566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
3639566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
36441ba4c6cSHeeho Park 
36541ba4c6cSHeeho Park     /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable)
36641ba4c6cSHeeho Park        for inner iteration and Cauchy direction calculation
36741ba4c6cSHeeho Park     */
36841ba4c6cSHeeho Park     if (bs > 1 && neP->auto_scale_multiphase) {
3699566063dSJacob Faibussowitsch       PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms));
37041ba4c6cSHeeho Park       for (j = 0; j < bs; j++) {
37141ba4c6cSHeeho Park         if (neP->auto_scale_max > 1.0) {
372ad540459SPierre Jolivet           if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max;
37341ba4c6cSHeeho Park         }
3749566063dSJacob Faibussowitsch         PetscCall(VecStrideSet(W, j, inorms[j]));
3759566063dSJacob Faibussowitsch         PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j]));
3769566063dSJacob Faibussowitsch         PetscCall(VecStrideScale(X, j, 1.0 / inorms[j]));
37741ba4c6cSHeeho Park       }
3789566063dSJacob Faibussowitsch       PetscCall(VecNorm(X, NORM_2, &xnorm));
37941ba4c6cSHeeho Park       if (i == 0) {
38041ba4c6cSHeeho Park         delta = neP->delta0 * xnorm;
38141ba4c6cSHeeho Park       } else {
38241ba4c6cSHeeho Park         delta = neP->delta * xnorm;
38341ba4c6cSHeeho Park       }
38441ba4c6cSHeeho Park       deltaM = neP->deltaM * xnorm;
385*f3fa974cSJacob Faibussowitsch       PetscCall(MatDiagonalScale(jac, NULL, W));
38641ba4c6cSHeeho Park     }
38741ba4c6cSHeeho Park 
38841ba4c6cSHeeho Park     /* calculating GradF of minimization function */
3899566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */
3909566063dSJacob Faibussowitsch     PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */
39141ba4c6cSHeeho Park 
39241ba4c6cSHeeho Park     inner_count        = 0;
39341ba4c6cSHeeho Park     neP->rho_satisfied = PETSC_FALSE;
39441ba4c6cSHeeho Park     while (1) {
39541ba4c6cSHeeho Park       if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */
3969566063dSJacob Faibussowitsch         PetscCall(VecCopy(YNtmp, Y));
39741ba4c6cSHeeho Park       } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */
3989566063dSJacob Faibussowitsch         PetscCall(MatMult(jac, GradF, W));
3999566063dSJacob Faibussowitsch         PetscCall(VecDotRealPart(W, W, &gTBg));     /* completes GradF^T J^T J GradF */
4009566063dSJacob Faibussowitsch         PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */
40141ba4c6cSHeeho Park         if (gTBg <= 0.0) {
40241ba4c6cSHeeho Park           auk = PETSC_MAX_REAL;
40341ba4c6cSHeeho Park         } else {
40441ba4c6cSHeeho Park           auk = PetscSqr(gfnorm) / gTBg;
40541ba4c6cSHeeho Park         }
40641ba4c6cSHeeho Park         auk = PetscMin(delta / gfnorm, auk);
4079566063dSJacob Faibussowitsch         PetscCall(VecCopy(GradF, YCtmp));           /* this could be improved */
4089566063dSJacob Faibussowitsch         PetscCall(VecScale(YCtmp, auk));            /* YCtmp, Cauchy solution*/
4099566063dSJacob Faibussowitsch         PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */
41041ba4c6cSHeeho Park         if (ycnorm >= delta) {                      /* see if the Cauchy solution meets the criteria */
4119566063dSJacob Faibussowitsch           PetscCall(VecCopy(YCtmp, Y));
4129566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm));
41341ba4c6cSHeeho Park         } else {                                  /* take ratio, tau, of Cauchy and Newton direction and step */
4149566063dSJacob Faibussowitsch           PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */
4159566063dSJacob Faibussowitsch           PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */
41641ba4c6cSHeeho Park           c0 = PetscSqr(c0);
4179566063dSJacob Faibussowitsch           PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1));
41841ba4c6cSHeeho Park           c1 = 2.0 * c1;
4199566063dSJacob Faibussowitsch           PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */
42041ba4c6cSHeeho Park           c2      = PetscSqr(c2) - PetscSqr(delta);
42141ba4c6cSHeeho Park           tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */
42241ba4c6cSHeeho Park           tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0);
42341ba4c6cSHeeho Park           tau     = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */
4249566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm));
4259566063dSJacob Faibussowitsch           PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp));
4269566063dSJacob Faibussowitsch           PetscCall(VecAXPY(W, -tau, YCtmp));
4279566063dSJacob Faibussowitsch           PetscCall(VecCopy(W, Y)); /* this could be improved */
42841ba4c6cSHeeho Park         }
42941ba4c6cSHeeho Park       } else {
43041ba4c6cSHeeho Park         /* if Cauchy is disabled, only use Newton direction */
43141ba4c6cSHeeho Park         auk = delta / ynnorm;
4329566063dSJacob Faibussowitsch         PetscCall(VecScale(YNtmp, auk));
4339566063dSJacob Faibussowitsch         PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/
43441ba4c6cSHeeho Park       }
43541ba4c6cSHeeho Park 
4369566063dSJacob Faibussowitsch       PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm  */
43741ba4c6cSHeeho Park       f0 = 0.5 * PetscSqr(fnorm);            /* minimizing function f(X) */
4389566063dSJacob Faibussowitsch       PetscCall(MatMult(jac, Y, W));
4399566063dSJacob Faibussowitsch       PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */
4409566063dSJacob Faibussowitsch       PetscCall(VecDotRealPart(GradF, Y, &gTy));
44141ba4c6cSHeeho Park       mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/
44241ba4c6cSHeeho Park 
44341ba4c6cSHeeho Park       /* scale back solution update */
44441ba4c6cSHeeho Park       if (bs > 1 && neP->auto_scale_multiphase) {
44541ba4c6cSHeeho Park         for (j = 0; j < bs; j++) {
4469566063dSJacob Faibussowitsch           PetscCall(VecStrideScale(Y, j, inorms[j]));
44741ba4c6cSHeeho Park           if (inner_count == 0) {
44841ba4c6cSHeeho Park             /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */
44941ba4c6cSHeeho Park             /* need to scale back X to match Y and provide proper update to the external code */
4509566063dSJacob Faibussowitsch             PetscCall(VecStrideScale(X, j, inorms[j]));
45141ba4c6cSHeeho Park           }
45241ba4c6cSHeeho Park         }
4539566063dSJacob Faibussowitsch         if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */
4549566063dSJacob Faibussowitsch         PetscCall(VecNorm(Y, NORM_2, &temp_ynorm));
45541ba4c6cSHeeho Park       } else {
45641ba4c6cSHeeho Park         temp_xnorm = xnorm;
45741ba4c6cSHeeho Park         temp_ynorm = ynorm;
45841ba4c6cSHeeho Park       }
45941ba4c6cSHeeho Park       inner_count++;
46041ba4c6cSHeeho Park 
46141ba4c6cSHeeho Park       /* Evaluate the solution to meet the improvement ratio criteria */
4629566063dSJacob Faibussowitsch       PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y));
4639566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(W, -1.0, Y, X));
4649566063dSJacob Faibussowitsch       PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w));
4659566063dSJacob Faibussowitsch       if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X));
4669566063dSJacob Faibussowitsch       PetscCall(VecCopy(Y, snes->vec_sol_update));
4679566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, W, G)); /*  F(X-Y) = G */
4689566063dSJacob Faibussowitsch       PetscCall(VecNorm(G, NORM_2, &gnorm));      /* gnorm <- || g || */
46941ba4c6cSHeeho Park       SNESCheckFunctionNorm(snes, gnorm);
47041ba4c6cSHeeho Park       g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */
47141ba4c6cSHeeho Park       if (f0 == mp) rho = 0.0;
47241ba4c6cSHeeho Park       else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */
47341ba4c6cSHeeho Park 
47441ba4c6cSHeeho Park       if (rho < neP->eta2) {
47541ba4c6cSHeeho Park         delta *= neP->t1; /* shrink the region */
47641ba4c6cSHeeho Park       } else if (rho > neP->eta3) {
47741ba4c6cSHeeho Park         delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */
47841ba4c6cSHeeho Park       }
47941ba4c6cSHeeho Park 
48041ba4c6cSHeeho Park       neP->delta = delta;
48141ba4c6cSHeeho Park       if (rho >= neP->eta1) {
48241ba4c6cSHeeho Park         /* unscale delta and xnorm before going to the next outer iteration */
48341ba4c6cSHeeho Park         if (bs > 1 && neP->auto_scale_multiphase) {
48441ba4c6cSHeeho Park           neP->delta = delta / xnorm;
48541ba4c6cSHeeho Park           xnorm      = temp_xnorm;
48641ba4c6cSHeeho Park           ynorm      = temp_ynorm;
48741ba4c6cSHeeho Park         }
48841ba4c6cSHeeho Park         neP->rho_satisfied = PETSC_TRUE;
48941ba4c6cSHeeho Park         break; /* the improvement ratio is satisfactory */
49041ba4c6cSHeeho Park       }
4919566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
49241ba4c6cSHeeho Park 
49341ba4c6cSHeeho Park       /* check to see if progress is hopeless */
49441ba4c6cSHeeho Park       neP->itflag = PETSC_FALSE;
49541ba4c6cSHeeho Park       /* both delta, ynorm, and xnorm are either scaled or unscaled */
4969566063dSJacob Faibussowitsch       PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP));
49741ba4c6cSHeeho Park       if (!reason) {
49841ba4c6cSHeeho Park         /* temp_xnorm, temp_ynorm is always unscaled */
49941ba4c6cSHeeho Park         /* also the inner iteration already calculated the Jacobian and solved the matrix */
50041ba4c6cSHeeho Park         /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */
5019566063dSJacob Faibussowitsch         PetscCall((*snes->ops->converged)(snes, snes->iter + 1, temp_xnorm, temp_ynorm, fnorm, &reason, snes->cnvP));
50241ba4c6cSHeeho Park       }
50341ba4c6cSHeeho Park       /* if multiphase state changes, break out inner iteration */
50441ba4c6cSHeeho Park       if (reason == SNES_BREAKOUT_INNER_ITER) {
50541ba4c6cSHeeho Park         if (bs > 1 && neP->auto_scale_multiphase) {
50641ba4c6cSHeeho Park           /* unscale delta and xnorm before going to the next outer iteration */
50741ba4c6cSHeeho Park           neP->delta = delta / xnorm;
50841ba4c6cSHeeho Park           xnorm      = temp_xnorm;
50941ba4c6cSHeeho Park           ynorm      = temp_ynorm;
51041ba4c6cSHeeho Park         }
51141ba4c6cSHeeho Park         reason = SNES_CONVERGED_ITERATING;
51241ba4c6cSHeeho Park         break;
51341ba4c6cSHeeho Park       }
51441ba4c6cSHeeho Park       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
51541ba4c6cSHeeho Park       if (reason) {
51641ba4c6cSHeeho Park         if (reason < 0) {
51741ba4c6cSHeeho Park           /* We're not progressing, so return with the current iterate */
5189566063dSJacob Faibussowitsch           PetscCall(SNESMonitor(snes, i + 1, fnorm));
51941ba4c6cSHeeho Park           breakout = PETSC_TRUE;
52041ba4c6cSHeeho Park           break;
52141ba4c6cSHeeho Park         } else if (reason > 0) {
52241ba4c6cSHeeho Park           /* We're converged, so return with the current iterate and update solution */
5239566063dSJacob Faibussowitsch           PetscCall(SNESMonitor(snes, i + 1, fnorm));
52441ba4c6cSHeeho Park           breakout = PETSC_FALSE;
52541ba4c6cSHeeho Park           break;
52641ba4c6cSHeeho Park         }
52741ba4c6cSHeeho Park       }
52841ba4c6cSHeeho Park       snes->numFailures++;
52941ba4c6cSHeeho Park     }
53041ba4c6cSHeeho Park     if (!breakout) {
53141ba4c6cSHeeho Park       /* Update function and solution vectors */
53241ba4c6cSHeeho Park       fnorm = gnorm;
5339566063dSJacob Faibussowitsch       PetscCall(VecCopy(G, F));
5349566063dSJacob Faibussowitsch       PetscCall(VecCopy(W, X));
53541ba4c6cSHeeho Park       /* Monitor convergence */
5369566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
53741ba4c6cSHeeho Park       snes->iter  = i + 1;
53841ba4c6cSHeeho Park       snes->norm  = fnorm;
53941ba4c6cSHeeho Park       snes->xnorm = xnorm;
54041ba4c6cSHeeho Park       snes->ynorm = ynorm;
5419566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5429566063dSJacob Faibussowitsch       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
5439566063dSJacob Faibussowitsch       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
54441ba4c6cSHeeho Park       /* Test for convergence, xnorm = || X || */
54541ba4c6cSHeeho Park       neP->itflag = PETSC_TRUE;
5469566063dSJacob Faibussowitsch       if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
547dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
54841ba4c6cSHeeho Park       if (reason) break;
54941ba4c6cSHeeho Park     } else break;
55041ba4c6cSHeeho Park   }
55141ba4c6cSHeeho Park 
5529566063dSJacob Faibussowitsch   /* PetscCall(PetscFree(inorms)); */
55341ba4c6cSHeeho Park   if (i == maxits) {
5549566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
55541ba4c6cSHeeho Park     if (!reason) reason = SNES_DIVERGED_MAX_IT;
55641ba4c6cSHeeho Park   }
5579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
55841ba4c6cSHeeho Park   snes->reason = reason;
5599566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
56041ba4c6cSHeeho Park   if (convtest != SNESTRDC_KSPConverged_Private) {
5619566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
5629566063dSJacob Faibussowitsch     PetscCall(PetscFree(ctx));
5639566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
56441ba4c6cSHeeho Park   }
56541ba4c6cSHeeho Park   PetscFunctionReturn(0);
56641ba4c6cSHeeho Park }
56741ba4c6cSHeeho Park 
568d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes)
569d71ae5a4SJacob Faibussowitsch {
57041ba4c6cSHeeho Park   PetscFunctionBegin;
5719566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 6));
5729566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
57341ba4c6cSHeeho Park   PetscFunctionReturn(0);
57441ba4c6cSHeeho Park }
57541ba4c6cSHeeho Park 
576d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONTRDC(SNES snes)
577d71ae5a4SJacob Faibussowitsch {
57841ba4c6cSHeeho Park   PetscFunctionBegin;
57941ba4c6cSHeeho Park   PetscFunctionReturn(0);
58041ba4c6cSHeeho Park }
58141ba4c6cSHeeho Park 
582d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes)
583d71ae5a4SJacob Faibussowitsch {
58441ba4c6cSHeeho Park   PetscFunctionBegin;
5859566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONTRDC(snes));
5869566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
58741ba4c6cSHeeho Park   PetscFunctionReturn(0);
58841ba4c6cSHeeho Park }
58941ba4c6cSHeeho Park 
590d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject)
591d71ae5a4SJacob Faibussowitsch {
59241ba4c6cSHeeho Park   SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data;
59341ba4c6cSHeeho Park 
59441ba4c6cSHeeho Park   PetscFunctionBegin;
595d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
5969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
5979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
5989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
5999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
6009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
6019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
6029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
6039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
6049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL));
6059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL));
6069566063dSJacob 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));
607d0609cedSBarry Smith   PetscOptionsHeadEnd();
60841ba4c6cSHeeho Park   PetscFunctionReturn(0);
60941ba4c6cSHeeho Park }
61041ba4c6cSHeeho Park 
611d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer)
612d71ae5a4SJacob Faibussowitsch {
61341ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
61441ba4c6cSHeeho Park   PetscBool        iascii;
61541ba4c6cSHeeho Park 
61641ba4c6cSHeeho Park   PetscFunctionBegin;
6179566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
61841ba4c6cSHeeho Park   if (iascii) {
61963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol));
6209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
6219566063dSJacob 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));
62241ba4c6cSHeeho Park   }
62341ba4c6cSHeeho Park   PetscFunctionReturn(0);
62441ba4c6cSHeeho Park }
625f6dfbefdSBarry Smith 
62641ba4c6cSHeeho Park /*MC
62741ba4c6cSHeeho Park       SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
62841ba4c6cSHeeho Park 
629f6dfbefdSBarry Smith    Options Database Keys:
63041ba4c6cSHeeho Park +   -snes_trdc_tol <tol> - trust region tolerance
63141ba4c6cSHeeho Park .   -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
63241ba4c6cSHeeho Park .   -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
63341ba4c6cSHeeho Park .   -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
63441ba4c6cSHeeho Park .   -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
63541ba4c6cSHeeho Park .   -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
63641ba4c6cSHeeho Park .   -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5)
63741ba4c6cSHeeho Park .   -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1)
63841ba4c6cSHeeho Park .   -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor
63941ba4c6cSHeeho Park .   -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm
64041ba4c6cSHeeho Park -   -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region
64141ba4c6cSHeeho Park 
642f6dfbefdSBarry Smith     Reference:
643f6dfbefdSBarry Smith .   - * "Linear and Nonlinear Solvers for Simulating Multiphase Flow
64441ba4c6cSHeeho Park     within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond,
64541ba4c6cSHeeho Park     Albert J. Valocchi, Tara LaForce.
64641ba4c6cSHeeho Park 
64741ba4c6cSHeeho Park    Level: intermediate
64841ba4c6cSHeeho Park 
649f6dfbefdSBarry Smith .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`,
650f6dfbefdSBarry Smith           `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`,
651f6dfbefdSBarry Smith           `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()`
65241ba4c6cSHeeho Park M*/
653d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes)
654d71ae5a4SJacob Faibussowitsch {
65541ba4c6cSHeeho Park   SNES_NEWTONTRDC *neP;
65641ba4c6cSHeeho Park 
65741ba4c6cSHeeho Park   PetscFunctionBegin;
65841ba4c6cSHeeho Park   snes->ops->setup          = SNESSetUp_NEWTONTRDC;
65941ba4c6cSHeeho Park   snes->ops->solve          = SNESSolve_NEWTONTRDC;
66041ba4c6cSHeeho Park   snes->ops->destroy        = SNESDestroy_NEWTONTRDC;
66141ba4c6cSHeeho Park   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC;
66241ba4c6cSHeeho Park   snes->ops->view           = SNESView_NEWTONTRDC;
66341ba4c6cSHeeho Park   snes->ops->reset          = SNESReset_NEWTONTRDC;
66441ba4c6cSHeeho Park 
66541ba4c6cSHeeho Park   snes->usesksp = PETSC_TRUE;
66641ba4c6cSHeeho Park   snes->usesnpc = PETSC_FALSE;
66741ba4c6cSHeeho Park 
66841ba4c6cSHeeho Park   snes->alwayscomputesfinalresidual = PETSC_TRUE;
66941ba4c6cSHeeho Park 
6704dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
67141ba4c6cSHeeho Park   snes->data                 = (void *)neP;
67241ba4c6cSHeeho Park   neP->delta                 = 0.0;
67341ba4c6cSHeeho Park   neP->delta0                = 0.1;
67441ba4c6cSHeeho Park   neP->eta1                  = 0.001;
67541ba4c6cSHeeho Park   neP->eta2                  = 0.25;
67641ba4c6cSHeeho Park   neP->eta3                  = 0.75;
67741ba4c6cSHeeho Park   neP->t1                    = 0.25;
67841ba4c6cSHeeho Park   neP->t2                    = 2.0;
67941ba4c6cSHeeho Park   neP->deltaM                = 0.5;
68041ba4c6cSHeeho Park   neP->sigma                 = 0.0001;
68141ba4c6cSHeeho Park   neP->itflag                = PETSC_FALSE;
68241ba4c6cSHeeho Park   neP->rnorm0                = 0.0;
68341ba4c6cSHeeho Park   neP->ttol                  = 0.0;
68441ba4c6cSHeeho Park   neP->use_cauchy            = PETSC_TRUE;
68541ba4c6cSHeeho Park   neP->auto_scale_multiphase = PETSC_FALSE;
68641ba4c6cSHeeho Park   neP->auto_scale_max        = -1.0;
68741ba4c6cSHeeho Park   neP->rho_satisfied         = PETSC_FALSE;
68841ba4c6cSHeeho Park   snes->deltatol             = 1.e-12;
68941ba4c6cSHeeho Park 
69041ba4c6cSHeeho Park   /* for multiphase (multivariable) scaling */
69141ba4c6cSHeeho Park   /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13
69241ba4c6cSHeeho Park      on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now.
6939566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(snes->work[0],&neP->bs));
6949566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(neP->bs,&neP->inorms));
69541ba4c6cSHeeho Park   */
69641ba4c6cSHeeho Park 
69741ba4c6cSHeeho Park   PetscFunctionReturn(0);
69841ba4c6cSHeeho Park }
699