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(TaoLineSearch ls, PetscOptionItems *PetscOptionsObject)
 32: {
 33:   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;

 35:   PetscOptionsHeadBegin(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:   PetscOptionsHeadEnd();
 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");
 56:     if (armP->nondescending) PetscViewerASCIIPrintf(pv, " (nondescending)");
 57:     if (ls->bounded) PetscViewerASCIIPrintf(pv, " (projected)");
 58:     PetscViewerASCIIPrintf(pv, ": alpha=%g beta=%g ", (double)armP->alpha, (double)armP->beta);
 59:     PetscViewerASCIIPrintf(pv, "sigma=%g ", (double)armP->sigma);
 60:     PetscViewerASCIIPrintf(pv, "memsize=%" PetscInt_FMT "\n", armP->memorySize);
 61:   }
 62:   return 0;
 63: }

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

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

 77:    Output parameters:
 78: +  f - merit function evaluated at new iterate, X + step*S
 79: .  G - gradient of merit function evaluated at new iterate, X + step*S
 80: .  X - new iterate
 81: -  step - final step length

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

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

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

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

135:   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) return 0;

137:   /* Check to see of the memory has been allocated.  If not, allocate
138:      the historical array and populate it with the initial function
139:      values. */
140:   if (!armP->memory) PetscMalloc1(armP->memorySize, &armP->memory);

142:   if (!armP->memorySetup) {
143:     for (i = 0; i < armP->memorySize; i++) armP->memory[i] = armP->alpha * (*f);

145:     armP->current       = 0;
146:     armP->lastReference = armP->memory[0];
147:     armP->memorySetup   = PETSC_TRUE;
148:   }

150:   /* Calculate reference value (MAX) */
151:   ref = armP->memory[0];
152:   idx = 0;

154:   for (i = 1; i < armP->memorySize; i++) {
155:     if (armP->memory[i] > ref) {
156:       ref = armP->memory[i];
157:       idx = i;
158:     }
159:   }

161:   if (armP->referencePolicy == REFERENCE_AVE) {
162:     ref = 0;
163:     for (i = 0; i < armP->memorySize; i++) ref += armP->memory[i];
164:     ref = ref / armP->memorySize;
165:     ref = PetscMax(ref, armP->memory[armP->current]);
166:   } else if (armP->referencePolicy == REFERENCE_MEAN) {
167:     ref = PetscMin(ref, 0.5 * (armP->lastReference + armP->memory[armP->current]));
168:   }
169:   VecDot(g, s, &gdx);

171:   if (PetscIsInfOrNanReal(gdx)) {
172:     PetscInfo(ls, "Initial Line Search step * g is Inf or Nan (%g)\n", (double)gdx);
173:     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
174:     return 0;
175:   }
176:   if (gdx >= 0.0) {
177:     PetscInfo(ls, "Initial Line Search step is not descent direction (g's=%g)\n", (double)gdx);
178:     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
179:     return 0;
180:   }

182:   if (armP->nondescending) {
183:     fact = armP->sigma;
184:   } else {
185:     fact = armP->sigma * gdx;
186:   }
187:   ls->step = ls->initstep;
188:   while (ls->step >= ls->stepmin && (ls->nfeval + ls->nfgeval) < ls->max_funcs) {
189:     /* Calculate iterate */
190:     ++its;
191:     VecWAXPY(armP->work, ls->step, s, x);
192:     if (ls->bounded) VecMedian(ls->lower, armP->work, ls->upper, armP->work);

194:     /* Calculate function at new iterate */
195:     if (ls->hasobjective) {
196:       TaoLineSearchComputeObjective(ls, armP->work, f);
197:       g_computed = PETSC_FALSE;
198:     } else if (ls->usegts) {
199:       TaoLineSearchComputeObjectiveAndGTS(ls, armP->work, f, &gdx);
200:       g_computed = PETSC_FALSE;
201:     } else {
202:       TaoLineSearchComputeObjectiveAndGradient(ls, armP->work, f, g);
203:       g_computed = PETSC_TRUE;
204:     }
205:     if (ls->step == ls->initstep) ls->f_fullstep = *f;

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

209:     if (PetscIsInfOrNanReal(*f)) {
210:       ls->step *= armP->beta_inf;
211:     } else {
212:       /* Check descent condition */
213:       if (armP->nondescending && *f <= ref - ls->step * fact * ref) break;
214:       if (!armP->nondescending && *f <= ref + ls->step * fact) break;

216:       ls->step *= armP->beta;
217:     }
218:   }

220:   /* Check termination */
221:   if (PetscIsInfOrNanReal(*f)) {
222:     PetscInfo(ls, "Function is inf or nan.\n");
223:     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
224:   } else if (ls->step < ls->stepmin) {
225:     PetscInfo(ls, "Step length is below tolerance.\n");
226:     ls->reason = TAOLINESEARCH_HALTED_RTOL;
227:   } else if ((ls->nfeval + ls->nfgeval) >= ls->max_funcs) {
228:     PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum allowed (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs);
229:     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
230:   }
231:   if (ls->reason) return 0;

233:   /* Successful termination, update memory */
234:   ls->reason          = TAOLINESEARCH_SUCCESS;
235:   armP->lastReference = ref;
236:   if (armP->replacementPolicy == REPLACE_FIFO) {
237:     armP->memory[armP->current++] = *f;
238:     if (armP->current >= armP->memorySize) armP->current = 0;
239:   } else {
240:     armP->current     = idx;
241:     armP->memory[idx] = *f;
242:   }

244:   /* Update iterate and compute gradient */
245:   VecCopy(armP->work, x);
246:   if (!g_computed) TaoLineSearchComputeGradient(ls, x, g);
247:   PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval, (double)ls->step);
248:   return 0;
249: }

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

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

257:    Level: developer

259: .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`

261: .keywords: Tao, linesearch
262: M*/
263: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
264: {
265:   TaoLineSearch_ARMIJO *armP;

268:   PetscNew(&armP);

270:   armP->memory            = NULL;
271:   armP->alpha             = 1.0;
272:   armP->beta              = 0.5;
273:   armP->beta_inf          = 0.5;
274:   armP->sigma             = 1e-4;
275:   armP->memorySize        = 1;
276:   armP->referencePolicy   = REFERENCE_MAX;
277:   armP->replacementPolicy = REPLACE_MRU;
278:   armP->nondescending     = PETSC_FALSE;
279:   ls->data                = (void *)armP;
280:   ls->initstep            = 1.0;
281:   ls->ops->setup          = NULL;
282:   ls->ops->monitor        = NULL;
283:   ls->ops->apply          = TaoLineSearchApply_Armijo;
284:   ls->ops->view           = TaoLineSearchView_Armijo;
285:   ls->ops->destroy        = TaoLineSearchDestroy_Armijo;
286:   ls->ops->reset          = TaoLineSearchReset_Armijo;
287:   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
288:   return 0;
289: }