xref: /petsc/src/tao/linesearch/interface/taolinesearch.c (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
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
22db781477SPatrick Sanan .seealso: `TaoLineSearch`, `TaoLineSearchView`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()`
23fe2efc57SMark @*/
24fe2efc57SMark PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A,PetscObject obj,const char name[])
25fe2efc57SMark {
26fe2efc57SMark   PetscFunctionBegin;
27fe2efc57SMark   PetscValidHeaderSpecific(A,TAOLINESEARCH_CLASSID,1);
289566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name));
29fe2efc57SMark   PetscFunctionReturn(0);
30fe2efc57SMark }
31fe2efc57SMark 
32fe2efc57SMark /*@C
330c52818fSSatish Balay   TaoLineSearchView - Prints information about the TaoLineSearch
340c52818fSSatish Balay 
350c52818fSSatish Balay   Collective on TaoLineSearch
360c52818fSSatish Balay 
370c52818fSSatish Balay   InputParameters:
38441846f8SBarry Smith + ls - the Tao context
390c52818fSSatish Balay - viewer - visualization context
400c52818fSSatish Balay 
410c52818fSSatish Balay   Options Database Key:
420c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search
430c52818fSSatish Balay 
440c52818fSSatish Balay   Notes:
450c52818fSSatish Balay   The available visualization contexts include
460c52818fSSatish Balay +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
470c52818fSSatish Balay -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
480c52818fSSatish Balay          output where only the first processor opens
490c52818fSSatish Balay          the file.  All other processors send their
500c52818fSSatish Balay          data to the first processor to print.
510c52818fSSatish Balay 
520c52818fSSatish Balay   Level: beginner
530c52818fSSatish Balay 
54db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`
550c52818fSSatish Balay @*/
560c52818fSSatish Balay 
570c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
580c52818fSSatish Balay {
590c52818fSSatish Balay   PetscBool               isascii, isstring;
60dedfbcbeSJed Brown   TaoLineSearchType       type;
61f06e3bfaSBarry Smith 
620c52818fSSatish Balay   PetscFunctionBegin;
630c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
640c52818fSSatish Balay   if (!viewer) {
659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer));
660c52818fSSatish Balay   }
670c52818fSSatish Balay   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
680c52818fSSatish Balay   PetscCheckSameComm(ls,1,viewer,2);
690c52818fSSatish Balay 
709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
720c52818fSSatish Balay   if (isascii) {
739566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer));
7463b15415SAlp Dener     if (ls->ops->view) {
759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
769566063dSJacob Faibussowitsch       PetscCall((*ls->ops->view)(ls,viewer));
779566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
780c52818fSSatish Balay     }
799566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
8063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%" PetscInt_FMT "\n",ls->max_funcs));
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol));
8263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%" PetscInt_FMT "\n",ls->nfeval));
8363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%" PetscInt_FMT "\n",ls->ngeval));
8463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%" PetscInt_FMT "\n",ls->nfgeval));
850c52818fSSatish Balay 
860c52818fSSatish Balay     if (ls->bounded) {
879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"using variable bounds\n"));
880c52818fSSatish Balay     }
899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason));
909566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
910c52818fSSatish Balay   } else if (isstring) {
929566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchGetType(ls,&type));
939566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer," %-3.3s",type));
940c52818fSSatish Balay   }
950c52818fSSatish Balay   PetscFunctionReturn(0);
960c52818fSSatish Balay }
970c52818fSSatish Balay 
980c52818fSSatish Balay /*@C
990c52818fSSatish Balay   TaoLineSearchCreate - Creates a TAO Line Search object.  Algorithms in TAO that use
1000c52818fSSatish Balay   line-searches will automatically create one.
1010c52818fSSatish Balay 
102d083f849SBarry Smith   Collective
1030c52818fSSatish Balay 
1040c52818fSSatish Balay   Input Parameter:
1050c52818fSSatish Balay . comm - MPI communicator
1060c52818fSSatish Balay 
1070c52818fSSatish Balay   Output Parameter:
1080c52818fSSatish Balay . newls - the new TaoLineSearch context
1090c52818fSSatish Balay 
1100c52818fSSatish Balay   Available methods include:
111147403d9SBarry Smith + more-thuente - the More-Thuente method
112147403d9SBarry Smith . gpcg - the GPCG method
1130c52818fSSatish Balay - unit - Do not perform any line search
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 
120db781477SPatrick Sanan .seealso: `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()`
1210c52818fSSatish Balay @*/
1220c52818fSSatish Balay 
1230c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
1240c52818fSSatish Balay {
1250c52818fSSatish Balay   TaoLineSearch  ls;
1260c52818fSSatish Balay 
1270c52818fSSatish Balay   PetscFunctionBegin;
1280c52818fSSatish Balay   PetscValidPointer(newls,2);
1299566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchInitializePackage());
1300c52818fSSatish Balay 
1319566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView));
1320c52818fSSatish Balay   ls->max_funcs = 30;
1330c52818fSSatish Balay   ls->ftol = 0.0001;
1340c52818fSSatish Balay   ls->gtol = 0.9;
1356f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1366f4723b1SBarry Smith   ls->rtol = 1.0e-5;
1376f4723b1SBarry Smith #else
1380c52818fSSatish Balay   ls->rtol = 1.0e-10;
1396f4723b1SBarry Smith #endif
1400c52818fSSatish Balay   ls->stepmin = 1.0e-20;
1410c52818fSSatish Balay   ls->stepmax = 1.0e+20;
1420c52818fSSatish Balay   ls->step = 1.0;
143a39c8e28SStefano Zampini   ls->initstep = 1.0;
1440c52818fSSatish Balay   *newls = ls;
1450c52818fSSatish Balay   PetscFunctionReturn(0);
1460c52818fSSatish Balay }
1470c52818fSSatish Balay 
1480c52818fSSatish Balay /*@
1490c52818fSSatish Balay   TaoLineSearchSetUp - Sets up the internal data structures for the later use
1500c52818fSSatish Balay   of a Tao solver
1510c52818fSSatish Balay 
1520c52818fSSatish Balay   Collective on ls
1530c52818fSSatish Balay 
1540c52818fSSatish Balay   Input Parameters:
1550c52818fSSatish Balay . ls - the TaoLineSearch context
1560c52818fSSatish Balay 
1570c52818fSSatish Balay   Notes:
1580c52818fSSatish Balay   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
1590c52818fSSatish Balay   automatically be called in TaoLineSearchSolve().  However, if the user
1600c52818fSSatish Balay   desires to call it explicitly, it should come after TaoLineSearchCreate()
1610c52818fSSatish Balay   but before TaoLineSearchApply().
1620c52818fSSatish Balay 
1630c52818fSSatish Balay   Level: developer
1640c52818fSSatish Balay 
165db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()`
1660c52818fSSatish Balay @*/
1670c52818fSSatish Balay 
1680c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
1690c52818fSSatish Balay {
1708caf6e8cSBarry Smith   const char     *default_type=TAOLINESEARCHMT;
1710c52818fSSatish Balay   PetscBool      flg;
172f06e3bfaSBarry Smith 
1730c52818fSSatish Balay   PetscFunctionBegin;
1740c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1750c52818fSSatish Balay   if (ls->setupcalled) PetscFunctionReturn(0);
1760c52818fSSatish Balay   if (!((PetscObject)ls)->type_name) {
1779566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(ls,default_type));
1780c52818fSSatish Balay   }
1790c52818fSSatish Balay   if (ls->ops->setup) {
1809566063dSJacob Faibussowitsch     PetscCall((*ls->ops->setup)(ls));
1810c52818fSSatish Balay   }
1820c52818fSSatish Balay   if (ls->usetaoroutines) {
1839566063dSJacob Faibussowitsch     PetscCall(TaoIsObjectiveDefined(ls->tao,&flg));
1840c52818fSSatish Balay     ls->hasobjective = flg;
1859566063dSJacob Faibussowitsch     PetscCall(TaoIsGradientDefined(ls->tao,&flg));
1860c52818fSSatish Balay     ls->hasgradient = flg;
1879566063dSJacob Faibussowitsch     PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao,&flg));
1880c52818fSSatish Balay     ls->hasobjectiveandgradient = flg;
1890c52818fSSatish Balay   } else {
1900c52818fSSatish Balay     if (ls->ops->computeobjective) {
1910c52818fSSatish Balay       ls->hasobjective = PETSC_TRUE;
1920c52818fSSatish Balay     } else {
1930c52818fSSatish Balay       ls->hasobjective = PETSC_FALSE;
1940c52818fSSatish Balay     }
1950c52818fSSatish Balay     if (ls->ops->computegradient) {
1960c52818fSSatish Balay       ls->hasgradient = PETSC_TRUE;
1970c52818fSSatish Balay     } else {
1980c52818fSSatish Balay       ls->hasgradient = PETSC_FALSE;
1990c52818fSSatish Balay     }
2000c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
2010c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_TRUE;
2020c52818fSSatish Balay     } else {
2030c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_FALSE;
2040c52818fSSatish Balay     }
2050c52818fSSatish Balay   }
2060c52818fSSatish Balay   ls->setupcalled = PETSC_TRUE;
2070c52818fSSatish Balay   PetscFunctionReturn(0);
2080c52818fSSatish Balay }
2090c52818fSSatish Balay 
2100c52818fSSatish Balay /*@
2110c52818fSSatish Balay   TaoLineSearchReset - Some line searches may carry state information
2120c52818fSSatish Balay   from one TaoLineSearchApply() to the next.  This function resets this
2130c52818fSSatish Balay   state information.
2140c52818fSSatish Balay 
2150c52818fSSatish Balay   Collective on TaoLineSearch
2160c52818fSSatish Balay 
2170c52818fSSatish Balay   Input Parameter:
2180c52818fSSatish Balay . ls - the TaoLineSearch context
2190c52818fSSatish Balay 
2200c52818fSSatish Balay   Level: developer
2210c52818fSSatish Balay 
222db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()`
2230c52818fSSatish Balay @*/
2240c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
2250c52818fSSatish Balay {
2260c52818fSSatish Balay   PetscFunctionBegin;
2270c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
2280c52818fSSatish Balay   if (ls->ops->reset) {
2299566063dSJacob Faibussowitsch     PetscCall((*ls->ops->reset)(ls));
2300c52818fSSatish Balay   }
2310c52818fSSatish Balay   PetscFunctionReturn(0);
2320c52818fSSatish Balay }
233f06e3bfaSBarry Smith 
2340c52818fSSatish Balay /*@
2350c52818fSSatish Balay   TaoLineSearchDestroy - Destroys the TAO context that was created with
2360c52818fSSatish Balay   TaoLineSearchCreate()
2370c52818fSSatish Balay 
2380c52818fSSatish Balay   Collective on TaoLineSearch
2390c52818fSSatish Balay 
2407a7aea1fSJed Brown   Input Parameter:
2410c52818fSSatish Balay . ls - the TaoLineSearch context
2420c52818fSSatish Balay 
2430c52818fSSatish Balay   Level: beginner
2440c52818fSSatish Balay 
2450c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
2460c52818fSSatish Balay @*/
2470c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
2480c52818fSSatish Balay {
2490c52818fSSatish Balay   PetscFunctionBegin;
2500c52818fSSatish Balay   if (!*ls) PetscFunctionReturn(0);
2510c52818fSSatish Balay   PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1);
25283c8fe1dSLisandro Dalcin   if (--((PetscObject)*ls)->refct > 0) {*ls = NULL; PetscFunctionReturn(0);}
2539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->stepdirection));
2549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->start_x));
2559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->upper));
2569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->lower));
2570c52818fSSatish Balay   if ((*ls)->ops->destroy) {
2589566063dSJacob Faibussowitsch     PetscCall((*(*ls)->ops->destroy)(*ls));
2590c52818fSSatish Balay   }
2602a0dac07SAlp Dener   if ((*ls)->usemonitor) {
2619566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&(*ls)->viewer));
2622a0dac07SAlp Dener   }
2639566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(ls));
2640c52818fSSatish Balay   PetscFunctionReturn(0);
2650c52818fSSatish Balay }
2660c52818fSSatish Balay 
2670c52818fSSatish Balay /*@
2680c52818fSSatish Balay   TaoLineSearchApply - Performs a line-search in a given step direction.  Criteria for acceptable step length depends on the line-search algorithm chosen
2690c52818fSSatish Balay 
2700c52818fSSatish Balay   Collective on TaoLineSearch
2710c52818fSSatish Balay 
2720c52818fSSatish Balay   Input Parameters:
273441846f8SBarry Smith + ls - the Tao context
2740c52818fSSatish Balay - s - search direction
2750c52818fSSatish Balay 
2766b867d5aSJose E. Roman   Input/Output Parameters:
2776b867d5aSJose E. Roman 
2780c52818fSSatish Balay   Output Parameters:
279f1a722f8SMatthew G. Knepley + x - On input the current solution, on output x contains the new solution determined by the line search
280f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution
281f1a722f8SMatthew G. Knepley . g - On input the gradient evaluated at x, on output contains the gradient at new solution
282f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x)
2830c52818fSSatish Balay - reason - reason why the line-search stopped
2840c52818fSSatish Balay 
28597bb3fdcSJose E. Roman   Notes:
2860c52818fSSatish Balay   reason will be set to one of:
2870c52818fSSatish Balay 
2880c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
2890c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
2900c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
2910c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
2920c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
2930c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
2940c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
2950c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
2960c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
2970c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search
2980c52818fSSatish Balay 
2990c52818fSSatish Balay   The algorithm developer must set up the TaoLineSearch with calls to
300441846f8SBarry Smith   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()
3010c52818fSSatish Balay 
3020c52818fSSatish Balay   You may or may not need to follow this with a call to
3030c52818fSSatish Balay   TaoAddLineSearchCounts(), depending on whether you want these
3040c52818fSSatish Balay   evaluations to count toward the total function/gradient evaluations.
3050c52818fSSatish Balay 
3060c52818fSSatish Balay   Level: beginner
3070c52818fSSatish Balay 
308db781477SPatrick Sanan   .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()`
3090c52818fSSatish Balay  @*/
3100c52818fSSatish Balay 
311e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
3120c52818fSSatish Balay {
3130c52818fSSatish Balay   PetscInt       low1,low2,low3,high1,high2,high3;
3140c52818fSSatish Balay 
3150c52818fSSatish Balay   PetscFunctionBegin;
3160c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
3170c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
3183f6ba705SLisandro Dalcin   PetscValidRealPointer(f,3);
3190c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
3200c52818fSSatish Balay   PetscValidHeaderSpecific(s,VEC_CLASSID,5);
3210c52818fSSatish Balay   PetscValidPointer(reason,7);
3220c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
3230c52818fSSatish Balay   PetscCheckSameTypeAndComm(x,2,g,4);
3240c52818fSSatish Balay   PetscCheckSameTypeAndComm(x,2,s,5);
3259566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &low1, &high1));
3269566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(g, &low2, &high2));
3279566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(s, &low3, &high3));
3283c859ba3SBarry Smith   PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible vector local lengths");
3290c52818fSSatish Balay 
33097ab8e72SStefano Zampini   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
3319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
3329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->stepdirection));
333050fc7a3SBarry Smith   ls->stepdirection = s;
3340c52818fSSatish Balay 
3359566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetUp(ls));
3363c859ba3SBarry Smith   PetscCheck(ls->ops->apply,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine");
3370c52818fSSatish Balay   ls->nfeval = 0;
3380c52818fSSatish Balay   ls->ngeval = 0;
3390c52818fSSatish Balay   ls->nfgeval = 0;
3400c52818fSSatish Balay   /* Check parameter values */
3410c52818fSSatish Balay   if (ls->ftol < 0.0) {
3429566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol));
3430c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3440c52818fSSatish Balay   }
3450c52818fSSatish Balay   if (ls->rtol < 0.0) {
3469566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol));
3470c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3480c52818fSSatish Balay   }
3490c52818fSSatish Balay   if (ls->gtol < 0.0) {
3509566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol));
3510c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3520c52818fSSatish Balay   }
3530c52818fSSatish Balay   if (ls->stepmin < 0.0) {
3549566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin));
3550c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3560c52818fSSatish Balay   }
3570c52818fSSatish Balay   if (ls->stepmax < ls->stepmin) {
3589566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax));
3590c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3600c52818fSSatish Balay   }
3610c52818fSSatish Balay   if (ls->max_funcs < 0) {
3629566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n",ls->max_funcs));
3630c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3640c52818fSSatish Balay   }
3650c52818fSSatish Balay   if (PetscIsInfOrNanReal(*f)) {
3669566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f));
3670c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_INFORNAN;
3680c52818fSSatish Balay   }
3690c52818fSSatish Balay 
3709566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)x));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->start_x));
3720c52818fSSatish Balay   ls->start_x = x;
373050fc7a3SBarry Smith 
3749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0));
3759566063dSJacob Faibussowitsch   PetscCall((*ls->ops->apply)(ls,x,f,g,s));
3769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0));
3770c52818fSSatish Balay   *reason = ls->reason;
3780c52818fSSatish Balay   ls->new_f = *f;
3790c52818fSSatish Balay 
38097ab8e72SStefano Zampini   if (steplength) *steplength = ls->step;
3810c52818fSSatish Balay 
3829566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view"));
3830c52818fSSatish Balay   PetscFunctionReturn(0);
3840c52818fSSatish Balay }
3850c52818fSSatish Balay 
3860c52818fSSatish Balay /*@C
3870c52818fSSatish Balay    TaoLineSearchSetType - Sets the algorithm used in a line search
3880c52818fSSatish Balay 
3890c52818fSSatish Balay    Collective on TaoLineSearch
3900c52818fSSatish Balay 
3910c52818fSSatish Balay    Input Parameters:
3920c52818fSSatish Balay +  ls - the TaoLineSearch context
393820824deSAlp Dener -  type - the TaoLineSearchType selection
3940c52818fSSatish Balay 
3950c52818fSSatish Balay   Available methods include:
396820824deSAlp Dener +  more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition
397820824deSAlp Dener .  armijo - simple backtracking line search enforcing only the sufficient decrease condition
398820824deSAlp Dener -  unit - do not perform a line search and always accept unit step length
3990c52818fSSatish Balay 
4000c52818fSSatish Balay   Options Database Keys:
401820824deSAlp Dener .  -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime
4020c52818fSSatish Balay 
4030c52818fSSatish Balay   Level: beginner
4040c52818fSSatish Balay 
405db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchGetType()`, `TaoLineSearchApply()`
4060c52818fSSatish Balay 
4070c52818fSSatish Balay @*/
4080c52818fSSatish Balay 
409dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
4100c52818fSSatish Balay {
4110c52818fSSatish Balay   PetscErrorCode (*r)(TaoLineSearch);
4120c52818fSSatish Balay   PetscBool      flg;
4130c52818fSSatish Balay 
4140c52818fSSatish Balay   PetscFunctionBegin;
4150c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
4160c52818fSSatish Balay   PetscValidCharPointer(type,2);
4179566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg));
4180c52818fSSatish Balay   if (flg) PetscFunctionReturn(0);
4190c52818fSSatish Balay 
4209566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r));
4213c859ba3SBarry Smith   PetscCheck(r,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type);
4220c52818fSSatish Balay   if (ls->ops->destroy) {
4239566063dSJacob Faibussowitsch     PetscCall((*(ls)->ops->destroy)(ls));
4240c52818fSSatish Balay   }
4250c52818fSSatish Balay   ls->max_funcs = 30;
4260c52818fSSatish Balay   ls->ftol = 0.0001;
4270c52818fSSatish Balay   ls->gtol = 0.9;
4286f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
4296f4723b1SBarry Smith   ls->rtol = 1.0e-5;
4306f4723b1SBarry Smith #else
4310c52818fSSatish Balay   ls->rtol = 1.0e-10;
4326f4723b1SBarry Smith #endif
4330c52818fSSatish Balay   ls->stepmin = 1.0e-20;
4340c52818fSSatish Balay   ls->stepmax = 1.0e+20;
4350c52818fSSatish Balay 
4360c52818fSSatish Balay   ls->nfeval = 0;
4370c52818fSSatish Balay   ls->ngeval = 0;
4380c52818fSSatish Balay   ls->nfgeval = 0;
43983c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
44083c8fe1dSLisandro Dalcin   ls->ops->apply = NULL;
44183c8fe1dSLisandro Dalcin   ls->ops->view = NULL;
44283c8fe1dSLisandro Dalcin   ls->ops->setfromoptions = NULL;
44383c8fe1dSLisandro Dalcin   ls->ops->destroy = NULL;
4440c52818fSSatish Balay   ls->setupcalled = PETSC_FALSE;
4459566063dSJacob Faibussowitsch   PetscCall((*r)(ls));
4469566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type));
4470c52818fSSatish Balay   PetscFunctionReturn(0);
4480c52818fSSatish Balay }
4490c52818fSSatish Balay 
4502a0dac07SAlp Dener /*@C
4512a0dac07SAlp Dener   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
4522a0dac07SAlp Dener   iteration number, step length, and function value before calling the implementation
4532a0dac07SAlp Dener   specific monitor.
4542a0dac07SAlp Dener 
4552a0dac07SAlp Dener    Input Parameters:
4562a0dac07SAlp Dener +  ls - the TaoLineSearch context
4572a0dac07SAlp Dener .  its - the current iterate number (>=0)
4582a0dac07SAlp Dener .  f - the current objective function value
4592a0dac07SAlp Dener -  step - the step length
4602a0dac07SAlp Dener 
4612a0dac07SAlp Dener    Options Database Key:
4622a0dac07SAlp Dener .  -tao_ls_monitor - Use the default monitor, which prints statistics to standard output
4632a0dac07SAlp Dener 
4642a0dac07SAlp Dener    Level: developer
4652a0dac07SAlp Dener 
4662a0dac07SAlp Dener @*/
4672a0dac07SAlp Dener PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
4682a0dac07SAlp Dener {
4692a0dac07SAlp Dener   PetscInt       tabs;
4702a0dac07SAlp Dener 
4712a0dac07SAlp Dener   PetscFunctionBegin;
4722a0dac07SAlp Dener   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
4732a0dac07SAlp Dener   if (ls->usemonitor) {
4749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs));
4759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel));
47663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its));
4779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f));
4789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step));
4792a0dac07SAlp Dener     if (ls->ops->monitor && its > 0) {
4809566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3));
4819566063dSJacob Faibussowitsch       PetscCall((*ls->ops->monitor)(ls));
4822a0dac07SAlp Dener     }
4839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs));
4842a0dac07SAlp Dener   }
4852a0dac07SAlp Dener   PetscFunctionReturn(0);
4862a0dac07SAlp Dener }
4872a0dac07SAlp Dener 
4880c52818fSSatish Balay /*@
4890c52818fSSatish Balay   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
4900c52818fSSatish Balay   options.
4910c52818fSSatish Balay 
4920c52818fSSatish Balay   Collective on TaoLineSearch
4930c52818fSSatish Balay 
49401d2d390SJose E. Roman   Input Parameter:
4950c52818fSSatish Balay . ls - the TaoLineSearch context
4960c52818fSSatish Balay 
4970c52818fSSatish Balay   Options Database Keys:
4980c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
4990c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease
5000c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition
5010c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step
502a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - initial steplength allowed
5030c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed
5040c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed
5050c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
5060c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output
5070c52818fSSatish Balay 
5080c52818fSSatish Balay   Level: beginner
5090c52818fSSatish Balay @*/
5100c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
5110c52818fSSatish Balay {
5128caf6e8cSBarry Smith   const char     *default_type=TAOLINESEARCHMT;
5132a0dac07SAlp Dener   char           type[256],monfilename[PETSC_MAX_PATH_LEN];
5142a0dac07SAlp Dener   PetscViewer    monviewer;
5150c52818fSSatish Balay   PetscBool      flg;
516f06e3bfaSBarry Smith 
5170c52818fSSatish Balay   PetscFunctionBegin;
5180c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
519d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)ls);
5200c52818fSSatish Balay   if (((PetscObject)ls)->type_name) {
5210c52818fSSatish Balay     default_type = ((PetscObject)ls)->type_name;
5220c52818fSSatish Balay   }
5230c52818fSSatish Balay   /* Check for type from options */
5249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg));
5250c52818fSSatish Balay   if (flg) {
5269566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(ls,type));
5270c52818fSSatish Balay   } else if (!((PetscObject)ls)->type_name) {
5289566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(ls,default_type));
5290c52818fSSatish Balay   }
5300c52818fSSatish Balay 
5319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL));
5329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL));
5339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL));
5349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL));
5359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL));
5369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL));
537a39c8e28SStefano Zampini   PetscCall(PetscOptionsReal("-tao_ls_stepinit","initial step","",ls->initstep,&ls->initstep,NULL));
5389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
5392a0dac07SAlp Dener   if (flg) {
5409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer));
5412a0dac07SAlp Dener     ls->viewer = monviewer;
5422a0dac07SAlp Dener     ls->usemonitor = PETSC_TRUE;
5432a0dac07SAlp Dener   }
5440c52818fSSatish Balay   if (ls->ops->setfromoptions) {
5459566063dSJacob Faibussowitsch     PetscCall((*ls->ops->setfromoptions)(PetscOptionsObject,ls));
5460c52818fSSatish Balay   }
547d0609cedSBarry Smith   PetscOptionsEnd();
5480c52818fSSatish Balay   PetscFunctionReturn(0);
5490c52818fSSatish Balay }
5500c52818fSSatish Balay 
5510c52818fSSatish Balay /*@C
5520c52818fSSatish Balay   TaoLineSearchGetType - Gets the current line search algorithm
5530c52818fSSatish Balay 
5540c52818fSSatish Balay   Not Collective
5550c52818fSSatish Balay 
5560c52818fSSatish Balay   Input Parameter:
5570c52818fSSatish Balay . ls - the TaoLineSearch context
5580c52818fSSatish Balay 
559fd292e60Sprj-   Output Parameter:
5600c52818fSSatish Balay . type - the line search algorithm in effect
5610c52818fSSatish Balay 
5620c52818fSSatish Balay   Level: developer
5630c52818fSSatish Balay 
5640c52818fSSatish Balay @*/
565dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
5660c52818fSSatish Balay {
5670c52818fSSatish Balay   PetscFunctionBegin;
5680c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
5690c52818fSSatish Balay   PetscValidPointer(type,2);
5700c52818fSSatish Balay   *type = ((PetscObject)ls)->type_name;
5710c52818fSSatish Balay   PetscFunctionReturn(0);
5720c52818fSSatish Balay }
5730c52818fSSatish Balay 
5740c52818fSSatish Balay /*@
5750c52818fSSatish Balay   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
5760c52818fSSatish Balay   routines used by the line search in last application (not cumulative).
5770c52818fSSatish Balay 
5780c52818fSSatish Balay   Not Collective
5790c52818fSSatish Balay 
5800c52818fSSatish Balay   Input Parameter:
5810c52818fSSatish Balay . ls - the TaoLineSearch context
5820c52818fSSatish Balay 
5830c52818fSSatish Balay   Output Parameters:
5840c52818fSSatish Balay + nfeval   - number of function evaluations
5850c52818fSSatish Balay . ngeval   - number of gradient evaluations
5860c52818fSSatish Balay - nfgeval  - number of function/gradient evaluations
5870c52818fSSatish Balay 
5880c52818fSSatish Balay   Level: intermediate
5890c52818fSSatish Balay 
5900c52818fSSatish Balay   Note:
591441846f8SBarry Smith   If the line search is using the Tao objective and gradient
592441846f8SBarry Smith   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
5930c52818fSSatish Balay   is already counting the number of evaluations.
5940c52818fSSatish Balay 
5950c52818fSSatish Balay @*/
5960c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
5970c52818fSSatish Balay {
5980c52818fSSatish Balay   PetscFunctionBegin;
5990c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6000c52818fSSatish Balay   *nfeval = ls->nfeval;
6010c52818fSSatish Balay   *ngeval = ls->ngeval;
6020c52818fSSatish Balay   *nfgeval = ls->nfgeval;
6030c52818fSSatish Balay   PetscFunctionReturn(0);
6040c52818fSSatish Balay }
6050c52818fSSatish Balay 
6060c52818fSSatish Balay /*@
607441846f8SBarry Smith   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
608441846f8SBarry Smith   Tao evaluation routines.
6090c52818fSSatish Balay 
6100c52818fSSatish Balay   Not Collective
6110c52818fSSatish Balay 
6120c52818fSSatish Balay   Input Parameter:
6130c52818fSSatish Balay . ls - the TaoLineSearch context
6140c52818fSSatish Balay 
6150c52818fSSatish Balay   Output Parameter:
616441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
6170c52818fSSatish Balay         otherwise PETSC_FALSE
6180c52818fSSatish Balay 
6190c52818fSSatish Balay   Level: developer
6200c52818fSSatish Balay @*/
621441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
6220c52818fSSatish Balay {
6230c52818fSSatish Balay   PetscFunctionBegin;
6240c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
625f06e3bfaSBarry Smith   *flg = ls->usetaoroutines;
6260c52818fSSatish Balay   PetscFunctionReturn(0);
6270c52818fSSatish Balay }
6280c52818fSSatish Balay 
6290c52818fSSatish Balay /*@C
6300c52818fSSatish Balay   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
6310c52818fSSatish Balay 
6320c52818fSSatish Balay   Logically Collective on TaoLineSearch
6330c52818fSSatish Balay 
634d8d19677SJose E. Roman   Input Parameters:
6350c52818fSSatish Balay + ls - the TaoLineSearch context
6360c52818fSSatish Balay . func - the objective function evaluation routine
6370c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6380c52818fSSatish Balay 
6390c52818fSSatish Balay   Calling sequence of func:
6400c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);
6410c52818fSSatish Balay 
6420c52818fSSatish Balay + x - input vector
6430c52818fSSatish Balay . f - function value
6440c52818fSSatish Balay - ctx (optional) user-defined context
6450c52818fSSatish Balay 
6460c52818fSSatish Balay   Level: beginner
6470c52818fSSatish Balay 
6480c52818fSSatish Balay   Note:
6490c52818fSSatish Balay   Use this routine only if you want the line search objective
650441846f8SBarry Smith   evaluation routine to be different from the Tao's objective
6510c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
6520c52818fSSatish Balay   the line search gradient and/or function/gradient routine.
6530c52818fSSatish Balay 
6540c52818fSSatish Balay   Note:
6550c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
6560c52818fSSatish Balay   line search, application programmers should be wary of overriding the
6570c52818fSSatish Balay   default objective routine.
6580c52818fSSatish Balay 
659db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
6600c52818fSSatish Balay @*/
6610c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
6620c52818fSSatish Balay {
6630c52818fSSatish Balay   PetscFunctionBegin;
6640c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
6650c52818fSSatish Balay 
6660c52818fSSatish Balay   ls->ops->computeobjective=func;
6670c52818fSSatish Balay   if (ctx) ls->userctx_func=ctx;
6680c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
6690c52818fSSatish Balay   PetscFunctionReturn(0);
6700c52818fSSatish Balay }
6710c52818fSSatish Balay 
6720c52818fSSatish Balay /*@C
6730c52818fSSatish Balay   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
6740c52818fSSatish Balay 
6750c52818fSSatish Balay   Logically Collective on TaoLineSearch
6760c52818fSSatish Balay 
677d8d19677SJose E. Roman   Input Parameters:
6780c52818fSSatish Balay + ls - the TaoLineSearch context
6790c52818fSSatish Balay . func - the gradient evaluation routine
6800c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6810c52818fSSatish Balay 
6820c52818fSSatish Balay   Calling sequence of func:
6830c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);
6840c52818fSSatish Balay 
6850c52818fSSatish Balay + x - input vector
6860c52818fSSatish Balay . g - gradient vector
6870c52818fSSatish Balay - ctx (optional) user-defined context
6880c52818fSSatish Balay 
6890c52818fSSatish Balay   Level: beginner
6900c52818fSSatish Balay 
6910c52818fSSatish Balay   Note:
6920c52818fSSatish Balay   Use this routine only if you want the line search gradient
693441846f8SBarry Smith   evaluation routine to be different from the Tao's gradient
6940c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
6950c52818fSSatish Balay   the line search function and/or function/gradient routine.
6960c52818fSSatish Balay 
6970c52818fSSatish Balay   Note:
6980c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own gradient routine for the
6990c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7000c52818fSSatish Balay   default gradient routine.
7010c52818fSSatish Balay 
702db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
7030c52818fSSatish Balay @*/
7040c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
7050c52818fSSatish Balay {
7060c52818fSSatish Balay   PetscFunctionBegin;
7070c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
7080c52818fSSatish Balay   ls->ops->computegradient=func;
7090c52818fSSatish Balay   if (ctx) ls->userctx_grad=ctx;
7100c52818fSSatish Balay   ls->usetaoroutines=PETSC_FALSE;
7110c52818fSSatish Balay   PetscFunctionReturn(0);
7120c52818fSSatish Balay }
7130c52818fSSatish Balay 
7140c52818fSSatish Balay /*@C
7150c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
7160c52818fSSatish Balay 
7170c52818fSSatish Balay   Logically Collective on TaoLineSearch
7180c52818fSSatish Balay 
719d8d19677SJose E. Roman   Input Parameters:
7200c52818fSSatish Balay + ls - the TaoLineSearch context
7210c52818fSSatish Balay . func - the objective and gradient evaluation routine
7220c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7230c52818fSSatish Balay 
7240c52818fSSatish Balay   Calling sequence of func:
7250c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);
7260c52818fSSatish Balay 
7270c52818fSSatish Balay + x - input vector
7280c52818fSSatish Balay . f - function value
7290c52818fSSatish Balay . g - gradient vector
7300c52818fSSatish Balay - ctx (optional) user-defined context
7310c52818fSSatish Balay 
7320c52818fSSatish Balay   Level: beginner
7330c52818fSSatish Balay 
7340c52818fSSatish Balay   Note:
7350c52818fSSatish Balay   Use this routine only if you want the line search objective and gradient
736441846f8SBarry Smith   evaluation routines to be different from the Tao's objective
7370c52818fSSatish Balay   and gradient evaluation routines.
7380c52818fSSatish Balay 
7390c52818fSSatish Balay   Note:
7400c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
7410c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7420c52818fSSatish Balay   default objective routine.
7430c52818fSSatish Balay 
744db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
7450c52818fSSatish Balay @*/
7460c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
7470c52818fSSatish Balay {
7480c52818fSSatish Balay   PetscFunctionBegin;
7490c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
7500c52818fSSatish Balay   ls->ops->computeobjectiveandgradient=func;
7510c52818fSSatish Balay   if (ctx) ls->userctx_funcgrad=ctx;
7520c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
7530c52818fSSatish Balay   PetscFunctionReturn(0);
7540c52818fSSatish Balay }
7550c52818fSSatish Balay 
7560c52818fSSatish Balay /*@C
7570c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
7580c52818fSSatish Balay   (gradient'*stepdirection) evaluation routine for the line search.
7590c52818fSSatish Balay   Sometimes it is more efficient to compute the inner product of the gradient
7600c52818fSSatish Balay   and the step direction than it is to compute the gradient, and this is all
7610c52818fSSatish Balay   the line search typically needs of the gradient.
7620c52818fSSatish Balay 
7630c52818fSSatish Balay   Logically Collective on TaoLineSearch
7640c52818fSSatish Balay 
765d8d19677SJose E. Roman   Input Parameters:
7660c52818fSSatish Balay + ls - the TaoLineSearch context
7670c52818fSSatish Balay . func - the objective and gradient evaluation routine
7680c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7690c52818fSSatish Balay 
7700c52818fSSatish Balay   Calling sequence of func:
7710c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);
7720c52818fSSatish Balay 
7730c52818fSSatish Balay + x - input vector
7740c52818fSSatish Balay . s - step direction
7750c52818fSSatish Balay . f - function value
7760c52818fSSatish Balay . gts - inner product of gradient and step direction vectors
7770c52818fSSatish Balay - ctx (optional) user-defined context
7780c52818fSSatish Balay 
7790c52818fSSatish Balay   Note: The gradient will still need to be computed at the end of the line
7800c52818fSSatish Balay   search, so you will still need to set a line search gradient evaluation
7810c52818fSSatish Balay   routine
7820c52818fSSatish Balay 
7830c52818fSSatish Balay   Note: Bounded line searches (those used in bounded optimization algorithms)
7840c52818fSSatish Balay   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
7850c52818fSSatish Balay   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()
7860c52818fSSatish Balay 
7870c52818fSSatish Balay   Level: advanced
7880c52818fSSatish Balay 
7890c52818fSSatish Balay   Note:
7900c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
7910c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7920c52818fSSatish Balay   default objective routine.
7930c52818fSSatish Balay 
794db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()`
7950c52818fSSatish Balay @*/
7960c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
7970c52818fSSatish Balay {
7980c52818fSSatish Balay   PetscFunctionBegin;
7990c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8000c52818fSSatish Balay   ls->ops->computeobjectiveandgts=func;
8010c52818fSSatish Balay   if (ctx) ls->userctx_funcgts=ctx;
8020c52818fSSatish Balay   ls->usegts = PETSC_TRUE;
8030c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
8040c52818fSSatish Balay   PetscFunctionReturn(0);
8050c52818fSSatish Balay }
8060c52818fSSatish Balay 
8070c52818fSSatish Balay /*@C
808441846f8SBarry Smith   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
809441846f8SBarry Smith   objective and gradient evaluation routines from the given Tao object.
8100c52818fSSatish Balay 
8110c52818fSSatish Balay   Logically Collective on TaoLineSearch
8120c52818fSSatish Balay 
813d8d19677SJose E. Roman   Input Parameters:
8140c52818fSSatish Balay + ls - the TaoLineSearch context
815441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines
8160c52818fSSatish Balay 
8170c52818fSSatish Balay   Level: developer
8180c52818fSSatish Balay 
819db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`
8200c52818fSSatish Balay @*/
821441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
8220c52818fSSatish Balay {
8230c52818fSSatish Balay   PetscFunctionBegin;
8240c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
825064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(ts,TAO_CLASSID,2);
826441846f8SBarry Smith   ls->tao = ts;
8270c52818fSSatish Balay   ls->usetaoroutines = PETSC_TRUE;
8280c52818fSSatish Balay   PetscFunctionReturn(0);
8290c52818fSSatish Balay }
8300c52818fSSatish Balay 
8310c52818fSSatish Balay /*@
8320c52818fSSatish Balay   TaoLineSearchComputeObjective - Computes the objective function value at a given point
8330c52818fSSatish Balay 
8340c52818fSSatish Balay   Collective on TaoLineSearch
8350c52818fSSatish Balay 
8360c52818fSSatish Balay   Input Parameters:
8370c52818fSSatish Balay + ls - the TaoLineSearch context
8380c52818fSSatish Balay - x - input vector
8390c52818fSSatish Balay 
8400c52818fSSatish Balay   Output Parameter:
8410c52818fSSatish Balay . f - Objective value at X
8420c52818fSSatish Balay 
84395452b02SPatrick Sanan   Notes:
84495452b02SPatrick Sanan     TaoLineSearchComputeObjective() is typically used within line searches
8450c52818fSSatish Balay   so most users would not generally call this routine themselves.
8460c52818fSSatish Balay 
8470c52818fSSatish Balay   Level: developer
8480c52818fSSatish Balay 
849db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
8500c52818fSSatish Balay @*/
8510c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
8520c52818fSSatish Balay {
8530c52818fSSatish Balay   Vec            gdummy;
8540c52818fSSatish Balay   PetscReal      gts;
855f06e3bfaSBarry Smith 
8560c52818fSSatish Balay   PetscFunctionBegin;
8570c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
8580c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
859dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f,3);
8600c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
8610c52818fSSatish Balay   if (ls->usetaoroutines) {
8629566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(ls->tao,x,f));
8630c52818fSSatish Balay   } else {
8643c859ba3SBarry Smith     PetscCheck(ls->ops->computeobjective || ls->ops->computeobjectiveandgradient || ls->ops->computeobjectiveandgts,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
8659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
8660c52818fSSatish Balay     PetscStackPush("TaoLineSearch user objective routine");
8670c52818fSSatish Balay     if (ls->ops->computeobjective) {
8689566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computeobjective)(ls,x,f,ls->userctx_func));
8690c52818fSSatish Balay     } else if (ls->ops->computeobjectiveandgradient) {
8709566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x,&gdummy));
8719566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad));
8729566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&gdummy));
8730c52818fSSatish Balay     } else {
8749566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts));
8750c52818fSSatish Balay     }
8760c52818fSSatish Balay     PetscStackPop;
8779566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
8780c52818fSSatish Balay   }
8790c52818fSSatish Balay   ls->nfeval++;
8800c52818fSSatish Balay   PetscFunctionReturn(0);
8810c52818fSSatish Balay }
8820c52818fSSatish Balay 
8830c52818fSSatish Balay /*@
8840c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
8850c52818fSSatish Balay 
886441846f8SBarry Smith   Collective on Tao
8870c52818fSSatish Balay 
8880c52818fSSatish Balay   Input Parameters:
8890c52818fSSatish Balay + ls - the TaoLineSearch context
8900c52818fSSatish Balay - x - input vector
8910c52818fSSatish Balay 
892d8d19677SJose E. Roman   Output Parameters:
8930c52818fSSatish Balay + f - Objective value at X
8940c52818fSSatish Balay - g - Gradient vector at X
8950c52818fSSatish Balay 
89695452b02SPatrick Sanan   Notes:
89795452b02SPatrick Sanan     TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
8980c52818fSSatish Balay   so most users would not generally call this routine themselves.
8990c52818fSSatish Balay 
9000c52818fSSatish Balay   Level: developer
9010c52818fSSatish Balay 
902db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
9030c52818fSSatish Balay @*/
9040c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
9050c52818fSSatish Balay {
9060c52818fSSatish Balay   PetscFunctionBegin;
9070c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
9080c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
909dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f,3);
9100c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
9110c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
9120c52818fSSatish Balay   PetscCheckSameComm(ls,1,g,4);
9130c52818fSSatish Balay   if (ls->usetaoroutines) {
9149566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(ls->tao,x,f,g));
9150c52818fSSatish Balay   } else {
9169566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
9170c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
918362febeeSStefano Zampini       PetscStackPush("TaoLineSearch user objective/gradient routine");
9199566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad));
9200c52818fSSatish Balay       PetscStackPop;
921362febeeSStefano Zampini     } else {
9223c859ba3SBarry Smith       PetscCheck(ls->ops->computeobjective,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
923362febeeSStefano Zampini       PetscStackPush("TaoLineSearch user objective routine");
9249566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computeobjective)(ls,x,f,ls->userctx_func));
925362febeeSStefano Zampini       PetscStackPop;
9263c859ba3SBarry Smith       PetscCheck(ls->ops->computegradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");
927362febeeSStefano Zampini       PetscStackPush("TaoLineSearch user gradient routine");
9289566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computegradient)(ls,x,g,ls->userctx_grad));
929362febeeSStefano Zampini       PetscStackPop;
930362febeeSStefano Zampini     }
9319566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
9329566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f)));
933f06e3bfaSBarry Smith   }
934fbe0838dSJason Sarich   ls->nfgeval++;
9350c52818fSSatish Balay   PetscFunctionReturn(0);
9360c52818fSSatish Balay }
9370c52818fSSatish Balay 
9380c52818fSSatish Balay /*@
9390c52818fSSatish Balay   TaoLineSearchComputeGradient - Computes the gradient of the objective function
9400c52818fSSatish Balay 
9410c52818fSSatish Balay   Collective on TaoLineSearch
9420c52818fSSatish Balay 
9430c52818fSSatish Balay   Input Parameters:
9440c52818fSSatish Balay + ls - the TaoLineSearch context
9450c52818fSSatish Balay - x - input vector
9460c52818fSSatish Balay 
9470c52818fSSatish Balay   Output Parameter:
9480c52818fSSatish Balay . g - gradient vector
9490c52818fSSatish Balay 
95095452b02SPatrick Sanan   Notes:
95195452b02SPatrick Sanan     TaoComputeGradient() is typically used within line searches
9520c52818fSSatish Balay   so most users would not generally call this routine themselves.
9530c52818fSSatish Balay 
9540c52818fSSatish Balay   Level: developer
9550c52818fSSatish Balay 
956db781477SPatrick Sanan .seealso: `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()`
9570c52818fSSatish Balay @*/
9580c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
9590c52818fSSatish Balay {
9600c52818fSSatish Balay   PetscReal      fdummy;
961f06e3bfaSBarry Smith 
9620c52818fSSatish Balay   PetscFunctionBegin;
9630c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
9640c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
9650c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,3);
9660c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
9670c52818fSSatish Balay   PetscCheckSameComm(ls,1,g,3);
9680c52818fSSatish Balay   if (ls->usetaoroutines) {
9699566063dSJacob Faibussowitsch     PetscCall(TaoComputeGradient(ls->tao,x,g));
9700c52818fSSatish Balay   } else {
9719566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
9720c52818fSSatish Balay     PetscStackPush("TaoLineSearch user gradient routine");
9730c52818fSSatish Balay     if (ls->ops->computegradient) {
9749566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computegradient)(ls,x,g,ls->userctx_grad));
9750c52818fSSatish Balay     } else {
9763c859ba3SBarry Smith       PetscCheck(ls->ops->computeobjectiveandgradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
9779566063dSJacob Faibussowitsch       PetscCall((*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad));
9780c52818fSSatish Balay     }
9790c52818fSSatish Balay     PetscStackPop;
9809566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
9810c52818fSSatish Balay   }
9820c52818fSSatish Balay   ls->ngeval++;
9830c52818fSSatish Balay   PetscFunctionReturn(0);
9840c52818fSSatish Balay }
9850c52818fSSatish Balay 
9860c52818fSSatish Balay /*@
9870c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point
9880c52818fSSatish Balay 
989441846f8SBarry Smith   Collective on Tao
9900c52818fSSatish Balay 
9910c52818fSSatish Balay   Input Parameters:
9920c52818fSSatish Balay + ls - the TaoLineSearch context
9930c52818fSSatish Balay - x - input vector
9940c52818fSSatish Balay 
995d8d19677SJose E. Roman   Output Parameters:
9960c52818fSSatish Balay + f - Objective value at X
9970c52818fSSatish Balay - gts - inner product of gradient and step direction at X
9980c52818fSSatish Balay 
99995452b02SPatrick Sanan   Notes:
100095452b02SPatrick Sanan     TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
10010c52818fSSatish Balay   so most users would not generally call this routine themselves.
10020c52818fSSatish Balay 
10030c52818fSSatish Balay   Level: developer
10040c52818fSSatish Balay 
1005db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
10060c52818fSSatish Balay @*/
10070c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
10080c52818fSSatish Balay {
10090c52818fSSatish Balay   PetscFunctionBegin;
10100c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
10110c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1012dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f,3);
1013dadcf809SJacob Faibussowitsch   PetscValidRealPointer(gts,4);
10140c52818fSSatish Balay   PetscCheckSameComm(ls,1,x,2);
10153c859ba3SBarry Smith   PetscCheck(ls->ops->computeobjectiveandgts,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
10169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0));
10170c52818fSSatish Balay   PetscStackPush("TaoLineSearch user objective/gts routine");
10189566063dSJacob Faibussowitsch   PetscCall((*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts));
10190c52818fSSatish Balay   PetscStackPop;
10209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0));
10219566063dSJacob Faibussowitsch   PetscCall(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f)));
10220c52818fSSatish Balay   ls->nfeval++;
10230c52818fSSatish Balay   PetscFunctionReturn(0);
10240c52818fSSatish Balay }
10250c52818fSSatish Balay 
10260c52818fSSatish Balay /*@
10270c52818fSSatish Balay   TaoLineSearchGetSolution - Returns the solution to the line search
10280c52818fSSatish Balay 
10290c52818fSSatish Balay   Collective on TaoLineSearch
10300c52818fSSatish Balay 
10310c52818fSSatish Balay   Input Parameter:
10320c52818fSSatish Balay . ls - the TaoLineSearch context
10330c52818fSSatish Balay 
1034d8d19677SJose E. Roman   Output Parameters:
10350c52818fSSatish Balay + x - the new solution
10360c52818fSSatish Balay . f - the objective function value at x
10370c52818fSSatish Balay . g - the gradient at x
10380c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search
10390c52818fSSatish Balay - reason - the reason why the line search terminated
10400c52818fSSatish Balay 
10410c52818fSSatish Balay   reason will be set to one of:
10420c52818fSSatish Balay 
10430c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
10440c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
10450c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
10460c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
10470c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
10480c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
10490c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
10500c52818fSSatish Balay 
10510c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
10520c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
10530c52818fSSatish Balay 
1054a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search
10550c52818fSSatish Balay 
10560c52818fSSatish Balay   Level: developer
10570c52818fSSatish Balay 
10580c52818fSSatish Balay @*/
1059e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
10600c52818fSSatish Balay {
10610c52818fSSatish Balay   PetscFunctionBegin;
10620c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
10630c52818fSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1064dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f,3);
10650c52818fSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
10660c52818fSSatish Balay   PetscValidIntPointer(reason,6);
10670c52818fSSatish Balay   if (ls->new_x) {
10689566063dSJacob Faibussowitsch     PetscCall(VecCopy(ls->new_x,x));
10690c52818fSSatish Balay   }
10700c52818fSSatish Balay   *f = ls->new_f;
10710c52818fSSatish Balay   if (ls->new_g) {
10729566063dSJacob Faibussowitsch     PetscCall(VecCopy(ls->new_g,g));
10730c52818fSSatish Balay   }
107497ab8e72SStefano Zampini   if (steplength) *steplength = ls->step;
10750c52818fSSatish Balay   *reason = ls->reason;
10760c52818fSSatish Balay   PetscFunctionReturn(0);
10770c52818fSSatish Balay }
10780c52818fSSatish Balay 
10790c52818fSSatish Balay /*@
10800c52818fSSatish Balay   TaoLineSearchGetStartingVector - Gets a the initial point of the line
10810c52818fSSatish Balay   search.
10820c52818fSSatish Balay 
10830c52818fSSatish Balay   Not Collective
10840c52818fSSatish Balay 
10850c52818fSSatish Balay   Input Parameter:
10860c52818fSSatish Balay . ls - the TaoLineSearch context
10870c52818fSSatish Balay 
10880c52818fSSatish Balay   Output Parameter:
10890c52818fSSatish Balay . x - The initial point of the line search
10900c52818fSSatish Balay 
10910c52818fSSatish Balay   Level: intermediate
10920c52818fSSatish Balay @*/
10930c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
10940c52818fSSatish Balay {
10950c52818fSSatish Balay   PetscFunctionBegin;
10960c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
109797ab8e72SStefano Zampini   if (x) *x = ls->start_x;
10980c52818fSSatish Balay   PetscFunctionReturn(0);
10990c52818fSSatish Balay }
11000c52818fSSatish Balay 
11010c52818fSSatish Balay /*@
11020c52818fSSatish Balay   TaoLineSearchGetStepDirection - Gets the step direction of the line
11030c52818fSSatish Balay   search.
11040c52818fSSatish Balay 
11050c52818fSSatish Balay   Not Collective
11060c52818fSSatish Balay 
11070c52818fSSatish Balay   Input Parameter:
11080c52818fSSatish Balay . ls - the TaoLineSearch context
11090c52818fSSatish Balay 
11100c52818fSSatish Balay   Output Parameter:
11110c52818fSSatish Balay . s - the step direction of the line search
11120c52818fSSatish Balay 
11130c52818fSSatish Balay   Level: advanced
11140c52818fSSatish Balay @*/
11150c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
11160c52818fSSatish Balay {
11170c52818fSSatish Balay   PetscFunctionBegin;
11180c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
111997ab8e72SStefano Zampini   if (s) *s = ls->stepdirection;
11200c52818fSSatish Balay   PetscFunctionReturn(0);
11210c52818fSSatish Balay }
11220c52818fSSatish Balay 
11230c52818fSSatish Balay /*@
11240c52818fSSatish Balay   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.
11250c52818fSSatish Balay 
11260c52818fSSatish Balay   Not Collective
11270c52818fSSatish Balay 
11280c52818fSSatish Balay   Input Parameter:
11290c52818fSSatish Balay . ls - the TaoLineSearch context
11300c52818fSSatish Balay 
11310c52818fSSatish Balay   Output Parameter:
11320c52818fSSatish Balay . f - the objective value at the full step length
11330c52818fSSatish Balay 
11340c52818fSSatish Balay   Level: developer
11350c52818fSSatish Balay @*/
11360c52818fSSatish Balay 
11370c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
11380c52818fSSatish Balay {
11390c52818fSSatish Balay   PetscFunctionBegin;
11400c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
11410c52818fSSatish Balay   *f_fullstep = ls->f_fullstep;
11420c52818fSSatish Balay   PetscFunctionReturn(0);
11430c52818fSSatish Balay }
11440c52818fSSatish Balay 
11450c52818fSSatish Balay /*@
11460c52818fSSatish Balay   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.
11470c52818fSSatish Balay 
1148441846f8SBarry Smith   Logically Collective on Tao
11490c52818fSSatish Balay 
11500c52818fSSatish Balay   Input Parameters:
11510c52818fSSatish Balay + ls - the TaoLineSearch context
11520c52818fSSatish Balay . xl  - vector of lower bounds
11530c52818fSSatish Balay - xu  - vector of upper bounds
11540c52818fSSatish Balay 
11550c52818fSSatish Balay   Note: If the variable bounds are not set with this routine, then
1156e270355aSBarry Smith   PETSC_NINFINITY and PETSC_INFINITY are assumed
11570c52818fSSatish Balay 
11580c52818fSSatish Balay   Level: beginner
11590c52818fSSatish Balay 
1160db781477SPatrick Sanan .seealso: `TaoSetVariableBounds()`, `TaoLineSearchCreate()`
11610c52818fSSatish Balay @*/
11620c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
11630c52818fSSatish Balay {
11640c52818fSSatish Balay   PetscFunctionBegin;
11650c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1166*76be6f4fSStefano Zampini   if (xl) PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1167*76be6f4fSStefano Zampini   if (xu) PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
11689566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
11699566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
11709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->lower));
11719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->upper));
11720c52818fSSatish Balay   ls->lower = xl;
11730c52818fSSatish Balay   ls->upper = xu;
1174*76be6f4fSStefano Zampini   ls->bounded = (PetscBool)(xl || xu);
11750c52818fSSatish Balay   PetscFunctionReturn(0);
11760c52818fSSatish Balay }
11770c52818fSSatish Balay 
11780c52818fSSatish Balay /*@
11790c52818fSSatish Balay   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
11800c52818fSSatish Balay   search.  If this value is not set then 1.0 is assumed.
11810c52818fSSatish Balay 
11820c52818fSSatish Balay   Logically Collective on TaoLineSearch
11830c52818fSSatish Balay 
11840c52818fSSatish Balay   Input Parameters:
11850c52818fSSatish Balay + ls - the TaoLineSearch context
11860c52818fSSatish Balay - s - the initial step size
11870c52818fSSatish Balay 
11880c52818fSSatish Balay   Level: intermediate
11890c52818fSSatish Balay 
1190db781477SPatrick Sanan .seealso: `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()`
11910c52818fSSatish Balay @*/
11920c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
11930c52818fSSatish Balay {
11940c52818fSSatish Balay   PetscFunctionBegin;
11950c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
1196a39c8e28SStefano Zampini   PetscValidLogicalCollectiveReal(ls,s,2);
11970c52818fSSatish Balay   ls->initstep = s;
11980c52818fSSatish Balay   PetscFunctionReturn(0);
11990c52818fSSatish Balay }
12000c52818fSSatish Balay 
12010c52818fSSatish Balay /*@
12020c52818fSSatish Balay   TaoLineSearchGetStepLength - Get the current step length
12030c52818fSSatish Balay 
12040c52818fSSatish Balay   Not Collective
12050c52818fSSatish Balay 
12060c52818fSSatish Balay   Input Parameters:
12070c52818fSSatish Balay . ls - the TaoLineSearch context
12080c52818fSSatish Balay 
12097a7aea1fSJed Brown   Output Parameters:
12100c52818fSSatish Balay . s - the current step length
12110c52818fSSatish Balay 
12120c52818fSSatish Balay   Level: beginner
12130c52818fSSatish Balay 
1214db781477SPatrick Sanan .seealso: `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()`
12150c52818fSSatish Balay @*/
12160c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
12170c52818fSSatish Balay {
12180c52818fSSatish Balay   PetscFunctionBegin;
12190c52818fSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
12200c52818fSSatish Balay   *s = ls->step;
12210c52818fSSatish Balay   PetscFunctionReturn(0);
12220c52818fSSatish Balay }
12230c52818fSSatish Balay 
12240ebee16dSLisandro Dalcin /*@C
12250c52818fSSatish Balay    TaoLineSearchRegister - Adds a line-search algorithm to the registry
12260c52818fSSatish Balay 
12270c52818fSSatish Balay    Not collective
12280c52818fSSatish Balay 
12290c52818fSSatish Balay    Input Parameters:
12300c52818fSSatish Balay +  sname - name of a new user-defined solver
12310c52818fSSatish Balay -  func - routine to Create method context
12320c52818fSSatish Balay 
12330c52818fSSatish Balay    Notes:
12340c52818fSSatish Balay    TaoLineSearchRegister() may be called multiple times to add several user-defined solvers.
12350c52818fSSatish Balay 
12360c52818fSSatish Balay    Sample usage:
12370c52818fSSatish Balay .vb
12380c52818fSSatish Balay    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
12390c52818fSSatish Balay .ve
12400c52818fSSatish Balay 
12410c52818fSSatish Balay    Then, your solver can be chosen with the procedural interface via
12420c52818fSSatish Balay $     TaoLineSearchSetType(ls,"my_linesearch")
12430c52818fSSatish Balay    or at runtime via the option
12440c52818fSSatish Balay $     -tao_ls_type my_linesearch
12450c52818fSSatish Balay 
12460c52818fSSatish Balay    Level: developer
12470c52818fSSatish Balay 
12480ebee16dSLisandro Dalcin @*/
12490c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
12500c52818fSSatish Balay {
12510c52818fSSatish Balay   PetscFunctionBegin;
12529566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchInitializePackage());
12539566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func));
12540c52818fSSatish Balay   PetscFunctionReturn(0);
12550c52818fSSatish Balay }
12560c52818fSSatish Balay 
12570c52818fSSatish Balay /*@C
12580c52818fSSatish Balay    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
12590c52818fSSatish Balay    for all TaoLineSearch options in the database.
12600c52818fSSatish Balay 
12610c52818fSSatish Balay    Collective on TaoLineSearch
12620c52818fSSatish Balay 
12630c52818fSSatish Balay    Input Parameters:
12640c52818fSSatish Balay +  ls - the TaoLineSearch solver context
12650c52818fSSatish Balay -  prefix - the prefix string to prepend to all line search requests
12660c52818fSSatish Balay 
12670c52818fSSatish Balay    Notes:
12680c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
12690c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
12700c52818fSSatish Balay 
12710c52818fSSatish Balay    Level: advanced
12720c52818fSSatish Balay 
1273db781477SPatrick Sanan .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
12740c52818fSSatish Balay @*/
12750c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
12760c52818fSSatish Balay {
1277f06e3bfaSBarry Smith   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
12780c52818fSSatish Balay }
12790c52818fSSatish Balay 
12800c52818fSSatish Balay /*@C
12810c52818fSSatish Balay   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
12820c52818fSSatish Balay   TaoLineSearch options in the database
12830c52818fSSatish Balay 
12840c52818fSSatish Balay   Not Collective
12850c52818fSSatish Balay 
12860c52818fSSatish Balay   Input Parameters:
12870c52818fSSatish Balay . ls - the TaoLineSearch context
12880c52818fSSatish Balay 
12890c52818fSSatish Balay   Output Parameters:
12900c52818fSSatish Balay . prefix - pointer to the prefix string used is returned
12910c52818fSSatish Balay 
129295452b02SPatrick Sanan   Notes:
129395452b02SPatrick Sanan     On the fortran side, the user should pass in a string 'prefix' of
12940c52818fSSatish Balay   sufficient length to hold the prefix.
12950c52818fSSatish Balay 
12960c52818fSSatish Balay   Level: advanced
12970c52818fSSatish Balay 
1298db781477SPatrick Sanan .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()`
12990c52818fSSatish Balay @*/
13000c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
13010c52818fSSatish Balay {
1302f06e3bfaSBarry Smith   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
13030c52818fSSatish Balay }
13040c52818fSSatish Balay 
13050c52818fSSatish Balay /*@C
13060c52818fSSatish Balay    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
13070c52818fSSatish Balay    TaoLineSearch options in the database.
13080c52818fSSatish Balay 
13090c52818fSSatish Balay    Logically Collective on TaoLineSearch
13100c52818fSSatish Balay 
13110c52818fSSatish Balay    Input Parameters:
13120c52818fSSatish Balay +  ls - the TaoLineSearch context
13130c52818fSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
13140c52818fSSatish Balay 
13150c52818fSSatish Balay    Notes:
13160c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
13170c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
13180c52818fSSatish Balay 
13190c52818fSSatish Balay    For example, to distinguish between the runtime options for two
13200c52818fSSatish Balay    different line searches, one could call
13210c52818fSSatish Balay .vb
13220c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
13230c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
13240c52818fSSatish Balay .ve
13250c52818fSSatish Balay 
13260c52818fSSatish Balay    This would enable use of different options for each system, such as
13270c52818fSSatish Balay .vb
13280c52818fSSatish Balay       -sys1_tao_ls_type mt
13290c52818fSSatish Balay       -sys2_tao_ls_type armijo
13300c52818fSSatish Balay .ve
13310c52818fSSatish Balay 
13320c52818fSSatish Balay    Level: advanced
13330c52818fSSatish Balay 
1334db781477SPatrick Sanan .seealso: `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
13350c52818fSSatish Balay @*/
13360c52818fSSatish Balay 
13370c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
13380c52818fSSatish Balay {
1339f06e3bfaSBarry Smith   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
13400c52818fSSatish Balay }
1341