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