xref: /petsc/src/snes/impls/tr/tr.c (revision 38f2d2fdb3b6f522a3102c6eb796cebecf3224c0)
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;
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;
20329f5518SBarry Smith   ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr);
21329f5518SBarry Smith   if (*reason) {
2271f87433Sdalcinl     ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr);
239b04593eSLois Curfman McInnes   }
24a935fc98SLois Curfman McInnes   /* Determine norm of solution */
2578b31e54SBarry Smith   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
26064f8208SBarry Smith   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
27064f8208SBarry Smith   if (nrm >= neP->delta) {
28ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr);
29329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
30df60cc22SBarry Smith   }
313a40ed3dSBarry Smith   PetscFunctionReturn(0);
32df60cc22SBarry Smith }
3382bf6240SBarry Smith 
34df60cc22SBarry Smith /*
354b27c08aSLois Curfman McInnes    SNESSolve_TR - Implements Newton's Method with a very simple trust
36ddbbdb52SLois Curfman McInnes    region approach for solving systems of nonlinear equations.
374800dd8cSBarry Smith 
38ddbbdb52SLois Curfman McInnes 
394800dd8cSBarry Smith */
404a2ae208SSatish Balay #undef __FUNCT__
414b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_TR"
426849ba73SBarry Smith static PetscErrorCode SNESSolve_TR(SNES snes)
434800dd8cSBarry Smith {
444b27c08aSLois Curfman McInnes   SNES_TR             *neP = (SNES_TR*)snes->data;
456b5873e3SBarry Smith   Vec                 X,F,Y,G,TMP,Ytmp;
466849ba73SBarry Smith   PetscErrorCode      ierr;
47a7cc72afSBarry Smith   PetscInt            maxits,i,lits;
48112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
49064f8208SBarry Smith   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1;
5079f36460SBarry Smith   PetscScalar         cnorm;
51df60cc22SBarry Smith   KSP                 ksp;
523f149594SLisandro Dalcin   SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
53a7cc72afSBarry Smith   PetscTruth          conv,breakout = PETSC_FALSE;
544800dd8cSBarry Smith 
553a40ed3dSBarry Smith   PetscFunctionBegin;
56fbe28522SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
57fbe28522SBarry Smith   X		= snes->vec_sol;	/* solution vector */
5839e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
59fbe28522SBarry Smith   Y		= snes->work[0];	/* work vectors */
60fbe28522SBarry Smith   G		= snes->work[1];
616b5873e3SBarry Smith   Ytmp          = snes->work[2];
624800dd8cSBarry Smith 
634c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
644c49b128SBarry Smith   snes->iter = 0;
654c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
664800dd8cSBarry Smith 
675334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
68cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
690f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
70fbe28522SBarry Smith   snes->norm = fnorm;
710f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
724800dd8cSBarry Smith   delta = neP->delta0*fnorm;
734800dd8cSBarry Smith   neP->delta = delta;
740462333dSBarry Smith   SNESLogConvHistory(snes,fnorm,0);
7594a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
76b37302c6SLois Curfman McInnes 
77d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
78d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
79d034289bSBarry Smith 
803f149594SLisandro Dalcin   /* XXX Sould we use snes->ops->converged like in SNESLS ?*/
813f149594SLisandro Dalcin   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
823f149594SLisandro Dalcin   if (snes->ops->converged) {
833f149594SLisandro Dalcin     ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
843f149594SLisandro Dalcin   }
853f149594SLisandro Dalcin 
86a935fc98SLois Curfman McInnes   /* Set the stopping criteria to use the More' trick. */
877c4af3dcSLois Curfman McInnes   ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr);
887c4af3dcSLois Curfman McInnes   if (!conv) {
893f149594SLisandro Dalcin     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
904b27c08aSLois Curfman McInnes     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr);
91ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
927c4af3dcSLois Curfman McInnes   }
93df60cc22SBarry Smith 
944800dd8cSBarry Smith   for (i=0; i<maxits; i++) {
9599a96b7cSMatthew Knepley 
9699a96b7cSMatthew Knepley     /* Call general purpose update function */
97e7788613SBarry Smith     if (snes->ops->update) {
98e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
9999a96b7cSMatthew Knepley     }
10099a96b7cSMatthew Knepley 
1015334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
10294b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
1032deea200SLois Curfman McInnes 
1042deea200SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1053f149594SLisandro Dalcin     ierr = SNES_KSPSolve(snes,snes->ksp,F,Ytmp);CHKERRQ(ierr);
1063f149594SLisandro Dalcin     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
107329f5518SBarry Smith     snes->linear_its += lits;
108ae15b995SBarry Smith     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
109064f8208SBarry Smith     ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
110064f8208SBarry Smith     norm1 = nrm;
1116b5873e3SBarry Smith     while(1) {
112393d2d9aSLois Curfman McInnes       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
113064f8208SBarry Smith       nrm = norm1;
1146b5873e3SBarry Smith 
1152deea200SLois Curfman McInnes       /* Scale Y if need be and predict new value of F norm */
116064f8208SBarry Smith       if (nrm >= delta) {
117064f8208SBarry Smith         nrm = delta/nrm;
118064f8208SBarry Smith         gpnorm = (1.0 - nrm)*fnorm;
119064f8208SBarry Smith         cnorm = nrm;
120ae15b995SBarry Smith         ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr);
1212dcb1b2aSMatthew Knepley         ierr = VecScale(Y,cnorm);CHKERRQ(ierr);
122064f8208SBarry Smith         nrm = gpnorm;
123fbe28522SBarry Smith         ynorm = delta;
124fbe28522SBarry Smith       } else {
125fbe28522SBarry Smith         gpnorm = 0.0;
126ae15b995SBarry Smith         ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
127064f8208SBarry Smith         ynorm = nrm;
128fbe28522SBarry Smith       }
12979f36460SBarry Smith       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
13081b6cf68SLois Curfman McInnes       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
1315334005bSBarry Smith       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
132cddf8d76SBarry Smith       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
1334800dd8cSBarry Smith       if (fnorm == gpnorm) rho = 0.0;
134fbe28522SBarry Smith       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
1354800dd8cSBarry Smith 
1364800dd8cSBarry Smith       /* Update size of trust region */
1374800dd8cSBarry Smith       if      (rho < neP->mu)  delta *= neP->delta1;
1384800dd8cSBarry Smith       else if (rho < neP->eta) delta *= neP->delta2;
1394800dd8cSBarry Smith       else                     delta *= neP->delta3;
140ae15b995SBarry Smith       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
141ae15b995SBarry Smith       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
1424800dd8cSBarry Smith       neP->delta = delta;
1434800dd8cSBarry Smith       if (rho > neP->sigma) break;
144ae15b995SBarry Smith       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
1456b5873e3SBarry Smith       /* check to see if progress is hopeless */
146ef87ac0dSKris Buschelman       neP->itflag = PETSC_FALSE;
1473f149594SLisandro Dalcin       if (snes->ops->converged) {
148e7788613SBarry Smith 	ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
1493f149594SLisandro Dalcin       }
150184914b5SBarry Smith       if (reason) {
15152392280SLois Curfman McInnes         /* We're not progressing, so return with the current iterate */
1525ed2d596SBarry Smith         SNESMonitor(snes,i+1,fnorm);
153a7cc72afSBarry Smith         breakout = PETSC_TRUE;
154454a90a3SBarry Smith         break;
15552392280SLois Curfman McInnes       }
15650ffb88aSMatthew Knepley       snes->numFailures++;
1574800dd8cSBarry Smith     }
1581acabf8cSLois Curfman McInnes     if (!breakout) {
1594800dd8cSBarry Smith       fnorm = gnorm;
1600f5bd95cSBarry Smith       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
161c509a162SBarry Smith       snes->iter = i+1;
162fbe28522SBarry Smith       snes->norm = fnorm;
1630f5bd95cSBarry Smith       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16439e2f89bSBarry Smith       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
16539e2f89bSBarry Smith       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
1660462333dSBarry Smith       SNESLogConvHistory(snes,fnorm,lits);
16794a424c1SBarry Smith       SNESMonitor(snes,i+1,fnorm);
1684800dd8cSBarry Smith       /* Test for convergence */
169ef87ac0dSKris Buschelman       neP->itflag = PETSC_TRUE;
1703f149594SLisandro Dalcin       if (snes->ops->converged) {
1713f149594SLisandro Dalcin 	ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
172e7788613SBarry Smith 	ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
17338442cffSBarry Smith       }
1743f149594SLisandro Dalcin       if (reason) break;
1751acabf8cSLois Curfman McInnes     } else {
1761acabf8cSLois Curfman McInnes       break;
1771acabf8cSLois Curfman McInnes     }
17838442cffSBarry Smith   }
17938442cffSBarry Smith   /* Verify solution is in corect location */
180e15f7bb5SBarry Smith   if (X != snes->vec_sol) {
181393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
182e15f7bb5SBarry Smith   }
183e15f7bb5SBarry Smith   if (F != snes->vec_func) {
184e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
185e15f7bb5SBarry Smith   }
18639e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
18739e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
18852392280SLois Curfman McInnes   if (i == maxits) {
189ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
190184914b5SBarry Smith     reason = SNES_DIVERGED_MAX_IT;
19152392280SLois Curfman McInnes   }
1920f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
193184914b5SBarry Smith   snes->reason = reason;
1940f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1964800dd8cSBarry Smith }
1974800dd8cSBarry Smith /*------------------------------------------------------------*/
1984a2ae208SSatish Balay #undef __FUNCT__
1994b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_TR"
2006849ba73SBarry Smith static PetscErrorCode SNESSetUp_TR(SNES snes)
2014800dd8cSBarry Smith {
202dfbe8321SBarry Smith   PetscErrorCode ierr;
2033a40ed3dSBarry Smith 
2043a40ed3dSBarry Smith   PetscFunctionBegin;
205410397dcSLisandro Dalcin   if (!snes->work) {
20681b6cf68SLois Curfman McInnes     snes->nwork = 4;
207d7e8b826SBarry Smith     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
208efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
209410397dcSLisandro Dalcin   }
21081b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
2124800dd8cSBarry Smith }
2134800dd8cSBarry Smith /*------------------------------------------------------------*/
2144a2ae208SSatish Balay #undef __FUNCT__
2154b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_TR"
2166849ba73SBarry Smith static PetscErrorCode SNESDestroy_TR(SNES snes)
2174800dd8cSBarry Smith {
218dfbe8321SBarry Smith   PetscErrorCode ierr;
2195baf8537SBarry Smith 
2203a40ed3dSBarry Smith   PetscFunctionBegin;
2215baf8537SBarry Smith   if (snes->nwork) {
2224b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
2235baf8537SBarry Smith   }
224606d414cSSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2253a40ed3dSBarry Smith   PetscFunctionReturn(0);
2264800dd8cSBarry Smith }
2274800dd8cSBarry Smith /*------------------------------------------------------------*/
2284800dd8cSBarry Smith 
2294a2ae208SSatish Balay #undef __FUNCT__
2304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_TR"
2316849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
2324800dd8cSBarry Smith {
2334b27c08aSLois Curfman McInnes   SNES_TR *ctx = (SNES_TR *)snes->data;
234dfbe8321SBarry Smith   PetscErrorCode ierr;
2354800dd8cSBarry Smith 
2363a40ed3dSBarry Smith   PetscFunctionBegin;
237b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
23887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
2394b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
2404b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
2414b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
2424b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
2434b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
2444b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
2454b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
246b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2473a40ed3dSBarry Smith   PetscFunctionReturn(0);
2484800dd8cSBarry Smith }
2494800dd8cSBarry Smith 
2504a2ae208SSatish Balay #undef __FUNCT__
2514b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_TR"
2526849ba73SBarry Smith static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
253a935fc98SLois Curfman McInnes {
2544b27c08aSLois Curfman McInnes   SNES_TR *tr = (SNES_TR *)snes->data;
255dfbe8321SBarry Smith   PetscErrorCode ierr;
25632077d6dSBarry Smith   PetscTruth iascii;
257a935fc98SLois Curfman McInnes 
2583a40ed3dSBarry Smith   PetscFunctionBegin;
25932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
26032077d6dSBarry Smith   if (iascii) {
261a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
262a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
2635cd90555SBarry Smith   } else {
26479a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
26519bcc07fSBarry Smith   }
2663a40ed3dSBarry Smith   PetscFunctionReturn(0);
267a935fc98SLois Curfman McInnes }
268a935fc98SLois Curfman McInnes 
26952392280SLois Curfman McInnes /* ---------------------------------------------------------------- */
2704a2ae208SSatish Balay #undef __FUNCT__
2714b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESConverged_TR"
2728ac28fe4SSatish Balay /*@C
2734b27c08aSLois Curfman McInnes    SNESConverged_TR - Monitors the convergence of the trust region
2744b27c08aSLois Curfman McInnes    method SNESTR for solving systems of nonlinear equations (default).
27552392280SLois Curfman McInnes 
276c7afd0dbSLois Curfman McInnes    Collective on SNES
277c7afd0dbSLois Curfman McInnes 
27852392280SLois Curfman McInnes    Input Parameters:
279c7afd0dbSLois Curfman McInnes +  snes - the SNES context
28052392280SLois Curfman McInnes .  xnorm - 2-norm of current iterate
28152392280SLois Curfman McInnes .  pnorm - 2-norm of current step
28252392280SLois Curfman McInnes .  fnorm - 2-norm of function
283c7afd0dbSLois Curfman McInnes -  dummy - unused context
284fee21e36SBarry Smith 
285184914b5SBarry Smith    Output Parameter:
286184914b5SBarry Smith .   reason - one of
28770441072SBarry Smith $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
288184914b5SBarry Smith $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
289184914b5SBarry Smith $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
290184914b5SBarry Smith $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
291184914b5SBarry Smith $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
292184914b5SBarry Smith $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
293184914b5SBarry Smith $  SNES_CONVERGED_ITERATING       - (otherwise)
29452392280SLois Curfman McInnes 
29552392280SLois Curfman McInnes    where
296c7afd0dbSLois Curfman McInnes +    delta    - trust region paramenter
297c7afd0dbSLois Curfman McInnes .    deltatol - trust region size tolerance,
298c7afd0dbSLois Curfman McInnes                 set with SNESSetTrustRegionTolerance()
299c7afd0dbSLois Curfman McInnes .    maxf - maximum number of function evaluations,
300c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
301c7afd0dbSLois Curfman McInnes .    nfct - number of function evaluations,
30270441072SBarry Smith .    abstol - absolute function norm tolerance,
303c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
304c7afd0dbSLois Curfman McInnes -    xtol - relative function norm tolerance,
305c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
30652392280SLois Curfman McInnes 
30715091d37SBarry Smith    Level: intermediate
30815091d37SBarry Smith 
30952392280SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
31052392280SLois Curfman McInnes 
31152392280SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
31252392280SLois Curfman McInnes @*/
31306ee9f85SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_TR(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
31452392280SLois Curfman McInnes {
3153f149594SLisandro Dalcin   SNES_TR        *neP;
316dfbe8321SBarry Smith   PetscErrorCode ierr;
317ddd12767SBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3193f149594SLisandro Dalcin   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3203f149594SLisandro Dalcin   PetscValidType(snes,1);
3213f149594SLisandro Dalcin   PetscValidPointer(reason,6);
3223f149594SLisandro Dalcin   neP = (SNES_TR *)snes->data;
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) {
3303f149594SLisandro Dalcin     ierr = SNESDefaultConverged(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 
383*38f2d2fdSLisandro Dalcin   ierr			= PetscNewLog(snes,SNES_TR,&neP);CHKERRQ(ierr);
384fbe28522SBarry Smith   snes->data	        = (void*)neP;
385fbe28522SBarry Smith   neP->mu		= 0.25;
386fbe28522SBarry Smith   neP->eta		= 0.75;
387fbe28522SBarry Smith   neP->delta		= 0.0;
388fbe28522SBarry Smith   neP->delta0		= 0.2;
389fbe28522SBarry Smith   neP->delta1		= 0.3;
390fbe28522SBarry Smith   neP->delta2		= 0.75;
391fbe28522SBarry Smith   neP->delta3		= 2.0;
392fbe28522SBarry Smith   neP->sigma		= 0.0001;
393ef87ac0dSKris Buschelman   neP->itflag		= PETSC_FALSE;
39452392280SLois Curfman McInnes   neP->rnorm0		= 0;
39552392280SLois Curfman McInnes   neP->ttol		= 0;
3963a40ed3dSBarry Smith   PetscFunctionReturn(0);
3974800dd8cSBarry Smith }
398fb2e594dSBarry Smith EXTERN_C_END
39982bf6240SBarry Smith 
400