Actual source code: umtr.c
1: /*$Id: umtr.c,v 1.114 2001/09/04 20:33:25 buschelm Exp $*/
3: #include src/snes/impls/umtr/umtr.h
5: /*
6: SNESSolve_UM_TR - Implements Newton's Method with a trust region approach
7: for solving unconstrained minimization problems.
9: The basic algorithm is taken from MINPACK-2 (dstrn).
11: SNESSolve_UM_TR computes a local minimizer of a twice differentiable function
12: f by applying a trust region variant of Newton's method. At each stage
13: of the algorithm, we us the prconditioned conjugate gradient method to
14: determine an approximate minimizer of the quadratic equation
16: q(s) = g^T * s + .5 * s^T * H * s
18: subject to the Euclidean norm trust region constraint
20: || D * s || <= delta,
22: where delta is the trust region radius and D is a scaling matrix.
23: Here g is the gradient and H is the Hessian matrix.
25: Note: SNESSolve_UM_TR MUST use the iterative solver KSPQCG; thus, we
26: set KSPQCG in this routine regardless of what the user may have
27: previously specified.
28: */
29: static int SNESSolve_UM_TR(SNES snes,int *outits)
30: {
31: SNES_UM_TR *neP = (SNES_UM_TR*)snes->data;
32: int maxits,i,nlconv,ierr,qits;
33: PetscTruth newton;
34: PetscReal xnorm,max_val,ftrial,delta,ltsnrm,quadratic;
35: PetscReal zero = 0.0,two = 2.0,four = 4.0;
36: PetscScalar one = 1.0;
37: Vec X,G,S,Xtrial;
38: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
39: SLES sles;
40: KSP ksp;
41: SNESConvergedReason reason;
42: KSPConvergedReason kreason;
46: snes->reason = SNES_CONVERGED_ITERATING;
48: nlconv = 0;
49: maxits = snes->max_its; /* maximum number of iterations */
50: X = snes->vec_sol; /* solution vector */
51: G = snes->vec_func; /* gradient vector */
52: S = snes->work[0]; /* work vectors */
53: Xtrial = snes->work[1];
54: delta = neP->delta; /* trust region radius */
56: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
57: PetscObjectTakeAccess(snes);
58: snes->iter = 0;
59: PetscObjectGrantAccess(snes);
60: SNESComputeMinimizationFunction(snes,X,&snes->fc); /* f(X) */
61: SNESComputeGradient(snes,X,G); /* G(X) <- gradient */
62: PetscObjectTakeAccess(snes);
63: VecNorm(G,NORM_2,&snes->norm); /* &snes->norm = || G || */
64: PetscObjectGrantAccess(snes);
65: SNESLogConvHistory(snes,snes->norm,0);
66: SNESMonitor(snes,0,snes->norm);
68: SNESGetSLES(snes,&sles);
69: SLESGetKSP(sles,&ksp);
70: KSPSetType(ksp,KSPQCG);
71: PetscLogInfo(snes,"SNESSolve_UM_TR: setting KSPType = KSPQCGn");
73: for (i=0; i<maxits && !nlconv; i++) {
74: PetscObjectTakeAccess(snes);
75: snes->iter = i+1;
76: PetscObjectGrantAccess(snes);
77: newton = PETSC_FALSE;
78: neP->success = 0;
79: snes->numFailures = 0;
80: SNESComputeHessian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
81: SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);
83: if (!i) { /* Initialize delta */
84: if (delta <= 0) {
85: if (xnorm > zero) delta = neP->factor1*xnorm;
86: else delta = neP->delta0;
87: MatNorm(snes->jacobian,NORM_1,&max_val);
88: if (ierr == PETSC_ERR_SUP) {
89: PetscLogInfo(snes,"SNESSolve_UM_TR: Initial delta computed without matrix norm infon");
90: } else if (ierr == 0) {
91: if (PetscAbsReal(max_val)<1.e-14)SETERRQ(PETSC_ERR_PLIB,"Hessian norm is too small");
92: delta = PetscMax(delta,snes->norm/max_val);
93: } else {
94:
95: }
96: } else {
97: delta = neP->delta0;
98: }
99: }
100: do {
101: /* Minimize the quadratic to compute the step s */
102: KSPQCGSetTrustRegionRadius(ksp,delta);
104: SLESSolve(snes->sles,G,S,&qits);
105: snes->linear_its += qits;
106: KSPQCGGetTrialStepNorm(ksp,<snrm);
107: KSPQCGGetQuadratic(ksp,&quadratic);
108: KSPGetConvergedReason(ksp,&kreason);
109: if ((int)kreason < 0) SETERRQ(PETSC_ERR_PLIB,"Failure in SLESSolve");
110: if (kreason != KSP_CONVERGED_QCG_NEG_CURVE && kreason != KSP_CONVERGED_QCG_CONSTRAINED) {
111: newton = PETSC_TRUE;
112: }
113: PetscLogInfo(snes,"SNESSolve_UM_TR: %d: ltsnrm=%g, delta=%g, q=%g, qits=%dn",
114: i,ltsnrm,delta,quadratic,qits);
116: VecWAXPY(&one,X,S,Xtrial); /* Xtrial <- X + S */
117: VecNorm(Xtrial,NORM_2,&xnorm);
119: /* ftrial = f(Xtrial) */
120: SNESComputeMinimizationFunction(snes,Xtrial,&ftrial);
122: /* Compute the function reduction and the step size */
123: neP->prered = -quadratic;
124: neP->actred = snes->fc - ftrial;
126: /* Adjust delta for the first Newton step */
127: if (!i && (newton)) delta = PetscMin(delta,ltsnrm);
129: if (neP->actred < neP->eta1 * neP->prered) { /* Unsuccessful step */
131: PetscLogInfo(snes,"SNESSolve_UM_TR: Rejecting stepn");
132: snes->numFailures += 1;
134: /* If iterate is Newton step, reduce delta to current step length */
135: if (newton) {
136: delta = ltsnrm;
137: newton = PETSC_FALSE;
138: }
139: delta /= four;
141: } else { /* Successful iteration; adjust trust radius */
143: neP->success = 1;
144: PetscLogInfo(snes,"SNESSolve_UM_TR: Accepting stepn");
145: if (newton) {
146: delta = sqrt(ltsnrm*delta);
147: if (neP->actred < neP->eta2 * neP->prered) delta /= two;
148: } else if (neP->actred < neP->eta2 * neP->prered)
149: delta /= delta;
150: else if ((neP->actred >= neP->eta3 * neP->prered) &&
151: (neP->actred < neP->eta4 * neP->prered))
152: delta *= two;
153: else if (neP->actred >= neP->eta4 * neP->prered)
154: delta *= four;
155: else neP->sflag = 1;
156: }
158: neP->delta = delta;
159: (*snes->converged)(snes,xnorm,snes->norm,ftrial,&reason,snes->cnvP);
160: if (reason) nlconv = 1;
161: } while (!neP->success && !nlconv);
163: /* Question: If (!neP->success && break), then last step was rejected,
164: but convergence was detected. Should this update really happen? */
165: snes->fc = ftrial;
166: VecCopy(Xtrial,X);
167: snes->vec_sol_always = X;
168: /* Note: At last iteration, the gradient evaluation is unnecessary */
169: SNESComputeGradient(snes,X,G);
170: PetscObjectTakeAccess(snes);
171: VecNorm(G,NORM_2,&snes->norm);
172: PetscObjectGrantAccess(snes);
173: snes->vec_func_always = G;
175: SNESLogConvHistory(snes,snes->norm,qits);
176: SNESMonitor(snes,i+1,snes->norm);
177: }
178: /* Verify solution is in corect location */
179: if (X != snes->vec_sol) {
180: VecCopy(X,snes->vec_sol);
181: snes->vec_sol_always = snes->vec_sol;
182: snes->vec_func_always = snes->vec_func;
183: }
184: if (i == maxits) {
185: PetscLogInfo(snes,"SNESSolve_UM_TR: Maximum number of iterations reached: %dn",maxits);
186: i--;
187: reason = SNES_DIVERGED_MAX_IT;
188: }
189: PetscObjectTakeAccess(snes);
190: snes->reason = reason;
191: PetscObjectGrantAccess(snes);
192: *outits = i; /* not i+1, since update for i happens in loop above */
193: return(0);
194: }
195: /*------------------------------------------------------------*/
196: static int SNESSetUp_UM_TR(SNES snes)
197: {
198: int ierr;
199: PetscTruth ilu,bjacobi;
200: SLES sles;
201: PC pc;
204: snes->nwork = 4;
205: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
206: PetscLogObjectParents(snes,snes->nwork,snes->work);
207: snes->vec_sol_update_always = snes->work[3];
209: /*
210: If PC was set by default to ILU then change it to Jacobi
211: */
212: SNESGetSLES(snes,&sles);
213: SLESGetPC(sles,&pc);
214: PetscTypeCompare((PetscObject)pc,PCILU,&ilu);
215: if (ilu) {
216: PCSetType(pc,PCJACOBI);
217: } else {
218: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&bjacobi);
219: if (bjacobi) {
220: /* cannot do this; since PC may not have been setup yet */
221: /* PCBJacobiGetSubSLES(pc,0,0,&subsles);
222: SLESGetPC(*subsles,&subpc);
223: PetscTypeCompare((PetscObject)subpc,PCILU,&ilu);
224: if (ilu) {
225: PCSetType(pc,PCJACOBI);
226: } */
227: /* don't really want to do this, since user may have selected BJacobi plus something
228: that is symmetric on each processor; really only want to catch the default ILU */
229: PCSetType(pc,PCJACOBI);
230: }
231: }
232: return(0);
233: }
234: /*------------------------------------------------------------*/
235: static int SNESDestroy_UM_TR(SNES snes)
236: {
237: int ierr;
240: if (snes->nwork) {
241: VecDestroyVecs(snes->work,snes->nwork);
242: }
243: PetscFree(snes->data);
244: return(0);
245: }
246: /*------------------------------------------------------------*/
247: /*@C
248: SNESConverged_UM_TR - Monitors the convergence of the SNESSolve_UM_TR()
249: routine (default).
251: Collective on SNES
253: Input Parameters:
254: + snes - the SNES context
255: . xnorm - 2-norm of current iterate
256: . gnorm - 2-norm of current gradient
257: . f - objective function value
258: - dummy - unused dummy context
260: Output Parameter:
261: . reason - one of
262: $ SNES_CONVERGED_FNORM_ABS (f < fmin),
263: $ SNES_CONVERGED_TR_REDUCTION (abs(ared) <= rtol*abs(f) && pred <= rtol*abs(f)),
264: $ SNES_CONVERGED_TR_DELTA (delta <= deltatol*xnorm),
265: $ SNES_DIVERGED_TR_REDUCTION (abs(ared) <= epsmch && pred <= epsmch),
266: $ SNES_DIVERGED_FUNCTION_COUNT (nfunc > max_func),
267: $ SNES_DIVERGED_FNORM_NAN (f = NaN),
268: $ SNES_CONVERGED_ITERATING (otherwise).
270: where
271: + ared - actual reduction
272: . delta - trust region paramenter
273: . deltatol - trust region size tolerance,
274: set with SNESSetTrustRegionTolerance()
275: . epsmch - machine epsilon
276: . fmin - lower bound on function value,
277: set with SNESSetMinimizationFunctionTolerance()
278: . nfunc - number of function evaluations
279: . maxfunc - maximum number of function evaluations,
280: set with SNESSetTolerances()
281: . pred - predicted reduction
282: - rtol - relative function tolerance,
283: set with SNESSetTolerances()
285: Level: intermediate
287: @*/
288: int SNESConverged_UM_TR(SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *dummy)
289: {
290: SNES_UM_TR *neP = (SNES_UM_TR*)snes->data;
291: PetscReal rtol = snes->rtol,delta = neP->delta,ared = neP->actred,pred = neP->prered;
292: PetscReal epsmch = 1.0e-14; /* This must be fixed */
295: if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
296: SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
297: } else if (f != f) {
298: PetscLogInfo(snes,"SNESConverged_UM_TR:Failed to converged, function is NaNn");
299: *reason = SNES_DIVERGED_FNORM_NAN;
300: } else if ((!neP->success || neP->sflag) && (delta <= snes->deltatol * xnorm)) {
301: neP->sflag = 0;
302: PetscLogInfo(snes,"SNESConverged_UM_TR: Trust region param satisfies tolerance: %g<=%g*%gn",
303: delta,snes->deltatol,xnorm);
304: *reason = SNES_CONVERGED_TR_DELTA;
305: } else if ((PetscAbsReal(ared) <= PetscAbsReal(f) * rtol) && (pred) <= rtol*PetscAbsReal(f)) {
306: PetscLogInfo(snes,"SNESConverged_UM_TR:Actual (%g) and predicted (%g) reductions<%g*%gn",
307: PetscAbsReal(ared),pred,rtol,PetscAbsReal(f));
308: *reason = SNES_CONVERGED_TR_REDUCTION;
309: } else if (f < snes->fmin) {
310: PetscLogInfo(snes,"SNESConverged_UM_TR:Function value (%g)<f_{minimum} (%g)n",f,snes->fmin);
311: *reason = SNES_CONVERGED_FNORM_ABS ;
312: } else if ((PetscAbsReal(ared) <= epsmch) && pred <= epsmch) {
313: PetscLogInfo(snes,"SNESConverged_UM_TR:Actual (%g) and predicted (%g) reductions<epsmch (%g)n",
314: PetscAbsReal(ared),pred,epsmch);
315: *reason = SNES_DIVERGED_TR_REDUCTION;
316: } else if (snes->nfuncs > snes->max_funcs) {
317: PetscLogInfo(snes,"SNESConverged_UM_TR:Exceeded maximum number of function evaluations:%d>%dn",
318: snes->nfuncs,snes->max_funcs);
319: *reason = SNES_DIVERGED_FUNCTION_COUNT;
320: } else {
321: *reason = SNES_CONVERGED_ITERATING;
322: }
323: return(0);
324: }
325: /*------------------------------------------------------------*/
326: static int SNESSetFromOptions_UM_TR(SNES snes)
327: {
328: SNES_UM_TR *ctx = (SNES_UM_TR *)snes->data;
329: int ierr;
330: SLES sles;
331: PC pc;
332: PetscTruth ismatshell,nopcset;
333: Mat pmat;
336: PetscOptionsHead("SNES trust region options for minimization");
337: PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);
338: PetscOptionsReal("-snes_um_eta1","eta1","None",ctx->eta1,&ctx->eta1,0);
339: PetscOptionsReal("-snes_um_eta2","step unsuccessful if reduction < eta1 * predicted reduction","None",ctx->eta2,&ctx->eta2,0);
340: PetscOptionsReal("-snes_um_eta3","eta3","None",ctx->eta3,&ctx->eta3,0);
341: PetscOptionsReal("-snes_um_eta4","eta4","None",ctx->eta4,&ctx->eta4,0);
342: PetscOptionsReal("-snes_um_delta0","delta0","None",ctx->delta,&ctx->delta,0);
343: PetscOptionsReal("-snes_um_factor1","factor1","None",ctx->factor1,&ctx->factor1,0);
344: PetscOptionsTail();
346: /* if preconditioner has not been set yet, and not using a matrix shell then
347: set preconditioner to Jacobi. This is to prevent PCSetFromOptions() from
348: setting a default of ILU or block Jacobi-ILU which won't work since TR
349: requires a symmetric preconditioner
350: */
351: SNESGetSLES(snes,&sles);
352: SLESGetPC(sles,&pc);
353: PetscTypeCompare((PetscObject)pc,0,&nopcset);
354: if (nopcset) {
355: PCGetOperators(pc,PETSC_NULL,&pmat,PETSC_NULL);
356: if (pmat) {
357: PetscTypeCompare((PetscObject)pmat,MATSHELL,&ismatshell);
358: if (!ismatshell) {
359: PCSetType(pc,PCJACOBI);
360: }
361: }
362: }
364: return(0);
365: }
367: /*------------------------------------------------------------*/
368: static int SNESView_UM_TR(SNES snes,PetscViewer viewer)
369: {
370: SNES_UM_TR *tr = (SNES_UM_TR *)snes->data;
371: int ierr;
372: PetscTruth isascii;
375: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
376: if (isascii) {
377: PetscViewerASCIIPrintf(viewer," eta1=%g, eta1=%g, eta3=%g, eta4=%gn",tr->eta1,tr->eta2,tr->eta3,tr->eta4);
378: PetscViewerASCIIPrintf(viewer," delta0=%g, factor1=%gn",tr->delta0,tr->factor1);
379: } else {
380: SETERRQ1(1,"Viewer type %s not supported for SNES UM TR",((PetscObject)viewer)->type_name);
381: }
382: return(0);
383: }
384: /*------------------------------------------------------------*/
385: EXTERN_C_BEGIN
386: int SNESCreate_UM_TR(SNES snes)
387: {
388: SNES_UM_TR *neP;
389: int ierr;
392: if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
393: SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
394: }
395: snes->setup = SNESSetUp_UM_TR;
396: snes->solve = SNESSolve_UM_TR;
397: snes->destroy = SNESDestroy_UM_TR;
398: snes->converged = SNESConverged_UM_TR;
399: snes->setfromoptions = SNESSetFromOptions_UM_TR;
400: snes->view = SNESView_UM_TR;
402: snes->nwork = 0;
404: ierr = PetscNew(SNES_UM_TR,&neP);
405: PetscLogObjectMemory(snes,sizeof(SNES_UM_TR));
406: snes->data = (void*)neP;
407: neP->delta0 = 1.0e-6;
408: neP->delta = 0.0;
409: neP->eta1 = 1.0e-4;
410: neP->eta2 = 0.25;
411: neP->eta3 = 0.50;
412: neP->eta4 = 0.90;
413: neP->factor1 = 1.0e-8;
414: neP->actred = 0.0;
415: neP->prered = 0.0;
416: neP->success = 0;
417: neP->sflag = 0;
418:
419: return(0);
420: }
421: EXTERN_C_END