Actual source code: armijo.c

petsc-3.7.7 2017-09-25
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:   PetscFree(armP->memory);
 35:   armP->memorySetup = PETSC_FALSE;
 36:   return(0);
 37: }

 41: static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
 42: {
 43:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
 44:   PetscErrorCode       ierr;

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

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

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

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

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

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

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


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

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

158:   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
159:     return(0);
160:   }

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

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

174:     armP->current = 0;
175:     armP->lastReference = armP->memory[0];
176:     armP->memorySetup=PETSC_TRUE;
177:   }

179:   /* Calculate reference value (MAX) */
180:   ref = armP->memory[0];
181:   idx = 0;

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

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

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

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

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

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

252:       ls->step *= armP->beta;
253:     }
254:   }

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

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

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

295: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
296: {
297:   TaoLineSearch_ARMIJO *armP;
298:   PetscErrorCode       ierr;

302:   PetscNewLog(ls,&armP);

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