xref: /petsc/src/snes/impls/tr/tr.c (revision 3312a946b8f7725d48c4ba71ff091ede799e5ff5)
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);
621c6b2ff8SBarry Smith     *reason = SNES_DIVERGED_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 
707cb011f5SBarry Smith /*@C
71c9368356SGlenn Hammond    SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined.
72c9368356SGlenn Hammond        Allows the user a chance to change or override the decision of the line search routine.
73c9368356SGlenn Hammond 
74c9368356SGlenn Hammond    Logically Collective on snes
75c9368356SGlenn Hammond 
76c9368356SGlenn Hammond    Input Parameters:
77c9368356SGlenn Hammond +  snes - the nonlinear solver object
78c9368356SGlenn Hammond .  func - [optional] function evaluation routine, see SNESNewtonTRPreCheck()  for the calling sequence
79c9368356SGlenn Hammond -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
80c9368356SGlenn Hammond 
81c9368356SGlenn Hammond    Level: intermediate
82c9368356SGlenn Hammond 
83c9368356SGlenn Hammond    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver.
84c9368356SGlenn Hammond 
85*3312a946SBarry Smith .seealso: SNESNewtonTRPreCheck(), SNESNewtonTRGetPreCheck(), SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck()
86c9368356SGlenn Hammond @*/
87c9368356SGlenn Hammond PetscErrorCode  SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx)
88c9368356SGlenn Hammond {
89c9368356SGlenn Hammond   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
90c9368356SGlenn Hammond 
91c9368356SGlenn Hammond   PetscFunctionBegin;
92c9368356SGlenn Hammond   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
93c9368356SGlenn Hammond   if (func) tr->precheck    = func;
94c9368356SGlenn Hammond   if (ctx)  tr->precheckctx = ctx;
95c9368356SGlenn Hammond   PetscFunctionReturn(0);
96c9368356SGlenn Hammond }
97c9368356SGlenn Hammond 
98c9368356SGlenn Hammond /*@C
99c9368356SGlenn Hammond    SNESNewtonTRGetPreCheck - Gets the pre-check function
100c9368356SGlenn Hammond 
101c9368356SGlenn Hammond    Not collective
102c9368356SGlenn Hammond 
103c9368356SGlenn Hammond    Input Parameter:
104c9368356SGlenn Hammond .  snes - the nonlinear solver context
105c9368356SGlenn Hammond 
106c9368356SGlenn Hammond    Output Parameters:
107c9368356SGlenn Hammond +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck()
108c9368356SGlenn Hammond -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
109c9368356SGlenn Hammond 
110c9368356SGlenn Hammond    Level: intermediate
111c9368356SGlenn Hammond 
112c9368356SGlenn Hammond .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRPreCheck()
113c9368356SGlenn Hammond @*/
114c9368356SGlenn Hammond PetscErrorCode  SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx)
115c9368356SGlenn Hammond {
116c9368356SGlenn Hammond   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
117c9368356SGlenn Hammond 
118c9368356SGlenn Hammond   PetscFunctionBegin;
119c9368356SGlenn Hammond   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
120c9368356SGlenn Hammond   if (func) *func = tr->precheck;
121c9368356SGlenn Hammond   if (ctx)  *ctx  = tr->precheckctx;
122c9368356SGlenn Hammond   PetscFunctionReturn(0);
123c9368356SGlenn Hammond }
124c9368356SGlenn Hammond 
125c9368356SGlenn Hammond /*@C
1267cb011f5SBarry Smith    SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
1277cb011f5SBarry Smith        function evaluation. Allows the user a chance to change or override the decision of the line search routine
1287cb011f5SBarry Smith 
1297cb011f5SBarry Smith    Logically Collective on snes
1307cb011f5SBarry Smith 
1317cb011f5SBarry Smith    Input Parameters:
1327cb011f5SBarry Smith +  snes - the nonlinear solver object
1337cb011f5SBarry Smith .  func - [optional] function evaluation routine, see SNESNewtonTRPostCheck()  for the calling sequence
1347cb011f5SBarry Smith -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
1357cb011f5SBarry Smith 
1367cb011f5SBarry Smith    Level: intermediate
1377cb011f5SBarry Smith 
138c9368356SGlenn Hammond    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in
1397cb011f5SBarry Smith    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
1407cb011f5SBarry Smith 
1417cb011f5SBarry Smith .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck()
1427cb011f5SBarry Smith @*/
143c9368356SGlenn Hammond PetscErrorCode  SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
1447cb011f5SBarry Smith {
1457cb011f5SBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
1467cb011f5SBarry Smith 
1477cb011f5SBarry Smith   PetscFunctionBegin;
1487cb011f5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1497cb011f5SBarry Smith   if (func) tr->postcheck    = func;
1507cb011f5SBarry Smith   if (ctx)  tr->postcheckctx = ctx;
1517cb011f5SBarry Smith   PetscFunctionReturn(0);
1527cb011f5SBarry Smith }
1537cb011f5SBarry Smith 
1547cb011f5SBarry Smith /*@C
1557cb011f5SBarry Smith    SNESNewtonTRGetPostCheck - Gets the post-check function
1567cb011f5SBarry Smith 
1577cb011f5SBarry Smith    Not collective
1587cb011f5SBarry Smith 
1597cb011f5SBarry Smith    Input Parameter:
1607cb011f5SBarry Smith .  snes - the nonlinear solver context
1617cb011f5SBarry Smith 
1627cb011f5SBarry Smith    Output Parameters:
1637cb011f5SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck()
1647cb011f5SBarry Smith -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
1657cb011f5SBarry Smith 
1667cb011f5SBarry Smith    Level: intermediate
1677cb011f5SBarry Smith 
1687cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck()
1697cb011f5SBarry Smith @*/
170c9368356SGlenn Hammond PetscErrorCode  SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
1717cb011f5SBarry Smith {
1727cb011f5SBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
1737cb011f5SBarry Smith 
1747cb011f5SBarry Smith   PetscFunctionBegin;
1757cb011f5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1767cb011f5SBarry Smith   if (func) *func = tr->postcheck;
1777cb011f5SBarry Smith   if (ctx)  *ctx  = tr->postcheckctx;
1787cb011f5SBarry Smith   PetscFunctionReturn(0);
1797cb011f5SBarry Smith }
1807cb011f5SBarry Smith 
1817cb011f5SBarry Smith /*@C
182c9368356SGlenn Hammond    SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR
183c9368356SGlenn Hammond 
184c9368356SGlenn Hammond    Logically Collective on snes
185c9368356SGlenn Hammond 
186c9368356SGlenn Hammond    Input Parameters:
187c9368356SGlenn Hammond +  snes - the solver
188c9368356SGlenn Hammond .  X - The last solution
189c9368356SGlenn Hammond -  Y - The step direction
190c9368356SGlenn Hammond 
191c9368356SGlenn Hammond    Output Parameters:
192c9368356SGlenn Hammond .  changed_Y - Indicator that the step direction Y has been changed.
193c9368356SGlenn Hammond 
194c9368356SGlenn Hammond    Level: developer
195c9368356SGlenn Hammond 
196c9368356SGlenn Hammond .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck()
197c9368356SGlenn Hammond @*/
198c9368356SGlenn Hammond static PetscErrorCode SNESNewtonTRPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y)
199c9368356SGlenn Hammond {
200c9368356SGlenn Hammond   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
201c9368356SGlenn Hammond   PetscErrorCode ierr;
202c9368356SGlenn Hammond 
203c9368356SGlenn Hammond   PetscFunctionBegin;
204c9368356SGlenn Hammond   *changed_Y = PETSC_FALSE;
205c9368356SGlenn Hammond   if (tr->precheck) {
206c9368356SGlenn Hammond     ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr);
207c9368356SGlenn Hammond     PetscValidLogicalCollectiveBool(snes,*changed_Y,4);
208c9368356SGlenn Hammond   }
209c9368356SGlenn Hammond   PetscFunctionReturn(0);
210c9368356SGlenn Hammond }
211c9368356SGlenn Hammond 
212c9368356SGlenn Hammond /*@C
2137cb011f5SBarry Smith    SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation
2147cb011f5SBarry Smith 
2157cb011f5SBarry Smith    Logically Collective on snes
2167cb011f5SBarry Smith 
2177cb011f5SBarry Smith    Input Parameters:
2187cb011f5SBarry Smith +  snes - the solver.  X - The last solution
2197cb011f5SBarry Smith .  Y - The full step direction
220*3312a946SBarry Smith -  W - The updated solution, W = X - Y
2217cb011f5SBarry Smith 
2227cb011f5SBarry Smith    Output Parameters:
223*3312a946SBarry Smith +  changed_Y - indicator if step has been changed
224*3312a946SBarry Smith -  changed_W - Indicator if the new candidate solution W has been changed.
2257cb011f5SBarry Smith 
2267cb011f5SBarry Smith    Notes:
227*3312a946SBarry Smith      If Y is changed then W is recomputed as X - Y
2287cb011f5SBarry Smith 
2297cb011f5SBarry Smith    Level: developer
2307cb011f5SBarry Smith 
2317cb011f5SBarry Smith .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck()
2327cb011f5SBarry Smith @*/
233c9368356SGlenn Hammond static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
2347cb011f5SBarry Smith {
2357cb011f5SBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
2367cb011f5SBarry Smith   PetscErrorCode ierr;
2377cb011f5SBarry Smith 
2387cb011f5SBarry Smith   PetscFunctionBegin;
239c9368356SGlenn Hammond   *changed_Y = PETSC_FALSE;
2407cb011f5SBarry Smith   *changed_W = PETSC_FALSE;
2417cb011f5SBarry Smith   if (tr->postcheck) {
242c9368356SGlenn Hammond     ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr);
243c9368356SGlenn Hammond     PetscValidLogicalCollectiveBool(snes,*changed_Y,5);
2447cb011f5SBarry Smith     PetscValidLogicalCollectiveBool(snes,*changed_W,6);
2457cb011f5SBarry Smith   }
2467cb011f5SBarry Smith   PetscFunctionReturn(0);
2477cb011f5SBarry Smith }
24885385478SLisandro Dalcin 
249df60cc22SBarry Smith /*
25004d7464bSBarry Smith    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
251ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
2524800dd8cSBarry Smith 
253ddbbdb52SLois Curfman McInnes 
2544800dd8cSBarry Smith */
25504d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
2564800dd8cSBarry Smith {
25704d7464bSBarry Smith   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
258c9368356SGlenn Hammond   Vec                      X,F,Y,G,Ytmp,W;
2596849ba73SBarry Smith   PetscErrorCode           ierr;
260a7cc72afSBarry Smith   PetscInt                 maxits,i,lits;
26185385478SLisandro Dalcin   PetscReal                rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
26279f36460SBarry Smith   PetscScalar              cnorm;
263df60cc22SBarry Smith   KSP                      ksp;
2643f149594SLisandro Dalcin   SNESConvergedReason      reason = SNES_CONVERGED_ITERATING;
265df8705c3SBarry Smith   PetscBool                breakout = PETSC_FALSE;
266df8705c3SBarry Smith   SNES_TR_KSPConverged_Ctx *ctx;
267df8705c3SBarry Smith   PetscErrorCode           (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
2684800dd8cSBarry Smith 
2693a40ed3dSBarry Smith   PetscFunctionBegin;
2706c4ed002SBarry 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);
271c579b300SPatrick Farrell 
272fbe28522SBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
273fbe28522SBarry Smith   X      = snes->vec_sol;               /* solution vector */
27439e2f89bSBarry Smith   F      = snes->vec_func;              /* residual vector */
275fbe28522SBarry Smith   Y      = snes->work[0];               /* work vectors */
276fbe28522SBarry Smith   G      = snes->work[1];
2776b5873e3SBarry Smith   Ytmp   = snes->work[2];
278c9368356SGlenn Hammond   W      = snes->work[3];
2794800dd8cSBarry Smith 
280e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
2814c49b128SBarry Smith   snes->iter = 0;
282e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
2834800dd8cSBarry Smith 
284df8705c3SBarry Smith   /* Set the linear stopping criteria to use the More' trick. */
285df8705c3SBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
286df8705c3SBarry Smith   ierr = KSPGetConvergenceTest(ksp,&convtest,NULL,NULL);CHKERRQ(ierr);
287df8705c3SBarry Smith   if (convtest != SNESTR_KSPConverged_Private) {
288df8705c3SBarry Smith     ierr                  = PetscNew(&ctx);CHKERRQ(ierr);
289df8705c3SBarry Smith     ctx->snes             = snes;
290df8705c3SBarry Smith     ierr                  = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
291df8705c3SBarry Smith     ierr                  = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr);
292df8705c3SBarry Smith     ierr                  = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr);
293df8705c3SBarry Smith   }
294df8705c3SBarry Smith 
295e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
2965334005bSBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
2971aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
298e4ed7901SPeter Brune 
299cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
300422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
3018b84755bSBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
302e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
303fbe28522SBarry Smith   snes->norm = fnorm;
304e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
3058b84755bSBarry Smith   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;
3064800dd8cSBarry Smith   neP->delta = delta;
307a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
3087a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
309b37302c6SLois Curfman McInnes 
31085385478SLisandro Dalcin   /* test convergence */
31185385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
31285385478SLisandro Dalcin   if (snes->reason) PetscFunctionReturn(0);
3133f149594SLisandro Dalcin 
314df60cc22SBarry Smith 
3154800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
31699a96b7cSMatthew Knepley 
31799a96b7cSMatthew Knepley     /* Call general purpose update function */
318e7788613SBarry Smith     if (snes->ops->update) {
319e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
32099a96b7cSMatthew Knepley     }
32199a96b7cSMatthew Knepley 
32285385478SLisandro Dalcin     /* Solve J Y = F, where J is Jacobian matrix */
323d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
32407b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
32523ee1639SBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
326d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
3273f149594SLisandro Dalcin     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
3281aa26658SKarl Rupp 
329329f5518SBarry Smith     snes->linear_its += lits;
3301aa26658SKarl Rupp 
331ae15b995SBarry Smith     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
332064f8208SBarry Smith     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
333064f8208SBarry Smith     norm1 = nrm;
3346b5873e3SBarry Smith     while (1) {
335c9368356SGlenn Hammond       PetscBool changed_y;
3367cb011f5SBarry Smith       PetscBool changed_w;
337393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
338064f8208SBarry Smith       nrm  = norm1;
3396b5873e3SBarry Smith 
3402deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
341064f8208SBarry Smith       if (nrm >= delta) {
342064f8208SBarry Smith         nrm    = delta/nrm;
343064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
344064f8208SBarry Smith         cnorm  = nrm;
34557622a8eSBarry Smith         ierr   = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr);
3462dcb1b2aSMatthew Knepley         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
347064f8208SBarry Smith         nrm    = gpnorm;
348fbe28522SBarry Smith         ynorm  = delta;
349fbe28522SBarry Smith       } else {
350fbe28522SBarry Smith         gpnorm = 0.0;
351ae15b995SBarry Smith         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
352064f8208SBarry Smith         ynorm  = nrm;
353fbe28522SBarry Smith       }
354c9368356SGlenn Hammond       /* PreCheck() allows for updates to Y prior to W <- X - Y */
355c9368356SGlenn Hammond       ierr = SNESNewtonTRPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr);
356c9368356SGlenn Hammond       ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);         /* W <- X - Y */
357c9368356SGlenn Hammond       ierr = SNESNewtonTRPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
358c9368356SGlenn Hammond       if (changed_y) ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);
359c1cd17f6SGlenn Hammond       ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
360c9368356SGlenn Hammond       ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /*  F(X) */
361cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
3627551ef16SBarry Smith       SNESCheckFunctionNorm(snes,gnorm);
3634800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
364fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
3654800dd8cSBarry Smith 
3664800dd8cSBarry Smith       /* Update size of trust region */
3674800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
3684800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
3694800dd8cSBarry Smith       else                     delta *= neP->delta3;
37057622a8eSBarry Smith       ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr);
37157622a8eSBarry Smith       ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr);
3721aa26658SKarl Rupp 
3734800dd8cSBarry Smith       neP->delta = delta;
3744800dd8cSBarry Smith       if (rho > neP->sigma) break;
375ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
3766b5873e3SBarry Smith       /* check to see if progress is hopeless */
377ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
378470880abSPatrick Sanan       ierr        = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
37985385478SLisandro Dalcin       if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);}
3807cb011f5SBarry Smith       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
381184914b5SBarry Smith       if (reason) {
38252392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
3837a03ce2fSLisandro Dalcin         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
384a7cc72afSBarry Smith         breakout = PETSC_TRUE;
385454a90a3SBarry Smith         break;
38652392280SLois Curfman McInnes       }
38750ffb88aSMatthew Knepley       snes->numFailures++;
3884800dd8cSBarry Smith     }
3891acabf8cSLois Curfman McInnes     if (!breakout) {
39085385478SLisandro Dalcin       /* Update function and solution vectors */
3914800dd8cSBarry Smith       fnorm = gnorm;
39285385478SLisandro Dalcin       ierr  = VecCopy(G,F);CHKERRQ(ierr);
393c9368356SGlenn Hammond       ierr  = VecCopy(W,X);CHKERRQ(ierr);
39485385478SLisandro Dalcin       /* Monitor convergence */
395e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
396c509a162SBarry Smith       snes->iter = i+1;
397fbe28522SBarry Smith       snes->norm = fnorm;
398c1e67a49SFande Kong       snes->xnorm = xnorm;
399c1e67a49SFande Kong       snes->ynorm = ynorm;
400e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
401a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
4027a03ce2fSLisandro Dalcin       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
40385385478SLisandro Dalcin       /* Test for convergence, xnorm = || X || */
404ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
405e2a6519dSDmitry Karpeev       if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
406e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
4073f149594SLisandro Dalcin       if (reason) break;
4081aa26658SKarl Rupp     } else break;
40938442cffSBarry Smith   }
41052392280SLois Curfman McInnes   if (i == maxits) {
411ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
41285385478SLisandro Dalcin     if (!reason) reason = SNES_DIVERGED_MAX_IT;
41352392280SLois Curfman McInnes   }
414e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
415184914b5SBarry Smith   snes->reason = reason;
416e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
4173a40ed3dSBarry Smith   PetscFunctionReturn(0);
4184800dd8cSBarry Smith }
4194800dd8cSBarry Smith /*------------------------------------------------------------*/
42004d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
4214800dd8cSBarry Smith {
422dfbe8321SBarry Smith   PetscErrorCode ierr;
4233a40ed3dSBarry Smith 
4243a40ed3dSBarry Smith   PetscFunctionBegin;
425c9368356SGlenn Hammond   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
4266cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4273a40ed3dSBarry Smith   PetscFunctionReturn(0);
4284800dd8cSBarry Smith }
4296b8b9a38SLisandro Dalcin 
43004d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes)
4316b8b9a38SLisandro Dalcin {
4326b8b9a38SLisandro Dalcin 
4336b8b9a38SLisandro Dalcin   PetscFunctionBegin;
4346b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
4356b8b9a38SLisandro Dalcin }
4366b8b9a38SLisandro Dalcin 
43704d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
4384800dd8cSBarry Smith {
439dfbe8321SBarry Smith   PetscErrorCode ierr;
4405baf8537SBarry Smith 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
44204d7464bSBarry Smith   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
443606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4443a40ed3dSBarry Smith   PetscFunctionReturn(0);
4454800dd8cSBarry Smith }
4464800dd8cSBarry Smith /*------------------------------------------------------------*/
4474800dd8cSBarry Smith 
4484416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes)
4494800dd8cSBarry Smith {
45004d7464bSBarry Smith   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
451dfbe8321SBarry Smith   PetscErrorCode ierr;
4524800dd8cSBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
454e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
45594ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
45694ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr);
45794ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr);
45894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr);
45994ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
46094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr);
46194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr);
46294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr);
463b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
4654800dd8cSBarry Smith }
4664800dd8cSBarry Smith 
46704d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
468a935fc98SLois Curfman McInnes {
46904d7464bSBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
470dfbe8321SBarry Smith   PetscErrorCode ierr;
471ace3abfcSBarry Smith   PetscBool      iascii;
472a935fc98SLois Curfman McInnes 
4733a40ed3dSBarry Smith   PetscFunctionBegin;
474251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
47532077d6dSBarry Smith   if (iascii) {
47654fe5c21SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr);
47757622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr);
47857622a8eSBarry 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);
47919bcc07fSBarry Smith   }
4803a40ed3dSBarry Smith   PetscFunctionReturn(0);
481a935fc98SLois Curfman McInnes }
48252392280SLois Curfman McInnes /* ------------------------------------------------------------ */
483ebe3b25bSBarry Smith /*MC
48404d7464bSBarry Smith       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
485ebe3b25bSBarry Smith 
486ebe3b25bSBarry Smith    Options Database:
4878e434772SBarry Smith +    -snes_trtol <tol> - trust region tolerance
4888e434772SBarry Smith .    -snes_tr_mu <mu> - trust region parameter
4898e434772SBarry Smith .    -snes_tr_eta <eta> - trust region parameter
4908e434772SBarry Smith .    -snes_tr_sigma <sigma> - trust region parameter
4918b84755bSBarry Smith .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
4928e434772SBarry Smith .    -snes_tr_delta1 <delta1> - trust region parameter
4938e434772SBarry Smith .    -snes_tr_delta2 <delta2> - trust region parameter
4948e434772SBarry Smith -    -snes_tr_delta3 <delta3> - trust region parameter
495ebe3b25bSBarry Smith 
496ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
497ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
498ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
499ebe3b25bSBarry Smith 
500ee3001cbSBarry Smith    Level: intermediate
501ee3001cbSBarry Smith 
50204d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
503ebe3b25bSBarry Smith 
504ebe3b25bSBarry Smith M*/
5058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
5064800dd8cSBarry Smith {
50704d7464bSBarry Smith   SNES_NEWTONTR  *neP;
508dfbe8321SBarry Smith   PetscErrorCode ierr;
5094800dd8cSBarry Smith 
5103a40ed3dSBarry Smith   PetscFunctionBegin;
51104d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONTR;
51204d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONTR;
51304d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONTR;
51404d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
51504d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONTR;
51604d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONTR;
517fbe28522SBarry Smith 
51842f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
519efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
52042f4f86dSBarry Smith 
5214fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5224fc747eaSLawrence Mitchell 
523b00a9115SJed Brown   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
524fbe28522SBarry Smith   snes->data  = (void*)neP;
525fbe28522SBarry Smith   neP->mu     = 0.25;
526fbe28522SBarry Smith   neP->eta    = 0.75;
527fbe28522SBarry Smith   neP->delta  = 0.0;
528fbe28522SBarry Smith   neP->delta0 = 0.2;
529fbe28522SBarry Smith   neP->delta1 = 0.3;
530fbe28522SBarry Smith   neP->delta2 = 0.75;
531fbe28522SBarry Smith   neP->delta3 = 2.0;
532fbe28522SBarry Smith   neP->sigma  = 0.0001;
533ef87ac0dSKris Buschelman   neP->itflag = PETSC_FALSE;
53475567043SBarry Smith   neP->rnorm0 = 0.0;
53575567043SBarry Smith   neP->ttol   = 0.0;
5363a40ed3dSBarry Smith   PetscFunctionReturn(0);
5374800dd8cSBarry Smith }
53882bf6240SBarry Smith 
539