Actual source code: owarmijo.c
petsc-3.10.5 2019-03-28
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, its=0;
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: }
150:
151: TaoLineSearchMonitor(ls, 0, *f, 0.0);
153: /* Check linesearch parameters */
154: if (armP->alpha < 1) {
155: PetscInfo1(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
156: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
157: } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
158: PetscInfo1(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta);
159: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
160: } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
161: PetscInfo1(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
162: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
163: } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
164: PetscInfo1(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
165: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
166: } else if (armP->memorySize < 1) {
167: PetscInfo1(ls,"OWArmijo line search error: memory_size (%D) < 1\n", armP->memorySize);
168: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
169: } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
170: PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n");
171: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
172: } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
173: PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n");
174: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
175: } else if (PetscIsInfOrNanReal(*f)) {
176: PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n");
177: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
178: }
180: if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) return(0);
182: /* Check to see of the memory has been allocated. If not, allocate
183: the historical array and populate it with the initial function
184: values. */
185: if (!armP->memory) {
186: PetscMalloc1(armP->memorySize, &armP->memory );
187: }
189: if (!armP->memorySetup) {
190: for (i = 0; i < armP->memorySize; i++) {
191: armP->memory[i] = armP->alpha*(*f);
192: }
193: armP->current = 0;
194: armP->lastReference = armP->memory[0];
195: armP->memorySetup=PETSC_TRUE;
196: }
198: /* Calculate reference value (MAX) */
199: ref = armP->memory[0];
200: idx = 0;
202: for (i = 1; i < armP->memorySize; i++) {
203: if (armP->memory[i] > ref) {
204: ref = armP->memory[i];
205: idx = i;
206: }
207: }
209: if (armP->referencePolicy == REFERENCE_AVE) {
210: ref = 0;
211: for (i = 0; i < armP->memorySize; i++) {
212: ref += armP->memory[i];
213: }
214: ref = ref / armP->memorySize;
215: ref = PetscMax(ref, armP->memory[armP->current]);
216: } else if (armP->referencePolicy == REFERENCE_MEAN) {
217: ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
218: }
220: if (armP->nondescending) {
221: fact = armP->sigma;
222: }
224: VecDuplicate(g,&g_old);
225: VecCopy(g,g_old);
227: ls->step = ls->initstep;
228: while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) {
229: /* Calculate iterate */
230: ++its;
231: VecCopy(x,armP->work);
232: VecAXPY(armP->work,ls->step,s);
234: partgdx=0.0;
235: ProjWork_OWLQN(armP->work,x,g_old,&partgdx);
236: MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm);
238: /* Check the condition of gdx */
239: if (PetscIsInfOrNanReal(gdx)) {
240: PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
241: ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
242: return(0);
243: }
244: if (gdx >= 0.0) {
245: PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
246: ls->reason = TAOLINESEARCH_FAILED_ASCENT;
247: return(0);
248: }
250: /* Calculate function at new iterate */
251: TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
252: g_computed=PETSC_TRUE;
253:
254: TaoLineSearchMonitor(ls, its, *f, ls->step);
256: if (ls->step == ls->initstep) {
257: ls->f_fullstep = *f;
258: }
260: if (PetscIsInfOrNanReal(*f)) {
261: ls->step *= armP->beta_inf;
262: } else {
263: /* Check descent condition */
264: if (armP->nondescending && *f <= ref - ls->step*fact*ref) break;
265: if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break;
266: ls->step *= armP->beta;
267: }
268: }
269: VecDestroy(&g_old);
271: /* Check termination */
272: if (PetscIsInfOrNanReal(*f)) {
273: PetscInfo(ls, "Function is inf or nan.\n");
274: ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER;
275: } else if (ls->step < owlqn_minstep) {
276: PetscInfo(ls, "Step length is below tolerance.\n");
277: ls->reason = TAOLINESEARCH_HALTED_RTOL;
278: } else if (ls->nfeval >= ls->max_funcs) {
279: PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs);
280: ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
281: }
282: if (ls->reason) return(0);
284: /* Successful termination, update memory */
285: ls->reason = TAOLINESEARCH_SUCCESS;
286: armP->lastReference = ref;
287: if (armP->replacementPolicy == REPLACE_FIFO) {
288: armP->memory[armP->current++] = *f;
289: if (armP->current >= armP->memorySize) {
290: armP->current = 0;
291: }
292: } else {
293: armP->current = idx;
294: armP->memory[idx] = *f;
295: }
297: /* Update iterate and compute gradient */
298: VecCopy(armP->work,x);
299: if (!g_computed) {
300: TaoLineSearchComputeGradient(ls, x, g);
301: }
302: PetscInfo2(ls, "%D function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step);
303: return(0);
304: }
306: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls)
307: {
308: TaoLineSearch_OWARMIJO *armP;
309: PetscErrorCode ierr;
313: PetscNewLog(ls,&armP);
315: armP->memory = NULL;
316: armP->alpha = 1.0;
317: armP->beta = 0.25;
318: armP->beta_inf = 0.25;
319: armP->sigma = 1e-4;
320: armP->memorySize = 1;
321: armP->referencePolicy = REFERENCE_MAX;
322: armP->replacementPolicy = REPLACE_MRU;
323: armP->nondescending=PETSC_FALSE;
324: ls->data = (void*)armP;
325: ls->initstep=0.1;
326: ls->ops->setup=0;
327: ls->ops->reset=0;
328: ls->ops->apply=TaoLineSearchApply_OWArmijo;
329: ls->ops->view = TaoLineSearchView_OWArmijo;
330: ls->ops->destroy = TaoLineSearchDestroy_OWArmijo;
331: ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo;
332: return(0);
333: }