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