xref: /petsc/src/snes/impls/tr/tr.c (revision 410397dcfd5ca393ac811969750485b0d82fb031)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
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__
104b27c08aSLois Curfman McInnes #define __FUNCT__ "SNES_TR_KSPConverged_Private"
1177431f27SBarry Smith PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal 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;
154b27c08aSLois Curfman McInnes   SNES_TR             *neP = (SNES_TR*)snes->data;
16df60cc22SBarry Smith   Vec                 x;
17064f8208SBarry Smith   PetscReal           nrm;
18dfbe8321SBarry Smith   PetscErrorCode      ierr;
19df60cc22SBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
2198120f82SLois Curfman McInnes   if (snes->ksp_ewconv) {
2290faec13SBarry Smith     if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker convergence 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) {
28ae15b995SBarry Smith     ierr = PetscInfo2(snes,"regular convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr);
299b04593eSLois Curfman McInnes   }
30df60cc22SBarry Smith 
31a935fc98SLois Curfman McInnes   /* Determine norm of solution */
3278b31e54SBarry Smith   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
33064f8208SBarry Smith   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
34064f8208SBarry Smith   if (nrm >= neP->delta) {
35ae15b995SBarry Smith     ierr = PetscInfo2(snes,"KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr);
36ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr);
37329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
38df60cc22SBarry Smith   }
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
40df60cc22SBarry Smith }
4182bf6240SBarry Smith 
42df60cc22SBarry Smith /*
434b27c08aSLois Curfman McInnes    SNESSolve_TR - Implements Newton's Method with a very simple trust
44ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
454800dd8cSBarry Smith 
46ddbbdb52SLois Curfman McInnes 
474800dd8cSBarry Smith */
484a2ae208SSatish Balay #undef __FUNCT__
494b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR"
506849ba73SBarry Smith static PetscErrorCode SNESSolve_TR(SNES snes)
514800dd8cSBarry Smith {
524b27c08aSLois Curfman McInnes   SNES_TR             *neP = (SNES_TR*)snes->data;
536b5873e3SBarry Smith   Vec                 X,F,Y,G,TMP,Ytmp;
546849ba73SBarry Smith   PetscErrorCode      ierr;
55a7cc72afSBarry Smith   PetscInt            maxits,i,lits;
56112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
57064f8208SBarry Smith   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1;
5879f36460SBarry Smith   PetscScalar         cnorm;
59df60cc22SBarry Smith   KSP                 ksp;
60184914b5SBarry Smith   SNESConvergedReason reason;
61a7cc72afSBarry Smith   PetscTruth          conv,breakout = PETSC_FALSE;
624800dd8cSBarry Smith 
633a40ed3dSBarry Smith   PetscFunctionBegin;
64fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
65fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
6639e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
67fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
68fbe28522SBarry Smith   G		= snes->work[1];
696b5873e3SBarry Smith   Ytmp          = snes->work[2];
704800dd8cSBarry Smith 
714c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
724c49b128SBarry Smith   snes->iter = 0;
734c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
747e6560a3SBarry Smith   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
754800dd8cSBarry Smith 
765334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
77cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
780f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
79fbe28522SBarry Smith   snes->norm = fnorm;
800f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
814800dd8cSBarry Smith   delta = neP->delta0*fnorm;
824800dd8cSBarry Smith   neP->delta = delta;
830462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
8494a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
8594b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
864800dd8cSBarry Smith 
8770441072SBarry Smith  if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
88b37302c6SLois Curfman McInnes 
89d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
90d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
91d034289bSBarry Smith 
92a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
937c4af3dcSLois Curfman McInnes   ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr);
947c4af3dcSLois Curfman McInnes   if (!conv) {
954b27c08aSLois Curfman McInnes     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr);
96ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
977c4af3dcSLois Curfman McInnes   }
98df60cc22SBarry Smith 
994800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
10099a96b7cSMatthew Knepley 
10199a96b7cSMatthew Knepley     /* Call general purpose update function */
102e7788613SBarry Smith     if (snes->ops->update) {
103e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
10499a96b7cSMatthew Knepley     }
10599a96b7cSMatthew Knepley 
1065334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
10794b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1082deea200SLois Curfman McInnes 
1092deea200SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
11023ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
111c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
112329f5518SBarry Smith     snes->linear_its += lits;
113ae15b995SBarry Smith     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
114064f8208SBarry Smith     ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
115064f8208SBarry Smith     norm1 = nrm;
1166b5873e3SBarry Smith     while(1) {
117393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
118064f8208SBarry Smith       nrm = norm1;
1196b5873e3SBarry Smith 
1202deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
121064f8208SBarry Smith       if (nrm >= delta) {
122064f8208SBarry Smith         nrm = delta/nrm;
123064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
124064f8208SBarry Smith         cnorm = nrm;
125ae15b995SBarry Smith         ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr);
1262dcb1b2aSMatthew Knepley         ierr = VecScale(Y,cnorm);CHKERRQ(ierr);
127064f8208SBarry Smith         nrm = gpnorm;
128fbe28522SBarry Smith         ynorm = delta;
129fbe28522SBarry Smith       } else {
130fbe28522SBarry Smith         gpnorm = 0.0;
131ae15b995SBarry Smith         ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
132064f8208SBarry Smith         ynorm = nrm;
133fbe28522SBarry Smith       }
13479f36460SBarry Smith       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
13581b6cf68SLois Curfman McInnes       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
1365334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
137cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1384800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
139fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1404800dd8cSBarry Smith 
1414800dd8cSBarry Smith       /* Update size of trust region */
1424800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1434800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1444800dd8cSBarry Smith       else                     delta *= neP->delta3;
145ae15b995SBarry Smith       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
146ae15b995SBarry Smith       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
1474800dd8cSBarry Smith       neP->delta = delta;
1484800dd8cSBarry Smith       if (rho > neP->sigma) break;
149ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
1506b5873e3SBarry Smith       /* check to see if progress is hopeless */
151ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
152e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
153184914b5SBarry Smith       if (reason) {
15452392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
1555ed2d596SBarry Smith         SNESMonitor(snes,i+1,fnorm);
156a7cc72afSBarry Smith         breakout = PETSC_TRUE;
157454a90a3SBarry Smith         break;
15852392280SLois Curfman McInnes       }
15950ffb88aSMatthew Knepley       snes->numFailures++;
1604800dd8cSBarry Smith     }
1611acabf8cSLois Curfman McInnes     if (!breakout) {
1624800dd8cSBarry Smith       fnorm = gnorm;
1630f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
164c509a162SBarry Smith       snes->iter = i+1;
165fbe28522SBarry Smith       snes->norm = fnorm;
1660f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16739e2f89bSBarry Smith       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
16839e2f89bSBarry Smith       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
169329f5518SBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);		/* xnorm = || X || */
1700462333dSBarry Smith       SNESLogConvHistory(snes,fnorm,lits);
17194a424c1SBarry Smith       SNESMonitor(snes,i+1,fnorm);
1724800dd8cSBarry Smith 
1734800dd8cSBarry Smith       /* Test for convergence */
174ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
175e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
176184914b5SBarry Smith       if (reason) {
17738442cffSBarry Smith         break;
17838442cffSBarry Smith       }
1791acabf8cSLois Curfman McInnes     } else {
1801acabf8cSLois Curfman McInnes       break;
1811acabf8cSLois Curfman McInnes     }
18238442cffSBarry Smith   }
18338442cffSBarry Smith   /* Verify solution is in corect location */
184e15f7bb5SBarry Smith   if (X != snes->vec_sol) {
185393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
186e15f7bb5SBarry Smith   }
187e15f7bb5SBarry Smith   if (F != snes->vec_func) {
188e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
189e15f7bb5SBarry Smith   }
19039e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
19139e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
19252392280SLois Curfman McInnes   if (i == maxits) {
193ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
194184914b5SBarry Smith     reason = SNES_DIVERGED_MAX_IT;
19552392280SLois Curfman McInnes   }
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__
2034b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR"
2046849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes)
2054800dd8cSBarry Smith {
206dfbe8321SBarry Smith   PetscErrorCode ierr;
2073a40ed3dSBarry Smith 
2083a40ed3dSBarry Smith   PetscFunctionBegin;
209*410397dcSLisandro Dalcin   if (!snes->work) {
21081b6cf68SLois Curfman McInnes     snes->nwork = 4;
211d7e8b826SBarry Smith     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
212efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
213*410397dcSLisandro Dalcin   }
21481b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2153a40ed3dSBarry Smith   PetscFunctionReturn(0);
2164800dd8cSBarry Smith }
2174800dd8cSBarry Smith /*------------------------------------------------------------*/
2184a2ae208SSatish Balay #undef __FUNCT__
2194b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR"
2206849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes)
2214800dd8cSBarry Smith {
222dfbe8321SBarry Smith   PetscErrorCode ierr;
2235baf8537SBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
2255baf8537SBarry Smith   if (snes->nwork) {
2264b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
2275baf8537SBarry Smith   }
228606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2293a40ed3dSBarry Smith   PetscFunctionReturn(0);
2304800dd8cSBarry Smith }
2314800dd8cSBarry Smith /*------------------------------------------------------------*/
2324800dd8cSBarry Smith 
2334a2ae208SSatish Balay #undef __FUNCT__
2344b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR"
2356849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
2364800dd8cSBarry Smith {
2374b27c08aSLois Curfman McInnes   SNES_TR *ctx = (SNES_TR *)snes->data;
238dfbe8321SBarry Smith   PetscErrorCode ierr;
2394800dd8cSBarry Smith 
2403a40ed3dSBarry Smith   PetscFunctionBegin;
241b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
24287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
2434b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
2444b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
2454b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
2464b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
2474b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
2484b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
2494b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
250b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2524800dd8cSBarry Smith }
2534800dd8cSBarry Smith 
2544a2ae208SSatish Balay #undef __FUNCT__
2554b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR"
2566849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
257a935fc98SLois Curfman McInnes {
2584b27c08aSLois Curfman McInnes   SNES_TR *tr = (SNES_TR *)snes->data;
259dfbe8321SBarry Smith   PetscErrorCode ierr;
26032077d6dSBarry Smith   PetscTruth iascii;
261a935fc98SLois Curfman McInnes 
2623a40ed3dSBarry Smith   PetscFunctionBegin;
26332077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
26432077d6dSBarry Smith   if (iascii) {
265a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
266a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
2675cd90555SBarry Smith   } else {
26879a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
26919bcc07fSBarry Smith   }
2703a40ed3dSBarry Smith   PetscFunctionReturn(0);
271a935fc98SLois Curfman McInnes }
272a935fc98SLois Curfman McInnes 
27352392280SLois Curfman McInnes /* ---------------------------------------------------------------- */
2744a2ae208SSatish Balay #undef __FUNCT__
2754b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR"
2768ac28fe4SSatish Balay /*@C
2774b27c08aSLois Curfman McInnes    SNESConverged_TR - Monitors the convergence of the trust region
2784b27c08aSLois Curfman McInnes    method SNESTR for solving systems of nonlinear equations (default).
27952392280SLois Curfman McInnes 
280c7afd0dbSLois Curfman McInnes    Collective on SNES
281c7afd0dbSLois Curfman McInnes 
28252392280SLois Curfman McInnes    Input Parameters:
283c7afd0dbSLois Curfman McInnes +  snes - the SNES context
28452392280SLois Curfman McInnes .  xnorm - 2-norm of current iterate
28552392280SLois Curfman McInnes .  pnorm - 2-norm of current step
28652392280SLois Curfman McInnes .  fnorm - 2-norm of function
287c7afd0dbSLois Curfman McInnes -  dummy - unused context
288fee21e36SBarry Smith 
289184914b5SBarry Smith    Output Parameter:
290184914b5SBarry Smith .   reason - one of
29170441072SBarry Smith $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
292184914b5SBarry Smith $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
293184914b5SBarry Smith $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
294184914b5SBarry Smith $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
295184914b5SBarry Smith $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
296184914b5SBarry Smith $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
297184914b5SBarry Smith $  SNES_CONVERGED_ITERATING       - (otherwise)
29852392280SLois Curfman McInnes 
29952392280SLois Curfman McInnes    where
300c7afd0dbSLois Curfman McInnes +    delta    - trust region paramenter
301c7afd0dbSLois Curfman McInnes .    deltatol - trust region size tolerance,
302c7afd0dbSLois Curfman McInnes                 set with SNESSetTrustRegionTolerance()
303c7afd0dbSLois Curfman McInnes .    maxf - maximum number of function evaluations,
304c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
305c7afd0dbSLois Curfman McInnes .    nfct - number of function evaluations,
30670441072SBarry Smith .    abstol - absolute function norm tolerance,
307c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
308c7afd0dbSLois Curfman McInnes -    xtol - relative function norm tolerance,
309c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
31052392280SLois Curfman McInnes 
31115091d37SBarry Smith    Level: intermediate
31215091d37SBarry Smith 
31352392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
31452392280SLois Curfman McInnes 
31552392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
31652392280SLois Curfman McInnes @*/
31706ee9f85SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_TR(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
31852392280SLois Curfman McInnes {
3194b27c08aSLois Curfman McInnes   SNES_TR *neP = (SNES_TR *)snes->data;
320dfbe8321SBarry Smith   PetscErrorCode ierr;
321ddd12767SBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
323d252947aSBarry Smith   if (fnorm != fnorm) {
324ae15b995SBarry Smith     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
325184914b5SBarry Smith     *reason = SNES_DIVERGED_FNORM_NAN;
326184914b5SBarry Smith   } else if (neP->delta < xnorm * snes->deltatol) {
327ae15b995SBarry Smith     ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,snes->deltatol);CHKERRQ(ierr);
328184914b5SBarry Smith     *reason = SNES_CONVERGED_TR_DELTA;
329184914b5SBarry Smith   } else if (neP->itflag) {
33006ee9f85SBarry Smith     ierr = SNESConverged_LS(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
33143e71028SBarry Smith   } else if (snes->nfuncs >= snes->max_funcs) {
332ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
333184914b5SBarry Smith     *reason = SNES_DIVERGED_FUNCTION_COUNT;
334184914b5SBarry Smith   } else {
335184914b5SBarry Smith     *reason = SNES_CONVERGED_ITERATING;
33652392280SLois Curfman McInnes   }
3373a40ed3dSBarry Smith   PetscFunctionReturn(0);
33852392280SLois Curfman McInnes }
33952392280SLois Curfman McInnes /* ------------------------------------------------------------ */
340ebe3b25bSBarry Smith /*MC
341ebe3b25bSBarry Smith       SNESTR - Newton based nonlinear solver that uses a trust region
342ebe3b25bSBarry Smith 
343ebe3b25bSBarry Smith    Options Database:
344ebe3b25bSBarry Smith +    -snes_trtol <tol> Trust region tolerance
345ebe3b25bSBarry Smith .    -snes_tr_mu <mu>
346ebe3b25bSBarry Smith .    -snes_tr_eta <eta>
347ebe3b25bSBarry Smith .    -snes_tr_sigma <sigma>
348ebe3b25bSBarry Smith .    -snes_tr_delta0 <delta0>
349ebe3b25bSBarry Smith .    -snes_tr_delta1 <delta1>
350ebe3b25bSBarry Smith .    -snes_tr_delta2 <delta2>
351ebe3b25bSBarry Smith -    -snes_tr_delta3 <delta3>
352ebe3b25bSBarry Smith 
353ebe3b25bSBarry Smith    The basic algorithm is taken from "The Minpack Project", by More',
354ebe3b25bSBarry Smith    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
355ebe3b25bSBarry Smith    of Mathematical Software", Wayne Cowell, editor.
356ebe3b25bSBarry Smith 
357ebe3b25bSBarry Smith    This is intended as a model implementation, since it does not
358ebe3b25bSBarry Smith    necessarily have many of the bells and whistles of other
359ebe3b25bSBarry Smith    implementations.
360ebe3b25bSBarry Smith 
361ee3001cbSBarry Smith    Level: intermediate
362ee3001cbSBarry Smith 
363ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
364ebe3b25bSBarry Smith 
365ebe3b25bSBarry Smith M*/
366fb2e594dSBarry Smith EXTERN_C_BEGIN
3674a2ae208SSatish Balay #undef __FUNCT__
3684b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_TR"
36963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_TR(SNES snes)
3704800dd8cSBarry Smith {
3714b27c08aSLois Curfman McInnes   SNES_TR        *neP;
372dfbe8321SBarry Smith   PetscErrorCode ierr;
3734800dd8cSBarry Smith 
3743a40ed3dSBarry Smith   PetscFunctionBegin;
375e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_TR;
376e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_TR;
377e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_TR;
378e7788613SBarry Smith   snes->ops->converged	     = SNESConverged_TR;
379e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_TR;
380e7788613SBarry Smith   snes->ops->view            = SNESView_TR;
3815baf8537SBarry Smith   snes->nwork                = 0;
382fbe28522SBarry Smith 
3834b27c08aSLois Curfman McInnes   ierr			= PetscNew(SNES_TR,&neP);CHKERRQ(ierr);
38452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_TR));CHKERRQ(ierr);
385fbe28522SBarry Smith   snes->data	        = (void*)neP;
386fbe28522SBarry Smith   neP->mu		= 0.25;
387fbe28522SBarry Smith   neP->eta		= 0.75;
388fbe28522SBarry Smith   neP->delta		= 0.0;
389fbe28522SBarry Smith   neP->delta0		= 0.2;
390fbe28522SBarry Smith   neP->delta1		= 0.3;
391fbe28522SBarry Smith   neP->delta2		= 0.75;
392fbe28522SBarry Smith   neP->delta3		= 2.0;
393fbe28522SBarry Smith   neP->sigma		= 0.0001;
394ef87ac0dSKris Buschelman   neP->itflag		= PETSC_FALSE;
39552392280SLois Curfman McInnes   neP->rnorm0		= 0;
39652392280SLois Curfman McInnes   neP->ttol		= 0;
3973a40ed3dSBarry Smith   PetscFunctionReturn(0);
3984800dd8cSBarry Smith }
399fb2e594dSBarry Smith EXTERN_C_END
40082bf6240SBarry Smith 
401