xref: /petsc/src/tao/linesearch/interface/taolinesearch.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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 
14*c3339decSBarry Smith    Collective
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 @*/
24d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A, PetscObject obj, const char name[])
25d71ae5a4SJacob Faibussowitsch {
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 
35*c3339decSBarry Smith   Collective
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 
57d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
58d71ae5a4SJacob Faibussowitsch {
590c52818fSSatish Balay   PetscBool         isascii, isstring;
60dedfbcbeSJed Brown   TaoLineSearchType type;
61f06e3bfaSBarry Smith 
620c52818fSSatish Balay   PetscFunctionBegin;
630c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
6448a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer));
650c52818fSSatish Balay   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
660c52818fSSatish Balay   PetscCheckSameComm(ls, 1, viewer, 2);
670c52818fSSatish Balay 
689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
699566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
700c52818fSSatish Balay   if (isascii) {
719566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer));
729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
73dbbe0bcdSBarry Smith     PetscTryTypeMethod(ls, view, viewer);
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
7663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "maximum function evaluations=%" PetscInt_FMT "\n", ls->max_funcs));
779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "tolerances: ftol=%g, rtol=%g, gtol=%g\n", (double)ls->ftol, (double)ls->rtol, (double)ls->gtol));
7863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT "\n", ls->nfeval));
7963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT "\n", ls->ngeval));
8063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT "\n", ls->nfgeval));
810c52818fSSatish Balay 
8248a46eb9SPierre Jolivet     if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(viewer, "using variable bounds\n"));
839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Termination reason: %d\n", (int)ls->reason));
849566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
850c52818fSSatish Balay   } else if (isstring) {
869566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchGetType(ls, &type));
879566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
880c52818fSSatish Balay   }
890c52818fSSatish Balay   PetscFunctionReturn(0);
900c52818fSSatish Balay }
910c52818fSSatish Balay 
920c52818fSSatish Balay /*@C
930c52818fSSatish Balay   TaoLineSearchCreate - Creates a TAO Line Search object.  Algorithms in TAO that use
940c52818fSSatish Balay   line-searches will automatically create one.
950c52818fSSatish Balay 
96d083f849SBarry Smith   Collective
970c52818fSSatish Balay 
980c52818fSSatish Balay   Input Parameter:
990c52818fSSatish Balay . comm - MPI communicator
1000c52818fSSatish Balay 
1010c52818fSSatish Balay   Output Parameter:
1020c52818fSSatish Balay . newls - the new TaoLineSearch context
1030c52818fSSatish Balay 
1040c52818fSSatish Balay   Available methods include:
105147403d9SBarry Smith + more-thuente - the More-Thuente method
106147403d9SBarry Smith . gpcg - the GPCG method
1070c52818fSSatish Balay - unit - Do not perform any line search
1080c52818fSSatish Balay 
1090c52818fSSatish Balay    Options Database Keys:
1100c52818fSSatish Balay .   -tao_ls_type - select which method TAO should use
1110c52818fSSatish Balay 
1120c52818fSSatish Balay    Level: beginner
1130c52818fSSatish Balay 
114db781477SPatrick Sanan .seealso: `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()`
1150c52818fSSatish Balay @*/
1160c52818fSSatish Balay 
117d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
118d71ae5a4SJacob Faibussowitsch {
1190c52818fSSatish Balay   TaoLineSearch ls;
1200c52818fSSatish Balay 
1210c52818fSSatish Balay   PetscFunctionBegin;
1220c52818fSSatish Balay   PetscValidPointer(newls, 2);
1239566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchInitializePackage());
1240c52818fSSatish Balay 
1259566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(ls, TAOLINESEARCH_CLASSID, "TaoLineSearch", "Linesearch", "Tao", comm, TaoLineSearchDestroy, TaoLineSearchView));
1260c52818fSSatish Balay   ls->max_funcs = 30;
1270c52818fSSatish Balay   ls->ftol      = 0.0001;
1280c52818fSSatish Balay   ls->gtol      = 0.9;
1296f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1306f4723b1SBarry Smith   ls->rtol = 1.0e-5;
1316f4723b1SBarry Smith #else
1320c52818fSSatish Balay   ls->rtol = 1.0e-10;
1336f4723b1SBarry Smith #endif
1340c52818fSSatish Balay   ls->stepmin  = 1.0e-20;
1350c52818fSSatish Balay   ls->stepmax  = 1.0e+20;
1360c52818fSSatish Balay   ls->step     = 1.0;
137a39c8e28SStefano Zampini   ls->initstep = 1.0;
1380c52818fSSatish Balay   *newls       = ls;
1390c52818fSSatish Balay   PetscFunctionReturn(0);
1400c52818fSSatish Balay }
1410c52818fSSatish Balay 
1420c52818fSSatish Balay /*@
1430c52818fSSatish Balay   TaoLineSearchSetUp - Sets up the internal data structures for the later use
1440c52818fSSatish Balay   of a Tao solver
1450c52818fSSatish Balay 
146*c3339decSBarry Smith   Collective
1470c52818fSSatish Balay 
1480c52818fSSatish Balay   Input Parameters:
1490c52818fSSatish Balay . ls - the TaoLineSearch context
1500c52818fSSatish Balay 
1510c52818fSSatish Balay   Notes:
1520c52818fSSatish Balay   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
1530c52818fSSatish Balay   automatically be called in TaoLineSearchSolve().  However, if the user
1540c52818fSSatish Balay   desires to call it explicitly, it should come after TaoLineSearchCreate()
1550c52818fSSatish Balay   but before TaoLineSearchApply().
1560c52818fSSatish Balay 
1570c52818fSSatish Balay   Level: developer
1580c52818fSSatish Balay 
159db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()`
1600c52818fSSatish Balay @*/
1610c52818fSSatish Balay 
162d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
163d71ae5a4SJacob Faibussowitsch {
1648caf6e8cSBarry Smith   const char *default_type = TAOLINESEARCHMT;
1650c52818fSSatish Balay   PetscBool   flg;
166f06e3bfaSBarry Smith 
1670c52818fSSatish Balay   PetscFunctionBegin;
1680c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1690c52818fSSatish Balay   if (ls->setupcalled) PetscFunctionReturn(0);
17048a46eb9SPierre Jolivet   if (!((PetscObject)ls)->type_name) PetscCall(TaoLineSearchSetType(ls, default_type));
171dbbe0bcdSBarry Smith   PetscTryTypeMethod(ls, setup);
1720c52818fSSatish Balay   if (ls->usetaoroutines) {
1739566063dSJacob Faibussowitsch     PetscCall(TaoIsObjectiveDefined(ls->tao, &flg));
1740c52818fSSatish Balay     ls->hasobjective = flg;
1759566063dSJacob Faibussowitsch     PetscCall(TaoIsGradientDefined(ls->tao, &flg));
1760c52818fSSatish Balay     ls->hasgradient = flg;
1779566063dSJacob Faibussowitsch     PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao, &flg));
1780c52818fSSatish Balay     ls->hasobjectiveandgradient = flg;
1790c52818fSSatish Balay   } else {
1800c52818fSSatish Balay     if (ls->ops->computeobjective) {
1810c52818fSSatish Balay       ls->hasobjective = PETSC_TRUE;
1820c52818fSSatish Balay     } else {
1830c52818fSSatish Balay       ls->hasobjective = PETSC_FALSE;
1840c52818fSSatish Balay     }
1850c52818fSSatish Balay     if (ls->ops->computegradient) {
1860c52818fSSatish Balay       ls->hasgradient = PETSC_TRUE;
1870c52818fSSatish Balay     } else {
1880c52818fSSatish Balay       ls->hasgradient = PETSC_FALSE;
1890c52818fSSatish Balay     }
1900c52818fSSatish Balay     if (ls->ops->computeobjectiveandgradient) {
1910c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_TRUE;
1920c52818fSSatish Balay     } else {
1930c52818fSSatish Balay       ls->hasobjectiveandgradient = PETSC_FALSE;
1940c52818fSSatish Balay     }
1950c52818fSSatish Balay   }
1960c52818fSSatish Balay   ls->setupcalled = PETSC_TRUE;
1970c52818fSSatish Balay   PetscFunctionReturn(0);
1980c52818fSSatish Balay }
1990c52818fSSatish Balay 
2000c52818fSSatish Balay /*@
2010c52818fSSatish Balay   TaoLineSearchReset - Some line searches may carry state information
2020c52818fSSatish Balay   from one TaoLineSearchApply() to the next.  This function resets this
2030c52818fSSatish Balay   state information.
2040c52818fSSatish Balay 
205*c3339decSBarry Smith   Collective
2060c52818fSSatish Balay 
2070c52818fSSatish Balay   Input Parameter:
2080c52818fSSatish Balay . ls - the TaoLineSearch context
2090c52818fSSatish Balay 
2100c52818fSSatish Balay   Level: developer
2110c52818fSSatish Balay 
212db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()`
2130c52818fSSatish Balay @*/
214d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
215d71ae5a4SJacob Faibussowitsch {
2160c52818fSSatish Balay   PetscFunctionBegin;
2170c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
218dbbe0bcdSBarry Smith   PetscTryTypeMethod(ls, reset);
2190c52818fSSatish Balay   PetscFunctionReturn(0);
2200c52818fSSatish Balay }
221f06e3bfaSBarry Smith 
2220c52818fSSatish Balay /*@
2230c52818fSSatish Balay   TaoLineSearchDestroy - Destroys the TAO context that was created with
2240c52818fSSatish Balay   TaoLineSearchCreate()
2250c52818fSSatish Balay 
226*c3339decSBarry Smith   Collective
2270c52818fSSatish Balay 
2287a7aea1fSJed Brown   Input Parameter:
2290c52818fSSatish Balay . ls - the TaoLineSearch context
2300c52818fSSatish Balay 
2310c52818fSSatish Balay   Level: beginner
2320c52818fSSatish Balay 
2330c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
2340c52818fSSatish Balay @*/
235d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
236d71ae5a4SJacob Faibussowitsch {
2370c52818fSSatish Balay   PetscFunctionBegin;
2380c52818fSSatish Balay   if (!*ls) PetscFunctionReturn(0);
2390c52818fSSatish Balay   PetscValidHeaderSpecific(*ls, TAOLINESEARCH_CLASSID, 1);
2409371c9d4SSatish Balay   if (--((PetscObject)*ls)->refct > 0) {
2419371c9d4SSatish Balay     *ls = NULL;
2429371c9d4SSatish Balay     PetscFunctionReturn(0);
2439371c9d4SSatish Balay   }
2449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->stepdirection));
2459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->start_x));
2469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->upper));
2479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ls)->lower));
24848a46eb9SPierre Jolivet   if ((*ls)->ops->destroy) PetscCall((*(*ls)->ops->destroy)(*ls));
24948a46eb9SPierre Jolivet   if ((*ls)->usemonitor) PetscCall(PetscViewerDestroy(&(*ls)->viewer));
2509566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(ls));
2510c52818fSSatish Balay   PetscFunctionReturn(0);
2520c52818fSSatish Balay }
2530c52818fSSatish Balay 
2540c52818fSSatish Balay /*@
2550c52818fSSatish Balay   TaoLineSearchApply - Performs a line-search in a given step direction.  Criteria for acceptable step length depends on the line-search algorithm chosen
2560c52818fSSatish Balay 
257*c3339decSBarry Smith   Collective
2580c52818fSSatish Balay 
2590c52818fSSatish Balay   Input Parameters:
260441846f8SBarry Smith + ls - the Tao context
2610c52818fSSatish Balay - s - search direction
2620c52818fSSatish Balay 
2636b867d5aSJose E. Roman   Input/Output Parameters:
2646b867d5aSJose E. Roman 
2650c52818fSSatish Balay   Output Parameters:
266f1a722f8SMatthew G. Knepley + x - On input the current solution, on output x contains the new solution determined by the line search
267f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution
268f1a722f8SMatthew G. Knepley . g - On input the gradient evaluated at x, on output contains the gradient at new solution
269f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x)
2700c52818fSSatish Balay - reason - reason why the line-search stopped
2710c52818fSSatish Balay 
27297bb3fdcSJose E. Roman   Notes:
2730c52818fSSatish Balay   reason will be set to one of:
2740c52818fSSatish Balay 
2750c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
2760c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
2770c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
2780c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
2790c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
2800c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
2810c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
2820c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
2830c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
2840c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search
2850c52818fSSatish Balay 
2860c52818fSSatish Balay   The algorithm developer must set up the TaoLineSearch with calls to
287441846f8SBarry Smith   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()
2880c52818fSSatish Balay 
2890c52818fSSatish Balay   You may or may not need to follow this with a call to
2900c52818fSSatish Balay   TaoAddLineSearchCounts(), depending on whether you want these
2910c52818fSSatish Balay   evaluations to count toward the total function/gradient evaluations.
2920c52818fSSatish Balay 
2930c52818fSSatish Balay   Level: beginner
2940c52818fSSatish Balay 
295db781477SPatrick Sanan   .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()`
2960c52818fSSatish Balay  @*/
2970c52818fSSatish Balay 
298d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
299d71ae5a4SJacob Faibussowitsch {
3000c52818fSSatish Balay   PetscInt low1, low2, low3, high1, high2, high3;
3010c52818fSSatish Balay 
3020c52818fSSatish Balay   PetscFunctionBegin;
3030c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
3040c52818fSSatish Balay   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
3053f6ba705SLisandro Dalcin   PetscValidRealPointer(f, 3);
3060c52818fSSatish Balay   PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
3070c52818fSSatish Balay   PetscValidHeaderSpecific(s, VEC_CLASSID, 5);
3080c52818fSSatish Balay   PetscValidPointer(reason, 7);
3090c52818fSSatish Balay   PetscCheckSameComm(ls, 1, x, 2);
3100c52818fSSatish Balay   PetscCheckSameTypeAndComm(x, 2, g, 4);
3110c52818fSSatish Balay   PetscCheckSameTypeAndComm(x, 2, s, 5);
3129566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &low1, &high1));
3139566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(g, &low2, &high2));
3149566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(s, &low3, &high3));
3153c859ba3SBarry Smith   PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible vector local lengths");
3160c52818fSSatish Balay 
31797ab8e72SStefano Zampini   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
3189566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
3199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->stepdirection));
320050fc7a3SBarry Smith   ls->stepdirection = s;
3210c52818fSSatish Balay 
3229566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetUp(ls));
3230c52818fSSatish Balay   ls->nfeval  = 0;
3240c52818fSSatish Balay   ls->ngeval  = 0;
3250c52818fSSatish Balay   ls->nfgeval = 0;
3260c52818fSSatish Balay   /* Check parameter values */
3270c52818fSSatish Balay   if (ls->ftol < 0.0) {
3289566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: ftol (%g) < 0\n", (double)ls->ftol));
3290c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3300c52818fSSatish Balay   }
3310c52818fSSatish Balay   if (ls->rtol < 0.0) {
3329566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: rtol (%g) < 0\n", (double)ls->rtol));
3330c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3340c52818fSSatish Balay   }
3350c52818fSSatish Balay   if (ls->gtol < 0.0) {
3369566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: gtol (%g) < 0\n", (double)ls->gtol));
3370c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3380c52818fSSatish Balay   }
3390c52818fSSatish Balay   if (ls->stepmin < 0.0) {
3409566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) < 0\n", (double)ls->stepmin));
3410c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3420c52818fSSatish Balay   }
3430c52818fSSatish Balay   if (ls->stepmax < ls->stepmin) {
3449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n", (double)ls->stepmin, (double)ls->stepmax));
3450c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3460c52818fSSatish Balay   }
3470c52818fSSatish Balay   if (ls->max_funcs < 0) {
3489566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n", ls->max_funcs));
3490c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
3500c52818fSSatish Balay   }
3510c52818fSSatish Balay   if (PetscIsInfOrNanReal(*f)) {
3529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Initial Line Search Function Value is Inf or Nan (%g)\n", (double)*f));
3530c52818fSSatish Balay     *reason = TAOLINESEARCH_FAILED_INFORNAN;
3540c52818fSSatish Balay   }
3550c52818fSSatish Balay 
3569566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)x));
3579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->start_x));
3580c52818fSSatish Balay   ls->start_x = x;
359050fc7a3SBarry Smith 
3609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply, ls, 0, 0, 0));
361dbbe0bcdSBarry Smith   PetscUseTypeMethod(ls, apply, x, f, g, s);
3629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0, 0, 0));
3630c52818fSSatish Balay   *reason   = ls->reason;
3640c52818fSSatish Balay   ls->new_f = *f;
3650c52818fSSatish Balay 
36697ab8e72SStefano Zampini   if (steplength) *steplength = ls->step;
3670c52818fSSatish Balay 
3689566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchViewFromOptions(ls, NULL, "-tao_ls_view"));
3690c52818fSSatish Balay   PetscFunctionReturn(0);
3700c52818fSSatish Balay }
3710c52818fSSatish Balay 
3720c52818fSSatish Balay /*@C
3730c52818fSSatish Balay    TaoLineSearchSetType - Sets the algorithm used in a line search
3740c52818fSSatish Balay 
375*c3339decSBarry Smith    Collective
3760c52818fSSatish Balay 
3770c52818fSSatish Balay    Input Parameters:
3780c52818fSSatish Balay +  ls - the TaoLineSearch context
379820824deSAlp Dener -  type - the TaoLineSearchType selection
3800c52818fSSatish Balay 
3810c52818fSSatish Balay   Available methods include:
382820824deSAlp Dener +  more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition
383820824deSAlp Dener .  armijo - simple backtracking line search enforcing only the sufficient decrease condition
384820824deSAlp Dener -  unit - do not perform a line search and always accept unit step length
3850c52818fSSatish Balay 
3860c52818fSSatish Balay   Options Database Keys:
387820824deSAlp Dener .  -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime
3880c52818fSSatish Balay 
3890c52818fSSatish Balay   Level: beginner
3900c52818fSSatish Balay 
391db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchGetType()`, `TaoLineSearchApply()`
3920c52818fSSatish Balay 
3930c52818fSSatish Balay @*/
3940c52818fSSatish Balay 
395d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
396d71ae5a4SJacob Faibussowitsch {
3970c52818fSSatish Balay   PetscErrorCode (*r)(TaoLineSearch);
3980c52818fSSatish Balay   PetscBool flg;
3990c52818fSSatish Balay 
4000c52818fSSatish Balay   PetscFunctionBegin;
4010c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
4020c52818fSSatish Balay   PetscValidCharPointer(type, 2);
4039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg));
4040c52818fSSatish Balay   if (flg) PetscFunctionReturn(0);
4050c52818fSSatish Balay 
4069566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TaoLineSearchList, type, (void (**)(void)) & r));
4073c859ba3SBarry Smith   PetscCheck(r, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested TaoLineSearch type %s", type);
408dbbe0bcdSBarry Smith   PetscTryTypeMethod(ls, destroy);
4090c52818fSSatish Balay   ls->max_funcs = 30;
4100c52818fSSatish Balay   ls->ftol      = 0.0001;
4110c52818fSSatish Balay   ls->gtol      = 0.9;
4126f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
4136f4723b1SBarry Smith   ls->rtol = 1.0e-5;
4146f4723b1SBarry Smith #else
4150c52818fSSatish Balay   ls->rtol = 1.0e-10;
4166f4723b1SBarry Smith #endif
4170c52818fSSatish Balay   ls->stepmin = 1.0e-20;
4180c52818fSSatish Balay   ls->stepmax = 1.0e+20;
4190c52818fSSatish Balay 
4200c52818fSSatish Balay   ls->nfeval              = 0;
4210c52818fSSatish Balay   ls->ngeval              = 0;
4220c52818fSSatish Balay   ls->nfgeval             = 0;
42383c8fe1dSLisandro Dalcin   ls->ops->setup          = NULL;
42483c8fe1dSLisandro Dalcin   ls->ops->apply          = NULL;
42583c8fe1dSLisandro Dalcin   ls->ops->view           = NULL;
42683c8fe1dSLisandro Dalcin   ls->ops->setfromoptions = NULL;
42783c8fe1dSLisandro Dalcin   ls->ops->destroy        = NULL;
4280c52818fSSatish Balay   ls->setupcalled         = PETSC_FALSE;
4299566063dSJacob Faibussowitsch   PetscCall((*r)(ls));
4309566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type));
4310c52818fSSatish Balay   PetscFunctionReturn(0);
4320c52818fSSatish Balay }
4330c52818fSSatish Balay 
4342a0dac07SAlp Dener /*@C
4352a0dac07SAlp Dener   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
4362a0dac07SAlp Dener   iteration number, step length, and function value before calling the implementation
4372a0dac07SAlp Dener   specific monitor.
4382a0dac07SAlp Dener 
4392a0dac07SAlp Dener    Input Parameters:
4402a0dac07SAlp Dener +  ls - the TaoLineSearch context
4412a0dac07SAlp Dener .  its - the current iterate number (>=0)
4422a0dac07SAlp Dener .  f - the current objective function value
4432a0dac07SAlp Dener -  step - the step length
4442a0dac07SAlp Dener 
4452a0dac07SAlp Dener    Options Database Key:
4462a0dac07SAlp Dener .  -tao_ls_monitor - Use the default monitor, which prints statistics to standard output
4472a0dac07SAlp Dener 
4482a0dac07SAlp Dener    Level: developer
4492a0dac07SAlp Dener 
4502a0dac07SAlp Dener @*/
451d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
452d71ae5a4SJacob Faibussowitsch {
4532a0dac07SAlp Dener   PetscInt tabs;
4542a0dac07SAlp Dener 
4552a0dac07SAlp Dener   PetscFunctionBegin;
4562a0dac07SAlp Dener   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
4572a0dac07SAlp Dener   if (ls->usemonitor) {
4589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs));
4599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel));
46063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its));
4619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f));
4629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step));
4632a0dac07SAlp Dener     if (ls->ops->monitor && its > 0) {
4649566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3));
465dbbe0bcdSBarry Smith       PetscUseTypeMethod(ls, monitor);
4662a0dac07SAlp Dener     }
4679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs));
4682a0dac07SAlp Dener   }
4692a0dac07SAlp Dener   PetscFunctionReturn(0);
4702a0dac07SAlp Dener }
4712a0dac07SAlp Dener 
4720c52818fSSatish Balay /*@
4730c52818fSSatish Balay   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
4740c52818fSSatish Balay   options.
4750c52818fSSatish Balay 
476*c3339decSBarry Smith   Collective
4770c52818fSSatish Balay 
47801d2d390SJose E. Roman   Input Parameter:
4790c52818fSSatish Balay . ls - the TaoLineSearch context
4800c52818fSSatish Balay 
4810c52818fSSatish Balay   Options Database Keys:
4820c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
4830c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease
4840c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition
4850c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step
486a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - initial steplength allowed
4870c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed
4880c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed
4890c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
4900c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output
4910c52818fSSatish Balay 
4920c52818fSSatish Balay   Level: beginner
4930c52818fSSatish Balay @*/
494d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
495d71ae5a4SJacob Faibussowitsch {
4968caf6e8cSBarry Smith   const char *default_type = TAOLINESEARCHMT;
4972a0dac07SAlp Dener   char        type[256], monfilename[PETSC_MAX_PATH_LEN];
4982a0dac07SAlp Dener   PetscViewer monviewer;
4990c52818fSSatish Balay   PetscBool   flg;
500f06e3bfaSBarry Smith 
5010c52818fSSatish Balay   PetscFunctionBegin;
5020c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
503d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)ls);
504ad540459SPierre Jolivet   if (((PetscObject)ls)->type_name) default_type = ((PetscObject)ls)->type_name;
5050c52818fSSatish Balay   /* Check for type from options */
5069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-tao_ls_type", "Tao Line Search type", "TaoLineSearchSetType", TaoLineSearchList, default_type, type, 256, &flg));
5070c52818fSSatish Balay   if (flg) {
5089566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(ls, type));
5090c52818fSSatish Balay   } else if (!((PetscObject)ls)->type_name) {
5109566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetType(ls, default_type));
5110c52818fSSatish Balay   }
5120c52818fSSatish Balay 
5139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_ls_max_funcs", "max function evals in line search", "", ls->max_funcs, &ls->max_funcs, NULL));
5149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_ftol", "tol for sufficient decrease", "", ls->ftol, &ls->ftol, NULL));
5159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_gtol", "tol for curvature condition", "", ls->gtol, &ls->gtol, NULL));
5169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_rtol", "relative tol for acceptable step", "", ls->rtol, &ls->rtol, NULL));
5179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_stepmin", "lower bound for step", "", ls->stepmin, &ls->stepmin, NULL));
5189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ls_stepmax", "upper bound for step", "", ls->stepmax, &ls->stepmax, NULL));
519a39c8e28SStefano Zampini   PetscCall(PetscOptionsReal("-tao_ls_stepinit", "initial step", "", ls->initstep, &ls->initstep, NULL));
5209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-tao_ls_monitor", "enable the basic monitor", "TaoLineSearchSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg));
5212a0dac07SAlp Dener   if (flg) {
5229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls), monfilename, &monviewer));
5232a0dac07SAlp Dener     ls->viewer     = monviewer;
5242a0dac07SAlp Dener     ls->usemonitor = PETSC_TRUE;
5252a0dac07SAlp Dener   }
526dbbe0bcdSBarry Smith   PetscTryTypeMethod(ls, setfromoptions, PetscOptionsObject);
527d0609cedSBarry Smith   PetscOptionsEnd();
5280c52818fSSatish Balay   PetscFunctionReturn(0);
5290c52818fSSatish Balay }
5300c52818fSSatish Balay 
5310c52818fSSatish Balay /*@C
5320c52818fSSatish Balay   TaoLineSearchGetType - Gets the current line search algorithm
5330c52818fSSatish Balay 
5340c52818fSSatish Balay   Not Collective
5350c52818fSSatish Balay 
5360c52818fSSatish Balay   Input Parameter:
5370c52818fSSatish Balay . ls - the TaoLineSearch context
5380c52818fSSatish Balay 
539fd292e60Sprj-   Output Parameter:
5400c52818fSSatish Balay . type - the line search algorithm in effect
5410c52818fSSatish Balay 
5420c52818fSSatish Balay   Level: developer
5430c52818fSSatish Balay 
5440c52818fSSatish Balay @*/
545d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
546d71ae5a4SJacob Faibussowitsch {
5470c52818fSSatish Balay   PetscFunctionBegin;
5480c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
5490c52818fSSatish Balay   PetscValidPointer(type, 2);
5500c52818fSSatish Balay   *type = ((PetscObject)ls)->type_name;
5510c52818fSSatish Balay   PetscFunctionReturn(0);
5520c52818fSSatish Balay }
5530c52818fSSatish Balay 
5540c52818fSSatish Balay /*@
5550c52818fSSatish Balay   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
5560c52818fSSatish Balay   routines used by the line search in last application (not cumulative).
5570c52818fSSatish Balay 
5580c52818fSSatish Balay   Not Collective
5590c52818fSSatish Balay 
5600c52818fSSatish Balay   Input Parameter:
5610c52818fSSatish Balay . ls - the TaoLineSearch context
5620c52818fSSatish Balay 
5630c52818fSSatish Balay   Output Parameters:
5640c52818fSSatish Balay + nfeval   - number of function evaluations
5650c52818fSSatish Balay . ngeval   - number of gradient evaluations
5660c52818fSSatish Balay - nfgeval  - number of function/gradient evaluations
5670c52818fSSatish Balay 
5680c52818fSSatish Balay   Level: intermediate
5690c52818fSSatish Balay 
5700c52818fSSatish Balay   Note:
571441846f8SBarry Smith   If the line search is using the Tao objective and gradient
572441846f8SBarry Smith   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
5730c52818fSSatish Balay   is already counting the number of evaluations.
5740c52818fSSatish Balay 
5750c52818fSSatish Balay @*/
576d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
577d71ae5a4SJacob Faibussowitsch {
5780c52818fSSatish Balay   PetscFunctionBegin;
5790c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
5800c52818fSSatish Balay   *nfeval  = ls->nfeval;
5810c52818fSSatish Balay   *ngeval  = ls->ngeval;
5820c52818fSSatish Balay   *nfgeval = ls->nfgeval;
5830c52818fSSatish Balay   PetscFunctionReturn(0);
5840c52818fSSatish Balay }
5850c52818fSSatish Balay 
5860c52818fSSatish Balay /*@
587441846f8SBarry Smith   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
588441846f8SBarry Smith   Tao evaluation routines.
5890c52818fSSatish Balay 
5900c52818fSSatish Balay   Not Collective
5910c52818fSSatish Balay 
5920c52818fSSatish Balay   Input Parameter:
5930c52818fSSatish Balay . ls - the TaoLineSearch context
5940c52818fSSatish Balay 
5950c52818fSSatish Balay   Output Parameter:
596441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
5970c52818fSSatish Balay         otherwise PETSC_FALSE
5980c52818fSSatish Balay 
5990c52818fSSatish Balay   Level: developer
6000c52818fSSatish Balay @*/
601d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
602d71ae5a4SJacob Faibussowitsch {
6030c52818fSSatish Balay   PetscFunctionBegin;
6040c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
605f06e3bfaSBarry Smith   *flg = ls->usetaoroutines;
6060c52818fSSatish Balay   PetscFunctionReturn(0);
6070c52818fSSatish Balay }
6080c52818fSSatish Balay 
6090c52818fSSatish Balay /*@C
6100c52818fSSatish Balay   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
6110c52818fSSatish Balay 
612*c3339decSBarry Smith   Logically Collective
6130c52818fSSatish Balay 
614d8d19677SJose E. Roman   Input Parameters:
6150c52818fSSatish Balay + ls - the TaoLineSearch context
6160c52818fSSatish Balay . func - the objective function evaluation routine
6170c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6180c52818fSSatish Balay 
6190c52818fSSatish Balay   Calling sequence of func:
6200c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);
6210c52818fSSatish Balay 
6220c52818fSSatish Balay + x - input vector
6230c52818fSSatish Balay . f - function value
6240c52818fSSatish Balay - ctx (optional) user-defined context
6250c52818fSSatish Balay 
6260c52818fSSatish Balay   Level: beginner
6270c52818fSSatish Balay 
6280c52818fSSatish Balay   Note:
6290c52818fSSatish Balay   Use this routine only if you want the line search objective
630441846f8SBarry Smith   evaluation routine to be different from the Tao's objective
6310c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
6320c52818fSSatish Balay   the line search gradient and/or function/gradient routine.
6330c52818fSSatish Balay 
6340c52818fSSatish Balay   Note:
6350c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
6360c52818fSSatish Balay   line search, application programmers should be wary of overriding the
6370c52818fSSatish Balay   default objective routine.
6380c52818fSSatish Balay 
639db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
6400c52818fSSatish Balay @*/
641d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *, void *), void *ctx)
642d71ae5a4SJacob Faibussowitsch {
6430c52818fSSatish Balay   PetscFunctionBegin;
6440c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
6450c52818fSSatish Balay 
6460c52818fSSatish Balay   ls->ops->computeobjective = func;
6470c52818fSSatish Balay   if (ctx) ls->userctx_func = ctx;
6480c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
6490c52818fSSatish Balay   PetscFunctionReturn(0);
6500c52818fSSatish Balay }
6510c52818fSSatish Balay 
6520c52818fSSatish Balay /*@C
6530c52818fSSatish Balay   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
6540c52818fSSatish Balay 
655*c3339decSBarry Smith   Logically Collective
6560c52818fSSatish Balay 
657d8d19677SJose E. Roman   Input Parameters:
6580c52818fSSatish Balay + ls - the TaoLineSearch context
6590c52818fSSatish Balay . func - the gradient evaluation routine
6600c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
6610c52818fSSatish Balay 
6620c52818fSSatish Balay   Calling sequence of func:
6630c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);
6640c52818fSSatish Balay 
6650c52818fSSatish Balay + x - input vector
6660c52818fSSatish Balay . g - gradient vector
6670c52818fSSatish Balay - ctx (optional) user-defined context
6680c52818fSSatish Balay 
6690c52818fSSatish Balay   Level: beginner
6700c52818fSSatish Balay 
6710c52818fSSatish Balay   Note:
6720c52818fSSatish Balay   Use this routine only if you want the line search gradient
673441846f8SBarry Smith   evaluation routine to be different from the Tao's gradient
6740c52818fSSatish Balay   evaluation routine. If you use this routine you must also set
6750c52818fSSatish Balay   the line search function and/or function/gradient routine.
6760c52818fSSatish Balay 
6770c52818fSSatish Balay   Note:
6780c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own gradient routine for the
6790c52818fSSatish Balay   line search, application programmers should be wary of overriding the
6800c52818fSSatish Balay   default gradient routine.
6810c52818fSSatish Balay 
682db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
6830c52818fSSatish Balay @*/
684d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec g, void *), void *ctx)
685d71ae5a4SJacob Faibussowitsch {
6860c52818fSSatish Balay   PetscFunctionBegin;
6870c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
6880c52818fSSatish Balay   ls->ops->computegradient = func;
6890c52818fSSatish Balay   if (ctx) ls->userctx_grad = ctx;
6900c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
6910c52818fSSatish Balay   PetscFunctionReturn(0);
6920c52818fSSatish Balay }
6930c52818fSSatish Balay 
6940c52818fSSatish Balay /*@C
6950c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
6960c52818fSSatish Balay 
697*c3339decSBarry Smith   Logically Collective
6980c52818fSSatish Balay 
699d8d19677SJose E. Roman   Input Parameters:
7000c52818fSSatish Balay + ls - the TaoLineSearch context
7010c52818fSSatish Balay . func - the objective and gradient evaluation routine
7020c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7030c52818fSSatish Balay 
7040c52818fSSatish Balay   Calling sequence of func:
7050c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);
7060c52818fSSatish Balay 
7070c52818fSSatish Balay + x - input vector
7080c52818fSSatish Balay . f - function value
7090c52818fSSatish Balay . g - gradient vector
7100c52818fSSatish Balay - ctx (optional) user-defined context
7110c52818fSSatish Balay 
7120c52818fSSatish Balay   Level: beginner
7130c52818fSSatish Balay 
7140c52818fSSatish Balay   Note:
7150c52818fSSatish Balay   Use this routine only if you want the line search objective and gradient
716441846f8SBarry Smith   evaluation routines to be different from the Tao's objective
7170c52818fSSatish Balay   and gradient evaluation routines.
7180c52818fSSatish Balay 
7190c52818fSSatish Balay   Note:
7200c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
7210c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7220c52818fSSatish Balay   default objective routine.
7230c52818fSSatish Balay 
724db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
7250c52818fSSatish Balay @*/
726d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void *), void *ctx)
727d71ae5a4SJacob Faibussowitsch {
7280c52818fSSatish Balay   PetscFunctionBegin;
7290c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
7300c52818fSSatish Balay   ls->ops->computeobjectiveandgradient = func;
7310c52818fSSatish Balay   if (ctx) ls->userctx_funcgrad = ctx;
7320c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
7330c52818fSSatish Balay   PetscFunctionReturn(0);
7340c52818fSSatish Balay }
7350c52818fSSatish Balay 
7360c52818fSSatish Balay /*@C
7370c52818fSSatish Balay   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
7380c52818fSSatish Balay   (gradient'*stepdirection) evaluation routine for the line search.
7390c52818fSSatish Balay   Sometimes it is more efficient to compute the inner product of the gradient
7400c52818fSSatish Balay   and the step direction than it is to compute the gradient, and this is all
7410c52818fSSatish Balay   the line search typically needs of the gradient.
7420c52818fSSatish Balay 
743*c3339decSBarry Smith   Logically Collective
7440c52818fSSatish Balay 
745d8d19677SJose E. Roman   Input Parameters:
7460c52818fSSatish Balay + ls - the TaoLineSearch context
7470c52818fSSatish Balay . func - the objective and gradient evaluation routine
7480c52818fSSatish Balay - ctx - the (optional) user-defined context for private data
7490c52818fSSatish Balay 
7500c52818fSSatish Balay   Calling sequence of func:
7510c52818fSSatish Balay $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);
7520c52818fSSatish Balay 
7530c52818fSSatish Balay + x - input vector
7540c52818fSSatish Balay . s - step direction
7550c52818fSSatish Balay . f - function value
7560c52818fSSatish Balay . gts - inner product of gradient and step direction vectors
7570c52818fSSatish Balay - ctx (optional) user-defined context
7580c52818fSSatish Balay 
7590c52818fSSatish Balay   Note: The gradient will still need to be computed at the end of the line
7600c52818fSSatish Balay   search, so you will still need to set a line search gradient evaluation
7610c52818fSSatish Balay   routine
7620c52818fSSatish Balay 
7630c52818fSSatish Balay   Note: Bounded line searches (those used in bounded optimization algorithms)
7640c52818fSSatish Balay   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
7650c52818fSSatish Balay   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()
7660c52818fSSatish Balay 
7670c52818fSSatish Balay   Level: advanced
7680c52818fSSatish Balay 
7690c52818fSSatish Balay   Note:
7700c52818fSSatish Balay   Some algorithms (lcl, gpcg) set their own objective routine for the
7710c52818fSSatish Balay   line search, application programmers should be wary of overriding the
7720c52818fSSatish Balay   default objective routine.
7730c52818fSSatish Balay 
774db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()`
7750c52818fSSatish Balay @*/
776d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void *), void *ctx)
777d71ae5a4SJacob Faibussowitsch {
7780c52818fSSatish Balay   PetscFunctionBegin;
7790c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
7800c52818fSSatish Balay   ls->ops->computeobjectiveandgts = func;
7810c52818fSSatish Balay   if (ctx) ls->userctx_funcgts = ctx;
7820c52818fSSatish Balay   ls->usegts         = PETSC_TRUE;
7830c52818fSSatish Balay   ls->usetaoroutines = PETSC_FALSE;
7840c52818fSSatish Balay   PetscFunctionReturn(0);
7850c52818fSSatish Balay }
7860c52818fSSatish Balay 
7870c52818fSSatish Balay /*@C
788441846f8SBarry Smith   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
789441846f8SBarry Smith   objective and gradient evaluation routines from the given Tao object.
7900c52818fSSatish Balay 
791*c3339decSBarry Smith   Logically Collective
7920c52818fSSatish Balay 
793d8d19677SJose E. Roman   Input Parameters:
7940c52818fSSatish Balay + ls - the TaoLineSearch context
795441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines
7960c52818fSSatish Balay 
7970c52818fSSatish Balay   Level: developer
7980c52818fSSatish Balay 
799db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`
8000c52818fSSatish Balay @*/
801d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
802d71ae5a4SJacob Faibussowitsch {
8030c52818fSSatish Balay   PetscFunctionBegin;
8040c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
805064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(ts, TAO_CLASSID, 2);
806441846f8SBarry Smith   ls->tao            = ts;
8070c52818fSSatish Balay   ls->usetaoroutines = PETSC_TRUE;
8080c52818fSSatish Balay   PetscFunctionReturn(0);
8090c52818fSSatish Balay }
8100c52818fSSatish Balay 
8110c52818fSSatish Balay /*@
8120c52818fSSatish Balay   TaoLineSearchComputeObjective - Computes the objective function value at a given point
8130c52818fSSatish Balay 
814*c3339decSBarry Smith   Collective
8150c52818fSSatish Balay 
8160c52818fSSatish Balay   Input Parameters:
8170c52818fSSatish Balay + ls - the TaoLineSearch context
8180c52818fSSatish Balay - x - input vector
8190c52818fSSatish Balay 
8200c52818fSSatish Balay   Output Parameter:
8210c52818fSSatish Balay . f - Objective value at X
8220c52818fSSatish Balay 
82395452b02SPatrick Sanan   Notes:
82495452b02SPatrick Sanan     TaoLineSearchComputeObjective() is typically used within line searches
8250c52818fSSatish Balay   so most users would not generally call this routine themselves.
8260c52818fSSatish Balay 
8270c52818fSSatish Balay   Level: developer
8280c52818fSSatish Balay 
829db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
8300c52818fSSatish Balay @*/
831d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
832d71ae5a4SJacob Faibussowitsch {
8330c52818fSSatish Balay   Vec       gdummy;
8340c52818fSSatish Balay   PetscReal gts;
835f06e3bfaSBarry Smith 
8360c52818fSSatish Balay   PetscFunctionBegin;
8370c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
8380c52818fSSatish Balay   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
839dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f, 3);
8400c52818fSSatish Balay   PetscCheckSameComm(ls, 1, x, 2);
8410c52818fSSatish Balay   if (ls->usetaoroutines) {
8429566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(ls->tao, x, f));
8430c52818fSSatish Balay   } else {
8443c859ba3SBarry 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");
8459566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
846792fecdfSBarry Smith     if (ls->ops->computeobjective) PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
847ef1023bdSBarry Smith     else if (ls->ops->computeobjectiveandgradient) {
8489566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &gdummy));
849792fecdfSBarry Smith       PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgradient)(ls, x, f, gdummy, ls->userctx_funcgrad));
8509566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&gdummy));
851792fecdfSBarry Smith     } else PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, &gts, ls->userctx_funcgts));
8529566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
8530c52818fSSatish Balay   }
8540c52818fSSatish Balay   ls->nfeval++;
8550c52818fSSatish Balay   PetscFunctionReturn(0);
8560c52818fSSatish Balay }
8570c52818fSSatish Balay 
8580c52818fSSatish Balay /*@
8590c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
8600c52818fSSatish Balay 
861*c3339decSBarry Smith   Collective
8620c52818fSSatish Balay 
8630c52818fSSatish Balay   Input Parameters:
8640c52818fSSatish Balay + ls - the TaoLineSearch context
8650c52818fSSatish Balay - x - input vector
8660c52818fSSatish Balay 
867d8d19677SJose E. Roman   Output Parameters:
8680c52818fSSatish Balay + f - Objective value at X
8690c52818fSSatish Balay - g - Gradient vector at X
8700c52818fSSatish Balay 
87195452b02SPatrick Sanan   Notes:
87295452b02SPatrick Sanan     TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
8730c52818fSSatish Balay   so most users would not generally call this routine themselves.
8740c52818fSSatish Balay 
8750c52818fSSatish Balay   Level: developer
8760c52818fSSatish Balay 
877db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
8780c52818fSSatish Balay @*/
879d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
880d71ae5a4SJacob Faibussowitsch {
8810c52818fSSatish Balay   PetscFunctionBegin;
8820c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
8830c52818fSSatish Balay   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
884dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f, 3);
8850c52818fSSatish Balay   PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
8860c52818fSSatish Balay   PetscCheckSameComm(ls, 1, x, 2);
8870c52818fSSatish Balay   PetscCheckSameComm(ls, 1, g, 4);
8880c52818fSSatish Balay   if (ls->usetaoroutines) {
8899566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(ls->tao, x, f, g));
8900c52818fSSatish Balay   } else {
8919566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
892792fecdfSBarry Smith     if (ls->ops->computeobjectiveandgradient) PetscCallBack("TaoLineSearch callback objective/gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, f, g, ls->userctx_funcgrad));
893ef1023bdSBarry Smith     else {
894792fecdfSBarry Smith       PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
895792fecdfSBarry Smith       PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
896362febeeSStefano Zampini     }
8979566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
8989566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
899f06e3bfaSBarry Smith   }
900fbe0838dSJason Sarich   ls->nfgeval++;
9010c52818fSSatish Balay   PetscFunctionReturn(0);
9020c52818fSSatish Balay }
9030c52818fSSatish Balay 
9040c52818fSSatish Balay /*@
9050c52818fSSatish Balay   TaoLineSearchComputeGradient - Computes the gradient of the objective function
9060c52818fSSatish Balay 
907*c3339decSBarry Smith   Collective
9080c52818fSSatish Balay 
9090c52818fSSatish Balay   Input Parameters:
9100c52818fSSatish Balay + ls - the TaoLineSearch context
9110c52818fSSatish Balay - x - input vector
9120c52818fSSatish Balay 
9130c52818fSSatish Balay   Output Parameter:
9140c52818fSSatish Balay . g - gradient vector
9150c52818fSSatish Balay 
91695452b02SPatrick Sanan   Notes:
91795452b02SPatrick Sanan     TaoComputeGradient() is typically used within line searches
9180c52818fSSatish Balay   so most users would not generally call this routine themselves.
9190c52818fSSatish Balay 
9200c52818fSSatish Balay   Level: developer
9210c52818fSSatish Balay 
922db781477SPatrick Sanan .seealso: `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()`
9230c52818fSSatish Balay @*/
924d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
925d71ae5a4SJacob Faibussowitsch {
9260c52818fSSatish Balay   PetscReal fdummy;
927f06e3bfaSBarry Smith 
9280c52818fSSatish Balay   PetscFunctionBegin;
9290c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
9300c52818fSSatish Balay   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
9310c52818fSSatish Balay   PetscValidHeaderSpecific(g, VEC_CLASSID, 3);
9320c52818fSSatish Balay   PetscCheckSameComm(ls, 1, x, 2);
9330c52818fSSatish Balay   PetscCheckSameComm(ls, 1, g, 3);
9340c52818fSSatish Balay   if (ls->usetaoroutines) {
9359566063dSJacob Faibussowitsch     PetscCall(TaoComputeGradient(ls->tao, x, g));
9360c52818fSSatish Balay   } else {
9379566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
938792fecdfSBarry Smith     if (ls->ops->computegradient) PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
93948a46eb9SPierre Jolivet     else PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, &fdummy, g, ls->userctx_funcgrad));
9409566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
9410c52818fSSatish Balay   }
9420c52818fSSatish Balay   ls->ngeval++;
9430c52818fSSatish Balay   PetscFunctionReturn(0);
9440c52818fSSatish Balay }
9450c52818fSSatish Balay 
9460c52818fSSatish Balay /*@
9470c52818fSSatish Balay   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point
9480c52818fSSatish Balay 
949*c3339decSBarry Smith   Collective
9500c52818fSSatish Balay 
9510c52818fSSatish Balay   Input Parameters:
9520c52818fSSatish Balay + ls - the TaoLineSearch context
9530c52818fSSatish Balay - x - input vector
9540c52818fSSatish Balay 
955d8d19677SJose E. Roman   Output Parameters:
9560c52818fSSatish Balay + f - Objective value at X
9570c52818fSSatish Balay - gts - inner product of gradient and step direction at X
9580c52818fSSatish Balay 
95995452b02SPatrick Sanan   Notes:
96095452b02SPatrick Sanan     TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
9610c52818fSSatish Balay   so most users would not generally call this routine themselves.
9620c52818fSSatish Balay 
9630c52818fSSatish Balay   Level: developer
9640c52818fSSatish Balay 
965db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
9660c52818fSSatish Balay @*/
967d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
968d71ae5a4SJacob Faibussowitsch {
9690c52818fSSatish Balay   PetscFunctionBegin;
9700c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
9710c52818fSSatish Balay   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
972dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f, 3);
973dadcf809SJacob Faibussowitsch   PetscValidRealPointer(gts, 4);
9740c52818fSSatish Balay   PetscCheckSameComm(ls, 1, x, 2);
9759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
976792fecdfSBarry Smith   PetscCallBack("TaoLineSearch callback objective/gts", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, gts, ls->userctx_funcgts));
9779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
9789566063dSJacob Faibussowitsch   PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
9790c52818fSSatish Balay   ls->nfeval++;
9800c52818fSSatish Balay   PetscFunctionReturn(0);
9810c52818fSSatish Balay }
9820c52818fSSatish Balay 
9830c52818fSSatish Balay /*@
9840c52818fSSatish Balay   TaoLineSearchGetSolution - Returns the solution to the line search
9850c52818fSSatish Balay 
986*c3339decSBarry Smith   Collective
9870c52818fSSatish Balay 
9880c52818fSSatish Balay   Input Parameter:
9890c52818fSSatish Balay . ls - the TaoLineSearch context
9900c52818fSSatish Balay 
991d8d19677SJose E. Roman   Output Parameters:
9920c52818fSSatish Balay + x - the new solution
9930c52818fSSatish Balay . f - the objective function value at x
9940c52818fSSatish Balay . g - the gradient at x
9950c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search
9960c52818fSSatish Balay - reason - the reason why the line search terminated
9970c52818fSSatish Balay 
9980c52818fSSatish Balay   reason will be set to one of:
9990c52818fSSatish Balay 
10000c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
10010c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
10020c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
10030c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
10040c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
10050c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
10060c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
10070c52818fSSatish Balay 
10080c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
10090c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason
10100c52818fSSatish Balay 
1011a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search
10120c52818fSSatish Balay 
10130c52818fSSatish Balay   Level: developer
10140c52818fSSatish Balay 
10150c52818fSSatish Balay @*/
1016d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1017d71ae5a4SJacob Faibussowitsch {
10180c52818fSSatish Balay   PetscFunctionBegin;
10190c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
10200c52818fSSatish Balay   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
1021dadcf809SJacob Faibussowitsch   PetscValidRealPointer(f, 3);
10220c52818fSSatish Balay   PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
10230c52818fSSatish Balay   PetscValidIntPointer(reason, 6);
10241baa6e33SBarry Smith   if (ls->new_x) PetscCall(VecCopy(ls->new_x, x));
10250c52818fSSatish Balay   *f = ls->new_f;
10261baa6e33SBarry Smith   if (ls->new_g) PetscCall(VecCopy(ls->new_g, g));
102797ab8e72SStefano Zampini   if (steplength) *steplength = ls->step;
10280c52818fSSatish Balay   *reason = ls->reason;
10290c52818fSSatish Balay   PetscFunctionReturn(0);
10300c52818fSSatish Balay }
10310c52818fSSatish Balay 
10320c52818fSSatish Balay /*@
10330c52818fSSatish Balay   TaoLineSearchGetStartingVector - Gets a the initial point of the line
10340c52818fSSatish Balay   search.
10350c52818fSSatish Balay 
10360c52818fSSatish Balay   Not Collective
10370c52818fSSatish Balay 
10380c52818fSSatish Balay   Input Parameter:
10390c52818fSSatish Balay . ls - the TaoLineSearch context
10400c52818fSSatish Balay 
10410c52818fSSatish Balay   Output Parameter:
10420c52818fSSatish Balay . x - The initial point of the line search
10430c52818fSSatish Balay 
10440c52818fSSatish Balay   Level: intermediate
10450c52818fSSatish Balay @*/
1046d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1047d71ae5a4SJacob Faibussowitsch {
10480c52818fSSatish Balay   PetscFunctionBegin;
10490c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
105097ab8e72SStefano Zampini   if (x) *x = ls->start_x;
10510c52818fSSatish Balay   PetscFunctionReturn(0);
10520c52818fSSatish Balay }
10530c52818fSSatish Balay 
10540c52818fSSatish Balay /*@
10550c52818fSSatish Balay   TaoLineSearchGetStepDirection - Gets the step direction of the line
10560c52818fSSatish Balay   search.
10570c52818fSSatish Balay 
10580c52818fSSatish Balay   Not Collective
10590c52818fSSatish Balay 
10600c52818fSSatish Balay   Input Parameter:
10610c52818fSSatish Balay . ls - the TaoLineSearch context
10620c52818fSSatish Balay 
10630c52818fSSatish Balay   Output Parameter:
10640c52818fSSatish Balay . s - the step direction of the line search
10650c52818fSSatish Balay 
10660c52818fSSatish Balay   Level: advanced
10670c52818fSSatish Balay @*/
1068d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1069d71ae5a4SJacob Faibussowitsch {
10700c52818fSSatish Balay   PetscFunctionBegin;
10710c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
107297ab8e72SStefano Zampini   if (s) *s = ls->stepdirection;
10730c52818fSSatish Balay   PetscFunctionReturn(0);
10740c52818fSSatish Balay }
10750c52818fSSatish Balay 
10760c52818fSSatish Balay /*@
10770c52818fSSatish Balay   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.
10780c52818fSSatish Balay 
10790c52818fSSatish Balay   Not Collective
10800c52818fSSatish Balay 
10810c52818fSSatish Balay   Input Parameter:
10820c52818fSSatish Balay . ls - the TaoLineSearch context
10830c52818fSSatish Balay 
10840c52818fSSatish Balay   Output Parameter:
10850c52818fSSatish Balay . f - the objective value at the full step length
10860c52818fSSatish Balay 
10870c52818fSSatish Balay   Level: developer
10880c52818fSSatish Balay @*/
10890c52818fSSatish Balay 
1090d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1091d71ae5a4SJacob Faibussowitsch {
10920c52818fSSatish Balay   PetscFunctionBegin;
10930c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
10940c52818fSSatish Balay   *f_fullstep = ls->f_fullstep;
10950c52818fSSatish Balay   PetscFunctionReturn(0);
10960c52818fSSatish Balay }
10970c52818fSSatish Balay 
10980c52818fSSatish Balay /*@
10990c52818fSSatish Balay   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.
11000c52818fSSatish Balay 
1101*c3339decSBarry Smith   Logically Collective
11020c52818fSSatish Balay 
11030c52818fSSatish Balay   Input Parameters:
11040c52818fSSatish Balay + ls - the TaoLineSearch context
11050c52818fSSatish Balay . xl  - vector of lower bounds
11060c52818fSSatish Balay - xu  - vector of upper bounds
11070c52818fSSatish Balay 
11080c52818fSSatish Balay   Note: If the variable bounds are not set with this routine, then
1109e270355aSBarry Smith   PETSC_NINFINITY and PETSC_INFINITY are assumed
11100c52818fSSatish Balay 
11110c52818fSSatish Balay   Level: beginner
11120c52818fSSatish Balay 
1113db781477SPatrick Sanan .seealso: `TaoSetVariableBounds()`, `TaoLineSearchCreate()`
11140c52818fSSatish Balay @*/
1115d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls, Vec xl, Vec xu)
1116d71ae5a4SJacob Faibussowitsch {
11170c52818fSSatish Balay   PetscFunctionBegin;
11180c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
111976be6f4fSStefano Zampini   if (xl) PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
112076be6f4fSStefano Zampini   if (xu) PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
11219566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
11229566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
11239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->lower));
11249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ls->upper));
11250c52818fSSatish Balay   ls->lower   = xl;
11260c52818fSSatish Balay   ls->upper   = xu;
112776be6f4fSStefano Zampini   ls->bounded = (PetscBool)(xl || xu);
11280c52818fSSatish Balay   PetscFunctionReturn(0);
11290c52818fSSatish Balay }
11300c52818fSSatish Balay 
11310c52818fSSatish Balay /*@
11320c52818fSSatish Balay   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
11330c52818fSSatish Balay   search.  If this value is not set then 1.0 is assumed.
11340c52818fSSatish Balay 
1135*c3339decSBarry Smith   Logically Collective
11360c52818fSSatish Balay 
11370c52818fSSatish Balay   Input Parameters:
11380c52818fSSatish Balay + ls - the TaoLineSearch context
11390c52818fSSatish Balay - s - the initial step size
11400c52818fSSatish Balay 
11410c52818fSSatish Balay   Level: intermediate
11420c52818fSSatish Balay 
1143db781477SPatrick Sanan .seealso: `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()`
11440c52818fSSatish Balay @*/
1145d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls, PetscReal s)
1146d71ae5a4SJacob Faibussowitsch {
11470c52818fSSatish Balay   PetscFunctionBegin;
11480c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1149a39c8e28SStefano Zampini   PetscValidLogicalCollectiveReal(ls, s, 2);
11500c52818fSSatish Balay   ls->initstep = s;
11510c52818fSSatish Balay   PetscFunctionReturn(0);
11520c52818fSSatish Balay }
11530c52818fSSatish Balay 
11540c52818fSSatish Balay /*@
11550c52818fSSatish Balay   TaoLineSearchGetStepLength - Get the current step length
11560c52818fSSatish Balay 
11570c52818fSSatish Balay   Not Collective
11580c52818fSSatish Balay 
11590c52818fSSatish Balay   Input Parameters:
11600c52818fSSatish Balay . ls - the TaoLineSearch context
11610c52818fSSatish Balay 
11627a7aea1fSJed Brown   Output Parameters:
11630c52818fSSatish Balay . s - the current step length
11640c52818fSSatish Balay 
11650c52818fSSatish Balay   Level: beginner
11660c52818fSSatish Balay 
1167db781477SPatrick Sanan .seealso: `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()`
11680c52818fSSatish Balay @*/
1169d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls, PetscReal *s)
1170d71ae5a4SJacob Faibussowitsch {
11710c52818fSSatish Balay   PetscFunctionBegin;
11720c52818fSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
11730c52818fSSatish Balay   *s = ls->step;
11740c52818fSSatish Balay   PetscFunctionReturn(0);
11750c52818fSSatish Balay }
11760c52818fSSatish Balay 
11770ebee16dSLisandro Dalcin /*@C
11780c52818fSSatish Balay    TaoLineSearchRegister - Adds a line-search algorithm to the registry
11790c52818fSSatish Balay 
11800c52818fSSatish Balay    Not collective
11810c52818fSSatish Balay 
11820c52818fSSatish Balay    Input Parameters:
11830c52818fSSatish Balay +  sname - name of a new user-defined solver
11840c52818fSSatish Balay -  func - routine to Create method context
11850c52818fSSatish Balay 
11860c52818fSSatish Balay    Notes:
11870c52818fSSatish Balay    TaoLineSearchRegister() may be called multiple times to add several user-defined solvers.
11880c52818fSSatish Balay 
11890c52818fSSatish Balay    Sample usage:
11900c52818fSSatish Balay .vb
11910c52818fSSatish Balay    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
11920c52818fSSatish Balay .ve
11930c52818fSSatish Balay 
11940c52818fSSatish Balay    Then, your solver can be chosen with the procedural interface via
11950c52818fSSatish Balay $     TaoLineSearchSetType(ls,"my_linesearch")
11960c52818fSSatish Balay    or at runtime via the option
11970c52818fSSatish Balay $     -tao_ls_type my_linesearch
11980c52818fSSatish Balay 
11990c52818fSSatish Balay    Level: developer
12000c52818fSSatish Balay 
12010ebee16dSLisandro Dalcin @*/
1202d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1203d71ae5a4SJacob Faibussowitsch {
12040c52818fSSatish Balay   PetscFunctionBegin;
12059566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchInitializePackage());
12069566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func));
12070c52818fSSatish Balay   PetscFunctionReturn(0);
12080c52818fSSatish Balay }
12090c52818fSSatish Balay 
12100c52818fSSatish Balay /*@C
12110c52818fSSatish Balay    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
12120c52818fSSatish Balay    for all TaoLineSearch options in the database.
12130c52818fSSatish Balay 
1214*c3339decSBarry Smith    Collective
12150c52818fSSatish Balay 
12160c52818fSSatish Balay    Input Parameters:
12170c52818fSSatish Balay +  ls - the TaoLineSearch solver context
12180c52818fSSatish Balay -  prefix - the prefix string to prepend to all line search requests
12190c52818fSSatish Balay 
12200c52818fSSatish Balay    Notes:
12210c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
12220c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
12230c52818fSSatish Balay 
12240c52818fSSatish Balay    Level: advanced
12250c52818fSSatish Balay 
1226db781477SPatrick Sanan .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
12270c52818fSSatish Balay @*/
1228d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1229d71ae5a4SJacob Faibussowitsch {
1230f06e3bfaSBarry Smith   return PetscObjectAppendOptionsPrefix((PetscObject)ls, p);
12310c52818fSSatish Balay }
12320c52818fSSatish Balay 
12330c52818fSSatish Balay /*@C
12340c52818fSSatish Balay   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
12350c52818fSSatish Balay   TaoLineSearch options in the database
12360c52818fSSatish Balay 
12370c52818fSSatish Balay   Not Collective
12380c52818fSSatish Balay 
12390c52818fSSatish Balay   Input Parameters:
12400c52818fSSatish Balay . ls - the TaoLineSearch context
12410c52818fSSatish Balay 
12420c52818fSSatish Balay   Output Parameters:
12430c52818fSSatish Balay . prefix - pointer to the prefix string used is returned
12440c52818fSSatish Balay 
124595452b02SPatrick Sanan   Notes:
124695452b02SPatrick Sanan     On the fortran side, the user should pass in a string 'prefix' of
12470c52818fSSatish Balay   sufficient length to hold the prefix.
12480c52818fSSatish Balay 
12490c52818fSSatish Balay   Level: advanced
12500c52818fSSatish Balay 
1251db781477SPatrick Sanan .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()`
12520c52818fSSatish Balay @*/
1253d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1254d71ae5a4SJacob Faibussowitsch {
1255f06e3bfaSBarry Smith   return PetscObjectGetOptionsPrefix((PetscObject)ls, p);
12560c52818fSSatish Balay }
12570c52818fSSatish Balay 
12580c52818fSSatish Balay /*@C
12590c52818fSSatish Balay    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
12600c52818fSSatish Balay    TaoLineSearch options in the database.
12610c52818fSSatish Balay 
1262*c3339decSBarry Smith    Logically Collective
12630c52818fSSatish Balay 
12640c52818fSSatish Balay    Input Parameters:
12650c52818fSSatish Balay +  ls - the TaoLineSearch context
12660c52818fSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
12670c52818fSSatish Balay 
12680c52818fSSatish Balay    Notes:
12690c52818fSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
12700c52818fSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
12710c52818fSSatish Balay 
12720c52818fSSatish Balay    For example, to distinguish between the runtime options for two
12730c52818fSSatish Balay    different line searches, one could call
12740c52818fSSatish Balay .vb
12750c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
12760c52818fSSatish Balay       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
12770c52818fSSatish Balay .ve
12780c52818fSSatish Balay 
12790c52818fSSatish Balay    This would enable use of different options for each system, such as
12800c52818fSSatish Balay .vb
12810c52818fSSatish Balay       -sys1_tao_ls_type mt
12820c52818fSSatish Balay       -sys2_tao_ls_type armijo
12830c52818fSSatish Balay .ve
12840c52818fSSatish Balay 
12850c52818fSSatish Balay    Level: advanced
12860c52818fSSatish Balay 
1287db781477SPatrick Sanan .seealso: `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
12880c52818fSSatish Balay @*/
12890c52818fSSatish Balay 
1290d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1291d71ae5a4SJacob Faibussowitsch {
1292f06e3bfaSBarry Smith   return PetscObjectSetOptionsPrefix((PetscObject)ls, p);
12930c52818fSSatish Balay }
1294