Actual source code: taocg.c

petsc-dev 2014-02-02
Report Typos and Errors
  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"};

 15: static PetscErrorCode TaoSolve_CG(Tao tao)
 16: {
 17:   TAO_CG                         *cgP = (TAO_CG*)tao->data;
 18:   PetscErrorCode                 ierr;
 19:   TaoTerminationReason     reason = TAO_CONTINUE_ITERATING;
 20:   TaoLineSearchTerminationReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
 21:   PetscReal                      step=1.0,f,gnorm,gnorm2,delta,gd,ginner,beta;
 22:   PetscReal                      gd_old,gnorm2_old,f_old;
 23:   PetscInt                       iter=0;

 26:   if (tao->XL || tao->XU || tao->ops->computebounds) {
 27:     PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by cg algorithm\n");
 28:   }

 30:   /*  Check convergence criteria */
 31:   TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);
 32:   VecNorm(tao->gradient,NORM_2,&gnorm);
 33:   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");

 35:   TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);
 36:   if (reason != TAO_CONTINUE_ITERATING) return(0);

 38:   /*  Set initial direction to -gradient */
 39:   VecCopy(tao->gradient, tao->stepdirection);
 40:   VecScale(tao->stepdirection, -1.0);
 41:   gnorm2 = gnorm*gnorm;

 43:   /*  Set initial scaling for the function */
 44:   if (f != 0.0) {
 45:     delta = 2.0*PetscAbsScalar(f) / gnorm2;
 46:     delta = PetscMax(delta,cgP->delta_min);
 47:     delta = PetscMin(delta,cgP->delta_max);
 48:   } else {
 49:     delta = 2.0 / gnorm2;
 50:     delta = PetscMax(delta,cgP->delta_min);
 51:     delta = PetscMin(delta,cgP->delta_max);
 52:   }
 53:   /*  Set counter for gradient and reset steps */
 54:   cgP->ngradsteps = 0;
 55:   cgP->nresetsteps = 0;

 57:   while (1) {
 58:     /*  Save the current gradient information */
 59:     f_old = f;
 60:     gnorm2_old = gnorm2;
 61:     VecCopy(tao->solution, cgP->X_old);
 62:     VecCopy(tao->gradient, cgP->G_old);
 63:     VecDot(tao->gradient, tao->stepdirection, &gd);
 64:     if ((gd >= 0) || PetscIsInfOrNanReal(gd)) {
 65:       ++cgP->ngradsteps;
 66:       if (f != 0.0) {
 67:         delta = 2.0*PetscAbsScalar(f) / gnorm2;
 68:         delta = PetscMax(delta,cgP->delta_min);
 69:         delta = PetscMin(delta,cgP->delta_max);
 70:       } else {
 71:         delta = 2.0 / gnorm2;
 72:         delta = PetscMax(delta,cgP->delta_min);
 73:         delta = PetscMin(delta,cgP->delta_max);
 74:       }

 76:       VecCopy(tao->gradient, tao->stepdirection);
 77:       VecScale(tao->stepdirection, -1.0);
 78:     }

 80:     /*  Search direction for improving point */
 81:     TaoLineSearchSetInitialStepLength(tao->linesearch,delta);
 82:     TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);
 83:     TaoAddLineSearchCounts(tao);
 84:     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
 85:       /*  Linesearch failed */
 86:       /*  Reset factors and use scaled gradient step */
 87:       ++cgP->nresetsteps;
 88:       f = f_old;
 89:       gnorm2 = gnorm2_old;
 90:       VecCopy(cgP->X_old, tao->solution);
 91:       VecCopy(cgP->G_old, tao->gradient);

 93:       if (f != 0.0) {
 94:         delta = 2.0*PetscAbsScalar(f) / gnorm2;
 95:         delta = PetscMax(delta,cgP->delta_min);
 96:         delta = PetscMin(delta,cgP->delta_max);
 97:       } else {
 98:         delta = 2.0 / gnorm2;
 99:         delta = PetscMax(delta,cgP->delta_min);
100:         delta = PetscMin(delta,cgP->delta_max);
101:       }

103:       VecCopy(tao->gradient, tao->stepdirection);
104:       VecScale(tao->stepdirection, -1.0);

106:       TaoLineSearchSetInitialStepLength(tao->linesearch,delta);
107:       TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);
108:       TaoAddLineSearchCounts(tao);

110:       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
111:         /*  Linesearch failed again */
112:         /*  switch to unscaled gradient */
113:         f = f_old;
114:         gnorm2 = gnorm2_old;
115:         VecCopy(cgP->X_old, tao->solution);
116:         VecCopy(cgP->G_old, tao->gradient);
117:         delta = 1.0;
118:         VecCopy(tao->solution, tao->stepdirection);
119:         VecScale(tao->stepdirection, -1.0);

121:         TaoLineSearchSetInitialStepLength(tao->linesearch,delta);
122:         TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);
123:         TaoAddLineSearchCounts(tao);
124:         if (ls_status != TAOLINESEARCH_SUCCESS &&
125:             ls_status != TAOLINESEARCH_SUCCESS_USER) {

127:           /*  Line search failed for last time -- give up */
128:           f = f_old;
129:           gnorm2 = gnorm2_old;
130:           VecCopy(cgP->X_old, tao->solution);
131:           VecCopy(cgP->G_old, tao->gradient);
132:           step = 0.0;
133:           reason = TAO_DIVERGED_LS_FAILURE;
134:           tao->reason = TAO_DIVERGED_LS_FAILURE;
135:         }
136:       }
137:     }

139:     /*  Check for bad value */
140:     VecNorm(tao->gradient,NORM_2,&gnorm);
141:     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User-provided compute function generated Inf or NaN");

143:     /*  Check for termination */
144:     gnorm2 =gnorm * gnorm;
145:     iter++;
146:     TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);
147:     if (reason != TAO_CONTINUE_ITERATING) {
148:       break;
149:     }

151:     /*  Check for restart condition */
152:     VecDot(tao->gradient, cgP->G_old, &ginner);
153:     if (PetscAbsScalar(ginner) >= cgP->eta * gnorm2) {
154:       /*  Gradients far from orthognal; use steepest descent direction */
155:       beta = 0.0;
156:     } else {
157:       /*  Gradients close to orthogonal; use conjugate gradient formula */
158:       switch (cgP->cg_type) {
159:       case CG_FletcherReeves:
160:         beta = gnorm2 / gnorm2_old;
161:         break;

163:       case CG_PolakRibiere:
164:         beta = (gnorm2 - ginner) / gnorm2_old;
165:         break;

167:       case CG_PolakRibierePlus:
168:         beta = PetscMax((gnorm2-ginner)/gnorm2_old, 0.0);
169:         break;

171:       case CG_HestenesStiefel:
172:         VecDot(tao->gradient, tao->stepdirection, &gd);
173:         VecDot(cgP->G_old, tao->stepdirection, &gd_old);
174:         beta = (gnorm2 - ginner) / (gd - gd_old);
175:         break;

177:       case CG_DaiYuan:
178:         VecDot(tao->gradient, tao->stepdirection, &gd);
179:         VecDot(cgP->G_old, tao->stepdirection, &gd_old);
180:         beta = gnorm2 / (gd - gd_old);
181:         break;

183:       default:
184:         beta = 0.0;
185:         break;
186:       }
187:     }

189:     /*  Compute the direction d=-g + beta*d */
190:     VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);

192:     /*  update initial steplength choice */
193:     delta = 1.0;
194:     delta = PetscMax(delta, cgP->delta_min);
195:     delta = PetscMin(delta, cgP->delta_max);
196:   }
197:   return(0);
198: }

202: static PetscErrorCode TaoSetUp_CG(Tao tao)
203: {
204:   TAO_CG         *cgP = (TAO_CG*)tao->data;

208:   if (!tao->gradient) {VecDuplicate(tao->solution,&tao->gradient);}
209:   if (!tao->stepdirection) {VecDuplicate(tao->solution,&tao->stepdirection); }
210:   if (!cgP->X_old) {VecDuplicate(tao->solution,&cgP->X_old);}
211:   if (!cgP->G_old) {VecDuplicate(tao->gradient,&cgP->G_old); }
212:    return(0);
213: }

217: static PetscErrorCode TaoDestroy_CG(Tao tao)
218: {
219:   TAO_CG         *cgP = (TAO_CG*) tao->data;

223:   if (tao->setupcalled) {
224:     VecDestroy(&cgP->X_old);
225:     VecDestroy(&cgP->G_old);
226:   }
227:   TaoLineSearchDestroy(&tao->linesearch);
228:   PetscFree(tao->data);
229:   return(0);
230: }

234: static PetscErrorCode TaoSetFromOptions_CG(Tao tao)
235: {
236:    TAO_CG         *cgP = (TAO_CG*)tao->data;

240:    TaoLineSearchSetFromOptions(tao->linesearch);
241:    PetscOptionsHead("Nonlinear Conjugate Gradient method for unconstrained optimization");
242:    PetscOptionsReal("-tao_cg_eta","restart tolerance", "", cgP->eta,&cgP->eta,0);
243:    PetscOptionsEList("-tao_cg_type","cg formula", "", CG_Table, CG_Types, CG_Table[cgP->cg_type], &cgP->cg_type, 0);
244:    PetscOptionsReal("-tao_cg_delta_min","minimum delta value", "", cgP->delta_min,&cgP->delta_min,0);
245:    PetscOptionsReal("-tao_cg_delta_max","maximum delta value", "", cgP->delta_max,&cgP->delta_max,0);
246:    PetscOptionsTail();
247:    return(0);
248: }

252: static PetscErrorCode TaoView_CG(Tao tao, PetscViewer viewer)
253: {
254:   PetscBool      isascii;
255:   TAO_CG         *cgP = (TAO_CG*)tao->data;

259:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
260:   if (isascii) {
261:     PetscViewerASCIIPushTab(viewer);
262:     PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cgP->cg_type]);
263:     PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", cgP->ngradsteps);
264:     ierr= PetscViewerASCIIPrintf(viewer, "Reset steps: %D\n", cgP->nresetsteps);
265:     PetscViewerASCIIPopTab(viewer);
266:   }
267:   return(0);
268: }


271: EXTERN_C_BEGIN
274: PetscErrorCode TaoCreate_CG(Tao tao)
275: {
276:   TAO_CG         *cgP;
277:   const char     *morethuente_type = TAOLINESEARCH_MT;

281:   tao->ops->setup = TaoSetUp_CG;
282:   tao->ops->solve = TaoSolve_CG;
283:   tao->ops->view = TaoView_CG;
284:   tao->ops->setfromoptions = TaoSetFromOptions_CG;
285:   tao->ops->destroy = TaoDestroy_CG;

287:   tao->max_it = 2000;
288:   tao->max_funcs = 4000;
289:   tao->fatol = 1e-4;
290:   tao->frtol = 1e-4;

292:   /*  Note: nondefault values should be used for nonlinear conjugate gradient  */
293:   /*  method.  In particular, gtol should be less that 0.5; the value used in  */
294:   /*  Nocedal and Wright is 0.10.  We use the default values for the  */
295:   /*  linesearch because it seems to work better. */
296:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
297:   TaoLineSearchSetType(tao->linesearch, morethuente_type);
298:   TaoLineSearchUseTaoRoutines(tao->linesearch, tao);

300:   PetscNewLog(tao,&cgP);
301:   tao->data = (void*)cgP;
302:   cgP->eta = 0.1;
303:   cgP->delta_min = 1e-7;
304:   cgP->delta_max = 100;
305:   cgP->cg_type = CG_PolakRibierePlus;
306:   return(0);
307: }
308: EXTERN_C_END