xref: /petsc/src/snes/impls/tr/tr.c (revision 6cab3a1b8d4fca83e1d7b61474c608460de73d5f)
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;
984800dd8cSBarry Smith 
993a40ed3dSBarry Smith   PetscFunctionBegin;
100fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
101fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
10239e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
103fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
104fbe28522SBarry Smith   G		= snes->work[1];
1056b5873e3SBarry Smith   Ytmp          = snes->work[2];
1064800dd8cSBarry Smith 
1074c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1084c49b128SBarry Smith   snes->iter = 0;
1094c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1104800dd8cSBarry Smith 
1115334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
112cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
1130f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
114fbe28522SBarry Smith   snes->norm = fnorm;
1150f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1164800dd8cSBarry Smith   delta = neP->delta0*fnorm;
1174800dd8cSBarry Smith   neP->delta = delta;
1180462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
1197a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
120b37302c6SLois Curfman McInnes 
121d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
122d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
12385385478SLisandro Dalcin   /* test convergence */
12485385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
12585385478SLisandro Dalcin   if (snes->reason) PetscFunctionReturn(0);
1263f149594SLisandro Dalcin 
127a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
128acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv,PETSC_NULL);CHKERRQ(ierr);
1297c4af3dcSLois Curfman McInnes   if (!conv) {
130971273eeSBarry Smith     SNES_TR_KSPConverged_Ctx *ctx;
1313f149594SLisandro Dalcin     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
132971273eeSBarry Smith     ierr = PetscNew(SNES_TR_KSPConverged_Ctx,&ctx);CHKERRQ(ierr);
133971273eeSBarry Smith     ctx->snes = snes;
134971273eeSBarry Smith     ierr = KSPDefaultConvergedCreate(&ctx->ctx);CHKERRQ(ierr);
135971273eeSBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr);
136ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
1377c4af3dcSLois Curfman McInnes   }
138df60cc22SBarry Smith 
1394800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
14099a96b7cSMatthew Knepley 
14199a96b7cSMatthew Knepley     /* Call general purpose update function */
142e7788613SBarry Smith     if (snes->ops->update) {
143e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
14499a96b7cSMatthew Knepley     }
14599a96b7cSMatthew Knepley 
14685385478SLisandro Dalcin     /* Solve J Y = F, where J is Jacobian matrix */
1475334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
14894b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1493f149594SLisandro Dalcin     ierr = SNES_KSPSolve(snes,snes->ksp,F,Ytmp);CHKERRQ(ierr);
1503f149594SLisandro Dalcin     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
151329f5518SBarry Smith     snes->linear_its += lits;
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;
164ae15b995SBarry Smith         ierr = PetscInfo1(snes,"Scaling direction by %G\n",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       }
17379f36460SBarry Smith       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
17485385478SLisandro Dalcin       ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr);
1755334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
176cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1774800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
178fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1794800dd8cSBarry Smith 
1804800dd8cSBarry Smith       /* Update size of trust region */
1814800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1824800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1834800dd8cSBarry Smith       else                     delta *= neP->delta3;
184ae15b995SBarry Smith       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
185ae15b995SBarry Smith       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
1864800dd8cSBarry Smith       neP->delta = delta;
1874800dd8cSBarry Smith       if (rho > neP->sigma) break;
188ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
1896b5873e3SBarry Smith       /* check to see if progress is hopeless */
190ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
19185385478SLisandro Dalcin       ierr = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
19285385478SLisandro Dalcin       if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); }
193184914b5SBarry Smith       if (reason) {
19452392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
1957a03ce2fSLisandro Dalcin         ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
196a7cc72afSBarry Smith         breakout = PETSC_TRUE;
197454a90a3SBarry Smith         break;
19852392280SLois Curfman McInnes       }
19950ffb88aSMatthew Knepley       snes->numFailures++;
2004800dd8cSBarry Smith     }
2011acabf8cSLois Curfman McInnes     if (!breakout) {
20285385478SLisandro Dalcin       /* Update function and solution vectors */
2034800dd8cSBarry Smith       fnorm = gnorm;
20485385478SLisandro Dalcin       ierr = VecCopy(G,F);CHKERRQ(ierr);
20585385478SLisandro Dalcin       ierr = VecCopy(Y,X);CHKERRQ(ierr);
20685385478SLisandro Dalcin       /* Monitor convergence */
2070f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
208c509a162SBarry Smith       snes->iter = i+1;
209fbe28522SBarry Smith       snes->norm = fnorm;
2100f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
21185385478SLisandro Dalcin       SNESLogConvHistory(snes,snes->norm,lits);
2127a03ce2fSLisandro Dalcin       ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
21385385478SLisandro Dalcin       /* Test for convergence, xnorm = || X || */
214ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
21585385478SLisandro Dalcin       if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
216e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
2173f149594SLisandro Dalcin       if (reason) break;
2181acabf8cSLois Curfman McInnes     } else {
2191acabf8cSLois Curfman McInnes       break;
2201acabf8cSLois Curfman McInnes     }
22138442cffSBarry Smith   }
22252392280SLois Curfman McInnes   if (i == maxits) {
223ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
22485385478SLisandro Dalcin     if (!reason) reason = SNES_DIVERGED_MAX_IT;
22552392280SLois Curfman McInnes   }
2260f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
227184914b5SBarry Smith   snes->reason = reason;
2280f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2293a40ed3dSBarry Smith   PetscFunctionReturn(0);
2304800dd8cSBarry Smith }
2314800dd8cSBarry Smith /*------------------------------------------------------------*/
2324a2ae208SSatish Balay #undef __FUNCT__
2334b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR"
2346849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes)
2354800dd8cSBarry Smith {
236dfbe8321SBarry Smith   PetscErrorCode ierr;
2373a40ed3dSBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
23958c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
240*6cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
2413a40ed3dSBarry Smith   PetscFunctionReturn(0);
2424800dd8cSBarry Smith }
2436b8b9a38SLisandro Dalcin 
2446b8b9a38SLisandro Dalcin #undef __FUNCT__
2456b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_TR"
2466b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_TR(SNES snes)
2476b8b9a38SLisandro Dalcin {
2486b8b9a38SLisandro Dalcin 
2496b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2506b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
2516b8b9a38SLisandro Dalcin }
2526b8b9a38SLisandro Dalcin 
2534a2ae208SSatish Balay #undef __FUNCT__
2544b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR"
2556849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes)
2564800dd8cSBarry Smith {
257dfbe8321SBarry Smith   PetscErrorCode ierr;
2585baf8537SBarry Smith 
2593a40ed3dSBarry Smith   PetscFunctionBegin;
2606b8b9a38SLisandro Dalcin   ierr = SNESReset_TR(snes);CHKERRQ(ierr);
261606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2623a40ed3dSBarry Smith   PetscFunctionReturn(0);
2634800dd8cSBarry Smith }
2644800dd8cSBarry Smith /*------------------------------------------------------------*/
2654800dd8cSBarry Smith 
2664a2ae208SSatish Balay #undef __FUNCT__
2674b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR"
2686849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
2694800dd8cSBarry Smith {
2704b27c08aSLois Curfman McInnes   SNES_TR *ctx = (SNES_TR *)snes->data;
271dfbe8321SBarry Smith   PetscErrorCode ierr;
2724800dd8cSBarry Smith 
2733a40ed3dSBarry Smith   PetscFunctionBegin;
274b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
27587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
2764b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
2774b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
2784b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
2794b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
2804b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
2814b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
2824b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
283b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2843a40ed3dSBarry Smith   PetscFunctionReturn(0);
2854800dd8cSBarry Smith }
2864800dd8cSBarry Smith 
2874a2ae208SSatish Balay #undef __FUNCT__
2884b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR"
2896849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
290a935fc98SLois Curfman McInnes {
2914b27c08aSLois Curfman McInnes   SNES_TR *tr = (SNES_TR *)snes->data;
292dfbe8321SBarry Smith   PetscErrorCode ierr;
293ace3abfcSBarry Smith   PetscBool  iascii;
294a935fc98SLois Curfman McInnes 
2953a40ed3dSBarry Smith   PetscFunctionBegin;
2962692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29732077d6dSBarry Smith   if (iascii) {
298a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
299a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
30019bcc07fSBarry Smith   }
3013a40ed3dSBarry Smith   PetscFunctionReturn(0);
302a935fc98SLois Curfman McInnes }
30352392280SLois Curfman McInnes /* ------------------------------------------------------------ */
304ebe3b25bSBarry Smith /*MC
305ebe3b25bSBarry Smith       SNESTR - Newton based nonlinear solver that uses a trust region
306ebe3b25bSBarry Smith 
307ebe3b25bSBarry Smith    Options Database:
308ebe3b25bSBarry Smith +    -snes_trtol <tol> Trust region tolerance
309ebe3b25bSBarry Smith .    -snes_tr_mu <mu>
310ebe3b25bSBarry Smith .    -snes_tr_eta <eta>
311ebe3b25bSBarry Smith .    -snes_tr_sigma <sigma>
312ebe3b25bSBarry Smith .    -snes_tr_delta0 <delta0>
313ebe3b25bSBarry Smith .    -snes_tr_delta1 <delta1>
314ebe3b25bSBarry Smith .    -snes_tr_delta2 <delta2>
315ebe3b25bSBarry Smith -    -snes_tr_delta3 <delta3>
316ebe3b25bSBarry Smith 
317ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
318ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
319ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
320ebe3b25bSBarry Smith 
321ebe3b25bSBarry Smith    This is intended as a model implementation, since it does not
322ebe3b25bSBarry Smith    necessarily have many of the bells and whistles of other
323ebe3b25bSBarry Smith    implementations.
324ebe3b25bSBarry Smith 
325ee3001cbSBarry Smith    Level: intermediate
326ee3001cbSBarry Smith 
327ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
328ebe3b25bSBarry Smith 
329ebe3b25bSBarry Smith M*/
330fb2e594dSBarry Smith EXTERN_C_BEGIN
3314a2ae208SSatish Balay #undef __FUNCT__
3324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR"
3337087cfbeSBarry Smith PetscErrorCode  SNESCreate_TR(SNES snes)
3344800dd8cSBarry Smith {
3354b27c08aSLois Curfman McInnes   SNES_TR        *neP;
336dfbe8321SBarry Smith   PetscErrorCode ierr;
3374800dd8cSBarry Smith 
3383a40ed3dSBarry Smith   PetscFunctionBegin;
339e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_TR;
340e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_TR;
341e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_TR;
342e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_TR;
343e7788613SBarry Smith   snes->ops->view            = SNESView_TR;
3446b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_TR;
345fbe28522SBarry Smith 
34642f4f86dSBarry Smith   snes->usesksp             = PETSC_TRUE;
34742f4f86dSBarry Smith   snes->usespc              = PETSC_FALSE;
34842f4f86dSBarry Smith 
34938f2d2fdSLisandro Dalcin   ierr			= PetscNewLog(snes,SNES_TR,&neP);CHKERRQ(ierr);
350fbe28522SBarry Smith   snes->data	        = (void*)neP;
351fbe28522SBarry Smith   neP->mu		= 0.25;
352fbe28522SBarry Smith   neP->eta		= 0.75;
353fbe28522SBarry Smith   neP->delta		= 0.0;
354fbe28522SBarry Smith   neP->delta0		= 0.2;
355fbe28522SBarry Smith   neP->delta1		= 0.3;
356fbe28522SBarry Smith   neP->delta2		= 0.75;
357fbe28522SBarry Smith   neP->delta3		= 2.0;
358fbe28522SBarry Smith   neP->sigma		= 0.0001;
359ef87ac0dSKris Buschelman   neP->itflag		= PETSC_FALSE;
36075567043SBarry Smith   neP->rnorm0		= 0.0;
36175567043SBarry Smith   neP->ttol		= 0.0;
3623a40ed3dSBarry Smith   PetscFunctionReturn(0);
3634800dd8cSBarry Smith }
364fb2e594dSBarry Smith EXTERN_C_END
36582bf6240SBarry Smith 
366