Actual source code: tr.c
1:
2: #include src/snes/impls/tr/tr.h
4: /*
5: This convergence test determines if the two norm of the
6: solution lies outside the trust region, if so it halts.
7: */
10: PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
11: {
12: SNES snes = (SNES) ctx;
13: SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
14: SNES_TR *neP = (SNES_TR*)snes->data;
15: Vec x;
16: PetscReal nrm;
17: PetscErrorCode ierr;
20: if (snes->ksp_ewconv) {
21: if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker onvergence context not created");
22: if (!n) {SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);}
23: kctx->lresid_last = rnorm;
24: }
25: KSPDefaultConverged(ksp,n,rnorm,reason,ctx);
26: if (*reason) {
27: PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: regular convergence test KSP iterations=%D, rnorm=%g\n",n,rnorm);
28: }
30: /* Determine norm of solution */
31: KSPBuildSolution(ksp,0,&x);
32: VecNorm(x,NORM_2,&nrm);
33: if (nrm >= neP->delta) {
34: PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: KSP iterations=%D, rnorm=%g\n",n,rnorm);
35: PetscLogInfo(snes,"SNES_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,nrm);
36: *reason = KSP_CONVERGED_STEP_LENGTH;
37: }
38: return(0);
39: }
41: /*
42: SNESSolve_TR - Implements Newton's Method with a very simple trust
43: region approach for solving systems of nonlinear equations.
45:
46: */
49: static PetscErrorCode SNESSolve_TR(SNES snes)
50: {
51: SNES_TR *neP = (SNES_TR*)snes->data;
52: Vec X,F,Y,G,TMP,Ytmp;
53: PetscErrorCode ierr;
54: PetscInt maxits,i,lits;
55: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
56: PetscReal rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1;
57: PetscScalar mone = -1.0,cnorm;
58: KSP ksp;
59: SNESConvergedReason reason;
60: PetscTruth conv,breakout = PETSC_FALSE;
63: maxits = snes->max_its; /* maximum number of iterations */
64: X = snes->vec_sol; /* solution vector */
65: F = snes->vec_func; /* residual vector */
66: Y = snes->work[0]; /* work vectors */
67: G = snes->work[1];
68: Ytmp = snes->work[2];
70: PetscObjectTakeAccess(snes);
71: snes->iter = 0;
72: PetscObjectGrantAccess(snes);
73: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
75: SNESComputeFunction(snes,X,F); /* F(X) */
76: VecNorm(F,NORM_2,&fnorm); /* fnorm <- || F || */
77: PetscObjectTakeAccess(snes);
78: snes->norm = fnorm;
79: PetscObjectGrantAccess(snes);
80: delta = neP->delta0*fnorm;
81: neP->delta = delta;
82: SNESLogConvHistory(snes,fnorm,0);
83: SNESMonitor(snes,0,fnorm);
84: SNESGetKSP(snes,&ksp);
86: if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}
88: /* set parameter for default relative tolerance convergence test */
89: snes->ttol = fnorm*snes->rtol;
91: /* Set the stopping criteria to use the More' trick. */
92: PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);
93: if (!conv) {
94: KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);
95: PetscLogInfo(snes,"SNESSolve_TR: Using Krylov convergence test SNES_TR_KSPConverged_Private\n");
96: }
97:
98: for (i=0; i<maxits; i++) {
99: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
100: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
102: /* Solve J Y = F, where J is Jacobian matrix */
103: KSPSolve(snes->ksp,F,Ytmp);
104: KSPGetIterationNumber(ksp,&lits);
105: snes->linear_its += lits;
106: PetscLogInfo(snes,"SNESSolve_TR: iter=%D, linear solve iterations=%D\n",snes->iter,lits);
107: VecNorm(Ytmp,NORM_2,&nrm);
108: norm1 = nrm;
109: while(1) {
110: VecCopy(Ytmp,Y);
111: nrm = norm1;
113: /* Scale Y if need be and predict new value of F norm */
114: if (nrm >= delta) {
115: nrm = delta/nrm;
116: gpnorm = (1.0 - nrm)*fnorm;
117: cnorm = nrm;
118: PetscLogInfo(snes,"SNESSolve_TR: Scaling direction by %g\n",nrm);
119: VecScale(&cnorm,Y);
120: nrm = gpnorm;
121: ynorm = delta;
122: } else {
123: gpnorm = 0.0;
124: PetscLogInfo(snes,"SNESSolve_TR: Direction is in Trust Region\n");
125: ynorm = nrm;
126: }
127: VecAYPX(&mone,X,Y); /* Y <- X - Y */
128: VecCopy(X,snes->vec_sol_update_always);
129: SNESComputeFunction(snes,Y,G); /* F(X) */
130: VecNorm(G,NORM_2,&gnorm); /* gnorm <- || g || */
131: if (fnorm == gpnorm) rho = 0.0;
132: else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
134: /* Update size of trust region */
135: if (rho < neP->mu) delta *= neP->delta1;
136: else if (rho < neP->eta) delta *= neP->delta2;
137: else delta *= neP->delta3;
138: PetscLogInfo(snes,"SNESSolve_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm);
139: PetscLogInfo(snes,"SNESSolve_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta);
140: neP->delta = delta;
141: if (rho > neP->sigma) break;
142: PetscLogInfo(snes,"SNESSolve_TR: Trying again in smaller region\n");
143: /* check to see if progress is hopeless */
144: neP->itflag = PETSC_FALSE;
145: (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);
146: if (reason) {
147: /* We're not progressing, so return with the current iterate */
148: SNESMonitor(snes,i+1,fnorm);
149: breakout = PETSC_TRUE;
150: break;
151: }
152: snes->numFailures++;
153: }
154: if (!breakout) {
155: fnorm = gnorm;
156: PetscObjectTakeAccess(snes);
157: snes->iter = i+1;
158: snes->norm = fnorm;
159: PetscObjectGrantAccess(snes);
160: TMP = F; F = G; snes->vec_func_always = F; G = TMP;
161: TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
162: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
163: SNESLogConvHistory(snes,fnorm,lits);
164: SNESMonitor(snes,i+1,fnorm);
166: /* Test for convergence */
167: neP->itflag = PETSC_TRUE;
168: (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);
169: if (reason) {
170: break;
171: }
172: } else {
173: break;
174: }
175: }
176: /* Verify solution is in corect location */
177: if (X != snes->vec_sol) {
178: VecCopy(X,snes->vec_sol);
179: }
180: if (F != snes->vec_func) {
181: VecCopy(F,snes->vec_func);
182: }
183: snes->vec_sol_always = snes->vec_sol;
184: snes->vec_func_always = snes->vec_func;
185: if (i == maxits) {
186: PetscLogInfo(snes,"SNESSolve_TR: Maximum number of iterations has been reached: %D\n",maxits);
187: reason = SNES_DIVERGED_MAX_IT;
188: }
189: PetscObjectTakeAccess(snes);
190: snes->reason = reason;
191: PetscObjectGrantAccess(snes);
192: return(0);
193: }
194: /*------------------------------------------------------------*/
197: static PetscErrorCode SNESSetUp_TR(SNES snes)
198: {
202: snes->nwork = 4;
203: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
204: PetscLogObjectParents(snes,snes->nwork,snes->work);
205: snes->vec_sol_update_always = snes->work[3];
206: return(0);
207: }
208: /*------------------------------------------------------------*/
211: static PetscErrorCode SNESDestroy_TR(SNES snes)
212: {
216: if (snes->nwork) {
217: VecDestroyVecs(snes->work,snes->nwork);
218: }
219: PetscFree(snes->data);
220: return(0);
221: }
222: /*------------------------------------------------------------*/
226: static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
227: {
228: SNES_TR *ctx = (SNES_TR *)snes->data;
232: PetscOptionsHead("SNES trust region options for nonlinear equations");
233: PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);
234: PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);
235: PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);
236: PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);
237: PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);
238: PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);
239: PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);
240: PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);
241: PetscOptionsTail();
242: return(0);
243: }
247: static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
248: {
249: SNES_TR *tr = (SNES_TR *)snes->data;
251: PetscTruth iascii;
254: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
255: if (iascii) {
256: PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);
257: PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);
258: } else {
259: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
260: }
261: return(0);
262: }
264: /* ---------------------------------------------------------------- */
267: /*@C
268: SNESConverged_TR - Monitors the convergence of the trust region
269: method SNESTR for solving systems of nonlinear equations (default).
271: Collective on SNES
273: Input Parameters:
274: + snes - the SNES context
275: . xnorm - 2-norm of current iterate
276: . pnorm - 2-norm of current step
277: . fnorm - 2-norm of function
278: - dummy - unused context
280: Output Parameter:
281: . reason - one of
282: $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol),
283: $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm),
284: $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0),
285: $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf),
286: $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN),
287: $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol),
288: $ SNES_CONVERGED_ITERATING - (otherwise)
290: where
291: + delta - trust region paramenter
292: . deltatol - trust region size tolerance,
293: set with SNESSetTrustRegionTolerance()
294: . maxf - maximum number of function evaluations,
295: set with SNESSetTolerances()
296: . nfct - number of function evaluations,
297: . abstol - absolute function norm tolerance,
298: set with SNESSetTolerances()
299: - xtol - relative function norm tolerance,
300: set with SNESSetTolerances()
302: Level: intermediate
304: .keywords: SNES, nonlinear, default, converged, convergence
306: .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
307: @*/
308: PetscErrorCode SNESConverged_TR(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
309: {
310: SNES_TR *neP = (SNES_TR *)snes->data;
314: if (fnorm != fnorm) {
315: PetscLogInfo(snes,"SNESConverged_TR:Failed to converged, function norm is NaN\n");
316: *reason = SNES_DIVERGED_FNORM_NAN;
317: } else if (neP->delta < xnorm * snes->deltatol) {
318: PetscLogInfo(snes,"SNESConverged_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol);
319: *reason = SNES_CONVERGED_TR_DELTA;
320: } else if (neP->itflag) {
321: SNESConverged_LS(snes,xnorm,pnorm,fnorm,reason,dummy);
322: } else if (snes->nfuncs > snes->max_funcs) {
323: PetscLogInfo(snes,"SNESConverged_TR: Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
324: *reason = SNES_DIVERGED_FUNCTION_COUNT;
325: } else {
326: *reason = SNES_CONVERGED_ITERATING;
327: }
328: return(0);
329: }
330: /* ------------------------------------------------------------ */
331: /*MC
332: SNESTR - Newton based nonlinear solver that uses a trust region
334: Options Database:
335: + -snes_trtol <tol> Trust region tolerance
336: . -snes_tr_mu <mu>
337: . -snes_tr_eta <eta>
338: . -snes_tr_sigma <sigma>
339: . -snes_tr_delta0 <delta0>
340: . -snes_tr_delta1 <delta1>
341: . -snes_tr_delta2 <delta2>
342: - -snes_tr_delta3 <delta3>
344: The basic algorithm is taken from "The Minpack Project", by More',
345: Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
346: of Mathematical Software", Wayne Cowell, editor.
348: This is intended as a model implementation, since it does not
349: necessarily have many of the bells and whistles of other
350: implementations.
352: Level: intermediate
354: .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
356: M*/
360: PetscErrorCode SNESCreate_TR(SNES snes)
361: {
362: SNES_TR *neP;
366: snes->setup = SNESSetUp_TR;
367: snes->solve = SNESSolve_TR;
368: snes->destroy = SNESDestroy_TR;
369: snes->converged = SNESConverged_TR;
370: snes->setfromoptions = SNESSetFromOptions_TR;
371: snes->view = SNESView_TR;
372: snes->nwork = 0;
373:
374: ierr = PetscNew(SNES_TR,&neP);
375: PetscLogObjectMemory(snes,sizeof(SNES_TR));
376: snes->data = (void*)neP;
377: neP->mu = 0.25;
378: neP->eta = 0.75;
379: neP->delta = 0.0;
380: neP->delta0 = 0.2;
381: neP->delta1 = 0.3;
382: neP->delta2 = 0.75;
383: neP->delta3 = 2.0;
384: neP->sigma = 0.0001;
385: neP->itflag = PETSC_FALSE;
386: neP->rnorm0 = 0;
387: neP->ttol = 0;
388: return(0);
389: }