Actual source code: owarmijo.c
petsc-3.7.3 2016-08-01
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(PetscOptionItems *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: ierr=PetscViewerASCIIPrintf(pv," OWArmijo linesearch",armP->alpha);
89: if (armP->nondescending) {
90: PetscViewerASCIIPrintf(pv, " (nondescending)");
91: }
92: ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
93: ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
94: ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
95: }
96: return(0);
97: }
101: /* @ TaoApply_OWArmijo - This routine performs a linesearch. It
102: backtracks until the (nonmonotone) OWArmijo conditions are satisfied.
104: Input Parameters:
105: + tao - TAO_SOLVER context
106: . X - current iterate (on output X contains new iterate, X + step*S)
107: . S - search direction
108: . f - merit function evaluated at X
109: . G - gradient of merit function evaluated at X
110: . W - work vector
111: - step - initial estimate of step length
113: Output parameters:
114: + f - merit function evaluated at new iterate, X + step*S
115: . G - gradient of merit function evaluated at new iterate, X + step*S
116: . X - new iterate
117: - step - final step length
119: Info is set to one of:
120: . 0 - the line search succeeds; the sufficient decrease
121: condition and the directional derivative condition hold
123: negative number if an input parameter is invalid
124: - -1 - step < 0
126: positive number > 1 if the line search otherwise terminates
127: + 1 - Step is at the lower bound, stepmin.
128: @ */
129: static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
130: {
131: TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
132: PetscErrorCode ierr;
133: PetscInt i;
134: PetscReal fact, ref, gdx;
135: PetscInt idx;
136: PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
137: Vec g_old;
138: PetscReal owlqn_minstep=0.005;
139: PetscReal partgdx;
140: MPI_Comm comm;
143: PetscObjectGetComm((PetscObject)ls,&comm);
144: fact = 0.0;
145: ls->nfeval=0;
146: ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
147: if (!armP->work) {
148: VecDuplicate(x,&armP->work);
149: armP->x = x;
150: PetscObjectReference((PetscObject)armP->x);
151: } else if (x != armP->x) {
152: VecDestroy(&armP->work);
153: VecDuplicate(x,&armP->work);
154: PetscObjectDereference((PetscObject)armP->x);
155: armP->x = x;
156: PetscObjectReference((PetscObject)armP->x);
157: }
159: /* Check linesearch parameters */
160: if (armP->alpha < 1) {
161: PetscInfo1(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
162: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
163: } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
164: PetscInfo1(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta);
165: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
166: } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
167: PetscInfo1(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
168: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
169: } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
170: PetscInfo1(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
171: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
172: } else if (armP->memorySize < 1) {
173: PetscInfo1(ls,"OWArmijo line search error: memory_size (%D) < 1\n", armP->memorySize);
174: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
175: } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
176: PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n");
177: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
178: } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
179: PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n");
180: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
181: } else if (PetscIsInfOrNanReal(*f)) {
182: PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n");
183: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
184: }
186: if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) return(0);
188: /* Check to see of the memory has been allocated. If not, allocate
189: the historical array and populate it with the initial function
190: values. */
191: if (!armP->memory) {
192: PetscMalloc1(armP->memorySize, &armP->memory );
193: }
195: if (!armP->memorySetup) {
196: for (i = 0; i < armP->memorySize; i++) {
197: armP->memory[i] = armP->alpha*(*f);
198: }
199: armP->current = 0;
200: armP->lastReference = armP->memory[0];
201: armP->memorySetup=PETSC_TRUE;
202: }
204: /* Calculate reference value (MAX) */
205: ref = armP->memory[0];
206: idx = 0;
208: for (i = 1; i < armP->memorySize; i++) {
209: if (armP->memory[i] > ref) {
210: ref = armP->memory[i];
211: idx = i;
212: }
213: }
215: if (armP->referencePolicy == REFERENCE_AVE) {
216: ref = 0;
217: for (i = 0; i < armP->memorySize; i++) {
218: ref += armP->memory[i];
219: }
220: ref = ref / armP->memorySize;
221: ref = PetscMax(ref, armP->memory[armP->current]);
222: } else if (armP->referencePolicy == REFERENCE_MEAN) {
223: ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
224: }
226: if (armP->nondescending) {
227: fact = armP->sigma;
228: }
230: VecDuplicate(g,&g_old);
231: VecCopy(g,g_old);
233: ls->step = ls->initstep;
234: while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) {
235: /* Calculate iterate */
236: VecCopy(x,armP->work);
237: VecAXPY(armP->work,ls->step,s);
239: partgdx=0.0;
240: ProjWork_OWLQN(armP->work,x,g_old,&partgdx);
241: MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm);
243: /* Check the condition of gdx */
244: if (PetscIsInfOrNanReal(gdx)) {
245: PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
246: ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
247: return(0);
248: }
249: if (gdx >= 0.0) {
250: PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
251: ls->reason = TAOLINESEARCH_FAILED_ASCENT;
252: return(0);
253: }
255: /* Calculate function at new iterate */
256: TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
257: g_computed=PETSC_TRUE;
259: if (ls->step == ls->initstep) {
260: ls->f_fullstep = *f;
261: }
263: if (PetscIsInfOrNanReal(*f)) {
264: ls->step *= armP->beta_inf;
265: } else {
266: /* Check descent condition */
267: if (armP->nondescending && *f <= ref - ls->step*fact*ref) break;
268: if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break;
269: ls->step *= armP->beta;
270: }
271: }
272: VecDestroy(&g_old);
274: /* Check termination */
275: if (PetscIsInfOrNanReal(*f)) {
276: PetscInfo(ls, "Function is inf or nan.\n");
277: ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER;
278: } else if (ls->step < owlqn_minstep) {
279: PetscInfo(ls, "Step length is below tolerance.\n");
280: ls->reason = TAOLINESEARCH_HALTED_RTOL;
281: } else if (ls->nfeval >= ls->max_funcs) {
282: PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs);
283: ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
284: }
285: if (ls->reason) return(0);
287: /* Successful termination, update memory */
288: ls->reason = TAOLINESEARCH_SUCCESS;
289: armP->lastReference = ref;
290: if (armP->replacementPolicy == REPLACE_FIFO) {
291: armP->memory[armP->current++] = *f;
292: if (armP->current >= armP->memorySize) {
293: armP->current = 0;
294: }
295: } else {
296: armP->current = idx;
297: armP->memory[idx] = *f;
298: }
300: /* Update iterate and compute gradient */
301: VecCopy(armP->work,x);
302: if (!g_computed) {
303: TaoLineSearchComputeGradient(ls, x, g);
304: }
305: PetscInfo2(ls, "%D function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step);
306: return(0);
307: }
311: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls)
312: {
313: TaoLineSearch_OWARMIJO *armP;
314: PetscErrorCode ierr;
318: PetscNewLog(ls,&armP);
320: armP->memory = NULL;
321: armP->alpha = 1.0;
322: armP->beta = 0.25;
323: armP->beta_inf = 0.25;
324: armP->sigma = 1e-4;
325: armP->memorySize = 1;
326: armP->referencePolicy = REFERENCE_MAX;
327: armP->replacementPolicy = REPLACE_MRU;
328: armP->nondescending=PETSC_FALSE;
329: ls->data = (void*)armP;
330: ls->initstep=0.1;
331: ls->ops->setup=0;
332: ls->ops->reset=0;
333: ls->ops->apply=TaoLineSearchApply_OWArmijo;
334: ls->ops->view = TaoLineSearchView_OWArmijo;
335: ls->ops->destroy = TaoLineSearchDestroy_OWArmijo;
336: ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo;
337: return(0);
338: }