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