Actual source code: armijo.c

  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

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

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

 24: static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
 25: {
 26:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 27:   PetscErrorCode       ierr;

 30:   PetscFree(armP->memory);
 31:   armP->memorySetup = PETSC_FALSE;
 32:   return(0);
 33: }

 35: static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
 36: {
 37:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 38:   PetscErrorCode       ierr;

 41:   PetscOptionsHead(PetscOptionsObject,"Armijo linesearch options");
 42:   PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL);
 43:   PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL);
 44:   PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL);
 45:   PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL);
 46:   PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL);
 47:   PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL);
 48:   PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL);
 49:   PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending,NULL);
 50:   PetscOptionsTail();
 51:   return(0);
 52: }

 54: static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
 55: {
 56:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 57:   PetscBool            isascii;
 58:   PetscErrorCode       ierr;

 61:   PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
 62:   if (isascii) {
 63:     ierr=PetscViewerASCIIPrintf(pv,"  Armijo linesearch",armP->alpha);
 64:     if (armP->nondescending) {
 65:       PetscViewerASCIIPrintf(pv, " (nondescending)");
 66:     }
 67:     if (ls->bounded) {
 68:       PetscViewerASCIIPrintf(pv," (projected)");
 69:     }
 70:     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
 71:     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
 72:     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
 73:   }
 74:   return(0);
 75: }

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

 80:    Input Parameters:
 81: +  tao - Tao context
 82: .  X - current iterate (on output X contains new iterate, X + step*S)
 83: .  S - search direction
 84: .  f - merit function evaluated at X
 85: .  G - gradient of merit function evaluated at X
 86: .  W - work vector
 87: -  step - initial estimate of step length

 89:    Output parameters:
 90: +  f - merit function evaluated at new iterate, X + step*S
 91: .  G - gradient of merit function evaluated at new iterate, X + step*S
 92: .  X - new iterate
 93: -  step - final step length

 95: @ */
 96: static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
 97: {
 98:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 99:   PetscErrorCode       ierr;
100:   PetscInt             i,its=0;
101:   PetscReal            fact, ref, gdx;
102:   PetscInt             idx;
103:   PetscBool            g_computed=PETSC_FALSE; /* to prevent extra gradient computation */

106:   TaoLineSearchMonitor(ls, 0, *f, 0.0);

108:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
109:   if (!armP->work) {
110:     VecDuplicate(x,&armP->work);
111:     armP->x = x;
112:     PetscObjectReference((PetscObject)armP->x);
113:   } else if (x != armP->x) {
114:     /* If x has changed, then recreate work */
115:     VecDestroy(&armP->work);
116:     VecDuplicate(x,&armP->work);
117:     PetscObjectDereference((PetscObject)armP->x);
118:     armP->x = x;
119:     PetscObjectReference((PetscObject)armP->x);
120:   }

122:   /* Check linesearch parameters */
123:   if (armP->alpha < 1) {
124:     PetscInfo1(ls,"Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
125:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
126:   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
127:     PetscInfo1(ls,"Armijo line search error: beta (%g) invalid\n", (double)armP->beta);
128:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
129:   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
130:     PetscInfo1(ls,"Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
131:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
132:   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
133:     PetscInfo1(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
134:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
135:   } else if (armP->memorySize < 1) {
136:     PetscInfo1(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize);
137:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
138:   } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
139:     PetscInfo(ls,"Armijo line search error: reference_policy invalid\n");
140:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
141:   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
142:     PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n");
143:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
144:   } else if (PetscIsInfOrNanReal(*f)) {
145:     PetscInfo(ls,"Armijo line search error: initial function inf or nan\n");
146:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
147:   }

149:   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
150:     return(0);
151:   }

153:   /* Check to see of the memory has been allocated.  If not, allocate
154:      the historical array and populate it with the initial function
155:      values. */
156:   if (!armP->memory) {
157:     PetscMalloc1(armP->memorySize, &armP->memory);
158:   }

160:   if (!armP->memorySetup) {
161:     for (i = 0; i < armP->memorySize; i++) {
162:       armP->memory[i] = armP->alpha*(*f);
163:     }

165:     armP->current = 0;
166:     armP->lastReference = armP->memory[0];
167:     armP->memorySetup=PETSC_TRUE;
168:   }

170:   /* Calculate reference value (MAX) */
171:   ref = armP->memory[0];
172:   idx = 0;

174:   for (i = 1; i < armP->memorySize; i++) {
175:     if (armP->memory[i] > ref) {
176:       ref = armP->memory[i];
177:       idx = i;
178:     }
179:   }

181:   if (armP->referencePolicy == REFERENCE_AVE) {
182:     ref = 0;
183:     for (i = 0; i < armP->memorySize; i++) {
184:       ref += armP->memory[i];
185:     }
186:     ref = ref / armP->memorySize;
187:     ref = PetscMax(ref, armP->memory[armP->current]);
188:   } else if (armP->referencePolicy == REFERENCE_MEAN) {
189:     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
190:   }
191:   VecDot(g,s,&gdx);

193:   if (PetscIsInfOrNanReal(gdx)) {
194:     PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
195:     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
196:     return(0);
197:   }
198:   if (gdx >= 0.0) {
199:     PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
200:     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
201:     return(0);
202:   }

204:   if (armP->nondescending) {
205:     fact = armP->sigma;
206:   } else {
207:     fact = armP->sigma * gdx;
208:   }
209:   ls->step = ls->initstep;
210:   while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) {
211:     /* Calculate iterate */
212:     ++its;
213:     VecCopy(x,armP->work);
214:     VecAXPY(armP->work,ls->step,s);
215:     if (ls->bounded) {
216:       VecMedian(ls->lower,armP->work,ls->upper,armP->work);
217:     }

219:     /* Calculate function at new iterate */
220:     if (ls->hasobjective) {
221:       TaoLineSearchComputeObjective(ls,armP->work,f);
222:       g_computed=PETSC_FALSE;
223:     } else if (ls->usegts) {
224:       TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx);
225:       g_computed=PETSC_FALSE;
226:     } else {
227:       TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
228:       g_computed=PETSC_TRUE;
229:     }
230:     if (ls->step == ls->initstep) {
231:       ls->f_fullstep = *f;
232:     }

234:     TaoLineSearchMonitor(ls, its, *f, ls->step);

236:     if (PetscIsInfOrNanReal(*f)) {
237:       ls->step *= armP->beta_inf;
238:     } else {
239:       /* Check descent condition */
240:       if (armP->nondescending && *f <= ref - ls->step*fact*ref)
241:         break;
242:       if (!armP->nondescending && *f <= ref + ls->step*fact) {
243:         break;
244:       }

246:       ls->step *= armP->beta;
247:     }
248:   }

250:   /* Check termination */
251:   if (PetscIsInfOrNanReal(*f)) {
252:     PetscInfo(ls, "Function is inf or nan.\n");
253:     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
254:   } else if (ls->step < ls->stepmin) {
255:     PetscInfo(ls, "Step length is below tolerance.\n");
256:     ls->reason = TAOLINESEARCH_HALTED_RTOL;
257:   } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
258:     PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs);
259:     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
260:   }
261:   if (ls->reason) {
262:     return(0);
263:   }

265:   /* Successful termination, update memory */
266:   ls->reason = TAOLINESEARCH_SUCCESS;
267:   armP->lastReference = ref;
268:   if (armP->replacementPolicy == REPLACE_FIFO) {
269:     armP->memory[armP->current++] = *f;
270:     if (armP->current >= armP->memorySize) {
271:       armP->current = 0;
272:     }
273:   } else {
274:     armP->current = idx;
275:     armP->memory[idx] = *f;
276:   }

278:   /* Update iterate and compute gradient */
279:   VecCopy(armP->work,x);
280:   if (!g_computed) {
281:     TaoLineSearchComputeGradient(ls, x, g);
282:   }
283:   PetscInfo2(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step);
284:   return(0);
285: }

287: /*MC
288:    TAOLINESEARCHARMIJO - Backtracking line-search that satisfies only the (nonmonotone) Armijo condition
289:    (i.e., sufficient decrease).

291:    Armijo line-search type can be selected with "-tao_ls_type armijo".

293:    Level: developer

295: .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchApply()

297: .keywords: Tao, linesearch
298: M*/
299: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
300: {
301:   TaoLineSearch_ARMIJO *armP;
302:   PetscErrorCode       ierr;

306:   PetscNewLog(ls,&armP);

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