xref: /petsc/src/snes/impls/tr/tr.c (revision 1c6b2ff8df46f8a4f5ce7a45c746a507d77a9320)
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);
62*1c6b2ff8SBarry 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 
7085385478SLisandro Dalcin 
71df60cc22SBarry Smith /*
7204d7464bSBarry Smith    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
73ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
744800dd8cSBarry Smith 
75ddbbdb52SLois Curfman McInnes 
764800dd8cSBarry Smith */
7704d7464bSBarry Smith static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
784800dd8cSBarry Smith {
7904d7464bSBarry Smith   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
8085385478SLisandro Dalcin   Vec                      X,F,Y,G,Ytmp;
816849ba73SBarry Smith   PetscErrorCode           ierr;
82a7cc72afSBarry Smith   PetscInt                 maxits,i,lits;
8385385478SLisandro Dalcin   PetscReal                rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
8479f36460SBarry Smith   PetscScalar              cnorm;
85df60cc22SBarry Smith   KSP                      ksp;
863f149594SLisandro Dalcin   SNESConvergedReason      reason = SNES_CONVERGED_ITERATING;
87df8705c3SBarry Smith   PetscBool                breakout = PETSC_FALSE;
88df8705c3SBarry Smith   SNES_TR_KSPConverged_Ctx *ctx;
89df8705c3SBarry Smith   PetscErrorCode           (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
904800dd8cSBarry Smith 
913a40ed3dSBarry Smith   PetscFunctionBegin;
926c4ed002SBarry 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);
93c579b300SPatrick Farrell 
94fbe28522SBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
95fbe28522SBarry Smith   X      = snes->vec_sol;               /* solution vector */
9639e2f89bSBarry Smith   F      = snes->vec_func;              /* residual vector */
97fbe28522SBarry Smith   Y      = snes->work[0];               /* work vectors */
98fbe28522SBarry Smith   G      = snes->work[1];
996b5873e3SBarry Smith   Ytmp   = snes->work[2];
1004800dd8cSBarry Smith 
101e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1024c49b128SBarry Smith   snes->iter = 0;
103e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
1044800dd8cSBarry Smith 
105df8705c3SBarry Smith   /* Set the linear stopping criteria to use the More' trick. */
106df8705c3SBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
107df8705c3SBarry Smith   ierr = KSPGetConvergenceTest(ksp,&convtest,NULL,NULL);CHKERRQ(ierr);
108df8705c3SBarry Smith   if (convtest != SNESTR_KSPConverged_Private) {
109df8705c3SBarry Smith     ierr                  = PetscNew(&ctx);CHKERRQ(ierr);
110df8705c3SBarry Smith     ctx->snes             = snes;
111df8705c3SBarry Smith     ierr                  = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
112df8705c3SBarry Smith     ierr                  = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr);
113df8705c3SBarry Smith     ierr                  = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr);
114df8705c3SBarry Smith   }
115df8705c3SBarry Smith 
116e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
1175334005bSBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
1181aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
119e4ed7901SPeter Brune 
120cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
121422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
1228b84755bSBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
123e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
124fbe28522SBarry Smith   snes->norm = fnorm;
125e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
1268b84755bSBarry Smith   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;
1274800dd8cSBarry Smith   neP->delta = delta;
128a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
1297a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
130b37302c6SLois Curfman McInnes 
13185385478SLisandro Dalcin   /* test convergence */
13285385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
13385385478SLisandro Dalcin   if (snes->reason) PetscFunctionReturn(0);
1343f149594SLisandro Dalcin 
135df60cc22SBarry Smith 
1364800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
13799a96b7cSMatthew Knepley 
13899a96b7cSMatthew Knepley     /* Call general purpose update function */
139e7788613SBarry Smith     if (snes->ops->update) {
140e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
14199a96b7cSMatthew Knepley     }
14299a96b7cSMatthew Knepley 
14385385478SLisandro Dalcin     /* Solve J Y = F, where J is Jacobian matrix */
144d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
14507b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
14623ee1639SBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
147d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
1483f149594SLisandro Dalcin     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1491aa26658SKarl Rupp 
150329f5518SBarry Smith     snes->linear_its += lits;
1511aa26658SKarl Rupp 
152ae15b995SBarry Smith     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
153064f8208SBarry Smith     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
154064f8208SBarry Smith     norm1 = nrm;
1556b5873e3SBarry Smith     while (1) {
156393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
157064f8208SBarry Smith       nrm  = norm1;
1586b5873e3SBarry Smith 
1592deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
160064f8208SBarry Smith       if (nrm >= delta) {
161064f8208SBarry Smith         nrm    = delta/nrm;
162064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
163064f8208SBarry Smith         cnorm  = nrm;
16457622a8eSBarry Smith         ierr   = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr);
1652dcb1b2aSMatthew Knepley         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
166064f8208SBarry Smith         nrm    = gpnorm;
167fbe28522SBarry Smith         ynorm  = delta;
168fbe28522SBarry Smith       } else {
169fbe28522SBarry Smith         gpnorm = 0.0;
170ae15b995SBarry Smith         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
171064f8208SBarry Smith         ynorm  = nrm;
172fbe28522SBarry Smith       }
1731592aa04SStefano Zampini       ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
17479f36460SBarry Smith       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
1755334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
176cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1777551ef16SBarry Smith       SNESCheckFunctionNorm(snes,gnorm);
1784800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
179fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1804800dd8cSBarry Smith 
1814800dd8cSBarry Smith       /* Update size of trust region */
1824800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1834800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1844800dd8cSBarry Smith       else                     delta *= neP->delta3;
18557622a8eSBarry Smith       ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr);
18657622a8eSBarry Smith       ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr);
1871aa26658SKarl Rupp 
1884800dd8cSBarry Smith       neP->delta = delta;
1894800dd8cSBarry Smith       if (rho > neP->sigma) break;
190ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
1916b5873e3SBarry Smith       /* check to see if progress is hopeless */
192ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
193470880abSPatrick Sanan       ierr        = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
19485385478SLisandro Dalcin       if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); }
195184914b5SBarry Smith       if (reason) {
19652392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
1977a03ce2fSLisandro Dalcin         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
198a7cc72afSBarry Smith         breakout = PETSC_TRUE;
199454a90a3SBarry Smith         break;
20052392280SLois Curfman McInnes       }
20150ffb88aSMatthew Knepley       snes->numFailures++;
2024800dd8cSBarry Smith     }
2031acabf8cSLois Curfman McInnes     if (!breakout) {
20485385478SLisandro Dalcin       /* Update function and solution vectors */
2054800dd8cSBarry Smith       fnorm = gnorm;
20685385478SLisandro Dalcin       ierr  = VecCopy(G,F);CHKERRQ(ierr);
20785385478SLisandro Dalcin       ierr  = VecCopy(Y,X);CHKERRQ(ierr);
20885385478SLisandro Dalcin       /* Monitor convergence */
209e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
210c509a162SBarry Smith       snes->iter = i+1;
211fbe28522SBarry Smith       snes->norm = fnorm;
212c1e67a49SFande Kong       snes->xnorm = xnorm;
213c1e67a49SFande Kong       snes->ynorm = ynorm;
214e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
215a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
2167a03ce2fSLisandro Dalcin       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
21785385478SLisandro Dalcin       /* Test for convergence, xnorm = || X || */
218ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
219e2a6519dSDmitry Karpeev       if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
220e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
2213f149594SLisandro Dalcin       if (reason) break;
2221aa26658SKarl Rupp     } else break;
22338442cffSBarry Smith   }
22452392280SLois Curfman McInnes   if (i == maxits) {
225ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
22685385478SLisandro Dalcin     if (!reason) reason = SNES_DIVERGED_MAX_IT;
22752392280SLois Curfman McInnes   }
228e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
229184914b5SBarry Smith   snes->reason = reason;
230e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
2313a40ed3dSBarry Smith   PetscFunctionReturn(0);
2324800dd8cSBarry Smith }
2334800dd8cSBarry Smith /*------------------------------------------------------------*/
23404d7464bSBarry Smith static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
2354800dd8cSBarry Smith {
236dfbe8321SBarry Smith   PetscErrorCode ierr;
2373a40ed3dSBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
239fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
2406cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
2413a40ed3dSBarry Smith   PetscFunctionReturn(0);
2424800dd8cSBarry Smith }
2436b8b9a38SLisandro Dalcin 
24404d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONTR(SNES snes)
2456b8b9a38SLisandro Dalcin {
2466b8b9a38SLisandro Dalcin 
2476b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2486b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
2496b8b9a38SLisandro Dalcin }
2506b8b9a38SLisandro Dalcin 
25104d7464bSBarry Smith static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
2524800dd8cSBarry Smith {
253dfbe8321SBarry Smith   PetscErrorCode ierr;
2545baf8537SBarry Smith 
2553a40ed3dSBarry Smith   PetscFunctionBegin;
25604d7464bSBarry Smith   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
257606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2583a40ed3dSBarry Smith   PetscFunctionReturn(0);
2594800dd8cSBarry Smith }
2604800dd8cSBarry Smith /*------------------------------------------------------------*/
2614800dd8cSBarry Smith 
2624416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes)
2634800dd8cSBarry Smith {
26404d7464bSBarry Smith   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
265dfbe8321SBarry Smith   PetscErrorCode ierr;
2664800dd8cSBarry Smith 
2673a40ed3dSBarry Smith   PetscFunctionBegin;
268e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
26994ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
27094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr);
27194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr);
27294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr);
27394ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
27494ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr);
27594ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr);
27694ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr);
277b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2783a40ed3dSBarry Smith   PetscFunctionReturn(0);
2794800dd8cSBarry Smith }
2804800dd8cSBarry Smith 
28104d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
282a935fc98SLois Curfman McInnes {
28304d7464bSBarry Smith   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
284dfbe8321SBarry Smith   PetscErrorCode ierr;
285ace3abfcSBarry Smith   PetscBool      iascii;
286a935fc98SLois Curfman McInnes 
2873a40ed3dSBarry Smith   PetscFunctionBegin;
288251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
28932077d6dSBarry Smith   if (iascii) {
29054fe5c21SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr);
29157622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr);
29257622a8eSBarry 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);
29319bcc07fSBarry Smith   }
2943a40ed3dSBarry Smith   PetscFunctionReturn(0);
295a935fc98SLois Curfman McInnes }
29652392280SLois Curfman McInnes /* ------------------------------------------------------------ */
297ebe3b25bSBarry Smith /*MC
29804d7464bSBarry Smith       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
299ebe3b25bSBarry Smith 
300ebe3b25bSBarry Smith    Options Database:
3018e434772SBarry Smith +    -snes_trtol <tol> - trust region tolerance
3028e434772SBarry Smith .    -snes_tr_mu <mu> - trust region parameter
3038e434772SBarry Smith .    -snes_tr_eta <eta> - trust region parameter
3048e434772SBarry Smith .    -snes_tr_sigma <sigma> - trust region parameter
3058b84755bSBarry Smith .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
3068e434772SBarry Smith .    -snes_tr_delta1 <delta1> - trust region parameter
3078e434772SBarry Smith .    -snes_tr_delta2 <delta2> - trust region parameter
3088e434772SBarry Smith -    -snes_tr_delta3 <delta3> - trust region parameter
309ebe3b25bSBarry Smith 
310ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
311ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
312ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
313ebe3b25bSBarry Smith 
314ee3001cbSBarry Smith    Level: intermediate
315ee3001cbSBarry Smith 
31604d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
317ebe3b25bSBarry Smith 
318ebe3b25bSBarry Smith M*/
3198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
3204800dd8cSBarry Smith {
32104d7464bSBarry Smith   SNES_NEWTONTR  *neP;
322dfbe8321SBarry Smith   PetscErrorCode ierr;
3234800dd8cSBarry Smith 
3243a40ed3dSBarry Smith   PetscFunctionBegin;
32504d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONTR;
32604d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONTR;
32704d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONTR;
32804d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
32904d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONTR;
33004d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONTR;
331fbe28522SBarry Smith 
33242f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
333efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
33442f4f86dSBarry Smith 
3354fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
3364fc747eaSLawrence Mitchell 
337b00a9115SJed Brown   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
338fbe28522SBarry Smith   snes->data  = (void*)neP;
339fbe28522SBarry Smith   neP->mu     = 0.25;
340fbe28522SBarry Smith   neP->eta    = 0.75;
341fbe28522SBarry Smith   neP->delta  = 0.0;
342fbe28522SBarry Smith   neP->delta0 = 0.2;
343fbe28522SBarry Smith   neP->delta1 = 0.3;
344fbe28522SBarry Smith   neP->delta2 = 0.75;
345fbe28522SBarry Smith   neP->delta3 = 2.0;
346fbe28522SBarry Smith   neP->sigma  = 0.0001;
347ef87ac0dSKris Buschelman   neP->itflag = PETSC_FALSE;
34875567043SBarry Smith   neP->rnorm0 = 0.0;
34975567043SBarry Smith   neP->ttol   = 0.0;
3503a40ed3dSBarry Smith   PetscFunctionReturn(0);
3514800dd8cSBarry Smith }
35282bf6240SBarry Smith 
353