Actual source code: armijo.c

petsc-3.9.4 2018-09-11
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

 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;
101:   PetscReal            fact, ref, gdx;
102:   PetscInt             idx;
103:   PetscBool            g_computed=PETSC_FALSE; /* to prevent extra gradient computation */


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

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

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

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

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

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

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

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

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

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

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

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

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

242:       ls->step *= armP->beta;
243:     }
244:   }

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

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

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

283: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
284: {
285:   TaoLineSearch_ARMIJO *armP;
286:   PetscErrorCode       ierr;

290:   PetscNewLog(ls,&armP);

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