Actual source code: owarmijo.c
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: PetscInt low,high,low1,high1,low2,high2,i;
18: VecGetOwnershipRange(w,&low,&high);
19: VecGetOwnershipRange(x,&low1,&high1);
20: VecGetOwnershipRange(gv,&low2,&high2);
22: *gdx=0.0;
23: VecGetArray(w,&wptr);
24: VecGetArrayRead(x,&xptr);
25: VecGetArrayRead(gv,&gptr);
27: for (i=0;i<high-low;i++) {
28: if (xptr[i]*wptr[i]<0.0) wptr[i]=0.0;
29: *gdx = *gdx + gptr[i]*(wptr[i]-xptr[i]);
30: }
31: VecRestoreArray(w,&wptr);
32: VecRestoreArrayRead(x,&xptr);
33: VecRestoreArrayRead(gv,&gptr);
34: return 0;
35: }
37: static PetscErrorCode TaoLineSearchDestroy_OWArmijo(TaoLineSearch ls)
38: {
39: TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
41: PetscFree(armP->memory);
42: if (armP->x) {
43: PetscObjectDereference((PetscObject)armP->x);
44: }
45: VecDestroy(&armP->work);
46: PetscFree(ls->data);
47: return 0;
48: }
50: static PetscErrorCode TaoLineSearchSetFromOptions_OWArmijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
51: {
52: TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
54: PetscOptionsHead(PetscOptionsObject,"OWArmijo linesearch options");
55: PetscOptionsReal("-tao_ls_OWArmijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL);
56: PetscOptionsReal("-tao_ls_OWArmijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL);
57: PetscOptionsReal("-tao_ls_OWArmijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL);
58: PetscOptionsReal("-tao_ls_OWArmijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL);
59: PetscOptionsInt("-tao_ls_OWArmijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL);
60: PetscOptionsInt("-tao_ls_OWArmijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL);
61: PetscOptionsInt("-tao_ls_OWArmijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL);
62: PetscOptionsBool("-tao_ls_OWArmijo_nondescending","Use nondescending OWArmijo algorithm","",armP->nondescending,&armP->nondescending,NULL);
63: PetscOptionsTail();
64: return 0;
65: }
67: static PetscErrorCode TaoLineSearchView_OWArmijo(TaoLineSearch ls, PetscViewer pv)
68: {
69: TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
70: PetscBool isascii;
72: PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
73: if (isascii) {
74: PetscViewerASCIIPrintf(pv," OWArmijo linesearch",armP->alpha);
75: if (armP->nondescending) {
76: PetscViewerASCIIPrintf(pv, " (nondescending)");
77: }
78: PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
79: PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
80: PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
81: }
82: return 0;
83: }
85: /* @ TaoApply_OWArmijo - This routine performs a linesearch. It
86: backtracks until the (nonmonotone) OWArmijo conditions are satisfied.
88: Input Parameters:
89: + tao - TAO_SOLVER context
90: . X - current iterate (on output X contains new iterate, X + step*S)
91: . S - search direction
92: . f - merit function evaluated at X
93: . G - gradient of merit function evaluated at X
94: . W - work vector
95: - step - initial estimate of step length
97: Output parameters:
98: + f - merit function evaluated at new iterate, X + step*S
99: . G - gradient of merit function evaluated at new iterate, X + step*S
100: . X - new iterate
101: - step - final step length
103: Info is set to one of:
104: . 0 - the line search succeeds; the sufficient decrease
105: condition and the directional derivative condition hold
107: negative number if an input parameter is invalid
108: - -1 - step < 0
110: positive number > 1 if the line search otherwise terminates
111: + 1 - Step is at the lower bound, stepmin.
112: @ */
113: static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
114: {
115: TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
116: PetscInt i, its=0;
117: PetscReal fact, ref, gdx;
118: PetscInt idx;
119: PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
120: Vec g_old;
121: PetscReal owlqn_minstep=0.005;
122: PetscReal partgdx;
123: MPI_Comm comm;
125: PetscObjectGetComm((PetscObject)ls,&comm);
126: fact = 0.0;
127: ls->nfeval=0;
128: ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
129: if (!armP->work) {
130: VecDuplicate(x,&armP->work);
131: armP->x = x;
132: PetscObjectReference((PetscObject)armP->x);
133: } else if (x != armP->x) {
134: VecDestroy(&armP->work);
135: VecDuplicate(x,&armP->work);
136: PetscObjectDereference((PetscObject)armP->x);
137: armP->x = x;
138: PetscObjectReference((PetscObject)armP->x);
139: }
141: TaoLineSearchMonitor(ls, 0, *f, 0.0);
143: /* Check linesearch parameters */
144: if (armP->alpha < 1) {
145: PetscInfo(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha);
146: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
147: } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
148: PetscInfo(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta);
149: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
150: } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
151: PetscInfo(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);
152: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
153: } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
154: PetscInfo(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
155: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
156: } else if (armP->memorySize < 1) {
157: PetscInfo(ls,"OWArmijo line search error: memory_size (%D) < 1\n", armP->memorySize);
158: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
159: } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
160: PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n");
161: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
162: } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
163: PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n");
164: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
165: } else if (PetscIsInfOrNanReal(*f)) {
166: PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n");
167: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
168: }
170: if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) return 0;
172: /* Check to see of the memory has been allocated. If not, allocate
173: the historical array and populate it with the initial function
174: values. */
175: if (!armP->memory) {
176: PetscMalloc1(armP->memorySize, &armP->memory);
177: }
179: if (!armP->memorySetup) {
180: for (i = 0; i < armP->memorySize; i++) {
181: armP->memory[i] = armP->alpha*(*f);
182: }
183: armP->current = 0;
184: armP->lastReference = armP->memory[0];
185: armP->memorySetup=PETSC_TRUE;
186: }
188: /* Calculate reference value (MAX) */
189: ref = armP->memory[0];
190: idx = 0;
192: for (i = 1; i < armP->memorySize; i++) {
193: if (armP->memory[i] > ref) {
194: ref = armP->memory[i];
195: idx = i;
196: }
197: }
199: if (armP->referencePolicy == REFERENCE_AVE) {
200: ref = 0;
201: for (i = 0; i < armP->memorySize; i++) {
202: ref += armP->memory[i];
203: }
204: ref = ref / armP->memorySize;
205: ref = PetscMax(ref, armP->memory[armP->current]);
206: } else if (armP->referencePolicy == REFERENCE_MEAN) {
207: ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
208: }
210: if (armP->nondescending) {
211: fact = armP->sigma;
212: }
214: VecDuplicate(g,&g_old);
215: VecCopy(g,g_old);
217: ls->step = ls->initstep;
218: while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) {
219: /* Calculate iterate */
220: ++its;
221: VecCopy(x,armP->work);
222: VecAXPY(armP->work,ls->step,s);
224: partgdx=0.0;
225: ProjWork_OWLQN(armP->work,x,g_old,&partgdx);
226: MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm);
228: /* Check the condition of gdx */
229: if (PetscIsInfOrNanReal(gdx)) {
230: PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
231: ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
232: return 0;
233: }
234: if (gdx >= 0.0) {
235: PetscInfo(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
236: ls->reason = TAOLINESEARCH_FAILED_ASCENT;
237: return 0;
238: }
240: /* Calculate function at new iterate */
241: TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
242: g_computed=PETSC_TRUE;
244: TaoLineSearchMonitor(ls, its, *f, ls->step);
246: if (ls->step == ls->initstep) {
247: ls->f_fullstep = *f;
248: }
250: if (PetscIsInfOrNanReal(*f)) {
251: ls->step *= armP->beta_inf;
252: } else {
253: /* Check descent condition */
254: if (armP->nondescending && *f <= ref - ls->step*fact*ref) break;
255: if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break;
256: ls->step *= armP->beta;
257: }
258: }
259: VecDestroy(&g_old);
261: /* Check termination */
262: if (PetscIsInfOrNanReal(*f)) {
263: PetscInfo(ls, "Function is inf or nan.\n");
264: ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER;
265: } else if (ls->step < owlqn_minstep) {
266: PetscInfo(ls, "Step length is below tolerance.\n");
267: ls->reason = TAOLINESEARCH_HALTED_RTOL;
268: } else if (ls->nfeval >= ls->max_funcs) {
269: PetscInfo(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs);
270: ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
271: }
272: if (ls->reason) return 0;
274: /* Successful termination, update memory */
275: ls->reason = TAOLINESEARCH_SUCCESS;
276: armP->lastReference = ref;
277: if (armP->replacementPolicy == REPLACE_FIFO) {
278: armP->memory[armP->current++] = *f;
279: if (armP->current >= armP->memorySize) {
280: armP->current = 0;
281: }
282: } else {
283: armP->current = idx;
284: armP->memory[idx] = *f;
285: }
287: /* Update iterate and compute gradient */
288: VecCopy(armP->work,x);
289: if (!g_computed) {
290: TaoLineSearchComputeGradient(ls, x, g);
291: }
292: PetscInfo(ls, "%D function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step);
293: return 0;
294: }
296: /*MC
297: TAOLINESEARCHOWARMIJO - Special line-search type for the Orthant-Wise Limited Quasi-Newton (TAOOWLQN) algorithm.
298: Should not be used with any other algorithm.
300: Level: developer
302: .keywords: Tao, linesearch
303: M*/
304: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls)
305: {
306: TaoLineSearch_OWARMIJO *armP;
309: PetscNewLog(ls,&armP);
311: armP->memory = NULL;
312: armP->alpha = 1.0;
313: armP->beta = 0.25;
314: armP->beta_inf = 0.25;
315: armP->sigma = 1e-4;
316: armP->memorySize = 1;
317: armP->referencePolicy = REFERENCE_MAX;
318: armP->replacementPolicy = REPLACE_MRU;
319: armP->nondescending=PETSC_FALSE;
320: ls->data = (void*)armP;
321: ls->initstep = 0.1;
322: ls->ops->monitor = NULL;
323: ls->ops->setup = NULL;
324: ls->ops->reset = NULL;
325: ls->ops->apply = TaoLineSearchApply_OWArmijo;
326: ls->ops->view = TaoLineSearchView_OWArmijo;
327: ls->ops->destroy = TaoLineSearchDestroy_OWArmijo;
328: ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo;
329: return 0;
330: }