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;

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

 22: static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
 23: {
 24:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;

 26:   PetscFree(armP->memory);
 27:   armP->memorySetup = PETSC_FALSE;
 28:   return 0;
 29: }

 31: static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
 32: {
 33:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;

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

 48: static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
 49: {
 50:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 51:   PetscBool            isascii;

 53:   PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
 54:   if (isascii) {
 55:     PetscViewerASCIIPrintf(pv,"  Armijo linesearch",armP->alpha);
 56:     if (armP->nondescending) {
 57:       PetscViewerASCIIPrintf(pv, " (nondescending)");
 58:     }
 59:     if (ls->bounded) {
 60:       PetscViewerASCIIPrintf(pv," (projected)");
 61:     }
 62:     PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
 63:     PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
 64:     PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
 65:   }
 66:   return 0;
 67: }

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

 72:    Input Parameters:
 73: +  tao - Tao context
 74: .  X - current iterate (on output X contains new iterate, X + step*S)
 75: .  S - search direction
 76: .  f - merit function evaluated at X
 77: .  G - gradient of merit function evaluated at X
 78: .  W - work vector
 79: -  step - initial estimate of step length

 81:    Output parameters:
 82: +  f - merit function evaluated at new iterate, X + step*S
 83: .  G - gradient of merit function evaluated at new iterate, X + step*S
 84: .  X - new iterate
 85: -  step - final step length

 87: @ */
 88: static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
 89: {
 90:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 91:   PetscInt             i,its=0;
 92:   PetscReal            fact, ref, gdx;
 93:   PetscInt             idx;
 94:   PetscBool            g_computed=PETSC_FALSE; /* to prevent extra gradient computation */

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

 98:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
 99:   if (!armP->work) {
100:     VecDuplicate(x,&armP->work);
101:     armP->x = x;
102:     PetscObjectReference((PetscObject)armP->x);
103:   } else if (x != armP->x) {
104:     /* If x has changed, then recreate work */
105:     VecDestroy(&armP->work);
106:     VecDuplicate(x,&armP->work);
107:     PetscObjectDereference((PetscObject)armP->x);
108:     armP->x = x;
109:     PetscObjectReference((PetscObject)armP->x);
110:   }

112:   /* Check linesearch parameters */
113:   if (armP->alpha < 1) {
114:     PetscInfo(ls,"Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
115:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
116:   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
117:     PetscInfo(ls,"Armijo line search error: beta (%g) invalid\n", (double)armP->beta);
118:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
119:   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
120:     PetscInfo(ls,"Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
121:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
122:   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
123:     PetscInfo(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
124:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
125:   } else if (armP->memorySize < 1) {
126:     PetscInfo(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize);
127:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
128:   } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
129:     PetscInfo(ls,"Armijo line search error: reference_policy invalid\n");
130:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
131:   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
132:     PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n");
133:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
134:   } else if (PetscIsInfOrNanReal(*f)) {
135:     PetscInfo(ls,"Armijo line search error: initial function inf or nan\n");
136:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
137:   }

139:   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
140:     return 0;
141:   }

143:   /* Check to see of the memory has been allocated.  If not, allocate
144:      the historical array and populate it with the initial function
145:      values. */
146:   if (!armP->memory) {
147:     PetscMalloc1(armP->memorySize, &armP->memory);
148:   }

150:   if (!armP->memorySetup) {
151:     for (i = 0; i < armP->memorySize; i++) {
152:       armP->memory[i] = armP->alpha*(*f);
153:     }

155:     armP->current = 0;
156:     armP->lastReference = armP->memory[0];
157:     armP->memorySetup=PETSC_TRUE;
158:   }

160:   /* Calculate reference value (MAX) */
161:   ref = armP->memory[0];
162:   idx = 0;

164:   for (i = 1; i < armP->memorySize; i++) {
165:     if (armP->memory[i] > ref) {
166:       ref = armP->memory[i];
167:       idx = i;
168:     }
169:   }

171:   if (armP->referencePolicy == REFERENCE_AVE) {
172:     ref = 0;
173:     for (i = 0; i < armP->memorySize; i++) {
174:       ref += armP->memory[i];
175:     }
176:     ref = ref / armP->memorySize;
177:     ref = PetscMax(ref, armP->memory[armP->current]);
178:   } else if (armP->referencePolicy == REFERENCE_MEAN) {
179:     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
180:   }
181:   VecDot(g,s,&gdx);

183:   if (PetscIsInfOrNanReal(gdx)) {
184:     PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
185:     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
186:     return 0;
187:   }
188:   if (gdx >= 0.0) {
189:     PetscInfo(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
190:     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
191:     return 0;
192:   }

194:   if (armP->nondescending) {
195:     fact = armP->sigma;
196:   } else {
197:     fact = armP->sigma * gdx;
198:   }
199:   ls->step = ls->initstep;
200:   while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) {
201:     /* Calculate iterate */
202:     ++its;
203:     VecCopy(x,armP->work);
204:     VecAXPY(armP->work,ls->step,s);
205:     if (ls->bounded) {
206:       VecMedian(ls->lower,armP->work,ls->upper,armP->work);
207:     }

209:     /* Calculate function at new iterate */
210:     if (ls->hasobjective) {
211:       TaoLineSearchComputeObjective(ls,armP->work,f);
212:       g_computed=PETSC_FALSE;
213:     } else if (ls->usegts) {
214:       TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx);
215:       g_computed=PETSC_FALSE;
216:     } else {
217:       TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
218:       g_computed=PETSC_TRUE;
219:     }
220:     if (ls->step == ls->initstep) {
221:       ls->f_fullstep = *f;
222:     }

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

226:     if (PetscIsInfOrNanReal(*f)) {
227:       ls->step *= armP->beta_inf;
228:     } else {
229:       /* Check descent condition */
230:       if (armP->nondescending && *f <= ref - ls->step*fact*ref)
231:         break;
232:       if (!armP->nondescending && *f <= ref + ls->step*fact) {
233:         break;
234:       }

236:       ls->step *= armP->beta;
237:     }
238:   }

240:   /* Check termination */
241:   if (PetscIsInfOrNanReal(*f)) {
242:     PetscInfo(ls, "Function is inf or nan.\n");
243:     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
244:   } else if (ls->step < ls->stepmin) {
245:     PetscInfo(ls, "Step length is below tolerance.\n");
246:     ls->reason = TAOLINESEARCH_HALTED_RTOL;
247:   } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
248:     PetscInfo(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs);
249:     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
250:   }
251:   if (ls->reason) {
252:     return 0;
253:   }

255:   /* Successful termination, update memory */
256:   ls->reason = TAOLINESEARCH_SUCCESS;
257:   armP->lastReference = ref;
258:   if (armP->replacementPolicy == REPLACE_FIFO) {
259:     armP->memory[armP->current++] = *f;
260:     if (armP->current >= armP->memorySize) {
261:       armP->current = 0;
262:     }
263:   } else {
264:     armP->current = idx;
265:     armP->memory[idx] = *f;
266:   }

268:   /* Update iterate and compute gradient */
269:   VecCopy(armP->work,x);
270:   if (!g_computed) {
271:     TaoLineSearchComputeGradient(ls, x, g);
272:   }
273:   PetscInfo(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step);
274:   return 0;
275: }

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

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

283:    Level: developer

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

287: .keywords: Tao, linesearch
288: M*/
289: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
290: {
291:   TaoLineSearch_ARMIJO *armP;

294:   PetscNewLog(ls,&armP);

296:   armP->memory = NULL;
297:   armP->alpha = 1.0;
298:   armP->beta = 0.5;
299:   armP->beta_inf = 0.5;
300:   armP->sigma = 1e-4;
301:   armP->memorySize = 1;
302:   armP->referencePolicy = REFERENCE_MAX;
303:   armP->replacementPolicy = REPLACE_MRU;
304:   armP->nondescending=PETSC_FALSE;
305:   ls->data = (void*)armP;
306:   ls->initstep = 1.0;
307:   ls->ops->setup = NULL;
308:   ls->ops->monitor = NULL;
309:   ls->ops->apply = TaoLineSearchApply_Armijo;
310:   ls->ops->view = TaoLineSearchView_Armijo;
311:   ls->ops->destroy = TaoLineSearchDestroy_Armijo;
312:   ls->ops->reset = TaoLineSearchReset_Armijo;
313:   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
314:   return 0;
315: }