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