xref: /petsc/src/tao/linesearch/impls/armijo/armijo.c (revision 2a0dac0797dfdf4a2942fb34b57c2a93ae7a8c7a)
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 static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls)
12 {
13   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
14   PetscErrorCode       ierr;
15 
16   PetscFunctionBegin;
17   ierr = PetscFree(armP->memory);CHKERRQ(ierr);
18   ierr = VecDestroy(&armP->x);CHKERRQ(ierr);
19   ierr = VecDestroy(&armP->work);CHKERRQ(ierr);
20   ierr = PetscFree(ls->data);CHKERRQ(ierr);
21   PetscFunctionReturn(0);
22 }
23 
24 static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
25 {
26   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
27   PetscErrorCode       ierr;
28 
29   PetscFunctionBegin;
30   ierr = PetscFree(armP->memory);CHKERRQ(ierr);
31   armP->memorySetup = PETSC_FALSE;
32   PetscFunctionReturn(0);
33 }
34 
35 static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
36 {
37   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
38   PetscErrorCode       ierr;
39 
40   PetscFunctionBegin;
41   ierr = PetscOptionsHead(PetscOptionsObject,"Armijo linesearch options");CHKERRQ(ierr);
42   ierr = PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL);CHKERRQ(ierr);
43   ierr = PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL);CHKERRQ(ierr);
44   ierr = PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL);CHKERRQ(ierr);
45   ierr = PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL);CHKERRQ(ierr);
46   ierr = PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL);CHKERRQ(ierr);
47   ierr = PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL);CHKERRQ(ierr);
48   ierr = PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL);CHKERRQ(ierr);
49   ierr = PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending,NULL);CHKERRQ(ierr);
50   ierr = PetscOptionsTail();CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
54 static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
55 {
56   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
57   PetscBool            isascii;
58   PetscErrorCode       ierr;
59 
60   PetscFunctionBegin;
61   ierr = PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
62   if (isascii) {
63     ierr=PetscViewerASCIIPrintf(pv,"  Armijo linesearch",armP->alpha);CHKERRQ(ierr);
64     if (armP->nondescending) {
65       ierr = PetscViewerASCIIPrintf(pv, " (nondescending)");CHKERRQ(ierr);
66     }
67     if (ls->bounded) {
68       ierr = PetscViewerASCIIPrintf(pv," (projected)");CHKERRQ(ierr);
69     }
70     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);CHKERRQ(ierr);
71     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);CHKERRQ(ierr);
72     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);CHKERRQ(ierr);
73   }
74   PetscFunctionReturn(0);
75 }
76 
77 /* @ TaoApply_Armijo - This routine performs a linesearch. It
78    backtracks until the (nonmonotone) Armijo conditions are satisfied.
79 
80    Input Parameters:
81 +  tao - Tao context
82 .  X - current iterate (on output X contains new iterate, X + step*S)
83 .  S - search direction
84 .  f - merit function evaluated at X
85 .  G - gradient of merit function evaluated at X
86 .  W - work vector
87 -  step - initial estimate of step length
88 
89    Output parameters:
90 +  f - merit function evaluated at new iterate, X + step*S
91 .  G - gradient of merit function evaluated at new iterate, X + step*S
92 .  X - new iterate
93 -  step - final step length
94 
95 @ */
96 static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
97 {
98   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
99   PetscErrorCode       ierr;
100   PetscInt             i, its=0;
101   PetscReal            fact, ref, gdx;
102   PetscInt             idx;
103   PetscBool            g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
104 
105   PetscFunctionBegin;
106 
107   ierr = TaoLineSearchMonitor(ls, 0, *f, 0.0);CHKERRQ(ierr);
108 
109   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
110   if (!armP->work) {
111     ierr = VecDuplicate(x,&armP->work);CHKERRQ(ierr);
112     armP->x = x;
113     ierr = PetscObjectReference((PetscObject)armP->x);CHKERRQ(ierr);
114   } else if (x != armP->x) {
115     /* If x has changed, then recreate work */
116     ierr = VecDestroy(&armP->work);CHKERRQ(ierr);
117     ierr = VecDuplicate(x,&armP->work);CHKERRQ(ierr);
118     ierr = PetscObjectDereference((PetscObject)armP->x);CHKERRQ(ierr);
119     armP->x = x;
120     ierr = PetscObjectReference((PetscObject)armP->x);CHKERRQ(ierr);
121   }
122 
123   /* Check linesearch parameters */
124   if (armP->alpha < 1) {
125     ierr = PetscInfo1(ls,"Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha);CHKERRQ(ierr);
126     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
127   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
128     ierr = PetscInfo1(ls,"Armijo line search error: beta (%g) invalid\n", (double)armP->beta);CHKERRQ(ierr);
129     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
130   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
131     ierr = PetscInfo1(ls,"Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);CHKERRQ(ierr);
132     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
133   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
134     ierr = PetscInfo1(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma);CHKERRQ(ierr);
135     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
136   } else if (armP->memorySize < 1) {
137     ierr = PetscInfo1(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize);CHKERRQ(ierr);
138     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
139   } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
140     ierr = PetscInfo(ls,"Armijo line search error: reference_policy invalid\n");CHKERRQ(ierr);
141     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
142   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
143     ierr = PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n");CHKERRQ(ierr);
144     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
145   } else if (PetscIsInfOrNanReal(*f)) {
146     ierr = PetscInfo(ls,"Armijo line search error: initial function inf or nan\n");CHKERRQ(ierr);
147     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
148   }
149 
150   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
151     PetscFunctionReturn(0);
152   }
153 
154   /* Check to see of the memory has been allocated.  If not, allocate
155      the historical array and populate it with the initial function
156      values. */
157   if (!armP->memory) {
158     ierr = PetscMalloc1(armP->memorySize, &armP->memory );CHKERRQ(ierr);
159   }
160 
161   if (!armP->memorySetup) {
162     for (i = 0; i < armP->memorySize; i++) {
163       armP->memory[i] = armP->alpha*(*f);
164     }
165 
166     armP->current = 0;
167     armP->lastReference = armP->memory[0];
168     armP->memorySetup=PETSC_TRUE;
169   }
170 
171   /* Calculate reference value (MAX) */
172   ref = armP->memory[0];
173   idx = 0;
174 
175   for (i = 1; i < armP->memorySize; i++) {
176     if (armP->memory[i] > ref) {
177       ref = armP->memory[i];
178       idx = i;
179     }
180   }
181 
182   if (armP->referencePolicy == REFERENCE_AVE) {
183     ref = 0;
184     for (i = 0; i < armP->memorySize; i++) {
185       ref += armP->memory[i];
186     }
187     ref = ref / armP->memorySize;
188     ref = PetscMax(ref, armP->memory[armP->current]);
189   } else if (armP->referencePolicy == REFERENCE_MEAN) {
190     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
191   }
192   ierr = VecDot(g,s,&gdx);CHKERRQ(ierr);
193 
194   if (PetscIsInfOrNanReal(gdx)) {
195     ierr = PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);CHKERRQ(ierr);
196     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
197     PetscFunctionReturn(0);
198   }
199   if (gdx >= 0.0) {
200     ierr = PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);CHKERRQ(ierr);
201     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
202     PetscFunctionReturn(0);
203   }
204 
205   if (armP->nondescending) {
206     fact = armP->sigma;
207   } else {
208     fact = armP->sigma * gdx;
209   }
210   ls->step = ls->initstep;
211   while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) {
212     /* Calculate iterate */
213     ++its;
214     ierr = VecCopy(x,armP->work);CHKERRQ(ierr);
215     ierr = VecAXPY(armP->work,ls->step,s);CHKERRQ(ierr);
216     if (ls->bounded) {
217       ierr = VecMedian(ls->lower,armP->work,ls->upper,armP->work);CHKERRQ(ierr);
218     }
219 
220     /* Calculate function at new iterate */
221     if (ls->hasobjective) {
222       ierr = TaoLineSearchComputeObjective(ls,armP->work,f);CHKERRQ(ierr);
223       g_computed=PETSC_FALSE;
224     } else if (ls->usegts) {
225       ierr = TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx);CHKERRQ(ierr);
226       g_computed=PETSC_FALSE;
227     } else {
228       ierr = TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);CHKERRQ(ierr);
229       g_computed=PETSC_TRUE;
230     }
231     if (ls->step == ls->initstep) {
232       ls->f_fullstep = *f;
233     }
234 
235     /* call monitor */
236     ierr = TaoLineSearchMonitor(ls, its, *f, ls->step);CHKERRQ(ierr);
237 
238     if (PetscIsInfOrNanReal(*f)) {
239       ls->step *= armP->beta_inf;
240     } else {
241       /* Check descent condition */
242       if (armP->nondescending && *f <= ref - ls->step*fact*ref)
243         break;
244       if (!armP->nondescending && *f <= ref + ls->step*fact) {
245         break;
246       }
247 
248       ls->step *= armP->beta;
249     }
250   }
251 
252   /* Check termination */
253   if (PetscIsInfOrNanReal(*f)) {
254     ierr = PetscInfo(ls, "Function is inf or nan.\n");CHKERRQ(ierr);
255     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
256   } else if (ls->step < ls->stepmin) {
257     ierr = PetscInfo(ls, "Step length is below tolerance.\n");CHKERRQ(ierr);
258     ls->reason = TAOLINESEARCH_HALTED_RTOL;
259   } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
260     ierr = PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs);CHKERRQ(ierr);
261     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
262   }
263   if (ls->reason) {
264     PetscFunctionReturn(0);
265   }
266 
267   /* Successful termination, update memory */
268   ls->reason = TAOLINESEARCH_SUCCESS;
269   armP->lastReference = ref;
270   if (armP->replacementPolicy == REPLACE_FIFO) {
271     armP->memory[armP->current++] = *f;
272     if (armP->current >= armP->memorySize) {
273       armP->current = 0;
274     }
275   } else {
276     armP->current = idx;
277     armP->memory[idx] = *f;
278   }
279 
280   /* Update iterate and compute gradient */
281   ierr = VecCopy(armP->work,x);CHKERRQ(ierr);
282   if (!g_computed) {
283     ierr = TaoLineSearchComputeGradient(ls, x, g);CHKERRQ(ierr);
284   }
285   ierr = PetscInfo2(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step);CHKERRQ(ierr);
286   PetscFunctionReturn(0);
287 }
288 
289 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
290 {
291   TaoLineSearch_ARMIJO *armP;
292   PetscErrorCode       ierr;
293 
294   PetscFunctionBegin;
295   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
296   ierr = PetscNewLog(ls,&armP);CHKERRQ(ierr);
297 
298   armP->memory = NULL;
299   armP->alpha = 1.0;
300   armP->beta = 0.5;
301   armP->beta_inf = 0.5;
302   armP->sigma = 1e-4;
303   armP->memorySize = 1;
304   armP->referencePolicy = REFERENCE_MAX;
305   armP->replacementPolicy = REPLACE_MRU;
306   armP->nondescending=PETSC_FALSE;
307   ls->data = (void*)armP;
308   ls->initstep=1.0;
309   ls->ops->setup=0;
310   ls->ops->apply=TaoLineSearchApply_Armijo;
311   ls->ops->view = TaoLineSearchView_Armijo;
312   ls->ops->destroy = TaoLineSearchDestroy_Armijo;
313   ls->ops->reset = TaoLineSearchReset_Armijo;
314   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
315   PetscFunctionReturn(0);
316 }
317 
318