xref: /petsc/src/snes/impls/tr/tr.c (revision e4ed7901070ce1ca36eb7d68dd07557542e66e31)
14800dd8cSBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h>                /*I   "petscsnes.h"   I*/
34800dd8cSBarry Smith 
4971273eeSBarry Smith typedef struct {
5971273eeSBarry Smith   void *ctx;
6971273eeSBarry Smith   SNES snes;
7971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx;
8971273eeSBarry Smith 
9fbe28522SBarry Smith /*
10df60cc22SBarry Smith    This convergence test determines if the two norm of the
11df60cc22SBarry Smith    solution lies outside the trust region, if so it halts.
12df60cc22SBarry Smith */
134a2ae208SSatish Balay #undef __FUNCT__
144b27c08aSLois Curfman McInnes #define __FUNCT__ "SNES_TR_KSPConverged_Private"
15971273eeSBarry Smith PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx)
16df60cc22SBarry Smith {
17971273eeSBarry Smith   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx;
18971273eeSBarry Smith   SNES                     snes = ctx->snes;
194b27c08aSLois Curfman McInnes   SNES_TR                  *neP = (SNES_TR*)snes->data;
20df60cc22SBarry Smith   Vec                      x;
21064f8208SBarry Smith   PetscReal                nrm;
22dfbe8321SBarry Smith   PetscErrorCode           ierr;
23df60cc22SBarry Smith 
243a40ed3dSBarry Smith   PetscFunctionBegin;
25971273eeSBarry Smith   ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx->ctx);CHKERRQ(ierr);
26329f5518SBarry Smith   if (*reason) {
2771f87433Sdalcinl     ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr);
289b04593eSLois Curfman McInnes   }
29a935fc98SLois Curfman McInnes   /* Determine norm of solution */
3078b31e54SBarry Smith   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
31064f8208SBarry Smith   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
32064f8208SBarry Smith   if (nrm >= neP->delta) {
33ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr);
34329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
35df60cc22SBarry Smith   }
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
37df60cc22SBarry Smith }
3882bf6240SBarry Smith 
39971273eeSBarry Smith #undef __FUNCT__
40971273eeSBarry Smith #define __FUNCT__ "SNES_TR_KSPConverged_Destroy"
41971273eeSBarry Smith PetscErrorCode SNES_TR_KSPConverged_Destroy(void *cctx)
42971273eeSBarry Smith {
43971273eeSBarry Smith   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx;
44971273eeSBarry Smith   PetscErrorCode           ierr;
45971273eeSBarry Smith 
46971273eeSBarry Smith   PetscFunctionBegin;
47971273eeSBarry Smith   ierr = KSPDefaultConvergedDestroy(ctx->ctx);CHKERRQ(ierr);
48971273eeSBarry Smith   ierr = PetscFree(ctx);CHKERRQ(ierr);
49971273eeSBarry Smith   PetscFunctionReturn(0);
50971273eeSBarry Smith }
51971273eeSBarry Smith 
5285385478SLisandro Dalcin /* ---------------------------------------------------------------- */
5385385478SLisandro Dalcin #undef __FUNCT__
5485385478SLisandro Dalcin #define __FUNCT__ "SNES_TR_Converged_Private"
5585385478SLisandro Dalcin /*
5685385478SLisandro Dalcin    SNES_TR_Converged_Private -test convergence JUST for
5785385478SLisandro Dalcin    the trust region tolerance.
5885385478SLisandro Dalcin 
5985385478SLisandro Dalcin */
6085385478SLisandro Dalcin static PetscErrorCode SNES_TR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
6185385478SLisandro Dalcin {
6285385478SLisandro Dalcin   SNES_TR        *neP = (SNES_TR *)snes->data;
6385385478SLisandro Dalcin   PetscErrorCode ierr;
6485385478SLisandro Dalcin 
6585385478SLisandro Dalcin   PetscFunctionBegin;
6685385478SLisandro Dalcin   *reason = SNES_CONVERGED_ITERATING;
6785385478SLisandro Dalcin   if (neP->delta < xnorm * snes->deltatol) {
6885385478SLisandro Dalcin     ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,snes->deltatol);CHKERRQ(ierr);
6985385478SLisandro Dalcin     *reason = SNES_CONVERGED_TR_DELTA;
7085385478SLisandro Dalcin   } else if (snes->nfuncs >= snes->max_funcs) {
7185385478SLisandro Dalcin     ierr = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr);
7285385478SLisandro Dalcin     *reason = SNES_DIVERGED_FUNCTION_COUNT;
7385385478SLisandro Dalcin   }
7485385478SLisandro Dalcin   PetscFunctionReturn(0);
7585385478SLisandro Dalcin }
7685385478SLisandro Dalcin 
7785385478SLisandro Dalcin 
78df60cc22SBarry Smith /*
794b27c08aSLois Curfman McInnes    SNESSolve_TR - Implements Newton's Method with a very simple trust
80ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
814800dd8cSBarry Smith 
82ddbbdb52SLois Curfman McInnes 
834800dd8cSBarry Smith */
844a2ae208SSatish Balay #undef __FUNCT__
854b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR"
866849ba73SBarry Smith static PetscErrorCode SNESSolve_TR(SNES snes)
874800dd8cSBarry Smith {
884b27c08aSLois Curfman McInnes   SNES_TR             *neP = (SNES_TR*)snes->data;
8985385478SLisandro Dalcin   Vec                 X,F,Y,G,Ytmp;
906849ba73SBarry Smith   PetscErrorCode      ierr;
91a7cc72afSBarry Smith   PetscInt            maxits,i,lits;
92112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
9385385478SLisandro Dalcin   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
9479f36460SBarry Smith   PetscScalar         cnorm;
95df60cc22SBarry Smith   KSP                 ksp;
963f149594SLisandro Dalcin   SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
97ace3abfcSBarry Smith   PetscBool           conv = PETSC_FALSE,breakout = PETSC_FALSE;
98*e4ed7901SPeter Brune   PetscBool          domainerror;
994800dd8cSBarry Smith 
1003a40ed3dSBarry Smith   PetscFunctionBegin;
101fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
102fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
10339e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
104fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
105fbe28522SBarry Smith   G		= snes->work[1];
1066b5873e3SBarry Smith   Ytmp          = snes->work[2];
1074800dd8cSBarry Smith 
1084c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1094c49b128SBarry Smith   snes->iter = 0;
1104c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1114800dd8cSBarry Smith 
112*e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
1135334005bSBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
114*e4ed7901SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
115*e4ed7901SPeter Brune     if (domainerror) {
116*e4ed7901SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
117*e4ed7901SPeter Brune       PetscFunctionReturn(0);
118*e4ed7901SPeter Brune     }
119*e4ed7901SPeter Brune   } else {
120*e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
121*e4ed7901SPeter Brune   }
122*e4ed7901SPeter Brune 
123*e4ed7901SPeter Brune   if (!snes->norm_init_set) {
124cddf8d76SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
125*e4ed7901SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
126*e4ed7901SPeter Brune   } else {
127*e4ed7901SPeter Brune     fnorm = snes->norm_init;
128*e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
129*e4ed7901SPeter Brune   }
130*e4ed7901SPeter Brune 
1310f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
132fbe28522SBarry Smith   snes->norm = fnorm;
1330f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1344800dd8cSBarry Smith   delta = neP->delta0*fnorm;
1354800dd8cSBarry Smith   neP->delta = delta;
1360462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
1377a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
138b37302c6SLois Curfman McInnes 
139d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
140d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
14185385478SLisandro Dalcin   /* test convergence */
14285385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
14385385478SLisandro Dalcin   if (snes->reason) PetscFunctionReturn(0);
1443f149594SLisandro Dalcin 
145a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
146acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv,PETSC_NULL);CHKERRQ(ierr);
1477c4af3dcSLois Curfman McInnes   if (!conv) {
148971273eeSBarry Smith     SNES_TR_KSPConverged_Ctx *ctx;
1493f149594SLisandro Dalcin     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
150971273eeSBarry Smith     ierr = PetscNew(SNES_TR_KSPConverged_Ctx,&ctx);CHKERRQ(ierr);
151971273eeSBarry Smith     ctx->snes = snes;
152971273eeSBarry Smith     ierr = KSPDefaultConvergedCreate(&ctx->ctx);CHKERRQ(ierr);
153971273eeSBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr);
154ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
1557c4af3dcSLois Curfman McInnes   }
156df60cc22SBarry Smith 
1574800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
15899a96b7cSMatthew Knepley 
15999a96b7cSMatthew Knepley     /* Call general purpose update function */
160e7788613SBarry Smith     if (snes->ops->update) {
161e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
16299a96b7cSMatthew Knepley     }
16399a96b7cSMatthew Knepley 
16485385478SLisandro Dalcin     /* Solve J Y = F, where J is Jacobian matrix */
1655334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
16694b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1673f149594SLisandro Dalcin     ierr = SNES_KSPSolve(snes,snes->ksp,F,Ytmp);CHKERRQ(ierr);
1683f149594SLisandro Dalcin     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
169329f5518SBarry Smith     snes->linear_its += lits;
170ae15b995SBarry Smith     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
171064f8208SBarry Smith     ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
172064f8208SBarry Smith     norm1 = nrm;
1736b5873e3SBarry Smith     while(1) {
174393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
175064f8208SBarry Smith       nrm = norm1;
1766b5873e3SBarry Smith 
1772deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
178064f8208SBarry Smith       if (nrm >= delta) {
179064f8208SBarry Smith         nrm = delta/nrm;
180064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
181064f8208SBarry Smith         cnorm = nrm;
182ae15b995SBarry Smith         ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr);
1832dcb1b2aSMatthew Knepley         ierr = VecScale(Y,cnorm);CHKERRQ(ierr);
184064f8208SBarry Smith         nrm = gpnorm;
185fbe28522SBarry Smith         ynorm = delta;
186fbe28522SBarry Smith       } else {
187fbe28522SBarry Smith         gpnorm = 0.0;
188ae15b995SBarry Smith         ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
189064f8208SBarry Smith         ynorm = nrm;
190fbe28522SBarry Smith       }
19179f36460SBarry Smith       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
19285385478SLisandro Dalcin       ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr);
1935334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
194cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1954800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
196fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1974800dd8cSBarry Smith 
1984800dd8cSBarry Smith       /* Update size of trust region */
1994800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
2004800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
2014800dd8cSBarry Smith       else                     delta *= neP->delta3;
202ae15b995SBarry Smith       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
203ae15b995SBarry Smith       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
2044800dd8cSBarry Smith       neP->delta = delta;
2054800dd8cSBarry Smith       if (rho > neP->sigma) break;
206ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
2076b5873e3SBarry Smith       /* check to see if progress is hopeless */
208ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
20985385478SLisandro Dalcin       ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
21085385478SLisandro Dalcin       if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); }
211184914b5SBarry Smith       if (reason) {
21252392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
2137a03ce2fSLisandro Dalcin         ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
214a7cc72afSBarry Smith         breakout = PETSC_TRUE;
215454a90a3SBarry Smith         break;
21652392280SLois Curfman McInnes       }
21750ffb88aSMatthew Knepley       snes->numFailures++;
2184800dd8cSBarry Smith     }
2191acabf8cSLois Curfman McInnes     if (!breakout) {
22085385478SLisandro Dalcin       /* Update function and solution vectors */
2214800dd8cSBarry Smith       fnorm = gnorm;
22285385478SLisandro Dalcin       ierr = VecCopy(G,F);CHKERRQ(ierr);
22385385478SLisandro Dalcin       ierr = VecCopy(Y,X);CHKERRQ(ierr);
22485385478SLisandro Dalcin       /* Monitor convergence */
2250f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
226c509a162SBarry Smith       snes->iter = i+1;
227fbe28522SBarry Smith       snes->norm = fnorm;
2280f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
22985385478SLisandro Dalcin       SNESLogConvHistory(snes,snes->norm,lits);
2307a03ce2fSLisandro Dalcin       ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
23185385478SLisandro Dalcin       /* Test for convergence, xnorm = || X || */
232ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
23385385478SLisandro Dalcin       if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
234e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
2353f149594SLisandro Dalcin       if (reason) break;
2361acabf8cSLois Curfman McInnes     } else {
2371acabf8cSLois Curfman McInnes       break;
2381acabf8cSLois Curfman McInnes     }
23938442cffSBarry Smith   }
24052392280SLois Curfman McInnes   if (i == maxits) {
241ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
24285385478SLisandro Dalcin     if (!reason) reason = SNES_DIVERGED_MAX_IT;
24352392280SLois Curfman McInnes   }
2440f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
245184914b5SBarry Smith   snes->reason = reason;
2460f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2473a40ed3dSBarry Smith   PetscFunctionReturn(0);
2484800dd8cSBarry Smith }
2494800dd8cSBarry Smith /*------------------------------------------------------------*/
2504a2ae208SSatish Balay #undef __FUNCT__
2514b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR"
2526849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes)
2534800dd8cSBarry Smith {
254dfbe8321SBarry Smith   PetscErrorCode ierr;
2553a40ed3dSBarry Smith 
2563a40ed3dSBarry Smith   PetscFunctionBegin;
25758c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
2586cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
2593a40ed3dSBarry Smith   PetscFunctionReturn(0);
2604800dd8cSBarry Smith }
2616b8b9a38SLisandro Dalcin 
2626b8b9a38SLisandro Dalcin #undef __FUNCT__
2636b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_TR"
2646b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_TR(SNES snes)
2656b8b9a38SLisandro Dalcin {
2666b8b9a38SLisandro Dalcin 
2676b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2686b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
2696b8b9a38SLisandro Dalcin }
2706b8b9a38SLisandro Dalcin 
2714a2ae208SSatish Balay #undef __FUNCT__
2724b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR"
2736849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes)
2744800dd8cSBarry Smith {
275dfbe8321SBarry Smith   PetscErrorCode ierr;
2765baf8537SBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
2786b8b9a38SLisandro Dalcin   ierr = SNESReset_TR(snes);CHKERRQ(ierr);
279606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
2814800dd8cSBarry Smith }
2824800dd8cSBarry Smith /*------------------------------------------------------------*/
2834800dd8cSBarry Smith 
2844a2ae208SSatish Balay #undef __FUNCT__
2854b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR"
2866849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
2874800dd8cSBarry Smith {
2884b27c08aSLois Curfman McInnes   SNES_TR *ctx = (SNES_TR *)snes->data;
289dfbe8321SBarry Smith   PetscErrorCode ierr;
2904800dd8cSBarry Smith 
2913a40ed3dSBarry Smith   PetscFunctionBegin;
292b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
29387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
2944b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
2954b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
2964b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
2974b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
2984b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
2994b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
3004b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
301b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3023a40ed3dSBarry Smith   PetscFunctionReturn(0);
3034800dd8cSBarry Smith }
3044800dd8cSBarry Smith 
3054a2ae208SSatish Balay #undef __FUNCT__
3064b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR"
3076849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
308a935fc98SLois Curfman McInnes {
3094b27c08aSLois Curfman McInnes   SNES_TR *tr = (SNES_TR *)snes->data;
310dfbe8321SBarry Smith   PetscErrorCode ierr;
311ace3abfcSBarry Smith   PetscBool  iascii;
312a935fc98SLois Curfman McInnes 
3133a40ed3dSBarry Smith   PetscFunctionBegin;
3142692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
31532077d6dSBarry Smith   if (iascii) {
316a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
317a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
31819bcc07fSBarry Smith   }
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
320a935fc98SLois Curfman McInnes }
32152392280SLois Curfman McInnes /* ------------------------------------------------------------ */
322ebe3b25bSBarry Smith /*MC
323ebe3b25bSBarry Smith       SNESTR - Newton based nonlinear solver that uses a trust region
324ebe3b25bSBarry Smith 
325ebe3b25bSBarry Smith    Options Database:
326ebe3b25bSBarry Smith +    -snes_trtol <tol> Trust region tolerance
327ebe3b25bSBarry Smith .    -snes_tr_mu <mu>
328ebe3b25bSBarry Smith .    -snes_tr_eta <eta>
329ebe3b25bSBarry Smith .    -snes_tr_sigma <sigma>
330ebe3b25bSBarry Smith .    -snes_tr_delta0 <delta0>
331ebe3b25bSBarry Smith .    -snes_tr_delta1 <delta1>
332ebe3b25bSBarry Smith .    -snes_tr_delta2 <delta2>
333ebe3b25bSBarry Smith -    -snes_tr_delta3 <delta3>
334ebe3b25bSBarry Smith 
335ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
336ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
337ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
338ebe3b25bSBarry Smith 
339ebe3b25bSBarry Smith    This is intended as a model implementation, since it does not
340ebe3b25bSBarry Smith    necessarily have many of the bells and whistles of other
341ebe3b25bSBarry Smith    implementations.
342ebe3b25bSBarry Smith 
343ee3001cbSBarry Smith    Level: intermediate
344ee3001cbSBarry Smith 
345ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
346ebe3b25bSBarry Smith 
347ebe3b25bSBarry Smith M*/
348fb2e594dSBarry Smith EXTERN_C_BEGIN
3494a2ae208SSatish Balay #undef __FUNCT__
3504b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR"
3517087cfbeSBarry Smith PetscErrorCode  SNESCreate_TR(SNES snes)
3524800dd8cSBarry Smith {
3534b27c08aSLois Curfman McInnes   SNES_TR        *neP;
354dfbe8321SBarry Smith   PetscErrorCode ierr;
3554800dd8cSBarry Smith 
3563a40ed3dSBarry Smith   PetscFunctionBegin;
357e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_TR;
358e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_TR;
359e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_TR;
360e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_TR;
361e7788613SBarry Smith   snes->ops->view            = SNESView_TR;
3626b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_TR;
363fbe28522SBarry Smith 
36442f4f86dSBarry Smith   snes->usesksp             = PETSC_TRUE;
36542f4f86dSBarry Smith   snes->usespc              = PETSC_FALSE;
36642f4f86dSBarry Smith 
36738f2d2fdSLisandro Dalcin   ierr			= PetscNewLog(snes,SNES_TR,&neP);CHKERRQ(ierr);
368fbe28522SBarry Smith   snes->data	        = (void*)neP;
369fbe28522SBarry Smith   neP->mu		= 0.25;
370fbe28522SBarry Smith   neP->eta		= 0.75;
371fbe28522SBarry Smith   neP->delta		= 0.0;
372fbe28522SBarry Smith   neP->delta0		= 0.2;
373fbe28522SBarry Smith   neP->delta1		= 0.3;
374fbe28522SBarry Smith   neP->delta2		= 0.75;
375fbe28522SBarry Smith   neP->delta3		= 2.0;
376fbe28522SBarry Smith   neP->sigma		= 0.0001;
377ef87ac0dSKris Buschelman   neP->itflag		= PETSC_FALSE;
37875567043SBarry Smith   neP->rnorm0		= 0.0;
37975567043SBarry Smith   neP->ttol		= 0.0;
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
3814800dd8cSBarry Smith }
382fb2e594dSBarry Smith EXTERN_C_END
38382bf6240SBarry Smith 
384