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