xref: /petsc/src/snes/impls/tr/tr.c (revision 7cb011f58c32087171e6c65e6b5293ee2291f1f1)
14800dd8cSBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h>                /*I   "petscsnes.h"   I*/
34800dd8cSBarry Smith 
4971273eeSBarry Smith typedef struct {
5971273eeSBarry Smith   SNES           snes;
6df8705c3SBarry Smith   /*  Information on the regular SNES convergence test; which may have been user provided */
7df8705c3SBarry Smith   PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
8df8705c3SBarry Smith   PetscErrorCode (*convdestroy)(void*);
9df8705c3SBarry Smith   void           *convctx;
10971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx;
11971273eeSBarry Smith 
12470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx)
13df60cc22SBarry Smith {
14971273eeSBarry Smith   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
15971273eeSBarry Smith   SNES                     snes = ctx->snes;
1604d7464bSBarry Smith   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
17df60cc22SBarry Smith   Vec                      x;
18064f8208SBarry Smith   PetscReal                nrm;
19dfbe8321SBarry Smith   PetscErrorCode           ierr;
20df60cc22SBarry Smith 
213a40ed3dSBarry Smith   PetscFunctionBegin;
22df8705c3SBarry Smith   ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr);
23329f5518SBarry Smith   if (*reason) {
24df8705c3SBarry Smith     ierr = PetscInfo2(snes,"Default or user provided convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr);
259b04593eSLois Curfman McInnes   }
26a935fc98SLois Curfman McInnes   /* Determine norm of solution */
2778b31e54SBarry Smith   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
28064f8208SBarry Smith   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
29064f8208SBarry Smith   if (nrm >= neP->delta) {
3057622a8eSBarry Smith     ierr    = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr);
31329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
32df60cc22SBarry Smith   }
333a40ed3dSBarry Smith   PetscFunctionReturn(0);
34df60cc22SBarry Smith }
3582bf6240SBarry Smith 
36470880abSPatrick Sanan static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx)
37971273eeSBarry Smith {
38971273eeSBarry Smith   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
39971273eeSBarry Smith   PetscErrorCode           ierr;
40971273eeSBarry Smith 
41971273eeSBarry Smith   PetscFunctionBegin;
42df8705c3SBarry Smith   ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr);
43971273eeSBarry Smith   ierr = PetscFree(ctx);CHKERRQ(ierr);
44971273eeSBarry Smith   PetscFunctionReturn(0);
45971273eeSBarry Smith }
46971273eeSBarry Smith 
4785385478SLisandro Dalcin /* ---------------------------------------------------------------- */
4885385478SLisandro Dalcin /*
49470880abSPatrick Sanan    SNESTR_Converged_Private -test convergence JUST for
5085385478SLisandro Dalcin    the trust region tolerance.
5185385478SLisandro Dalcin 
5285385478SLisandro Dalcin */
53470880abSPatrick Sanan static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
5485385478SLisandro Dalcin {
5504d7464bSBarry Smith   SNES_NEWTONTR  *neP = (SNES_NEWTONTR*)snes->data;
5685385478SLisandro Dalcin   PetscErrorCode ierr;
5785385478SLisandro Dalcin 
5885385478SLisandro Dalcin   PetscFunctionBegin;
5985385478SLisandro Dalcin   *reason = SNES_CONVERGED_ITERATING;
6085385478SLisandro Dalcin   if (neP->delta < xnorm * snes->deltatol) {
6157622a8eSBarry Smith     ierr    = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr);
6285385478SLisandro Dalcin     *reason = SNES_CONVERGED_TR_DELTA;
63e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
6485385478SLisandro Dalcin     ierr    = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr);
6585385478SLisandro Dalcin     *reason = SNES_DIVERGED_FUNCTION_COUNT;
6685385478SLisandro Dalcin   }
6785385478SLisandro Dalcin   PetscFunctionReturn(0);
6885385478SLisandro Dalcin }
6985385478SLisandro Dalcin 
70*7cb011f5SBarry Smith /*@C
71*7cb011f5SBarry Smith    SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
72*7cb011f5SBarry Smith        function evaluation. Allows the user a chance to change or override the decision of the line search routine
73*7cb011f5SBarry Smith 
74*7cb011f5SBarry Smith    Logically Collective on snes
75*7cb011f5SBarry Smith 
76*7cb011f5SBarry Smith    Input Parameters:
77*7cb011f5SBarry Smith +  snes - the nonlinear solver object
78*7cb011f5SBarry Smith .  func - [optional] function evaluation routine, see SNESNewtonTRPostCheck()  for the calling sequence
79*7cb011f5SBarry Smith -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
80*7cb011f5SBarry Smith 
81*7cb011f5SBarry Smith    Level: intermediate
82*7cb011f5SBarry Smith 
83*7cb011f5SBarry Smith    Note: This function is called BEFORE the function evalaulation within the SNESNEWTONTR solver while the function set in
84*7cb011f5SBarry Smith    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
85*7cb011f5SBarry Smith 
86*7cb011f5SBarry Smith .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck()
87*7cb011f5SBarry Smith @*/
88*7cb011f5SBarry Smith PetscErrorCode  SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,void*),void *ctx)
89*7cb011f5SBarry Smith {
90*7cb011f5SBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
91*7cb011f5SBarry Smith 
92*7cb011f5SBarry Smith   PetscFunctionBegin;
93*7cb011f5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
94*7cb011f5SBarry Smith   if (func) tr->postcheck    = func;
95*7cb011f5SBarry Smith   if (ctx)  tr->postcheckctx = ctx;
96*7cb011f5SBarry Smith   PetscFunctionReturn(0);
97*7cb011f5SBarry Smith }
98*7cb011f5SBarry Smith 
99*7cb011f5SBarry Smith /*@C
100*7cb011f5SBarry Smith    SNESNewtonTRGetPostCheck - Gets the post-check function
101*7cb011f5SBarry Smith 
102*7cb011f5SBarry Smith    Not collective
103*7cb011f5SBarry Smith 
104*7cb011f5SBarry Smith    Input Parameter:
105*7cb011f5SBarry Smith .  snes - the nonlinear solver context
106*7cb011f5SBarry Smith 
107*7cb011f5SBarry Smith    Output Parameters:
108*7cb011f5SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck()
109*7cb011f5SBarry Smith -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
110*7cb011f5SBarry Smith 
111*7cb011f5SBarry Smith    Level: intermediate
112*7cb011f5SBarry Smith 
113*7cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck()
114*7cb011f5SBarry Smith @*/
115*7cb011f5SBarry Smith PetscErrorCode  SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,void*),void **ctx)
116*7cb011f5SBarry Smith {
117*7cb011f5SBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
118*7cb011f5SBarry Smith 
119*7cb011f5SBarry Smith   PetscFunctionBegin;
120*7cb011f5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
121*7cb011f5SBarry Smith   if (func) *func = tr->postcheck;
122*7cb011f5SBarry Smith   if (ctx)  *ctx  = tr->postcheckctx;
123*7cb011f5SBarry Smith   PetscFunctionReturn(0);
124*7cb011f5SBarry Smith }
125*7cb011f5SBarry Smith 
126*7cb011f5SBarry Smith /*@C
127*7cb011f5SBarry Smith    SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation
128*7cb011f5SBarry Smith 
129*7cb011f5SBarry Smith    Logically Collective on snes
130*7cb011f5SBarry Smith 
131*7cb011f5SBarry Smith    Input Parameters:
132*7cb011f5SBarry Smith +  snes - the solver.  X - The last solution
133*7cb011f5SBarry Smith .  Y - The full step direction
134*7cb011f5SBarry Smith -  W - The updated solution, W = X + tao*Y for some tao
135*7cb011f5SBarry Smith 
136*7cb011f5SBarry Smith    Output Parameters:
137*7cb011f5SBarry Smith .  changed_W - Indicator if the new candidate solution W has been changed.
138*7cb011f5SBarry Smith 
139*7cb011f5SBarry Smith    Notes:
140*7cb011f5SBarry Smith      If Y is changed then W is recomputed as X + tao*Y where tao is the current update value in SNESNEWTONTR
141*7cb011f5SBarry Smith 
142*7cb011f5SBarry Smith    Level: developer
143*7cb011f5SBarry Smith 
144*7cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck()
145*7cb011f5SBarry Smith @*/
146*7cb011f5SBarry Smith static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_W)
147*7cb011f5SBarry Smith {
148*7cb011f5SBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
149*7cb011f5SBarry Smith   PetscErrorCode ierr;
150*7cb011f5SBarry Smith 
151*7cb011f5SBarry Smith   PetscFunctionBegin;
152*7cb011f5SBarry Smith   *changed_W = PETSC_FALSE;
153*7cb011f5SBarry Smith   if (tr->postcheck) {
154*7cb011f5SBarry Smith     ierr = (*tr->postcheck)(snes,X,Y,W,changed_W,tr->postcheckctx);CHKERRQ(ierr);
155*7cb011f5SBarry Smith     PetscValidLogicalCollectiveBool(snes,*changed_W,6);
156*7cb011f5SBarry Smith   }
157*7cb011f5SBarry Smith   PetscFunctionReturn(0);
158*7cb011f5SBarry Smith }
15985385478SLisandro Dalcin 
160df60cc22SBarry Smith /*
16104d7464bSBarry Smith    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
162ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
1634800dd8cSBarry Smith 
164ddbbdb52SLois Curfman McInnes 
1654800dd8cSBarry Smith */
16604d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
1674800dd8cSBarry Smith {
16804d7464bSBarry Smith   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
16985385478SLisandro Dalcin   Vec                      X,F,Y,G,Ytmp;
1706849ba73SBarry Smith   PetscErrorCode           ierr;
171a7cc72afSBarry Smith   PetscInt                 maxits,i,lits;
17285385478SLisandro Dalcin   PetscReal                rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
17379f36460SBarry Smith   PetscScalar              cnorm;
174df60cc22SBarry Smith   KSP                      ksp;
1753f149594SLisandro Dalcin   SNESConvergedReason      reason = SNES_CONVERGED_ITERATING;
176df8705c3SBarry Smith   PetscBool                breakout = PETSC_FALSE;
177df8705c3SBarry Smith   SNES_TR_KSPConverged_Ctx *ctx;
178df8705c3SBarry Smith   PetscErrorCode           (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
1794800dd8cSBarry Smith 
1803a40ed3dSBarry Smith   PetscFunctionBegin;
1816c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
182c579b300SPatrick Farrell 
183fbe28522SBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
184fbe28522SBarry Smith   X      = snes->vec_sol;               /* solution vector */
18539e2f89bSBarry Smith   F      = snes->vec_func;              /* residual vector */
186fbe28522SBarry Smith   Y      = snes->work[0];               /* work vectors */
187fbe28522SBarry Smith   G      = snes->work[1];
1886b5873e3SBarry Smith   Ytmp   = snes->work[2];
1894800dd8cSBarry Smith 
190e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1914c49b128SBarry Smith   snes->iter = 0;
192e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
1934800dd8cSBarry Smith 
194df8705c3SBarry Smith   /* Set the linear stopping criteria to use the More' trick. */
195df8705c3SBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
196df8705c3SBarry Smith   ierr = KSPGetConvergenceTest(ksp,&convtest,NULL,NULL);CHKERRQ(ierr);
197df8705c3SBarry Smith   if (convtest != SNESTR_KSPConverged_Private) {
198df8705c3SBarry Smith     ierr                  = PetscNew(&ctx);CHKERRQ(ierr);
199df8705c3SBarry Smith     ctx->snes             = snes;
200df8705c3SBarry Smith     ierr                  = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
201df8705c3SBarry Smith     ierr                  = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr);
202df8705c3SBarry Smith     ierr                  = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr);
203df8705c3SBarry Smith   }
204df8705c3SBarry Smith 
205e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
2065334005bSBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
2071aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
208e4ed7901SPeter Brune 
209cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
210422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
2118b84755bSBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
212e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
213fbe28522SBarry Smith   snes->norm = fnorm;
214e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
2158b84755bSBarry Smith   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;
2164800dd8cSBarry Smith   neP->delta = delta;
217a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
2187a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
219b37302c6SLois Curfman McInnes 
22085385478SLisandro Dalcin   /* test convergence */
22185385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
22285385478SLisandro Dalcin   if (snes->reason) PetscFunctionReturn(0);
2233f149594SLisandro Dalcin 
224df60cc22SBarry Smith 
2254800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
22699a96b7cSMatthew Knepley 
22799a96b7cSMatthew Knepley     /* Call general purpose update function */
228e7788613SBarry Smith     if (snes->ops->update) {
229e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
23099a96b7cSMatthew Knepley     }
23199a96b7cSMatthew Knepley 
23285385478SLisandro Dalcin     /* Solve J Y = F, where J is Jacobian matrix */
233d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
23407b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
23523ee1639SBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
236d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
2373f149594SLisandro Dalcin     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2381aa26658SKarl Rupp 
239329f5518SBarry Smith     snes->linear_its += lits;
2401aa26658SKarl Rupp 
241ae15b995SBarry Smith     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
242064f8208SBarry Smith     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
243064f8208SBarry Smith     norm1 = nrm;
2446b5873e3SBarry Smith     while (1) {
245*7cb011f5SBarry Smith       PetscBool changed_w;
246393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
247064f8208SBarry Smith       nrm  = norm1;
2486b5873e3SBarry Smith 
2492deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
250064f8208SBarry Smith       if (nrm >= delta) {
251064f8208SBarry Smith         nrm    = delta/nrm;
252064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
253064f8208SBarry Smith         cnorm  = nrm;
25457622a8eSBarry Smith         ierr   = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr);
2552dcb1b2aSMatthew Knepley         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
256064f8208SBarry Smith         nrm    = gpnorm;
257fbe28522SBarry Smith         ynorm  = delta;
258fbe28522SBarry Smith       } else {
259fbe28522SBarry Smith         gpnorm = 0.0;
260ae15b995SBarry Smith         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
261064f8208SBarry Smith         ynorm  = nrm;
262fbe28522SBarry Smith       }
2631592aa04SStefano Zampini       ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
26479f36460SBarry Smith       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
265*7cb011f5SBarry Smith       ierr = SNESNewtonTRPostCheck(snes,X,Ytmp,Y,&changed_w);CHKERRQ(ierr);
2665334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
267cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
2687551ef16SBarry Smith       SNESCheckFunctionNorm(snes,gnorm);
2694800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
270fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
2714800dd8cSBarry Smith 
2724800dd8cSBarry Smith       /* Update size of trust region */
2734800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
2744800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
2754800dd8cSBarry Smith       else                     delta *= neP->delta3;
27657622a8eSBarry Smith       ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr);
27757622a8eSBarry Smith       ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr);
2781aa26658SKarl Rupp 
2794800dd8cSBarry Smith       neP->delta = delta;
2804800dd8cSBarry Smith       if (rho > neP->sigma) break;
281ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
2826b5873e3SBarry Smith       /* check to see if progress is hopeless */
283ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
284470880abSPatrick Sanan       ierr        = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
28585385478SLisandro Dalcin       if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);}
286*7cb011f5SBarry Smith       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
287184914b5SBarry Smith       if (reason) {
28852392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
2897a03ce2fSLisandro Dalcin         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
290a7cc72afSBarry Smith         breakout = PETSC_TRUE;
291454a90a3SBarry Smith         break;
29252392280SLois Curfman McInnes       }
29350ffb88aSMatthew Knepley       snes->numFailures++;
2944800dd8cSBarry Smith     }
2951acabf8cSLois Curfman McInnes     if (!breakout) {
29685385478SLisandro Dalcin       /* Update function and solution vectors */
2974800dd8cSBarry Smith       fnorm = gnorm;
29885385478SLisandro Dalcin       ierr  = VecCopy(G,F);CHKERRQ(ierr);
29985385478SLisandro Dalcin       ierr  = VecCopy(Y,X);CHKERRQ(ierr);
30085385478SLisandro Dalcin       /* Monitor convergence */
301e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
302c509a162SBarry Smith       snes->iter = i+1;
303fbe28522SBarry Smith       snes->norm = fnorm;
304c1e67a49SFande Kong       snes->xnorm = xnorm;
305c1e67a49SFande Kong       snes->ynorm = ynorm;
306e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
307a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
3087a03ce2fSLisandro Dalcin       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
30985385478SLisandro Dalcin       /* Test for convergence, xnorm = || X || */
310ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
311e2a6519dSDmitry Karpeev       if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
312e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
3133f149594SLisandro Dalcin       if (reason) break;
3141aa26658SKarl Rupp     } else break;
31538442cffSBarry Smith   }
31652392280SLois Curfman McInnes   if (i == maxits) {
317ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
31885385478SLisandro Dalcin     if (!reason) reason = SNES_DIVERGED_MAX_IT;
31952392280SLois Curfman McInnes   }
320e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
321184914b5SBarry Smith   snes->reason = reason;
322e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
3244800dd8cSBarry Smith }
3254800dd8cSBarry Smith /*------------------------------------------------------------*/
32604d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
3274800dd8cSBarry Smith {
328dfbe8321SBarry Smith   PetscErrorCode ierr;
3293a40ed3dSBarry Smith 
3303a40ed3dSBarry Smith   PetscFunctionBegin;
331fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
3326cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
3344800dd8cSBarry Smith }
3356b8b9a38SLisandro Dalcin 
33604d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes)
3376b8b9a38SLisandro Dalcin {
3386b8b9a38SLisandro Dalcin 
3396b8b9a38SLisandro Dalcin   PetscFunctionBegin;
3406b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
3416b8b9a38SLisandro Dalcin }
3426b8b9a38SLisandro Dalcin 
34304d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
3444800dd8cSBarry Smith {
345dfbe8321SBarry Smith   PetscErrorCode ierr;
3465baf8537SBarry Smith 
3473a40ed3dSBarry Smith   PetscFunctionBegin;
34804d7464bSBarry Smith   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
349606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3503a40ed3dSBarry Smith   PetscFunctionReturn(0);
3514800dd8cSBarry Smith }
3524800dd8cSBarry Smith /*------------------------------------------------------------*/
3534800dd8cSBarry Smith 
3544416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes)
3554800dd8cSBarry Smith {
35604d7464bSBarry Smith   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
357dfbe8321SBarry Smith   PetscErrorCode ierr;
3584800dd8cSBarry Smith 
3593a40ed3dSBarry Smith   PetscFunctionBegin;
360e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
36194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
36294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr);
36394ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr);
36494ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr);
36594ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
36694ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr);
36794ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr);
36894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr);
369b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3703a40ed3dSBarry Smith   PetscFunctionReturn(0);
3714800dd8cSBarry Smith }
3724800dd8cSBarry Smith 
37304d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
374a935fc98SLois Curfman McInnes {
37504d7464bSBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
376dfbe8321SBarry Smith   PetscErrorCode ierr;
377ace3abfcSBarry Smith   PetscBool      iascii;
378a935fc98SLois Curfman McInnes 
3793a40ed3dSBarry Smith   PetscFunctionBegin;
380251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
38132077d6dSBarry Smith   if (iascii) {
38254fe5c21SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr);
38357622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr);
38457622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);CHKERRQ(ierr);
38519bcc07fSBarry Smith   }
3863a40ed3dSBarry Smith   PetscFunctionReturn(0);
387a935fc98SLois Curfman McInnes }
38852392280SLois Curfman McInnes /* ------------------------------------------------------------ */
389ebe3b25bSBarry Smith /*MC
39004d7464bSBarry Smith       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
391ebe3b25bSBarry Smith 
392ebe3b25bSBarry Smith    Options Database:
3938e434772SBarry Smith +    -snes_trtol <tol> - trust region tolerance
3948e434772SBarry Smith .    -snes_tr_mu <mu> - trust region parameter
3958e434772SBarry Smith .    -snes_tr_eta <eta> - trust region parameter
3968e434772SBarry Smith .    -snes_tr_sigma <sigma> - trust region parameter
3978b84755bSBarry Smith .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
3988e434772SBarry Smith .    -snes_tr_delta1 <delta1> - trust region parameter
3998e434772SBarry Smith .    -snes_tr_delta2 <delta2> - trust region parameter
4008e434772SBarry Smith -    -snes_tr_delta3 <delta3> - trust region parameter
401ebe3b25bSBarry Smith 
402ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
403ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
404ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
405ebe3b25bSBarry Smith 
406ee3001cbSBarry Smith    Level: intermediate
407ee3001cbSBarry Smith 
40804d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
409ebe3b25bSBarry Smith 
410ebe3b25bSBarry Smith M*/
4118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
4124800dd8cSBarry Smith {
41304d7464bSBarry Smith   SNES_NEWTONTR  *neP;
414dfbe8321SBarry Smith   PetscErrorCode ierr;
4154800dd8cSBarry Smith 
4163a40ed3dSBarry Smith   PetscFunctionBegin;
41704d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONTR;
41804d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONTR;
41904d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONTR;
42004d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
42104d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONTR;
42204d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONTR;
423fbe28522SBarry Smith 
42442f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
425efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
42642f4f86dSBarry Smith 
4274fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4284fc747eaSLawrence Mitchell 
429b00a9115SJed Brown   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
430fbe28522SBarry Smith   snes->data  = (void*)neP;
431fbe28522SBarry Smith   neP->mu     = 0.25;
432fbe28522SBarry Smith   neP->eta    = 0.75;
433fbe28522SBarry Smith   neP->delta  = 0.0;
434fbe28522SBarry Smith   neP->delta0 = 0.2;
435fbe28522SBarry Smith   neP->delta1 = 0.3;
436fbe28522SBarry Smith   neP->delta2 = 0.75;
437fbe28522SBarry Smith   neP->delta3 = 2.0;
438fbe28522SBarry Smith   neP->sigma  = 0.0001;
439ef87ac0dSKris Buschelman   neP->itflag = PETSC_FALSE;
44075567043SBarry Smith   neP->rnorm0 = 0.0;
44175567043SBarry Smith   neP->ttol   = 0.0;
4423a40ed3dSBarry Smith   PetscFunctionReturn(0);
4434800dd8cSBarry Smith }
44482bf6240SBarry Smith 
445