xref: /petsc/src/snes/impls/ntrdc/ntrdc.c (revision db2b9530cea7203cdddcdfcaf19eada576ef3459)
1 
2 #include <../src/snes/impls/ntrdc/ntrdcimpl.h>                /*I   "petscsnes.h"   I*/
3 
4 typedef struct {
5   SNES           snes;
6   /*  Information on the regular SNES convergence test; which may have been user provided
7       Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho
8       Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private
9  */
10 
11   PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
12   PetscErrorCode (*convdestroy)(void*);
13   void           *convctx;
14 } SNES_TRDC_KSPConverged_Ctx;
15 
16 static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx)
17 {
18   SNES_TRDC_KSPConverged_Ctx  *ctx = (SNES_TRDC_KSPConverged_Ctx*)cctx;
19   SNES                        snes = ctx->snes;
20   SNES_NEWTONTRDC             *neP = (SNES_NEWTONTRDC*)snes->data;
21   Vec                         x;
22   PetscReal                   nrm;
23   PetscErrorCode              ierr;
24 
25   PetscFunctionBegin;
26   ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr);
27   if (*reason) {
28     ierr = PetscInfo(snes,"Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr);
29   }
30   /* Determine norm of solution */
31   ierr = KSPBuildSolution(ksp,NULL,&x);CHKERRQ(ierr);
32   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
33   if (nrm >= neP->delta) {
34     ierr    = PetscInfo(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr);
35     *reason = KSP_CONVERGED_STEP_LENGTH;
36   }
37   PetscFunctionReturn(0);
38 }
39 
40 static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx)
41 {
42   SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx*)cctx;
43   PetscErrorCode             ierr;
44 
45   PetscFunctionBegin;
46   ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr);
47   ierr = PetscFree(ctx);CHKERRQ(ierr);
48 
49   PetscFunctionReturn(0);
50 }
51 
52 /* ---------------------------------------------------------------- */
53 /*
54    SNESTRDC_Converged_Private -test convergence JUST for
55    the trust region tolerance.
56 
57 */
58 static PetscErrorCode SNESTRDC_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
59 {
60   SNES_NEWTONTRDC  *neP = (SNES_NEWTONTRDC*)snes->data;
61   PetscErrorCode   ierr;
62 
63   PetscFunctionBegin;
64   *reason = SNES_CONVERGED_ITERATING;
65   if (neP->delta < xnorm * snes->deltatol) {
66     ierr    = PetscInfo(snes,"Diverged due to too small a trust region %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr);
67     *reason = SNES_DIVERGED_TR_DELTA;
68   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
69     ierr    = PetscInfo(snes,"Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n",snes->max_funcs);CHKERRQ(ierr);
70     *reason = SNES_DIVERGED_FUNCTION_COUNT;
71   }
72   PetscFunctionReturn(0);
73 }
74 
75 /*@
76   SNESNewtonTRDCGetRhoFlag - Get whether the solution update is within the trust-region.
77 
78   Input Parameters:
79 . snes - the nonlinear solver object
80 
81   Output Parameters:
82 . rho_flag: PETSC_TRUE if the solution update is in the trust-region; otherwise, PETSC_FALSE
83 
84   Level: developer
85 
86 @*/
87 PetscErrorCode  SNESNewtonTRDCGetRhoFlag(SNES snes,PetscBool *rho_flag)
88 {
89   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
90 
91   PetscFunctionBegin;
92   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
93   PetscValidBoolPointer(rho_flag,2);
94   *rho_flag = tr->rho_satisfied;
95   PetscFunctionReturn(0);
96 }
97 
98 /*@C
99    SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined.
100        Allows the user a chance to change or override the trust region decision.
101 
102    Logically Collective on snes
103 
104    Input Parameters:
105 +  snes - the nonlinear solver object
106 .  func - [optional] function evaluation routine, see SNESNewtonTRDCPreCheck()  for the calling sequence
107 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
108 
109    Level: intermediate
110 
111    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver.
112 
113 .seealso: SNESNewtonTRDCPreCheck(), SNESNewtonTRDCGetPreCheck(), SNESNewtonTRDCSetPostCheck(), SNESNewtonTRDCGetPostCheck()
114 @*/
115 PetscErrorCode  SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx)
116 {
117   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
121   if (func) tr->precheck    = func;
122   if (ctx)  tr->precheckctx = ctx;
123   PetscFunctionReturn(0);
124 }
125 
126 /*@C
127    SNESNewtonTRDCGetPreCheck - Gets the pre-check function
128 
129    Not collective
130 
131    Input Parameter:
132 .  snes - the nonlinear solver context
133 
134    Output Parameters:
135 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPreCheck()
136 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
137 
138    Level: intermediate
139 
140 .seealso: SNESNewtonTRDCSetPreCheck(), SNESNewtonTRDCPreCheck()
141 @*/
142 PetscErrorCode  SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx)
143 {
144   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
145 
146   PetscFunctionBegin;
147   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
148   if (func) *func = tr->precheck;
149   if (ctx)  *ctx  = tr->precheckctx;
150   PetscFunctionReturn(0);
151 }
152 
153 /*@C
154    SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
155        function evaluation. Allows the user a chance to change or override the decision of the line search routine
156 
157    Logically Collective on snes
158 
159    Input Parameters:
160 +  snes - the nonlinear solver object
161 .  func - [optional] function evaluation routine, see SNESNewtonTRDCPostCheck()  for the calling sequence
162 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
163 
164    Level: intermediate
165 
166    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver while the function set in
167    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
168 
169 .seealso: SNESNewtonTRDCPostCheck(), SNESNewtonTRDCGetPostCheck()
170 @*/
171 PetscErrorCode  SNESNewtonTRDCSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
172 {
173   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
174 
175   PetscFunctionBegin;
176   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
177   if (func) tr->postcheck    = func;
178   if (ctx)  tr->postcheckctx = ctx;
179   PetscFunctionReturn(0);
180 }
181 
182 /*@C
183    SNESNewtonTRDCGetPostCheck - Gets the post-check function
184 
185    Not collective
186 
187    Input Parameter:
188 .  snes - the nonlinear solver context
189 
190    Output Parameters:
191 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck()
192 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
193 
194    Level: intermediate
195 
196 .seealso: SNESNewtonTRDCSetPostCheck(), SNESNewtonTRDCPostCheck()
197 @*/
198 PetscErrorCode  SNESNewtonTRDCGetPostCheck(SNES snes,PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
199 {
200   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
201 
202   PetscFunctionBegin;
203   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
204   if (func) *func = tr->postcheck;
205   if (ctx)  *ctx  = tr->postcheckctx;
206   PetscFunctionReturn(0);
207 }
208 
209 /*@C
210    SNESNewtonTRDCPreCheck - Called before the step has been determined in SNESNEWTONTRDC
211 
212    Logically Collective on snes
213 
214    Input Parameters:
215 +  snes - the solver
216 .  X - The last solution
217 -  Y - The step direction
218 
219    Output Parameters:
220 .  changed_Y - Indicator that the step direction Y has been changed.
221 
222    Level: developer
223 
224 .seealso: SNESNewtonTRDCSetPreCheck(), SNESNewtonTRDCGetPreCheck()
225 @*/
226 static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y)
227 {
228   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
229   PetscErrorCode   ierr;
230 
231   PetscFunctionBegin;
232   *changed_Y = PETSC_FALSE;
233   if (tr->precheck) {
234     ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr);
235     PetscValidLogicalCollectiveBool(snes,*changed_Y,4);
236   }
237   PetscFunctionReturn(0);
238 }
239 
240 /*@C
241    SNESNewtonTRDCPostCheck - Called after the step has been determined in SNESNEWTONTRDC but before the function evaluation
242 
243    Logically Collective on snes
244 
245    Input Parameters:
246 +  snes - the solver
247 .  X - The last solution
248 .  Y - The full step direction
249 -  W - The updated solution, W = X - Y
250 
251    Output Parameters:
252 +  changed_Y - indicator if step has been changed
253 -  changed_W - Indicator if the new candidate solution W has been changed.
254 
255    Notes:
256      If Y is changed then W is recomputed as X - Y
257 
258    Level: developer
259 
260 .seealso: SNESNewtonTRDCSetPostCheck(), SNESNewtonTRDCGetPostCheck()
261 @*/
262 static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
263 {
264   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
265   PetscErrorCode   ierr;
266 
267   PetscFunctionBegin;
268   *changed_Y = PETSC_FALSE;
269   *changed_W = PETSC_FALSE;
270   if (tr->postcheck) {
271     ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr);
272     PetscValidLogicalCollectiveBool(snes,*changed_Y,5);
273     PetscValidLogicalCollectiveBool(snes,*changed_W,6);
274   }
275   PetscFunctionReturn(0);
276 }
277 
278 /*
279    SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
280    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
281    nonlinear equations
282 
283 */
284 static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes)
285 {
286   SNES_NEWTONTRDC            *neP = (SNES_NEWTONTRDC*)snes->data;
287   Vec                        X,F,Y,G,W,GradF,YNtmp;
288   Vec                        YCtmp;
289   Mat                        jac;
290   PetscErrorCode             ierr;
291   PetscInt                   maxits,i,j,lits,inner_count,bs;
292   PetscReal                  rho,fnorm,gnorm,xnorm=0,delta,ynorm,temp_xnorm,temp_ynorm;  /* TRDC inner iteration */
293   PetscReal                  inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */
294   PetscReal                  deltaM,ynnorm,f0,mp,gTy,g,yTHy;  /* rho calculation */
295   PetscReal                  auk,gfnorm,ycnorm,c0,c1,c2,tau,tau_pos,tau_neg,gTBg;  /* Cauchy Point */
296   KSP                        ksp;
297   SNESConvergedReason        reason = SNES_CONVERGED_ITERATING;
298   PetscBool                  breakout = PETSC_FALSE;
299   SNES_TRDC_KSPConverged_Ctx *ctx;
300   PetscErrorCode             (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),(*convdestroy)(void*);
301   void                       *convctx;
302 
303   PetscFunctionBegin;
304   maxits = snes->max_its;               /* maximum number of iterations */
305   X      = snes->vec_sol;               /* solution vector */
306   F      = snes->vec_func;              /* residual vector */
307   Y      = snes->work[0];               /* update vector */
308   G      = snes->work[1];               /* updated residual */
309   W      = snes->work[2];               /* temporary vector */
310   GradF  = snes->work[3];               /* grad f = J^T F */
311   YNtmp  = snes->work[4];               /* Newton solution */
312   YCtmp  = snes->work[5];               /* Cauchy solution */
313 
314   PetscCheckFalse(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);
315 
316   ierr = VecGetBlockSize(YNtmp,&bs);CHKERRQ(ierr);
317 
318   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
319   snes->iter = 0;
320   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
321 
322   /* Set the linear stopping criteria to use the More' trick. From tr.c */
323   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
324   ierr = KSPGetConvergenceTest(ksp,&convtest,&convctx,&convdestroy);CHKERRQ(ierr);
325   if (convtest != SNESTRDC_KSPConverged_Private) {
326     ierr                  = PetscNew(&ctx);CHKERRQ(ierr);
327     ctx->snes             = snes;
328     ierr                  = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
329     ierr                  = KSPSetConvergenceTest(ksp,SNESTRDC_KSPConverged_Private,ctx,SNESTRDC_KSPConverged_Destroy);CHKERRQ(ierr);
330     ierr                  = PetscInfo(snes,"Using Krylov convergence test SNESTRDC_KSPConverged_Private\n");CHKERRQ(ierr);
331   }
332 
333   if (!snes->vec_func_init_set) {
334     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
335   } else snes->vec_func_init_set = PETSC_FALSE;
336 
337   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
338   SNESCheckFunctionNorm(snes,fnorm);
339   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);             /* xnorm <- || X || */
340 
341   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
342   snes->norm = fnorm;
343   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
344   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;  /* initial trust region size scaled by xnorm */
345   deltaM     = xnorm ? neP->deltaM*xnorm : neP->deltaM;  /* maximum trust region size scaled by xnorm */
346   neP->delta = delta;
347   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
348   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
349 
350   neP->rho_satisfied = PETSC_FALSE;
351 
352   /* test convergence */
353   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
354   if (snes->reason) PetscFunctionReturn(0);
355 
356   for (i=0; i<maxits; i++) {
357     PetscBool changed_y;
358     PetscBool changed_w;
359 
360     /* dogleg method */
361     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
362     SNESCheckJacobianDomainerror(snes);
363     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian);CHKERRQ(ierr);
364     ierr = KSPSolve(snes->ksp,F,YNtmp);CHKERRQ(ierr);   /* Quasi Newton Solution */
365     SNESCheckKSPSolve(snes);  /* this is necessary but old tr.c did not have it*/
366     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
367     ierr = SNESGetJacobian(snes,&jac,NULL,NULL,NULL);CHKERRQ(ierr);
368 
369     /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable)
370        for inner iteration and Cauchy direction calculation
371     */
372     if (bs > 1 && neP->auto_scale_multiphase) {
373       ierr = VecStrideNormAll(YNtmp,NORM_INFINITY,inorms);CHKERRQ(ierr);
374       for (j=0; j<bs; j++) {
375         if (neP->auto_scale_max > 1.0) {
376           if (inorms[j] < 1.0/neP->auto_scale_max) {
377             inorms[j] = 1.0/neP->auto_scale_max;
378           }
379         }
380         ierr = VecStrideSet(W,j,inorms[j]);CHKERRQ(ierr);
381         ierr = VecStrideScale(YNtmp,j,1.0/inorms[j]);
382         ierr = VecStrideScale(X,j,1.0/inorms[j]);
383       }
384       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);
385       if (i == 0) {
386         delta = neP->delta0*xnorm;
387       } else {
388         delta = neP->delta*xnorm;
389       }
390       deltaM = neP->deltaM*xnorm;
391       ierr = MatDiagonalScale(jac,PETSC_NULL,W);CHKERRQ(ierr);
392     }
393 
394     /* calculating GradF of minimization function */
395     ierr = MatMultTranspose(jac,F,GradF);CHKERRQ(ierr);  /* grad f = J^T F */
396     ierr = VecNorm(YNtmp,NORM_2,&ynnorm);CHKERRQ(ierr);  /* ynnorm <- || Y_newton || */
397 
398     inner_count = 0;
399     neP->rho_satisfied = PETSC_FALSE;
400     while (1) {
401       if (ynnorm <= delta) {  /* see if the Newton solution is within the trust region */
402         ierr = VecCopy(YNtmp,Y);CHKERRQ(ierr);
403       } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */
404         ierr = MatMult(jac,GradF,W);CHKERRQ(ierr);
405         ierr = VecDotRealPart(W,W,&gTBg);CHKERRQ(ierr);  /* completes GradF^T J^T J GradF */
406         ierr = VecNorm(GradF,NORM_2,&gfnorm);CHKERRQ(ierr);  /* grad f norm <- || grad f || */
407         if (gTBg <= 0.0) {
408           auk = PETSC_MAX_REAL;
409         } else {
410           auk = PetscSqr(gfnorm)/gTBg;
411         }
412         auk  = PetscMin(delta/gfnorm,auk);
413         ierr = VecCopy(GradF,YCtmp);CHKERRQ(ierr);  /* this could be improved */
414         ierr = VecScale(YCtmp,auk);CHKERRQ(ierr);  /* YCtmp, Cauchy solution*/
415         ierr = VecNorm(YCtmp,NORM_2,&ycnorm);CHKERRQ(ierr);  /* ycnorm <- || Y_cauchy || */
416         if (ycnorm >= delta) {  /* see if the Cauchy solution meets the criteria */
417             ierr = VecCopy(YCtmp,Y);CHKERRQ(ierr);
418             ierr = PetscInfo(snes,"DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n",(double)delta,(double)ynnorm,(double)ycnorm);CHKERRQ(ierr);
419         } else {  /* take ratio, tau, of Cauchy and Newton direction and step */
420           ierr    = VecAXPY(YNtmp,-1.0,YCtmp);CHKERRQ(ierr);  /* YCtmp = A, YNtmp = B */
421           ierr    = VecNorm(YNtmp,NORM_2,&c0);CHKERRQ(ierr);  /* this could be improved */
422           c0      = PetscSqr(c0);
423           ierr    = VecDotRealPart(YCtmp,YNtmp,&c1);CHKERRQ(ierr);
424           c1      = 2.0*c1;
425           ierr    = VecNorm(YCtmp,NORM_2,&c2);CHKERRQ(ierr);  /* this could be improved */
426           c2      = PetscSqr(c2) - PetscSqr(delta);
427           tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4.*c0*c2))/(2.*c0); /* quadratic formula */
428           tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4.*c0*c2))/(2.*c0);
429           tau     = PetscMax(tau_pos, tau_neg);  /* can tau_neg > tau_pos? I don't think so, but just in case. */
430           ierr    = PetscInfo(snes,"DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n",(double)tau,(double)ynnorm,(double)ycnorm);CHKERRQ(ierr);
431           ierr    = VecWAXPY(W,tau,YNtmp,YCtmp);CHKERRQ(ierr);
432           ierr    = VecAXPY(W,-tau,YCtmp);CHKERRQ(ierr);
433           ierr    = VecCopy(W, Y);CHKERRQ(ierr); /* this could be improved */
434         }
435       } else {
436         /* if Cauchy is disabled, only use Newton direction */
437         auk = delta/ynnorm;
438         ierr = VecScale(YNtmp,auk);CHKERRQ(ierr);
439         ierr = VecCopy(YNtmp,Y);CHKERRQ(ierr); /* this could be improved (many VecCopy, VecNorm)*/
440       }
441 
442       ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);  /* compute the final ynorm  */
443       f0 = 0.5*PetscSqr(fnorm);  /* minimizing function f(X) */
444       ierr = MatMult(jac,Y,W);CHKERRQ(ierr);
445       ierr = VecDotRealPart(W,W,&yTHy);CHKERRQ(ierr);  /* completes GradY^T J^T J GradY */
446       ierr = VecDotRealPart(GradF,Y,&gTy);CHKERRQ(ierr);
447       mp = f0 - gTy + 0.5*yTHy;  /* quadratic model to satisfy, -gTy because our update is X-Y*/
448 
449       /* scale back solution update */
450       if (bs > 1 && neP->auto_scale_multiphase) {
451         for (j=0; j<bs; j++) {
452           ierr = VecStrideScale(Y,j,inorms[j]);
453           if (inner_count == 0) {
454             /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */
455             /* need to scale back X to match Y and provide proper update to the external code */
456             ierr = VecStrideScale(X,j,inorms[j]);
457           }
458         }
459         if (inner_count == 0) {ierr = VecNorm(X,NORM_2,&temp_xnorm);CHKERRQ(ierr);}  /* only in the first iteration */
460         ierr = VecNorm(Y,NORM_2,&temp_ynorm);CHKERRQ(ierr);
461       } else {
462         temp_xnorm = xnorm;
463         temp_ynorm = ynorm;
464       }
465       inner_count++;
466 
467       /* Evaluate the solution to meet the improvement ratio criteria */
468       ierr = SNESNewtonTRDCPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr);
469       ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);
470       ierr = SNESNewtonTRDCPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
471       if (changed_y) {ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);}
472       ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
473       ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /*  F(X-Y) = G */
474       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
475       SNESCheckFunctionNorm(snes,gnorm);
476       g = 0.5*PetscSqr(gnorm); /* minimizing function g(W) */
477       if (f0 == mp) rho = 0.0;
478       else rho = (f0 - g)/(f0 - mp);  /* actual improvement over predicted improvement */
479 
480       if (rho < neP->eta2) {
481         delta *= neP->t1;  /* shrink the region */
482       } else if (rho > neP->eta3) {
483         delta = PetscMin(neP->t2*delta,deltaM); /* expand the region, but not greater than deltaM */
484       }
485 
486       neP->delta = delta;
487       if (rho >= neP->eta1) {
488         /* unscale delta and xnorm before going to the next outer iteration */
489         if (bs > 1 && neP->auto_scale_multiphase) {
490           neP->delta = delta/xnorm;
491           xnorm      = temp_xnorm;
492           ynorm      = temp_ynorm;
493         }
494         neP->rho_satisfied = PETSC_TRUE;
495         break;  /* the improvement ratio is satisfactory */
496       }
497       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
498 
499       /* check to see if progress is hopeless */
500       neP->itflag = PETSC_FALSE;
501       /* both delta, ynorm, and xnorm are either scaled or unscaled */
502       ierr        = SNESTRDC_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
503       if (!reason) {
504          /* temp_xnorm, temp_ynorm is always unscaled */
505          /* also the inner iteration already calculated the Jacobian and solved the matrix */
506          /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */
507          ierr = (*snes->ops->converged)(snes,snes->iter+1,temp_xnorm,temp_ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
508       }
509       /* if multiphase state changes, break out inner iteration */
510       if (reason == SNES_BREAKOUT_INNER_ITER) {
511         if (bs > 1 && neP->auto_scale_multiphase) {
512           /* unscale delta and xnorm before going to the next outer iteration */
513           neP->delta = delta/xnorm;
514           xnorm      = temp_xnorm;
515           ynorm      = temp_ynorm;
516         }
517         reason = SNES_CONVERGED_ITERATING;
518         break;
519       }
520       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
521       if (reason) {
522         if (reason < 0) {
523             /* We're not progressing, so return with the current iterate */
524             ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
525             breakout = PETSC_TRUE;
526             break;
527         } else if (reason > 0) {
528             /* We're converged, so return with the current iterate and update solution */
529             ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
530             breakout = PETSC_FALSE;
531             break;
532         }
533       }
534       snes->numFailures++;
535     }
536     if (!breakout) {
537       /* Update function and solution vectors */
538       fnorm       = gnorm;
539       ierr        = VecCopy(G,F);CHKERRQ(ierr);
540       ierr        = VecCopy(W,X);CHKERRQ(ierr);
541       /* Monitor convergence */
542       ierr        = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
543       snes->iter  = i+1;
544       snes->norm  = fnorm;
545       snes->xnorm = xnorm;
546       snes->ynorm = ynorm;
547       ierr        = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
548       ierr        = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
549       ierr        = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
550       /* Test for convergence, xnorm = || X || */
551       neP->itflag = PETSC_TRUE;
552       if (snes->ops->converged != SNESConvergedSkip) {ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);}
553       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
554       if (reason) break;
555     } else break;
556   }
557 
558   /* ierr         = PetscFree(inorms);CHKERRQ(ierr); */
559   if (i == maxits) {
560     ierr = PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",maxits);CHKERRQ(ierr);
561     if (!reason) reason = SNES_DIVERGED_MAX_IT;
562   }
563   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
564   snes->reason = reason;
565   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
566   if (convtest != SNESTRDC_KSPConverged_Private) {
567     ierr       = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
568     ierr       = PetscFree(ctx);CHKERRQ(ierr);
569     ierr       = KSPSetConvergenceTest(ksp,convtest,convctx,convdestroy);CHKERRQ(ierr);
570   }
571   PetscFunctionReturn(0);
572 }
573 
574 /*------------------------------------------------------------*/
575 static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes)
576 {
577   PetscErrorCode ierr;
578 
579   PetscFunctionBegin;
580   ierr = SNESSetWorkVecs(snes,6);CHKERRQ(ierr);
581   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
582   PetscFunctionReturn(0);
583 }
584 
585 PetscErrorCode SNESReset_NEWTONTRDC(SNES snes)
586 {
587 
588   PetscFunctionBegin;
589   PetscFunctionReturn(0);
590 }
591 
592 static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes)
593 {
594   PetscErrorCode ierr;
595 
596   PetscFunctionBegin;
597   ierr = SNESReset_NEWTONTRDC(snes);CHKERRQ(ierr);
598   ierr = PetscFree(snes->data);CHKERRQ(ierr);
599   PetscFunctionReturn(0);
600 }
601 /*------------------------------------------------------------*/
602 
603 static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(PetscOptionItems *PetscOptionsObject,SNES snes)
604 {
605   SNES_NEWTONTRDC  *ctx = (SNES_NEWTONTRDC*)snes->data;
606   PetscErrorCode   ierr;
607 
608   PetscFunctionBegin;
609   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
610   ierr = PetscOptionsReal("-snes_trdc_tol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
611   ierr = PetscOptionsReal("-snes_trdc_eta1","eta1","None",ctx->eta1,&ctx->eta1,NULL);CHKERRQ(ierr);
612   ierr = PetscOptionsReal("-snes_trdc_eta2","eta2","None",ctx->eta2,&ctx->eta2,NULL);CHKERRQ(ierr);
613   ierr = PetscOptionsReal("-snes_trdc_eta3","eta3","None",ctx->eta3,&ctx->eta3,NULL);CHKERRQ(ierr);
614   ierr = PetscOptionsReal("-snes_trdc_t1","t1","None",ctx->t1,&ctx->t1,NULL);CHKERRQ(ierr);
615   ierr = PetscOptionsReal("-snes_trdc_t2","t2","None",ctx->t2,&ctx->t2,NULL);CHKERRQ(ierr);
616   ierr = PetscOptionsReal("-snes_trdc_deltaM","deltaM","None",ctx->deltaM,&ctx->deltaM,NULL);CHKERRQ(ierr);
617   ierr = PetscOptionsReal("-snes_trdc_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
618   ierr = PetscOptionsReal("-snes_trdc_auto_scale_max","auto_scale_max","None",ctx->auto_scale_max,&ctx->auto_scale_max,NULL);CHKERRQ(ierr);
619   ierr = PetscOptionsBool("-snes_trdc_use_cauchy","use_cauchy","use Cauchy step and direction",ctx->use_cauchy,&ctx->use_cauchy,NULL);CHKERRQ(ierr);
620   ierr = PetscOptionsBool("-snes_trdc_auto_scale_multiphase","auto_scale_multiphase","Auto scaling for proper cauchy direction",ctx->auto_scale_multiphase,&ctx->auto_scale_multiphase,NULL);CHKERRQ(ierr);
621   ierr = PetscOptionsTail();CHKERRQ(ierr);
622   PetscFunctionReturn(0);
623 }
624 
625 static PetscErrorCode SNESView_NEWTONTRDC(SNES snes,PetscViewer viewer)
626 {
627   SNES_NEWTONTRDC  *tr = (SNES_NEWTONTRDC*)snes->data;
628   PetscErrorCode   ierr;
629   PetscBool        iascii;
630 
631   PetscFunctionBegin;
632   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
633   if (iascii) {
634     ierr = PetscViewerASCIIPrintf(viewer,"  Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr);
635     ierr = PetscViewerASCIIPrintf(viewer,"  eta1=%g, eta2=%g, eta3=%g\n",(double)tr->eta1,(double)tr->eta2,(double)tr->eta3);CHKERRQ(ierr);
636     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, t1=%g, t2=%g, deltaM=%g\n",(double)tr->delta0,(double)tr->t1,(double)tr->t2,(double)tr->deltaM);CHKERRQ(ierr);
637   }
638   PetscFunctionReturn(0);
639 }
640 /* ------------------------------------------------------------ */
641 /*MC
642       SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
643 
644    Options Database:
645 +   -snes_trdc_tol <tol> - trust region tolerance
646 .   -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
647 .   -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
648 .   -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
649 .   -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
650 .   -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
651 .   -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5)
652 .   -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1)
653 .   -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor
654 .   -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm
655 -   -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region
656 
657     Notes:
658     The algorithm is taken from "Linear and Nonlinear Solvers for Simulating Multiphase Flow
659     within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond,
660     Albert J. Valocchi, Tara LaForce.
661 
662    Level: intermediate
663 
664 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance(), SNESNEWTONTRDC
665 
666 M*/
667 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes)
668 {
669   SNES_NEWTONTRDC  *neP;
670   PetscErrorCode   ierr;
671 
672   PetscFunctionBegin;
673   snes->ops->setup          = SNESSetUp_NEWTONTRDC;
674   snes->ops->solve          = SNESSolve_NEWTONTRDC;
675   snes->ops->destroy        = SNESDestroy_NEWTONTRDC;
676   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC;
677   snes->ops->view           = SNESView_NEWTONTRDC;
678   snes->ops->reset          = SNESReset_NEWTONTRDC;
679 
680   snes->usesksp = PETSC_TRUE;
681   snes->usesnpc = PETSC_FALSE;
682 
683   snes->alwayscomputesfinalresidual = PETSC_TRUE;
684 
685   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
686   snes->data  = (void*)neP;
687   neP->delta  = 0.0;
688   neP->delta0 = 0.1;
689   neP->eta1   = 0.001;
690   neP->eta2   = 0.25;
691   neP->eta3   = 0.75;
692   neP->t1     = 0.25;
693   neP->t2     = 2.0;
694   neP->deltaM = 0.5;
695   neP->sigma  = 0.0001;
696   neP->itflag = PETSC_FALSE;
697   neP->rnorm0 = 0.0;
698   neP->ttol   = 0.0;
699   neP->use_cauchy            = PETSC_TRUE;
700   neP->auto_scale_multiphase = PETSC_FALSE;
701   neP->auto_scale_max        = -1.0;
702   neP->rho_satisfied         = PETSC_FALSE;
703   snes->deltatol             = 1.e-12;
704 
705   /* for multiphase (multivariable) scaling */
706   /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13
707      on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now.
708   ierr = VecGetBlockSize(snes->work[0],&neP->bs);CHKERRQ(ierr);
709   ierr = PetscCalloc1(neP->bs,&neP->inorms);CHKERRQ(ierr);
710   */
711 
712   PetscFunctionReturn(0);
713 }
714