xref: /petsc/src/snes/impls/tr/tr.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
14800dd8cSBarry Smith 
2e090d566SSatish Balay #include "src/snes/impls/tr/tr.h"                /*I   "petscsnes.h"   I*/
34800dd8cSBarry Smith 
4fbe28522SBarry Smith /*
5df60cc22SBarry Smith    This convergence test determines if the two norm of the
6df60cc22SBarry Smith    solution lies outside the trust region, if so it halts.
7df60cc22SBarry Smith */
84a2ae208SSatish Balay #undef __FUNCT__
94b27c08aSLois Curfman McInnes #define __FUNCT__ "SNES_TR_KSPConverged_Private"
10dfbe8321SBarry Smith PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,int n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
11df60cc22SBarry Smith {
12df60cc22SBarry Smith   SNES                snes = (SNES) ctx;
1398120f82SLois Curfman McInnes   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
144b27c08aSLois Curfman McInnes   SNES_TR             *neP = (SNES_TR*)snes->data;
15df60cc22SBarry Smith   Vec                 x;
16064f8208SBarry Smith   PetscReal           nrm;
17dfbe8321SBarry Smith   PetscErrorCode ierr;
18df60cc22SBarry Smith 
193a40ed3dSBarry Smith   PetscFunctionBegin;
2098120f82SLois Curfman McInnes   if (snes->ksp_ewconv) {
2129bbc08cSBarry Smith     if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker onvergence context not created");
22c50d6201SSatish Balay     if (!n) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);}
2398120f82SLois Curfman McInnes     kctx->lresid_last = rnorm;
24df60cc22SBarry Smith   }
25329f5518SBarry Smith   ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr);
26329f5518SBarry Smith   if (*reason) {
274b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: regular convergence test KSP iterations=%d, rnorm=%g\n",n,rnorm);
289b04593eSLois Curfman McInnes   }
29df60cc22SBarry Smith 
30a935fc98SLois Curfman McInnes   /* Determine norm of solution */
3178b31e54SBarry Smith   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
32064f8208SBarry Smith   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
33064f8208SBarry Smith   if (nrm >= neP->delta) {
344b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm);
354b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,nrm);
36329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
37df60cc22SBarry Smith   }
383a40ed3dSBarry Smith   PetscFunctionReturn(0);
39df60cc22SBarry Smith }
4082bf6240SBarry Smith 
41df60cc22SBarry Smith /*
424b27c08aSLois Curfman McInnes    SNESSolve_TR - Implements Newton's Method with a very simple trust
43ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
444800dd8cSBarry Smith 
45ddbbdb52SLois Curfman McInnes 
464800dd8cSBarry Smith */
474a2ae208SSatish Balay #undef __FUNCT__
484b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR"
49*6849ba73SBarry Smith static PetscErrorCode SNESSolve_TR(SNES snes)
504800dd8cSBarry Smith {
514b27c08aSLois Curfman McInnes   SNES_TR             *neP = (SNES_TR*)snes->data;
526b5873e3SBarry Smith   Vec                 X,F,Y,G,TMP,Ytmp;
53*6849ba73SBarry Smith   PetscErrorCode ierr;
54*6849ba73SBarry Smith   int                 maxits,i,lits,breakout = 0;
55112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
56064f8208SBarry Smith   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1;
57ea709b57SSatish Balay   PetscScalar         mone = -1.0,cnorm;
58df60cc22SBarry Smith   KSP                 ksp;
59184914b5SBarry Smith   SNESConvergedReason reason;
607c4af3dcSLois Curfman McInnes   PetscTruth          conv;
614800dd8cSBarry Smith 
623a40ed3dSBarry Smith   PetscFunctionBegin;
63fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
64fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
6539e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
66fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
67fbe28522SBarry Smith   G		= snes->work[1];
686b5873e3SBarry Smith   Ytmp          = snes->work[2];
694800dd8cSBarry Smith 
704c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
714c49b128SBarry Smith   snes->iter = 0;
724c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
737e6560a3SBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
744800dd8cSBarry Smith 
755334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
76cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
770f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
78fbe28522SBarry Smith   snes->norm = fnorm;
790f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
804800dd8cSBarry Smith   delta = neP->delta0*fnorm;
814800dd8cSBarry Smith   neP->delta = delta;
820462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
8394a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
8494b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
854800dd8cSBarry Smith 
86c9780b6fSBarry Smith  if (fnorm < snes->atol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
87b37302c6SLois Curfman McInnes 
88d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
89d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
90d034289bSBarry Smith 
91a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
927c4af3dcSLois Curfman McInnes   ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr);
937c4af3dcSLois Curfman McInnes   if (!conv) {
944b27c08aSLois Curfman McInnes     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr);
954b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_TR: Using Krylov convergence test SNES_TR_KSPConverged_Private\n");
967c4af3dcSLois Curfman McInnes   }
97df60cc22SBarry Smith 
984800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
995334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
10094b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1012deea200SLois Curfman McInnes 
1022deea200SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
10323ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
104c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
105329f5518SBarry Smith     snes->linear_its += lits;
1064b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
107064f8208SBarry Smith     ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
108064f8208SBarry Smith     norm1 = nrm;
1096b5873e3SBarry Smith     while(1) {
110393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
111064f8208SBarry Smith       nrm = norm1;
1126b5873e3SBarry Smith 
1132deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
114064f8208SBarry Smith       if (nrm >= delta) {
115064f8208SBarry Smith         nrm = delta/nrm;
116064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
117064f8208SBarry Smith         cnorm = nrm;
1184b27c08aSLois Curfman McInnes         PetscLogInfo(snes,"SNESSolve_TR: Scaling direction by %g\n",nrm);
119393d2d9aSLois Curfman McInnes         ierr = VecScale(&cnorm,Y);CHKERRQ(ierr);
120064f8208SBarry Smith         nrm = gpnorm;
121fbe28522SBarry Smith         ynorm = delta;
122fbe28522SBarry Smith       } else {
123fbe28522SBarry Smith         gpnorm = 0.0;
1244b27c08aSLois Curfman McInnes         PetscLogInfo(snes,"SNESSolve_TR: Direction is in Trust Region\n");
125064f8208SBarry Smith         ynorm = nrm;
126fbe28522SBarry Smith       }
1272deea200SLois Curfman McInnes       ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr);            /* Y <- X - Y */
12881b6cf68SLois Curfman McInnes       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
1295334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
130cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1314800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
132fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1334800dd8cSBarry Smith 
1344800dd8cSBarry Smith       /* Update size of trust region */
1354800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1364800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1374800dd8cSBarry Smith       else                     delta *= neP->delta3;
1384b27c08aSLois Curfman McInnes       PetscLogInfo(snes,"SNESSolve_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm);
1394b27c08aSLois Curfman McInnes       PetscLogInfo(snes,"SNESSolve_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta);
1404800dd8cSBarry Smith       neP->delta = delta;
1414800dd8cSBarry Smith       if (rho > neP->sigma) break;
1424b27c08aSLois Curfman McInnes       PetscLogInfo(snes,"SNESSolve_TR: Trying again in smaller region\n");
1436b5873e3SBarry Smith       /* check to see if progress is hopeless */
14452392280SLois Curfman McInnes       neP->itflag = 0;
145184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
146184914b5SBarry Smith       if (reason) {
14752392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
1485ed2d596SBarry Smith         SNESMonitor(snes,i+1,fnorm);
149454a90a3SBarry Smith         breakout = 1;
150454a90a3SBarry Smith         break;
15152392280SLois Curfman McInnes       }
15250ffb88aSMatthew Knepley       snes->numFailures++;
1534800dd8cSBarry Smith     }
1541acabf8cSLois Curfman McInnes     if (!breakout) {
1554800dd8cSBarry Smith       fnorm = gnorm;
1560f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
157c509a162SBarry Smith       snes->iter = i+1;
158fbe28522SBarry Smith       snes->norm = fnorm;
1590f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16039e2f89bSBarry Smith       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
16139e2f89bSBarry Smith       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
162329f5518SBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);		/* xnorm = || X || */
1630462333dSBarry Smith       SNESLogConvHistory(snes,fnorm,lits);
16494a424c1SBarry Smith       SNESMonitor(snes,i+1,fnorm);
1654800dd8cSBarry Smith 
1664800dd8cSBarry Smith       /* Test for convergence */
16752392280SLois Curfman McInnes       neP->itflag = 1;
168184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
169184914b5SBarry Smith       if (reason) {
17038442cffSBarry Smith         break;
17138442cffSBarry Smith       }
1721acabf8cSLois Curfman McInnes     } else {
1731acabf8cSLois Curfman McInnes       break;
1741acabf8cSLois Curfman McInnes     }
17538442cffSBarry Smith   }
17638442cffSBarry Smith   /* Verify solution is in corect location */
177e15f7bb5SBarry Smith   if (X != snes->vec_sol) {
178393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
179e15f7bb5SBarry Smith   }
180e15f7bb5SBarry Smith   if (F != snes->vec_func) {
181e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
182e15f7bb5SBarry Smith   }
18339e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
18439e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
18552392280SLois Curfman McInnes   if (i == maxits) {
1864b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_TR: Maximum number of iterations has been reached: %d\n",maxits);
187184914b5SBarry Smith     reason = SNES_DIVERGED_MAX_IT;
18852392280SLois Curfman McInnes   }
1890f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
190184914b5SBarry Smith   snes->reason = reason;
1910f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
1934800dd8cSBarry Smith }
1944800dd8cSBarry Smith /*------------------------------------------------------------*/
1954a2ae208SSatish Balay #undef __FUNCT__
1964b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR"
197*6849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes)
1984800dd8cSBarry Smith {
199dfbe8321SBarry Smith   PetscErrorCode ierr;
2003a40ed3dSBarry Smith 
2013a40ed3dSBarry Smith   PetscFunctionBegin;
20281b6cf68SLois Curfman McInnes   snes->nwork = 4;
203d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
204b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
20581b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2063a40ed3dSBarry Smith   PetscFunctionReturn(0);
2074800dd8cSBarry Smith }
2084800dd8cSBarry Smith /*------------------------------------------------------------*/
2094a2ae208SSatish Balay #undef __FUNCT__
2104b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR"
211*6849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes)
2124800dd8cSBarry Smith {
213dfbe8321SBarry Smith   PetscErrorCode ierr;
2145baf8537SBarry Smith 
2153a40ed3dSBarry Smith   PetscFunctionBegin;
2165baf8537SBarry Smith   if (snes->nwork) {
2174b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
2185baf8537SBarry Smith   }
219606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2203a40ed3dSBarry Smith   PetscFunctionReturn(0);
2214800dd8cSBarry Smith }
2224800dd8cSBarry Smith /*------------------------------------------------------------*/
2234800dd8cSBarry Smith 
2244a2ae208SSatish Balay #undef __FUNCT__
2254b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR"
226*6849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
2274800dd8cSBarry Smith {
2284b27c08aSLois Curfman McInnes   SNES_TR *ctx = (SNES_TR *)snes->data;
229dfbe8321SBarry Smith   PetscErrorCode ierr;
2304800dd8cSBarry Smith 
2313a40ed3dSBarry Smith   PetscFunctionBegin;
232b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
23387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
2344b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
2354b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
2364b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
2374b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
2384b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
2394b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
2404b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
241b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2423a40ed3dSBarry Smith   PetscFunctionReturn(0);
2434800dd8cSBarry Smith }
2444800dd8cSBarry Smith 
2454a2ae208SSatish Balay #undef __FUNCT__
2464b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR"
247*6849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
248a935fc98SLois Curfman McInnes {
2494b27c08aSLois Curfman McInnes   SNES_TR *tr = (SNES_TR *)snes->data;
250dfbe8321SBarry Smith   PetscErrorCode ierr;
25132077d6dSBarry Smith   PetscTruth iascii;
252a935fc98SLois Curfman McInnes 
2533a40ed3dSBarry Smith   PetscFunctionBegin;
25432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
25532077d6dSBarry Smith   if (iascii) {
256b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
257b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
2585cd90555SBarry Smith   } else {
25929bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
26019bcc07fSBarry Smith   }
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
262a935fc98SLois Curfman McInnes }
263a935fc98SLois Curfman McInnes 
26452392280SLois Curfman McInnes /* ---------------------------------------------------------------- */
2654a2ae208SSatish Balay #undef __FUNCT__
2664b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR"
2678ac28fe4SSatish Balay /*@C
2684b27c08aSLois Curfman McInnes    SNESConverged_TR - Monitors the convergence of the trust region
2694b27c08aSLois Curfman McInnes    method SNESTR for solving systems of nonlinear equations (default).
27052392280SLois Curfman McInnes 
271c7afd0dbSLois Curfman McInnes    Collective on SNES
272c7afd0dbSLois Curfman McInnes 
27352392280SLois Curfman McInnes    Input Parameters:
274c7afd0dbSLois Curfman McInnes +  snes - the SNES context
27552392280SLois Curfman McInnes .  xnorm - 2-norm of current iterate
27652392280SLois Curfman McInnes .  pnorm - 2-norm of current step
27752392280SLois Curfman McInnes .  fnorm - 2-norm of function
278c7afd0dbSLois Curfman McInnes -  dummy - unused context
279fee21e36SBarry Smith 
280184914b5SBarry Smith    Output Parameter:
281184914b5SBarry Smith .   reason - one of
282184914b5SBarry Smith $  SNES_CONVERGED_FNORM_ABS       - (fnorm < atol),
283184914b5SBarry Smith $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
284184914b5SBarry Smith $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
285184914b5SBarry Smith $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
286184914b5SBarry Smith $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
287184914b5SBarry Smith $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
288184914b5SBarry Smith $  SNES_CONVERGED_ITERATING       - (otherwise)
28952392280SLois Curfman McInnes 
29052392280SLois Curfman McInnes    where
291c7afd0dbSLois Curfman McInnes +    delta    - trust region paramenter
292c7afd0dbSLois Curfman McInnes .    deltatol - trust region size tolerance,
293c7afd0dbSLois Curfman McInnes                 set with SNESSetTrustRegionTolerance()
294c7afd0dbSLois Curfman McInnes .    maxf - maximum number of function evaluations,
295c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
296c7afd0dbSLois Curfman McInnes .    nfct - number of function evaluations,
297c7afd0dbSLois Curfman McInnes .    atol - absolute function norm tolerance,
298c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
299c7afd0dbSLois Curfman McInnes -    xtol - relative function norm tolerance,
300c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
30152392280SLois Curfman McInnes 
30215091d37SBarry Smith    Level: intermediate
30315091d37SBarry Smith 
30452392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
30552392280SLois Curfman McInnes 
30652392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
30752392280SLois Curfman McInnes @*/
308dfbe8321SBarry Smith PetscErrorCode SNESConverged_TR(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
30952392280SLois Curfman McInnes {
3104b27c08aSLois Curfman McInnes   SNES_TR *neP = (SNES_TR *)snes->data;
311dfbe8321SBarry Smith   PetscErrorCode ierr;
312ddd12767SBarry Smith 
3133a40ed3dSBarry Smith   PetscFunctionBegin;
314d252947aSBarry Smith   if (fnorm != fnorm) {
3154b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESConverged_TR:Failed to converged, function norm is NaN\n");
316184914b5SBarry Smith     *reason = SNES_DIVERGED_FNORM_NAN;
317184914b5SBarry Smith   } else if (neP->delta < xnorm * snes->deltatol) {
3184b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESConverged_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol);
319184914b5SBarry Smith     *reason = SNES_CONVERGED_TR_DELTA;
320184914b5SBarry Smith   } else if (neP->itflag) {
3214b27c08aSLois Curfman McInnes     ierr = SNESConverged_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
3221acabf8cSLois Curfman McInnes   } else if (snes->nfuncs > snes->max_funcs) {
3234b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESConverged_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs);
324184914b5SBarry Smith     *reason = SNES_DIVERGED_FUNCTION_COUNT;
325184914b5SBarry Smith   } else {
326184914b5SBarry Smith     *reason = SNES_CONVERGED_ITERATING;
32752392280SLois Curfman McInnes   }
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
32952392280SLois Curfman McInnes }
33052392280SLois Curfman McInnes /* ------------------------------------------------------------ */
331ebe3b25bSBarry Smith /*MC
332ebe3b25bSBarry Smith       SNESTR - Newton based nonlinear solver that uses a trust region
333ebe3b25bSBarry Smith 
334ebe3b25bSBarry Smith    Options Database:
335ebe3b25bSBarry Smith +    -snes_trtol <tol> Trust region tolerance
336ebe3b25bSBarry Smith .    -snes_tr_mu <mu>
337ebe3b25bSBarry Smith .    -snes_tr_eta <eta>
338ebe3b25bSBarry Smith .    -snes_tr_sigma <sigma>
339ebe3b25bSBarry Smith .    -snes_tr_delta0 <delta0>
340ebe3b25bSBarry Smith .    -snes_tr_delta1 <delta1>
341ebe3b25bSBarry Smith .    -snes_tr_delta2 <delta2>
342ebe3b25bSBarry Smith -    -snes_tr_delta3 <delta3>
343ebe3b25bSBarry Smith 
344ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
345ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
346ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
347ebe3b25bSBarry Smith 
348ebe3b25bSBarry Smith    This is intended as a model implementation, since it does not
349ebe3b25bSBarry Smith    necessarily have many of the bells and whistles of other
350ebe3b25bSBarry Smith    implementations.
351ebe3b25bSBarry Smith 
352ee3001cbSBarry Smith    Level: intermediate
353ee3001cbSBarry Smith 
354ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
355ebe3b25bSBarry Smith 
356ebe3b25bSBarry Smith M*/
357fb2e594dSBarry Smith EXTERN_C_BEGIN
3584a2ae208SSatish Balay #undef __FUNCT__
3594b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR"
360dfbe8321SBarry Smith PetscErrorCode SNESCreate_TR(SNES snes)
3614800dd8cSBarry Smith {
3624b27c08aSLois Curfman McInnes   SNES_TR *neP;
363dfbe8321SBarry Smith   PetscErrorCode ierr;
3644800dd8cSBarry Smith 
3653a40ed3dSBarry Smith   PetscFunctionBegin;
3664b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_TR;
3674b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_TR;
3684b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_TR;
3694b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_TR;
3704b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_TR;
3714b27c08aSLois Curfman McInnes   snes->view            = SNESView_TR;
3725baf8537SBarry Smith   snes->nwork           = 0;
373fbe28522SBarry Smith 
3744b27c08aSLois Curfman McInnes   ierr			= PetscNew(SNES_TR,&neP);CHKERRQ(ierr);
3754b27c08aSLois Curfman McInnes   PetscLogObjectMemory(snes,sizeof(SNES_TR));
376fbe28522SBarry Smith   snes->data	        = (void*)neP;
377fbe28522SBarry Smith   neP->mu		= 0.25;
378fbe28522SBarry Smith   neP->eta		= 0.75;
379fbe28522SBarry Smith   neP->delta		= 0.0;
380fbe28522SBarry Smith   neP->delta0		= 0.2;
381fbe28522SBarry Smith   neP->delta1		= 0.3;
382fbe28522SBarry Smith   neP->delta2		= 0.75;
383fbe28522SBarry Smith   neP->delta3		= 2.0;
384fbe28522SBarry Smith   neP->sigma		= 0.0001;
385fbe28522SBarry Smith   neP->itflag		= 0;
38652392280SLois Curfman McInnes   neP->rnorm0		= 0;
38752392280SLois Curfman McInnes   neP->ttol		= 0;
3883a40ed3dSBarry Smith   PetscFunctionReturn(0);
3894800dd8cSBarry Smith }
390fb2e594dSBarry Smith EXTERN_C_END
39182bf6240SBarry Smith 
392