xref: /petsc/src/tao/linesearch/interface/taolinesearch.c (revision 5b9dfbb6d981a8e1bcd44e9e88084b3910bf6321)
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   const char     *default_type=TAOLINESEARCHMT;
532   char           type[256],monfilename[PETSC_MAX_PATH_LEN];
533   PetscViewer    monviewer;
534   PetscBool      flg;
535 
536   PetscFunctionBegin;
537   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
538   PetscObjectOptionsBegin((PetscObject)ls);
539   if (((PetscObject)ls)->type_name) {
540     default_type = ((PetscObject)ls)->type_name;
541   }
542   /* Check for type from options */
543   PetscCall(PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg));
544   if (flg) {
545     PetscCall(TaoLineSearchSetType(ls,type));
546   } else if (!((PetscObject)ls)->type_name) {
547     PetscCall(TaoLineSearchSetType(ls,default_type));
548   }
549 
550   PetscCall(PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL));
551   PetscCall(PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL));
552   PetscCall(PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL));
553   PetscCall(PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL));
554   PetscCall(PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL));
555   PetscCall(PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL));
556   PetscCall(PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
557   if (flg) {
558     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer));
559     ls->viewer = monviewer;
560     ls->usemonitor = PETSC_TRUE;
561   }
562   if (ls->ops->setfromoptions) {
563     PetscCall((*ls->ops->setfromoptions)(PetscOptionsObject,ls));
564   }
565   PetscOptionsEnd();
566   PetscFunctionReturn(0);
567 }
568 
569 /*@C
570   TaoLineSearchGetType - Gets the current line search algorithm
571 
572   Not Collective
573 
574   Input Parameter:
575 . ls - the TaoLineSearch context
576 
577   Output Parameter:
578 . type - the line search algorithm in effect
579 
580   Level: developer
581 
582 @*/
583 PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
584 {
585   PetscFunctionBegin;
586   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
587   PetscValidPointer(type,2);
588   *type = ((PetscObject)ls)->type_name;
589   PetscFunctionReturn(0);
590 }
591 
592 /*@
593   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
594   routines used by the line search in last application (not cumulative).
595 
596   Not Collective
597 
598   Input Parameter:
599 . ls - the TaoLineSearch context
600 
601   Output Parameters:
602 + nfeval   - number of function evaluations
603 . ngeval   - number of gradient evaluations
604 - nfgeval  - number of function/gradient evaluations
605 
606   Level: intermediate
607 
608   Note:
609   If the line search is using the Tao objective and gradient
610   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
611   is already counting the number of evaluations.
612 
613 @*/
614 PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
615 {
616   PetscFunctionBegin;
617   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
618   *nfeval = ls->nfeval;
619   *ngeval = ls->ngeval;
620   *nfgeval = ls->nfgeval;
621   PetscFunctionReturn(0);
622 }
623 
624 /*@
625   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
626   Tao evaluation routines.
627 
628   Not Collective
629 
630   Input Parameter:
631 . ls - the TaoLineSearch context
632 
633   Output Parameter:
634 . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
635         otherwise PETSC_FALSE
636 
637   Level: developer
638 @*/
639 PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
640 {
641   PetscFunctionBegin;
642   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
643   *flg = ls->usetaoroutines;
644   PetscFunctionReturn(0);
645 }
646 
647 /*@C
648   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
649 
650   Logically Collective on TaoLineSearch
651 
652   Input Parameters:
653 + ls - the TaoLineSearch context
654 . func - the objective function evaluation routine
655 - ctx - the (optional) user-defined context for private data
656 
657   Calling sequence of func:
658 $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);
659 
660 + x - input vector
661 . f - function value
662 - ctx (optional) user-defined context
663 
664   Level: beginner
665 
666   Note:
667   Use this routine only if you want the line search objective
668   evaluation routine to be different from the Tao's objective
669   evaluation routine. If you use this routine you must also set
670   the line search gradient and/or function/gradient routine.
671 
672   Note:
673   Some algorithms (lcl, gpcg) set their own objective routine for the
674   line search, application programmers should be wary of overriding the
675   default objective routine.
676 
677 .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
678 @*/
679 PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
680 {
681   PetscFunctionBegin;
682   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
683 
684   ls->ops->computeobjective=func;
685   if (ctx) ls->userctx_func=ctx;
686   ls->usetaoroutines=PETSC_FALSE;
687   PetscFunctionReturn(0);
688 }
689 
690 /*@C
691   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
692 
693   Logically Collective on TaoLineSearch
694 
695   Input Parameters:
696 + ls - the TaoLineSearch context
697 . func - the gradient evaluation routine
698 - ctx - the (optional) user-defined context for private data
699 
700   Calling sequence of func:
701 $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);
702 
703 + x - input vector
704 . g - gradient vector
705 - ctx (optional) user-defined context
706 
707   Level: beginner
708 
709   Note:
710   Use this routine only if you want the line search gradient
711   evaluation routine to be different from the Tao's gradient
712   evaluation routine. If you use this routine you must also set
713   the line search function and/or function/gradient routine.
714 
715   Note:
716   Some algorithms (lcl, gpcg) set their own gradient routine for the
717   line search, application programmers should be wary of overriding the
718   default gradient routine.
719 
720 .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
721 @*/
722 PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
723 {
724   PetscFunctionBegin;
725   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
726   ls->ops->computegradient=func;
727   if (ctx) ls->userctx_grad=ctx;
728   ls->usetaoroutines=PETSC_FALSE;
729   PetscFunctionReturn(0);
730 }
731 
732 /*@C
733   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
734 
735   Logically Collective on TaoLineSearch
736 
737   Input Parameters:
738 + ls - the TaoLineSearch context
739 . func - the objective and gradient evaluation routine
740 - ctx - the (optional) user-defined context for private data
741 
742   Calling sequence of func:
743 $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);
744 
745 + x - input vector
746 . f - function value
747 . g - gradient vector
748 - ctx (optional) user-defined context
749 
750   Level: beginner
751 
752   Note:
753   Use this routine only if you want the line search objective and gradient
754   evaluation routines to be different from the Tao's objective
755   and gradient evaluation routines.
756 
757   Note:
758   Some algorithms (lcl, gpcg) set their own objective routine for the
759   line search, application programmers should be wary of overriding the
760   default objective routine.
761 
762 .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines()
763 @*/
764 PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
765 {
766   PetscFunctionBegin;
767   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
768   ls->ops->computeobjectiveandgradient=func;
769   if (ctx) ls->userctx_funcgrad=ctx;
770   ls->usetaoroutines = PETSC_FALSE;
771   PetscFunctionReturn(0);
772 }
773 
774 /*@C
775   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
776   (gradient'*stepdirection) evaluation routine for the line search.
777   Sometimes it is more efficient to compute the inner product of the gradient
778   and the step direction than it is to compute the gradient, and this is all
779   the line search typically needs of the gradient.
780 
781   Logically Collective on TaoLineSearch
782 
783   Input Parameters:
784 + ls - the TaoLineSearch context
785 . func - the objective and gradient evaluation routine
786 - ctx - the (optional) user-defined context for private data
787 
788   Calling sequence of func:
789 $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);
790 
791 + x - input vector
792 . s - step direction
793 . f - function value
794 . gts - inner product of gradient and step direction vectors
795 - ctx (optional) user-defined context
796 
797   Note: The gradient will still need to be computed at the end of the line
798   search, so you will still need to set a line search gradient evaluation
799   routine
800 
801   Note: Bounded line searches (those used in bounded optimization algorithms)
802   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
803   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()
804 
805   Level: advanced
806 
807   Note:
808   Some algorithms (lcl, gpcg) set their own objective routine for the
809   line search, application programmers should be wary of overriding the
810   default objective routine.
811 
812 .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
813 @*/
814 PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
815 {
816   PetscFunctionBegin;
817   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
818   ls->ops->computeobjectiveandgts=func;
819   if (ctx) ls->userctx_funcgts=ctx;
820   ls->usegts = PETSC_TRUE;
821   ls->usetaoroutines = PETSC_FALSE;
822   PetscFunctionReturn(0);
823 }
824 
825 /*@C
826   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
827   objective and gradient evaluation routines from the given Tao object.
828 
829   Logically Collective on TaoLineSearch
830 
831   Input Parameters:
832 + ls - the TaoLineSearch context
833 - ts - the Tao context with defined objective/gradient evaluation routines
834 
835   Level: developer
836 
837 .seealso: TaoLineSearchCreate()
838 @*/
839 PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
840 {
841   PetscFunctionBegin;
842   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
843   PetscValidHeaderSpecific(ts,TAO_CLASSID,2);
844   ls->tao = ts;
845   ls->usetaoroutines = PETSC_TRUE;
846   PetscFunctionReturn(0);
847 }
848 
849 /*@
850   TaoLineSearchComputeObjective - Computes the objective function value at a given point
851 
852   Collective on TaoLineSearch
853 
854   Input Parameters:
855 + ls - the TaoLineSearch context
856 - x - input vector
857 
858   Output Parameter:
859 . f - Objective value at X
860 
861   Notes:
862     TaoLineSearchComputeObjective() is typically used within line searches
863   so most users would not generally call this routine themselves.
864 
865   Level: developer
866 
867 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
868 @*/
869 PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
870 {
871   Vec            gdummy;
872   PetscReal      gts;
873 
874   PetscFunctionBegin;
875   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
876   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
877   PetscValidRealPointer(f,3);
878   PetscCheckSameComm(ls,1,x,2);
879   if (ls->usetaoroutines) {
880     PetscCall(TaoComputeObjective(ls->tao,x,f));
881   } else {
882     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");
883     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
884     PetscStackPush("TaoLineSearch user objective routine");
885     if (ls->ops->computeobjective) {
886       PetscCall((*ls->ops->computeobjective)(ls,x,f,ls->userctx_func));
887     } else if (ls->ops->computeobjectiveandgradient) {
888       PetscCall(VecDuplicate(x,&gdummy));
889       PetscCall((*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad));
890       PetscCall(VecDestroy(&gdummy));
891     } else {
892       PetscCall((*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts));
893     }
894     PetscStackPop;
895     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
896   }
897   ls->nfeval++;
898   PetscFunctionReturn(0);
899 }
900 
901 /*@
902   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
903 
904   Collective on Tao
905 
906   Input Parameters:
907 + ls - the TaoLineSearch context
908 - x - input vector
909 
910   Output Parameters:
911 + f - Objective value at X
912 - g - Gradient vector at X
913 
914   Notes:
915     TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
916   so most users would not generally call this routine themselves.
917 
918   Level: developer
919 
920 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
921 @*/
922 PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
923 {
924   PetscFunctionBegin;
925   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
926   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
927   PetscValidRealPointer(f,3);
928   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
929   PetscCheckSameComm(ls,1,x,2);
930   PetscCheckSameComm(ls,1,g,4);
931   if (ls->usetaoroutines) {
932     PetscCall(TaoComputeObjectiveAndGradient(ls->tao,x,f,g));
933   } else {
934     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
935     if (ls->ops->computeobjectiveandgradient) {
936       PetscStackPush("TaoLineSearch user objective/gradient routine");
937       PetscCall((*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad));
938       PetscStackPop;
939     } else {
940       PetscCheck(ls->ops->computeobjective,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
941       PetscStackPush("TaoLineSearch user objective routine");
942       PetscCall((*ls->ops->computeobjective)(ls,x,f,ls->userctx_func));
943       PetscStackPop;
944       PetscCheck(ls->ops->computegradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");
945       PetscStackPush("TaoLineSearch user gradient routine");
946       PetscCall((*ls->ops->computegradient)(ls,x,g,ls->userctx_grad));
947       PetscStackPop;
948     }
949     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
950     PetscCall(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f)));
951   }
952   ls->nfgeval++;
953   PetscFunctionReturn(0);
954 }
955 
956 /*@
957   TaoLineSearchComputeGradient - Computes the gradient of the objective function
958 
959   Collective on TaoLineSearch
960 
961   Input Parameters:
962 + ls - the TaoLineSearch context
963 - x - input vector
964 
965   Output Parameter:
966 . g - gradient vector
967 
968   Notes:
969     TaoComputeGradient() is typically used within line searches
970   so most users would not generally call this routine themselves.
971 
972   Level: developer
973 
974 .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
975 @*/
976 PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
977 {
978   PetscReal      fdummy;
979 
980   PetscFunctionBegin;
981   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
982   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
983   PetscValidHeaderSpecific(g,VEC_CLASSID,3);
984   PetscCheckSameComm(ls,1,x,2);
985   PetscCheckSameComm(ls,1,g,3);
986   if (ls->usetaoroutines) {
987     PetscCall(TaoComputeGradient(ls->tao,x,g));
988   } else {
989     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
990     PetscStackPush("TaoLineSearch user gradient routine");
991     if (ls->ops->computegradient) {
992       PetscCall((*ls->ops->computegradient)(ls,x,g,ls->userctx_grad));
993     } else {
994       PetscCheck(ls->ops->computeobjectiveandgradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
995       PetscCall((*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad));
996     }
997     PetscStackPop;
998     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
999   }
1000   ls->ngeval++;
1001   PetscFunctionReturn(0);
1002 }
1003 
1004 /*@
1005   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point
1006 
1007   Collective on Tao
1008 
1009   Input Parameters:
1010 + ls - the TaoLineSearch context
1011 - x - input vector
1012 
1013   Output Parameters:
1014 + f - Objective value at X
1015 - gts - inner product of gradient and step direction at X
1016 
1017   Notes:
1018     TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
1019   so most users would not generally call this routine themselves.
1020 
1021   Level: developer
1022 
1023 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
1024 @*/
1025 PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
1026 {
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1029   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1030   PetscValidRealPointer(f,3);
1031   PetscValidRealPointer(gts,4);
1032   PetscCheckSameComm(ls,1,x,2);
1033   PetscCheck(ls->ops->computeobjectiveandgts,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
1034   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
1035   PetscStackPush("TaoLineSearch user objective/gts routine");
1036   PetscCall((*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts));
1037   PetscStackPop;
1038   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
1039   PetscCall(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f)));
1040   ls->nfeval++;
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 /*@
1045   TaoLineSearchGetSolution - Returns the solution to the line search
1046 
1047   Collective on TaoLineSearch
1048 
1049   Input Parameter:
1050 . ls - the TaoLineSearch context
1051 
1052   Output Parameters:
1053 + x - the new solution
1054 . f - the objective function value at x
1055 . g - the gradient at x
1056 . steplength - the multiple of the step direction taken by the line search
1057 - reason - the reason why the line search terminated
1058 
1059   reason will be set to one of:
1060 
1061 + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
1062 . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
1063 . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
1064 . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
1065 . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
1066 . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
1067 . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
1068 
1069 . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
1070 . TAOLINESEARCH_HALTED_OTHER - any other reason
1071 
1072 - TAOLINESEARCH_SUCCESS - successful line search
1073 
1074   Level: developer
1075 
1076 @*/
1077 PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1078 {
1079   PetscFunctionBegin;
1080   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1081   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1082   PetscValidRealPointer(f,3);
1083   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
1084   PetscValidIntPointer(reason,6);
1085   if (ls->new_x) {
1086     PetscCall(VecCopy(ls->new_x,x));
1087   }
1088   *f = ls->new_f;
1089   if (ls->new_g) {
1090     PetscCall(VecCopy(ls->new_g,g));
1091   }
1092   if (steplength) *steplength = ls->step;
1093   *reason = ls->reason;
1094   PetscFunctionReturn(0);
1095 }
1096 
1097 /*@
1098   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1099   search.
1100 
1101   Not Collective
1102 
1103   Input Parameter:
1104 . ls - the TaoLineSearch context
1105 
1106   Output Parameter:
1107 . x - The initial point of the line search
1108 
1109   Level: intermediate
1110 @*/
1111 PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1112 {
1113   PetscFunctionBegin;
1114   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1115   if (x) *x = ls->start_x;
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 /*@
1120   TaoLineSearchGetStepDirection - Gets the step direction of the line
1121   search.
1122 
1123   Not Collective
1124 
1125   Input Parameter:
1126 . ls - the TaoLineSearch context
1127 
1128   Output Parameter:
1129 . s - the step direction of the line search
1130 
1131   Level: advanced
1132 @*/
1133 PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1134 {
1135   PetscFunctionBegin;
1136   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1137   if (s) *s = ls->stepdirection;
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 /*@
1142   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.
1143 
1144   Not Collective
1145 
1146   Input Parameter:
1147 . ls - the TaoLineSearch context
1148 
1149   Output Parameter:
1150 . f - the objective value at the full step length
1151 
1152   Level: developer
1153 @*/
1154 
1155 PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1156 {
1157   PetscFunctionBegin;
1158   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1159   *f_fullstep = ls->f_fullstep;
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 /*@
1164   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.
1165 
1166   Logically Collective on Tao
1167 
1168   Input Parameters:
1169 + ls - the TaoLineSearch context
1170 . xl  - vector of lower bounds
1171 - xu  - vector of upper bounds
1172 
1173   Note: If the variable bounds are not set with this routine, then
1174   PETSC_NINFINITY and PETSC_INFINITY are assumed
1175 
1176   Level: beginner
1177 
1178 .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1179 @*/
1180 PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1181 {
1182   PetscFunctionBegin;
1183   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1184   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1185   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
1186   PetscCall(PetscObjectReference((PetscObject)xl));
1187   PetscCall(PetscObjectReference((PetscObject)xu));
1188   PetscCall(VecDestroy(&ls->lower));
1189   PetscCall(VecDestroy(&ls->upper));
1190   ls->lower = xl;
1191   ls->upper = xu;
1192   ls->bounded = 1;
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 /*@
1197   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
1198   search.  If this value is not set then 1.0 is assumed.
1199 
1200   Logically Collective on TaoLineSearch
1201 
1202   Input Parameters:
1203 + ls - the TaoLineSearch context
1204 - s - the initial step size
1205 
1206   Level: intermediate
1207 
1208 .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1209 @*/
1210 PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1211 {
1212   PetscFunctionBegin;
1213   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1214   ls->initstep = s;
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 /*@
1219   TaoLineSearchGetStepLength - Get the current step length
1220 
1221   Not Collective
1222 
1223   Input Parameters:
1224 . ls - the TaoLineSearch context
1225 
1226   Output Parameters:
1227 . s - the current step length
1228 
1229   Level: beginner
1230 
1231 .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
1232 @*/
1233 PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
1234 {
1235   PetscFunctionBegin;
1236   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1237   *s = ls->step;
1238   PetscFunctionReturn(0);
1239 }
1240 
1241 /*@C
1242    TaoLineSearchRegister - Adds a line-search algorithm to the registry
1243 
1244    Not collective
1245 
1246    Input Parameters:
1247 +  sname - name of a new user-defined solver
1248 -  func - routine to Create method context
1249 
1250    Notes:
1251    TaoLineSearchRegister() may be called multiple times to add several user-defined solvers.
1252 
1253    Sample usage:
1254 .vb
1255    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1256 .ve
1257 
1258    Then, your solver can be chosen with the procedural interface via
1259 $     TaoLineSearchSetType(ls,"my_linesearch")
1260    or at runtime via the option
1261 $     -tao_ls_type my_linesearch
1262 
1263    Level: developer
1264 
1265 @*/
1266 PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1267 {
1268   PetscFunctionBegin;
1269   PetscCall(TaoLineSearchInitializePackage());
1270   PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func));
1271   PetscFunctionReturn(0);
1272 }
1273 
1274 /*@C
1275    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1276    for all TaoLineSearch options in the database.
1277 
1278    Collective on TaoLineSearch
1279 
1280    Input Parameters:
1281 +  ls - the TaoLineSearch solver context
1282 -  prefix - the prefix string to prepend to all line search requests
1283 
1284    Notes:
1285    A hyphen (-) must NOT be given at the beginning of the prefix name.
1286    The first character of all runtime options is AUTOMATICALLY the hyphen.
1287 
1288    Level: advanced
1289 
1290 .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1291 @*/
1292 PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1293 {
1294   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1295 }
1296 
1297 /*@C
1298   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1299   TaoLineSearch options in the database
1300 
1301   Not Collective
1302 
1303   Input Parameters:
1304 . ls - the TaoLineSearch context
1305 
1306   Output Parameters:
1307 . prefix - pointer to the prefix string used is returned
1308 
1309   Notes:
1310     On the fortran side, the user should pass in a string 'prefix' of
1311   sufficient length to hold the prefix.
1312 
1313   Level: advanced
1314 
1315 .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
1316 @*/
1317 PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1318 {
1319   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
1320 }
1321 
1322 /*@C
1323    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1324    TaoLineSearch options in the database.
1325 
1326    Logically Collective on TaoLineSearch
1327 
1328    Input Parameters:
1329 +  ls - the TaoLineSearch context
1330 -  prefix - the prefix string to prepend to all TAO option requests
1331 
1332    Notes:
1333    A hyphen (-) must NOT be given at the beginning of the prefix name.
1334    The first character of all runtime options is AUTOMATICALLY the hyphen.
1335 
1336    For example, to distinguish between the runtime options for two
1337    different line searches, one could call
1338 .vb
1339       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1340       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1341 .ve
1342 
1343    This would enable use of different options for each system, such as
1344 .vb
1345       -sys1_tao_ls_type mt
1346       -sys2_tao_ls_type armijo
1347 .ve
1348 
1349    Level: advanced
1350 
1351 .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1352 @*/
1353 
1354 PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1355 {
1356   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1357 }
1358