Actual source code: armijo.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: #include <petsc/private/taolinesearchimpl.h>
  2: #include <../src/tao/linesearch/impls/armijo/armijo.h>

  4: #define REPLACE_FIFO 1
  5: #define REPLACE_MRU  2

  7: #define REFERENCE_MAX  1
  8: #define REFERENCE_AVE  2
  9: #define REFERENCE_MEAN 3

 13: static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls)
 14: {
 15:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 16:   PetscErrorCode       ierr;

 19:   PetscFree(armP->memory);
 20:   VecDestroy(&armP->x);
 21:   VecDestroy(&armP->work);
 22:   PetscFree(ls->data);
 23:   return(0);
 24: }

 28: static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
 29: {
 30:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 31:   PetscErrorCode       ierr;

 34:   if (armP->memory != NULL) {
 35:     PetscFree(armP->memory);
 36:   }
 37:   armP->memorySetup = PETSC_FALSE;
 38:   return(0);
 39: }

 43: static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptions *PetscOptionsObject,TaoLineSearch ls)
 44: {
 45:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 46:   PetscErrorCode       ierr;

 49:   PetscOptionsHead(PetscOptionsObject,"Armijo linesearch options");
 50:   PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL);
 51:   PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL);
 52:   PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL);
 53:   PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL);
 54:   PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL);
 55:   PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL);
 56:   PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL);
 57:   PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending,NULL);
 58:   PetscOptionsTail();
 59:   return(0);
 60: }

 64: static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
 65: {
 66:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 67:   PetscBool            isascii;
 68:   PetscErrorCode       ierr;

 71:   PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
 72:   if (isascii) {
 73:     PetscViewerASCIIPrintf(pv,"  maxf=%D, ftol=%g, gtol=%g\n",ls->max_funcs, (double)ls->rtol, (double)ls->ftol);
 74:     ierr=PetscViewerASCIIPrintf(pv,"  Armijo linesearch",armP->alpha);
 75:     if (armP->nondescending) {
 76:       PetscViewerASCIIPrintf(pv, " (nondescending)");
 77:     }
 78:     if (ls->bounded) {
 79:       PetscViewerASCIIPrintf(pv," (projected)");
 80:     }
 81:     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
 82:     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
 83:     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
 84:   }
 85:   return(0);
 86: }

 90: /* @ TaoApply_Armijo - This routine performs a linesearch. It
 91:    backtracks until the (nonmonotone) Armijo conditions are satisfied.

 93:    Input Parameters:
 94: +  tao - Tao context
 95: .  X - current iterate (on output X contains new iterate, X + step*S)
 96: .  S - search direction
 97: .  f - merit function evaluated at X
 98: .  G - gradient of merit function evaluated at X
 99: .  W - work vector
100: -  step - initial estimate of step length

102:    Output parameters:
103: +  f - merit function evaluated at new iterate, X + step*S
104: .  G - gradient of merit function evaluated at new iterate, X + step*S
105: .  X - new iterate
106: -  step - final step length

108: @ */
109: static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
110: {
111:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
112:   PetscErrorCode       ierr;
113:   PetscInt             i;
114:   PetscReal            fact, ref, gdx;
115:   PetscInt             idx;
116:   PetscBool            g_computed=PETSC_FALSE; /* to prevent extra gradient computation */


120:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
121:   if (!armP->work) {
122:     VecDuplicate(x,&armP->work);
123:     armP->x = x;
124:     PetscObjectReference((PetscObject)armP->x);
125:   } else if (x != armP->x) {
126:     /* If x has changed, then recreate work */
127:     VecDestroy(&armP->work);
128:     VecDuplicate(x,&armP->work);
129:     PetscObjectDereference((PetscObject)armP->x);
130:     armP->x = x;
131:     PetscObjectReference((PetscObject)armP->x);
132:   }

134:   /* Check linesearch parameters */
135:   if (armP->alpha < 1) {
136:     PetscInfo1(ls,"Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
137:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
138:   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
139:     PetscInfo1(ls,"Armijo line search error: beta (%g) invalid\n", (double)armP->beta);
140:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
141:   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
142:     PetscInfo1(ls,"Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
143:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
144:   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
145:     PetscInfo1(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
146:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
147:   } else if (armP->memorySize < 1) {
148:     PetscInfo1(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize);
149:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
150:   } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
151:     PetscInfo(ls,"Armijo line search error: reference_policy invalid\n");
152:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
153:   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
154:     PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n");
155:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
156:   } else if (PetscIsInfOrNanReal(*f)) {
157:     PetscInfo(ls,"Armijo line search error: initial function inf or nan\n");
158:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
159:   }

161:   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
162:     return(0);
163:   }

165:   /* Check to see of the memory has been allocated.  If not, allocate
166:      the historical array and populate it with the initial function
167:      values. */
168:   if (!armP->memory) {
169:     PetscMalloc1(armP->memorySize, &armP->memory );
170:   }

172:   if (!armP->memorySetup) {
173:     for (i = 0; i < armP->memorySize; i++) {
174:       armP->memory[i] = armP->alpha*(*f);
175:     }

177:     armP->current = 0;
178:     armP->lastReference = armP->memory[0];
179:     armP->memorySetup=PETSC_TRUE;
180:   }

182:   /* Calculate reference value (MAX) */
183:   ref = armP->memory[0];
184:   idx = 0;

186:   for (i = 1; i < armP->memorySize; i++) {
187:     if (armP->memory[i] > ref) {
188:       ref = armP->memory[i];
189:       idx = i;
190:     }
191:   }

193:   if (armP->referencePolicy == REFERENCE_AVE) {
194:     ref = 0;
195:     for (i = 0; i < armP->memorySize; i++) {
196:       ref += armP->memory[i];
197:     }
198:     ref = ref / armP->memorySize;
199:     ref = PetscMax(ref, armP->memory[armP->current]);
200:   } else if (armP->referencePolicy == REFERENCE_MEAN) {
201:     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
202:   }
203:   VecDot(g,s,&gdx);

205:   if (PetscIsInfOrNanReal(gdx)) {
206:     PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
207:     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
208:     return(0);
209:   }
210:   if (gdx >= 0.0) {
211:     PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
212:     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
213:     return(0);
214:   }

216:   if (armP->nondescending) {
217:     fact = armP->sigma;
218:   } else {
219:     fact = armP->sigma * gdx;
220:   }
221:   ls->step = ls->initstep;
222:   while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) {
223:     /* Calculate iterate */
224:     VecCopy(x,armP->work);
225:     VecAXPY(armP->work,ls->step,s);
226:     if (ls->bounded) {
227:       VecMedian(ls->lower,armP->work,ls->upper,armP->work);
228:     }

230:     /* Calculate function at new iterate */
231:     if (ls->hasobjective) {
232:       TaoLineSearchComputeObjective(ls,armP->work,f);
233:       g_computed=PETSC_FALSE;
234:     } else if (ls->usegts) {
235:       TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx);
236:       g_computed=PETSC_FALSE;
237:     } else {
238:       TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
239:       g_computed=PETSC_TRUE;
240:     }
241:     if (ls->step == ls->initstep) {
242:       ls->f_fullstep = *f;
243:     }

245:     if (PetscIsInfOrNanReal(*f)) {
246:       ls->step *= armP->beta_inf;
247:     } else {
248:       /* Check descent condition */
249:       if (armP->nondescending && *f <= ref - ls->step*fact*ref)
250:         break;
251:       if (!armP->nondescending && *f <= ref + ls->step*fact) {
252:         break;
253:       }

255:       ls->step *= armP->beta;
256:     }
257:   }

259:   /* Check termination */
260:   if (PetscIsInfOrNanReal(*f)) {
261:     PetscInfo(ls, "Function is inf or nan.\n");
262:     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
263:   } else if (ls->step < ls->stepmin) {
264:     PetscInfo(ls, "Step length is below tolerance.\n");
265:     ls->reason = TAOLINESEARCH_HALTED_RTOL;
266:   } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
267:     PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs);
268:     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
269:   }
270:   if (ls->reason) {
271:     return(0);
272:   }

274:   /* Successful termination, update memory */
275:   ls->reason = TAOLINESEARCH_SUCCESS;
276:   armP->lastReference = ref;
277:   if (armP->replacementPolicy == REPLACE_FIFO) {
278:     armP->memory[armP->current++] = *f;
279:     if (armP->current >= armP->memorySize) {
280:       armP->current = 0;
281:     }
282:   } else {
283:     armP->current = idx;
284:     armP->memory[idx] = *f;
285:   }

287:   /* Update iterate and compute gradient */
288:   VecCopy(armP->work,x);
289:   if (!g_computed) {
290:     TaoLineSearchComputeGradient(ls, x, g);
291:   }
292:   PetscInfo2(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step);
293:   return(0);
294: }

298: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
299: {
300:   TaoLineSearch_ARMIJO *armP;
301:   PetscErrorCode       ierr;

305:   PetscNewLog(ls,&armP);

307:   armP->memory = NULL;
308:   armP->alpha = 1.0;
309:   armP->beta = 0.5;
310:   armP->beta_inf = 0.5;
311:   armP->sigma = 1e-4;
312:   armP->memorySize = 1;
313:   armP->referencePolicy = REFERENCE_MAX;
314:   armP->replacementPolicy = REPLACE_MRU;
315:   armP->nondescending=PETSC_FALSE;
316:   ls->data = (void*)armP;
317:   ls->initstep=1.0;
318:   ls->ops->setup=0;
319:   ls->ops->apply=TaoLineSearchApply_Armijo;
320:   ls->ops->view = TaoLineSearchView_Armijo;
321:   ls->ops->destroy = TaoLineSearchDestroy_Armijo;
322:   ls->ops->reset = TaoLineSearchReset_Armijo;
323:   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
324:   return(0);
325: }