xref: /petsc/src/snes/impls/tr/tr.c (revision 7c4af3dc9754d1d13541fb941ee023709650f9f9)
1*7c4af3dcSLois Curfman McInnes /*$Id: tr.c,v 1.123 2001/03/23 23:24:15 balay Exp curfman $*/
24800dd8cSBarry Smith 
3e090d566SSatish Balay #include "src/snes/impls/tr/tr.h"                /*I   "petscsnes.h"   I*/
44800dd8cSBarry Smith 
5fbe28522SBarry Smith /*
6df60cc22SBarry Smith    This convergence test determines if the two norm of the
7df60cc22SBarry Smith    solution lies outside the trust region, if so it halts.
8df60cc22SBarry Smith */
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "SNES_EQ_TR_KSPConverged_Private"
113acb151aSSatish Balay int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n,double rnorm,KSPConvergedReason *reason,void *ctx)
12df60cc22SBarry Smith {
13df60cc22SBarry Smith   SNES                snes = (SNES) ctx;
1498120f82SLois Curfman McInnes   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
156831982aSBarry Smith   SNES_EQ_TR          *neP = (SNES_EQ_TR*)snes->data;
16df60cc22SBarry Smith   Vec                 x;
1798120f82SLois Curfman McInnes   double              norm;
183acb151aSSatish Balay   int                 ierr;
19df60cc22SBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
2198120f82SLois Curfman McInnes   if (snes->ksp_ewconv) {
2229bbc08cSBarry Smith     if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker onvergence context not created");
23c50d6201SSatish Balay     if (!n) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);}
2498120f82SLois Curfman McInnes     kctx->lresid_last = rnorm;
25df60cc22SBarry Smith   }
26329f5518SBarry Smith   ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr);
27329f5518SBarry Smith   if (*reason) {
28b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: regular convergence test KSP iterations=%d, rnorm=%g\n",n,rnorm);
299b04593eSLois Curfman McInnes   }
30df60cc22SBarry Smith 
31a935fc98SLois Curfman McInnes   /* Determine norm of solution */
3278b31e54SBarry Smith   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
33cddf8d76SBarry Smith   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
34df60cc22SBarry Smith   if (norm >= neP->delta) {
35b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm);
36b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,norm);
37329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
38df60cc22SBarry Smith   }
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
40df60cc22SBarry Smith }
4182bf6240SBarry Smith 
42df60cc22SBarry Smith /*
43f63b844aSLois Curfman McInnes    SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust
44ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
454800dd8cSBarry Smith 
464800dd8cSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
474800dd8cSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
48ddbbdb52SLois Curfman McInnes    of Mathematical Software", Wayne Cowell, editor.
49ddbbdb52SLois Curfman McInnes 
504800dd8cSBarry Smith    This is intended as a model implementation, since it does not
514800dd8cSBarry Smith    necessarily have many of the bells and whistles of other
524800dd8cSBarry Smith    implementations.
534800dd8cSBarry Smith */
544a2ae208SSatish Balay #undef __FUNCT__
554a2ae208SSatish Balay #define __FUNCT__ "SNESSolve_EQ_TR"
56f63b844aSLois Curfman McInnes static int SNESSolve_EQ_TR(SNES snes,int *its)
574800dd8cSBarry Smith {
586831982aSBarry Smith   SNES_EQ_TR          *neP = (SNES_EQ_TR*)snes->data;
596b5873e3SBarry Smith   Vec                 X,F,Y,G,TMP,Ytmp;
600462333dSBarry Smith   int                 maxits,i,ierr,lits,breakout = 0;
61112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
620462333dSBarry Smith   double              rho,fnorm,gnorm,gpnorm,xnorm,delta,norm,ynorm,norm1;
635334005bSBarry Smith   Scalar              mone = -1.0,cnorm;
64df60cc22SBarry Smith   KSP                 ksp;
65df60cc22SBarry Smith   SLES                sles;
66184914b5SBarry Smith   SNESConvergedReason reason;
67*7c4af3dcSLois Curfman McInnes   PetscTruth          conv;
684800dd8cSBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionBegin;
70fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
71fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
7239e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
73fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
74fbe28522SBarry Smith   G		= snes->work[1];
756b5873e3SBarry Smith   Ytmp          = snes->work[2];
764800dd8cSBarry Smith 
774c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
784c49b128SBarry Smith   snes->iter = 0;
794c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
807e6560a3SBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
814800dd8cSBarry Smith 
825334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
83cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
840f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
85fbe28522SBarry Smith   snes->norm = fnorm;
860f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
874800dd8cSBarry Smith   delta = neP->delta0*fnorm;
884800dd8cSBarry Smith   neP->delta = delta;
890462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
9094a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
914800dd8cSBarry Smith 
92184914b5SBarry Smith  if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
93b37302c6SLois Curfman McInnes 
94d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
95d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
96d034289bSBarry Smith 
97a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
98*7c4af3dcSLois Curfman McInnes   ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr);
99*7c4af3dcSLois Curfman McInnes   if (!conv) {
10078b31e54SBarry Smith     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
10178b31e54SBarry Smith     ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1026831982aSBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr);
103*7c4af3dcSLois Curfman McInnes   }
104df60cc22SBarry Smith 
1054800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
1065334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1075334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1082deea200SLois Curfman McInnes 
1092deea200SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
11078b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr);
111329f5518SBarry Smith     snes->linear_its += lits;
112b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
113cddf8d76SBarry Smith     ierr = VecNorm(Ytmp,NORM_2,&norm);CHKERRQ(ierr);
1142deea200SLois Curfman McInnes     norm1 = norm;
1156b5873e3SBarry Smith     while(1) {
116393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
1172deea200SLois Curfman McInnes       norm = norm1;
1186b5873e3SBarry Smith 
1192deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
120eafb4bcbSBarry Smith       if (norm >= delta) {
121fbe28522SBarry Smith         norm = delta/norm;
122fbe28522SBarry Smith         gpnorm = (1.0 - norm)*fnorm;
123eafb4bcbSBarry Smith         cnorm = norm;
124b0a32e0cSBarry Smith         PetscLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",norm);
125393d2d9aSLois Curfman McInnes         ierr = VecScale(&cnorm,Y);CHKERRQ(ierr);
126eafb4bcbSBarry Smith         norm = gpnorm;
127fbe28522SBarry Smith         ynorm = delta;
128fbe28522SBarry Smith       } else {
129fbe28522SBarry Smith         gpnorm = 0.0;
130b0a32e0cSBarry Smith         PetscLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n");
131fbe28522SBarry Smith         ynorm = norm;
132fbe28522SBarry Smith       }
1332deea200SLois Curfman McInnes       ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr);            /* Y <- X - Y */
13481b6cf68SLois Curfman McInnes       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
1355334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
136cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1374800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
138fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1394800dd8cSBarry Smith 
1404800dd8cSBarry Smith       /* Update size of trust region */
1414800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1424800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1434800dd8cSBarry Smith       else                     delta *= neP->delta3;
144b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm);
145b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta);
1464800dd8cSBarry Smith       neP->delta = delta;
1474800dd8cSBarry Smith       if (rho > neP->sigma) break;
148b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n");
1496b5873e3SBarry Smith       /* check to see if progress is hopeless */
15052392280SLois Curfman McInnes       neP->itflag = 0;
151184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
152184914b5SBarry Smith       if (reason) {
15352392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
154454a90a3SBarry Smith         breakout = 1;
155454a90a3SBarry Smith         break;
15652392280SLois Curfman McInnes       }
15752392280SLois Curfman McInnes       snes->nfailures++;
1584800dd8cSBarry Smith     }
1591acabf8cSLois Curfman McInnes     if (!breakout) {
1604800dd8cSBarry Smith       fnorm = gnorm;
1610f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
162c509a162SBarry Smith       snes->iter = i+1;
163fbe28522SBarry Smith       snes->norm = fnorm;
1640f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16539e2f89bSBarry Smith       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
16639e2f89bSBarry Smith       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
167329f5518SBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);		/* xnorm = || X || */
1680462333dSBarry Smith       SNESLogConvHistory(snes,fnorm,lits);
16994a424c1SBarry Smith       SNESMonitor(snes,i+1,fnorm);
1704800dd8cSBarry Smith 
1714800dd8cSBarry Smith       /* Test for convergence */
17252392280SLois Curfman McInnes       neP->itflag = 1;
173184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
174184914b5SBarry Smith       if (reason) {
17538442cffSBarry Smith         break;
17638442cffSBarry Smith       }
1771acabf8cSLois Curfman McInnes     } else {
1781acabf8cSLois Curfman McInnes       break;
1791acabf8cSLois Curfman McInnes     }
18038442cffSBarry Smith   }
18138442cffSBarry Smith   /* Verify solution is in corect location */
182e15f7bb5SBarry Smith   if (X != snes->vec_sol) {
183393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
184e15f7bb5SBarry Smith   }
185e15f7bb5SBarry Smith   if (F != snes->vec_func) {
186e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
187e15f7bb5SBarry Smith   }
18839e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
18939e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
19052392280SLois Curfman McInnes   if (i == maxits) {
191b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits);
19252392280SLois Curfman McInnes     i--;
193184914b5SBarry Smith     reason = SNES_DIVERGED_MAX_IT;
19452392280SLois Curfman McInnes   }
19552392280SLois Curfman McInnes   *its = i+1;
1960f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
197184914b5SBarry Smith   snes->reason = reason;
1980f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1993a40ed3dSBarry Smith   PetscFunctionReturn(0);
2004800dd8cSBarry Smith }
2014800dd8cSBarry Smith /*------------------------------------------------------------*/
2024a2ae208SSatish Balay #undef __FUNCT__
2034a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp_EQ_TR"
204f63b844aSLois Curfman McInnes static int SNESSetUp_EQ_TR(SNES snes)
2054800dd8cSBarry Smith {
206fbe28522SBarry Smith   int ierr;
2073a40ed3dSBarry Smith 
2083a40ed3dSBarry Smith   PetscFunctionBegin;
20981b6cf68SLois Curfman McInnes   snes->nwork = 4;
210d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
211b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
21281b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2133a40ed3dSBarry Smith   PetscFunctionReturn(0);
2144800dd8cSBarry Smith }
2154800dd8cSBarry Smith /*------------------------------------------------------------*/
2164a2ae208SSatish Balay #undef __FUNCT__
2174a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy_EQ_TR"
218e1311b90SBarry Smith static int SNESDestroy_EQ_TR(SNES snes)
2194800dd8cSBarry Smith {
220393d2d9aSLois Curfman McInnes   int  ierr;
2215baf8537SBarry Smith 
2223a40ed3dSBarry Smith   PetscFunctionBegin;
2235baf8537SBarry Smith   if (snes->nwork) {
2244b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
2255baf8537SBarry Smith   }
226606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2273a40ed3dSBarry Smith   PetscFunctionReturn(0);
2284800dd8cSBarry Smith }
2294800dd8cSBarry Smith /*------------------------------------------------------------*/
2304800dd8cSBarry Smith 
2314a2ae208SSatish Balay #undef __FUNCT__
2324a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions_EQ_TR"
233f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_TR(SNES snes)
2344800dd8cSBarry Smith {
2356831982aSBarry Smith   SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data;
236f1af5d2fSBarry Smith   int        ierr;
2374800dd8cSBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
239b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
240b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
241b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
242b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
243b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
244b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
245b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
246b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
247b0a32e0cSBarry Smith     ierr = PetscOptionsDouble("-snes_eq_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
248b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2493a40ed3dSBarry Smith   PetscFunctionReturn(0);
2504800dd8cSBarry Smith }
2514800dd8cSBarry Smith 
2524a2ae208SSatish Balay #undef __FUNCT__
2534a2ae208SSatish Balay #define __FUNCT__ "SNESView_EQ_TR"
254b0a32e0cSBarry Smith static int SNESView_EQ_TR(SNES snes,PetscViewer viewer)
255a935fc98SLois Curfman McInnes {
2566831982aSBarry Smith   SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data;
25751695f54SLois Curfman McInnes   int        ierr;
2586831982aSBarry Smith   PetscTruth isascii;
259a935fc98SLois Curfman McInnes 
2603a40ed3dSBarry Smith   PetscFunctionBegin;
261b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
2620f5bd95cSBarry Smith   if (isascii) {
263b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
264b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
2655cd90555SBarry Smith   } else {
26629bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
26719bcc07fSBarry Smith   }
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
269a935fc98SLois Curfman McInnes }
270a935fc98SLois Curfman McInnes 
27152392280SLois Curfman McInnes /* ---------------------------------------------------------------- */
2724a2ae208SSatish Balay #undef __FUNCT__
2734a2ae208SSatish Balay #define __FUNCT__ "SNESConverged_EQ_TR"
2748ac28fe4SSatish Balay /*@C
275f525115eSLois Curfman McInnes    SNESConverged_EQ_TR - Monitors the convergence of the trust region
2766831982aSBarry Smith    method SNESEQTR for solving systems of nonlinear equations (default).
27752392280SLois Curfman McInnes 
278c7afd0dbSLois Curfman McInnes    Collective on SNES
279c7afd0dbSLois Curfman McInnes 
28052392280SLois Curfman McInnes    Input Parameters:
281c7afd0dbSLois Curfman McInnes +  snes - the SNES context
28252392280SLois Curfman McInnes .  xnorm - 2-norm of current iterate
28352392280SLois Curfman McInnes .  pnorm - 2-norm of current step
28452392280SLois Curfman McInnes .  fnorm - 2-norm of function
285c7afd0dbSLois Curfman McInnes -  dummy - unused context
286fee21e36SBarry Smith 
287184914b5SBarry Smith    Output Parameter:
288184914b5SBarry Smith .   reason - one of
289184914b5SBarry Smith $  SNES_CONVERGED_FNORM_ABS       - (fnorm < atol),
290184914b5SBarry Smith $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
291184914b5SBarry Smith $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
292184914b5SBarry Smith $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
293184914b5SBarry Smith $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
294184914b5SBarry Smith $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
295184914b5SBarry Smith $  SNES_CONVERGED_ITERATING       - (otherwise)
29652392280SLois Curfman McInnes 
29752392280SLois Curfman McInnes    where
298c7afd0dbSLois Curfman McInnes +    delta    - trust region paramenter
299c7afd0dbSLois Curfman McInnes .    deltatol - trust region size tolerance,
300c7afd0dbSLois Curfman McInnes                 set with SNESSetTrustRegionTolerance()
301c7afd0dbSLois Curfman McInnes .    maxf - maximum number of function evaluations,
302c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
303c7afd0dbSLois Curfman McInnes .    nfct - number of function evaluations,
304c7afd0dbSLois Curfman McInnes .    atol - absolute function norm tolerance,
305c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
306c7afd0dbSLois Curfman McInnes -    xtol - relative function norm tolerance,
307c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
30852392280SLois Curfman McInnes 
30915091d37SBarry Smith    Level: intermediate
31015091d37SBarry Smith 
31152392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
31252392280SLois Curfman McInnes 
31352392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
31452392280SLois Curfman McInnes @*/
315184914b5SBarry Smith int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy)
31652392280SLois Curfman McInnes {
3176831982aSBarry Smith   SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data;
318184914b5SBarry Smith   int        ierr;
319ddd12767SBarry Smith 
3203a40ed3dSBarry Smith   PetscFunctionBegin;
3213a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
32229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
3233a40ed3dSBarry Smith   }
32452392280SLois Curfman McInnes 
325d252947aSBarry Smith   if (fnorm != fnorm) {
326b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n");
327184914b5SBarry Smith     *reason = SNES_DIVERGED_FNORM_NAN;
328184914b5SBarry Smith   } else if (neP->delta < xnorm * snes->deltatol) {
329b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol);
330184914b5SBarry Smith     *reason = SNES_CONVERGED_TR_DELTA;
331184914b5SBarry Smith   } else if (neP->itflag) {
332184914b5SBarry Smith     ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
3331acabf8cSLois Curfman McInnes   } else if (snes->nfuncs > snes->max_funcs) {
334b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs);
335184914b5SBarry Smith     *reason = SNES_DIVERGED_FUNCTION_COUNT;
336184914b5SBarry Smith   } else {
337184914b5SBarry Smith     *reason = SNES_CONVERGED_ITERATING;
33852392280SLois Curfman McInnes   }
3393a40ed3dSBarry Smith   PetscFunctionReturn(0);
34052392280SLois Curfman McInnes }
34152392280SLois Curfman McInnes /* ------------------------------------------------------------ */
342fb2e594dSBarry Smith EXTERN_C_BEGIN
3434a2ae208SSatish Balay #undef __FUNCT__
3444a2ae208SSatish Balay #define __FUNCT__ "SNESCreate_EQ_TR"
345f63b844aSLois Curfman McInnes int SNESCreate_EQ_TR(SNES snes)
3464800dd8cSBarry Smith {
3476831982aSBarry Smith   SNES_EQ_TR *neP;
34882502324SSatish Balay   int ierr;
3494800dd8cSBarry Smith 
3503a40ed3dSBarry Smith   PetscFunctionBegin;
3513a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
35229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
3533a40ed3dSBarry Smith   }
354f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_TR;
355f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_TR;
356f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_TR;
357f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_TR;
358f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_TR;
359f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_TR;
3605baf8537SBarry Smith   snes->nwork           = 0;
361fbe28522SBarry Smith 
362b0a32e0cSBarry Smith   ierr			= PetscNew(SNES_EQ_TR,&neP);CHKERRQ(ierr);
363b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_EQ_TR));
364fbe28522SBarry Smith   snes->data	        = (void*)neP;
365fbe28522SBarry Smith   neP->mu		= 0.25;
366fbe28522SBarry Smith   neP->eta		= 0.75;
367fbe28522SBarry Smith   neP->delta		= 0.0;
368fbe28522SBarry Smith   neP->delta0		= 0.2;
369fbe28522SBarry Smith   neP->delta1		= 0.3;
370fbe28522SBarry Smith   neP->delta2		= 0.75;
371fbe28522SBarry Smith   neP->delta3		= 2.0;
372fbe28522SBarry Smith   neP->sigma		= 0.0001;
373fbe28522SBarry Smith   neP->itflag		= 0;
37452392280SLois Curfman McInnes   neP->rnorm0		= 0;
37552392280SLois Curfman McInnes   neP->ttol		= 0;
3763a40ed3dSBarry Smith   PetscFunctionReturn(0);
3774800dd8cSBarry Smith }
378fb2e594dSBarry Smith EXTERN_C_END
37982bf6240SBarry Smith 
380