Actual source code: tr.c
1: /*$Id: tr.c,v 1.123 2001/03/23 23:24:15 balay Exp $*/
3: #include src/snes/impls/tr/tr.h
5: /*
6: This convergence test determines if the two norm of the
7: solution lies outside the trust region, if so it halts.
8: */
9: int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n,double rnorm,KSPConvergedReason *reason,void *ctx)
10: {
11: SNES snes = (SNES) ctx;
12: SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
13: SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data;
14: Vec x;
15: double norm;
16: int ierr;
19: if (snes->ksp_ewconv) {
20: if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker onvergence context not created");
21: if (!n) {SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);}
22: kctx->lresid_last = rnorm;
23: }
24: KSPDefaultConverged(ksp,n,rnorm,reason,ctx);
25: if (*reason) {
26: PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: regular convergence test KSP iterations=%d, rnorm=%gn",n,rnorm);
27: }
29: /* Determine norm of solution */
30: KSPBuildSolution(ksp,0,&x);
31: VecNorm(x,NORM_2,&norm);
32: if (norm >= neP->delta) {
33: PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%gn",n,rnorm);
34: PetscLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%gn",neP->delta,norm);
35: *reason = KSP_CONVERGED_STEP_LENGTH;
36: }
37: return(0);
38: }
40: /*
41: SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust
42: region approach for solving systems of nonlinear equations.
44: The basic algorithm is taken from "The Minpack Project", by More',
45: Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
46: of Mathematical Software", Wayne Cowell, editor.
48: This is intended as a model implementation, since it does not
49: necessarily have many of the bells and whistles of other
50: implementations.
51: */
52: static int SNESSolve_EQ_TR(SNES snes,int *its)
53: {
54: SNES_EQ_TR *neP = (SNES_EQ_TR*)snes->data;
55: Vec X,F,Y,G,TMP,Ytmp;
56: int maxits,i,ierr,lits,breakout = 0;
57: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
58: double rho,fnorm,gnorm,gpnorm,xnorm,delta,norm,ynorm,norm1;
59: Scalar mone = -1.0,cnorm;
60: KSP ksp;
61: SLES sles;
62: SNESConvergedReason reason;
65: maxits = snes->max_its; /* maximum number of iterations */
66: X = snes->vec_sol; /* solution vector */
67: F = snes->vec_func; /* residual vector */
68: Y = snes->work[0]; /* work vectors */
69: G = snes->work[1];
70: Ytmp = snes->work[2];
72: PetscObjectTakeAccess(snes);
73: snes->iter = 0;
74: PetscObjectGrantAccess(snes);
75: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
77: SNESComputeFunction(snes,X,F); /* F(X) */
78: VecNorm(F,NORM_2,&fnorm); /* fnorm <- || F || */
79: PetscObjectTakeAccess(snes);
80: snes->norm = fnorm;
81: PetscObjectGrantAccess(snes);
82: delta = neP->delta0*fnorm;
83: neP->delta = delta;
84: SNESLogConvHistory(snes,fnorm,0);
85: SNESMonitor(snes,0,fnorm);
87: if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}
89: /* set parameter for default relative tolerance convergence test */
90: snes->ttol = fnorm*snes->rtol;
92: /* Set the stopping criteria to use the More' trick. */
93: SNESGetSLES(snes,&sles);
94: SLESGetKSP(sles,&ksp);
95: KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);
96:
97: for (i=0; i<maxits; i++) {
98: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
99: SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);
101: /* Solve J Y = F, where J is Jacobian matrix */
102: SLESSolve(snes->sles,F,Ytmp,&lits);
103: snes->linear_its += lits;
104: PetscLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%dn",snes->iter,lits);
105: VecNorm(Ytmp,NORM_2,&norm);
106: norm1 = norm;
107: while(1) {
108: VecCopy(Ytmp,Y);
109: norm = norm1;
111: /* Scale Y if need be and predict new value of F norm */
112: if (norm >= delta) {
113: norm = delta/norm;
114: gpnorm = (1.0 - norm)*fnorm;
115: cnorm = norm;
116: PetscLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %gn",norm);
117: VecScale(&cnorm,Y);
118: norm = gpnorm;
119: ynorm = delta;
120: } else {
121: gpnorm = 0.0;
122: PetscLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Regionn");
123: ynorm = norm;
124: }
125: VecAYPX(&mone,X,Y); /* Y <- X - Y */
126: VecCopy(X,snes->vec_sol_update_always);
127: SNESComputeFunction(snes,Y,G); /* F(X) */
128: VecNorm(G,NORM_2,&gnorm); /* gnorm <- || g || */
129: if (fnorm == gpnorm) rho = 0.0;
130: else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
132: /* Update size of trust region */
133: if (rho < neP->mu) delta *= neP->delta1;
134: else if (rho < neP->eta) delta *= neP->delta2;
135: else delta *= neP->delta3;
136: PetscLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%gn",fnorm,gnorm,ynorm);
137: PetscLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%gn",gpnorm,rho,delta);
138: neP->delta = delta;
139: if (rho > neP->sigma) break;
140: PetscLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller regionn");
141: /* check to see if progress is hopeless */
142: neP->itflag = 0;
143: (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);
144: if (reason) {
145: /* We're not progressing, so return with the current iterate */
146: breakout = 1;
147: break;
148: }
149: snes->nfailures++;
150: }
151: if (!breakout) {
152: fnorm = gnorm;
153: PetscObjectTakeAccess(snes);
154: snes->iter = i+1;
155: snes->norm = fnorm;
156: PetscObjectGrantAccess(snes);
157: TMP = F; F = G; snes->vec_func_always = F; G = TMP;
158: TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
159: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
160: SNESLogConvHistory(snes,fnorm,lits);
161: SNESMonitor(snes,i+1,fnorm);
163: /* Test for convergence */
164: neP->itflag = 1;
165: (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);
166: if (reason) {
167: break;
168: }
169: } else {
170: break;
171: }
172: }
173: /* Verify solution is in corect location */
174: if (X != snes->vec_sol) {
175: VecCopy(X,snes->vec_sol);
176: }
177: if (F != snes->vec_func) {
178: VecCopy(F,snes->vec_func);
179: }
180: snes->vec_sol_always = snes->vec_sol;
181: snes->vec_func_always = snes->vec_func;
182: if (i == maxits) {
183: PetscLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %dn",maxits);
184: i--;
185: reason = SNES_DIVERGED_MAX_IT;
186: }
187: *its = i+1;
188: PetscObjectTakeAccess(snes);
189: snes->reason = reason;
190: PetscObjectGrantAccess(snes);
191: return(0);
192: }
193: /*------------------------------------------------------------*/
194: static int SNESSetUp_EQ_TR(SNES snes)
195: {
199: snes->nwork = 4;
200: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
201: PetscLogObjectParents(snes,snes->nwork,snes->work);
202: snes->vec_sol_update_always = snes->work[3];
203: return(0);
204: }
205: /*------------------------------------------------------------*/
206: static int SNESDestroy_EQ_TR(SNES snes)
207: {
208: int ierr;
211: if (snes->nwork) {
212: VecDestroyVecs(snes->work,snes->nwork);
213: }
214: PetscFree(snes->data);
215: return(0);
216: }
217: /*------------------------------------------------------------*/
219: static int SNESSetFromOptions_EQ_TR(SNES snes)
220: {
221: SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data;
222: int ierr;
225: PetscOptionsHead("SNES trust region options for nonlinear equations");
226: PetscOptionsDouble("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);
227: PetscOptionsDouble("-snes_eq_tr_mu","mu","None",ctx->mu,&ctx->mu,0);
228: PetscOptionsDouble("-snes_eq_tr_eta","eta","None",ctx->eta,&ctx->eta,0);
229: PetscOptionsDouble("-snes_eq_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);
230: PetscOptionsDouble("-snes_eq_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);
231: PetscOptionsDouble("-snes_eq_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);
232: PetscOptionsDouble("-snes_eq_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);
233: PetscOptionsDouble("-snes_eq_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);
234: PetscOptionsTail();
235: return(0);
236: }
238: static int SNESView_EQ_TR(SNES snes,PetscViewer viewer)
239: {
240: SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data;
241: int ierr;
242: PetscTruth isascii;
245: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
246: if (isascii) {
247: PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%gn",tr->mu,tr->eta,tr->sigma);
248: PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%gn",tr->delta0,tr->delta1,tr->delta2,tr->delta3);
249: } else {
250: SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
251: }
252: return(0);
253: }
255: /* ---------------------------------------------------------------- */
256: /*@C
257: SNESConverged_EQ_TR - Monitors the convergence of the trust region
258: method SNESEQTR for solving systems of nonlinear equations (default).
260: Collective on SNES
262: Input Parameters:
263: + snes - the SNES context
264: . xnorm - 2-norm of current iterate
265: . pnorm - 2-norm of current step
266: . fnorm - 2-norm of function
267: - dummy - unused context
269: Output Parameter:
270: . reason - one of
271: $ SNES_CONVERGED_FNORM_ABS - (fnorm < atol),
272: $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm),
273: $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0),
274: $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf),
275: $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN),
276: $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol),
277: $ SNES_CONVERGED_ITERATING - (otherwise)
279: where
280: + delta - trust region paramenter
281: . deltatol - trust region size tolerance,
282: set with SNESSetTrustRegionTolerance()
283: . maxf - maximum number of function evaluations,
284: set with SNESSetTolerances()
285: . nfct - number of function evaluations,
286: . atol - absolute function norm tolerance,
287: set with SNESSetTolerances()
288: - xtol - relative function norm tolerance,
289: set with SNESSetTolerances()
291: Level: intermediate
293: .keywords: SNES, nonlinear, default, converged, convergence
295: .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
296: @*/
297: int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy)
298: {
299: SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data;
300: int ierr;
303: if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
304: SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
305: }
307: if (fnorm != fnorm) {
308: PetscLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaNn");
309: *reason = SNES_DIVERGED_FNORM_NAN;
310: } else if (neP->delta < xnorm * snes->deltatol) {
311: PetscLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%gn",neP->delta,xnorm,snes->deltatol);
312: *reason = SNES_CONVERGED_TR_DELTA;
313: } else if (neP->itflag) {
314: SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);
315: } else if (snes->nfuncs > snes->max_funcs) {
316: PetscLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %dn",snes->nfuncs,snes->max_funcs);
317: *reason = SNES_DIVERGED_FUNCTION_COUNT;
318: } else {
319: *reason = SNES_CONVERGED_ITERATING;
320: }
321: return(0);
322: }
323: /* ------------------------------------------------------------ */
324: EXTERN_C_BEGIN
325: int SNESCreate_EQ_TR(SNES snes)
326: {
327: SNES_EQ_TR *neP;
331: if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
332: SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
333: }
334: snes->setup = SNESSetUp_EQ_TR;
335: snes->solve = SNESSolve_EQ_TR;
336: snes->destroy = SNESDestroy_EQ_TR;
337: snes->converged = SNESConverged_EQ_TR;
338: snes->setfromoptions = SNESSetFromOptions_EQ_TR;
339: snes->view = SNESView_EQ_TR;
340: snes->nwork = 0;
341:
342: ierr = PetscNew(SNES_EQ_TR,&neP);
343: PetscLogObjectMemory(snes,sizeof(SNES_EQ_TR));
344: snes->data = (void*)neP;
345: neP->mu = 0.25;
346: neP->eta = 0.75;
347: neP->delta = 0.0;
348: neP->delta0 = 0.2;
349: neP->delta1 = 0.3;
350: neP->delta2 = 0.75;
351: neP->delta3 = 2.0;
352: neP->sigma = 0.0001;
353: neP->itflag = 0;
354: neP->rnorm0 = 0;
355: neP->ttol = 0;
356: return(0);
357: }
358: EXTERN_C_END