xref: /petsc/src/tao/linesearch/interface/taolinesearch.c (revision cedd07cade5cbdfdad435c8172b7ec8972d9cd8d)
1 #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/
2 #include <petsc/private/taolinesearchimpl.h>
3 
4 PetscFunctionList TaoLineSearchList = NULL;
5 
6 PetscClassId TAOLINESEARCH_CLASSID = 0;
7 
8 PetscLogEvent TAOLINESEARCH_Apply;
9 PetscLogEvent TAOLINESEARCH_Eval;
10 
11 /*@C
12    TaoLineSearchViewFromOptions - View from Options
13 
14    Collective
15 
16    Input Parameters:
17 +  A - the Tao context
18 .  obj - Optional object
19 -  name - command line option
20 
21    Level: intermediate
22 .seealso: `TaoLineSearch`, `TaoLineSearchView`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()`
23 @*/
24 PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A, PetscObject obj, const char name[])
25 {
26   PetscFunctionBegin;
27   PetscValidHeaderSpecific(A, TAOLINESEARCH_CLASSID, 1);
28   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
29   PetscFunctionReturn(0);
30 }
31 
32 /*@C
33   TaoLineSearchView - Prints information about the TaoLineSearch
34 
35   Collective
36 
37   InputParameters:
38 + ls - the Tao context
39 - viewer - visualization context
40 
41   Options Database Key:
42 . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search
43 
44   Notes:
45   The available visualization contexts include
46 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
47 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
48          output where only the first processor opens
49          the file.  All other processors send their
50          data to the first processor to print.
51 
52   Level: beginner
53 
54 .seealso: `PetscViewerASCIIOpen()`
55 @*/
56 
57 PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
58 {
59   PetscBool         isascii, isstring;
60   TaoLineSearchType type;
61 
62   PetscFunctionBegin;
63   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
64   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer));
65   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
66   PetscCheckSameComm(ls, 1, viewer, 2);
67 
68   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
69   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
70   if (isascii) {
71     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer));
72     PetscCall(PetscViewerASCIIPushTab(viewer));
73     PetscTryTypeMethod(ls, view, viewer);
74     PetscCall(PetscViewerASCIIPopTab(viewer));
75     PetscCall(PetscViewerASCIIPushTab(viewer));
76     PetscCall(PetscViewerASCIIPrintf(viewer, "maximum function evaluations=%" PetscInt_FMT "\n", ls->max_funcs));
77     PetscCall(PetscViewerASCIIPrintf(viewer, "tolerances: ftol=%g, rtol=%g, gtol=%g\n", (double)ls->ftol, (double)ls->rtol, (double)ls->gtol));
78     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT "\n", ls->nfeval));
79     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT "\n", ls->ngeval));
80     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT "\n", ls->nfgeval));
81 
82     if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(viewer, "using variable bounds\n"));
83     PetscCall(PetscViewerASCIIPrintf(viewer, "Termination reason: %d\n", (int)ls->reason));
84     PetscCall(PetscViewerASCIIPopTab(viewer));
85   } else if (isstring) {
86     PetscCall(TaoLineSearchGetType(ls, &type));
87     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
88   }
89   PetscFunctionReturn(0);
90 }
91 
92 /*@C
93   TaoLineSearchCreate - Creates a TAO Line Search object.  Algorithms in TAO that use
94   line-searches will automatically create one.
95 
96   Collective
97 
98   Input Parameter:
99 . comm - MPI communicator
100 
101   Output Parameter:
102 . newls - the new TaoLineSearch context
103 
104   Available methods include:
105 + more-thuente - the More-Thuente method
106 . gpcg - the GPCG method
107 - unit - Do not perform any line search
108 
109    Options Database Keys:
110 .   -tao_ls_type - select which method TAO should use
111 
112    Level: beginner
113 
114 .seealso: `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()`
115 @*/
116 
117 PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
118 {
119   TaoLineSearch ls;
120 
121   PetscFunctionBegin;
122   PetscValidPointer(newls, 2);
123   PetscCall(TaoLineSearchInitializePackage());
124 
125   PetscCall(PetscHeaderCreate(ls, TAOLINESEARCH_CLASSID, "TaoLineSearch", "Linesearch", "Tao", comm, TaoLineSearchDestroy, TaoLineSearchView));
126   ls->max_funcs = 30;
127   ls->ftol      = 0.0001;
128   ls->gtol      = 0.9;
129 #if defined(PETSC_USE_REAL_SINGLE)
130   ls->rtol = 1.0e-5;
131 #else
132   ls->rtol = 1.0e-10;
133 #endif
134   ls->stepmin  = 1.0e-20;
135   ls->stepmax  = 1.0e+20;
136   ls->step     = 1.0;
137   ls->initstep = 1.0;
138   *newls       = ls;
139   PetscFunctionReturn(0);
140 }
141 
142 /*@
143   TaoLineSearchSetUp - Sets up the internal data structures for the later use
144   of a Tao solver
145 
146   Collective
147 
148   Input Parameters:
149 . ls - the TaoLineSearch context
150 
151   Notes:
152   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
153   automatically be called in TaoLineSearchSolve().  However, if the user
154   desires to call it explicitly, it should come after TaoLineSearchCreate()
155   but before TaoLineSearchApply().
156 
157   Level: developer
158 
159 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()`
160 @*/
161 
162 PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
163 {
164   const char *default_type = TAOLINESEARCHMT;
165   PetscBool   flg;
166 
167   PetscFunctionBegin;
168   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
169   if (ls->setupcalled) PetscFunctionReturn(0);
170   if (!((PetscObject)ls)->type_name) PetscCall(TaoLineSearchSetType(ls, default_type));
171   PetscTryTypeMethod(ls, setup);
172   if (ls->usetaoroutines) {
173     PetscCall(TaoIsObjectiveDefined(ls->tao, &flg));
174     ls->hasobjective = flg;
175     PetscCall(TaoIsGradientDefined(ls->tao, &flg));
176     ls->hasgradient = flg;
177     PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao, &flg));
178     ls->hasobjectiveandgradient = flg;
179   } else {
180     if (ls->ops->computeobjective) {
181       ls->hasobjective = PETSC_TRUE;
182     } else {
183       ls->hasobjective = PETSC_FALSE;
184     }
185     if (ls->ops->computegradient) {
186       ls->hasgradient = PETSC_TRUE;
187     } else {
188       ls->hasgradient = PETSC_FALSE;
189     }
190     if (ls->ops->computeobjectiveandgradient) {
191       ls->hasobjectiveandgradient = PETSC_TRUE;
192     } else {
193       ls->hasobjectiveandgradient = PETSC_FALSE;
194     }
195   }
196   ls->setupcalled = PETSC_TRUE;
197   PetscFunctionReturn(0);
198 }
199 
200 /*@
201   TaoLineSearchReset - Some line searches may carry state information
202   from one TaoLineSearchApply() to the next.  This function resets this
203   state information.
204 
205   Collective
206 
207   Input Parameter:
208 . ls - the TaoLineSearch context
209 
210   Level: developer
211 
212 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()`
213 @*/
214 PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
215 {
216   PetscFunctionBegin;
217   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
218   PetscTryTypeMethod(ls, reset);
219   PetscFunctionReturn(0);
220 }
221 
222 /*@
223   TaoLineSearchDestroy - Destroys the TAO context that was created with
224   TaoLineSearchCreate()
225 
226   Collective
227 
228   Input Parameter:
229 . ls - the TaoLineSearch context
230 
231   Level: beginner
232 
233 .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
234 @*/
235 PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
236 {
237   PetscFunctionBegin;
238   if (!*ls) PetscFunctionReturn(0);
239   PetscValidHeaderSpecific(*ls, TAOLINESEARCH_CLASSID, 1);
240   if (--((PetscObject)*ls)->refct > 0) {
241     *ls = NULL;
242     PetscFunctionReturn(0);
243   }
244   PetscCall(VecDestroy(&(*ls)->stepdirection));
245   PetscCall(VecDestroy(&(*ls)->start_x));
246   PetscCall(VecDestroy(&(*ls)->upper));
247   PetscCall(VecDestroy(&(*ls)->lower));
248   if ((*ls)->ops->destroy) PetscCall((*(*ls)->ops->destroy)(*ls));
249   if ((*ls)->usemonitor) PetscCall(PetscViewerDestroy(&(*ls)->viewer));
250   PetscCall(PetscHeaderDestroy(ls));
251   PetscFunctionReturn(0);
252 }
253 
254 /*@
255   TaoLineSearchApply - Performs a line-search in a given step direction.  Criteria for acceptable step length depends on the line-search algorithm chosen
256 
257   Collective
258 
259   Input Parameters:
260 + ls - the Tao context
261 - s - search direction
262 
263   Input/Output Parameters:
264 
265   Output Parameters:
266 + x - On input the current solution, on output x contains the new solution determined by the line search
267 . f - On input the objective function value at current solution, on output contains the objective function value at new solution
268 . g - On input the gradient evaluated at x, on output contains the gradient at new solution
269 . steplength - scalar multiplier of s used ( x = x0 + steplength * x)
270 - reason - reason why the line-search stopped
271 
272   Notes:
273   reason will be set to one of:
274 
275 + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
276 . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
277 . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
278 . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
279 . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
280 . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
281 . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
282 . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
283 . TAOLINESEARCH_HALTED_OTHER - any other reason
284 - TAOLINESEARCH_SUCCESS - successful line search
285 
286   The algorithm developer must set up the TaoLineSearch with calls to
287   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()
288 
289   You may or may not need to follow this with a call to
290   TaoAddLineSearchCounts(), depending on whether you want these
291   evaluations to count toward the total function/gradient evaluations.
292 
293   Level: beginner
294 
295   .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()`
296  @*/
297 
298 PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
299 {
300   PetscInt low1, low2, low3, high1, high2, high3;
301 
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
304   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
305   PetscValidRealPointer(f, 3);
306   PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
307   PetscValidHeaderSpecific(s, VEC_CLASSID, 5);
308   PetscValidPointer(reason, 7);
309   PetscCheckSameComm(ls, 1, x, 2);
310   PetscCheckSameTypeAndComm(x, 2, g, 4);
311   PetscCheckSameTypeAndComm(x, 2, s, 5);
312   PetscCall(VecGetOwnershipRange(x, &low1, &high1));
313   PetscCall(VecGetOwnershipRange(g, &low2, &high2));
314   PetscCall(VecGetOwnershipRange(s, &low3, &high3));
315   PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible vector local lengths");
316 
317   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
318   PetscCall(PetscObjectReference((PetscObject)s));
319   PetscCall(VecDestroy(&ls->stepdirection));
320   ls->stepdirection = s;
321 
322   PetscCall(TaoLineSearchSetUp(ls));
323   ls->nfeval  = 0;
324   ls->ngeval  = 0;
325   ls->nfgeval = 0;
326   /* Check parameter values */
327   if (ls->ftol < 0.0) {
328     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: ftol (%g) < 0\n", (double)ls->ftol));
329     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
330   }
331   if (ls->rtol < 0.0) {
332     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: rtol (%g) < 0\n", (double)ls->rtol));
333     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
334   }
335   if (ls->gtol < 0.0) {
336     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: gtol (%g) < 0\n", (double)ls->gtol));
337     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
338   }
339   if (ls->stepmin < 0.0) {
340     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) < 0\n", (double)ls->stepmin));
341     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
342   }
343   if (ls->stepmax < ls->stepmin) {
344     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n", (double)ls->stepmin, (double)ls->stepmax));
345     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
346   }
347   if (ls->max_funcs < 0) {
348     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n", ls->max_funcs));
349     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
350   }
351   if (PetscIsInfOrNanReal(*f)) {
352     PetscCall(PetscInfo(ls, "Initial Line Search Function Value is Inf or Nan (%g)\n", (double)*f));
353     *reason = TAOLINESEARCH_FAILED_INFORNAN;
354   }
355 
356   PetscCall(PetscObjectReference((PetscObject)x));
357   PetscCall(VecDestroy(&ls->start_x));
358   ls->start_x = x;
359 
360   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply, ls, 0, 0, 0));
361   PetscUseTypeMethod(ls, apply, x, f, g, s);
362   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0, 0, 0));
363   *reason   = ls->reason;
364   ls->new_f = *f;
365 
366   if (steplength) *steplength = ls->step;
367 
368   PetscCall(TaoLineSearchViewFromOptions(ls, NULL, "-tao_ls_view"));
369   PetscFunctionReturn(0);
370 }
371 
372 /*@C
373    TaoLineSearchSetType - Sets the algorithm used in a line search
374 
375    Collective
376 
377    Input Parameters:
378 +  ls - the TaoLineSearch context
379 -  type - the TaoLineSearchType selection
380 
381   Available methods include:
382 +  more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition
383 .  armijo - simple backtracking line search enforcing only the sufficient decrease condition
384 -  unit - do not perform a line search and always accept unit step length
385 
386   Options Database Keys:
387 .  -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime
388 
389   Level: beginner
390 
391 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchGetType()`, `TaoLineSearchApply()`
392 
393 @*/
394 
395 PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
396 {
397   PetscErrorCode (*r)(TaoLineSearch);
398   PetscBool flg;
399 
400   PetscFunctionBegin;
401   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
402   PetscValidCharPointer(type, 2);
403   PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg));
404   if (flg) PetscFunctionReturn(0);
405 
406   PetscCall(PetscFunctionListFind(TaoLineSearchList, type, (void (**)(void)) & r));
407   PetscCheck(r, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested TaoLineSearch type %s", type);
408   PetscTryTypeMethod(ls, destroy);
409   ls->max_funcs = 30;
410   ls->ftol      = 0.0001;
411   ls->gtol      = 0.9;
412 #if defined(PETSC_USE_REAL_SINGLE)
413   ls->rtol = 1.0e-5;
414 #else
415   ls->rtol = 1.0e-10;
416 #endif
417   ls->stepmin = 1.0e-20;
418   ls->stepmax = 1.0e+20;
419 
420   ls->nfeval              = 0;
421   ls->ngeval              = 0;
422   ls->nfgeval             = 0;
423   ls->ops->setup          = NULL;
424   ls->ops->apply          = NULL;
425   ls->ops->view           = NULL;
426   ls->ops->setfromoptions = NULL;
427   ls->ops->destroy        = NULL;
428   ls->setupcalled         = PETSC_FALSE;
429   PetscCall((*r)(ls));
430   PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type));
431   PetscFunctionReturn(0);
432 }
433 
434 /*@C
435   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
436   iteration number, step length, and function value before calling the implementation
437   specific monitor.
438 
439    Input Parameters:
440 +  ls - the TaoLineSearch context
441 .  its - the current iterate number (>=0)
442 .  f - the current objective function value
443 -  step - the step length
444 
445    Options Database Key:
446 .  -tao_ls_monitor - Use the default monitor, which prints statistics to standard output
447 
448    Level: developer
449 
450 @*/
451 PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
452 {
453   PetscInt tabs;
454 
455   PetscFunctionBegin;
456   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
457   if (ls->usemonitor) {
458     PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs));
459     PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel));
460     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its));
461     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f));
462     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step));
463     if (ls->ops->monitor && its > 0) {
464       PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3));
465       PetscUseTypeMethod(ls, monitor);
466     }
467     PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs));
468   }
469   PetscFunctionReturn(0);
470 }
471 
472 /*@
473   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
474   options.
475 
476   Collective
477 
478   Input Parameter:
479 . ls - the TaoLineSearch context
480 
481   Options Database Keys:
482 + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
483 . -tao_ls_ftol <tol> - tolerance for sufficient decrease
484 . -tao_ls_gtol <tol> - tolerance for curvature condition
485 . -tao_ls_rtol <tol> - relative tolerance for acceptable step
486 . -tao_ls_stepinit <step> - initial steplength allowed
487 . -tao_ls_stepmin <step> - minimum steplength allowed
488 . -tao_ls_stepmax <step> - maximum steplength allowed
489 . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
490 - -tao_ls_view - display line-search results to standard output
491 
492   Level: beginner
493 @*/
494 PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
495 {
496   const char *default_type = TAOLINESEARCHMT;
497   char        type[256], monfilename[PETSC_MAX_PATH_LEN];
498   PetscViewer monviewer;
499   PetscBool   flg;
500 
501   PetscFunctionBegin;
502   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
503   PetscObjectOptionsBegin((PetscObject)ls);
504   if (((PetscObject)ls)->type_name) default_type = ((PetscObject)ls)->type_name;
505   /* Check for type from options */
506   PetscCall(PetscOptionsFList("-tao_ls_type", "Tao Line Search type", "TaoLineSearchSetType", TaoLineSearchList, default_type, type, 256, &flg));
507   if (flg) {
508     PetscCall(TaoLineSearchSetType(ls, type));
509   } else if (!((PetscObject)ls)->type_name) {
510     PetscCall(TaoLineSearchSetType(ls, default_type));
511   }
512 
513   PetscCall(PetscOptionsInt("-tao_ls_max_funcs", "max function evals in line search", "", ls->max_funcs, &ls->max_funcs, NULL));
514   PetscCall(PetscOptionsReal("-tao_ls_ftol", "tol for sufficient decrease", "", ls->ftol, &ls->ftol, NULL));
515   PetscCall(PetscOptionsReal("-tao_ls_gtol", "tol for curvature condition", "", ls->gtol, &ls->gtol, NULL));
516   PetscCall(PetscOptionsReal("-tao_ls_rtol", "relative tol for acceptable step", "", ls->rtol, &ls->rtol, NULL));
517   PetscCall(PetscOptionsReal("-tao_ls_stepmin", "lower bound for step", "", ls->stepmin, &ls->stepmin, NULL));
518   PetscCall(PetscOptionsReal("-tao_ls_stepmax", "upper bound for step", "", ls->stepmax, &ls->stepmax, NULL));
519   PetscCall(PetscOptionsReal("-tao_ls_stepinit", "initial step", "", ls->initstep, &ls->initstep, NULL));
520   PetscCall(PetscOptionsString("-tao_ls_monitor", "enable the basic monitor", "TaoLineSearchSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg));
521   if (flg) {
522     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls), monfilename, &monviewer));
523     ls->viewer     = monviewer;
524     ls->usemonitor = PETSC_TRUE;
525   }
526   PetscTryTypeMethod(ls, setfromoptions, PetscOptionsObject);
527   PetscOptionsEnd();
528   PetscFunctionReturn(0);
529 }
530 
531 /*@C
532   TaoLineSearchGetType - Gets the current line search algorithm
533 
534   Not Collective
535 
536   Input Parameter:
537 . ls - the TaoLineSearch context
538 
539   Output Parameter:
540 . type - the line search algorithm in effect
541 
542   Level: developer
543 
544 @*/
545 PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
546 {
547   PetscFunctionBegin;
548   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
549   PetscValidPointer(type, 2);
550   *type = ((PetscObject)ls)->type_name;
551   PetscFunctionReturn(0);
552 }
553 
554 /*@
555   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
556   routines used by the line search in last application (not cumulative).
557 
558   Not Collective
559 
560   Input Parameter:
561 . ls - the TaoLineSearch context
562 
563   Output Parameters:
564 + nfeval   - number of function evaluations
565 . ngeval   - number of gradient evaluations
566 - nfgeval  - number of function/gradient evaluations
567 
568   Level: intermediate
569 
570   Note:
571   If the line search is using the Tao objective and gradient
572   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
573   is already counting the number of evaluations.
574 
575 @*/
576 PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
577 {
578   PetscFunctionBegin;
579   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
580   *nfeval  = ls->nfeval;
581   *ngeval  = ls->ngeval;
582   *nfgeval = ls->nfgeval;
583   PetscFunctionReturn(0);
584 }
585 
586 /*@
587   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
588   Tao evaluation routines.
589 
590   Not Collective
591 
592   Input Parameter:
593 . ls - the TaoLineSearch context
594 
595   Output Parameter:
596 . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
597         otherwise PETSC_FALSE
598 
599   Level: developer
600 @*/
601 PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
602 {
603   PetscFunctionBegin;
604   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
605   *flg = ls->usetaoroutines;
606   PetscFunctionReturn(0);
607 }
608 
609 /*@C
610   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
611 
612   Logically Collective
613 
614   Input Parameters:
615 + ls - the TaoLineSearch context
616 . func - the objective function evaluation routine
617 - ctx - the (optional) user-defined context for private data
618 
619   Calling sequence of func:
620 $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);
621 
622 + x - input vector
623 . f - function value
624 - ctx (optional) user-defined context
625 
626   Level: beginner
627 
628   Note:
629   Use this routine only if you want the line search objective
630   evaluation routine to be different from the Tao's objective
631   evaluation routine. If you use this routine you must also set
632   the line search gradient and/or function/gradient routine.
633 
634   Note:
635   Some algorithms (lcl, gpcg) set their own objective routine for the
636   line search, application programmers should be wary of overriding the
637   default objective routine.
638 
639 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
640 @*/
641 PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *, void *), void *ctx)
642 {
643   PetscFunctionBegin;
644   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
645 
646   ls->ops->computeobjective = func;
647   if (ctx) ls->userctx_func = ctx;
648   ls->usetaoroutines = PETSC_FALSE;
649   PetscFunctionReturn(0);
650 }
651 
652 /*@C
653   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
654 
655   Logically Collective
656 
657   Input Parameters:
658 + ls - the TaoLineSearch context
659 . func - the gradient evaluation routine
660 - ctx - the (optional) user-defined context for private data
661 
662   Calling sequence of func:
663 $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);
664 
665 + x - input vector
666 . g - gradient vector
667 - ctx (optional) user-defined context
668 
669   Level: beginner
670 
671   Note:
672   Use this routine only if you want the line search gradient
673   evaluation routine to be different from the Tao's gradient
674   evaluation routine. If you use this routine you must also set
675   the line search function and/or function/gradient routine.
676 
677   Note:
678   Some algorithms (lcl, gpcg) set their own gradient routine for the
679   line search, application programmers should be wary of overriding the
680   default gradient routine.
681 
682 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
683 @*/
684 PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec g, void *), void *ctx)
685 {
686   PetscFunctionBegin;
687   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
688   ls->ops->computegradient = func;
689   if (ctx) ls->userctx_grad = ctx;
690   ls->usetaoroutines = PETSC_FALSE;
691   PetscFunctionReturn(0);
692 }
693 
694 /*@C
695   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
696 
697   Logically Collective
698 
699   Input Parameters:
700 + ls - the TaoLineSearch context
701 . func - the objective and gradient evaluation routine
702 - ctx - the (optional) user-defined context for private data
703 
704   Calling sequence of func:
705 $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);
706 
707 + x - input vector
708 . f - function value
709 . g - gradient vector
710 - ctx (optional) user-defined context
711 
712   Level: beginner
713 
714   Note:
715   Use this routine only if you want the line search objective and gradient
716   evaluation routines to be different from the Tao's objective
717   and gradient evaluation routines.
718 
719   Note:
720   Some algorithms (lcl, gpcg) set their own objective routine for the
721   line search, application programmers should be wary of overriding the
722   default objective routine.
723 
724 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
725 @*/
726 PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void *), void *ctx)
727 {
728   PetscFunctionBegin;
729   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
730   ls->ops->computeobjectiveandgradient = func;
731   if (ctx) ls->userctx_funcgrad = ctx;
732   ls->usetaoroutines = PETSC_FALSE;
733   PetscFunctionReturn(0);
734 }
735 
736 /*@C
737   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
738   (gradient'*stepdirection) evaluation routine for the line search.
739   Sometimes it is more efficient to compute the inner product of the gradient
740   and the step direction than it is to compute the gradient, and this is all
741   the line search typically needs of the gradient.
742 
743   Logically Collective
744 
745   Input Parameters:
746 + ls - the TaoLineSearch context
747 . func - the objective and gradient evaluation routine
748 - ctx - the (optional) user-defined context for private data
749 
750   Calling sequence of func:
751 $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);
752 
753 + x - input vector
754 . s - step direction
755 . f - function value
756 . gts - inner product of gradient and step direction vectors
757 - ctx (optional) user-defined context
758 
759   Note: The gradient will still need to be computed at the end of the line
760   search, so you will still need to set a line search gradient evaluation
761   routine
762 
763   Note: Bounded line searches (those used in bounded optimization algorithms)
764   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
765   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()
766 
767   Level: advanced
768 
769   Note:
770   Some algorithms (lcl, gpcg) set their own objective routine for the
771   line search, application programmers should be wary of overriding the
772   default objective routine.
773 
774 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()`
775 @*/
776 PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void *), void *ctx)
777 {
778   PetscFunctionBegin;
779   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
780   ls->ops->computeobjectiveandgts = func;
781   if (ctx) ls->userctx_funcgts = ctx;
782   ls->usegts         = PETSC_TRUE;
783   ls->usetaoroutines = PETSC_FALSE;
784   PetscFunctionReturn(0);
785 }
786 
787 /*@C
788   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
789   objective and gradient evaluation routines from the given Tao object.
790 
791   Logically Collective
792 
793   Input Parameters:
794 + ls - the TaoLineSearch context
795 - ts - the Tao context with defined objective/gradient evaluation routines
796 
797   Level: developer
798 
799 .seealso: `TaoLineSearchCreate()`
800 @*/
801 PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
802 {
803   PetscFunctionBegin;
804   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
805   PetscValidHeaderSpecific(ts, TAO_CLASSID, 2);
806   ls->tao            = ts;
807   ls->usetaoroutines = PETSC_TRUE;
808   PetscFunctionReturn(0);
809 }
810 
811 /*@
812   TaoLineSearchComputeObjective - Computes the objective function value at a given point
813 
814   Collective
815 
816   Input Parameters:
817 + ls - the TaoLineSearch context
818 - x - input vector
819 
820   Output Parameter:
821 . f - Objective value at X
822 
823   Notes:
824     TaoLineSearchComputeObjective() is typically used within line searches
825   so most users would not generally call this routine themselves.
826 
827   Level: developer
828 
829 .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
830 @*/
831 PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
832 {
833   Vec       gdummy;
834   PetscReal gts;
835 
836   PetscFunctionBegin;
837   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
838   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
839   PetscValidRealPointer(f, 3);
840   PetscCheckSameComm(ls, 1, x, 2);
841   if (ls->usetaoroutines) {
842     PetscCall(TaoComputeObjective(ls->tao, x, f));
843   } else {
844     PetscCheck(ls->ops->computeobjective || ls->ops->computeobjectiveandgradient || ls->ops->computeobjectiveandgts, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_WRONGSTATE, "Line Search does not have objective function set");
845     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
846     if (ls->ops->computeobjective) PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
847     else if (ls->ops->computeobjectiveandgradient) {
848       PetscCall(VecDuplicate(x, &gdummy));
849       PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgradient)(ls, x, f, gdummy, ls->userctx_funcgrad));
850       PetscCall(VecDestroy(&gdummy));
851     } else PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, &gts, ls->userctx_funcgts));
852     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
853   }
854   ls->nfeval++;
855   PetscFunctionReturn(0);
856 }
857 
858 /*@
859   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
860 
861   Collective
862 
863   Input Parameters:
864 + ls - the TaoLineSearch context
865 - x - input vector
866 
867   Output Parameters:
868 + f - Objective value at X
869 - g - Gradient vector at X
870 
871   Notes:
872     TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
873   so most users would not generally call this routine themselves.
874 
875   Level: developer
876 
877 .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
878 @*/
879 PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
880 {
881   PetscFunctionBegin;
882   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
883   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
884   PetscValidRealPointer(f, 3);
885   PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
886   PetscCheckSameComm(ls, 1, x, 2);
887   PetscCheckSameComm(ls, 1, g, 4);
888   if (ls->usetaoroutines) {
889     PetscCall(TaoComputeObjectiveAndGradient(ls->tao, x, f, g));
890   } else {
891     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
892     if (ls->ops->computeobjectiveandgradient) PetscCallBack("TaoLineSearch callback objective/gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, f, g, ls->userctx_funcgrad));
893     else {
894       PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
895       PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
896     }
897     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
898     PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
899   }
900   ls->nfgeval++;
901   PetscFunctionReturn(0);
902 }
903 
904 /*@
905   TaoLineSearchComputeGradient - Computes the gradient of the objective function
906 
907   Collective
908 
909   Input Parameters:
910 + ls - the TaoLineSearch context
911 - x - input vector
912 
913   Output Parameter:
914 . g - gradient vector
915 
916   Notes:
917     TaoComputeGradient() is typically used within line searches
918   so most users would not generally call this routine themselves.
919 
920   Level: developer
921 
922 .seealso: `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()`
923 @*/
924 PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
925 {
926   PetscReal fdummy;
927 
928   PetscFunctionBegin;
929   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
930   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
931   PetscValidHeaderSpecific(g, VEC_CLASSID, 3);
932   PetscCheckSameComm(ls, 1, x, 2);
933   PetscCheckSameComm(ls, 1, g, 3);
934   if (ls->usetaoroutines) {
935     PetscCall(TaoComputeGradient(ls->tao, x, g));
936   } else {
937     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
938     if (ls->ops->computegradient) PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
939     else PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, &fdummy, g, ls->userctx_funcgrad));
940     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
941   }
942   ls->ngeval++;
943   PetscFunctionReturn(0);
944 }
945 
946 /*@
947   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point
948 
949   Collective
950 
951   Input Parameters:
952 + ls - the TaoLineSearch context
953 - x - input vector
954 
955   Output Parameters:
956 + f - Objective value at X
957 - gts - inner product of gradient and step direction at X
958 
959   Notes:
960     TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
961   so most users would not generally call this routine themselves.
962 
963   Level: developer
964 
965 .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
966 @*/
967 PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
968 {
969   PetscFunctionBegin;
970   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
971   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
972   PetscValidRealPointer(f, 3);
973   PetscValidRealPointer(gts, 4);
974   PetscCheckSameComm(ls, 1, x, 2);
975   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
976   PetscCallBack("TaoLineSearch callback objective/gts", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, gts, ls->userctx_funcgts));
977   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
978   PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
979   ls->nfeval++;
980   PetscFunctionReturn(0);
981 }
982 
983 /*@
984   TaoLineSearchGetSolution - Returns the solution to the line search
985 
986   Collective
987 
988   Input Parameter:
989 . ls - the TaoLineSearch context
990 
991   Output Parameters:
992 + x - the new solution
993 . f - the objective function value at x
994 . g - the gradient at x
995 . steplength - the multiple of the step direction taken by the line search
996 - reason - the reason why the line search terminated
997 
998   reason will be set to one of:
999 
1000 + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
1001 . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
1002 . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
1003 . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
1004 . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
1005 . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
1006 . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
1007 
1008 . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
1009 . TAOLINESEARCH_HALTED_OTHER - any other reason
1010 
1011 - TAOLINESEARCH_SUCCESS - successful line search
1012 
1013   Level: developer
1014 
1015 @*/
1016 PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1017 {
1018   PetscFunctionBegin;
1019   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1020   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
1021   PetscValidRealPointer(f, 3);
1022   PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
1023   PetscValidIntPointer(reason, 6);
1024   if (ls->new_x) PetscCall(VecCopy(ls->new_x, x));
1025   *f = ls->new_f;
1026   if (ls->new_g) PetscCall(VecCopy(ls->new_g, g));
1027   if (steplength) *steplength = ls->step;
1028   *reason = ls->reason;
1029   PetscFunctionReturn(0);
1030 }
1031 
1032 /*@
1033   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1034   search.
1035 
1036   Not Collective
1037 
1038   Input Parameter:
1039 . ls - the TaoLineSearch context
1040 
1041   Output Parameter:
1042 . x - The initial point of the line search
1043 
1044   Level: intermediate
1045 @*/
1046 PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1047 {
1048   PetscFunctionBegin;
1049   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1050   if (x) *x = ls->start_x;
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 /*@
1055   TaoLineSearchGetStepDirection - Gets the step direction of the line
1056   search.
1057 
1058   Not Collective
1059 
1060   Input Parameter:
1061 . ls - the TaoLineSearch context
1062 
1063   Output Parameter:
1064 . s - the step direction of the line search
1065 
1066   Level: advanced
1067 @*/
1068 PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1069 {
1070   PetscFunctionBegin;
1071   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1072   if (s) *s = ls->stepdirection;
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 /*@
1077   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.
1078 
1079   Not Collective
1080 
1081   Input Parameter:
1082 . ls - the TaoLineSearch context
1083 
1084   Output Parameter:
1085 . f - the objective value at the full step length
1086 
1087   Level: developer
1088 @*/
1089 
1090 PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1091 {
1092   PetscFunctionBegin;
1093   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1094   *f_fullstep = ls->f_fullstep;
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 /*@
1099   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.
1100 
1101   Logically Collective
1102 
1103   Input Parameters:
1104 + ls - the TaoLineSearch context
1105 . xl  - vector of lower bounds
1106 - xu  - vector of upper bounds
1107 
1108   Note: If the variable bounds are not set with this routine, then
1109   PETSC_NINFINITY and PETSC_INFINITY are assumed
1110 
1111   Level: beginner
1112 
1113 .seealso: `TaoSetVariableBounds()`, `TaoLineSearchCreate()`
1114 @*/
1115 PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls, Vec xl, Vec xu)
1116 {
1117   PetscFunctionBegin;
1118   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1119   if (xl) PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
1120   if (xu) PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
1121   PetscCall(PetscObjectReference((PetscObject)xl));
1122   PetscCall(PetscObjectReference((PetscObject)xu));
1123   PetscCall(VecDestroy(&ls->lower));
1124   PetscCall(VecDestroy(&ls->upper));
1125   ls->lower   = xl;
1126   ls->upper   = xu;
1127   ls->bounded = (PetscBool)(xl || xu);
1128   PetscFunctionReturn(0);
1129 }
1130 
1131 /*@
1132   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
1133   search.  If this value is not set then 1.0 is assumed.
1134 
1135   Logically Collective
1136 
1137   Input Parameters:
1138 + ls - the TaoLineSearch context
1139 - s - the initial step size
1140 
1141   Level: intermediate
1142 
1143 .seealso: `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()`
1144 @*/
1145 PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls, PetscReal s)
1146 {
1147   PetscFunctionBegin;
1148   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1149   PetscValidLogicalCollectiveReal(ls, s, 2);
1150   ls->initstep = s;
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 /*@
1155   TaoLineSearchGetStepLength - Get the current step length
1156 
1157   Not Collective
1158 
1159   Input Parameters:
1160 . ls - the TaoLineSearch context
1161 
1162   Output Parameters:
1163 . s - the current step length
1164 
1165   Level: beginner
1166 
1167 .seealso: `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()`
1168 @*/
1169 PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls, PetscReal *s)
1170 {
1171   PetscFunctionBegin;
1172   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1173   *s = ls->step;
1174   PetscFunctionReturn(0);
1175 }
1176 
1177 /*@C
1178    TaoLineSearchRegister - Adds a line-search algorithm to the registry
1179 
1180    Not collective
1181 
1182    Input Parameters:
1183 +  sname - name of a new user-defined solver
1184 -  func - routine to Create method context
1185 
1186    Notes:
1187    TaoLineSearchRegister() may be called multiple times to add several user-defined solvers.
1188 
1189    Sample usage:
1190 .vb
1191    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1192 .ve
1193 
1194    Then, your solver can be chosen with the procedural interface via
1195 $     TaoLineSearchSetType(ls,"my_linesearch")
1196    or at runtime via the option
1197 $     -tao_ls_type my_linesearch
1198 
1199    Level: developer
1200 
1201 @*/
1202 PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1203 {
1204   PetscFunctionBegin;
1205   PetscCall(TaoLineSearchInitializePackage());
1206   PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func));
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 /*@C
1211    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1212    for all TaoLineSearch options in the database.
1213 
1214    Collective
1215 
1216    Input Parameters:
1217 +  ls - the TaoLineSearch solver context
1218 -  prefix - the prefix string to prepend to all line search requests
1219 
1220    Notes:
1221    A hyphen (-) must NOT be given at the beginning of the prefix name.
1222    The first character of all runtime options is AUTOMATICALLY the hyphen.
1223 
1224    Level: advanced
1225 
1226 .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
1227 @*/
1228 PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1229 {
1230   return PetscObjectAppendOptionsPrefix((PetscObject)ls, p);
1231 }
1232 
1233 /*@C
1234   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1235   TaoLineSearch options in the database
1236 
1237   Not Collective
1238 
1239   Input Parameters:
1240 . ls - the TaoLineSearch context
1241 
1242   Output Parameters:
1243 . prefix - pointer to the prefix string used is returned
1244 
1245   Notes:
1246     On the fortran side, the user should pass in a string 'prefix' of
1247   sufficient length to hold the prefix.
1248 
1249   Level: advanced
1250 
1251 .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()`
1252 @*/
1253 PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1254 {
1255   return PetscObjectGetOptionsPrefix((PetscObject)ls, p);
1256 }
1257 
1258 /*@C
1259    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1260    TaoLineSearch options in the database.
1261 
1262    Logically Collective
1263 
1264    Input Parameters:
1265 +  ls - the TaoLineSearch context
1266 -  prefix - the prefix string to prepend to all TAO option requests
1267 
1268    Notes:
1269    A hyphen (-) must NOT be given at the beginning of the prefix name.
1270    The first character of all runtime options is AUTOMATICALLY the hyphen.
1271 
1272    For example, to distinguish between the runtime options for two
1273    different line searches, one could call
1274 .vb
1275       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1276       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1277 .ve
1278 
1279    This would enable use of different options for each system, such as
1280 .vb
1281       -sys1_tao_ls_type mt
1282       -sys2_tao_ls_type armijo
1283 .ve
1284 
1285    Level: advanced
1286 
1287 .seealso: `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
1288 @*/
1289 
1290 PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1291 {
1292   return PetscObjectSetOptionsPrefix((PetscObject)ls, p);
1293 }
1294