Actual source code: owarmijo.c

petsc-3.9.4 2018-09-11
Report Typos and Errors

  2:  #include <petsc/private/taolinesearchimpl.h>
  3:  #include <../src/tao/linesearch/impls/owarmijo/owarmijo.h>

  5: #define REPLACE_FIFO 1
  6: #define REPLACE_MRU  2

  8: #define REFERENCE_MAX  1
  9: #define REFERENCE_AVE  2
 10: #define REFERENCE_MEAN 3

 12: static PetscErrorCode ProjWork_OWLQN(Vec w,Vec x,Vec gv,PetscReal *gdx)
 13: {
 14:   const PetscReal *xptr,*gptr;
 15:   PetscReal       *wptr;
 16:   PetscErrorCode  ierr;
 17:   PetscInt        low,high,low1,high1,low2,high2,i;

 20:   ierr=VecGetOwnershipRange(w,&low,&high);
 21:   ierr=VecGetOwnershipRange(x,&low1,&high1);
 22:   ierr=VecGetOwnershipRange(gv,&low2,&high2);

 24:   *gdx=0.0;
 25:   VecGetArray(w,&wptr);
 26:   VecGetArrayRead(x,&xptr);
 27:   VecGetArrayRead(gv,&gptr);

 29:   for (i=0;i<high-low;i++) {
 30:     if (xptr[i]*wptr[i]<0.0) wptr[i]=0.0;
 31:     *gdx = *gdx + gptr[i]*(wptr[i]-xptr[i]);
 32:   }
 33:   VecRestoreArray(w,&wptr);
 34:   VecRestoreArrayRead(x,&xptr);
 35:   VecRestoreArrayRead(gv,&gptr);
 36:   return(0);
 37: }

 39: static PetscErrorCode TaoLineSearchDestroy_OWArmijo(TaoLineSearch ls)
 40: {
 41:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
 42:   PetscErrorCode         ierr;

 45:   PetscFree(armP->memory);
 46:   if (armP->x) {
 47:     PetscObjectDereference((PetscObject)armP->x);
 48:   }
 49:   VecDestroy(&armP->work);
 50:   PetscFree(ls->data);
 51:   return(0);
 52: }

 54: static PetscErrorCode TaoLineSearchSetFromOptions_OWArmijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
 55: {
 56:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
 57:   PetscErrorCode         ierr;

 60:   PetscOptionsHead(PetscOptionsObject,"OWArmijo linesearch options");
 61:   PetscOptionsReal("-tao_ls_OWArmijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL);
 62:   PetscOptionsReal("-tao_ls_OWArmijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL);
 63:   PetscOptionsReal("-tao_ls_OWArmijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL);
 64:   PetscOptionsReal("-tao_ls_OWArmijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL);
 65:   PetscOptionsInt("-tao_ls_OWArmijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL);
 66:   PetscOptionsInt("-tao_ls_OWArmijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL);
 67:   PetscOptionsInt("-tao_ls_OWArmijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL);
 68:   PetscOptionsBool("-tao_ls_OWArmijo_nondescending","Use nondescending OWArmijo algorithm","",armP->nondescending,&armP->nondescending,NULL);
 69:   PetscOptionsTail();
 70:   return(0);
 71: }

 73: static PetscErrorCode TaoLineSearchView_OWArmijo(TaoLineSearch ls, PetscViewer pv)
 74: {
 75:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
 76:   PetscBool              isascii;
 77:   PetscErrorCode         ierr;

 80:   PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
 81:   if (isascii) {
 82:     ierr=PetscViewerASCIIPrintf(pv,"  OWArmijo linesearch",armP->alpha);
 83:     if (armP->nondescending) {
 84:       PetscViewerASCIIPrintf(pv, " (nondescending)");
 85:     }
 86:     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
 87:     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
 88:     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
 89:   }
 90:   return(0);
 91: }

 93: /* @ TaoApply_OWArmijo - This routine performs a linesearch. It
 94:    backtracks until the (nonmonotone) OWArmijo conditions are satisfied.

 96:    Input Parameters:
 97: +  tao - TAO_SOLVER context
 98: .  X - current iterate (on output X contains new iterate, X + step*S)
 99: .  S - search direction
100: .  f - merit function evaluated at X
101: .  G - gradient of merit function evaluated at X
102: .  W - work vector
103: -  step - initial estimate of step length

105:    Output parameters:
106: +  f - merit function evaluated at new iterate, X + step*S
107: .  G - gradient of merit function evaluated at new iterate, X + step*S
108: .  X - new iterate
109: -  step - final step length

111:    Info is set to one of:
112: .   0 - the line search succeeds; the sufficient decrease
113:    condition and the directional derivative condition hold

115:    negative number if an input parameter is invalid
116: -   -1 -  step < 0

118:    positive number > 1 if the line search otherwise terminates
119: +    1 -  Step is at the lower bound, stepmin.
120: @ */
121: static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
122: {
123:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
124:   PetscErrorCode         ierr;
125:   PetscInt               i;
126:   PetscReal              fact, ref, gdx;
127:   PetscInt               idx;
128:   PetscBool              g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
129:   Vec                    g_old;
130:   PetscReal              owlqn_minstep=0.005;
131:   PetscReal              partgdx;
132:   MPI_Comm               comm;

135:   PetscObjectGetComm((PetscObject)ls,&comm);
136:   fact = 0.0;
137:   ls->nfeval=0;
138:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
139:   if (!armP->work) {
140:     VecDuplicate(x,&armP->work);
141:     armP->x = x;
142:     PetscObjectReference((PetscObject)armP->x);
143:   } else if (x != armP->x) {
144:     VecDestroy(&armP->work);
145:     VecDuplicate(x,&armP->work);
146:     PetscObjectDereference((PetscObject)armP->x);
147:     armP->x = x;
148:     PetscObjectReference((PetscObject)armP->x);
149:   }

151:   /* Check linesearch parameters */
152:   if (armP->alpha < 1) {
153:     PetscInfo1(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
154:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
155:   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
156:     PetscInfo1(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta);
157:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
158:   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
159:     PetscInfo1(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
160:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
161:   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
162:     PetscInfo1(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
163:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
164:   } else if (armP->memorySize < 1) {
165:     PetscInfo1(ls,"OWArmijo line search error: memory_size (%D) < 1\n", armP->memorySize);
166:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
167:   }  else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
168:     PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n");
169:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
170:   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
171:     PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n");
172:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
173:   } else if (PetscIsInfOrNanReal(*f)) {
174:     PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n");
175:     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
176:   }

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

180:   /* Check to see of the memory has been allocated.  If not, allocate
181:      the historical array and populate it with the initial function
182:      values. */
183:   if (!armP->memory) {
184:     PetscMalloc1(armP->memorySize, &armP->memory );
185:   }

187:   if (!armP->memorySetup) {
188:     for (i = 0; i < armP->memorySize; i++) {
189:       armP->memory[i] = armP->alpha*(*f);
190:     }
191:     armP->current = 0;
192:     armP->lastReference = armP->memory[0];
193:     armP->memorySetup=PETSC_TRUE;
194:   }

196:   /* Calculate reference value (MAX) */
197:   ref = armP->memory[0];
198:   idx = 0;

200:   for (i = 1; i < armP->memorySize; i++) {
201:     if (armP->memory[i] > ref) {
202:       ref = armP->memory[i];
203:       idx = i;
204:     }
205:   }

207:   if (armP->referencePolicy == REFERENCE_AVE) {
208:     ref = 0;
209:     for (i = 0; i < armP->memorySize; i++) {
210:       ref += armP->memory[i];
211:     }
212:     ref = ref / armP->memorySize;
213:     ref = PetscMax(ref, armP->memory[armP->current]);
214:   } else if (armP->referencePolicy == REFERENCE_MEAN) {
215:     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
216:   }

218:   if (armP->nondescending) {
219:     fact = armP->sigma;
220:   }

222:   VecDuplicate(g,&g_old);
223:   VecCopy(g,g_old);

225:   ls->step = ls->initstep;
226:   while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) {
227:     /* Calculate iterate */
228:     VecCopy(x,armP->work);
229:     VecAXPY(armP->work,ls->step,s);

231:     partgdx=0.0;
232:     ProjWork_OWLQN(armP->work,x,g_old,&partgdx);
233:     MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm);

235:     /* Check the condition of gdx */
236:     if (PetscIsInfOrNanReal(gdx)) {
237:       PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
238:       ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
239:       return(0);
240:     }
241:     if (gdx >= 0.0) {
242:       PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
243:       ls->reason = TAOLINESEARCH_FAILED_ASCENT;
244:       return(0);
245:     }

247:     /* Calculate function at new iterate */
248:     TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
249:     g_computed=PETSC_TRUE;

251:     if (ls->step == ls->initstep) {
252:       ls->f_fullstep = *f;
253:     }

255:     if (PetscIsInfOrNanReal(*f)) {
256:       ls->step *= armP->beta_inf;
257:     } else {
258:       /* Check descent condition */
259:       if (armP->nondescending && *f <= ref - ls->step*fact*ref) break;
260:       if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break;
261:       ls->step *= armP->beta;
262:     }
263:   }
264:   VecDestroy(&g_old);

266:   /* Check termination */
267:   if (PetscIsInfOrNanReal(*f)) {
268:     PetscInfo(ls, "Function is inf or nan.\n");
269:     ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER;
270:   } else if (ls->step < owlqn_minstep) {
271:     PetscInfo(ls, "Step length is below tolerance.\n");
272:     ls->reason = TAOLINESEARCH_HALTED_RTOL;
273:   } else if (ls->nfeval >= ls->max_funcs) {
274:     PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs);
275:     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
276:   }
277:   if (ls->reason) return(0);

279:   /* Successful termination, update memory */
280:   ls->reason = TAOLINESEARCH_SUCCESS;
281:   armP->lastReference = ref;
282:   if (armP->replacementPolicy == REPLACE_FIFO) {
283:     armP->memory[armP->current++] = *f;
284:     if (armP->current >= armP->memorySize) {
285:       armP->current = 0;
286:     }
287:   } else {
288:     armP->current = idx;
289:     armP->memory[idx] = *f;
290:   }

292:   /* Update iterate and compute gradient */
293:   VecCopy(armP->work,x);
294:   if (!g_computed) {
295:     TaoLineSearchComputeGradient(ls, x, g);
296:   }
297:   PetscInfo2(ls, "%D function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step);
298:   return(0);
299: }

301: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls)
302: {
303:   TaoLineSearch_OWARMIJO *armP;
304:   PetscErrorCode         ierr;

308:   PetscNewLog(ls,&armP);

310:   armP->memory = NULL;
311:   armP->alpha = 1.0;
312:   armP->beta = 0.25;
313:   armP->beta_inf = 0.25;
314:   armP->sigma = 1e-4;
315:   armP->memorySize = 1;
316:   armP->referencePolicy = REFERENCE_MAX;
317:   armP->replacementPolicy = REPLACE_MRU;
318:   armP->nondescending=PETSC_FALSE;
319:   ls->data = (void*)armP;
320:   ls->initstep=0.1;
321:   ls->ops->setup=0;
322:   ls->ops->reset=0;
323:   ls->ops->apply=TaoLineSearchApply_OWArmijo;
324:   ls->ops->view = TaoLineSearchView_OWArmijo;
325:   ls->ops->destroy = TaoLineSearchDestroy_OWArmijo;
326:   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo;
327:   return(0);
328: }