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