xref: /petsc/src/snes/impls/tr/tr.c (revision 32077d6df18512e94702c86e7d562481ed0cebd0)
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"
104b27c08aSLois Curfman McInnes int 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;
173acb151aSSatish Balay   int                 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"
49c9780b6fSBarry Smith static int 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;
530462333dSBarry Smith   int                 maxits,i,ierr,lits,breakout = 0;
54112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
55064f8208SBarry Smith   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1;
56ea709b57SSatish Balay   PetscScalar         mone = -1.0,cnorm;
57df60cc22SBarry Smith   KSP                 ksp;
58184914b5SBarry Smith   SNESConvergedReason reason;
597c4af3dcSLois Curfman McInnes   PetscTruth          conv;
604800dd8cSBarry Smith 
613a40ed3dSBarry Smith   PetscFunctionBegin;
62fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
63fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
6439e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
65fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
66fbe28522SBarry Smith   G		= snes->work[1];
676b5873e3SBarry Smith   Ytmp          = snes->work[2];
684800dd8cSBarry Smith 
694c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
704c49b128SBarry Smith   snes->iter = 0;
714c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
727e6560a3SBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
734800dd8cSBarry Smith 
745334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
75cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
760f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
77fbe28522SBarry Smith   snes->norm = fnorm;
780f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
794800dd8cSBarry Smith   delta = neP->delta0*fnorm;
804800dd8cSBarry Smith   neP->delta = delta;
810462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
8294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
8394b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
844800dd8cSBarry Smith 
85c9780b6fSBarry Smith  if (fnorm < snes->atol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
86b37302c6SLois Curfman McInnes 
87d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
88d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
89d034289bSBarry Smith 
90a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
917c4af3dcSLois Curfman McInnes   ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr);
927c4af3dcSLois Curfman McInnes   if (!conv) {
934b27c08aSLois Curfman McInnes     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr);
944b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_TR: Using Krylov convergence test SNES_TR_KSPConverged_Private\n");
957c4af3dcSLois Curfman McInnes   }
96df60cc22SBarry Smith 
974800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
985334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
9994b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1002deea200SLois Curfman McInnes 
1012deea200SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
10223ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
103c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
104329f5518SBarry Smith     snes->linear_its += lits;
1054b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
106064f8208SBarry Smith     ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
107064f8208SBarry Smith     norm1 = nrm;
1086b5873e3SBarry Smith     while(1) {
109393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
110064f8208SBarry Smith       nrm = norm1;
1116b5873e3SBarry Smith 
1122deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
113064f8208SBarry Smith       if (nrm >= delta) {
114064f8208SBarry Smith         nrm = delta/nrm;
115064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
116064f8208SBarry Smith         cnorm = nrm;
1174b27c08aSLois Curfman McInnes         PetscLogInfo(snes,"SNESSolve_TR: Scaling direction by %g\n",nrm);
118393d2d9aSLois Curfman McInnes         ierr = VecScale(&cnorm,Y);CHKERRQ(ierr);
119064f8208SBarry Smith         nrm = gpnorm;
120fbe28522SBarry Smith         ynorm = delta;
121fbe28522SBarry Smith       } else {
122fbe28522SBarry Smith         gpnorm = 0.0;
1234b27c08aSLois Curfman McInnes         PetscLogInfo(snes,"SNESSolve_TR: Direction is in Trust Region\n");
124064f8208SBarry Smith         ynorm = nrm;
125fbe28522SBarry Smith       }
1262deea200SLois Curfman McInnes       ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr);            /* Y <- X - Y */
12781b6cf68SLois Curfman McInnes       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
1285334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
129cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1304800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
131fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1324800dd8cSBarry Smith 
1334800dd8cSBarry Smith       /* Update size of trust region */
1344800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1354800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1364800dd8cSBarry Smith       else                     delta *= neP->delta3;
1374b27c08aSLois Curfman McInnes       PetscLogInfo(snes,"SNESSolve_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm);
1384b27c08aSLois Curfman McInnes       PetscLogInfo(snes,"SNESSolve_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta);
1394800dd8cSBarry Smith       neP->delta = delta;
1404800dd8cSBarry Smith       if (rho > neP->sigma) break;
1414b27c08aSLois Curfman McInnes       PetscLogInfo(snes,"SNESSolve_TR: Trying again in smaller region\n");
1426b5873e3SBarry Smith       /* check to see if progress is hopeless */
14352392280SLois Curfman McInnes       neP->itflag = 0;
144184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
145184914b5SBarry Smith       if (reason) {
14652392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
1475ed2d596SBarry Smith         SNESMonitor(snes,i+1,fnorm);
148454a90a3SBarry Smith         breakout = 1;
149454a90a3SBarry Smith         break;
15052392280SLois Curfman McInnes       }
15150ffb88aSMatthew Knepley       snes->numFailures++;
1524800dd8cSBarry Smith     }
1531acabf8cSLois Curfman McInnes     if (!breakout) {
1544800dd8cSBarry Smith       fnorm = gnorm;
1550f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
156c509a162SBarry Smith       snes->iter = i+1;
157fbe28522SBarry Smith       snes->norm = fnorm;
1580f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15939e2f89bSBarry Smith       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
16039e2f89bSBarry Smith       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
161329f5518SBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);		/* xnorm = || X || */
1620462333dSBarry Smith       SNESLogConvHistory(snes,fnorm,lits);
16394a424c1SBarry Smith       SNESMonitor(snes,i+1,fnorm);
1644800dd8cSBarry Smith 
1654800dd8cSBarry Smith       /* Test for convergence */
16652392280SLois Curfman McInnes       neP->itflag = 1;
167184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
168184914b5SBarry Smith       if (reason) {
16938442cffSBarry Smith         break;
17038442cffSBarry Smith       }
1711acabf8cSLois Curfman McInnes     } else {
1721acabf8cSLois Curfman McInnes       break;
1731acabf8cSLois Curfman McInnes     }
17438442cffSBarry Smith   }
17538442cffSBarry Smith   /* Verify solution is in corect location */
176e15f7bb5SBarry Smith   if (X != snes->vec_sol) {
177393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
178e15f7bb5SBarry Smith   }
179e15f7bb5SBarry Smith   if (F != snes->vec_func) {
180e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
181e15f7bb5SBarry Smith   }
18239e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
18339e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
18452392280SLois Curfman McInnes   if (i == maxits) {
1854b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_TR: Maximum number of iterations has been reached: %d\n",maxits);
186184914b5SBarry Smith     reason = SNES_DIVERGED_MAX_IT;
18752392280SLois Curfman McInnes   }
1880f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
189184914b5SBarry Smith   snes->reason = reason;
1900f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1913a40ed3dSBarry Smith   PetscFunctionReturn(0);
1924800dd8cSBarry Smith }
1934800dd8cSBarry Smith /*------------------------------------------------------------*/
1944a2ae208SSatish Balay #undef __FUNCT__
1954b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR"
1964b27c08aSLois Curfman McInnes static int SNESSetUp_TR(SNES snes)
1974800dd8cSBarry Smith {
198fbe28522SBarry Smith   int ierr;
1993a40ed3dSBarry Smith 
2003a40ed3dSBarry Smith   PetscFunctionBegin;
20181b6cf68SLois Curfman McInnes   snes->nwork = 4;
202d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
203b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
20481b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2053a40ed3dSBarry Smith   PetscFunctionReturn(0);
2064800dd8cSBarry Smith }
2074800dd8cSBarry Smith /*------------------------------------------------------------*/
2084a2ae208SSatish Balay #undef __FUNCT__
2094b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR"
2104b27c08aSLois Curfman McInnes static int SNESDestroy_TR(SNES snes)
2114800dd8cSBarry Smith {
212393d2d9aSLois Curfman McInnes   int  ierr;
2135baf8537SBarry Smith 
2143a40ed3dSBarry Smith   PetscFunctionBegin;
2155baf8537SBarry Smith   if (snes->nwork) {
2164b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
2175baf8537SBarry Smith   }
218606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2193a40ed3dSBarry Smith   PetscFunctionReturn(0);
2204800dd8cSBarry Smith }
2214800dd8cSBarry Smith /*------------------------------------------------------------*/
2224800dd8cSBarry Smith 
2234a2ae208SSatish Balay #undef __FUNCT__
2244b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR"
2254b27c08aSLois Curfman McInnes static int SNESSetFromOptions_TR(SNES snes)
2264800dd8cSBarry Smith {
2274b27c08aSLois Curfman McInnes   SNES_TR *ctx = (SNES_TR *)snes->data;
228f1af5d2fSBarry Smith   int     ierr;
2294800dd8cSBarry Smith 
2303a40ed3dSBarry Smith   PetscFunctionBegin;
231b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
23287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
2334b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
2344b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
2354b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
2364b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
2374b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
2384b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
2394b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
240b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2413a40ed3dSBarry Smith   PetscFunctionReturn(0);
2424800dd8cSBarry Smith }
2434800dd8cSBarry Smith 
2444a2ae208SSatish Balay #undef __FUNCT__
2454b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR"
2464b27c08aSLois Curfman McInnes static int SNESView_TR(SNES snes,PetscViewer viewer)
247a935fc98SLois Curfman McInnes {
2484b27c08aSLois Curfman McInnes   SNES_TR *tr = (SNES_TR *)snes->data;
24951695f54SLois Curfman McInnes   int        ierr;
250*32077d6dSBarry Smith   PetscTruth iascii;
251a935fc98SLois Curfman McInnes 
2523a40ed3dSBarry Smith   PetscFunctionBegin;
253*32077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
254*32077d6dSBarry Smith   if (iascii) {
255b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
256b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
2575cd90555SBarry Smith   } else {
25829bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
25919bcc07fSBarry Smith   }
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
261a935fc98SLois Curfman McInnes }
262a935fc98SLois Curfman McInnes 
26352392280SLois Curfman McInnes /* ---------------------------------------------------------------- */
2644a2ae208SSatish Balay #undef __FUNCT__
2654b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR"
2668ac28fe4SSatish Balay /*@C
2674b27c08aSLois Curfman McInnes    SNESConverged_TR - Monitors the convergence of the trust region
2684b27c08aSLois Curfman McInnes    method SNESTR for solving systems of nonlinear equations (default).
26952392280SLois Curfman McInnes 
270c7afd0dbSLois Curfman McInnes    Collective on SNES
271c7afd0dbSLois Curfman McInnes 
27252392280SLois Curfman McInnes    Input Parameters:
273c7afd0dbSLois Curfman McInnes +  snes - the SNES context
27452392280SLois Curfman McInnes .  xnorm - 2-norm of current iterate
27552392280SLois Curfman McInnes .  pnorm - 2-norm of current step
27652392280SLois Curfman McInnes .  fnorm - 2-norm of function
277c7afd0dbSLois Curfman McInnes -  dummy - unused context
278fee21e36SBarry Smith 
279184914b5SBarry Smith    Output Parameter:
280184914b5SBarry Smith .   reason - one of
281184914b5SBarry Smith $  SNES_CONVERGED_FNORM_ABS       - (fnorm < atol),
282184914b5SBarry Smith $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
283184914b5SBarry Smith $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
284184914b5SBarry Smith $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
285184914b5SBarry Smith $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
286184914b5SBarry Smith $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
287184914b5SBarry Smith $  SNES_CONVERGED_ITERATING       - (otherwise)
28852392280SLois Curfman McInnes 
28952392280SLois Curfman McInnes    where
290c7afd0dbSLois Curfman McInnes +    delta    - trust region paramenter
291c7afd0dbSLois Curfman McInnes .    deltatol - trust region size tolerance,
292c7afd0dbSLois Curfman McInnes                 set with SNESSetTrustRegionTolerance()
293c7afd0dbSLois Curfman McInnes .    maxf - maximum number of function evaluations,
294c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
295c7afd0dbSLois Curfman McInnes .    nfct - number of function evaluations,
296c7afd0dbSLois Curfman McInnes .    atol - absolute function norm tolerance,
297c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
298c7afd0dbSLois Curfman McInnes -    xtol - relative function norm tolerance,
299c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
30052392280SLois Curfman McInnes 
30115091d37SBarry Smith    Level: intermediate
30215091d37SBarry Smith 
30352392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
30452392280SLois Curfman McInnes 
30552392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
30652392280SLois Curfman McInnes @*/
3074b27c08aSLois Curfman McInnes int SNESConverged_TR(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
30852392280SLois Curfman McInnes {
3094b27c08aSLois Curfman McInnes   SNES_TR *neP = (SNES_TR *)snes->data;
310184914b5SBarry Smith   int     ierr;
311ddd12767SBarry Smith 
3123a40ed3dSBarry Smith   PetscFunctionBegin;
313d252947aSBarry Smith   if (fnorm != fnorm) {
3144b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESConverged_TR:Failed to converged, function norm is NaN\n");
315184914b5SBarry Smith     *reason = SNES_DIVERGED_FNORM_NAN;
316184914b5SBarry Smith   } else if (neP->delta < xnorm * snes->deltatol) {
3174b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESConverged_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol);
318184914b5SBarry Smith     *reason = SNES_CONVERGED_TR_DELTA;
319184914b5SBarry Smith   } else if (neP->itflag) {
3204b27c08aSLois Curfman McInnes     ierr = SNESConverged_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
3211acabf8cSLois Curfman McInnes   } else if (snes->nfuncs > snes->max_funcs) {
3224b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESConverged_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs);
323184914b5SBarry Smith     *reason = SNES_DIVERGED_FUNCTION_COUNT;
324184914b5SBarry Smith   } else {
325184914b5SBarry Smith     *reason = SNES_CONVERGED_ITERATING;
32652392280SLois Curfman McInnes   }
3273a40ed3dSBarry Smith   PetscFunctionReturn(0);
32852392280SLois Curfman McInnes }
32952392280SLois Curfman McInnes /* ------------------------------------------------------------ */
330ebe3b25bSBarry Smith /*MC
331ebe3b25bSBarry Smith       SNESTR - Newton based nonlinear solver that uses a trust region
332ebe3b25bSBarry Smith 
333ebe3b25bSBarry Smith    Options Database:
334ebe3b25bSBarry Smith +    -snes_trtol <tol> Trust region tolerance
335ebe3b25bSBarry Smith .    -snes_tr_mu <mu>
336ebe3b25bSBarry Smith .    -snes_tr_eta <eta>
337ebe3b25bSBarry Smith .    -snes_tr_sigma <sigma>
338ebe3b25bSBarry Smith .    -snes_tr_delta0 <delta0>
339ebe3b25bSBarry Smith .    -snes_tr_delta1 <delta1>
340ebe3b25bSBarry Smith .    -snes_tr_delta2 <delta2>
341ebe3b25bSBarry Smith -    -snes_tr_delta3 <delta3>
342ebe3b25bSBarry Smith 
343ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
344ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
345ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
346ebe3b25bSBarry Smith 
347ebe3b25bSBarry Smith    This is intended as a model implementation, since it does not
348ebe3b25bSBarry Smith    necessarily have many of the bells and whistles of other
349ebe3b25bSBarry Smith    implementations.
350ebe3b25bSBarry Smith 
351ee3001cbSBarry Smith    Level: intermediate
352ee3001cbSBarry Smith 
353ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
354ebe3b25bSBarry Smith 
355ebe3b25bSBarry Smith M*/
356fb2e594dSBarry Smith EXTERN_C_BEGIN
3574a2ae208SSatish Balay #undef __FUNCT__
3584b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR"
3594b27c08aSLois Curfman McInnes int SNESCreate_TR(SNES snes)
3604800dd8cSBarry Smith {
3614b27c08aSLois Curfman McInnes   SNES_TR *neP;
36282502324SSatish Balay   int     ierr;
3634800dd8cSBarry Smith 
3643a40ed3dSBarry Smith   PetscFunctionBegin;
3654b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_TR;
3664b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_TR;
3674b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_TR;
3684b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_TR;
3694b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_TR;
3704b27c08aSLois Curfman McInnes   snes->view            = SNESView_TR;
3715baf8537SBarry Smith   snes->nwork           = 0;
372fbe28522SBarry Smith 
3734b27c08aSLois Curfman McInnes   ierr			= PetscNew(SNES_TR,&neP);CHKERRQ(ierr);
3744b27c08aSLois Curfman McInnes   PetscLogObjectMemory(snes,sizeof(SNES_TR));
375fbe28522SBarry Smith   snes->data	        = (void*)neP;
376fbe28522SBarry Smith   neP->mu		= 0.25;
377fbe28522SBarry Smith   neP->eta		= 0.75;
378fbe28522SBarry Smith   neP->delta		= 0.0;
379fbe28522SBarry Smith   neP->delta0		= 0.2;
380fbe28522SBarry Smith   neP->delta1		= 0.3;
381fbe28522SBarry Smith   neP->delta2		= 0.75;
382fbe28522SBarry Smith   neP->delta3		= 2.0;
383fbe28522SBarry Smith   neP->sigma		= 0.0001;
384fbe28522SBarry Smith   neP->itflag		= 0;
38552392280SLois Curfman McInnes   neP->rnorm0		= 0;
38652392280SLois Curfman McInnes   neP->ttol		= 0;
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
3884800dd8cSBarry Smith }
389fb2e594dSBarry Smith EXTERN_C_END
39082bf6240SBarry Smith 
391