xref: /petsc/src/tao/linesearch/impls/armijo/armijo.c (revision ef46b1a67e276116c83b5d4ce8efc2932ea4fc0a)
1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/armijo/armijo.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay #define REPLACE_FIFO 1
5a7e14dcfSSatish Balay #define REPLACE_MRU  2
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay #define REFERENCE_MAX  1
8a7e14dcfSSatish Balay #define REFERENCE_AVE  2
9a7e14dcfSSatish Balay #define REFERENCE_MEAN 3
10a7e14dcfSSatish Balay 
11a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls)
12a7e14dcfSSatish Balay {
138caf6e8cSBarry Smith   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(PetscFree(armP->memory));
179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&armP->x));
189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&armP->work));
199566063dSJacob Faibussowitsch   PetscCall(PetscFree(ls->data));
20a7e14dcfSSatish Balay   PetscFunctionReturn(0);
21a7e14dcfSSatish Balay }
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls)
24a7e14dcfSSatish Balay {
258caf6e8cSBarry Smith   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
26a7e14dcfSSatish Balay 
27a7e14dcfSSatish Balay   PetscFunctionBegin;
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(armP->memory));
29a7e14dcfSSatish Balay   armP->memorySetup = PETSC_FALSE;
30a7e14dcfSSatish Balay   PetscFunctionReturn(0);
31a7e14dcfSSatish Balay }
32a7e14dcfSSatish Balay 
334416b707SBarry Smith static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
34a7e14dcfSSatish Balay {
358caf6e8cSBarry Smith   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay   PetscFunctionBegin;
38d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Armijo linesearch options");
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL));
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending,NULL));
47d0609cedSBarry Smith   PetscOptionsHeadEnd();
48a7e14dcfSSatish Balay   PetscFunctionReturn(0);
49a7e14dcfSSatish Balay }
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv)
52a7e14dcfSSatish Balay {
538caf6e8cSBarry Smith   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
54a7e14dcfSSatish Balay   PetscBool            isascii;
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii));
58a7e14dcfSSatish Balay   if (isascii) {
599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"  Armijo linesearch",armP->alpha));
60a7e14dcfSSatish Balay     if (armP->nondescending) {
619566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv, " (nondescending)"));
62a7e14dcfSSatish Balay     }
63a7e14dcfSSatish Balay     if (ls->bounded) {
649566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv," (projected)"));
65a7e14dcfSSatish Balay     }
669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta));
679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma));
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize));
69a7e14dcfSSatish Balay   }
70a7e14dcfSSatish Balay   PetscFunctionReturn(0);
71a7e14dcfSSatish Balay }
72a7e14dcfSSatish Balay 
73a7e14dcfSSatish Balay /* @ TaoApply_Armijo - This routine performs a linesearch. It
74a7e14dcfSSatish Balay    backtracks until the (nonmonotone) Armijo conditions are satisfied.
75a7e14dcfSSatish Balay 
76a7e14dcfSSatish Balay    Input Parameters:
77441846f8SBarry Smith +  tao - Tao context
78a7e14dcfSSatish Balay .  X - current iterate (on output X contains new iterate, X + step*S)
79a7e14dcfSSatish Balay .  S - search direction
80a7e14dcfSSatish Balay .  f - merit function evaluated at X
81a7e14dcfSSatish Balay .  G - gradient of merit function evaluated at X
82a7e14dcfSSatish Balay .  W - work vector
83a7e14dcfSSatish Balay -  step - initial estimate of step length
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay    Output parameters:
86a7e14dcfSSatish Balay +  f - merit function evaluated at new iterate, X + step*S
87a7e14dcfSSatish Balay .  G - gradient of merit function evaluated at new iterate, X + step*S
88a7e14dcfSSatish Balay .  X - new iterate
89a7e14dcfSSatish Balay -  step - final step length
90a7e14dcfSSatish Balay 
91a7e14dcfSSatish Balay @ */
92a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
93a7e14dcfSSatish Balay {
948caf6e8cSBarry Smith   TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data;
952a0dac07SAlp Dener   PetscInt             i,its=0;
96a7e14dcfSSatish Balay   PetscReal            fact, ref, gdx;
97a7e14dcfSSatish Balay   PetscInt             idx;
98a7e14dcfSSatish Balay   PetscBool            g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
99a7e14dcfSSatish Balay 
100a7e14dcfSSatish Balay   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
1022a0dac07SAlp Dener 
103a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
104a7e14dcfSSatish Balay   if (!armP->work) {
1059566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&armP->work));
106a7e14dcfSSatish Balay     armP->x = x;
1079566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)armP->x));
10847a47007SBarry Smith   } else if (x != armP->x) {
109a7e14dcfSSatish Balay     /* If x has changed, then recreate work */
1109566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&armP->work));
1119566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&armP->work));
1129566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)armP->x));
113a7e14dcfSSatish Balay     armP->x = x;
1149566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)armP->x));
115a7e14dcfSSatish Balay   }
116a7e14dcfSSatish Balay 
117a7e14dcfSSatish Balay   /* Check linesearch parameters */
118a7e14dcfSSatish Balay   if (armP->alpha < 1) {
1199566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha));
120a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
12147a47007SBarry Smith   } else if ((armP->beta <= 0) || (armP->beta >= 1)) {
1229566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: beta (%g) invalid\n", (double)armP->beta));
123a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
12447a47007SBarry Smith   } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) {
1259566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf));
126a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
12747a47007SBarry Smith   } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) {
1289566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma));
129a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
13047a47007SBarry Smith   } else if (armP->memorySize < 1) {
1319566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize));
132a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
13347a47007SBarry Smith   } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) {
1349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: reference_policy invalid\n"));
135a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
13647a47007SBarry Smith   } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) {
1379566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n"));
138a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
13947a47007SBarry Smith   } else if (PetscIsInfOrNanReal(*f)) {
1409566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Armijo line search error: initial function inf or nan\n"));
141a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER;
142a7e14dcfSSatish Balay   }
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay   if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) {
145a7e14dcfSSatish Balay     PetscFunctionReturn(0);
146a7e14dcfSSatish Balay   }
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay   /* Check to see of the memory has been allocated.  If not, allocate
149a7e14dcfSSatish Balay      the historical array and populate it with the initial function
150a7e14dcfSSatish Balay      values. */
1516c23d075SBarry Smith   if (!armP->memory) {
1529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(armP->memorySize, &armP->memory));
153a7e14dcfSSatish Balay   }
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay   if (!armP->memorySetup) {
156a7e14dcfSSatish Balay     for (i = 0; i < armP->memorySize; i++) {
157a7e14dcfSSatish Balay       armP->memory[i] = armP->alpha*(*f);
158a7e14dcfSSatish Balay     }
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay     armP->current = 0;
161a7e14dcfSSatish Balay     armP->lastReference = armP->memory[0];
162a7e14dcfSSatish Balay     armP->memorySetup=PETSC_TRUE;
163a7e14dcfSSatish Balay   }
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay   /* Calculate reference value (MAX) */
166a7e14dcfSSatish Balay   ref = armP->memory[0];
167a7e14dcfSSatish Balay   idx = 0;
168a7e14dcfSSatish Balay 
169a7e14dcfSSatish Balay   for (i = 1; i < armP->memorySize; i++) {
170a7e14dcfSSatish Balay     if (armP->memory[i] > ref) {
171a7e14dcfSSatish Balay       ref = armP->memory[i];
172a7e14dcfSSatish Balay       idx = i;
173a7e14dcfSSatish Balay     }
174a7e14dcfSSatish Balay   }
175a7e14dcfSSatish Balay 
176a7e14dcfSSatish Balay   if (armP->referencePolicy == REFERENCE_AVE) {
177a7e14dcfSSatish Balay     ref = 0;
178a7e14dcfSSatish Balay     for (i = 0; i < armP->memorySize; i++) {
179a7e14dcfSSatish Balay       ref += armP->memory[i];
180a7e14dcfSSatish Balay     }
181a7e14dcfSSatish Balay     ref = ref / armP->memorySize;
182a7e14dcfSSatish Balay     ref = PetscMax(ref, armP->memory[armP->current]);
18347a47007SBarry Smith   } else if (armP->referencePolicy == REFERENCE_MEAN) {
184a7e14dcfSSatish Balay     ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current]));
185a7e14dcfSSatish Balay   }
1869566063dSJacob Faibussowitsch   PetscCall(VecDot(g,s,&gdx));
187a7e14dcfSSatish Balay 
188a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(gdx)) {
1899566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx));
190a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
191a7e14dcfSSatish Balay     PetscFunctionReturn(0);
192a7e14dcfSSatish Balay   }
193a7e14dcfSSatish Balay   if (gdx >= 0.0) {
1949566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx));
195a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
196a7e14dcfSSatish Balay     PetscFunctionReturn(0);
197a7e14dcfSSatish Balay   }
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay   if (armP->nondescending) {
200a7e14dcfSSatish Balay     fact = armP->sigma;
201a7e14dcfSSatish Balay   } else {
202a7e14dcfSSatish Balay     fact = armP->sigma * gdx;
203a7e14dcfSSatish Balay   }
204a7e14dcfSSatish Balay   ls->step = ls->initstep;
205a7e14dcfSSatish Balay   while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) {
206a7e14dcfSSatish Balay     /* Calculate iterate */
2072a0dac07SAlp Dener     ++its;
208*ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(armP->work,ls->step,s,x));
209a7e14dcfSSatish Balay     if (ls->bounded) {
2109566063dSJacob Faibussowitsch       PetscCall(VecMedian(ls->lower,armP->work,ls->upper,armP->work));
211a7e14dcfSSatish Balay     }
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay     /* Calculate function at new iterate */
214a7e14dcfSSatish Balay     if (ls->hasobjective) {
2159566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjective(ls,armP->work,f));
216a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
217a7e14dcfSSatish Balay     } else if (ls->usegts) {
2189566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx));
219a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
220a7e14dcfSSatish Balay     } else {
2219566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g));
222a7e14dcfSSatish Balay       g_computed=PETSC_TRUE;
223a7e14dcfSSatish Balay     }
224a7e14dcfSSatish Balay     if (ls->step == ls->initstep) {
225a7e14dcfSSatish Balay       ls->f_fullstep = *f;
226a7e14dcfSSatish Balay     }
227a7e14dcfSSatish Balay 
2289566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchMonitor(ls, its, *f, ls->step));
2292a0dac07SAlp Dener 
230a7e14dcfSSatish Balay     if (PetscIsInfOrNanReal(*f)) {
231a7e14dcfSSatish Balay       ls->step *= armP->beta_inf;
23247a47007SBarry Smith     } else {
233a7e14dcfSSatish Balay       /* Check descent condition */
234a7e14dcfSSatish Balay       if (armP->nondescending && *f <= ref - ls->step*fact*ref)
235a7e14dcfSSatish Balay         break;
236a7e14dcfSSatish Balay       if (!armP->nondescending && *f <= ref + ls->step*fact) {
237a7e14dcfSSatish Balay         break;
238a7e14dcfSSatish Balay       }
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay       ls->step *= armP->beta;
241a7e14dcfSSatish Balay     }
242a7e14dcfSSatish Balay   }
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay   /* Check termination */
245a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(*f)) {
2469566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Function is inf or nan.\n"));
247a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
248a7e14dcfSSatish Balay   } else if (ls->step < ls->stepmin) {
2499566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Step length is below tolerance.\n"));
250a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_RTOL;
251a7e14dcfSSatish Balay   } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
2529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs));
253a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
254a7e14dcfSSatish Balay   }
255a7e14dcfSSatish Balay   if (ls->reason) {
256a7e14dcfSSatish Balay     PetscFunctionReturn(0);
257a7e14dcfSSatish Balay   }
258a7e14dcfSSatish Balay 
259a7e14dcfSSatish Balay   /* Successful termination, update memory */
260a31dc708SJason Sarich   ls->reason = TAOLINESEARCH_SUCCESS;
261a7e14dcfSSatish Balay   armP->lastReference = ref;
262a7e14dcfSSatish Balay   if (armP->replacementPolicy == REPLACE_FIFO) {
263a7e14dcfSSatish Balay     armP->memory[armP->current++] = *f;
264a7e14dcfSSatish Balay     if (armP->current >= armP->memorySize) {
265a7e14dcfSSatish Balay       armP->current = 0;
266a7e14dcfSSatish Balay     }
267a7e14dcfSSatish Balay   } else {
268a7e14dcfSSatish Balay     armP->current = idx;
269a7e14dcfSSatish Balay     armP->memory[idx] = *f;
270a7e14dcfSSatish Balay   }
271a7e14dcfSSatish Balay 
272a7e14dcfSSatish Balay   /* Update iterate and compute gradient */
2739566063dSJacob Faibussowitsch   PetscCall(VecCopy(armP->work,x));
274a7e14dcfSSatish Balay   if (!g_computed) {
2759566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchComputeGradient(ls, x, g));
276a7e14dcfSSatish Balay   }
2779566063dSJacob Faibussowitsch   PetscCall(PetscInfo(ls, "%D function evals in line search, step = %g\n",ls->nfeval, (double)ls->step));
278a7e14dcfSSatish Balay   PetscFunctionReturn(0);
279a7e14dcfSSatish Balay }
280a7e14dcfSSatish Balay 
28190b6438dSAlp Dener /*MC
28290b6438dSAlp Dener    TAOLINESEARCHARMIJO - Backtracking line-search that satisfies only the (nonmonotone) Armijo condition
28390b6438dSAlp Dener    (i.e., sufficient decrease).
28490b6438dSAlp Dener 
28590b6438dSAlp Dener    Armijo line-search type can be selected with "-tao_ls_type armijo".
28690b6438dSAlp Dener 
28790b6438dSAlp Dener    Level: developer
28890b6438dSAlp Dener 
28990b6438dSAlp Dener .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchApply()
29090b6438dSAlp Dener 
29190b6438dSAlp Dener .keywords: Tao, linesearch
29290b6438dSAlp Dener M*/
293728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls)
294a7e14dcfSSatish Balay {
2958caf6e8cSBarry Smith   TaoLineSearch_ARMIJO *armP;
296a7e14dcfSSatish Balay 
297a7e14dcfSSatish Balay   PetscFunctionBegin;
298a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
2999566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ls,&armP));
300a7e14dcfSSatish Balay 
3016c23d075SBarry Smith   armP->memory = NULL;
302a7e14dcfSSatish Balay   armP->alpha = 1.0;
303a7e14dcfSSatish Balay   armP->beta = 0.5;
304a7e14dcfSSatish Balay   armP->beta_inf = 0.5;
305a7e14dcfSSatish Balay   armP->sigma = 1e-4;
306a7e14dcfSSatish Balay   armP->memorySize = 1;
307a7e14dcfSSatish Balay   armP->referencePolicy = REFERENCE_MAX;
308a7e14dcfSSatish Balay   armP->replacementPolicy = REPLACE_MRU;
309a7e14dcfSSatish Balay   armP->nondescending=PETSC_FALSE;
310a7e14dcfSSatish Balay   ls->data = (void*)armP;
311a7e14dcfSSatish Balay   ls->initstep = 1.0;
31283c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
3134664e3ffSStefano Zampini   ls->ops->monitor = NULL;
314a7e14dcfSSatish Balay   ls->ops->apply = TaoLineSearchApply_Armijo;
315a7e14dcfSSatish Balay   ls->ops->view = TaoLineSearchView_Armijo;
316a7e14dcfSSatish Balay   ls->ops->destroy = TaoLineSearchDestroy_Armijo;
317a7e14dcfSSatish Balay   ls->ops->reset = TaoLineSearchReset_Armijo;
318a7e14dcfSSatish Balay   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo;
319a7e14dcfSSatish Balay   PetscFunctionReturn(0);
320a7e14dcfSSatish Balay }
321