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