Actual source code: owarmijo.c

petsc-3.6.1 2015-08-06
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: }

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

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

 58: static PetscErrorCode TaoLineSearchSetFromOptions_OWArmijo(PetscOptions *PetscOptionsObject,TaoLineSearch ls)
 59: {
 60:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
 61:   PetscErrorCode         ierr;

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

 79: static PetscErrorCode TaoLineSearchView_OWArmijo(TaoLineSearch ls, PetscViewer pv)
 80: {
 81:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
 82:   PetscBool              isascii;
 83:   PetscErrorCode         ierr;

 86:   PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
 87:   if (isascii) {
 88:     PetscViewerASCIIPrintf(pv,"  maxf=%D, ftol=%g, gtol=%g\n",ls->max_funcs, (double)ls->rtol, (double)ls->ftol);
 89:     ierr=PetscViewerASCIIPrintf(pv,"  OWArmijo linesearch",armP->alpha);
 90:     if (armP->nondescending) {
 91:       PetscViewerASCIIPrintf(pv, " (nondescending)");
 92:     }
 93:     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
 94:     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
 95:     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
 96:   }
 97:   return(0);
 98: }

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

105:    Input Parameters:
106: +  tao - TAO_SOLVER context
107: .  X - current iterate (on output X contains new iterate, X + step*S)
108: .  S - search direction
109: .  f - merit function evaluated at X
110: .  G - gradient of merit function evaluated at X
111: .  W - work vector
112: -  step - initial estimate of step length

114:    Output parameters:
115: +  f - merit function evaluated at new iterate, X + step*S
116: .  G - gradient of merit function evaluated at new iterate, X + step*S
117: .  X - new iterate
118: -  step - final step length

120:    Info is set to one of:
121: .   0 - the line search succeeds; the sufficient decrease
122:    condition and the directional derivative condition hold

124:    negative number if an input parameter is invalid
125: -   -1 -  step < 0

127:    positive number > 1 if the line search otherwise terminates
128: +    1 -  Step is at the lower bound, stepmin.
129: @ */
130: static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
131: {
132:   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
133:   PetscErrorCode         ierr;
134:   PetscInt               i;
135:   PetscReal              fact, ref, gdx;
136:   PetscInt               idx;
137:   PetscBool              g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
138:   Vec                    g_old;
139:   PetscReal              owlqn_minstep=0.005;
140:   PetscReal              partgdx;
141:   MPI_Comm               comm;

144:   PetscObjectGetComm((PetscObject)ls,&comm);
145:   fact = 0.0;
146:   ls->nfeval=0;
147:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
148:   if (!armP->work) {
149:     VecDuplicate(x,&armP->work);
150:     armP->x = x;
151:     PetscObjectReference((PetscObject)armP->x);
152:   } else if (x != armP->x) {
153:     VecDestroy(&armP->work);
154:     VecDuplicate(x,&armP->work);
155:     PetscObjectDereference((PetscObject)armP->x);
156:     armP->x = x;
157:     PetscObjectReference((PetscObject)armP->x);
158:   }

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

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

189:   /* Check to see of the memory has been allocated.  If not, allocate
190:      the historical array and populate it with the initial function
191:      values. */
192:   if (!armP->memory) {
193:     PetscMalloc1(armP->memorySize, &armP->memory );
194:   }

196:   if (!armP->memorySetup) {
197:     for (i = 0; i < armP->memorySize; i++) {
198:       armP->memory[i] = armP->alpha*(*f);
199:     }
200:     armP->current = 0;
201:     armP->lastReference = armP->memory[0];
202:     armP->memorySetup=PETSC_TRUE;
203:   }

205:   /* Calculate reference value (MAX) */
206:   ref = armP->memory[0];
207:   idx = 0;

209:   for (i = 1; i < armP->memorySize; i++) {
210:     if (armP->memory[i] > ref) {
211:       ref = armP->memory[i];
212:       idx = i;
213:     }
214:   }

216:   if (armP->referencePolicy == REFERENCE_AVE) {
217:     ref = 0;
218:     for (i = 0; i < armP->memorySize; i++) {
219:       ref += armP->memory[i];
220:     }
221:     ref = ref / armP->memorySize;
222:     ref = PetscMax(ref, armP->memory[armP->current]);
223:   } else if (armP->referencePolicy == REFERENCE_MEAN) {
224:     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
225:   }

227:   if (armP->nondescending) {
228:     fact = armP->sigma;
229:   }

231:   VecDuplicate(g,&g_old);
232:   VecCopy(g,g_old);

234:   ls->step = ls->initstep;
235:   while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) {
236:     /* Calculate iterate */
237:     VecCopy(x,armP->work);
238:     VecAXPY(armP->work,ls->step,s);

240:     partgdx=0.0;
241:     ProjWork_OWLQN(armP->work,x,g_old,&partgdx);
242:     MPI_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm);

244:     /* Check the condition of gdx */
245:     if (PetscIsInfOrNanReal(gdx)) {
246:       PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
247:       ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
248:       return(0);
249:     }
250:     if (gdx >= 0.0) {
251:       PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
252:       ls->reason = TAOLINESEARCH_FAILED_ASCENT;
253:       return(0);
254:     }

256:     /* Calculate function at new iterate */
257:     TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
258:     g_computed=PETSC_TRUE;

260:     if (ls->step == ls->initstep) {
261:       ls->f_fullstep = *f;
262:     }

264:     if (PetscIsInfOrNanReal(*f)) {
265:       ls->step *= armP->beta_inf;
266:     } else {
267:       /* Check descent condition */
268:       if (armP->nondescending && *f <= ref - ls->step*fact*ref) break;
269:       if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break;
270:       ls->step *= armP->beta;
271:     }
272:   }
273:   VecDestroy(&g_old);

275:   /* Check termination */
276:   if (PetscIsInfOrNanReal(*f)) {
277:     PetscInfo(ls, "Function is inf or nan.\n");
278:     ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER;
279:   } else if (ls->step < owlqn_minstep) {
280:     PetscInfo(ls, "Step length is below tolerance.\n");
281:     ls->reason = TAOLINESEARCH_HALTED_RTOL;
282:   } else if (ls->nfeval >= ls->max_funcs) {
283:     PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs);
284:     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
285:   }
286:   if (ls->reason) return(0);

288:   /* Successful termination, update memory */
289:   ls->reason = TAOLINESEARCH_SUCCESS;
290:   armP->lastReference = ref;
291:   if (armP->replacementPolicy == REPLACE_FIFO) {
292:     armP->memory[armP->current++] = *f;
293:     if (armP->current >= armP->memorySize) {
294:       armP->current = 0;
295:     }
296:   } else {
297:     armP->current = idx;
298:     armP->memory[idx] = *f;
299:   }

301:   /* Update iterate and compute gradient */
302:   VecCopy(armP->work,x);
303:   if (!g_computed) {
304:     TaoLineSearchComputeGradient(ls, x, g);
305:   }
306:   PetscInfo2(ls, "%D function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step);
307:   return(0);
308: }

312: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls)
313: {
314:   TaoLineSearch_OWARMIJO *armP;
315:   PetscErrorCode         ierr;

319:   PetscNewLog(ls,&armP);

321:   armP->memory = NULL;
322:   armP->alpha = 1.0;
323:   armP->beta = 0.25;
324:   armP->beta_inf = 0.25;
325:   armP->sigma = 1e-4;
326:   armP->memorySize = 1;
327:   armP->referencePolicy = REFERENCE_MAX;
328:   armP->replacementPolicy = REPLACE_MRU;
329:   armP->nondescending=PETSC_FALSE;
330:   ls->data = (void*)armP;
331:   ls->initstep=0.1;
332:   ls->ops->setup=0;
333:   ls->ops->reset=0;
334:   ls->ops->apply=TaoLineSearchApply_OWArmijo;
335:   ls->ops->view = TaoLineSearchView_OWArmijo;
336:   ls->ops->destroy = TaoLineSearchDestroy_OWArmijo;
337:   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo;
338:   return(0);
339: }