Actual source code: armijo.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/taolinesearchimpl.h>
2: #include <../src/tao/linesearch/impls/armijo/armijo.h>
4: #define REPLACE_FIFO 1
5: #define REPLACE_MRU 2
7: #define REFERENCE_MAX 1
8: #define REFERENCE_AVE 2
9: #define REFERENCE_MEAN 3
13: static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls)
14: {
15: TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
16: PetscErrorCode ierr;
19: PetscFree(armP->memory);
20: VecDestroy(&armP->x);
21: VecDestroy(&armP->work);
22: PetscFree(ls->data);
23: return(0);
24: }
28: static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
29: {
30: TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
31: PetscErrorCode ierr;
34: if (armP->memory != NULL) {
35: PetscFree(armP->memory);
36: }
37: armP->memorySetup = PETSC_FALSE;
38: return(0);
39: }
43: static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(TaoLineSearch ls)
44: {
45: TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
46: PetscErrorCode ierr;
49: PetscOptionsHead("Armijo linesearch options");
50: PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha, 0);
51: PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf, 0);
52: PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta, 0);
53: PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma, 0);
54: PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize, 0);
55: PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy, 0);
56: PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy, 0);
57: PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending, 0);
58: PetscOptionsTail();
59: return(0);
60: }
64: static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
65: {
66: TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
67: PetscBool isascii;
68: PetscErrorCode ierr;
71: PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);
72: if (isascii) {
73: PetscViewerASCIIPrintf(pv," maxf=%D, ftol=%g, gtol=%g\n",ls->max_funcs, (double)ls->rtol, (double)ls->ftol);
74: ierr=PetscViewerASCIIPrintf(pv," Armijo linesearch",armP->alpha);
75: if (armP->nondescending) {
76: PetscViewerASCIIPrintf(pv, " (nondescending)");
77: }
78: if (ls->bounded) {
79: PetscViewerASCIIPrintf(pv," (projected)");
80: }
81: ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);
82: ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);
83: ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);
84: }
85: return(0);
86: }
90: /* @ TaoApply_Armijo - This routine performs a linesearch. It
91: backtracks until the (nonmonotone) Armijo conditions are satisfied.
93: Input Parameters:
94: + tao - Tao context
95: . X - current iterate (on output X contains new iterate, X + step*S)
96: . S - search direction
97: . f - merit function evaluated at X
98: . G - gradient of merit function evaluated at X
99: . W - work vector
100: - step - initial estimate of step length
102: Output parameters:
103: + f - merit function evaluated at new iterate, X + step*S
104: . G - gradient of merit function evaluated at new iterate, X + step*S
105: . X - new iterate
106: - step - final step length
108: Info is set to one of:
109: . 0 - the line search succeeds; the sufficient decrease
110: condition and the directional derivative condition hold
112: negative number if an input parameter is invalid
113: - -1 - step < 0
115: positive number > 1 if the line search otherwise terminates
116: + 1 - Step is at the lower bound, stepmin.
117: @ */
118: static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
119: {
120: TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
121: PetscErrorCode ierr;
122: PetscInt i;
123: PetscReal fact, ref, gdx;
124: PetscInt idx;
125: PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
129: ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
130: if (!armP->work) {
131: VecDuplicate(x,&armP->work);
132: armP->x = x;
133: PetscObjectReference((PetscObject)armP->x);
134: } else if (x != armP->x) {
135: /* If x has changed, then recreate work */
136: VecDestroy(&armP->work);
137: VecDuplicate(x,&armP->work);
138: PetscObjectDereference((PetscObject)armP->x);
139: armP->x = x;
140: PetscObjectReference((PetscObject)armP->x);
141: }
143: /* Check linesearch parameters */
144: if (armP->alpha < 1) {
145: PetscInfo1(ls,"Armijo 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: PetscInfo1(ls,"Armijo 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: PetscInfo1(ls,"Armijo 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: PetscInfo1(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma);
155: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
156: } else if (armP->memorySize < 1) {
157: PetscInfo1(ls,"Armijo 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,"Armijo 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,"Armijo line search error: replacement_policy invalid\n");
164: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
165: } else if (PetscIsInfOrNanReal(*f)) {
166: PetscInfo(ls,"Armijo line search error: initial function inf or nan\n");
167: ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
168: }
170: if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
171: return(0);
172: }
174: /* Check to see of the memory has been allocated. If not, allocate
175: the historical array and populate it with the initial function
176: values. */
177: if (!armP->memory) {
178: PetscMalloc1(armP->memorySize, &armP->memory );
179: }
181: if (!armP->memorySetup) {
182: for (i = 0; i < armP->memorySize; i++) {
183: armP->memory[i] = armP->alpha*(*f);
184: }
186: armP->current = 0;
187: armP->lastReference = armP->memory[0];
188: armP->memorySetup=PETSC_TRUE;
189: }
191: /* Calculate reference value (MAX) */
192: ref = armP->memory[0];
193: idx = 0;
195: for (i = 1; i < armP->memorySize; i++) {
196: if (armP->memory[i] > ref) {
197: ref = armP->memory[i];
198: idx = i;
199: }
200: }
202: if (armP->referencePolicy == REFERENCE_AVE) {
203: ref = 0;
204: for (i = 0; i < armP->memorySize; i++) {
205: ref += armP->memory[i];
206: }
207: ref = ref / armP->memorySize;
208: ref = PetscMax(ref, armP->memory[armP->current]);
209: } else if (armP->referencePolicy == REFERENCE_MEAN) {
210: ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
211: }
212: VecDot(g,s,&gdx);
214: if (PetscIsInfOrNanReal(gdx)) {
215: PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);
216: ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
217: return(0);
218: }
219: if (gdx >= 0.0) {
220: PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);
221: ls->reason = TAOLINESEARCH_FAILED_ASCENT;
222: return(0);
223: }
225: if (armP->nondescending) {
226: fact = armP->sigma;
227: } else {
228: fact = armP->sigma * gdx;
229: }
230: ls->step = ls->initstep;
231: while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) {
232: /* Calculate iterate */
233: VecCopy(x,armP->work);
234: VecAXPY(armP->work,ls->step,s);
235: if (ls->bounded) {
236: VecMedian(ls->lower,armP->work,ls->upper,armP->work);
237: }
239: /* Calculate function at new iterate */
240: if (ls->hasobjective) {
241: TaoLineSearchComputeObjective(ls,armP->work,f);
242: g_computed=PETSC_FALSE;
243: } else if (ls->usegts) {
244: TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx);
245: g_computed=PETSC_FALSE;
246: } else {
247: TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);
248: g_computed=PETSC_TRUE;
249: }
250: if (ls->step == ls->initstep) {
251: ls->f_fullstep = *f;
252: }
254: if (PetscIsInfOrNanReal(*f)) {
255: ls->step *= armP->beta_inf;
256: } else {
257: /* Check descent condition */
258: if (armP->nondescending && *f <= ref - ls->step*fact*ref)
259: break;
260: if (!armP->nondescending && *f <= ref + ls->step*fact) {
261: break;
262: }
264: ls->step *= armP->beta;
265: }
266: }
268: /* Check termination */
269: if (PetscIsInfOrNanReal(*f)) {
270: PetscInfo(ls, "Function is inf or nan.\n");
271: ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
272: } else if (ls->step < ls->stepmin) {
273: PetscInfo(ls, "Step length is below tolerance.\n");
274: ls->reason = TAOLINESEARCH_HALTED_RTOL;
275: } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
276: PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs);
277: ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
278: }
279: if (ls->reason) {
280: return(0);
281: }
283: /* Successful termination, update memory */
284: armP->lastReference = ref;
285: if (armP->replacementPolicy == REPLACE_FIFO) {
286: armP->memory[armP->current++] = *f;
287: if (armP->current >= armP->memorySize) {
288: armP->current = 0;
289: }
290: } else {
291: armP->current = idx;
292: armP->memory[idx] = *f;
293: }
295: /* Update iterate and compute gradient */
296: VecCopy(armP->work,x);
297: if (!g_computed) {
298: TaoLineSearchComputeGradient(ls, x, g);
299: }
300: PetscInfo2(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step);
301: return(0);
302: }
304: EXTERN_C_BEGIN
307: PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
308: {
309: TaoLineSearch_ARMIJO *armP;
310: PetscErrorCode ierr;
314: PetscNewLog(ls,&armP);
316: armP->memory = NULL;
317: armP->alpha = 1.0;
318: armP->beta = 0.5;
319: armP->beta_inf = 0.5;
320: armP->sigma = 1e-4;
321: armP->memorySize = 1;
322: armP->referencePolicy = REFERENCE_MAX;
323: armP->replacementPolicy = REPLACE_MRU;
324: armP->nondescending=PETSC_FALSE;
325: ls->data = (void*)armP;
326: ls->initstep=1.0;
327: ls->ops->setup=0;
328: ls->ops->apply=TaoLineSearchApply_Armijo;
329: ls->ops->view = TaoLineSearchView_Armijo;
330: ls->ops->destroy = TaoLineSearchDestroy_Armijo;
331: ls->ops->reset = TaoLineSearchReset_Armijo;
332: ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
333: return(0);
334: }
335: EXTERN_C_END