xref: /petsc/src/snes/impls/ntrdc/ntrdc.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
169371c9d4SSatish Balay static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) {
1741ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx  = (SNES_TRDC_KSPConverged_Ctx *)cctx;
1841ba4c6cSHeeho Park   SNES                        snes = ctx->snes;
1941ba4c6cSHeeho Park   SNES_NEWTONTRDC            *neP  = (SNES_NEWTONTRDC *)snes->data;
2041ba4c6cSHeeho Park   Vec                         x;
2141ba4c6cSHeeho Park   PetscReal                   nrm;
2241ba4c6cSHeeho Park 
2341ba4c6cSHeeho Park   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
25*48a46eb9SPierre Jolivet   if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm));
2641ba4c6cSHeeho Park   /* Determine norm of solution */
279566063dSJacob Faibussowitsch   PetscCall(KSPBuildSolution(ksp, NULL, &x));
289566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &nrm));
2941ba4c6cSHeeho Park   if (nrm >= neP->delta) {
309566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
3141ba4c6cSHeeho Park     *reason = KSP_CONVERGED_STEP_LENGTH;
3241ba4c6cSHeeho Park   }
3341ba4c6cSHeeho Park   PetscFunctionReturn(0);
3441ba4c6cSHeeho Park }
3541ba4c6cSHeeho Park 
369371c9d4SSatish Balay static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) {
3741ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx;
3841ba4c6cSHeeho Park 
3941ba4c6cSHeeho Park   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall((*ctx->convdestroy)(ctx->convctx));
419566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
4241ba4c6cSHeeho Park 
4341ba4c6cSHeeho Park   PetscFunctionReturn(0);
4441ba4c6cSHeeho Park }
4541ba4c6cSHeeho Park 
4641ba4c6cSHeeho Park /* ---------------------------------------------------------------- */
4741ba4c6cSHeeho Park /*
4841ba4c6cSHeeho Park    SNESTRDC_Converged_Private -test convergence JUST for
4941ba4c6cSHeeho Park    the trust region tolerance.
5041ba4c6cSHeeho Park 
5141ba4c6cSHeeho Park */
529371c9d4SSatish Balay static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) {
5341ba4c6cSHeeho Park   SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
5441ba4c6cSHeeho Park 
5541ba4c6cSHeeho Park   PetscFunctionBegin;
5641ba4c6cSHeeho Park   *reason = SNES_CONVERGED_ITERATING;
5741ba4c6cSHeeho Park   if (neP->delta < xnorm * snes->deltatol) {
589566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol));
5941ba4c6cSHeeho Park     *reason = SNES_DIVERGED_TR_DELTA;
6041ba4c6cSHeeho Park   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
619566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
6241ba4c6cSHeeho Park     *reason = SNES_DIVERGED_FUNCTION_COUNT;
6341ba4c6cSHeeho Park   }
6441ba4c6cSHeeho Park   PetscFunctionReturn(0);
6541ba4c6cSHeeho Park }
6641ba4c6cSHeeho Park 
6741ba4c6cSHeeho Park /*@
6841ba4c6cSHeeho Park   SNESNewtonTRDCGetRhoFlag - Get whether the solution update is within the trust-region.
6941ba4c6cSHeeho Park 
7041ba4c6cSHeeho Park   Input Parameters:
7141ba4c6cSHeeho Park . snes - the nonlinear solver object
7241ba4c6cSHeeho Park 
7341ba4c6cSHeeho Park   Output Parameters:
7441ba4c6cSHeeho Park . rho_flag: PETSC_TRUE if the solution update is in the trust-region; otherwise, PETSC_FALSE
7541ba4c6cSHeeho Park 
7641ba4c6cSHeeho Park   Level: developer
7741ba4c6cSHeeho Park 
7841ba4c6cSHeeho Park @*/
799371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag) {
8041ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
8141ba4c6cSHeeho Park 
8241ba4c6cSHeeho Park   PetscFunctionBegin;
8341ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
8441ba4c6cSHeeho Park   PetscValidBoolPointer(rho_flag, 2);
8541ba4c6cSHeeho Park   *rho_flag = tr->rho_satisfied;
8641ba4c6cSHeeho Park   PetscFunctionReturn(0);
8741ba4c6cSHeeho Park }
8841ba4c6cSHeeho Park 
8941ba4c6cSHeeho Park /*@C
9041ba4c6cSHeeho Park    SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined.
9141ba4c6cSHeeho Park        Allows the user a chance to change or override the trust region decision.
9241ba4c6cSHeeho Park 
9341ba4c6cSHeeho Park    Logically Collective on snes
9441ba4c6cSHeeho Park 
9541ba4c6cSHeeho Park    Input Parameters:
9641ba4c6cSHeeho Park +  snes - the nonlinear solver object
9741ba4c6cSHeeho Park .  func - [optional] function evaluation routine, see SNESNewtonTRDCPreCheck()  for the calling sequence
9841ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
9941ba4c6cSHeeho Park 
10041ba4c6cSHeeho Park    Level: intermediate
10141ba4c6cSHeeho Park 
10241ba4c6cSHeeho Park    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver.
10341ba4c6cSHeeho Park 
104db781477SPatrick Sanan .seealso: `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
10541ba4c6cSHeeho Park @*/
1069371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) {
10741ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
10841ba4c6cSHeeho Park 
10941ba4c6cSHeeho Park   PetscFunctionBegin;
11041ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
11141ba4c6cSHeeho Park   if (func) tr->precheck = func;
11241ba4c6cSHeeho Park   if (ctx) tr->precheckctx = ctx;
11341ba4c6cSHeeho Park   PetscFunctionReturn(0);
11441ba4c6cSHeeho Park }
11541ba4c6cSHeeho Park 
11641ba4c6cSHeeho Park /*@C
11741ba4c6cSHeeho Park    SNESNewtonTRDCGetPreCheck - Gets the pre-check function
11841ba4c6cSHeeho Park 
11941ba4c6cSHeeho Park    Not collective
12041ba4c6cSHeeho Park 
12141ba4c6cSHeeho Park    Input Parameter:
12241ba4c6cSHeeho Park .  snes - the nonlinear solver context
12341ba4c6cSHeeho Park 
12441ba4c6cSHeeho Park    Output Parameters:
12541ba4c6cSHeeho Park +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPreCheck()
12641ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
12741ba4c6cSHeeho Park 
12841ba4c6cSHeeho Park    Level: intermediate
12941ba4c6cSHeeho Park 
130db781477SPatrick Sanan .seealso: `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()`
13141ba4c6cSHeeho Park @*/
1329371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) {
13341ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
13441ba4c6cSHeeho Park 
13541ba4c6cSHeeho Park   PetscFunctionBegin;
13641ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
13741ba4c6cSHeeho Park   if (func) *func = tr->precheck;
13841ba4c6cSHeeho Park   if (ctx) *ctx = tr->precheckctx;
13941ba4c6cSHeeho Park   PetscFunctionReturn(0);
14041ba4c6cSHeeho Park }
14141ba4c6cSHeeho Park 
14241ba4c6cSHeeho Park /*@C
14341ba4c6cSHeeho Park    SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
14441ba4c6cSHeeho Park        function evaluation. Allows the user a chance to change or override the decision of the line search routine
14541ba4c6cSHeeho Park 
14641ba4c6cSHeeho Park    Logically Collective on snes
14741ba4c6cSHeeho Park 
14841ba4c6cSHeeho Park    Input Parameters:
14941ba4c6cSHeeho Park +  snes - the nonlinear solver object
15041ba4c6cSHeeho Park .  func - [optional] function evaluation routine, see SNESNewtonTRDCPostCheck()  for the calling sequence
15141ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
15241ba4c6cSHeeho Park 
15341ba4c6cSHeeho Park    Level: intermediate
15441ba4c6cSHeeho Park 
15541ba4c6cSHeeho Park    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver while the function set in
15641ba4c6cSHeeho Park    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
15741ba4c6cSHeeho Park 
158db781477SPatrick Sanan .seealso: `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
15941ba4c6cSHeeho Park @*/
1609371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) {
16141ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
16241ba4c6cSHeeho Park 
16341ba4c6cSHeeho Park   PetscFunctionBegin;
16441ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
16541ba4c6cSHeeho Park   if (func) tr->postcheck = func;
16641ba4c6cSHeeho Park   if (ctx) tr->postcheckctx = ctx;
16741ba4c6cSHeeho Park   PetscFunctionReturn(0);
16841ba4c6cSHeeho Park }
16941ba4c6cSHeeho Park 
17041ba4c6cSHeeho Park /*@C
17141ba4c6cSHeeho Park    SNESNewtonTRDCGetPostCheck - Gets the post-check function
17241ba4c6cSHeeho Park 
17341ba4c6cSHeeho Park    Not collective
17441ba4c6cSHeeho Park 
17541ba4c6cSHeeho Park    Input Parameter:
17641ba4c6cSHeeho Park .  snes - the nonlinear solver context
17741ba4c6cSHeeho Park 
17841ba4c6cSHeeho Park    Output Parameters:
17941ba4c6cSHeeho Park +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck()
18041ba4c6cSHeeho Park -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
18141ba4c6cSHeeho Park 
18241ba4c6cSHeeho Park    Level: intermediate
18341ba4c6cSHeeho Park 
184db781477SPatrick Sanan .seealso: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`
18541ba4c6cSHeeho Park @*/
1869371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) {
18741ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
18841ba4c6cSHeeho Park 
18941ba4c6cSHeeho Park   PetscFunctionBegin;
19041ba4c6cSHeeho Park   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
19141ba4c6cSHeeho Park   if (func) *func = tr->postcheck;
19241ba4c6cSHeeho Park   if (ctx) *ctx = tr->postcheckctx;
19341ba4c6cSHeeho Park   PetscFunctionReturn(0);
19441ba4c6cSHeeho Park }
19541ba4c6cSHeeho Park 
19641ba4c6cSHeeho Park /*@C
19741ba4c6cSHeeho Park    SNESNewtonTRDCPreCheck - Called before the step has been determined in SNESNEWTONTRDC
19841ba4c6cSHeeho Park 
19941ba4c6cSHeeho Park    Logically Collective on snes
20041ba4c6cSHeeho Park 
20141ba4c6cSHeeho Park    Input Parameters:
20241ba4c6cSHeeho Park +  snes - the solver
20341ba4c6cSHeeho Park .  X - The last solution
20441ba4c6cSHeeho Park -  Y - The step direction
20541ba4c6cSHeeho Park 
20641ba4c6cSHeeho Park    Output Parameters:
20741ba4c6cSHeeho Park .  changed_Y - Indicator that the step direction Y has been changed.
20841ba4c6cSHeeho Park 
20941ba4c6cSHeeho Park    Level: developer
21041ba4c6cSHeeho Park 
211db781477SPatrick Sanan .seealso: `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`
21241ba4c6cSHeeho Park @*/
2139371c9d4SSatish Balay static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) {
21441ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
21541ba4c6cSHeeho Park 
21641ba4c6cSHeeho Park   PetscFunctionBegin;
21741ba4c6cSHeeho Park   *changed_Y = PETSC_FALSE;
21841ba4c6cSHeeho Park   if (tr->precheck) {
2199566063dSJacob Faibussowitsch     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
22041ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_Y, 4);
22141ba4c6cSHeeho Park   }
22241ba4c6cSHeeho Park   PetscFunctionReturn(0);
22341ba4c6cSHeeho Park }
22441ba4c6cSHeeho Park 
22541ba4c6cSHeeho Park /*@C
22641ba4c6cSHeeho Park    SNESNewtonTRDCPostCheck - Called after the step has been determined in SNESNEWTONTRDC but before the function evaluation
22741ba4c6cSHeeho Park 
22841ba4c6cSHeeho Park    Logically Collective on snes
22941ba4c6cSHeeho Park 
23041ba4c6cSHeeho Park    Input Parameters:
231f1a722f8SMatthew G. Knepley +  snes - the solver
232f1a722f8SMatthew G. Knepley .  X - The last solution
23341ba4c6cSHeeho Park .  Y - The full step direction
23441ba4c6cSHeeho Park -  W - The updated solution, W = X - Y
23541ba4c6cSHeeho Park 
23641ba4c6cSHeeho Park    Output Parameters:
23741ba4c6cSHeeho Park +  changed_Y - indicator if step has been changed
23841ba4c6cSHeeho Park -  changed_W - Indicator if the new candidate solution W has been changed.
23941ba4c6cSHeeho Park 
24041ba4c6cSHeeho Park    Notes:
24141ba4c6cSHeeho Park      If Y is changed then W is recomputed as X - Y
24241ba4c6cSHeeho Park 
24341ba4c6cSHeeho Park    Level: developer
24441ba4c6cSHeeho Park 
245db781477SPatrick Sanan .seealso: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
24641ba4c6cSHeeho Park @*/
2479371c9d4SSatish Balay static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) {
24841ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
24941ba4c6cSHeeho Park 
25041ba4c6cSHeeho Park   PetscFunctionBegin;
25141ba4c6cSHeeho Park   *changed_Y = PETSC_FALSE;
25241ba4c6cSHeeho Park   *changed_W = PETSC_FALSE;
25341ba4c6cSHeeho Park   if (tr->postcheck) {
2549566063dSJacob Faibussowitsch     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
25541ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_Y, 5);
25641ba4c6cSHeeho Park     PetscValidLogicalCollectiveBool(snes, *changed_W, 6);
25741ba4c6cSHeeho Park   }
25841ba4c6cSHeeho Park   PetscFunctionReturn(0);
25941ba4c6cSHeeho Park }
26041ba4c6cSHeeho Park 
26141ba4c6cSHeeho Park /*
26241ba4c6cSHeeho Park    SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
26341ba4c6cSHeeho Park    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
26441ba4c6cSHeeho Park    nonlinear equations
26541ba4c6cSHeeho Park 
26641ba4c6cSHeeho Park */
2679371c9d4SSatish Balay static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) {
26841ba4c6cSHeeho Park   SNES_NEWTONTRDC            *neP = (SNES_NEWTONTRDC *)snes->data;
26941ba4c6cSHeeho Park   Vec                         X, F, Y, G, W, GradF, YNtmp;
27041ba4c6cSHeeho Park   Vec                         YCtmp;
27141ba4c6cSHeeho Park   Mat                         jac;
27241ba4c6cSHeeho Park   PetscInt                    maxits, i, j, lits, inner_count, bs;
27341ba4c6cSHeeho Park   PetscReal                   rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */
27441ba4c6cSHeeho Park   PetscReal                   inorms[99];                                                         /* need to make it dynamic eventually, fixed max block size of 99 for now */
27541ba4c6cSHeeho Park   PetscReal                   deltaM, ynnorm, f0, mp, gTy, g, yTHy;                               /* rho calculation */
27641ba4c6cSHeeho Park   PetscReal                   auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg;       /* Cauchy Point */
27741ba4c6cSHeeho Park   KSP                         ksp;
27841ba4c6cSHeeho Park   SNESConvergedReason         reason   = SNES_CONVERGED_ITERATING;
27941ba4c6cSHeeho Park   PetscBool                   breakout = PETSC_FALSE;
28041ba4c6cSHeeho Park   SNES_TRDC_KSPConverged_Ctx *ctx;
28141ba4c6cSHeeho Park   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
28241ba4c6cSHeeho Park   void *convctx;
28341ba4c6cSHeeho Park 
28441ba4c6cSHeeho Park   PetscFunctionBegin;
28541ba4c6cSHeeho Park   maxits = snes->max_its;  /* maximum number of iterations */
28641ba4c6cSHeeho Park   X      = snes->vec_sol;  /* solution vector */
28741ba4c6cSHeeho Park   F      = snes->vec_func; /* residual vector */
28841ba4c6cSHeeho Park   Y      = snes->work[0];  /* update vector */
28941ba4c6cSHeeho Park   G      = snes->work[1];  /* updated residual */
29041ba4c6cSHeeho Park   W      = snes->work[2];  /* temporary vector */
29141ba4c6cSHeeho Park   GradF  = snes->work[3];  /* grad f = J^T F */
29241ba4c6cSHeeho Park   YNtmp  = snes->work[4];  /* Newton solution */
29341ba4c6cSHeeho Park   YCtmp  = snes->work[5];  /* Cauchy solution */
29441ba4c6cSHeeho Park 
2950b121fc5SBarry 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);
29641ba4c6cSHeeho Park 
2979566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(YNtmp, &bs));
29841ba4c6cSHeeho Park 
2999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
30041ba4c6cSHeeho Park   snes->iter = 0;
3019566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
30241ba4c6cSHeeho Park 
30341ba4c6cSHeeho Park   /* Set the linear stopping criteria to use the More' trick. From tr.c */
3049566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
3059566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
30641ba4c6cSHeeho Park   if (convtest != SNESTRDC_KSPConverged_Private) {
3079566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
30841ba4c6cSHeeho Park     ctx->snes = snes;
3099566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
3109566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy));
3119566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n"));
31241ba4c6cSHeeho Park   }
31341ba4c6cSHeeho Park 
31441ba4c6cSHeeho Park   if (!snes->vec_func_init_set) {
3159566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
31641ba4c6cSHeeho Park   } else snes->vec_func_init_set = PETSC_FALSE;
31741ba4c6cSHeeho Park 
3189566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
31941ba4c6cSHeeho Park   SNESCheckFunctionNorm(snes, fnorm);
3209566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
32141ba4c6cSHeeho Park 
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
32341ba4c6cSHeeho Park   snes->norm = fnorm;
3249566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
32541ba4c6cSHeeho Park   delta      = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */
32641ba4c6cSHeeho Park   deltaM     = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */
32741ba4c6cSHeeho Park   neP->delta = delta;
3289566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3299566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
33041ba4c6cSHeeho Park 
33141ba4c6cSHeeho Park   neP->rho_satisfied = PETSC_FALSE;
33241ba4c6cSHeeho Park 
33341ba4c6cSHeeho Park   /* test convergence */
334dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
33541ba4c6cSHeeho Park   if (snes->reason) PetscFunctionReturn(0);
33641ba4c6cSHeeho Park 
33741ba4c6cSHeeho Park   for (i = 0; i < maxits; i++) {
33841ba4c6cSHeeho Park     PetscBool changed_y;
33941ba4c6cSHeeho Park     PetscBool changed_w;
34041ba4c6cSHeeho Park 
34141ba4c6cSHeeho Park     /* dogleg method */
3429566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
34341ba4c6cSHeeho Park     SNESCheckJacobianDomainerror(snes);
3449566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian));
3459566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */
34641ba4c6cSHeeho Park     SNESCheckKSPSolve(snes);                  /* this is necessary but old tr.c did not have it*/
3479566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
3489566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
34941ba4c6cSHeeho Park 
35041ba4c6cSHeeho Park     /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable)
35141ba4c6cSHeeho Park        for inner iteration and Cauchy direction calculation
35241ba4c6cSHeeho Park     */
35341ba4c6cSHeeho Park     if (bs > 1 && neP->auto_scale_multiphase) {
3549566063dSJacob Faibussowitsch       PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms));
35541ba4c6cSHeeho Park       for (j = 0; j < bs; j++) {
35641ba4c6cSHeeho Park         if (neP->auto_scale_max > 1.0) {
3579371c9d4SSatish Balay           if (inorms[j] < 1.0 / neP->auto_scale_max) { inorms[j] = 1.0 / neP->auto_scale_max; }
35841ba4c6cSHeeho Park         }
3599566063dSJacob Faibussowitsch         PetscCall(VecStrideSet(W, j, inorms[j]));
3609566063dSJacob Faibussowitsch         PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j]));
3619566063dSJacob Faibussowitsch         PetscCall(VecStrideScale(X, j, 1.0 / inorms[j]));
36241ba4c6cSHeeho Park       }
3639566063dSJacob Faibussowitsch       PetscCall(VecNorm(X, NORM_2, &xnorm));
36441ba4c6cSHeeho Park       if (i == 0) {
36541ba4c6cSHeeho Park         delta = neP->delta0 * xnorm;
36641ba4c6cSHeeho Park       } else {
36741ba4c6cSHeeho Park         delta = neP->delta * xnorm;
36841ba4c6cSHeeho Park       }
36941ba4c6cSHeeho Park       deltaM = neP->deltaM * xnorm;
3709566063dSJacob Faibussowitsch       PetscCall(MatDiagonalScale(jac, PETSC_NULL, W));
37141ba4c6cSHeeho Park     }
37241ba4c6cSHeeho Park 
37341ba4c6cSHeeho Park     /* calculating GradF of minimization function */
3749566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */
3759566063dSJacob Faibussowitsch     PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */
37641ba4c6cSHeeho Park 
37741ba4c6cSHeeho Park     inner_count        = 0;
37841ba4c6cSHeeho Park     neP->rho_satisfied = PETSC_FALSE;
37941ba4c6cSHeeho Park     while (1) {
38041ba4c6cSHeeho Park       if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */
3819566063dSJacob Faibussowitsch         PetscCall(VecCopy(YNtmp, Y));
38241ba4c6cSHeeho Park       } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */
3839566063dSJacob Faibussowitsch         PetscCall(MatMult(jac, GradF, W));
3849566063dSJacob Faibussowitsch         PetscCall(VecDotRealPart(W, W, &gTBg));     /* completes GradF^T J^T J GradF */
3859566063dSJacob Faibussowitsch         PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */
38641ba4c6cSHeeho Park         if (gTBg <= 0.0) {
38741ba4c6cSHeeho Park           auk = PETSC_MAX_REAL;
38841ba4c6cSHeeho Park         } else {
38941ba4c6cSHeeho Park           auk = PetscSqr(gfnorm) / gTBg;
39041ba4c6cSHeeho Park         }
39141ba4c6cSHeeho Park         auk = PetscMin(delta / gfnorm, auk);
3929566063dSJacob Faibussowitsch         PetscCall(VecCopy(GradF, YCtmp));           /* this could be improved */
3939566063dSJacob Faibussowitsch         PetscCall(VecScale(YCtmp, auk));            /* YCtmp, Cauchy solution*/
3949566063dSJacob Faibussowitsch         PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */
39541ba4c6cSHeeho Park         if (ycnorm >= delta) {                      /* see if the Cauchy solution meets the criteria */
3969566063dSJacob Faibussowitsch           PetscCall(VecCopy(YCtmp, Y));
3979566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm));
39841ba4c6cSHeeho Park         } else {                                  /* take ratio, tau, of Cauchy and Newton direction and step */
3999566063dSJacob Faibussowitsch           PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */
4009566063dSJacob Faibussowitsch           PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */
40141ba4c6cSHeeho Park           c0 = PetscSqr(c0);
4029566063dSJacob Faibussowitsch           PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1));
40341ba4c6cSHeeho Park           c1 = 2.0 * c1;
4049566063dSJacob Faibussowitsch           PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */
40541ba4c6cSHeeho Park           c2      = PetscSqr(c2) - PetscSqr(delta);
40641ba4c6cSHeeho Park           tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */
40741ba4c6cSHeeho Park           tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0);
40841ba4c6cSHeeho Park           tau     = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */
4099566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm));
4109566063dSJacob Faibussowitsch           PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp));
4119566063dSJacob Faibussowitsch           PetscCall(VecAXPY(W, -tau, YCtmp));
4129566063dSJacob Faibussowitsch           PetscCall(VecCopy(W, Y)); /* this could be improved */
41341ba4c6cSHeeho Park         }
41441ba4c6cSHeeho Park       } else {
41541ba4c6cSHeeho Park         /* if Cauchy is disabled, only use Newton direction */
41641ba4c6cSHeeho Park         auk = delta / ynnorm;
4179566063dSJacob Faibussowitsch         PetscCall(VecScale(YNtmp, auk));
4189566063dSJacob Faibussowitsch         PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/
41941ba4c6cSHeeho Park       }
42041ba4c6cSHeeho Park 
4219566063dSJacob Faibussowitsch       PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm  */
42241ba4c6cSHeeho Park       f0 = 0.5 * PetscSqr(fnorm);            /* minimizing function f(X) */
4239566063dSJacob Faibussowitsch       PetscCall(MatMult(jac, Y, W));
4249566063dSJacob Faibussowitsch       PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */
4259566063dSJacob Faibussowitsch       PetscCall(VecDotRealPart(GradF, Y, &gTy));
42641ba4c6cSHeeho Park       mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/
42741ba4c6cSHeeho Park 
42841ba4c6cSHeeho Park       /* scale back solution update */
42941ba4c6cSHeeho Park       if (bs > 1 && neP->auto_scale_multiphase) {
43041ba4c6cSHeeho Park         for (j = 0; j < bs; j++) {
4319566063dSJacob Faibussowitsch           PetscCall(VecStrideScale(Y, j, inorms[j]));
43241ba4c6cSHeeho Park           if (inner_count == 0) {
43341ba4c6cSHeeho Park             /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */
43441ba4c6cSHeeho Park             /* need to scale back X to match Y and provide proper update to the external code */
4359566063dSJacob Faibussowitsch             PetscCall(VecStrideScale(X, j, inorms[j]));
43641ba4c6cSHeeho Park           }
43741ba4c6cSHeeho Park         }
4389566063dSJacob Faibussowitsch         if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */
4399566063dSJacob Faibussowitsch         PetscCall(VecNorm(Y, NORM_2, &temp_ynorm));
44041ba4c6cSHeeho Park       } else {
44141ba4c6cSHeeho Park         temp_xnorm = xnorm;
44241ba4c6cSHeeho Park         temp_ynorm = ynorm;
44341ba4c6cSHeeho Park       }
44441ba4c6cSHeeho Park       inner_count++;
44541ba4c6cSHeeho Park 
44641ba4c6cSHeeho Park       /* Evaluate the solution to meet the improvement ratio criteria */
4479566063dSJacob Faibussowitsch       PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y));
4489566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(W, -1.0, Y, X));
4499566063dSJacob Faibussowitsch       PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w));
4509566063dSJacob Faibussowitsch       if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X));
4519566063dSJacob Faibussowitsch       PetscCall(VecCopy(Y, snes->vec_sol_update));
4529566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, W, G)); /*  F(X-Y) = G */
4539566063dSJacob Faibussowitsch       PetscCall(VecNorm(G, NORM_2, &gnorm));      /* gnorm <- || g || */
45441ba4c6cSHeeho Park       SNESCheckFunctionNorm(snes, gnorm);
45541ba4c6cSHeeho Park       g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */
45641ba4c6cSHeeho Park       if (f0 == mp) rho = 0.0;
45741ba4c6cSHeeho Park       else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */
45841ba4c6cSHeeho Park 
45941ba4c6cSHeeho Park       if (rho < neP->eta2) {
46041ba4c6cSHeeho Park         delta *= neP->t1; /* shrink the region */
46141ba4c6cSHeeho Park       } else if (rho > neP->eta3) {
46241ba4c6cSHeeho Park         delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */
46341ba4c6cSHeeho Park       }
46441ba4c6cSHeeho Park 
46541ba4c6cSHeeho Park       neP->delta = delta;
46641ba4c6cSHeeho Park       if (rho >= neP->eta1) {
46741ba4c6cSHeeho Park         /* unscale delta and xnorm before going to the next outer iteration */
46841ba4c6cSHeeho Park         if (bs > 1 && neP->auto_scale_multiphase) {
46941ba4c6cSHeeho Park           neP->delta = delta / xnorm;
47041ba4c6cSHeeho Park           xnorm      = temp_xnorm;
47141ba4c6cSHeeho Park           ynorm      = temp_ynorm;
47241ba4c6cSHeeho Park         }
47341ba4c6cSHeeho Park         neP->rho_satisfied = PETSC_TRUE;
47441ba4c6cSHeeho Park         break; /* the improvement ratio is satisfactory */
47541ba4c6cSHeeho Park       }
4769566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
47741ba4c6cSHeeho Park 
47841ba4c6cSHeeho Park       /* check to see if progress is hopeless */
47941ba4c6cSHeeho Park       neP->itflag = PETSC_FALSE;
48041ba4c6cSHeeho Park       /* both delta, ynorm, and xnorm are either scaled or unscaled */
4819566063dSJacob Faibussowitsch       PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP));
48241ba4c6cSHeeho Park       if (!reason) {
48341ba4c6cSHeeho Park         /* temp_xnorm, temp_ynorm is always unscaled */
48441ba4c6cSHeeho Park         /* also the inner iteration already calculated the Jacobian and solved the matrix */
48541ba4c6cSHeeho Park         /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */
4869566063dSJacob Faibussowitsch         PetscCall((*snes->ops->converged)(snes, snes->iter + 1, temp_xnorm, temp_ynorm, fnorm, &reason, snes->cnvP));
48741ba4c6cSHeeho Park       }
48841ba4c6cSHeeho Park       /* if multiphase state changes, break out inner iteration */
48941ba4c6cSHeeho Park       if (reason == SNES_BREAKOUT_INNER_ITER) {
49041ba4c6cSHeeho Park         if (bs > 1 && neP->auto_scale_multiphase) {
49141ba4c6cSHeeho Park           /* unscale delta and xnorm before going to the next outer iteration */
49241ba4c6cSHeeho Park           neP->delta = delta / xnorm;
49341ba4c6cSHeeho Park           xnorm      = temp_xnorm;
49441ba4c6cSHeeho Park           ynorm      = temp_ynorm;
49541ba4c6cSHeeho Park         }
49641ba4c6cSHeeho Park         reason = SNES_CONVERGED_ITERATING;
49741ba4c6cSHeeho Park         break;
49841ba4c6cSHeeho Park       }
49941ba4c6cSHeeho Park       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
50041ba4c6cSHeeho Park       if (reason) {
50141ba4c6cSHeeho Park         if (reason < 0) {
50241ba4c6cSHeeho Park           /* We're not progressing, so return with the current iterate */
5039566063dSJacob Faibussowitsch           PetscCall(SNESMonitor(snes, i + 1, fnorm));
50441ba4c6cSHeeho Park           breakout = PETSC_TRUE;
50541ba4c6cSHeeho Park           break;
50641ba4c6cSHeeho Park         } else if (reason > 0) {
50741ba4c6cSHeeho Park           /* We're converged, so return with the current iterate and update solution */
5089566063dSJacob Faibussowitsch           PetscCall(SNESMonitor(snes, i + 1, fnorm));
50941ba4c6cSHeeho Park           breakout = PETSC_FALSE;
51041ba4c6cSHeeho Park           break;
51141ba4c6cSHeeho Park         }
51241ba4c6cSHeeho Park       }
51341ba4c6cSHeeho Park       snes->numFailures++;
51441ba4c6cSHeeho Park     }
51541ba4c6cSHeeho Park     if (!breakout) {
51641ba4c6cSHeeho Park       /* Update function and solution vectors */
51741ba4c6cSHeeho Park       fnorm = gnorm;
5189566063dSJacob Faibussowitsch       PetscCall(VecCopy(G, F));
5199566063dSJacob Faibussowitsch       PetscCall(VecCopy(W, X));
52041ba4c6cSHeeho Park       /* Monitor convergence */
5219566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
52241ba4c6cSHeeho Park       snes->iter  = i + 1;
52341ba4c6cSHeeho Park       snes->norm  = fnorm;
52441ba4c6cSHeeho Park       snes->xnorm = xnorm;
52541ba4c6cSHeeho Park       snes->ynorm = ynorm;
5269566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5279566063dSJacob Faibussowitsch       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
5289566063dSJacob Faibussowitsch       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
52941ba4c6cSHeeho Park       /* Test for convergence, xnorm = || X || */
53041ba4c6cSHeeho Park       neP->itflag = PETSC_TRUE;
5319566063dSJacob Faibussowitsch       if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
532dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
53341ba4c6cSHeeho Park       if (reason) break;
53441ba4c6cSHeeho Park     } else break;
53541ba4c6cSHeeho Park   }
53641ba4c6cSHeeho Park 
5379566063dSJacob Faibussowitsch   /* PetscCall(PetscFree(inorms)); */
53841ba4c6cSHeeho Park   if (i == maxits) {
5399566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
54041ba4c6cSHeeho Park     if (!reason) reason = SNES_DIVERGED_MAX_IT;
54141ba4c6cSHeeho Park   }
5429566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
54341ba4c6cSHeeho Park   snes->reason = reason;
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
54541ba4c6cSHeeho Park   if (convtest != SNESTRDC_KSPConverged_Private) {
5469566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
5479566063dSJacob Faibussowitsch     PetscCall(PetscFree(ctx));
5489566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
54941ba4c6cSHeeho Park   }
55041ba4c6cSHeeho Park   PetscFunctionReturn(0);
55141ba4c6cSHeeho Park }
55241ba4c6cSHeeho Park 
55341ba4c6cSHeeho Park /*------------------------------------------------------------*/
5549371c9d4SSatish Balay static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) {
55541ba4c6cSHeeho Park   PetscFunctionBegin;
5569566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 6));
5579566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
55841ba4c6cSHeeho Park   PetscFunctionReturn(0);
55941ba4c6cSHeeho Park }
56041ba4c6cSHeeho Park 
5619371c9d4SSatish Balay PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) {
56241ba4c6cSHeeho Park   PetscFunctionBegin;
56341ba4c6cSHeeho Park   PetscFunctionReturn(0);
56441ba4c6cSHeeho Park }
56541ba4c6cSHeeho Park 
5669371c9d4SSatish Balay static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) {
56741ba4c6cSHeeho Park   PetscFunctionBegin;
5689566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONTRDC(snes));
5699566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
57041ba4c6cSHeeho Park   PetscFunctionReturn(0);
57141ba4c6cSHeeho Park }
57241ba4c6cSHeeho Park /*------------------------------------------------------------*/
57341ba4c6cSHeeho Park 
5749371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject) {
57541ba4c6cSHeeho Park   SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data;
57641ba4c6cSHeeho Park 
57741ba4c6cSHeeho Park   PetscFunctionBegin;
578d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
5799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
5809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
5819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
5829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
5839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
5849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
5859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
5869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
5879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL));
5889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL));
5899566063dSJacob 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));
590d0609cedSBarry Smith   PetscOptionsHeadEnd();
59141ba4c6cSHeeho Park   PetscFunctionReturn(0);
59241ba4c6cSHeeho Park }
59341ba4c6cSHeeho Park 
5949371c9d4SSatish Balay static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer) {
59541ba4c6cSHeeho Park   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
59641ba4c6cSHeeho Park   PetscBool        iascii;
59741ba4c6cSHeeho Park 
59841ba4c6cSHeeho Park   PetscFunctionBegin;
5999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
60041ba4c6cSHeeho Park   if (iascii) {
60163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol));
6029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
6039566063dSJacob 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));
60441ba4c6cSHeeho Park   }
60541ba4c6cSHeeho Park   PetscFunctionReturn(0);
60641ba4c6cSHeeho Park }
60741ba4c6cSHeeho Park /* ------------------------------------------------------------ */
60841ba4c6cSHeeho Park /*MC
60941ba4c6cSHeeho Park       SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
61041ba4c6cSHeeho Park 
61141ba4c6cSHeeho Park    Options Database:
61241ba4c6cSHeeho Park +   -snes_trdc_tol <tol> - trust region tolerance
61341ba4c6cSHeeho Park .   -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
61441ba4c6cSHeeho Park .   -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
61541ba4c6cSHeeho Park .   -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
61641ba4c6cSHeeho Park .   -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
61741ba4c6cSHeeho Park .   -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
61841ba4c6cSHeeho Park .   -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5)
61941ba4c6cSHeeho Park .   -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1)
62041ba4c6cSHeeho Park .   -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor
62141ba4c6cSHeeho Park .   -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm
62241ba4c6cSHeeho Park -   -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region
62341ba4c6cSHeeho Park 
62441ba4c6cSHeeho Park     Notes:
62541ba4c6cSHeeho Park     The algorithm is taken from "Linear and Nonlinear Solvers for Simulating Multiphase Flow
62641ba4c6cSHeeho Park     within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond,
62741ba4c6cSHeeho Park     Albert J. Valocchi, Tara LaForce.
62841ba4c6cSHeeho Park 
62941ba4c6cSHeeho Park    Level: intermediate
63041ba4c6cSHeeho Park 
631db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, `SNESNEWTONTRDC`
63241ba4c6cSHeeho Park 
63341ba4c6cSHeeho Park M*/
6349371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) {
63541ba4c6cSHeeho Park   SNES_NEWTONTRDC *neP;
63641ba4c6cSHeeho Park 
63741ba4c6cSHeeho Park   PetscFunctionBegin;
63841ba4c6cSHeeho Park   snes->ops->setup          = SNESSetUp_NEWTONTRDC;
63941ba4c6cSHeeho Park   snes->ops->solve          = SNESSolve_NEWTONTRDC;
64041ba4c6cSHeeho Park   snes->ops->destroy        = SNESDestroy_NEWTONTRDC;
64141ba4c6cSHeeho Park   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC;
64241ba4c6cSHeeho Park   snes->ops->view           = SNESView_NEWTONTRDC;
64341ba4c6cSHeeho Park   snes->ops->reset          = SNESReset_NEWTONTRDC;
64441ba4c6cSHeeho Park 
64541ba4c6cSHeeho Park   snes->usesksp = PETSC_TRUE;
64641ba4c6cSHeeho Park   snes->usesnpc = PETSC_FALSE;
64741ba4c6cSHeeho Park 
64841ba4c6cSHeeho Park   snes->alwayscomputesfinalresidual = PETSC_TRUE;
64941ba4c6cSHeeho Park 
6509566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &neP));
65141ba4c6cSHeeho Park   snes->data                 = (void *)neP;
65241ba4c6cSHeeho Park   neP->delta                 = 0.0;
65341ba4c6cSHeeho Park   neP->delta0                = 0.1;
65441ba4c6cSHeeho Park   neP->eta1                  = 0.001;
65541ba4c6cSHeeho Park   neP->eta2                  = 0.25;
65641ba4c6cSHeeho Park   neP->eta3                  = 0.75;
65741ba4c6cSHeeho Park   neP->t1                    = 0.25;
65841ba4c6cSHeeho Park   neP->t2                    = 2.0;
65941ba4c6cSHeeho Park   neP->deltaM                = 0.5;
66041ba4c6cSHeeho Park   neP->sigma                 = 0.0001;
66141ba4c6cSHeeho Park   neP->itflag                = PETSC_FALSE;
66241ba4c6cSHeeho Park   neP->rnorm0                = 0.0;
66341ba4c6cSHeeho Park   neP->ttol                  = 0.0;
66441ba4c6cSHeeho Park   neP->use_cauchy            = PETSC_TRUE;
66541ba4c6cSHeeho Park   neP->auto_scale_multiphase = PETSC_FALSE;
66641ba4c6cSHeeho Park   neP->auto_scale_max        = -1.0;
66741ba4c6cSHeeho Park   neP->rho_satisfied         = PETSC_FALSE;
66841ba4c6cSHeeho Park   snes->deltatol             = 1.e-12;
66941ba4c6cSHeeho Park 
67041ba4c6cSHeeho Park   /* for multiphase (multivariable) scaling */
67141ba4c6cSHeeho Park   /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13
67241ba4c6cSHeeho Park      on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now.
6739566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(snes->work[0],&neP->bs));
6749566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(neP->bs,&neP->inorms));
67541ba4c6cSHeeho Park   */
67641ba4c6cSHeeho Park 
67741ba4c6cSHeeho Park   PetscFunctionReturn(0);
67841ba4c6cSHeeho Park }
679