xref: /petsc/src/tao/linesearch/interface/taolinesearch.c (revision ba92ff593176f3ffed64b48a0721b2817410e47a)
1*ba92ff59SBarry Smith #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/
2aaa7dc30SBarry Smith #include <petsc-private/taolinesearchimpl.h>
30c52818fSSatish Balay 
40c52818fSSatish Balay PetscBool TaoLineSearchInitialized = PETSC_FALSE;
56c23d075SBarry Smith PetscFunctionList TaoLineSearchList = NULL;
60c52818fSSatish Balay 
70c52818fSSatish Balay PetscClassId TAOLINESEARCH_CLASSID=0;
80c52818fSSatish Balay PetscLogEvent TaoLineSearch_ApplyEvent = 0, TaoLineSearch_EvalEvent=0;
90c52818fSSatish Balay 
100c52818fSSatish Balay #undef __FUNCT__
110c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchView"
120c52818fSSatish Balay /*@C
130c52818fSSatish Balay   TaoLineSearchView - Prints information about the TaoLineSearch
140c52818fSSatish Balay 
150c52818fSSatish Balay   Collective on TaoLineSearch
160c52818fSSatish Balay 
170c52818fSSatish Balay   InputParameters:
180c52818fSSatish Balay + ls - the TaoSolver context
190c52818fSSatish Balay - viewer - visualization context
200c52818fSSatish Balay 
210c52818fSSatish Balay   Options Database Key:
220c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search
230c52818fSSatish Balay 
240c52818fSSatish Balay   Notes:
250c52818fSSatish Balay   The available visualization contexts include
260c52818fSSatish Balay +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
270c52818fSSatish Balay -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
280c52818fSSatish Balay          output where only the first processor opens
290c52818fSSatish Balay          the file.  All other processors send their
300c52818fSSatish Balay          data to the first processor to print.
310c52818fSSatish Balay 
320c52818fSSatish Balay   Level: beginner
330c52818fSSatish Balay 
340c52818fSSatish Balay .seealso: PetscViewerASCIIOpen()
350c52818fSSatish Balay @*/
360c52818fSSatish Balay 
370c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
380c52818fSSatish Balay {
390c52818fSSatish Balay   PetscErrorCode          ierr;
400c52818fSSatish Balay   PetscBool               isascii, isstring;
410c52818fSSatish Balay   const TaoLineSearchType type;
420c52818fSSatish Balay   PetscBool               destroyviewer = PETSC_FALSE;
43f06e3bfaSBarry Smith 
440c52818fSSatish Balay   PetscFunctionBegin;
450c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
460c52818fSSatish Balay   if (!viewer) {
470c52818fSSatish Balay     destroyviewer = PETSC_TRUE;
480c52818fSSatish Balay     ierr = PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);CHKERRQ(ierr);
490c52818fSSatish Balay   }
500c52818fSSatish Balay   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
510c52818fSSatish Balay   PetscCheckSameComm(ls,1,viewer,2);
520c52818fSSatish Balay 
530c52818fSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
540c52818fSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
550c52818fSSatish Balay   if (isascii) {
560c52818fSSatish Balay     if (((PetscObject)ls)->prefix) {
570c52818fSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:(%s)\n",((PetscObject)ls)->prefix);CHKERRQ(ierr);
580c52818fSSatish Balay     } else {
590c52818fSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:\n");CHKERRQ(ierr);
600c52818fSSatish Balay     }
610c52818fSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
620c52818fSSatish Balay     ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr);
630c52818fSSatish Balay     if (type) {
640c52818fSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"type: %s\n",type);CHKERRQ(ierr);
650c52818fSSatish Balay     } else {
660c52818fSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"type: not set yet\n");CHKERRQ(ierr);
670c52818fSSatish Balay     }
680c52818fSSatish Balay     if (ls->ops->view) {
690c52818fSSatish Balay       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
700c52818fSSatish Balay       ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr);
710c52818fSSatish Balay       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
720c52818fSSatish Balay     }
730c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr);
74f06e3bfaSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);CHKERRQ(ierr);
750c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr);
760c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr);
770c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr);
780c52818fSSatish Balay 
790c52818fSSatish Balay     if (ls->bounded) {
800c52818fSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr);
810c52818fSSatish Balay     }
820c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %D\n",(int)ls->reason);CHKERRQ(ierr);
830c52818fSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
840c52818fSSatish Balay 
850c52818fSSatish Balay   } else if (isstring) {
860c52818fSSatish Balay     ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr);
870c52818fSSatish Balay     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
880c52818fSSatish Balay   }
890c52818fSSatish Balay   if (destroyviewer) {
900c52818fSSatish Balay     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
910c52818fSSatish Balay   }
920c52818fSSatish Balay   PetscFunctionReturn(0);
930c52818fSSatish Balay }
940c52818fSSatish Balay 
950c52818fSSatish Balay #undef __FUNCT__
960c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchCreate"
970c52818fSSatish Balay /*@C
980c52818fSSatish Balay   TaoLineSearchCreate - Creates a TAO Line Search object.  Algorithms in TAO that use
990c52818fSSatish Balay   line-searches will automatically create one.
1000c52818fSSatish Balay 
1010c52818fSSatish Balay   Collective on MPI_Comm
1020c52818fSSatish Balay 
1030c52818fSSatish Balay   Input Parameter:
1040c52818fSSatish Balay . comm - MPI communicator
1050c52818fSSatish Balay 
1060c52818fSSatish Balay   Output Parameter:
1070c52818fSSatish Balay . newls - the new TaoLineSearch context
1080c52818fSSatish Balay 
1090c52818fSSatish Balay   Available methods include:
1100c52818fSSatish Balay + more-thuente
1110c52818fSSatish Balay . gpcg
1120c52818fSSatish Balay - unit - Do not perform any line search
1130c52818fSSatish Balay 
1140c52818fSSatish Balay 
1150c52818fSSatish Balay    Options Database Keys:
1160c52818fSSatish Balay .   -tao_ls_type - select which method TAO should use
1170c52818fSSatish Balay 
1180c52818fSSatish Balay    Level: beginner
1190c52818fSSatish Balay 
1200c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
1210c52818fSSatish Balay @*/
1220c52818fSSatish Balay 
1230c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
1240c52818fSSatish Balay {
1250c52818fSSatish Balay   PetscErrorCode ierr;
1260c52818fSSatish Balay   TaoLineSearch  ls;
1270c52818fSSatish Balay 
1280c52818fSSatish Balay   PetscFunctionBegin;
1290c52818fSSatish Balay   PetscValidPointer(newls,2);
1306c23d075SBarry Smith   *newls = NULL;
1310c52818fSSatish Balay 
1320c52818fSSatish Balay  #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1330c52818fSSatish Balay   ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr);
1340c52818fSSatish Balay  #endif
1350c52818fSSatish Balay 
136f06e3bfaSBarry Smith   ierr = PetscHeaderCreate(ls,_p_TaoLineSearch,struct _TaoLineSearchOps,TAOLINESEARCH_CLASSID,"TaoLineSearch",0,0,comm,TaoLineSearchDestroy,TaoLineSearchView);CHKERRQ(ierr);
1370c52818fSSatish Balay   ls->bounded = 0;
1380c52818fSSatish Balay   ls->max_funcs=30;
1390c52818fSSatish Balay   ls->ftol = 0.0001;
1400c52818fSSatish Balay   ls->gtol = 0.9;
1410c52818fSSatish Balay   ls->rtol = 1.0e-10;
1420c52818fSSatish Balay   ls->stepmin=1.0e-20;
1430c52818fSSatish Balay   ls->stepmax=1.0e+20;
1440c52818fSSatish Balay   ls->step=1.0;
1450c52818fSSatish Balay   ls->nfeval=0;
1460c52818fSSatish Balay   ls->ngeval=0;
1470c52818fSSatish Balay   ls->nfgeval=0;
1480c52818fSSatish Balay 
1490c52818fSSatish Balay   ls->ops->computeobjective=0;
1500c52818fSSatish Balay   ls->ops->computegradient=0;
1510c52818fSSatish Balay   ls->ops->computeobjectiveandgradient=0;
1520c52818fSSatish Balay   ls->ops->computeobjectiveandgts=0;
1530c52818fSSatish Balay   ls->ops->setup=0;
1540c52818fSSatish Balay   ls->ops->apply=0;
1550c52818fSSatish Balay   ls->ops->view=0;
1560c52818fSSatish Balay   ls->ops->setfromoptions=0;
1570c52818fSSatish Balay   ls->ops->reset=0;
1580c52818fSSatish Balay   ls->ops->destroy=0;
1590c52818fSSatish Balay   ls->setupcalled=PETSC_FALSE;
1600c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
1610c52818fSSatish Balay   *newls = ls;
1620c52818fSSatish Balay   PetscFunctionReturn(0);
1630c52818fSSatish Balay }
1640c52818fSSatish Balay 
1650c52818fSSatish Balay #undef __FUNCT__
1660c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetUp"
1670c52818fSSatish Balay /*@
1680c52818fSSatish Balay   TaoLineSearchSetUp - Sets up the internal data structures for the later use
1690c52818fSSatish Balay   of a Tao solver
1700c52818fSSatish Balay 
1710c52818fSSatish Balay   Collective on ls
1720c52818fSSatish Balay 
1730c52818fSSatish Balay   Input Parameters:
1740c52818fSSatish Balay . ls - the TaoLineSearch context
1750c52818fSSatish Balay 
1760c52818fSSatish Balay   Notes:
1770c52818fSSatish Balay   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
1780c52818fSSatish Balay   automatically be called in TaoLineSearchSolve().  However, if the user
1790c52818fSSatish Balay   desires to call it explicitly, it should come after TaoLineSearchCreate()
1800c52818fSSatish Balay   but before TaoLineSearchApply().
1810c52818fSSatish Balay 
1820c52818fSSatish Balay   Level: developer
1830c52818fSSatish Balay 
1840c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
1850c52818fSSatish Balay @*/
1860c52818fSSatish Balay 
1870c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
1880c52818fSSatish Balay {
1890c52818fSSatish Balay   PetscErrorCode ierr;
1900c52818fSSatish Balay   const char     *default_type=TAOLINESEARCH_MT;
1910c52818fSSatish Balay   PetscBool      flg;
192f06e3bfaSBarry Smith 
1930c52818fSSatish Balay   PetscFunctionBegin;
1940c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1950c52818fSSatish Balay   if (ls->setupcalled) PetscFunctionReturn(0);
1960c52818fSSatish Balay   if (!((PetscObject)ls)->type_name) {
1970c52818fSSatish Balay     ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr);
1980c52818fSSatish Balay   }
1990c52818fSSatish Balay   if (ls->ops->setup) {
2000c52818fSSatish Balay     ierr = (*ls->ops->setup)(ls);CHKERRQ(ierr);
2010c52818fSSatish Balay   }
2020c52818fSSatish Balay   if (ls->usetaoroutines) {
2030c52818fSSatish Balay     ierr = TaoIsObjectiveDefined(ls->taosolver,&flg);CHKERRQ(ierr);
2040c52818fSSatish Balay     ls->hasobjective = flg;
2050c52818fSSatish Balay     ierr = TaoIsGradientDefined(ls->taosolver,&flg);CHKERRQ(ierr);
2060c52818fSSatish Balay     ls->hasgradient = flg;
2070c52818fSSatish Balay     ierr = TaoIsObjectiveAndGradientDefined(ls->taosolver,&flg);CHKERRQ(ierr);
2080c52818fSSatish Balay     ls->hasobjectiveandgradient = flg;
2090c52818fSSatish Balay   } else {
2100c52818fSSatish Balay     if (ls->ops->computeobjective) {
2110c52818fSSatish Balay       ls->hasobjective = PETSC_TRUE;
2120c52818fSSatish Balay     } else {
2130c52818fSSatish Balay       ls->hasobjective = PETSC_FALSE;
2140c52818fSSatish Balay     }
2150c52818fSSatish Balay     if (ls->ops->computegradient) {
2160c52818fSSatish Balay       ls->hasgradient = PETSC_TRUE;
2170c52818fSSatish Balay     } else {
2180c52818fSSatish Balay       ls->hasgradient = PETSC_FALSE;
2190c52818fSSatish Balay     }
2200c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
2210c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_TRUE;
2220c52818fSSatish Balay     } else {
2230c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_FALSE;
2240c52818fSSatish Balay     }
2250c52818fSSatish Balay   }
2260c52818fSSatish Balay   ls->setupcalled = PETSC_TRUE;
2270c52818fSSatish Balay   PetscFunctionReturn(0);
2280c52818fSSatish Balay }
2290c52818fSSatish Balay 
2300c52818fSSatish Balay #undef __FUNCT__
2310c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchReset"
2320c52818fSSatish Balay /*@
2330c52818fSSatish Balay   TaoLineSearchReset - Some line searches may carry state information
2340c52818fSSatish Balay   from one TaoLineSearchApply() to the next.  This function resets this
2350c52818fSSatish Balay   state information.
2360c52818fSSatish Balay 
2370c52818fSSatish Balay   Collective on TaoLineSearch
2380c52818fSSatish Balay 
2390c52818fSSatish Balay   Input Parameter:
2400c52818fSSatish Balay . ls - the TaoLineSearch context
2410c52818fSSatish Balay 
2420c52818fSSatish Balay   Level: developer
2430c52818fSSatish Balay 
2440c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
2450c52818fSSatish Balay @*/
2460c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
2470c52818fSSatish Balay {
2480c52818fSSatish Balay   PetscErrorCode ierr;
249050fc7a3SBarry Smith 
2500c52818fSSatish Balay   PetscFunctionBegin;
2510c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
2520c52818fSSatish Balay   if (ls->ops->reset) {
2530c52818fSSatish Balay     ierr = (*ls->ops->reset)(ls);CHKERRQ(ierr);
2540c52818fSSatish Balay   }
2550c52818fSSatish Balay   PetscFunctionReturn(0);
2560c52818fSSatish Balay }
257f06e3bfaSBarry Smith 
2580c52818fSSatish Balay #undef __FUNCT__
2590c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchDestroy"
2600c52818fSSatish Balay /*@
2610c52818fSSatish Balay   TaoLineSearchDestroy - Destroys the TAO context that was created with
2620c52818fSSatish Balay   TaoLineSearchCreate()
2630c52818fSSatish Balay 
2640c52818fSSatish Balay   Collective on TaoLineSearch
2650c52818fSSatish Balay 
2660c52818fSSatish Balay   Input Parameter
2670c52818fSSatish Balay . ls - the TaoLineSearch context
2680c52818fSSatish Balay 
2690c52818fSSatish Balay   Level: beginner
2700c52818fSSatish Balay 
2710c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
2720c52818fSSatish Balay @*/
2730c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
2740c52818fSSatish Balay {
2750c52818fSSatish Balay   PetscErrorCode ierr;
276050fc7a3SBarry Smith 
2770c52818fSSatish Balay   PetscFunctionBegin;
2780c52818fSSatish Balay   if (!*ls) PetscFunctionReturn(0);
2790c52818fSSatish Balay   PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1);
2800c52818fSSatish Balay   if (--((PetscObject)*ls)->refct > 0) {*ls=0; PetscFunctionReturn(0);}
281050fc7a3SBarry Smith   ierr = VecDestroy(&(*ls)->stepdirection);CHKERRQ(ierr);
282050fc7a3SBarry Smith   ierr = VecDestroy(&(*ls)->start_x);CHKERRQ(ierr);
2830c52818fSSatish Balay   if ((*ls)->ops->destroy) {
2840c52818fSSatish Balay     ierr = (*(*ls)->ops->destroy)(*ls);CHKERRQ(ierr);
2850c52818fSSatish Balay   }
2860c52818fSSatish Balay   ierr = PetscHeaderDestroy(ls);CHKERRQ(ierr);
2870c52818fSSatish Balay   PetscFunctionReturn(0);
2880c52818fSSatish Balay }
2890c52818fSSatish Balay 
2900c52818fSSatish Balay #undef __FUNCT__
2910c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchApply"
2920c52818fSSatish Balay /*@
2930c52818fSSatish Balay   TaoLineSearchApply - Performs a line-search in a given step direction.  Criteria for acceptable step length depends on the line-search algorithm chosen
2940c52818fSSatish Balay 
2950c52818fSSatish Balay   Collective on TaoLineSearch
2960c52818fSSatish Balay 
2970c52818fSSatish Balay   Input Parameters:
2980c52818fSSatish Balay + ls - the TaoSolver context
2990c52818fSSatish Balay . x - The current solution (on output x contains the new solution determined by the line search)
3000c52818fSSatish Balay . f - objective function value at current solution (on output contains the objective function value at new solution)
3010c52818fSSatish Balay . g - gradient evaluated at x (on output contains the gradient at new solution)
3020c52818fSSatish Balay - s - search direction
3030c52818fSSatish Balay 
3040c52818fSSatish Balay   Output Parameters:
3050c52818fSSatish Balay + x - new solution
3060c52818fSSatish Balay . f - objective function value at x
3070c52818fSSatish Balay . g - gradient vector at x
3080c52818fSSatish Balay . steplength - scalar multiplier of s used ( x = x0 + steplength * x )
3090c52818fSSatish Balay - reason - reason why the line-search stopped
3100c52818fSSatish Balay 
3110c52818fSSatish Balay   reason will be set to one of:
3120c52818fSSatish Balay 
3130c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
3140c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
3150c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
3160c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
3170c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
3180c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
3190c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
3200c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
3210c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
3220c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search
3230c52818fSSatish Balay 
3240c52818fSSatish Balay   Note:
3250c52818fSSatish Balay   The algorithm developer must set up the TaoLineSearch with calls to
3260c52818fSSatish Balay   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoSolverRoutines()
3270c52818fSSatish Balay 
3280c52818fSSatish Balay   Note:
3290c52818fSSatish Balay   You may or may not need to follow this with a call to
3300c52818fSSatish Balay   TaoAddLineSearchCounts(), depending on whether you want these
3310c52818fSSatish Balay   evaluations to count toward the total function/gradient evaluations.
3320c52818fSSatish Balay 
3330c52818fSSatish Balay   Level: beginner
3340c52818fSSatish Balay 
3350c52818fSSatish Balay   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
3360c52818fSSatish Balay  @*/
3370c52818fSSatish Balay 
3380c52818fSSatish Balay PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchTerminationReason *reason)
3390c52818fSSatish Balay {
3400c52818fSSatish Balay   PetscErrorCode ierr;
3410c52818fSSatish Balay   PetscViewer    viewer;
3420c52818fSSatish Balay   PetscInt       low1,low2,low3,high1,high2,high3;
3430c52818fSSatish Balay   PetscBool      flg;
3440c52818fSSatish Balay   char           filename[PETSC_MAX_PATH_LEN];
3450c52818fSSatish Balay 
3460c52818fSSatish Balay   PetscFunctionBegin;
3470c52818fSSatish Balay   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
3480c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
3490c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
3500c52818fSSatish Balay   PetscValidScalarPointer(f,3);
3510c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
3520c52818fSSatish Balay   PetscValidHeaderSpecific(s,VEC_CLASSID,5);
3530c52818fSSatish Balay   PetscValidPointer(reason,7);
3540c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
3550c52818fSSatish Balay   PetscCheckSameTypeAndComm(x,2,g,4);
3560c52818fSSatish Balay   PetscCheckSameTypeAndComm(x,2,s,5);
3570c52818fSSatish Balay   ierr = VecGetOwnershipRange(x, &low1, &high1);CHKERRQ(ierr);
3580c52818fSSatish Balay   ierr = VecGetOwnershipRange(g, &low2, &high2);CHKERRQ(ierr);
3590c52818fSSatish Balay   ierr = VecGetOwnershipRange(s, &low3, &high3);CHKERRQ(ierr);
360f06e3bfaSBarry Smith   if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths");
3610c52818fSSatish Balay 
3620c52818fSSatish Balay   ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr);
363050fc7a3SBarry Smith   ierr = VecDestroy(&ls->stepdirection);CHKERRQ(ierr);
364050fc7a3SBarry Smith   ls->stepdirection = s;
3650c52818fSSatish Balay 
3660c52818fSSatish Balay   ierr = TaoLineSearchSetUp(ls);CHKERRQ(ierr);
367f06e3bfaSBarry Smith   if (!ls->ops->apply) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine");
3680c52818fSSatish Balay   ls->nfeval=0;
3690c52818fSSatish Balay   ls->ngeval=0;
3700c52818fSSatish Balay   ls->nfgeval=0;
3710c52818fSSatish Balay   /* Check parameter values */
3720c52818fSSatish Balay   if (ls->ftol < 0.0) {
373f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);CHKERRQ(ierr);
3740c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3750c52818fSSatish Balay   }
3760c52818fSSatish Balay   if (ls->rtol < 0.0) {
377f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);CHKERRQ(ierr);
3780c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3790c52818fSSatish Balay   }
3800c52818fSSatish Balay   if (ls->gtol < 0.0) {
381f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);CHKERRQ(ierr);
3820c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3830c52818fSSatish Balay   }
3840c52818fSSatish Balay   if (ls->stepmin < 0.0) {
385f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);CHKERRQ(ierr);
3860c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3870c52818fSSatish Balay   }
3880c52818fSSatish Balay   if (ls->stepmax < ls->stepmin) {
389f06e3bfaSBarry Smith     ierr = PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);CHKERRQ(ierr);
3900c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3910c52818fSSatish Balay   }
3920c52818fSSatish Balay   if (ls->max_funcs < 0) {
3930c52818fSSatish Balay     ierr = PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);CHKERRQ(ierr);
3940c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3950c52818fSSatish Balay   }
3960c52818fSSatish Balay   if (PetscIsInfOrNanReal(*f)) {
397f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);CHKERRQ(ierr);
3980c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_INFORNAN;
3990c52818fSSatish Balay   }
4000c52818fSSatish Balay 
4010c52818fSSatish Balay   ierr = PetscObjectReference((PetscObject)x);
4020c52818fSSatish Balay   ierr = VecDestroy(&ls->start_x);CHKERRQ(ierr);
4030c52818fSSatish Balay   ls->start_x = x;
404050fc7a3SBarry Smith 
4050c52818fSSatish Balay   ierr = PetscLogEventBegin(TaoLineSearch_ApplyEvent,ls,0,0,0);CHKERRQ(ierr);
4060c52818fSSatish Balay   ierr = (*ls->ops->apply)(ls,x,f,g,s);CHKERRQ(ierr);
4070c52818fSSatish Balay   ierr = PetscLogEventEnd(TaoLineSearch_ApplyEvent, ls, 0,0,0);CHKERRQ(ierr);
4080c52818fSSatish Balay   *reason=ls->reason;
4090c52818fSSatish Balay   ls->new_f = *f;
4100c52818fSSatish Balay 
4110c52818fSSatish Balay   if (steplength) {
4120c52818fSSatish Balay     *steplength=ls->step;
4130c52818fSSatish Balay   }
4140c52818fSSatish Balay 
4150c52818fSSatish Balay   ierr = PetscOptionsGetString(((PetscObject)ls)->prefix,"-tao_ls_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
4160c52818fSSatish Balay   if (ls->viewls && !PetscPreLoadingOn) {
4170c52818fSSatish Balay     ierr = PetscViewerASCIIOpen(((PetscObject)ls)->comm,filename,&viewer);CHKERRQ(ierr);
4180c52818fSSatish Balay     ierr = TaoLineSearchView(ls,viewer);CHKERRQ(ierr);
4190c52818fSSatish Balay     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
4200c52818fSSatish Balay   }
4210c52818fSSatish Balay   PetscFunctionReturn(0);
4220c52818fSSatish Balay }
4230c52818fSSatish Balay 
4240c52818fSSatish Balay #undef __FUNCT__
4250c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetType"
4260c52818fSSatish Balay /*@C
4270c52818fSSatish Balay    TaoLineSearchSetType - Sets the algorithm used in a line search
4280c52818fSSatish Balay 
4290c52818fSSatish Balay    Collective on TaoLineSearch
4300c52818fSSatish Balay 
4310c52818fSSatish Balay    Input Parameters:
4320c52818fSSatish Balay +  ls - the TaoLineSearch context
4330c52818fSSatish Balay -  type - a known method
4340c52818fSSatish Balay 
4350c52818fSSatish Balay   Available methods include:
4360c52818fSSatish Balay + more-thuente
4370c52818fSSatish Balay . gpcg
4380c52818fSSatish Balay - unit - Do not perform any line search
4390c52818fSSatish Balay 
4400c52818fSSatish Balay 
4410c52818fSSatish Balay   Options Database Keys:
4420c52818fSSatish Balay .   -tao_ls_type - select which method TAO should use
4430c52818fSSatish Balay 
4440c52818fSSatish Balay   Level: beginner
4450c52818fSSatish Balay 
4460c52818fSSatish Balay 
4470c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply()
4480c52818fSSatish Balay 
4490c52818fSSatish Balay @*/
4500c52818fSSatish Balay 
4510c52818fSSatish Balay PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, const TaoLineSearchType type)
4520c52818fSSatish Balay {
4530c52818fSSatish Balay   PetscErrorCode ierr;
4540c52818fSSatish Balay   PetscErrorCode (*r)(TaoLineSearch);
4550c52818fSSatish Balay   PetscBool      flg;
4560c52818fSSatish Balay 
4570c52818fSSatish Balay   PetscFunctionBegin;
4580c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
4590c52818fSSatish Balay   PetscValidCharPointer(type,2);
4600c52818fSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg);CHKERRQ(ierr);
4610c52818fSSatish Balay   if (flg) PetscFunctionReturn(0);
4620c52818fSSatish Balay 
4630c52818fSSatish Balay   ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);CHKERRQ(ierr);
4640c52818fSSatish Balay   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type);
4650c52818fSSatish Balay   if (ls->ops->destroy) {
4660c52818fSSatish Balay     ierr = (*(ls)->ops->destroy)(ls);CHKERRQ(ierr);
4670c52818fSSatish Balay   }
4680c52818fSSatish Balay   ls->max_funcs=30;
4690c52818fSSatish Balay   ls->ftol = 0.0001;
4700c52818fSSatish Balay   ls->gtol = 0.9;
4710c52818fSSatish Balay   ls->rtol = 1.0e-10;
4720c52818fSSatish Balay   ls->stepmin=1.0e-20;
4730c52818fSSatish Balay   ls->stepmax=1.0e+20;
4740c52818fSSatish Balay 
4750c52818fSSatish Balay   ls->nfeval=0;
4760c52818fSSatish Balay   ls->ngeval=0;
4770c52818fSSatish Balay   ls->nfgeval=0;
4780c52818fSSatish Balay   ls->ops->setup=0;
4790c52818fSSatish Balay   ls->ops->apply=0;
4800c52818fSSatish Balay   ls->ops->view=0;
4810c52818fSSatish Balay   ls->ops->setfromoptions=0;
4820c52818fSSatish Balay   ls->ops->destroy=0;
4830c52818fSSatish Balay   ls->setupcalled = PETSC_FALSE;
4840c52818fSSatish Balay   ierr = (*r)(ls);CHKERRQ(ierr);
4850c52818fSSatish Balay   ierr = PetscObjectChangeTypeName((PetscObject)ls, type);CHKERRQ(ierr);
4860c52818fSSatish Balay   PetscFunctionReturn(0);
4870c52818fSSatish Balay }
4880c52818fSSatish Balay 
4890c52818fSSatish Balay #undef __FUNCT__
4900c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetFromOptions"
4910c52818fSSatish Balay /*@
4920c52818fSSatish Balay   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
4930c52818fSSatish Balay   options.
4940c52818fSSatish Balay 
4950c52818fSSatish Balay   Collective on TaoLineSearch
4960c52818fSSatish Balay 
4970c52818fSSatish Balay   Input Paremeter:
4980c52818fSSatish Balay . ls - the TaoLineSearch context
4990c52818fSSatish Balay 
5000c52818fSSatish Balay   Options Database Keys:
5010c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
5020c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease
5030c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition
5040c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step
5050c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed
5060c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed
5070c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
5080c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output
5090c52818fSSatish Balay 
5100c52818fSSatish Balay   Level: beginner
5110c52818fSSatish Balay @*/
5120c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
5130c52818fSSatish Balay {
5140c52818fSSatish Balay   PetscErrorCode ierr;
5150c52818fSSatish Balay   const char     *default_type=TAOLINESEARCH_MT;
5160c52818fSSatish Balay   char           type[256];
5170c52818fSSatish Balay   PetscBool      flg;
518f06e3bfaSBarry Smith 
5190c52818fSSatish Balay   PetscFunctionBegin;
5200c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
5210c52818fSSatish Balay   ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr);
5220c52818fSSatish Balay   if (!TaoLineSearchInitialized) {
5230c52818fSSatish Balay     ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr);
5240c52818fSSatish Balay   }
5250c52818fSSatish Balay   if (((PetscObject)ls)->type_name) {
5260c52818fSSatish Balay     default_type = ((PetscObject)ls)->type_name;
5270c52818fSSatish Balay   }
5280c52818fSSatish Balay   /* Check for type from options */
5290c52818fSSatish Balay   ierr = PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);CHKERRQ(ierr);
5300c52818fSSatish Balay   if (flg) {
5310c52818fSSatish Balay     ierr = TaoLineSearchSetType(ls,type);CHKERRQ(ierr);
5320c52818fSSatish Balay   } else if (!((PetscObject)ls)->type_name) {
5330c52818fSSatish Balay     ierr = TaoLineSearchSetType(ls,default_type);
5340c52818fSSatish Balay   }
5350c52818fSSatish Balay 
536f06e3bfaSBarry Smith   ierr = PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,0);CHKERRQ(ierr);
537f06e3bfaSBarry Smith   ierr = PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,0);CHKERRQ(ierr);
538f06e3bfaSBarry Smith   ierr = PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,0);CHKERRQ(ierr);
539f06e3bfaSBarry Smith   ierr = PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,0);CHKERRQ(ierr);
540f06e3bfaSBarry Smith   ierr = PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,0);CHKERRQ(ierr);
541f06e3bfaSBarry Smith   ierr = PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,0);CHKERRQ(ierr);
5426c23d075SBarry Smith   ierr = PetscOptionsBool("-tao_ls_view","view TaoLineSearch info after each line search has completed","TaoLineSearchView",PETSC_FALSE,&ls->viewls,NULL);CHKERRQ(ierr);
5430c52818fSSatish Balay   if (ls->ops->setfromoptions) {
5440c52818fSSatish Balay     ierr = (*ls->ops->setfromoptions)(ls);CHKERRQ(ierr);
5450c52818fSSatish Balay   }
5460c52818fSSatish Balay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5470c52818fSSatish Balay   PetscFunctionReturn(0);
5480c52818fSSatish Balay }
5490c52818fSSatish Balay 
5500c52818fSSatish Balay #undef __FUNCT__
5510c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetType"
5520c52818fSSatish Balay /*@C
5530c52818fSSatish Balay   TaoLineSearchGetType - Gets the current line search algorithm
5540c52818fSSatish Balay 
5550c52818fSSatish Balay   Not Collective
5560c52818fSSatish Balay 
5570c52818fSSatish Balay   Input Parameter:
5580c52818fSSatish Balay . ls - the TaoLineSearch context
5590c52818fSSatish Balay 
5600c52818fSSatish Balay   Output Paramter:
5610c52818fSSatish Balay . type - the line search algorithm in effect
5620c52818fSSatish Balay 
5630c52818fSSatish Balay   Level: developer
5640c52818fSSatish Balay 
5650c52818fSSatish Balay @*/
5660c52818fSSatish Balay PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, const TaoLineSearchType *type)
5670c52818fSSatish Balay {
5680c52818fSSatish Balay   PetscFunctionBegin;
5690c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
5700c52818fSSatish Balay   PetscValidPointer(type,2);
5710c52818fSSatish Balay   *type = ((PetscObject)ls)->type_name;
5720c52818fSSatish Balay   PetscFunctionReturn(0);
5730c52818fSSatish Balay }
5740c52818fSSatish Balay 
5750c52818fSSatish Balay #undef __FUNCT__
5760c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetNumberFunctionEvaluations"
5770c52818fSSatish Balay /*@
5780c52818fSSatish Balay   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
5790c52818fSSatish Balay   routines used by the line search in last application (not cumulative).
5800c52818fSSatish Balay 
5810c52818fSSatish Balay   Not Collective
5820c52818fSSatish Balay 
5830c52818fSSatish Balay   Input Parameter:
5840c52818fSSatish Balay . ls - the TaoLineSearch context
5850c52818fSSatish Balay 
5860c52818fSSatish Balay   Output Parameters:
5870c52818fSSatish Balay + nfeval   - number of function evaluations
5880c52818fSSatish Balay . ngeval   - number of gradient evaluations
5890c52818fSSatish Balay - nfgeval  - number of function/gradient evaluations
5900c52818fSSatish Balay 
5910c52818fSSatish Balay   Level: intermediate
5920c52818fSSatish Balay 
5930c52818fSSatish Balay   Note:
5940c52818fSSatish Balay   If the line search is using the TaoSolver objective and gradient
5950c52818fSSatish Balay   routines directly (see TaoLineSearchUseTaoSolverRoutines()), then TAO
5960c52818fSSatish Balay   is already counting the number of evaluations.
5970c52818fSSatish Balay 
5980c52818fSSatish Balay @*/
5990c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
6000c52818fSSatish Balay {
6010c52818fSSatish Balay   PetscFunctionBegin;
6020c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6030c52818fSSatish Balay   *nfeval = ls->nfeval;
6040c52818fSSatish Balay   *ngeval = ls->ngeval;
6050c52818fSSatish Balay   *nfgeval = ls->nfgeval;
6060c52818fSSatish Balay   PetscFunctionReturn(0);
6070c52818fSSatish Balay }
6080c52818fSSatish Balay 
6090c52818fSSatish Balay #undef __FUNCT__
6100c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchIsUsingTaoSolverRoutines"
6110c52818fSSatish Balay /*@
6120c52818fSSatish Balay   TaoLineSearchIsUsingTaoSolverRoutines - Checks whether the line search is using
6130c52818fSSatish Balay   TaoSolver evaluation routines.
6140c52818fSSatish Balay 
6150c52818fSSatish Balay   Not Collective
6160c52818fSSatish Balay 
6170c52818fSSatish Balay   Input Parameter:
6180c52818fSSatish Balay . ls - the TaoLineSearch context
6190c52818fSSatish Balay 
6200c52818fSSatish Balay   Output Parameter:
6210c52818fSSatish Balay . flg - PETSC_TRUE if the line search is using TaoSolver evaluation routines,
6220c52818fSSatish Balay         otherwise PETSC_FALSE
6230c52818fSSatish Balay 
6240c52818fSSatish Balay   Level: developer
6250c52818fSSatish Balay @*/
6260c52818fSSatish Balay PetscErrorCode TaoLineSearchIsUsingTaoSolverRoutines(TaoLineSearch ls, PetscBool *flg)
6270c52818fSSatish Balay {
6280c52818fSSatish Balay   PetscFunctionBegin;
6290c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
630f06e3bfaSBarry Smith   *flg = ls->usetaoroutines;
6310c52818fSSatish Balay   PetscFunctionReturn(0);
6320c52818fSSatish Balay }
6330c52818fSSatish Balay 
6340c52818fSSatish Balay #undef __FUNCT__
6350c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetObjectiveRoutine"
6360c52818fSSatish Balay /*@C
6370c52818fSSatish Balay   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
6380c52818fSSatish Balay 
6390c52818fSSatish Balay   Logically Collective on TaoLineSearch
6400c52818fSSatish Balay 
6410c52818fSSatish Balay   Input Parameter:
6420c52818fSSatish Balay + ls - the TaoLineSearch context
6430c52818fSSatish Balay . func - the objective function evaluation routine
6440c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6450c52818fSSatish Balay 
6460c52818fSSatish Balay   Calling sequence of func:
6470c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);
6480c52818fSSatish Balay 
6490c52818fSSatish Balay + x - input vector
6500c52818fSSatish Balay . f - function value
6510c52818fSSatish Balay - ctx (optional) user-defined context
6520c52818fSSatish Balay 
6530c52818fSSatish Balay   Level: beginner
6540c52818fSSatish Balay 
6550c52818fSSatish Balay   Note:
6560c52818fSSatish Balay   Use this routine only if you want the line search objective
6570c52818fSSatish Balay   evaluation routine to be different from the TaoSolver's objective
6580c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
6590c52818fSSatish Balay   the line search gradient and/or function/gradient routine.
6600c52818fSSatish Balay 
6610c52818fSSatish Balay   Note:
6620c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
6630c52818fSSatish Balay   line search, application programmers should be wary of overriding the
6640c52818fSSatish Balay   default objective routine.
6650c52818fSSatish Balay 
6660c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoSolverRoutines()
6670c52818fSSatish Balay @*/
6680c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
6690c52818fSSatish Balay {
6700c52818fSSatish Balay   PetscFunctionBegin;
6710c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6720c52818fSSatish Balay 
6730c52818fSSatish Balay   ls->ops->computeobjective=func;
6740c52818fSSatish Balay   if (ctx) ls->userctx_func=ctx;
6750c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
6760c52818fSSatish Balay   PetscFunctionReturn(0);
6770c52818fSSatish Balay }
6780c52818fSSatish Balay 
6790c52818fSSatish Balay #undef __FUNCT__
6800c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetGradientRoutine"
6810c52818fSSatish Balay /*@C
6820c52818fSSatish Balay   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
6830c52818fSSatish Balay 
6840c52818fSSatish Balay   Logically Collective on TaoLineSearch
6850c52818fSSatish Balay 
6860c52818fSSatish Balay   Input Parameter:
6870c52818fSSatish Balay + ls - the TaoLineSearch context
6880c52818fSSatish Balay . func - the gradient evaluation routine
6890c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6900c52818fSSatish Balay 
6910c52818fSSatish Balay   Calling sequence of func:
6920c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);
6930c52818fSSatish Balay 
6940c52818fSSatish Balay + x - input vector
6950c52818fSSatish Balay . g - gradient vector
6960c52818fSSatish Balay - ctx (optional) user-defined context
6970c52818fSSatish Balay 
6980c52818fSSatish Balay   Level: beginner
6990c52818fSSatish Balay 
7000c52818fSSatish Balay   Note:
7010c52818fSSatish Balay   Use this routine only if you want the line search gradient
7020c52818fSSatish Balay   evaluation routine to be different from the TaoSolver's gradient
7030c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
7040c52818fSSatish Balay   the line search function and/or function/gradient routine.
7050c52818fSSatish Balay 
7060c52818fSSatish Balay   Note:
7070c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own gradient routine for the
7080c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7090c52818fSSatish Balay   default gradient routine.
7100c52818fSSatish Balay 
7110c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoSolverRoutines()
7120c52818fSSatish Balay @*/
7130c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
7140c52818fSSatish Balay {
7150c52818fSSatish Balay   PetscFunctionBegin;
7160c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
7170c52818fSSatish Balay   ls->ops->computegradient=func;
7180c52818fSSatish Balay   if (ctx) ls->userctx_grad=ctx;
7190c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
7200c52818fSSatish Balay   PetscFunctionReturn(0);
7210c52818fSSatish Balay }
7220c52818fSSatish Balay 
7230c52818fSSatish Balay #undef __FUNCT__
7240c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetObjectiveAndGradientRoutine"
7250c52818fSSatish Balay /*@C
7260c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
7270c52818fSSatish Balay 
7280c52818fSSatish Balay   Logically Collective on TaoLineSearch
7290c52818fSSatish Balay 
7300c52818fSSatish Balay   Input Parameter:
7310c52818fSSatish Balay + ls - the TaoLineSearch context
7320c52818fSSatish Balay . func - the objective and gradient evaluation routine
7330c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7340c52818fSSatish Balay 
7350c52818fSSatish Balay   Calling sequence of func:
7360c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);
7370c52818fSSatish Balay 
7380c52818fSSatish Balay + x - input vector
7390c52818fSSatish Balay . f - function value
7400c52818fSSatish Balay . g - gradient vector
7410c52818fSSatish Balay - ctx (optional) user-defined context
7420c52818fSSatish Balay 
7430c52818fSSatish Balay   Level: beginner
7440c52818fSSatish Balay 
7450c52818fSSatish Balay   Note:
7460c52818fSSatish Balay   Use this routine only if you want the line search objective and gradient
7470c52818fSSatish Balay   evaluation routines to be different from the TaoSolver's objective
7480c52818fSSatish Balay   and gradient evaluation routines.
7490c52818fSSatish Balay 
7500c52818fSSatish Balay   Note:
7510c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
7520c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7530c52818fSSatish Balay   default objective routine.
7540c52818fSSatish Balay 
7550c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoSolverRoutines()
7560c52818fSSatish Balay @*/
7570c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
7580c52818fSSatish Balay {
7590c52818fSSatish Balay   PetscFunctionBegin;
7600c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
7610c52818fSSatish Balay   ls->ops->computeobjectiveandgradient=func;
7620c52818fSSatish Balay   if (ctx) ls->userctx_funcgrad=ctx;
7630c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
7640c52818fSSatish Balay   PetscFunctionReturn(0);
7650c52818fSSatish Balay }
7660c52818fSSatish Balay 
7670c52818fSSatish Balay #undef __FUNCT__
7680c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetObjectiveAndGTSRoutine"
7690c52818fSSatish Balay /*@C
7700c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
7710c52818fSSatish Balay   (gradient'*stepdirection) evaluation routine for the line search.
7720c52818fSSatish Balay   Sometimes it is more efficient to compute the inner product of the gradient
7730c52818fSSatish Balay   and the step direction than it is to compute the gradient, and this is all
7740c52818fSSatish Balay   the line search typically needs of the gradient.
7750c52818fSSatish Balay 
7760c52818fSSatish Balay   Logically Collective on TaoLineSearch
7770c52818fSSatish Balay 
7780c52818fSSatish Balay   Input Parameter:
7790c52818fSSatish Balay + ls - the TaoLineSearch context
7800c52818fSSatish Balay . func - the objective and gradient evaluation routine
7810c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7820c52818fSSatish Balay 
7830c52818fSSatish Balay   Calling sequence of func:
7840c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);
7850c52818fSSatish Balay 
7860c52818fSSatish Balay + x - input vector
7870c52818fSSatish Balay . s - step direction
7880c52818fSSatish Balay . f - function value
7890c52818fSSatish Balay . gts - inner product of gradient and step direction vectors
7900c52818fSSatish Balay - ctx (optional) user-defined context
7910c52818fSSatish Balay 
7920c52818fSSatish Balay   Note: The gradient will still need to be computed at the end of the line
7930c52818fSSatish Balay   search, so you will still need to set a line search gradient evaluation
7940c52818fSSatish Balay   routine
7950c52818fSSatish Balay 
7960c52818fSSatish Balay   Note: Bounded line searches (those used in bounded optimization algorithms)
7970c52818fSSatish Balay   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
7980c52818fSSatish Balay   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()
7990c52818fSSatish Balay 
8000c52818fSSatish Balay   Level: advanced
8010c52818fSSatish Balay 
8020c52818fSSatish Balay   Note:
8030c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
8040c52818fSSatish Balay   line search, application programmers should be wary of overriding the
8050c52818fSSatish Balay   default objective routine.
8060c52818fSSatish Balay 
8070c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoSolverRoutines()
8080c52818fSSatish Balay @*/
8090c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
8100c52818fSSatish Balay {
8110c52818fSSatish Balay   PetscFunctionBegin;
8120c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8130c52818fSSatish Balay   ls->ops->computeobjectiveandgts=func;
8140c52818fSSatish Balay   if (ctx) ls->userctx_funcgts=ctx;
8150c52818fSSatish Balay   ls->usegts = PETSC_TRUE;
8160c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
8170c52818fSSatish Balay   PetscFunctionReturn(0);
8180c52818fSSatish Balay }
8190c52818fSSatish Balay 
8200c52818fSSatish Balay #undef __FUNCT__
8210c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchUseTaoSolverRoutines"
8220c52818fSSatish Balay /*@C
8230c52818fSSatish Balay   TaoLineSearchUseTaoSolverRoutines - Informs the TaoLineSearch to use the
8240c52818fSSatish Balay   objective and gradient evaluation routines from the given TaoSolver object.
8250c52818fSSatish Balay 
8260c52818fSSatish Balay   Logically Collective on TaoLineSearch
8270c52818fSSatish Balay 
8280c52818fSSatish Balay   Input Parameter:
8290c52818fSSatish Balay + ls - the TaoLineSearch context
8300c52818fSSatish Balay - ts - the TaoSolver context with defined objective/gradient evaluation routines
8310c52818fSSatish Balay 
8320c52818fSSatish Balay   Level: developer
8330c52818fSSatish Balay 
8340c52818fSSatish Balay .seealso: TaoLineSearchCreate()
8350c52818fSSatish Balay @*/
8360c52818fSSatish Balay PetscErrorCode TaoLineSearchUseTaoSolverRoutines(TaoLineSearch ls, TaoSolver ts)
8370c52818fSSatish Balay {
8380c52818fSSatish Balay   PetscFunctionBegin;
8390c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8400c52818fSSatish Balay   PetscValidHeaderSpecific(ts,TAOSOLVER_CLASSID,1);
8410c52818fSSatish Balay   ls->taosolver = ts;
8420c52818fSSatish Balay   ls->usetaoroutines=PETSC_TRUE;
8430c52818fSSatish Balay   PetscFunctionReturn(0);
8440c52818fSSatish Balay }
8450c52818fSSatish Balay 
8460c52818fSSatish Balay #undef __FUNCT__
8470c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeObjective"
8480c52818fSSatish Balay /*@
8490c52818fSSatish Balay   TaoLineSearchComputeObjective - Computes the objective function value at a given point
8500c52818fSSatish Balay 
8510c52818fSSatish Balay   Collective on TaoLineSearch
8520c52818fSSatish Balay 
8530c52818fSSatish Balay   Input Parameters:
8540c52818fSSatish Balay + ls - the TaoLineSearch context
8550c52818fSSatish Balay - x - input vector
8560c52818fSSatish Balay 
8570c52818fSSatish Balay   Output Parameter:
8580c52818fSSatish Balay . f - Objective value at X
8590c52818fSSatish Balay 
8600c52818fSSatish Balay   Notes: TaoLineSearchComputeObjective() is typically used within line searches
8610c52818fSSatish Balay   so most users would not generally call this routine themselves.
8620c52818fSSatish Balay 
8630c52818fSSatish Balay   Level: developer
8640c52818fSSatish Balay 
8650c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
8660c52818fSSatish Balay @*/
8670c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
8680c52818fSSatish Balay {
8690c52818fSSatish Balay   PetscErrorCode ierr;
8700c52818fSSatish Balay   Vec            gdummy;
8710c52818fSSatish Balay   PetscReal      gts;
872f06e3bfaSBarry Smith 
8730c52818fSSatish Balay   PetscFunctionBegin;
8740c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8750c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
8760c52818fSSatish Balay   PetscValidPointer(f,3);
8770c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
8780c52818fSSatish Balay   if (ls->usetaoroutines) {
8790c52818fSSatish Balay     ierr = TaoComputeObjective(ls->taosolver,x,f);CHKERRQ(ierr);
8800c52818fSSatish Balay   } else {
8810c52818fSSatish Balay     ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
882f06e3bfaSBarry Smith     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");
8830c52818fSSatish Balay     PetscStackPush("TaoLineSearch user objective routine");
8840c52818fSSatish Balay     if (ls->ops->computeobjective) {
8850c52818fSSatish Balay       ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr);
8860c52818fSSatish Balay     } else if (ls->ops->computeobjectiveandgradient) {
8870c52818fSSatish Balay       ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr);
8880c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr);
8890c52818fSSatish Balay       ierr = VecDestroy(&gdummy);CHKERRQ(ierr);
8900c52818fSSatish Balay     } else {
8910c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);CHKERRQ(ierr);
8920c52818fSSatish Balay     }
8930c52818fSSatish Balay     PetscStackPop;
8940c52818fSSatish Balay     ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
8950c52818fSSatish Balay   }
8960c52818fSSatish Balay   ls->nfeval++;
8970c52818fSSatish Balay   PetscFunctionReturn(0);
8980c52818fSSatish Balay }
8990c52818fSSatish Balay 
9000c52818fSSatish Balay #undef __FUNCT__
9010c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeObjectiveAndGradient"
9020c52818fSSatish Balay /*@
9030c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
9040c52818fSSatish Balay 
9050c52818fSSatish Balay   Collective on TaoSOlver
9060c52818fSSatish Balay 
9070c52818fSSatish Balay   Input Parameters:
9080c52818fSSatish Balay + ls - the TaoLineSearch context
9090c52818fSSatish Balay - x - input vector
9100c52818fSSatish Balay 
9110c52818fSSatish Balay   Output Parameter:
9120c52818fSSatish Balay + f - Objective value at X
9130c52818fSSatish Balay - g - Gradient vector at X
9140c52818fSSatish Balay 
9150c52818fSSatish Balay   Notes: TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
9160c52818fSSatish Balay   so most users would not generally call this routine themselves.
9170c52818fSSatish Balay 
9180c52818fSSatish Balay   Level: developer
9190c52818fSSatish Balay 
9200c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
9210c52818fSSatish Balay @*/
9220c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
9230c52818fSSatish Balay {
9240c52818fSSatish Balay   PetscErrorCode ierr;
925f06e3bfaSBarry Smith 
9260c52818fSSatish Balay   PetscFunctionBegin;
9270c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
9280c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
9290c52818fSSatish Balay   PetscValidPointer(f,3);
9300c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
9310c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
9320c52818fSSatish Balay   PetscCheckSameComm(ls,1,g,4);
9330c52818fSSatish Balay   if (ls->usetaoroutines) {
934f06e3bfaSBarry Smith       ierr = TaoComputeObjectiveAndGradient(ls->taosolver,x,f,g);CHKERRQ(ierr);
9350c52818fSSatish Balay   } else {
9360c52818fSSatish Balay     ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
937f06e3bfaSBarry Smith     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
938f06e3bfaSBarry Smith     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");
9390c52818fSSatish Balay 
9400c52818fSSatish Balay     PetscStackPush("TaoLineSearch user objective/gradient routine");
9410c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
9420c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr);
9430c52818fSSatish Balay     } else {
9440c52818fSSatish Balay       ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr);
9450c52818fSSatish Balay       ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr);
9460c52818fSSatish Balay     }
9470c52818fSSatish Balay     PetscStackPop;
9480c52818fSSatish Balay     ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
9490c52818fSSatish Balay     ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",*f);CHKERRQ(ierr);
9500c52818fSSatish Balay     ls->nfgeval++;
951f06e3bfaSBarry Smith   }
9520c52818fSSatish Balay   PetscFunctionReturn(0);
9530c52818fSSatish Balay }
9540c52818fSSatish Balay 
9550c52818fSSatish Balay #undef __FUNCT__
9560c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeGradient"
9570c52818fSSatish Balay /*@
9580c52818fSSatish Balay   TaoLineSearchComputeGradient - Computes the gradient of the objective function
9590c52818fSSatish Balay 
9600c52818fSSatish Balay   Collective on TaoLineSearch
9610c52818fSSatish Balay 
9620c52818fSSatish Balay   Input Parameters:
9630c52818fSSatish Balay + ls - the TaoLineSearch context
9640c52818fSSatish Balay - x - input vector
9650c52818fSSatish Balay 
9660c52818fSSatish Balay   Output Parameter:
9670c52818fSSatish Balay . g - gradient vector
9680c52818fSSatish Balay 
9690c52818fSSatish Balay   Notes: TaoComputeGradient() is typically used within line searches
9700c52818fSSatish Balay   so most users would not generally call this routine themselves.
9710c52818fSSatish Balay 
9720c52818fSSatish Balay   Level: developer
9730c52818fSSatish Balay 
9740c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
9750c52818fSSatish Balay @*/
9760c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
9770c52818fSSatish Balay {
9780c52818fSSatish Balay   PetscErrorCode ierr;
9790c52818fSSatish Balay   PetscReal      fdummy;
980f06e3bfaSBarry Smith 
9810c52818fSSatish Balay   PetscFunctionBegin;
9820c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
9830c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
9840c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,3);
9850c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
9860c52818fSSatish Balay   PetscCheckSameComm(ls,1,g,3);
9870c52818fSSatish Balay   if (ls->usetaoroutines) {
9880c52818fSSatish Balay     ierr = TaoComputeGradient(ls->taosolver,x,g);CHKERRQ(ierr);
9890c52818fSSatish Balay   } else {
9900c52818fSSatish Balay     ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
991f06e3bfaSBarry Smith     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
9920c52818fSSatish Balay     PetscStackPush("TaoLineSearch user gradient routine");
9930c52818fSSatish Balay     if (ls->ops->computegradient) {
9940c52818fSSatish Balay       ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr);
9950c52818fSSatish Balay     } else {
9960c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr);
9970c52818fSSatish Balay     }
9980c52818fSSatish Balay     PetscStackPop;
9990c52818fSSatish Balay     ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
10000c52818fSSatish Balay   }
10010c52818fSSatish Balay   ls->ngeval++;
10020c52818fSSatish Balay   PetscFunctionReturn(0);
10030c52818fSSatish Balay }
10040c52818fSSatish Balay 
10050c52818fSSatish Balay #undef __FUNCT__
10060c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeObjectiveAndGTS"
10070c52818fSSatish Balay /*@
10080c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point
10090c52818fSSatish Balay 
10100c52818fSSatish Balay   Collective on TaoSolver
10110c52818fSSatish Balay 
10120c52818fSSatish Balay   Input Parameters:
10130c52818fSSatish Balay + ls - the TaoLineSearch context
10140c52818fSSatish Balay - x - input vector
10150c52818fSSatish Balay 
10160c52818fSSatish Balay   Output Parameter:
10170c52818fSSatish Balay + f - Objective value at X
10180c52818fSSatish Balay - gts - inner product of gradient and step direction at X
10190c52818fSSatish Balay 
10200c52818fSSatish Balay   Notes: TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
10210c52818fSSatish Balay   so most users would not generally call this routine themselves.
10220c52818fSSatish Balay 
10230c52818fSSatish Balay   Level: developer
10240c52818fSSatish Balay 
10250c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
10260c52818fSSatish Balay @*/
10270c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
10280c52818fSSatish Balay {
10290c52818fSSatish Balay   PetscErrorCode ierr;
10300c52818fSSatish Balay   PetscFunctionBegin;
10310c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
10320c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
10330c52818fSSatish Balay   PetscValidPointer(f,3);
10340c52818fSSatish Balay   PetscValidPointer(gts,4);
10350c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
10360c52818fSSatish Balay   ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
1037f06e3bfaSBarry Smith   if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
10380c52818fSSatish Balay   PetscStackPush("TaoLineSearch user objective/gts routine");
10390c52818fSSatish Balay   ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr);
10400c52818fSSatish Balay   PetscStackPop;
10410c52818fSSatish Balay   ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr);
10420c52818fSSatish Balay   ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",*f);CHKERRQ(ierr);
10430c52818fSSatish Balay   ls->nfeval++;
10440c52818fSSatish Balay   PetscFunctionReturn(0);
10450c52818fSSatish Balay }
10460c52818fSSatish Balay 
10470c52818fSSatish Balay #undef __FUNCT__
10480c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetSolution"
10490c52818fSSatish Balay /*@
10500c52818fSSatish Balay   TaoLineSearchGetSolution - Returns the solution to the line search
10510c52818fSSatish Balay 
10520c52818fSSatish Balay   Collective on TaoLineSearch
10530c52818fSSatish Balay 
10540c52818fSSatish Balay   Input Parameter:
10550c52818fSSatish Balay . ls - the TaoLineSearch context
10560c52818fSSatish Balay 
10570c52818fSSatish Balay   Output Parameter:
10580c52818fSSatish Balay + x - the new solution
10590c52818fSSatish Balay . f - the objective function value at x
10600c52818fSSatish Balay . g - the gradient at x
10610c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search
10620c52818fSSatish Balay - reason - the reason why the line search terminated
10630c52818fSSatish Balay 
10640c52818fSSatish Balay   reason will be set to one of:
10650c52818fSSatish Balay 
10660c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
10670c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
10680c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
10690c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
10700c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
10710c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
10720c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
10730c52818fSSatish Balay 
10740c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
10750c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
10760c52818fSSatish Balay 
10770c52818fSSatish Balay + TAOLINESEARCH_SUCCESS - successful line search
10780c52818fSSatish Balay 
10790c52818fSSatish Balay   Level: developer
10800c52818fSSatish Balay 
10810c52818fSSatish Balay @*/
10820c52818fSSatish Balay PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchTerminationReason *reason)
10830c52818fSSatish Balay {
10840c52818fSSatish Balay   PetscErrorCode ierr;
1085f06e3bfaSBarry Smith 
10860c52818fSSatish Balay   PetscFunctionBegin;
10870c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
10880c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
10890c52818fSSatish Balay   PetscValidPointer(f,3);
10900c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
10910c52818fSSatish Balay   PetscValidIntPointer(reason,6);
10920c52818fSSatish Balay   if (ls->new_x) {
10930c52818fSSatish Balay     ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr);
10940c52818fSSatish Balay   }
10950c52818fSSatish Balay   *f = ls->new_f;
10960c52818fSSatish Balay   if (ls->new_g) {
10970c52818fSSatish Balay     ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr);
10980c52818fSSatish Balay   }
10990c52818fSSatish Balay   if (steplength) {
11000c52818fSSatish Balay     *steplength=ls->step;
11010c52818fSSatish Balay   }
11020c52818fSSatish Balay   *reason = ls->reason;
11030c52818fSSatish Balay   PetscFunctionReturn(0);
11040c52818fSSatish Balay }
11050c52818fSSatish Balay 
11060c52818fSSatish Balay #undef __FUNCT__
11070c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetStartingVector"
11080c52818fSSatish Balay /*@
11090c52818fSSatish Balay   TaoLineSearchGetStartingVector - Gets a the initial point of the line
11100c52818fSSatish Balay   search.
11110c52818fSSatish Balay 
11120c52818fSSatish Balay   Not Collective
11130c52818fSSatish Balay 
11140c52818fSSatish Balay   Input Parameter:
11150c52818fSSatish Balay . ls - the TaoLineSearch context
11160c52818fSSatish Balay 
11170c52818fSSatish Balay   Output Parameter:
11180c52818fSSatish Balay . x - The initial point of the line search
11190c52818fSSatish Balay 
11200c52818fSSatish Balay   Level: intermediate
11210c52818fSSatish Balay @*/
11220c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
11230c52818fSSatish Balay {
11240c52818fSSatish Balay   PetscFunctionBegin;
11250c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11260c52818fSSatish Balay   if (x) {
11270c52818fSSatish Balay     *x = ls->start_x;
11280c52818fSSatish Balay   }
11290c52818fSSatish Balay   PetscFunctionReturn(0);
11300c52818fSSatish Balay }
11310c52818fSSatish Balay 
11320c52818fSSatish Balay #undef __FUNCT__
11330c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetStepDirection"
11340c52818fSSatish Balay /*@
11350c52818fSSatish Balay   TaoLineSearchGetStepDirection - Gets the step direction of the line
11360c52818fSSatish Balay   search.
11370c52818fSSatish Balay 
11380c52818fSSatish Balay   Not Collective
11390c52818fSSatish Balay 
11400c52818fSSatish Balay   Input Parameter:
11410c52818fSSatish Balay . ls - the TaoLineSearch context
11420c52818fSSatish Balay 
11430c52818fSSatish Balay   Output Parameter:
11440c52818fSSatish Balay . s - the step direction of the line search
11450c52818fSSatish Balay 
11460c52818fSSatish Balay   Level: advanced
11470c52818fSSatish Balay @*/
11480c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
11490c52818fSSatish Balay {
11500c52818fSSatish Balay   PetscFunctionBegin;
11510c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11520c52818fSSatish Balay   if (s) {
11530c52818fSSatish Balay     *s = ls->stepdirection;
11540c52818fSSatish Balay   }
11550c52818fSSatish Balay   PetscFunctionReturn(0);
11560c52818fSSatish Balay 
11570c52818fSSatish Balay }
11580c52818fSSatish Balay 
11590c52818fSSatish Balay #undef __FUNCT__
11600c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetFullStepObjective"
11610c52818fSSatish Balay /*@
11620c52818fSSatish Balay   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.
11630c52818fSSatish Balay 
11640c52818fSSatish Balay   Not Collective
11650c52818fSSatish Balay 
11660c52818fSSatish Balay   Input Parameter:
11670c52818fSSatish Balay . ls - the TaoLineSearch context
11680c52818fSSatish Balay 
11690c52818fSSatish Balay   Output Parameter:
11700c52818fSSatish Balay . f - the objective value at the full step length
11710c52818fSSatish Balay 
11720c52818fSSatish Balay   Level: developer
11730c52818fSSatish Balay @*/
11740c52818fSSatish Balay 
11750c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
11760c52818fSSatish Balay {
11770c52818fSSatish Balay   PetscFunctionBegin;
11780c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11790c52818fSSatish Balay   *f_fullstep = ls->f_fullstep;
11800c52818fSSatish Balay   PetscFunctionReturn(0);
11810c52818fSSatish Balay }
11820c52818fSSatish Balay 
11830c52818fSSatish Balay #undef __FUNCT__
11840c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetVariableBounds"
11850c52818fSSatish Balay /*@
11860c52818fSSatish Balay   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.
11870c52818fSSatish Balay 
11880c52818fSSatish Balay   Logically Collective on TaoSolver
11890c52818fSSatish Balay 
11900c52818fSSatish Balay   Input Parameters:
11910c52818fSSatish Balay + ls - the TaoLineSearch context
11920c52818fSSatish Balay . xl  - vector of lower bounds
11930c52818fSSatish Balay - xu  - vector of upper bounds
11940c52818fSSatish Balay 
11950c52818fSSatish Balay   Note: If the variable bounds are not set with this routine, then
1196e270355aSBarry Smith   PETSC_NINFINITY and PETSC_INFINITY are assumed
11970c52818fSSatish Balay 
11980c52818fSSatish Balay   Level: beginner
11990c52818fSSatish Balay 
12000c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
12010c52818fSSatish Balay @*/
12020c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
12030c52818fSSatish Balay {
12040c52818fSSatish Balay   PetscFunctionBegin;
12050c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12060c52818fSSatish Balay   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
12070c52818fSSatish Balay   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
12080c52818fSSatish Balay   ls->lower = xl;
12090c52818fSSatish Balay   ls->upper = xu;
12100c52818fSSatish Balay   ls->bounded = 1;
12110c52818fSSatish Balay   PetscFunctionReturn(0);
12120c52818fSSatish Balay }
12130c52818fSSatish Balay 
12140c52818fSSatish Balay #undef __FUNCT__
12150c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetInitialStepLength"
12160c52818fSSatish Balay /*@
12170c52818fSSatish Balay   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
12180c52818fSSatish Balay   search.  If this value is not set then 1.0 is assumed.
12190c52818fSSatish Balay 
12200c52818fSSatish Balay   Logically Collective on TaoLineSearch
12210c52818fSSatish Balay 
12220c52818fSSatish Balay   Input Parameters:
12230c52818fSSatish Balay + ls - the TaoLineSearch context
12240c52818fSSatish Balay - s - the initial step size
12250c52818fSSatish Balay 
12260c52818fSSatish Balay   Level: intermediate
12270c52818fSSatish Balay 
12280c52818fSSatish Balay .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
12290c52818fSSatish Balay @*/
12300c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
12310c52818fSSatish Balay {
12320c52818fSSatish Balay   PetscFunctionBegin;
12330c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12340c52818fSSatish Balay   ls->initstep = s;
12350c52818fSSatish Balay   PetscFunctionReturn(0);
12360c52818fSSatish Balay }
12370c52818fSSatish Balay 
12380c52818fSSatish Balay #undef __FUNCT__
12390c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetStepLength"
12400c52818fSSatish Balay /*@
12410c52818fSSatish Balay   TaoLineSearchGetStepLength - Get the current step length
12420c52818fSSatish Balay 
12430c52818fSSatish Balay   Not Collective
12440c52818fSSatish Balay 
12450c52818fSSatish Balay   Input Parameters:
12460c52818fSSatish Balay . ls - the TaoLineSearch context
12470c52818fSSatish Balay 
12480c52818fSSatish Balay   Output Parameters
12490c52818fSSatish Balay . s - the current step length
12500c52818fSSatish Balay 
12510c52818fSSatish Balay   Level: beginner
12520c52818fSSatish Balay 
12530c52818fSSatish Balay .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
12540c52818fSSatish Balay @*/
12550c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
12560c52818fSSatish Balay {
12570c52818fSSatish Balay   PetscFunctionBegin;
12580c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12590c52818fSSatish Balay   *s = ls->step;
12600c52818fSSatish Balay   PetscFunctionReturn(0);
12610c52818fSSatish Balay }
12620c52818fSSatish Balay 
12630c52818fSSatish Balay #undef __FUNCT__
12640c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchRegister"
12650c52818fSSatish Balay /*MC
12660c52818fSSatish Balay    TaoLineSearchRegister - Adds a line-search algorithm to the registry
12670c52818fSSatish Balay 
12680c52818fSSatish Balay    Not collective
12690c52818fSSatish Balay 
12700c52818fSSatish Balay    Input Parameters:
12710c52818fSSatish Balay +  sname - name of a new user-defined solver
12720c52818fSSatish Balay -  func - routine to Create method context
12730c52818fSSatish Balay 
12740c52818fSSatish Balay    Notes:
12750c52818fSSatish Balay    TaoLineSearchRegister() may be called multiple times to add several user-defined solvers.
12760c52818fSSatish Balay 
12770c52818fSSatish Balay    Sample usage:
12780c52818fSSatish Balay .vb
12790c52818fSSatish Balay    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
12800c52818fSSatish Balay .ve
12810c52818fSSatish Balay 
12820c52818fSSatish Balay    Then, your solver can be chosen with the procedural interface via
12830c52818fSSatish Balay $     TaoLineSearchSetType(ls,"my_linesearch")
12840c52818fSSatish Balay    or at runtime via the option
12850c52818fSSatish Balay $     -tao_ls_type my_linesearch
12860c52818fSSatish Balay 
12870c52818fSSatish Balay    Level: developer
12880c52818fSSatish Balay 
12890c52818fSSatish Balay .seealso: TaoLineSearchRegisterDestroy()
12900c52818fSSatish Balay M*/
12910c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
12920c52818fSSatish Balay {
12930c52818fSSatish Balay   PetscErrorCode ierr;
12940c52818fSSatish Balay   PetscFunctionBegin;
12950c52818fSSatish Balay   ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr);
12960c52818fSSatish Balay   PetscFunctionReturn(0);
12970c52818fSSatish Balay }
12980c52818fSSatish Balay 
12990c52818fSSatish Balay #undef __FUNCT__
13000c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchRegisterDestroy"
13010c52818fSSatish Balay /*@C
13020c52818fSSatish Balay    TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were
13030c52818fSSatish Balay    registered by TaoLineSearchRegister().
13040c52818fSSatish Balay 
13050c52818fSSatish Balay    Not Collective
13060c52818fSSatish Balay 
13070c52818fSSatish Balay    Level: developer
13080c52818fSSatish Balay 
13090c52818fSSatish Balay .seealso: TaoLineSearchRegister()
13100c52818fSSatish Balay @*/
13110c52818fSSatish Balay PetscErrorCode TaoLineSearchRegisterDestroy(void)
13120c52818fSSatish Balay {
13130c52818fSSatish Balay   PetscErrorCode ierr;
13140c52818fSSatish Balay   PetscFunctionBegin;
13150c52818fSSatish Balay   ierr = PetscFunctionListDestroy(&TaoLineSearchList);CHKERRQ(ierr);
13160c52818fSSatish Balay   TaoLineSearchInitialized = PETSC_FALSE;
13170c52818fSSatish Balay   PetscFunctionReturn(0);
13180c52818fSSatish Balay }
13190c52818fSSatish Balay 
13200c52818fSSatish Balay #undef __FUNCT__
13210c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchAppendOptionsPrefix"
13220c52818fSSatish Balay /*@C
13230c52818fSSatish Balay    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
13240c52818fSSatish Balay    for all TaoLineSearch options in the database.
13250c52818fSSatish Balay 
13260c52818fSSatish Balay 
13270c52818fSSatish Balay    Collective on TaoLineSearch
13280c52818fSSatish Balay 
13290c52818fSSatish Balay    Input Parameters:
13300c52818fSSatish Balay +  ls - the TaoLineSearch solver context
13310c52818fSSatish Balay -  prefix - the prefix string to prepend to all line search requests
13320c52818fSSatish Balay 
13330c52818fSSatish Balay    Notes:
13340c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
13350c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
13360c52818fSSatish Balay 
13370c52818fSSatish Balay 
13380c52818fSSatish Balay    Level: advanced
13390c52818fSSatish Balay 
13400c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
13410c52818fSSatish Balay @*/
13420c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
13430c52818fSSatish Balay {
1344f06e3bfaSBarry Smith   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
13450c52818fSSatish Balay }
13460c52818fSSatish Balay 
13470c52818fSSatish Balay #undef __FUNCT__
13480c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetOptionsPrefix"
13490c52818fSSatish Balay /*@C
13500c52818fSSatish Balay   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
13510c52818fSSatish Balay   TaoLineSearch options in the database
13520c52818fSSatish Balay 
13530c52818fSSatish Balay   Not Collective
13540c52818fSSatish Balay 
13550c52818fSSatish Balay   Input Parameters:
13560c52818fSSatish Balay . ls - the TaoLineSearch context
13570c52818fSSatish Balay 
13580c52818fSSatish Balay   Output Parameters:
13590c52818fSSatish Balay . prefix - pointer to the prefix string used is returned
13600c52818fSSatish Balay 
13610c52818fSSatish Balay   Notes: On the fortran side, the user should pass in a string 'prefix' of
13620c52818fSSatish Balay   sufficient length to hold the prefix.
13630c52818fSSatish Balay 
13640c52818fSSatish Balay   Level: advanced
13650c52818fSSatish Balay 
13660c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
13670c52818fSSatish Balay @*/
13680c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
13690c52818fSSatish Balay {
1370f06e3bfaSBarry Smith   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
13710c52818fSSatish Balay }
13720c52818fSSatish Balay 
13730c52818fSSatish Balay #undef __FUNCT__
13740c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetOptionsPrefix"
13750c52818fSSatish Balay /*@C
13760c52818fSSatish Balay    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
13770c52818fSSatish Balay    TaoLineSearch options in the database.
13780c52818fSSatish Balay 
13790c52818fSSatish Balay 
13800c52818fSSatish Balay    Logically Collective on TaoLineSearch
13810c52818fSSatish Balay 
13820c52818fSSatish Balay    Input Parameters:
13830c52818fSSatish Balay +  ls - the TaoLineSearch context
13840c52818fSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
13850c52818fSSatish Balay 
13860c52818fSSatish Balay    Notes:
13870c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
13880c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
13890c52818fSSatish Balay 
13900c52818fSSatish Balay    For example, to distinguish between the runtime options for two
13910c52818fSSatish Balay    different line searches, one could call
13920c52818fSSatish Balay .vb
13930c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
13940c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
13950c52818fSSatish Balay .ve
13960c52818fSSatish Balay 
13970c52818fSSatish Balay    This would enable use of different options for each system, such as
13980c52818fSSatish Balay .vb
13990c52818fSSatish Balay       -sys1_tao_ls_type mt
14000c52818fSSatish Balay       -sys2_tao_ls_type armijo
14010c52818fSSatish Balay .ve
14020c52818fSSatish Balay 
14030c52818fSSatish Balay    Level: advanced
14040c52818fSSatish Balay 
14050c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
14060c52818fSSatish Balay @*/
14070c52818fSSatish Balay 
14080c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
14090c52818fSSatish Balay {
1410f06e3bfaSBarry Smith   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
14110c52818fSSatish Balay }
1412