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