xref: /petsc/src/tao/linesearch/impls/armijo/armijo.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 #include <petsc/private/taolinesearchimpl.h>
2 #include <../src/tao/linesearch/impls/armijo/armijo.h>
3 
4 #define REPLACE_FIFO 1
5 #define REPLACE_MRU  2
6 
7 #define REFERENCE_MAX  1
8 #define REFERENCE_AVE  2
9 #define REFERENCE_MEAN 3
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "TaoLineSearchDestroy_Armijo"
13 static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls)
14 {
15   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
16   PetscErrorCode       ierr;
17 
18   PetscFunctionBegin;
19   ierr = PetscFree(armP->memory);CHKERRQ(ierr);
20   ierr = VecDestroy(&armP->x);CHKERRQ(ierr);
21   ierr = VecDestroy(&armP->work);CHKERRQ(ierr);
22   ierr = PetscFree(ls->data);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "TaoLineSearchReset_Armijo"
28 static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
29 {
30   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
31   PetscErrorCode       ierr;
32 
33   PetscFunctionBegin;
34   if (armP->memory != NULL) {
35     ierr = PetscFree(armP->memory);CHKERRQ(ierr);
36   }
37   armP->memorySetup = PETSC_FALSE;
38   PetscFunctionReturn(0);
39 }
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "TaoLineSearchSetFromOptions_Armijo"
43 static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptions *PetscOptionsObject,TaoLineSearch ls)
44 {
45   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
46   PetscErrorCode       ierr;
47 
48   PetscFunctionBegin;
49   ierr = PetscOptionsHead(PetscOptionsObject,"Armijo linesearch options");CHKERRQ(ierr);
50   ierr = PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL);CHKERRQ(ierr);
51   ierr = PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL);CHKERRQ(ierr);
52   ierr = PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL);CHKERRQ(ierr);
53   ierr = PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL);CHKERRQ(ierr);
54   ierr = PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL);CHKERRQ(ierr);
55   ierr = PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL);CHKERRQ(ierr);
56   ierr = PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL);CHKERRQ(ierr);
57   ierr = PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending,NULL);CHKERRQ(ierr);
58   ierr = PetscOptionsTail();CHKERRQ(ierr);
59   PetscFunctionReturn(0);
60 }
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "TaoLineSearchView_Armijo"
64 static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
65 {
66   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
67   PetscBool            isascii;
68   PetscErrorCode       ierr;
69 
70   PetscFunctionBegin;
71   ierr = PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
72   if (isascii) {
73     ierr = PetscViewerASCIIPrintf(pv,"  maxf=%D, ftol=%g, gtol=%g\n",ls->max_funcs, (double)ls->rtol, (double)ls->ftol);CHKERRQ(ierr);
74     ierr=PetscViewerASCIIPrintf(pv,"  Armijo linesearch",armP->alpha);CHKERRQ(ierr);
75     if (armP->nondescending) {
76       ierr = PetscViewerASCIIPrintf(pv, " (nondescending)");CHKERRQ(ierr);
77     }
78     if (ls->bounded) {
79       ierr = PetscViewerASCIIPrintf(pv," (projected)");CHKERRQ(ierr);
80     }
81     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);CHKERRQ(ierr);
82     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);CHKERRQ(ierr);
83     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);CHKERRQ(ierr);
84   }
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "TaoLineSearchApply_Armijo"
90 /* @ TaoApply_Armijo - This routine performs a linesearch. It
91    backtracks until the (nonmonotone) Armijo conditions are satisfied.
92 
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
101 
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
107 
108    Info is set to one of:
109 .   0 - the line search succeeds; the sufficient decrease
110    condition and the directional derivative condition hold
111 
112    negative number if an input parameter is invalid
113 -   -1 -  step < 0
114 
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 */
126 
127   PetscFunctionBegin;
128 
129   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
130   if (!armP->work) {
131     ierr = VecDuplicate(x,&armP->work);CHKERRQ(ierr);
132     armP->x = x;
133     ierr = PetscObjectReference((PetscObject)armP->x);CHKERRQ(ierr);
134   } else if (x != armP->x) {
135     /* If x has changed, then recreate work */
136     ierr = VecDestroy(&armP->work);CHKERRQ(ierr);
137     ierr = VecDuplicate(x,&armP->work);CHKERRQ(ierr);
138     ierr = PetscObjectDereference((PetscObject)armP->x);CHKERRQ(ierr);
139     armP->x = x;
140     ierr = PetscObjectReference((PetscObject)armP->x);CHKERRQ(ierr);
141   }
142 
143   /* Check linesearch parameters */
144   if (armP->alpha < 1) {
145     ierr = PetscInfo1(ls,"Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha);CHKERRQ(ierr);
146     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
147   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
148     ierr = PetscInfo1(ls,"Armijo line search error: beta (%g) invalid\n", (double)armP->beta);CHKERRQ(ierr);
149     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
150   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
151     ierr = PetscInfo1(ls,"Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);CHKERRQ(ierr);
152     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
153   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
154     ierr = PetscInfo1(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma);CHKERRQ(ierr);
155     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
156   } else if (armP->memorySize < 1) {
157     ierr = PetscInfo1(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize);CHKERRQ(ierr);
158     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
159   } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
160     ierr = PetscInfo(ls,"Armijo line search error: reference_policy invalid\n");CHKERRQ(ierr);
161     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
162   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
163     ierr = PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n");CHKERRQ(ierr);
164     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
165   } else if (PetscIsInfOrNanReal(*f)) {
166     ierr = PetscInfo(ls,"Armijo line search error: initial function inf or nan\n");CHKERRQ(ierr);
167     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
168   }
169 
170   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
171     PetscFunctionReturn(0);
172   }
173 
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     ierr = PetscMalloc1(armP->memorySize, &armP->memory );CHKERRQ(ierr);
179   }
180 
181   if (!armP->memorySetup) {
182     for (i = 0; i < armP->memorySize; i++) {
183       armP->memory[i] = armP->alpha*(*f);
184     }
185 
186     armP->current = 0;
187     armP->lastReference = armP->memory[0];
188     armP->memorySetup=PETSC_TRUE;
189   }
190 
191   /* Calculate reference value (MAX) */
192   ref = armP->memory[0];
193   idx = 0;
194 
195   for (i = 1; i < armP->memorySize; i++) {
196     if (armP->memory[i] > ref) {
197       ref = armP->memory[i];
198       idx = i;
199     }
200   }
201 
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   ierr = VecDot(g,s,&gdx);CHKERRQ(ierr);
213 
214   if (PetscIsInfOrNanReal(gdx)) {
215     ierr = PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);CHKERRQ(ierr);
216     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
217     PetscFunctionReturn(0);
218   }
219   if (gdx >= 0.0) {
220     ierr = PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);CHKERRQ(ierr);
221     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
222     PetscFunctionReturn(0);
223   }
224 
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     ierr = VecCopy(x,armP->work);CHKERRQ(ierr);
234     ierr = VecAXPY(armP->work,ls->step,s);CHKERRQ(ierr);
235     if (ls->bounded) {
236       ierr = VecMedian(ls->lower,armP->work,ls->upper,armP->work);CHKERRQ(ierr);
237     }
238 
239     /* Calculate function at new iterate */
240     if (ls->hasobjective) {
241       ierr = TaoLineSearchComputeObjective(ls,armP->work,f);CHKERRQ(ierr);
242       g_computed=PETSC_FALSE;
243     } else if (ls->usegts) {
244       ierr = TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx);CHKERRQ(ierr);
245       g_computed=PETSC_FALSE;
246     } else {
247       ierr = TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);CHKERRQ(ierr);
248       g_computed=PETSC_TRUE;
249     }
250     if (ls->step == ls->initstep) {
251       ls->f_fullstep = *f;
252     }
253 
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       }
263 
264       ls->step *= armP->beta;
265     }
266   }
267 
268   /* Check termination */
269   if (PetscIsInfOrNanReal(*f)) {
270     ierr = PetscInfo(ls, "Function is inf or nan.\n");CHKERRQ(ierr);
271     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
272   } else if (ls->step < ls->stepmin) {
273     ierr = PetscInfo(ls, "Step length is below tolerance.\n");CHKERRQ(ierr);
274     ls->reason = TAOLINESEARCH_HALTED_RTOL;
275   } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
276     ierr = PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs);CHKERRQ(ierr);
277     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
278   }
279   if (ls->reason) {
280     PetscFunctionReturn(0);
281   }
282 
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   }
294 
295   /* Update iterate and compute gradient */
296   ierr = VecCopy(armP->work,x);CHKERRQ(ierr);
297   if (!g_computed) {
298     ierr = TaoLineSearchComputeGradient(ls, x, g);CHKERRQ(ierr);
299   }
300   ierr = PetscInfo2(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "TaoLineSearchCreate_Armijo"
306 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
307 {
308   TaoLineSearch_ARMIJO *armP;
309   PetscErrorCode       ierr;
310 
311   PetscFunctionBegin;
312   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
313   ierr = PetscNewLog(ls,&armP);CHKERRQ(ierr);
314 
315   armP->memory = NULL;
316   armP->alpha = 1.0;
317   armP->beta = 0.5;
318   armP->beta_inf = 0.5;
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=1.0;
326   ls->ops->setup=0;
327   ls->ops->apply=TaoLineSearchApply_Armijo;
328   ls->ops->view = TaoLineSearchView_Armijo;
329   ls->ops->destroy = TaoLineSearchDestroy_Armijo;
330   ls->ops->reset = TaoLineSearchReset_Armijo;
331   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
332   PetscFunctionReturn(0);
333 }
334 
335