xref: /petsc/src/tao/linesearch/impls/owarmijo/owarmijo.c (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
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(PetscOptionItems *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,"  OWArmijo linesearch",armP->alpha);CHKERRQ(ierr);
89     if (armP->nondescending) {
90       ierr = PetscViewerASCIIPrintf(pv, " (nondescending)");CHKERRQ(ierr);
91     }
92     ierr=PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta);CHKERRQ(ierr);
93     ierr=PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma);CHKERRQ(ierr);
94     ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);CHKERRQ(ierr);
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 #undef __FUNCT__
100 #define __FUNCT__ "TaoLineSearchApply_OWArmijo"
101 /* @ TaoApply_OWArmijo - This routine performs a linesearch. It
102    backtracks until the (nonmonotone) OWArmijo conditions are satisfied.
103 
104    Input Parameters:
105 +  tao - TAO_SOLVER context
106 .  X - current iterate (on output X contains new iterate, X + step*S)
107 .  S - search direction
108 .  f - merit function evaluated at X
109 .  G - gradient of merit function evaluated at X
110 .  W - work vector
111 -  step - initial estimate of step length
112 
113    Output parameters:
114 +  f - merit function evaluated at new iterate, X + step*S
115 .  G - gradient of merit function evaluated at new iterate, X + step*S
116 .  X - new iterate
117 -  step - final step length
118 
119    Info is set to one of:
120 .   0 - the line search succeeds; the sufficient decrease
121    condition and the directional derivative condition hold
122 
123    negative number if an input parameter is invalid
124 -   -1 -  step < 0
125 
126    positive number > 1 if the line search otherwise terminates
127 +    1 -  Step is at the lower bound, stepmin.
128 @ */
129 static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
130 {
131   TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data;
132   PetscErrorCode         ierr;
133   PetscInt               i;
134   PetscReal              fact, ref, gdx;
135   PetscInt               idx;
136   PetscBool              g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
137   Vec                    g_old;
138   PetscReal              owlqn_minstep=0.005;
139   PetscReal              partgdx;
140   MPI_Comm               comm;
141 
142   PetscFunctionBegin;
143   ierr = PetscObjectGetComm((PetscObject)ls,&comm);CHKERRQ(ierr);
144   fact = 0.0;
145   ls->nfeval=0;
146   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
147   if (!armP->work) {
148     ierr = VecDuplicate(x,&armP->work);CHKERRQ(ierr);
149     armP->x = x;
150     ierr = PetscObjectReference((PetscObject)armP->x);CHKERRQ(ierr);
151   } else if (x != armP->x) {
152     ierr = VecDestroy(&armP->work);CHKERRQ(ierr);
153     ierr = VecDuplicate(x,&armP->work);CHKERRQ(ierr);
154     ierr = PetscObjectDereference((PetscObject)armP->x);CHKERRQ(ierr);
155     armP->x = x;
156     ierr = PetscObjectReference((PetscObject)armP->x);CHKERRQ(ierr);
157   }
158 
159   /* Check linesearch parameters */
160   if (armP->alpha < 1) {
161     ierr = PetscInfo1(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha);CHKERRQ(ierr);
162     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
163   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
164     ierr = PetscInfo1(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta);CHKERRQ(ierr);
165     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
166   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
167     ierr = PetscInfo1(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf);CHKERRQ(ierr);
168     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
169   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
170     ierr = PetscInfo1(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma);CHKERRQ(ierr);
171     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
172   } else if (armP->memorySize < 1) {
173     ierr = PetscInfo1(ls,"OWArmijo line search error: memory_size (%D) < 1\n", armP->memorySize);CHKERRQ(ierr);
174     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
175   }  else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
176     ierr = PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n");CHKERRQ(ierr);
177     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
178   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
179     ierr = PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n");CHKERRQ(ierr);
180     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
181   } else if (PetscIsInfOrNanReal(*f)) {
182     ierr = PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n");CHKERRQ(ierr);
183     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
184   }
185 
186   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) PetscFunctionReturn(0);
187 
188   /* Check to see of the memory has been allocated.  If not, allocate
189      the historical array and populate it with the initial function
190      values. */
191   if (!armP->memory) {
192     ierr = PetscMalloc1(armP->memorySize, &armP->memory );CHKERRQ(ierr);
193   }
194 
195   if (!armP->memorySetup) {
196     for (i = 0; i < armP->memorySize; i++) {
197       armP->memory[i] = armP->alpha*(*f);
198     }
199     armP->current = 0;
200     armP->lastReference = armP->memory[0];
201     armP->memorySetup=PETSC_TRUE;
202   }
203 
204   /* Calculate reference value (MAX) */
205   ref = armP->memory[0];
206   idx = 0;
207 
208   for (i = 1; i < armP->memorySize; i++) {
209     if (armP->memory[i] > ref) {
210       ref = armP->memory[i];
211       idx = i;
212     }
213   }
214 
215   if (armP->referencePolicy == REFERENCE_AVE) {
216     ref = 0;
217     for (i = 0; i < armP->memorySize; i++) {
218       ref += armP->memory[i];
219     }
220     ref = ref / armP->memorySize;
221     ref = PetscMax(ref, armP->memory[armP->current]);
222   } else if (armP->referencePolicy == REFERENCE_MEAN) {
223     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
224   }
225 
226   if (armP->nondescending) {
227     fact = armP->sigma;
228   }
229 
230   ierr = VecDuplicate(g,&g_old);CHKERRQ(ierr);
231   ierr = VecCopy(g,g_old);CHKERRQ(ierr);
232 
233   ls->step = ls->initstep;
234   while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) {
235     /* Calculate iterate */
236     ierr = VecCopy(x,armP->work);CHKERRQ(ierr);
237     ierr = VecAXPY(armP->work,ls->step,s);CHKERRQ(ierr);
238 
239     partgdx=0.0;
240     ierr = ProjWork_OWLQN(armP->work,x,g_old,&partgdx);CHKERRQ(ierr);
241     ierr = MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr);
242 
243     /* Check the condition of gdx */
244     if (PetscIsInfOrNanReal(gdx)) {
245       ierr = PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx);CHKERRQ(ierr);
246       ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
247       PetscFunctionReturn(0);
248     }
249     if (gdx >= 0.0) {
250       ierr = PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx);CHKERRQ(ierr);
251       ls->reason = TAOLINESEARCH_FAILED_ASCENT;
252       PetscFunctionReturn(0);
253     }
254 
255     /* Calculate function at new iterate */
256     ierr = TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g);CHKERRQ(ierr);
257     g_computed=PETSC_TRUE;
258 
259     if (ls->step == ls->initstep) {
260       ls->f_fullstep = *f;
261     }
262 
263     if (PetscIsInfOrNanReal(*f)) {
264       ls->step *= armP->beta_inf;
265     } else {
266       /* Check descent condition */
267       if (armP->nondescending && *f <= ref - ls->step*fact*ref) break;
268       if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break;
269       ls->step *= armP->beta;
270     }
271   }
272   ierr = VecDestroy(&g_old);CHKERRQ(ierr);
273 
274   /* Check termination */
275   if (PetscIsInfOrNanReal(*f)) {
276     ierr = PetscInfo(ls, "Function is inf or nan.\n");CHKERRQ(ierr);
277     ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER;
278   } else if (ls->step < owlqn_minstep) {
279     ierr = PetscInfo(ls, "Step length is below tolerance.\n");CHKERRQ(ierr);
280     ls->reason = TAOLINESEARCH_HALTED_RTOL;
281   } else if (ls->nfeval >= ls->max_funcs) {
282     ierr = PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs);CHKERRQ(ierr);
283     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
284   }
285   if (ls->reason) PetscFunctionReturn(0);
286 
287   /* Successful termination, update memory */
288   ls->reason = TAOLINESEARCH_SUCCESS;
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