Actual source code: taocg.c
petsc-3.9.4 2018-09-11
1: #include <petsctaolinesearch.h>
2: #include <../src/tao/unconstrained/impls/cg/taocg.h>
4: #define CG_FletcherReeves 0
5: #define CG_PolakRibiere 1
6: #define CG_PolakRibierePlus 2
7: #define CG_HestenesStiefel 3
8: #define CG_DaiYuan 4
9: #define CG_Types 5
11: static const char *CG_Table[64] = {"fr", "pr", "prp", "hs", "dy"};
13: static PetscErrorCode TaoSolve_CG(Tao tao)
14: {
15: TAO_CG *cgP = (TAO_CG*)tao->data;
16: PetscErrorCode ierr;
17: TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
18: PetscReal step=1.0,f,gnorm,gnorm2,delta,gd,ginner,beta;
19: PetscReal gd_old,gnorm2_old,f_old;
22: if (tao->XL || tao->XU || tao->ops->computebounds) {
23: PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by cg algorithm\n");
24: }
26: /* Check convergence criteria */
27: TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);
28: VecNorm(tao->gradient,NORM_2,&gnorm);
29: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
30:
31: tao->reason = TAO_CONTINUE_ITERATING;
32: TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
33: TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);
34: (*tao->ops->convergencetest)(tao,tao->cnvP);
35: if (tao->reason != TAO_CONTINUE_ITERATING) return(0);
37: /* Set initial direction to -gradient */
38: VecCopy(tao->gradient, tao->stepdirection);
39: VecScale(tao->stepdirection, -1.0);
40: gnorm2 = gnorm*gnorm;
42: /* Set initial scaling for the function */
43: if (f != 0.0) {
44: delta = 2.0*PetscAbsScalar(f) / gnorm2;
45: delta = PetscMax(delta,cgP->delta_min);
46: delta = PetscMin(delta,cgP->delta_max);
47: } else {
48: delta = 2.0 / gnorm2;
49: delta = PetscMax(delta,cgP->delta_min);
50: delta = PetscMin(delta,cgP->delta_max);
51: }
52: /* Set counter for gradient and reset steps */
53: cgP->ngradsteps = 0;
54: cgP->nresetsteps = 0;
56: while (1) {
57: /* Save the current gradient information */
58: f_old = f;
59: gnorm2_old = gnorm2;
60: VecCopy(tao->solution, cgP->X_old);
61: VecCopy(tao->gradient, cgP->G_old);
62: VecDot(tao->gradient, tao->stepdirection, &gd);
63: if ((gd >= 0) || PetscIsInfOrNanReal(gd)) {
64: ++cgP->ngradsteps;
65: if (f != 0.0) {
66: delta = 2.0*PetscAbsScalar(f) / gnorm2;
67: delta = PetscMax(delta,cgP->delta_min);
68: delta = PetscMin(delta,cgP->delta_max);
69: } else {
70: delta = 2.0 / gnorm2;
71: delta = PetscMax(delta,cgP->delta_min);
72: delta = PetscMin(delta,cgP->delta_max);
73: }
75: VecCopy(tao->gradient, tao->stepdirection);
76: VecScale(tao->stepdirection, -1.0);
77: }
79: /* Search direction for improving point */
80: TaoLineSearchSetInitialStepLength(tao->linesearch,delta);
81: TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);
82: TaoAddLineSearchCounts(tao);
83: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
84: /* Linesearch failed */
85: /* Reset factors and use scaled gradient step */
86: ++cgP->nresetsteps;
87: f = f_old;
88: gnorm2 = gnorm2_old;
89: VecCopy(cgP->X_old, tao->solution);
90: VecCopy(cgP->G_old, tao->gradient);
92: if (f != 0.0) {
93: delta = 2.0*PetscAbsScalar(f) / gnorm2;
94: delta = PetscMax(delta,cgP->delta_min);
95: delta = PetscMin(delta,cgP->delta_max);
96: } else {
97: delta = 2.0 / gnorm2;
98: delta = PetscMax(delta,cgP->delta_min);
99: delta = PetscMin(delta,cgP->delta_max);
100: }
102: VecCopy(tao->gradient, tao->stepdirection);
103: VecScale(tao->stepdirection, -1.0);
105: TaoLineSearchSetInitialStepLength(tao->linesearch,delta);
106: TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);
107: TaoAddLineSearchCounts(tao);
109: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
110: /* Linesearch failed again */
111: /* switch to unscaled gradient */
112: f = f_old;
113: VecCopy(cgP->X_old, tao->solution);
114: VecCopy(cgP->G_old, tao->gradient);
115: delta = 1.0;
116: VecCopy(tao->solution, tao->stepdirection);
117: VecScale(tao->stepdirection, -1.0);
119: TaoLineSearchSetInitialStepLength(tao->linesearch,delta);
120: TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);
121: TaoAddLineSearchCounts(tao);
122: if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
124: /* Line search failed for last time -- give up */
125: f = f_old;
126: VecCopy(cgP->X_old, tao->solution);
127: VecCopy(cgP->G_old, tao->gradient);
128: step = 0.0;
129: tao->reason = TAO_DIVERGED_LS_FAILURE;
130: }
131: }
132: }
134: /* Check for bad value */
135: VecNorm(tao->gradient,NORM_2,&gnorm);
136: if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User-provided compute function generated Inf or NaN");
138: /* Check for termination */
139: gnorm2 =gnorm * gnorm;
140: tao->niter++;
141: TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);
142: TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);
143: (*tao->ops->convergencetest)(tao,tao->cnvP);
144: if (tao->reason != TAO_CONTINUE_ITERATING) {
145: break;
146: }
148: /* Check for restart condition */
149: VecDot(tao->gradient, cgP->G_old, &ginner);
150: if (PetscAbsScalar(ginner) >= cgP->eta * gnorm2) {
151: /* Gradients far from orthognal; use steepest descent direction */
152: beta = 0.0;
153: } else {
154: /* Gradients close to orthogonal; use conjugate gradient formula */
155: switch (cgP->cg_type) {
156: case CG_FletcherReeves:
157: beta = gnorm2 / gnorm2_old;
158: break;
160: case CG_PolakRibiere:
161: beta = (gnorm2 - ginner) / gnorm2_old;
162: break;
164: case CG_PolakRibierePlus:
165: beta = PetscMax((gnorm2-ginner)/gnorm2_old, 0.0);
166: break;
168: case CG_HestenesStiefel:
169: VecDot(tao->gradient, tao->stepdirection, &gd);
170: VecDot(cgP->G_old, tao->stepdirection, &gd_old);
171: beta = (gnorm2 - ginner) / (gd - gd_old);
172: break;
174: case CG_DaiYuan:
175: VecDot(tao->gradient, tao->stepdirection, &gd);
176: VecDot(cgP->G_old, tao->stepdirection, &gd_old);
177: beta = gnorm2 / (gd - gd_old);
178: break;
180: default:
181: beta = 0.0;
182: break;
183: }
184: }
186: /* Compute the direction d=-g + beta*d */
187: VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);
189: /* update initial steplength choice */
190: delta = 1.0;
191: delta = PetscMax(delta, cgP->delta_min);
192: delta = PetscMin(delta, cgP->delta_max);
193: }
194: return(0);
195: }
197: static PetscErrorCode TaoSetUp_CG(Tao tao)
198: {
199: TAO_CG *cgP = (TAO_CG*)tao->data;
203: if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient);}
204: if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection); }
205: if (!cgP->X_old) {VecDuplicate(tao->solution,&cgP->X_old);}
206: if (!cgP->G_old) {VecDuplicate(tao->gradient,&cgP->G_old); }
207: return(0);
208: }
210: static PetscErrorCode TaoDestroy_CG(Tao tao)
211: {
212: TAO_CG *cgP = (TAO_CG*) tao->data;
216: if (tao->setupcalled) {
217: VecDestroy(&cgP->X_old);
218: VecDestroy(&cgP->G_old);
219: }
220: TaoLineSearchDestroy(&tao->linesearch);
221: PetscFree(tao->data);
222: return(0);
223: }
225: static PetscErrorCode TaoSetFromOptions_CG(PetscOptionItems *PetscOptionsObject,Tao tao)
226: {
227: TAO_CG *cgP = (TAO_CG*)tao->data;
231: TaoLineSearchSetFromOptions(tao->linesearch);
232: PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");
233: PetscOptionsReal("-tao_cg_eta","restart tolerance", "", cgP->eta,&cgP->eta,NULL);
234: PetscOptionsEList("-tao_cg_type","cg formula", "", CG_Table, CG_Types, CG_Table[cgP->cg_type], &cgP->cg_type,NULL);
235: PetscOptionsReal("-tao_cg_delta_min","minimum delta value", "", cgP->delta_min,&cgP->delta_min,NULL);
236: PetscOptionsReal("-tao_cg_delta_max","maximum delta value", "", cgP->delta_max,&cgP->delta_max,NULL);
237: PetscOptionsTail();
238: return(0);
239: }
241: static PetscErrorCode TaoView_CG(Tao tao, PetscViewer viewer)
242: {
243: PetscBool isascii;
244: TAO_CG *cgP = (TAO_CG*)tao->data;
248: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
249: if (isascii) {
250: PetscViewerASCIIPushTab(viewer);
251: PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cgP->cg_type]);
252: PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", cgP->ngradsteps);
253: ierr= PetscViewerASCIIPrintf(viewer, "Reset steps: %D\n", cgP->nresetsteps);
254: PetscViewerASCIIPopTab(viewer);
255: }
256: return(0);
257: }
259: /*MC
260: TAOCG - Nonlinear conjugate gradient method is an extension of the
261: nonlinear conjugate gradient solver for nonlinear optimization.
263: Options Database Keys:
264: + -tao_cg_eta <r> - restart tolerance
265: . -tao_cg_type <taocg_type> - cg formula
266: . -tao_cg_delta_min <r> - minimum delta value
267: - -tao_cg_delta_max <r> - maximum delta value
269: Notes:
270: CG formulas are:
271: "fr" - Fletcher-Reeves
272: "pr" - Polak-Ribiere
273: "prp" - Polak-Ribiere-Plus
274: "hs" - Hestenes-Steifel
275: "dy" - Dai-Yuan
276: Level: beginner
277: M*/
280: PETSC_EXTERN PetscErrorCode TaoCreate_CG(Tao tao)
281: {
282: TAO_CG *cgP;
283: const char *morethuente_type = TAOLINESEARCHMT;
287: tao->ops->setup = TaoSetUp_CG;
288: tao->ops->solve = TaoSolve_CG;
289: tao->ops->view = TaoView_CG;
290: tao->ops->setfromoptions = TaoSetFromOptions_CG;
291: tao->ops->destroy = TaoDestroy_CG;
293: /* Override default settings (unless already changed) */
294: if (!tao->max_it_changed) tao->max_it = 2000;
295: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
297: /* Note: nondefault values should be used for nonlinear conjugate gradient */
298: /* method. In particular, gtol should be less that 0.5; the value used in */
299: /* Nocedal and Wright is 0.10. We use the default values for the */
300: /* linesearch because it seems to work better. */
301: TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
302: PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
303: TaoLineSearchSetType(tao->linesearch, morethuente_type);
304: TaoLineSearchUseTaoRoutines(tao->linesearch, tao);
305: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
307: PetscNewLog(tao,&cgP);
308: tao->data = (void*)cgP;
309: cgP->eta = 0.1;
310: cgP->delta_min = 1e-7;
311: cgP->delta_max = 100;
312: cgP->cg_type = CG_PolakRibierePlus;
313: return(0);
314: }