xref: /petsc/src/tao/linesearch/interface/taolinesearch.c (revision 6b867d5ac32ed0c728f185df9d084acdf26f70bf)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/
2af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
30c52818fSSatish Balay 
46c23d075SBarry Smith PetscFunctionList TaoLineSearchList = NULL;
50c52818fSSatish Balay 
60c52818fSSatish Balay PetscClassId TAOLINESEARCH_CLASSID=0;
70ebee16dSLisandro Dalcin 
80ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Apply;
90ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Eval;
100c52818fSSatish Balay 
110c52818fSSatish Balay /*@C
12fe2efc57SMark    TaoLineSearchViewFromOptions - View from Options
13fe2efc57SMark 
14fe2efc57SMark    Collective on TaoLineSearch
15fe2efc57SMark 
16fe2efc57SMark    Input Parameters:
17fe2efc57SMark +  A - the Tao context
18736c3998SJose E. Roman .  obj - Optional object
19736c3998SJose E. Roman -  name - command line option
20fe2efc57SMark 
21fe2efc57SMark    Level: intermediate
22fe2efc57SMark .seealso:  TaoLineSearch, TaoLineSearchView, PetscObjectViewFromOptions(), TaoLineSearchCreate()
23fe2efc57SMark @*/
24fe2efc57SMark PetscErrorCode  TaoLineSearchViewFromOptions(TaoLineSearch A,PetscObject obj,const char name[])
25fe2efc57SMark {
26fe2efc57SMark   PetscErrorCode ierr;
27fe2efc57SMark 
28fe2efc57SMark   PetscFunctionBegin;
29fe2efc57SMark   PetscValidHeaderSpecific(A,TAOLINESEARCH_CLASSID,1);
30fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
31fe2efc57SMark   PetscFunctionReturn(0);
32fe2efc57SMark }
33fe2efc57SMark 
34fe2efc57SMark /*@C
350c52818fSSatish Balay   TaoLineSearchView - Prints information about the TaoLineSearch
360c52818fSSatish Balay 
370c52818fSSatish Balay   Collective on TaoLineSearch
380c52818fSSatish Balay 
390c52818fSSatish Balay   InputParameters:
40441846f8SBarry Smith + ls - the Tao context
410c52818fSSatish Balay - viewer - visualization context
420c52818fSSatish Balay 
430c52818fSSatish Balay   Options Database Key:
440c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search
450c52818fSSatish Balay 
460c52818fSSatish Balay   Notes:
470c52818fSSatish Balay   The available visualization contexts include
480c52818fSSatish Balay +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
490c52818fSSatish Balay -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
500c52818fSSatish Balay          output where only the first processor opens
510c52818fSSatish Balay          the file.  All other processors send their
520c52818fSSatish Balay          data to the first processor to print.
530c52818fSSatish Balay 
540c52818fSSatish Balay   Level: beginner
550c52818fSSatish Balay 
560c52818fSSatish Balay .seealso: PetscViewerASCIIOpen()
570c52818fSSatish Balay @*/
580c52818fSSatish Balay 
590c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
600c52818fSSatish Balay {
610c52818fSSatish Balay   PetscErrorCode          ierr;
620c52818fSSatish Balay   PetscBool               isascii, isstring;
63dedfbcbeSJed Brown   TaoLineSearchType       type;
64f06e3bfaSBarry Smith 
650c52818fSSatish Balay   PetscFunctionBegin;
660c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
670c52818fSSatish Balay   if (!viewer) {
680c52818fSSatish Balay     ierr = PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);CHKERRQ(ierr);
690c52818fSSatish Balay   }
700c52818fSSatish Balay   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
710c52818fSSatish Balay   PetscCheckSameComm(ls,1,viewer,2);
720c52818fSSatish Balay 
730c52818fSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
740c52818fSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
750c52818fSSatish Balay   if (isascii) {
7663b15415SAlp Dener     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);CHKERRQ(ierr);
7763b15415SAlp Dener     if (ls->ops->view) {
7863b15415SAlp Dener       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7963b15415SAlp Dener       ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr);
8063b15415SAlp Dener       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
810c52818fSSatish Balay     }
820c52818fSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
830c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr);
84f06e3bfaSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);CHKERRQ(ierr);
850c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr);
860c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr);
870c52818fSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr);
880c52818fSSatish Balay 
890c52818fSSatish Balay     if (ls->bounded) {
900c52818fSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr);
910c52818fSSatish Balay     }
92de196d00SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);CHKERRQ(ierr);
930c52818fSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
940c52818fSSatish Balay   } else if (isstring) {
950c52818fSSatish Balay     ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr);
960c52818fSSatish Balay     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
970c52818fSSatish Balay   }
980c52818fSSatish Balay   PetscFunctionReturn(0);
990c52818fSSatish Balay }
1000c52818fSSatish Balay 
1010c52818fSSatish Balay /*@C
1020c52818fSSatish Balay   TaoLineSearchCreate - Creates a TAO Line Search object.  Algorithms in TAO that use
1030c52818fSSatish Balay   line-searches will automatically create one.
1040c52818fSSatish Balay 
105d083f849SBarry Smith   Collective
1060c52818fSSatish Balay 
1070c52818fSSatish Balay   Input Parameter:
1080c52818fSSatish Balay . comm - MPI communicator
1090c52818fSSatish Balay 
1100c52818fSSatish Balay   Output Parameter:
1110c52818fSSatish Balay . newls - the new TaoLineSearch context
1120c52818fSSatish Balay 
1130c52818fSSatish Balay   Available methods include:
1140c52818fSSatish Balay + more-thuente
1150c52818fSSatish Balay . gpcg
1160c52818fSSatish Balay - unit - Do not perform any line search
1170c52818fSSatish Balay 
1180c52818fSSatish Balay    Options Database Keys:
1190c52818fSSatish Balay .   -tao_ls_type - select which method TAO should use
1200c52818fSSatish Balay 
1210c52818fSSatish Balay    Level: beginner
1220c52818fSSatish Balay 
1230c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
1240c52818fSSatish Balay @*/
1250c52818fSSatish Balay 
1260c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
1270c52818fSSatish Balay {
1280c52818fSSatish Balay   PetscErrorCode ierr;
1290c52818fSSatish Balay   TaoLineSearch  ls;
1300c52818fSSatish Balay 
1310c52818fSSatish Balay   PetscFunctionBegin;
1320c52818fSSatish Balay   PetscValidPointer(newls,2);
1336c23d075SBarry Smith   *newls = NULL;
1340c52818fSSatish Balay 
1350c52818fSSatish Balay   ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr);
1360c52818fSSatish Balay 
13773107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);CHKERRQ(ierr);
1380c52818fSSatish Balay   ls->bounded = 0;
1390c52818fSSatish Balay   ls->max_funcs=30;
1400c52818fSSatish Balay   ls->ftol = 0.0001;
1410c52818fSSatish Balay   ls->gtol = 0.9;
1426f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1436f4723b1SBarry Smith   ls->rtol = 1.0e-5;
1446f4723b1SBarry Smith #else
1450c52818fSSatish Balay   ls->rtol = 1.0e-10;
1466f4723b1SBarry Smith #endif
1470c52818fSSatish Balay   ls->stepmin=1.0e-20;
1480c52818fSSatish Balay   ls->stepmax=1.0e+20;
1490c52818fSSatish Balay   ls->step=1.0;
1500c52818fSSatish Balay   ls->nfeval=0;
1510c52818fSSatish Balay   ls->ngeval=0;
1520c52818fSSatish Balay   ls->nfgeval=0;
1530c52818fSSatish Balay 
15483c8fe1dSLisandro Dalcin   ls->ops->computeobjective = NULL;
15583c8fe1dSLisandro Dalcin   ls->ops->computegradient = NULL;
15683c8fe1dSLisandro Dalcin   ls->ops->computeobjectiveandgradient = NULL;
15783c8fe1dSLisandro Dalcin   ls->ops->computeobjectiveandgts = NULL;
15883c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
15983c8fe1dSLisandro Dalcin   ls->ops->apply = NULL;
16083c8fe1dSLisandro Dalcin   ls->ops->view = NULL;
16183c8fe1dSLisandro Dalcin   ls->ops->setfromoptions = NULL;
16283c8fe1dSLisandro Dalcin   ls->ops->reset = NULL;
16383c8fe1dSLisandro Dalcin   ls->ops->destroy = NULL;
16483c8fe1dSLisandro Dalcin   ls->ops->monitor = NULL;
1652a0dac07SAlp Dener   ls->usemonitor=PETSC_FALSE;
1660c52818fSSatish Balay   ls->setupcalled=PETSC_FALSE;
1670c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
1680c52818fSSatish Balay   *newls = ls;
1690c52818fSSatish Balay   PetscFunctionReturn(0);
1700c52818fSSatish Balay }
1710c52818fSSatish Balay 
1720c52818fSSatish Balay /*@
1730c52818fSSatish Balay   TaoLineSearchSetUp - Sets up the internal data structures for the later use
1740c52818fSSatish Balay   of a Tao solver
1750c52818fSSatish Balay 
1760c52818fSSatish Balay   Collective on ls
1770c52818fSSatish Balay 
1780c52818fSSatish Balay   Input Parameters:
1790c52818fSSatish Balay . ls - the TaoLineSearch context
1800c52818fSSatish Balay 
1810c52818fSSatish Balay   Notes:
1820c52818fSSatish Balay   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
1830c52818fSSatish Balay   automatically be called in TaoLineSearchSolve().  However, if the user
1840c52818fSSatish Balay   desires to call it explicitly, it should come after TaoLineSearchCreate()
1850c52818fSSatish Balay   but before TaoLineSearchApply().
1860c52818fSSatish Balay 
1870c52818fSSatish Balay   Level: developer
1880c52818fSSatish Balay 
1890c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
1900c52818fSSatish Balay @*/
1910c52818fSSatish Balay 
1920c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
1930c52818fSSatish Balay {
1940c52818fSSatish Balay   PetscErrorCode ierr;
1958caf6e8cSBarry Smith   const char     *default_type=TAOLINESEARCHMT;
1960c52818fSSatish Balay   PetscBool      flg;
197f06e3bfaSBarry Smith 
1980c52818fSSatish Balay   PetscFunctionBegin;
1990c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
2000c52818fSSatish Balay   if (ls->setupcalled) PetscFunctionReturn(0);
2010c52818fSSatish Balay   if (!((PetscObject)ls)->type_name) {
2020c52818fSSatish Balay     ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr);
2030c52818fSSatish Balay   }
2040c52818fSSatish Balay   if (ls->ops->setup) {
2050c52818fSSatish Balay     ierr = (*ls->ops->setup)(ls);CHKERRQ(ierr);
2060c52818fSSatish Balay   }
2070c52818fSSatish Balay   if (ls->usetaoroutines) {
208441846f8SBarry Smith     ierr = TaoIsObjectiveDefined(ls->tao,&flg);CHKERRQ(ierr);
2090c52818fSSatish Balay     ls->hasobjective = flg;
210441846f8SBarry Smith     ierr = TaoIsGradientDefined(ls->tao,&flg);CHKERRQ(ierr);
2110c52818fSSatish Balay     ls->hasgradient = flg;
212441846f8SBarry Smith     ierr = TaoIsObjectiveAndGradientDefined(ls->tao,&flg);CHKERRQ(ierr);
2130c52818fSSatish Balay     ls->hasobjectiveandgradient = flg;
2140c52818fSSatish Balay   } else {
2150c52818fSSatish Balay     if (ls->ops->computeobjective) {
2160c52818fSSatish Balay       ls->hasobjective = PETSC_TRUE;
2170c52818fSSatish Balay     } else {
2180c52818fSSatish Balay       ls->hasobjective = PETSC_FALSE;
2190c52818fSSatish Balay     }
2200c52818fSSatish Balay     if (ls->ops->computegradient) {
2210c52818fSSatish Balay       ls->hasgradient = PETSC_TRUE;
2220c52818fSSatish Balay     } else {
2230c52818fSSatish Balay       ls->hasgradient = PETSC_FALSE;
2240c52818fSSatish Balay     }
2250c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
2260c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_TRUE;
2270c52818fSSatish Balay     } else {
2280c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_FALSE;
2290c52818fSSatish Balay     }
2300c52818fSSatish Balay   }
2310c52818fSSatish Balay   ls->setupcalled = PETSC_TRUE;
2320c52818fSSatish Balay   PetscFunctionReturn(0);
2330c52818fSSatish Balay }
2340c52818fSSatish Balay 
2350c52818fSSatish Balay /*@
2360c52818fSSatish Balay   TaoLineSearchReset - Some line searches may carry state information
2370c52818fSSatish Balay   from one TaoLineSearchApply() to the next.  This function resets this
2380c52818fSSatish Balay   state information.
2390c52818fSSatish Balay 
2400c52818fSSatish Balay   Collective on TaoLineSearch
2410c52818fSSatish Balay 
2420c52818fSSatish Balay   Input Parameter:
2430c52818fSSatish Balay . ls - the TaoLineSearch context
2440c52818fSSatish Balay 
2450c52818fSSatish Balay   Level: developer
2460c52818fSSatish Balay 
2470c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
2480c52818fSSatish Balay @*/
2490c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
2500c52818fSSatish Balay {
2510c52818fSSatish Balay   PetscErrorCode ierr;
252050fc7a3SBarry Smith 
2530c52818fSSatish Balay   PetscFunctionBegin;
2540c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
2550c52818fSSatish Balay   if (ls->ops->reset) {
2560c52818fSSatish Balay     ierr = (*ls->ops->reset)(ls);CHKERRQ(ierr);
2570c52818fSSatish Balay   }
2580c52818fSSatish Balay   PetscFunctionReturn(0);
2590c52818fSSatish Balay }
260f06e3bfaSBarry Smith 
2610c52818fSSatish Balay /*@
2620c52818fSSatish Balay   TaoLineSearchDestroy - Destroys the TAO context that was created with
2630c52818fSSatish Balay   TaoLineSearchCreate()
2640c52818fSSatish Balay 
2650c52818fSSatish Balay   Collective on TaoLineSearch
2660c52818fSSatish Balay 
2677a7aea1fSJed Brown   Input Parameter:
2680c52818fSSatish Balay . ls - the TaoLineSearch context
2690c52818fSSatish Balay 
2700c52818fSSatish Balay   Level: beginner
2710c52818fSSatish Balay 
2720c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
2730c52818fSSatish Balay @*/
2740c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
2750c52818fSSatish Balay {
2760c52818fSSatish Balay   PetscErrorCode ierr;
277050fc7a3SBarry Smith 
2780c52818fSSatish Balay   PetscFunctionBegin;
2790c52818fSSatish Balay   if (!*ls) PetscFunctionReturn(0);
2800c52818fSSatish Balay   PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1);
28183c8fe1dSLisandro Dalcin   if (--((PetscObject)*ls)->refct > 0) {*ls = NULL; PetscFunctionReturn(0);}
282050fc7a3SBarry Smith   ierr = VecDestroy(&(*ls)->stepdirection);CHKERRQ(ierr);
283050fc7a3SBarry Smith   ierr = VecDestroy(&(*ls)->start_x);CHKERRQ(ierr);
2840c52818fSSatish Balay   if ((*ls)->ops->destroy) {
2850c52818fSSatish Balay     ierr = (*(*ls)->ops->destroy)(*ls);CHKERRQ(ierr);
2860c52818fSSatish Balay   }
2872a0dac07SAlp Dener   if ((*ls)->usemonitor) {
2882a0dac07SAlp Dener     ierr = PetscViewerDestroy(&(*ls)->viewer);CHKERRQ(ierr);
2892a0dac07SAlp Dener   }
2900c52818fSSatish Balay   ierr = PetscHeaderDestroy(ls);CHKERRQ(ierr);
2910c52818fSSatish Balay   PetscFunctionReturn(0);
2920c52818fSSatish Balay }
2930c52818fSSatish Balay 
2940c52818fSSatish Balay /*@
2950c52818fSSatish Balay   TaoLineSearchApply - Performs a line-search in a given step direction.  Criteria for acceptable step length depends on the line-search algorithm chosen
2960c52818fSSatish Balay 
2970c52818fSSatish Balay   Collective on TaoLineSearch
2980c52818fSSatish Balay 
2990c52818fSSatish Balay   Input Parameters:
300441846f8SBarry Smith + ls - the Tao context
3010c52818fSSatish Balay - s - search direction
3020c52818fSSatish Balay 
303*6b867d5aSJose E. Roman   Input/Output Parameters:
304*6b867d5aSJose E. Roman + x - The current solution (on output x contains the new solution determined by the line search)
305*6b867d5aSJose E. Roman . f - objective function value at current solution (on output contains the objective function value at new solution)
306*6b867d5aSJose E. Roman - g - gradient evaluated at x (on output contains the gradient at new solution)
307*6b867d5aSJose E. Roman 
3080c52818fSSatish Balay   Output Parameters:
309*6b867d5aSJose E. Roman + steplength - scalar multiplier of s used ( x = x0 + steplength * x)
3100c52818fSSatish Balay - reason - reason why the line-search stopped
3110c52818fSSatish Balay 
3120c52818fSSatish Balay   reason will be set to one of:
3130c52818fSSatish Balay 
3140c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
3150c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
3160c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
3170c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
3180c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
3190c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
3200c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
3210c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
3220c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
3230c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search
3240c52818fSSatish Balay 
3250c52818fSSatish Balay   Note:
3260c52818fSSatish Balay   The algorithm developer must set up the TaoLineSearch with calls to
327441846f8SBarry Smith   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()
3280c52818fSSatish Balay 
3290c52818fSSatish Balay   Note:
3300c52818fSSatish Balay   You may or may not need to follow this with a call to
3310c52818fSSatish Balay   TaoAddLineSearchCounts(), depending on whether you want these
3320c52818fSSatish Balay   evaluations to count toward the total function/gradient evaluations.
3330c52818fSSatish Balay 
3340c52818fSSatish Balay   Level: beginner
3350c52818fSSatish Balay 
3360c52818fSSatish Balay   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
3370c52818fSSatish Balay  @*/
3380c52818fSSatish Balay 
339e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
3400c52818fSSatish Balay {
3410c52818fSSatish Balay   PetscErrorCode ierr;
3420c52818fSSatish Balay   PetscInt       low1,low2,low3,high1,high2,high3;
3430c52818fSSatish Balay 
3440c52818fSSatish Balay   PetscFunctionBegin;
3450c52818fSSatish Balay   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
3460c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
3470c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
3483f6ba705SLisandro Dalcin   PetscValidRealPointer(f,3);
3490c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
3500c52818fSSatish Balay   PetscValidHeaderSpecific(s,VEC_CLASSID,5);
3510c52818fSSatish Balay   PetscValidPointer(reason,7);
3520c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
3530c52818fSSatish Balay   PetscCheckSameTypeAndComm(x,2,g,4);
3540c52818fSSatish Balay   PetscCheckSameTypeAndComm(x,2,s,5);
3550c52818fSSatish Balay   ierr = VecGetOwnershipRange(x, &low1, &high1);CHKERRQ(ierr);
3560c52818fSSatish Balay   ierr = VecGetOwnershipRange(g, &low2, &high2);CHKERRQ(ierr);
3570c52818fSSatish Balay   ierr = VecGetOwnershipRange(s, &low3, &high3);CHKERRQ(ierr);
358691b26d3SBarry Smith   if (low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible vector local lengths");
3590c52818fSSatish Balay 
3600c52818fSSatish Balay   ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr);
361050fc7a3SBarry Smith   ierr = VecDestroy(&ls->stepdirection);CHKERRQ(ierr);
362050fc7a3SBarry Smith   ls->stepdirection = s;
3630c52818fSSatish Balay 
3640c52818fSSatish Balay   ierr = TaoLineSearchSetUp(ls);CHKERRQ(ierr);
365691b26d3SBarry Smith   if (!ls->ops->apply) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine");
3660c52818fSSatish Balay   ls->nfeval=0;
3670c52818fSSatish Balay   ls->ngeval=0;
3680c52818fSSatish Balay   ls->nfgeval=0;
3690c52818fSSatish Balay   /* Check parameter values */
3700c52818fSSatish Balay   if (ls->ftol < 0.0) {
371f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);CHKERRQ(ierr);
3720c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3730c52818fSSatish Balay   }
3740c52818fSSatish Balay   if (ls->rtol < 0.0) {
375f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);CHKERRQ(ierr);
3760c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3770c52818fSSatish Balay   }
3780c52818fSSatish Balay   if (ls->gtol < 0.0) {
379f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);CHKERRQ(ierr);
3800c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3810c52818fSSatish Balay   }
3820c52818fSSatish Balay   if (ls->stepmin < 0.0) {
383f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);CHKERRQ(ierr);
3840c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3850c52818fSSatish Balay   }
3860c52818fSSatish Balay   if (ls->stepmax < ls->stepmin) {
387f06e3bfaSBarry Smith     ierr = PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);CHKERRQ(ierr);
3880c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3890c52818fSSatish Balay   }
3900c52818fSSatish Balay   if (ls->max_funcs < 0) {
3910c52818fSSatish Balay     ierr = PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);CHKERRQ(ierr);
3920c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
3930c52818fSSatish Balay   }
3940c52818fSSatish Balay   if (PetscIsInfOrNanReal(*f)) {
395f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);CHKERRQ(ierr);
3960c52818fSSatish Balay     *reason=TAOLINESEARCH_FAILED_INFORNAN;
3970c52818fSSatish Balay   }
3980c52818fSSatish Balay 
399302440fdSBarry Smith   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
4000c52818fSSatish Balay   ierr = VecDestroy(&ls->start_x);CHKERRQ(ierr);
4010c52818fSSatish Balay   ls->start_x = x;
402050fc7a3SBarry Smith 
4030ebee16dSLisandro Dalcin   ierr = PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0);CHKERRQ(ierr);
4040c52818fSSatish Balay   ierr = (*ls->ops->apply)(ls,x,f,g,s);CHKERRQ(ierr);
4050ebee16dSLisandro Dalcin   ierr = PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0);CHKERRQ(ierr);
4060c52818fSSatish Balay   *reason=ls->reason;
4070c52818fSSatish Balay   ls->new_f = *f;
4080c52818fSSatish Balay 
4090c52818fSSatish Balay   if (steplength) {
4100c52818fSSatish Balay     *steplength=ls->step;
4110c52818fSSatish Balay   }
4120c52818fSSatish Balay 
413fbe0838dSJason Sarich   ierr = TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");CHKERRQ(ierr);
4140c52818fSSatish Balay   PetscFunctionReturn(0);
4150c52818fSSatish Balay }
4160c52818fSSatish Balay 
4170c52818fSSatish Balay /*@C
4180c52818fSSatish Balay    TaoLineSearchSetType - Sets the algorithm used in a line search
4190c52818fSSatish Balay 
4200c52818fSSatish Balay    Collective on TaoLineSearch
4210c52818fSSatish Balay 
4220c52818fSSatish Balay    Input Parameters:
4230c52818fSSatish Balay +  ls - the TaoLineSearch context
424820824deSAlp Dener -  type - the TaoLineSearchType selection
4250c52818fSSatish Balay 
4260c52818fSSatish Balay   Available methods include:
427820824deSAlp Dener +  more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition
428820824deSAlp Dener .  armijo - simple backtracking line search enforcing only the sufficient decrease condition
429820824deSAlp Dener -  unit - do not perform a line search and always accept unit step length
4300c52818fSSatish Balay 
4310c52818fSSatish Balay   Options Database Keys:
432820824deSAlp Dener .  -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime
4330c52818fSSatish Balay 
4340c52818fSSatish Balay   Level: beginner
4350c52818fSSatish Balay 
4360c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply()
4370c52818fSSatish Balay 
4380c52818fSSatish Balay @*/
4390c52818fSSatish Balay 
440dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
4410c52818fSSatish Balay {
4420c52818fSSatish Balay   PetscErrorCode ierr;
4430c52818fSSatish Balay   PetscErrorCode (*r)(TaoLineSearch);
4440c52818fSSatish Balay   PetscBool      flg;
4450c52818fSSatish Balay 
4460c52818fSSatish Balay   PetscFunctionBegin;
4470c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
4480c52818fSSatish Balay   PetscValidCharPointer(type,2);
4490c52818fSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg);CHKERRQ(ierr);
4500c52818fSSatish Balay   if (flg) PetscFunctionReturn(0);
4510c52818fSSatish Balay 
4520c52818fSSatish Balay   ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);CHKERRQ(ierr);
453691b26d3SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type);
4540c52818fSSatish Balay   if (ls->ops->destroy) {
4550c52818fSSatish Balay     ierr = (*(ls)->ops->destroy)(ls);CHKERRQ(ierr);
4560c52818fSSatish Balay   }
4570c52818fSSatish Balay   ls->max_funcs=30;
4580c52818fSSatish Balay   ls->ftol = 0.0001;
4590c52818fSSatish Balay   ls->gtol = 0.9;
4606f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
4616f4723b1SBarry Smith   ls->rtol = 1.0e-5;
4626f4723b1SBarry Smith #else
4630c52818fSSatish Balay   ls->rtol = 1.0e-10;
4646f4723b1SBarry Smith #endif
4650c52818fSSatish Balay   ls->stepmin=1.0e-20;
4660c52818fSSatish Balay   ls->stepmax=1.0e+20;
4670c52818fSSatish Balay 
4680c52818fSSatish Balay   ls->nfeval=0;
4690c52818fSSatish Balay   ls->ngeval=0;
4700c52818fSSatish Balay   ls->nfgeval=0;
47183c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
47283c8fe1dSLisandro Dalcin   ls->ops->apply = NULL;
47383c8fe1dSLisandro Dalcin   ls->ops->view = NULL;
47483c8fe1dSLisandro Dalcin   ls->ops->setfromoptions = NULL;
47583c8fe1dSLisandro Dalcin   ls->ops->destroy = NULL;
4760c52818fSSatish Balay   ls->setupcalled = PETSC_FALSE;
4770c52818fSSatish Balay   ierr = (*r)(ls);CHKERRQ(ierr);
4780c52818fSSatish Balay   ierr = PetscObjectChangeTypeName((PetscObject)ls, type);CHKERRQ(ierr);
4790c52818fSSatish Balay   PetscFunctionReturn(0);
4800c52818fSSatish Balay }
4810c52818fSSatish Balay 
4822a0dac07SAlp Dener /*@C
4832a0dac07SAlp Dener   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
4842a0dac07SAlp Dener   iteration number, step length, and function value before calling the implementation
4852a0dac07SAlp Dener   specific monitor.
4862a0dac07SAlp Dener 
4872a0dac07SAlp Dener    Input Parameters:
4882a0dac07SAlp Dener +  ls - the TaoLineSearch context
4892a0dac07SAlp Dener .  its - the current iterate number (>=0)
4902a0dac07SAlp Dener .  f - the current objective function value
4912a0dac07SAlp Dener -  step - the step length
4922a0dac07SAlp Dener 
4932a0dac07SAlp Dener    Options Database Key:
4942a0dac07SAlp Dener .  -tao_ls_monitor - Use the default monitor, which prints statistics to standard output
4952a0dac07SAlp Dener 
4962a0dac07SAlp Dener    Level: developer
4972a0dac07SAlp Dener 
4982a0dac07SAlp Dener @*/
4992a0dac07SAlp Dener PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
5002a0dac07SAlp Dener {
5012a0dac07SAlp Dener   PetscErrorCode ierr;
5022a0dac07SAlp Dener   PetscInt       tabs;
5032a0dac07SAlp Dener 
5042a0dac07SAlp Dener   PetscFunctionBegin;
5052a0dac07SAlp Dener   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
5062a0dac07SAlp Dener   if (ls->usemonitor) {
5072a0dac07SAlp Dener     ierr = PetscViewerASCIIGetTab(ls->viewer, &tabs);CHKERRQ(ierr);
5082a0dac07SAlp Dener     ierr = PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel);CHKERRQ(ierr);
5092a0dac07SAlp Dener     ierr = PetscViewerASCIIPrintf(ls->viewer, "%3D LS", its);CHKERRQ(ierr);
5102a0dac07SAlp Dener     ierr = PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f);CHKERRQ(ierr);
511ad362606SAlp Dener     ierr = PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step);CHKERRQ(ierr);
5122a0dac07SAlp Dener     if (ls->ops->monitor && its > 0) {
5132a0dac07SAlp Dener       ierr = PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3);CHKERRQ(ierr);
5142a0dac07SAlp Dener       ierr = (*ls->ops->monitor)(ls);CHKERRQ(ierr);
5152a0dac07SAlp Dener     }
5162a0dac07SAlp Dener     ierr = PetscViewerASCIISetTab(ls->viewer, tabs);CHKERRQ(ierr);
5172a0dac07SAlp Dener   }
5182a0dac07SAlp Dener   PetscFunctionReturn(0);
5192a0dac07SAlp Dener }
5202a0dac07SAlp Dener 
5210c52818fSSatish Balay /*@
5220c52818fSSatish Balay   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
5230c52818fSSatish Balay   options.
5240c52818fSSatish Balay 
5250c52818fSSatish Balay   Collective on TaoLineSearch
5260c52818fSSatish Balay 
52701d2d390SJose E. Roman   Input Parameter:
5280c52818fSSatish Balay . ls - the TaoLineSearch context
5290c52818fSSatish Balay 
5300c52818fSSatish Balay   Options Database Keys:
5310c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
5320c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease
5330c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition
5340c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step
5350c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed
5360c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed
5370c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
5380c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output
5390c52818fSSatish Balay 
5400c52818fSSatish Balay   Level: beginner
5410c52818fSSatish Balay @*/
5420c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
5430c52818fSSatish Balay {
5440c52818fSSatish Balay   PetscErrorCode ierr;
5458caf6e8cSBarry Smith   const char     *default_type=TAOLINESEARCHMT;
5462a0dac07SAlp Dener   char           type[256],monfilename[PETSC_MAX_PATH_LEN];
5472a0dac07SAlp Dener   PetscViewer    monviewer;
5480c52818fSSatish Balay   PetscBool      flg;
549f06e3bfaSBarry Smith 
5500c52818fSSatish Balay   PetscFunctionBegin;
5510c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
5520c52818fSSatish Balay   ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr);
5530c52818fSSatish Balay   if (((PetscObject)ls)->type_name) {
5540c52818fSSatish Balay     default_type = ((PetscObject)ls)->type_name;
5550c52818fSSatish Balay   }
5560c52818fSSatish Balay   /* Check for type from options */
5570c52818fSSatish Balay   ierr = PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);CHKERRQ(ierr);
5580c52818fSSatish Balay   if (flg) {
5590c52818fSSatish Balay     ierr = TaoLineSearchSetType(ls,type);CHKERRQ(ierr);
5600c52818fSSatish Balay   } else if (!((PetscObject)ls)->type_name) {
561302440fdSBarry Smith     ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr);
5620c52818fSSatish Balay   }
5630c52818fSSatish Balay 
56494ae4db5SBarry Smith   ierr = PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);CHKERRQ(ierr);
56594ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);CHKERRQ(ierr);
56694ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);CHKERRQ(ierr);
56794ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);CHKERRQ(ierr);
56894ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);CHKERRQ(ierr);
56994ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);CHKERRQ(ierr);
570589a23caSBarry Smith   ierr = PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,sizeof(monfilename),&flg);CHKERRQ(ierr);
5712a0dac07SAlp Dener   if (flg) {
5722a0dac07SAlp Dener     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer);CHKERRQ(ierr);
5732a0dac07SAlp Dener     ls->viewer = monviewer;
5742a0dac07SAlp Dener     ls->usemonitor = PETSC_TRUE;
5752a0dac07SAlp Dener   }
5760c52818fSSatish Balay   if (ls->ops->setfromoptions) {
5773a004c28SBarry Smith     ierr = (*ls->ops->setfromoptions)(PetscOptionsObject,ls);CHKERRQ(ierr);
5780c52818fSSatish Balay   }
5790c52818fSSatish Balay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5800c52818fSSatish Balay   PetscFunctionReturn(0);
5810c52818fSSatish Balay }
5820c52818fSSatish Balay 
5830c52818fSSatish Balay /*@C
5840c52818fSSatish Balay   TaoLineSearchGetType - Gets the current line search algorithm
5850c52818fSSatish Balay 
5860c52818fSSatish Balay   Not Collective
5870c52818fSSatish Balay 
5880c52818fSSatish Balay   Input Parameter:
5890c52818fSSatish Balay . ls - the TaoLineSearch context
5900c52818fSSatish Balay 
591fd292e60Sprj-   Output Parameter:
5920c52818fSSatish Balay . type - the line search algorithm in effect
5930c52818fSSatish Balay 
5940c52818fSSatish Balay   Level: developer
5950c52818fSSatish Balay 
5960c52818fSSatish Balay @*/
597dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
5980c52818fSSatish Balay {
5990c52818fSSatish Balay   PetscFunctionBegin;
6000c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6010c52818fSSatish Balay   PetscValidPointer(type,2);
6020c52818fSSatish Balay   *type = ((PetscObject)ls)->type_name;
6030c52818fSSatish Balay   PetscFunctionReturn(0);
6040c52818fSSatish Balay }
6050c52818fSSatish Balay 
6060c52818fSSatish Balay /*@
6070c52818fSSatish Balay   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
6080c52818fSSatish Balay   routines used by the line search in last application (not cumulative).
6090c52818fSSatish Balay 
6100c52818fSSatish Balay   Not Collective
6110c52818fSSatish Balay 
6120c52818fSSatish Balay   Input Parameter:
6130c52818fSSatish Balay . ls - the TaoLineSearch context
6140c52818fSSatish Balay 
6150c52818fSSatish Balay   Output Parameters:
6160c52818fSSatish Balay + nfeval   - number of function evaluations
6170c52818fSSatish Balay . ngeval   - number of gradient evaluations
6180c52818fSSatish Balay - nfgeval  - number of function/gradient evaluations
6190c52818fSSatish Balay 
6200c52818fSSatish Balay   Level: intermediate
6210c52818fSSatish Balay 
6220c52818fSSatish Balay   Note:
623441846f8SBarry Smith   If the line search is using the Tao objective and gradient
624441846f8SBarry Smith   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
6250c52818fSSatish Balay   is already counting the number of evaluations.
6260c52818fSSatish Balay 
6270c52818fSSatish Balay @*/
6280c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
6290c52818fSSatish Balay {
6300c52818fSSatish Balay   PetscFunctionBegin;
6310c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6320c52818fSSatish Balay   *nfeval = ls->nfeval;
6330c52818fSSatish Balay   *ngeval = ls->ngeval;
6340c52818fSSatish Balay   *nfgeval = ls->nfgeval;
6350c52818fSSatish Balay   PetscFunctionReturn(0);
6360c52818fSSatish Balay }
6370c52818fSSatish Balay 
6380c52818fSSatish Balay /*@
639441846f8SBarry Smith   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
640441846f8SBarry Smith   Tao evaluation routines.
6410c52818fSSatish Balay 
6420c52818fSSatish Balay   Not Collective
6430c52818fSSatish Balay 
6440c52818fSSatish Balay   Input Parameter:
6450c52818fSSatish Balay . ls - the TaoLineSearch context
6460c52818fSSatish Balay 
6470c52818fSSatish Balay   Output Parameter:
648441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
6490c52818fSSatish Balay         otherwise PETSC_FALSE
6500c52818fSSatish Balay 
6510c52818fSSatish Balay   Level: developer
6520c52818fSSatish Balay @*/
653441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
6540c52818fSSatish Balay {
6550c52818fSSatish Balay   PetscFunctionBegin;
6560c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
657f06e3bfaSBarry Smith   *flg = ls->usetaoroutines;
6580c52818fSSatish Balay   PetscFunctionReturn(0);
6590c52818fSSatish Balay }
6600c52818fSSatish Balay 
6610c52818fSSatish Balay /*@C
6620c52818fSSatish Balay   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
6630c52818fSSatish Balay 
6640c52818fSSatish Balay   Logically Collective on TaoLineSearch
6650c52818fSSatish Balay 
666d8d19677SJose E. Roman   Input Parameters:
6670c52818fSSatish Balay + ls - the TaoLineSearch context
6680c52818fSSatish Balay . func - the objective function evaluation routine
6690c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6700c52818fSSatish Balay 
6710c52818fSSatish Balay   Calling sequence of func:
6720c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);
6730c52818fSSatish Balay 
6740c52818fSSatish Balay + x - input vector
6750c52818fSSatish Balay . f - function value
6760c52818fSSatish Balay - ctx (optional) user-defined context
6770c52818fSSatish Balay 
6780c52818fSSatish Balay   Level: beginner
6790c52818fSSatish Balay 
6800c52818fSSatish Balay   Note:
6810c52818fSSatish Balay   Use this routine only if you want the line search objective
682441846f8SBarry Smith   evaluation routine to be different from the Tao's objective
6830c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
6840c52818fSSatish Balay   the line search gradient and/or function/gradient routine.
6850c52818fSSatish Balay 
6860c52818fSSatish Balay   Note:
6870c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
6880c52818fSSatish Balay   line search, application programmers should be wary of overriding the
6890c52818fSSatish Balay   default objective routine.
6900c52818fSSatish Balay 
691441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
6920c52818fSSatish Balay @*/
6930c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
6940c52818fSSatish Balay {
6950c52818fSSatish Balay   PetscFunctionBegin;
6960c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6970c52818fSSatish Balay 
6980c52818fSSatish Balay   ls->ops->computeobjective=func;
6990c52818fSSatish Balay   if (ctx) ls->userctx_func=ctx;
7000c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
7010c52818fSSatish Balay   PetscFunctionReturn(0);
7020c52818fSSatish Balay }
7030c52818fSSatish Balay 
7040c52818fSSatish Balay /*@C
7050c52818fSSatish Balay   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
7060c52818fSSatish Balay 
7070c52818fSSatish Balay   Logically Collective on TaoLineSearch
7080c52818fSSatish Balay 
709d8d19677SJose E. Roman   Input Parameters:
7100c52818fSSatish Balay + ls - the TaoLineSearch context
7110c52818fSSatish Balay . func - the gradient evaluation routine
7120c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7130c52818fSSatish Balay 
7140c52818fSSatish Balay   Calling sequence of func:
7150c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);
7160c52818fSSatish Balay 
7170c52818fSSatish Balay + x - input vector
7180c52818fSSatish Balay . g - gradient vector
7190c52818fSSatish Balay - ctx (optional) user-defined context
7200c52818fSSatish Balay 
7210c52818fSSatish Balay   Level: beginner
7220c52818fSSatish Balay 
7230c52818fSSatish Balay   Note:
7240c52818fSSatish Balay   Use this routine only if you want the line search gradient
725441846f8SBarry Smith   evaluation routine to be different from the Tao's gradient
7260c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
7270c52818fSSatish Balay   the line search function and/or function/gradient routine.
7280c52818fSSatish Balay 
7290c52818fSSatish Balay   Note:
7300c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own gradient routine for the
7310c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7320c52818fSSatish Balay   default gradient routine.
7330c52818fSSatish Balay 
734441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
7350c52818fSSatish Balay @*/
7360c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
7370c52818fSSatish Balay {
7380c52818fSSatish Balay   PetscFunctionBegin;
7390c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
7400c52818fSSatish Balay   ls->ops->computegradient=func;
7410c52818fSSatish Balay   if (ctx) ls->userctx_grad=ctx;
7420c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
7430c52818fSSatish Balay   PetscFunctionReturn(0);
7440c52818fSSatish Balay }
7450c52818fSSatish Balay 
7460c52818fSSatish Balay /*@C
7470c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
7480c52818fSSatish Balay 
7490c52818fSSatish Balay   Logically Collective on TaoLineSearch
7500c52818fSSatish Balay 
751d8d19677SJose E. Roman   Input Parameters:
7520c52818fSSatish Balay + ls - the TaoLineSearch context
7530c52818fSSatish Balay . func - the objective and gradient evaluation routine
7540c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7550c52818fSSatish Balay 
7560c52818fSSatish Balay   Calling sequence of func:
7570c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);
7580c52818fSSatish Balay 
7590c52818fSSatish Balay + x - input vector
7600c52818fSSatish Balay . f - function value
7610c52818fSSatish Balay . g - gradient vector
7620c52818fSSatish Balay - ctx (optional) user-defined context
7630c52818fSSatish Balay 
7640c52818fSSatish Balay   Level: beginner
7650c52818fSSatish Balay 
7660c52818fSSatish Balay   Note:
7670c52818fSSatish Balay   Use this routine only if you want the line search objective and gradient
768441846f8SBarry Smith   evaluation routines to be different from the Tao's objective
7690c52818fSSatish Balay   and gradient evaluation routines.
7700c52818fSSatish Balay 
7710c52818fSSatish Balay   Note:
7720c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
7730c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7740c52818fSSatish Balay   default objective routine.
7750c52818fSSatish Balay 
776441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines()
7770c52818fSSatish Balay @*/
7780c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
7790c52818fSSatish Balay {
7800c52818fSSatish Balay   PetscFunctionBegin;
7810c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
7820c52818fSSatish Balay   ls->ops->computeobjectiveandgradient=func;
7830c52818fSSatish Balay   if (ctx) ls->userctx_funcgrad=ctx;
7840c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
7850c52818fSSatish Balay   PetscFunctionReturn(0);
7860c52818fSSatish Balay }
7870c52818fSSatish Balay 
7880c52818fSSatish Balay /*@C
7890c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
7900c52818fSSatish Balay   (gradient'*stepdirection) evaluation routine for the line search.
7910c52818fSSatish Balay   Sometimes it is more efficient to compute the inner product of the gradient
7920c52818fSSatish Balay   and the step direction than it is to compute the gradient, and this is all
7930c52818fSSatish Balay   the line search typically needs of the gradient.
7940c52818fSSatish Balay 
7950c52818fSSatish Balay   Logically Collective on TaoLineSearch
7960c52818fSSatish Balay 
797d8d19677SJose E. Roman   Input Parameters:
7980c52818fSSatish Balay + ls - the TaoLineSearch context
7990c52818fSSatish Balay . func - the objective and gradient evaluation routine
8000c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
8010c52818fSSatish Balay 
8020c52818fSSatish Balay   Calling sequence of func:
8030c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);
8040c52818fSSatish Balay 
8050c52818fSSatish Balay + x - input vector
8060c52818fSSatish Balay . s - step direction
8070c52818fSSatish Balay . f - function value
8080c52818fSSatish Balay . gts - inner product of gradient and step direction vectors
8090c52818fSSatish Balay - ctx (optional) user-defined context
8100c52818fSSatish Balay 
8110c52818fSSatish Balay   Note: The gradient will still need to be computed at the end of the line
8120c52818fSSatish Balay   search, so you will still need to set a line search gradient evaluation
8130c52818fSSatish Balay   routine
8140c52818fSSatish Balay 
8150c52818fSSatish Balay   Note: Bounded line searches (those used in bounded optimization algorithms)
8160c52818fSSatish Balay   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
8170c52818fSSatish Balay   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()
8180c52818fSSatish Balay 
8190c52818fSSatish Balay   Level: advanced
8200c52818fSSatish Balay 
8210c52818fSSatish Balay   Note:
8220c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
8230c52818fSSatish Balay   line search, application programmers should be wary of overriding the
8240c52818fSSatish Balay   default objective routine.
8250c52818fSSatish Balay 
826441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
8270c52818fSSatish Balay @*/
8280c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
8290c52818fSSatish Balay {
8300c52818fSSatish Balay   PetscFunctionBegin;
8310c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8320c52818fSSatish Balay   ls->ops->computeobjectiveandgts=func;
8330c52818fSSatish Balay   if (ctx) ls->userctx_funcgts=ctx;
8340c52818fSSatish Balay   ls->usegts = PETSC_TRUE;
8350c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
8360c52818fSSatish Balay   PetscFunctionReturn(0);
8370c52818fSSatish Balay }
8380c52818fSSatish Balay 
8390c52818fSSatish Balay /*@C
840441846f8SBarry Smith   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
841441846f8SBarry Smith   objective and gradient evaluation routines from the given Tao object.
8420c52818fSSatish Balay 
8430c52818fSSatish Balay   Logically Collective on TaoLineSearch
8440c52818fSSatish Balay 
845d8d19677SJose E. Roman   Input Parameters:
8460c52818fSSatish Balay + ls - the TaoLineSearch context
847441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines
8480c52818fSSatish Balay 
8490c52818fSSatish Balay   Level: developer
8500c52818fSSatish Balay 
8510c52818fSSatish Balay .seealso: TaoLineSearchCreate()
8520c52818fSSatish Balay @*/
853441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
8540c52818fSSatish Balay {
8550c52818fSSatish Balay   PetscFunctionBegin;
8560c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
857064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(ts,TAO_CLASSID,2);
858441846f8SBarry Smith   ls->tao = ts;
8590c52818fSSatish Balay   ls->usetaoroutines=PETSC_TRUE;
8600c52818fSSatish Balay   PetscFunctionReturn(0);
8610c52818fSSatish Balay }
8620c52818fSSatish Balay 
8630c52818fSSatish Balay /*@
8640c52818fSSatish Balay   TaoLineSearchComputeObjective - Computes the objective function value at a given point
8650c52818fSSatish Balay 
8660c52818fSSatish Balay   Collective on TaoLineSearch
8670c52818fSSatish Balay 
8680c52818fSSatish Balay   Input Parameters:
8690c52818fSSatish Balay + ls - the TaoLineSearch context
8700c52818fSSatish Balay - x - input vector
8710c52818fSSatish Balay 
8720c52818fSSatish Balay   Output Parameter:
8730c52818fSSatish Balay . f - Objective value at X
8740c52818fSSatish Balay 
87595452b02SPatrick Sanan   Notes:
87695452b02SPatrick Sanan     TaoLineSearchComputeObjective() is typically used within line searches
8770c52818fSSatish Balay   so most users would not generally call this routine themselves.
8780c52818fSSatish Balay 
8790c52818fSSatish Balay   Level: developer
8800c52818fSSatish Balay 
8810c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
8820c52818fSSatish Balay @*/
8830c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
8840c52818fSSatish Balay {
8850c52818fSSatish Balay   PetscErrorCode ierr;
8860c52818fSSatish Balay   Vec            gdummy;
8870c52818fSSatish Balay   PetscReal      gts;
888f06e3bfaSBarry Smith 
8890c52818fSSatish Balay   PetscFunctionBegin;
8900c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8910c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
8920c52818fSSatish Balay   PetscValidPointer(f,3);
8930c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
8940c52818fSSatish Balay   if (ls->usetaoroutines) {
895441846f8SBarry Smith     ierr = TaoComputeObjective(ls->tao,x,f);CHKERRQ(ierr);
8960c52818fSSatish Balay   } else {
897691b26d3SBarry Smith     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient && !ls->ops->computeobjectiveandgts) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
8980ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
8990c52818fSSatish Balay     PetscStackPush("TaoLineSearch user objective routine");
9000c52818fSSatish Balay     if (ls->ops->computeobjective) {
9010c52818fSSatish Balay       ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr);
9020c52818fSSatish Balay     } else if (ls->ops->computeobjectiveandgradient) {
9030c52818fSSatish Balay       ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr);
9040c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr);
9050c52818fSSatish Balay       ierr = VecDestroy(&gdummy);CHKERRQ(ierr);
9060c52818fSSatish Balay     } else {
9070c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);CHKERRQ(ierr);
9080c52818fSSatish Balay     }
9090c52818fSSatish Balay     PetscStackPop;
9100ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
9110c52818fSSatish Balay   }
9120c52818fSSatish Balay   ls->nfeval++;
9130c52818fSSatish Balay   PetscFunctionReturn(0);
9140c52818fSSatish Balay }
9150c52818fSSatish Balay 
9160c52818fSSatish Balay /*@
9170c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
9180c52818fSSatish Balay 
919441846f8SBarry Smith   Collective on Tao
9200c52818fSSatish Balay 
9210c52818fSSatish Balay   Input Parameters:
9220c52818fSSatish Balay + ls - the TaoLineSearch context
9230c52818fSSatish Balay - x - input vector
9240c52818fSSatish Balay 
925d8d19677SJose E. Roman   Output Parameters:
9260c52818fSSatish Balay + f - Objective value at X
9270c52818fSSatish Balay - g - Gradient vector at X
9280c52818fSSatish Balay 
92995452b02SPatrick Sanan   Notes:
93095452b02SPatrick Sanan     TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
9310c52818fSSatish Balay   so most users would not generally call this routine themselves.
9320c52818fSSatish Balay 
9330c52818fSSatish Balay   Level: developer
9340c52818fSSatish Balay 
9350c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
9360c52818fSSatish Balay @*/
9370c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
9380c52818fSSatish Balay {
9390c52818fSSatish Balay   PetscErrorCode ierr;
940f06e3bfaSBarry Smith 
9410c52818fSSatish Balay   PetscFunctionBegin;
9420c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
9430c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
9440c52818fSSatish Balay   PetscValidPointer(f,3);
9450c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
9460c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
9470c52818fSSatish Balay   PetscCheckSameComm(ls,1,g,4);
9480c52818fSSatish Balay   if (ls->usetaoroutines) {
949441846f8SBarry Smith     ierr = TaoComputeObjectiveAndGradient(ls->tao,x,f,g);CHKERRQ(ierr);
9500c52818fSSatish Balay   } else {
951691b26d3SBarry Smith     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
952691b26d3SBarry Smith     if (!ls->ops->computegradient  && !ls->ops->computeobjectiveandgradient) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");
9530ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
9540c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
955362febeeSStefano Zampini       PetscStackPush("TaoLineSearch user objective/gradient routine");
9560c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr);
9570c52818fSSatish Balay       PetscStackPop;
958362febeeSStefano Zampini     } else {
959362febeeSStefano Zampini       PetscStackPush("TaoLineSearch user objective routine");
960362febeeSStefano Zampini       ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr);
961362febeeSStefano Zampini       PetscStackPop;
962362febeeSStefano Zampini       PetscStackPush("TaoLineSearch user gradient routine");
963362febeeSStefano Zampini       ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr);
964362febeeSStefano Zampini       PetscStackPop;
965362febeeSStefano Zampini     }
9660ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
967335036cbSBarry Smith     ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr);
968f06e3bfaSBarry Smith   }
969fbe0838dSJason Sarich   ls->nfgeval++;
9700c52818fSSatish Balay   PetscFunctionReturn(0);
9710c52818fSSatish Balay }
9720c52818fSSatish Balay 
9730c52818fSSatish Balay /*@
9740c52818fSSatish Balay   TaoLineSearchComputeGradient - Computes the gradient of the objective function
9750c52818fSSatish Balay 
9760c52818fSSatish Balay   Collective on TaoLineSearch
9770c52818fSSatish Balay 
9780c52818fSSatish Balay   Input Parameters:
9790c52818fSSatish Balay + ls - the TaoLineSearch context
9800c52818fSSatish Balay - x - input vector
9810c52818fSSatish Balay 
9820c52818fSSatish Balay   Output Parameter:
9830c52818fSSatish Balay . g - gradient vector
9840c52818fSSatish Balay 
98595452b02SPatrick Sanan   Notes:
98695452b02SPatrick Sanan     TaoComputeGradient() is typically used within line searches
9870c52818fSSatish Balay   so most users would not generally call this routine themselves.
9880c52818fSSatish Balay 
9890c52818fSSatish Balay   Level: developer
9900c52818fSSatish Balay 
9910c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
9920c52818fSSatish Balay @*/
9930c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
9940c52818fSSatish Balay {
9950c52818fSSatish Balay   PetscErrorCode ierr;
9960c52818fSSatish Balay   PetscReal      fdummy;
997f06e3bfaSBarry Smith 
9980c52818fSSatish Balay   PetscFunctionBegin;
9990c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
10000c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
10010c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,3);
10020c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
10030c52818fSSatish Balay   PetscCheckSameComm(ls,1,g,3);
10040c52818fSSatish Balay   if (ls->usetaoroutines) {
1005441846f8SBarry Smith     ierr = TaoComputeGradient(ls->tao,x,g);CHKERRQ(ierr);
10060c52818fSSatish Balay   } else {
1007691b26d3SBarry Smith     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
10080ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
10090c52818fSSatish Balay     PetscStackPush("TaoLineSearch user gradient routine");
10100c52818fSSatish Balay     if (ls->ops->computegradient) {
10110c52818fSSatish Balay       ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr);
10120c52818fSSatish Balay     } else {
10130c52818fSSatish Balay       ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr);
10140c52818fSSatish Balay     }
10150c52818fSSatish Balay     PetscStackPop;
10160ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
10170c52818fSSatish Balay   }
10180c52818fSSatish Balay   ls->ngeval++;
10190c52818fSSatish Balay   PetscFunctionReturn(0);
10200c52818fSSatish Balay }
10210c52818fSSatish Balay 
10220c52818fSSatish Balay /*@
10230c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point
10240c52818fSSatish Balay 
1025441846f8SBarry Smith   Collective on Tao
10260c52818fSSatish Balay 
10270c52818fSSatish Balay   Input Parameters:
10280c52818fSSatish Balay + ls - the TaoLineSearch context
10290c52818fSSatish Balay - x - input vector
10300c52818fSSatish Balay 
1031d8d19677SJose E. Roman   Output Parameters:
10320c52818fSSatish Balay + f - Objective value at X
10330c52818fSSatish Balay - gts - inner product of gradient and step direction at X
10340c52818fSSatish Balay 
103595452b02SPatrick Sanan   Notes:
103695452b02SPatrick Sanan     TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
10370c52818fSSatish Balay   so most users would not generally call this routine themselves.
10380c52818fSSatish Balay 
10390c52818fSSatish Balay   Level: developer
10400c52818fSSatish Balay 
10410c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
10420c52818fSSatish Balay @*/
10430c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
10440c52818fSSatish Balay {
10450c52818fSSatish Balay   PetscErrorCode ierr;
10460c52818fSSatish Balay   PetscFunctionBegin;
10470c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
10480c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
10490c52818fSSatish Balay   PetscValidPointer(f,3);
10500c52818fSSatish Balay   PetscValidPointer(gts,4);
10510c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
1052691b26d3SBarry Smith   if (!ls->ops->computeobjectiveandgts) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
10530ebee16dSLisandro Dalcin   ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
10540c52818fSSatish Balay   PetscStackPush("TaoLineSearch user objective/gts routine");
10550c52818fSSatish Balay   ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr);
10560c52818fSSatish Balay   PetscStackPop;
10570ebee16dSLisandro Dalcin   ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr);
1058335036cbSBarry Smith   ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr);
10590c52818fSSatish Balay   ls->nfeval++;
10600c52818fSSatish Balay   PetscFunctionReturn(0);
10610c52818fSSatish Balay }
10620c52818fSSatish Balay 
10630c52818fSSatish Balay /*@
10640c52818fSSatish Balay   TaoLineSearchGetSolution - Returns the solution to the line search
10650c52818fSSatish Balay 
10660c52818fSSatish Balay   Collective on TaoLineSearch
10670c52818fSSatish Balay 
10680c52818fSSatish Balay   Input Parameter:
10690c52818fSSatish Balay . ls - the TaoLineSearch context
10700c52818fSSatish Balay 
1071d8d19677SJose E. Roman   Output Parameters:
10720c52818fSSatish Balay + x - the new solution
10730c52818fSSatish Balay . f - the objective function value at x
10740c52818fSSatish Balay . g - the gradient at x
10750c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search
10760c52818fSSatish Balay - reason - the reason why the line search terminated
10770c52818fSSatish Balay 
10780c52818fSSatish Balay   reason will be set to one of:
10790c52818fSSatish Balay 
10800c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
10810c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
10820c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
10830c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
10840c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
10850c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
10860c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
10870c52818fSSatish Balay 
10880c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
10890c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
10900c52818fSSatish Balay 
1091a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search
10920c52818fSSatish Balay 
10930c52818fSSatish Balay   Level: developer
10940c52818fSSatish Balay 
10950c52818fSSatish Balay @*/
1096e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
10970c52818fSSatish Balay {
10980c52818fSSatish Balay   PetscErrorCode ierr;
1099f06e3bfaSBarry Smith 
11000c52818fSSatish Balay   PetscFunctionBegin;
11010c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11020c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
11030c52818fSSatish Balay   PetscValidPointer(f,3);
11040c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
11050c52818fSSatish Balay   PetscValidIntPointer(reason,6);
11060c52818fSSatish Balay   if (ls->new_x) {
11070c52818fSSatish Balay     ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr);
11080c52818fSSatish Balay   }
11090c52818fSSatish Balay   *f = ls->new_f;
11100c52818fSSatish Balay   if (ls->new_g) {
11110c52818fSSatish Balay     ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr);
11120c52818fSSatish Balay   }
11130c52818fSSatish Balay   if (steplength) {
11140c52818fSSatish Balay     *steplength=ls->step;
11150c52818fSSatish Balay   }
11160c52818fSSatish Balay   *reason = ls->reason;
11170c52818fSSatish Balay   PetscFunctionReturn(0);
11180c52818fSSatish Balay }
11190c52818fSSatish Balay 
11200c52818fSSatish Balay /*@
11210c52818fSSatish Balay   TaoLineSearchGetStartingVector - Gets a the initial point of the line
11220c52818fSSatish Balay   search.
11230c52818fSSatish Balay 
11240c52818fSSatish Balay   Not Collective
11250c52818fSSatish Balay 
11260c52818fSSatish Balay   Input Parameter:
11270c52818fSSatish Balay . ls - the TaoLineSearch context
11280c52818fSSatish Balay 
11290c52818fSSatish Balay   Output Parameter:
11300c52818fSSatish Balay . x - The initial point of the line search
11310c52818fSSatish Balay 
11320c52818fSSatish Balay   Level: intermediate
11330c52818fSSatish Balay @*/
11340c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
11350c52818fSSatish Balay {
11360c52818fSSatish Balay   PetscFunctionBegin;
11370c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11380c52818fSSatish Balay   if (x) {
11390c52818fSSatish Balay     *x = ls->start_x;
11400c52818fSSatish Balay   }
11410c52818fSSatish Balay   PetscFunctionReturn(0);
11420c52818fSSatish Balay }
11430c52818fSSatish Balay 
11440c52818fSSatish Balay /*@
11450c52818fSSatish Balay   TaoLineSearchGetStepDirection - Gets the step direction of the line
11460c52818fSSatish Balay   search.
11470c52818fSSatish Balay 
11480c52818fSSatish Balay   Not Collective
11490c52818fSSatish Balay 
11500c52818fSSatish Balay   Input Parameter:
11510c52818fSSatish Balay . ls - the TaoLineSearch context
11520c52818fSSatish Balay 
11530c52818fSSatish Balay   Output Parameter:
11540c52818fSSatish Balay . s - the step direction of the line search
11550c52818fSSatish Balay 
11560c52818fSSatish Balay   Level: advanced
11570c52818fSSatish Balay @*/
11580c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
11590c52818fSSatish Balay {
11600c52818fSSatish Balay   PetscFunctionBegin;
11610c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11620c52818fSSatish Balay   if (s) {
11630c52818fSSatish Balay     *s = ls->stepdirection;
11640c52818fSSatish Balay   }
11650c52818fSSatish Balay   PetscFunctionReturn(0);
11660c52818fSSatish Balay 
11670c52818fSSatish Balay }
11680c52818fSSatish Balay 
11690c52818fSSatish Balay /*@
11700c52818fSSatish Balay   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.
11710c52818fSSatish Balay 
11720c52818fSSatish Balay   Not Collective
11730c52818fSSatish Balay 
11740c52818fSSatish Balay   Input Parameter:
11750c52818fSSatish Balay . ls - the TaoLineSearch context
11760c52818fSSatish Balay 
11770c52818fSSatish Balay   Output Parameter:
11780c52818fSSatish Balay . f - the objective value at the full step length
11790c52818fSSatish Balay 
11800c52818fSSatish Balay   Level: developer
11810c52818fSSatish Balay @*/
11820c52818fSSatish Balay 
11830c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
11840c52818fSSatish Balay {
11850c52818fSSatish Balay   PetscFunctionBegin;
11860c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11870c52818fSSatish Balay   *f_fullstep = ls->f_fullstep;
11880c52818fSSatish Balay   PetscFunctionReturn(0);
11890c52818fSSatish Balay }
11900c52818fSSatish Balay 
11910c52818fSSatish Balay /*@
11920c52818fSSatish Balay   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.
11930c52818fSSatish Balay 
1194441846f8SBarry Smith   Logically Collective on Tao
11950c52818fSSatish Balay 
11960c52818fSSatish Balay   Input Parameters:
11970c52818fSSatish Balay + ls - the TaoLineSearch context
11980c52818fSSatish Balay . xl  - vector of lower bounds
11990c52818fSSatish Balay - xu  - vector of upper bounds
12000c52818fSSatish Balay 
12010c52818fSSatish Balay   Note: If the variable bounds are not set with this routine, then
1202e270355aSBarry Smith   PETSC_NINFINITY and PETSC_INFINITY are assumed
12030c52818fSSatish Balay 
12040c52818fSSatish Balay   Level: beginner
12050c52818fSSatish Balay 
12060c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
12070c52818fSSatish Balay @*/
12080c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
12090c52818fSSatish Balay {
12100c52818fSSatish Balay   PetscFunctionBegin;
12110c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12120c52818fSSatish Balay   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
12130c52818fSSatish Balay   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
12140c52818fSSatish Balay   ls->lower = xl;
12150c52818fSSatish Balay   ls->upper = xu;
12160c52818fSSatish Balay   ls->bounded = 1;
12170c52818fSSatish Balay   PetscFunctionReturn(0);
12180c52818fSSatish Balay }
12190c52818fSSatish Balay 
12200c52818fSSatish Balay /*@
12210c52818fSSatish Balay   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
12220c52818fSSatish Balay   search.  If this value is not set then 1.0 is assumed.
12230c52818fSSatish Balay 
12240c52818fSSatish Balay   Logically Collective on TaoLineSearch
12250c52818fSSatish Balay 
12260c52818fSSatish Balay   Input Parameters:
12270c52818fSSatish Balay + ls - the TaoLineSearch context
12280c52818fSSatish Balay - s - the initial step size
12290c52818fSSatish Balay 
12300c52818fSSatish Balay   Level: intermediate
12310c52818fSSatish Balay 
12320c52818fSSatish Balay .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
12330c52818fSSatish Balay @*/
12340c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
12350c52818fSSatish Balay {
12360c52818fSSatish Balay   PetscFunctionBegin;
12370c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12380c52818fSSatish Balay   ls->initstep = s;
12390c52818fSSatish Balay   PetscFunctionReturn(0);
12400c52818fSSatish Balay }
12410c52818fSSatish Balay 
12420c52818fSSatish Balay /*@
12430c52818fSSatish Balay   TaoLineSearchGetStepLength - Get the current step length
12440c52818fSSatish Balay 
12450c52818fSSatish Balay   Not Collective
12460c52818fSSatish Balay 
12470c52818fSSatish Balay   Input Parameters:
12480c52818fSSatish Balay . ls - the TaoLineSearch context
12490c52818fSSatish Balay 
12507a7aea1fSJed Brown   Output Parameters:
12510c52818fSSatish Balay . s - the current step length
12520c52818fSSatish Balay 
12530c52818fSSatish Balay   Level: beginner
12540c52818fSSatish Balay 
12550c52818fSSatish Balay .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
12560c52818fSSatish Balay @*/
12570c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
12580c52818fSSatish Balay {
12590c52818fSSatish Balay   PetscFunctionBegin;
12600c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12610c52818fSSatish Balay   *s = ls->step;
12620c52818fSSatish Balay   PetscFunctionReturn(0);
12630c52818fSSatish Balay }
12640c52818fSSatish Balay 
12650ebee16dSLisandro Dalcin /*@C
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 
12890ebee16dSLisandro Dalcin @*/
12900c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
12910c52818fSSatish Balay {
12920c52818fSSatish Balay   PetscErrorCode ierr;
12930c52818fSSatish Balay   PetscFunctionBegin;
12941d36bdfdSBarry Smith   ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr);
12950c52818fSSatish Balay   ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr);
12960c52818fSSatish Balay   PetscFunctionReturn(0);
12970c52818fSSatish Balay }
12980c52818fSSatish Balay 
12990c52818fSSatish Balay /*@C
13000c52818fSSatish Balay    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
13010c52818fSSatish Balay    for all TaoLineSearch options in the database.
13020c52818fSSatish Balay 
13030c52818fSSatish Balay    Collective on TaoLineSearch
13040c52818fSSatish Balay 
13050c52818fSSatish Balay    Input Parameters:
13060c52818fSSatish Balay +  ls - the TaoLineSearch solver context
13070c52818fSSatish Balay -  prefix - the prefix string to prepend to all line search requests
13080c52818fSSatish Balay 
13090c52818fSSatish Balay    Notes:
13100c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
13110c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
13120c52818fSSatish Balay 
13130c52818fSSatish Balay    Level: advanced
13140c52818fSSatish Balay 
13150c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
13160c52818fSSatish Balay @*/
13170c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
13180c52818fSSatish Balay {
1319f06e3bfaSBarry Smith   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
13200c52818fSSatish Balay }
13210c52818fSSatish Balay 
13220c52818fSSatish Balay /*@C
13230c52818fSSatish Balay   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
13240c52818fSSatish Balay   TaoLineSearch options in the database
13250c52818fSSatish Balay 
13260c52818fSSatish Balay   Not Collective
13270c52818fSSatish Balay 
13280c52818fSSatish Balay   Input Parameters:
13290c52818fSSatish Balay . ls - the TaoLineSearch context
13300c52818fSSatish Balay 
13310c52818fSSatish Balay   Output Parameters:
13320c52818fSSatish Balay . prefix - pointer to the prefix string used is returned
13330c52818fSSatish Balay 
133495452b02SPatrick Sanan   Notes:
133595452b02SPatrick Sanan     On the fortran side, the user should pass in a string 'prefix' of
13360c52818fSSatish Balay   sufficient length to hold the prefix.
13370c52818fSSatish Balay 
13380c52818fSSatish Balay   Level: advanced
13390c52818fSSatish Balay 
13400c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
13410c52818fSSatish Balay @*/
13420c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
13430c52818fSSatish Balay {
1344f06e3bfaSBarry Smith   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
13450c52818fSSatish Balay }
13460c52818fSSatish Balay 
13470c52818fSSatish Balay /*@C
13480c52818fSSatish Balay    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
13490c52818fSSatish Balay    TaoLineSearch options in the database.
13500c52818fSSatish Balay 
13510c52818fSSatish Balay    Logically Collective on TaoLineSearch
13520c52818fSSatish Balay 
13530c52818fSSatish Balay    Input Parameters:
13540c52818fSSatish Balay +  ls - the TaoLineSearch context
13550c52818fSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
13560c52818fSSatish Balay 
13570c52818fSSatish Balay    Notes:
13580c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
13590c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
13600c52818fSSatish Balay 
13610c52818fSSatish Balay    For example, to distinguish between the runtime options for two
13620c52818fSSatish Balay    different line searches, one could call
13630c52818fSSatish Balay .vb
13640c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
13650c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
13660c52818fSSatish Balay .ve
13670c52818fSSatish Balay 
13680c52818fSSatish Balay    This would enable use of different options for each system, such as
13690c52818fSSatish Balay .vb
13700c52818fSSatish Balay       -sys1_tao_ls_type mt
13710c52818fSSatish Balay       -sys2_tao_ls_type armijo
13720c52818fSSatish Balay .ve
13730c52818fSSatish Balay 
13740c52818fSSatish Balay    Level: advanced
13750c52818fSSatish Balay 
13760c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
13770c52818fSSatish Balay @*/
13780c52818fSSatish Balay 
13790c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
13800c52818fSSatish Balay {
1381f06e3bfaSBarry Smith   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
13820c52818fSSatish Balay }
1383