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 11*1404853cSMatthew G. Knepley const char *const TaoLineSearchConvergedReasons_Shifted[] = {"FAILED_ASCENT", "FAILED_BADPARAMETER", "FAILED_INFORNAN", "CONTINUE_ITERATING", "SUCCESS", "SUCCESS_USER", "HALTED_OTHER", "HALTED_MAXFCN", "HALTED_UPPERBOUND", "HALTED_LOWERBOUND", "HALTED_RTOL", "HALTED_USER", "TaoLineSearchConvergedReason", "TAOLINESEARCH_", NULL}; 12*1404853cSMatthew G. Knepley const char *const *TaoLineSearchConvergedReasons = TaoLineSearchConvergedReasons_Shifted + 3; 13*1404853cSMatthew G. Knepley 14ffeef943SBarry Smith /*@ 1520f4b53cSBarry Smith TaoLineSearchViewFromOptions - View a `TaoLineSearch` object based on values in the options database 16fe2efc57SMark 17c3339decSBarry Smith Collective 18fe2efc57SMark 19fe2efc57SMark Input Parameters: 2020f4b53cSBarry Smith + A - the `Tao` context 21736c3998SJose E. Roman . obj - Optional object 22736c3998SJose E. Roman - name - command line option 23fe2efc57SMark 24fe2efc57SMark Level: intermediate 2520f4b53cSBarry Smith 2620f4b53cSBarry Smith Note: 2720f4b53cSBarry Smith See `PetscObjectViewFromOptions()` for available viewer options 2820f4b53cSBarry Smith 291cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchView()`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()` 30fe2efc57SMark @*/ 31d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A, PetscObject obj, const char name[]) 32d71ae5a4SJacob Faibussowitsch { 33fe2efc57SMark PetscFunctionBegin; 34fe2efc57SMark PetscValidHeaderSpecific(A, TAOLINESEARCH_CLASSID, 1); 359566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37fe2efc57SMark } 38fe2efc57SMark 39ffeef943SBarry Smith /*@ 4020f4b53cSBarry Smith TaoLineSearchView - Prints information about the `TaoLineSearch` 410c52818fSSatish Balay 42c3339decSBarry Smith Collective 430c52818fSSatish Balay 440c52818fSSatish Balay Input Parameters: 4520f4b53cSBarry Smith + ls - the `TaoLineSearch` context 460c52818fSSatish Balay - viewer - visualization context 470c52818fSSatish Balay 480c52818fSSatish Balay Options Database Key: 4920f4b53cSBarry Smith . -tao_ls_view - Calls `TaoLineSearchView()` at the end of each line search 5020f4b53cSBarry Smith 5120f4b53cSBarry Smith Level: beginner 520c52818fSSatish Balay 530c52818fSSatish Balay Notes: 540c52818fSSatish Balay The available visualization contexts include 5520f4b53cSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 5620f4b53cSBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 570c52818fSSatish Balay output where only the first processor opens 580c52818fSSatish Balay the file. All other processors send their 590c52818fSSatish Balay data to the first processor to print. 600c52818fSSatish Balay 611cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `PetscViewerASCIIOpen()`, `TaoLineSearchViewFromOptions()` 620c52818fSSatish Balay @*/ 63d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 64d71ae5a4SJacob Faibussowitsch { 650c52818fSSatish Balay PetscBool isascii, isstring; 66dedfbcbeSJed Brown TaoLineSearchType type; 67f06e3bfaSBarry Smith 680c52818fSSatish Balay PetscFunctionBegin; 690c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 7048a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer)); 710c52818fSSatish Balay PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 720c52818fSSatish Balay PetscCheckSameComm(ls, 1, viewer, 2); 730c52818fSSatish Balay 749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 760c52818fSSatish Balay if (isascii) { 779566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer)); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 79dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, view, viewer); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 8263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "maximum function evaluations=%" PetscInt_FMT "\n", ls->max_funcs)); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "tolerances: ftol=%g, rtol=%g, gtol=%g\n", (double)ls->ftol, (double)ls->rtol, (double)ls->gtol)); 8463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT "\n", ls->nfeval)); 8563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT "\n", ls->ngeval)); 8663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT "\n", ls->nfgeval)); 870c52818fSSatish Balay 8848a46eb9SPierre Jolivet if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(viewer, "using variable bounds\n")); 89*1404853cSMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "Termination reason: %s\n", TaoLineSearchConvergedReasons[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 } 953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 960c52818fSSatish Balay } 970c52818fSSatish Balay 980c52818fSSatish Balay /*@C 9920f4b53cSBarry Smith TaoLineSearchCreate - Creates a `TaoLineSearch` object. Algorithms in `Tao` that use 10020f4b53cSBarry Smith line-searches will automatically create one so this all is rarely needed 1010c52818fSSatish Balay 102d083f849SBarry Smith Collective 1030c52818fSSatish Balay 1040c52818fSSatish Balay Input Parameter: 1050c52818fSSatish Balay . comm - MPI communicator 1060c52818fSSatish Balay 1070c52818fSSatish Balay Output Parameter: 10820f4b53cSBarry Smith . newls - the new `TaoLineSearch` context 1090c52818fSSatish Balay 11020f4b53cSBarry Smith Options Database Key: 11120f4b53cSBarry Smith . -tao_ls_type - select which method `Tao` should use 1120c52818fSSatish Balay 11320f4b53cSBarry Smith Level: developer 1140c52818fSSatish Balay 1151cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()` 1160c52818fSSatish Balay @*/ 117d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 118d71ae5a4SJacob Faibussowitsch { 1190c52818fSSatish Balay TaoLineSearch ls; 1200c52818fSSatish Balay 1210c52818fSSatish Balay PetscFunctionBegin; 1224f572ea9SToby Isaac PetscAssertPointer(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; 1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1400c52818fSSatish Balay } 1410c52818fSSatish Balay 1420c52818fSSatish Balay /*@ 1430c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 14420f4b53cSBarry Smith of a `TaoLineSearch` 1450c52818fSSatish Balay 146c3339decSBarry Smith Collective 1470c52818fSSatish Balay 14820f4b53cSBarry Smith Input Parameter: 14920f4b53cSBarry Smith . ls - the `TaoLineSearch` context 1500c52818fSSatish Balay 1510c52818fSSatish Balay Level: developer 1520c52818fSSatish Balay 15320f4b53cSBarry Smith Note: 15420f4b53cSBarry Smith The user will not need to explicitly call `TaoLineSearchSetUp()`, as it will 15520f4b53cSBarry Smith automatically be called in `TaoLineSearchSolve()`. However, if the user 15620f4b53cSBarry Smith desires to call it explicitly, it should come after `TaoLineSearchCreate()` 15720f4b53cSBarry Smith but before `TaoLineSearchApply()`. 1580c52818fSSatish Balay 1591cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()` 16020f4b53cSBarry Smith @*/ 161d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 162d71ae5a4SJacob Faibussowitsch { 1638caf6e8cSBarry Smith const char *default_type = TAOLINESEARCHMT; 1640c52818fSSatish Balay PetscBool flg; 165f06e3bfaSBarry Smith 1660c52818fSSatish Balay PetscFunctionBegin; 1670c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 1683ba16761SJacob Faibussowitsch if (ls->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 16948a46eb9SPierre Jolivet if (!((PetscObject)ls)->type_name) PetscCall(TaoLineSearchSetType(ls, default_type)); 170dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, setup); 1710c52818fSSatish Balay if (ls->usetaoroutines) { 1729566063dSJacob Faibussowitsch PetscCall(TaoIsObjectiveDefined(ls->tao, &flg)); 1730c52818fSSatish Balay ls->hasobjective = flg; 1749566063dSJacob Faibussowitsch PetscCall(TaoIsGradientDefined(ls->tao, &flg)); 1750c52818fSSatish Balay ls->hasgradient = flg; 1769566063dSJacob Faibussowitsch PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao, &flg)); 1770c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 1780c52818fSSatish Balay } else { 1790c52818fSSatish Balay if (ls->ops->computeobjective) { 1800c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 1810c52818fSSatish Balay } else { 1820c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 1830c52818fSSatish Balay } 1840c52818fSSatish Balay if (ls->ops->computegradient) { 1850c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 1860c52818fSSatish Balay } else { 1870c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 1880c52818fSSatish Balay } 1890c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 1900c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 1910c52818fSSatish Balay } else { 1920c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 1930c52818fSSatish Balay } 1940c52818fSSatish Balay } 1950c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1970c52818fSSatish Balay } 1980c52818fSSatish Balay 1990c52818fSSatish Balay /*@ 2000c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 20120f4b53cSBarry Smith from one `TaoLineSearchApply()` to the next. This function resets this 2020c52818fSSatish Balay state information. 2030c52818fSSatish Balay 204c3339decSBarry Smith Collective 2050c52818fSSatish Balay 2060c52818fSSatish Balay Input Parameter: 20720f4b53cSBarry Smith . ls - the `TaoLineSearch` context 2080c52818fSSatish Balay 2090c52818fSSatish Balay Level: developer 2100c52818fSSatish Balay 2111cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()` 2120c52818fSSatish Balay @*/ 213d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 214d71ae5a4SJacob Faibussowitsch { 2150c52818fSSatish Balay PetscFunctionBegin; 2160c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 217dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, reset); 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2190c52818fSSatish Balay } 220f06e3bfaSBarry Smith 2210c52818fSSatish Balay /*@ 22220f4b53cSBarry Smith TaoLineSearchDestroy - Destroys the `TaoLineSearch` context that was created with 22320f4b53cSBarry Smith `TaoLineSearchCreate()` 2240c52818fSSatish Balay 225c3339decSBarry Smith Collective 2260c52818fSSatish Balay 2277a7aea1fSJed Brown Input Parameter: 22820f4b53cSBarry Smith . ls - the `TaoLineSearch` context 2290c52818fSSatish Balay 23020f4b53cSBarry Smith Level: developer 2310c52818fSSatish Balay 232e056e8ceSJacob Faibussowitsch .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSolve()` 2330c52818fSSatish Balay @*/ 234d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 235d71ae5a4SJacob Faibussowitsch { 2360c52818fSSatish Balay PetscFunctionBegin; 2373ba16761SJacob Faibussowitsch if (!*ls) PetscFunctionReturn(PETSC_SUCCESS); 2380c52818fSSatish Balay PetscValidHeaderSpecific(*ls, TAOLINESEARCH_CLASSID, 1); 2399371c9d4SSatish Balay if (--((PetscObject)*ls)->refct > 0) { 2409371c9d4SSatish Balay *ls = NULL; 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2429371c9d4SSatish Balay } 2439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->stepdirection)); 2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->start_x)); 2459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->upper)); 2469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->lower)); 247213acdd3SPierre Jolivet PetscTryTypeMethod(*ls, destroy); 24848a46eb9SPierre Jolivet if ((*ls)->usemonitor) PetscCall(PetscViewerDestroy(&(*ls)->viewer)); 2499566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(ls)); 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2510c52818fSSatish Balay } 2520c52818fSSatish Balay 2530c52818fSSatish Balay /*@ 25420f4b53cSBarry Smith TaoLineSearchApply - Performs a line-search in a given step direction. 25520f4b53cSBarry Smith Criteria for acceptable step length depends on the line-search algorithm chosen 2560c52818fSSatish Balay 257c3339decSBarry Smith Collective 2580c52818fSSatish Balay 2590c52818fSSatish Balay Input Parameters: 26020f4b53cSBarry Smith + ls - the `TaoLineSearch` context 2610c52818fSSatish Balay - s - search direction 2620c52818fSSatish Balay 2630c52818fSSatish Balay Output Parameters: 26420f4b53cSBarry Smith + x - On input the current solution, on output `x` contains the new solution determined by the line search 265f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution 26620f4b53cSBarry Smith . g - On input the gradient evaluated at `x`, on output contains the gradient at new solution 267f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x) 26820f4b53cSBarry Smith - reason - `TaoLineSearchConvergedReason` reason why the line-search stopped 26920f4b53cSBarry Smith 27020f4b53cSBarry Smith Level: advanced 2710c52818fSSatish Balay 27297bb3fdcSJose E. Roman Notes: 27320f4b53cSBarry Smith The algorithm developer must set up the `TaoLineSearch` with calls to 27420f4b53cSBarry Smith `TaoLineSearchSetObjectiveRoutine()` and `TaoLineSearchSetGradientRoutine()`, 27520f4b53cSBarry Smith `TaoLineSearchSetObjectiveAndGradientRoutine()`, or `TaoLineSearchUseTaoRoutines()`. 27620f4b53cSBarry Smith The latter is done automatically by default and thus requires no user input. 2770c52818fSSatish Balay 2780c52818fSSatish Balay You may or may not need to follow this with a call to 27920f4b53cSBarry Smith `TaoAddLineSearchCounts()`, depending on whether you want these 2800c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 2810c52818fSSatish Balay 2821cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearchConvergedReason`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, 28320f4b53cSBarry Smith `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()` 2840c52818fSSatish Balay @*/ 285d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 286d71ae5a4SJacob Faibussowitsch { 2870c52818fSSatish Balay PetscInt low1, low2, low3, high1, high2, high3; 2880c52818fSSatish Balay 2890c52818fSSatish Balay PetscFunctionBegin; 2900c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 2910c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 2924f572ea9SToby Isaac PetscAssertPointer(f, 3); 2930c52818fSSatish Balay PetscValidHeaderSpecific(g, VEC_CLASSID, 4); 2940c52818fSSatish Balay PetscValidHeaderSpecific(s, VEC_CLASSID, 5); 2954f572ea9SToby Isaac PetscAssertPointer(reason, 7); 2960c52818fSSatish Balay PetscCheckSameComm(ls, 1, x, 2); 2970c52818fSSatish Balay PetscCheckSameTypeAndComm(x, 2, g, 4); 2980c52818fSSatish Balay PetscCheckSameTypeAndComm(x, 2, s, 5); 2999566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low1, &high1)); 3009566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(g, &low2, &high2)); 3019566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(s, &low3, &high3)); 3023c859ba3SBarry Smith PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible vector local lengths"); 3030c52818fSSatish Balay 30497ab8e72SStefano Zampini *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s)); 3069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->stepdirection)); 307050fc7a3SBarry Smith ls->stepdirection = s; 3080c52818fSSatish Balay 3099566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetUp(ls)); 3100c52818fSSatish Balay ls->nfeval = 0; 3110c52818fSSatish Balay ls->ngeval = 0; 3120c52818fSSatish Balay ls->nfgeval = 0; 3130c52818fSSatish Balay /* Check parameter values */ 3140c52818fSSatish Balay if (ls->ftol < 0.0) { 3159566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: ftol (%g) < 0\n", (double)ls->ftol)); 3160c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3170c52818fSSatish Balay } 3180c52818fSSatish Balay if (ls->rtol < 0.0) { 3199566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: rtol (%g) < 0\n", (double)ls->rtol)); 3200c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3210c52818fSSatish Balay } 3220c52818fSSatish Balay if (ls->gtol < 0.0) { 3239566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: gtol (%g) < 0\n", (double)ls->gtol)); 3240c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3250c52818fSSatish Balay } 3260c52818fSSatish Balay if (ls->stepmin < 0.0) { 3279566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) < 0\n", (double)ls->stepmin)); 3280c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3290c52818fSSatish Balay } 3300c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 3319566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n", (double)ls->stepmin, (double)ls->stepmax)); 3320c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3330c52818fSSatish Balay } 3340c52818fSSatish Balay if (ls->max_funcs < 0) { 3359566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n", ls->max_funcs)); 3360c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3370c52818fSSatish Balay } 3380c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 3399566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Initial Line Search Function Value is Inf or Nan (%g)\n", (double)*f)); 3400c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_INFORNAN; 3410c52818fSSatish Balay } 3420c52818fSSatish Balay 3439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)x)); 3449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->start_x)); 3450c52818fSSatish Balay ls->start_x = x; 346050fc7a3SBarry Smith 3479566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply, ls, 0, 0, 0)); 348dbbe0bcdSBarry Smith PetscUseTypeMethod(ls, apply, x, f, g, s); 3499566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0, 0, 0)); 3500c52818fSSatish Balay *reason = ls->reason; 3510c52818fSSatish Balay ls->new_f = *f; 3520c52818fSSatish Balay 35397ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 3540c52818fSSatish Balay 3559566063dSJacob Faibussowitsch PetscCall(TaoLineSearchViewFromOptions(ls, NULL, "-tao_ls_view")); 3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3570c52818fSSatish Balay } 3580c52818fSSatish Balay 359cc4c1da9SBarry Smith /*@ 3600c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 3610c52818fSSatish Balay 362c3339decSBarry Smith Collective 3630c52818fSSatish Balay 3640c52818fSSatish Balay Input Parameters: 36520f4b53cSBarry Smith + ls - the `TaoLineSearch` context 36620f4b53cSBarry Smith - type - the `TaoLineSearchType` selection 3670c52818fSSatish Balay 36820f4b53cSBarry Smith Options Database Key: 36920f4b53cSBarry Smith . -tao_ls_type <type> - select which method Tao should use at runtime 3700c52818fSSatish Balay 3710c52818fSSatish Balay Level: beginner 3720c52818fSSatish Balay 3731cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchCreate()`, `TaoLineSearchGetType()`, 37420f4b53cSBarry Smith `TaoLineSearchApply()` 3750c52818fSSatish Balay @*/ 376d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 377d71ae5a4SJacob Faibussowitsch { 3780c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 3790c52818fSSatish Balay PetscBool flg; 3800c52818fSSatish Balay 3810c52818fSSatish Balay PetscFunctionBegin; 3820c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 3834f572ea9SToby Isaac PetscAssertPointer(type, 2); 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg)); 3853ba16761SJacob Faibussowitsch if (flg) PetscFunctionReturn(PETSC_SUCCESS); 3860c52818fSSatish Balay 387835f2295SStefano Zampini PetscCall(PetscFunctionListFind(TaoLineSearchList, type, &r)); 3883c859ba3SBarry Smith PetscCheck(r, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested TaoLineSearch type %s", type); 389dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, destroy); 3900c52818fSSatish Balay ls->max_funcs = 30; 3910c52818fSSatish Balay ls->ftol = 0.0001; 3920c52818fSSatish Balay ls->gtol = 0.9; 3936f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 3946f4723b1SBarry Smith ls->rtol = 1.0e-5; 3956f4723b1SBarry Smith #else 3960c52818fSSatish Balay ls->rtol = 1.0e-10; 3976f4723b1SBarry Smith #endif 3980c52818fSSatish Balay ls->stepmin = 1.0e-20; 3990c52818fSSatish Balay ls->stepmax = 1.0e+20; 4000c52818fSSatish Balay 4010c52818fSSatish Balay ls->nfeval = 0; 4020c52818fSSatish Balay ls->ngeval = 0; 4030c52818fSSatish Balay ls->nfgeval = 0; 40483c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 40583c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 40683c8fe1dSLisandro Dalcin ls->ops->view = NULL; 40783c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 40883c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 4090c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4109566063dSJacob Faibussowitsch PetscCall((*r)(ls)); 4119566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type)); 4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4130c52818fSSatish Balay } 4140c52818fSSatish Balay 4152a0dac07SAlp Dener /*@C 41646091a0eSPierre Jolivet TaoLineSearchMonitor - Monitor the line search steps. This routine will output the 4172a0dac07SAlp Dener iteration number, step length, and function value before calling the implementation 4182a0dac07SAlp Dener specific monitor. 4192a0dac07SAlp Dener 4202a0dac07SAlp Dener Input Parameters: 42120f4b53cSBarry Smith + ls - the `TaoLineSearch` context 4222a0dac07SAlp Dener . its - the current iterate number (>=0) 4232a0dac07SAlp Dener . f - the current objective function value 4242a0dac07SAlp Dener - step - the step length 4252a0dac07SAlp Dener 4262a0dac07SAlp Dener Options Database Key: 4272a0dac07SAlp Dener . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output 4282a0dac07SAlp Dener 4292a0dac07SAlp Dener Level: developer 4302a0dac07SAlp Dener 43120f4b53cSBarry Smith .seealso: `TaoLineSearch` 4322a0dac07SAlp Dener @*/ 433d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step) 434d71ae5a4SJacob Faibussowitsch { 4352a0dac07SAlp Dener PetscInt tabs; 4362a0dac07SAlp Dener 4372a0dac07SAlp Dener PetscFunctionBegin; 4382a0dac07SAlp Dener PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 4392a0dac07SAlp Dener if (ls->usemonitor) { 4409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs)); 4419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel)); 44263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its)); 4439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f)); 4449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step)); 4452a0dac07SAlp Dener if (ls->ops->monitor && its > 0) { 4469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3)); 447dbbe0bcdSBarry Smith PetscUseTypeMethod(ls, monitor); 4482a0dac07SAlp Dener } 4499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs)); 4502a0dac07SAlp Dener } 4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4522a0dac07SAlp Dener } 4532a0dac07SAlp Dener 4540c52818fSSatish Balay /*@ 45520f4b53cSBarry Smith TaoLineSearchSetFromOptions - Sets various `TaoLineSearch` parameters from user 4560c52818fSSatish Balay options. 4570c52818fSSatish Balay 458c3339decSBarry Smith Collective 4590c52818fSSatish Balay 46001d2d390SJose E. Roman Input Parameter: 46120f4b53cSBarry Smith . ls - the `TaoLineSearch` context 4620c52818fSSatish Balay 4630c52818fSSatish Balay Options Database Keys: 46420f4b53cSBarry Smith + -tao_ls_type <type> - The algorithm that `TaoLineSearch` uses (more-thuente, gpcg, unit) 4650c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 4660c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 4670c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 468a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - initial steplength allowed 4690c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 4700c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 4710c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 4720c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 4730c52818fSSatish Balay 4740c52818fSSatish Balay Level: beginner 47520f4b53cSBarry Smith 47620f4b53cSBarry Smith .seealso: `TaoLineSearch` 4770c52818fSSatish Balay @*/ 478d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 479d71ae5a4SJacob Faibussowitsch { 4808caf6e8cSBarry Smith const char *default_type = TAOLINESEARCHMT; 4812a0dac07SAlp Dener char type[256], monfilename[PETSC_MAX_PATH_LEN]; 4822a0dac07SAlp Dener PetscViewer monviewer; 4830c52818fSSatish Balay PetscBool flg; 484f06e3bfaSBarry Smith 4850c52818fSSatish Balay PetscFunctionBegin; 4860c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 487d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)ls); 488ad540459SPierre Jolivet if (((PetscObject)ls)->type_name) default_type = ((PetscObject)ls)->type_name; 4890c52818fSSatish Balay /* Check for type from options */ 4909566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-tao_ls_type", "Tao Line Search type", "TaoLineSearchSetType", TaoLineSearchList, default_type, type, 256, &flg)); 4910c52818fSSatish Balay if (flg) { 4929566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls, type)); 4930c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 4949566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls, default_type)); 4950c52818fSSatish Balay } 4960c52818fSSatish Balay 4979566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_max_funcs", "max function evals in line search", "", ls->max_funcs, &ls->max_funcs, NULL)); 4989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_ftol", "tol for sufficient decrease", "", ls->ftol, &ls->ftol, NULL)); 4999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_gtol", "tol for curvature condition", "", ls->gtol, &ls->gtol, NULL)); 5009566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_rtol", "relative tol for acceptable step", "", ls->rtol, &ls->rtol, NULL)); 5019566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_stepmin", "lower bound for step", "", ls->stepmin, &ls->stepmin, NULL)); 5029566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_stepmax", "upper bound for step", "", ls->stepmax, &ls->stepmax, NULL)); 503a39c8e28SStefano Zampini PetscCall(PetscOptionsReal("-tao_ls_stepinit", "initial step", "", ls->initstep, &ls->initstep, NULL)); 5049566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-tao_ls_monitor", "enable the basic monitor", "TaoLineSearchSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 5052a0dac07SAlp Dener if (flg) { 5069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls), monfilename, &monviewer)); 5072a0dac07SAlp Dener ls->viewer = monviewer; 5082a0dac07SAlp Dener ls->usemonitor = PETSC_TRUE; 5092a0dac07SAlp Dener } 510dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, setfromoptions, PetscOptionsObject); 511d0609cedSBarry Smith PetscOptionsEnd(); 5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5130c52818fSSatish Balay } 5140c52818fSSatish Balay 515cc4c1da9SBarry Smith /*@ 5160c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5170c52818fSSatish Balay 5180c52818fSSatish Balay Not Collective 5190c52818fSSatish Balay 5200c52818fSSatish Balay Input Parameter: 52120f4b53cSBarry Smith . ls - the `TaoLineSearch` context 5220c52818fSSatish Balay 523fd292e60Sprj- Output Parameter: 5240c52818fSSatish Balay . type - the line search algorithm in effect 5250c52818fSSatish Balay 5260c52818fSSatish Balay Level: developer 5270c52818fSSatish Balay 52820f4b53cSBarry Smith .seealso: `TaoLineSearch` 5290c52818fSSatish Balay @*/ 530d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 531d71ae5a4SJacob Faibussowitsch { 5320c52818fSSatish Balay PetscFunctionBegin; 5330c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 5344f572ea9SToby Isaac PetscAssertPointer(type, 2); 5350c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 5363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5370c52818fSSatish Balay } 5380c52818fSSatish Balay 5390c52818fSSatish Balay /*@ 5400c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 5410c52818fSSatish Balay routines used by the line search in last application (not cumulative). 5420c52818fSSatish Balay 5430c52818fSSatish Balay Not Collective 5440c52818fSSatish Balay 5450c52818fSSatish Balay Input Parameter: 54620f4b53cSBarry Smith . ls - the `TaoLineSearch` context 5470c52818fSSatish Balay 5480c52818fSSatish Balay Output Parameters: 5490c52818fSSatish Balay + nfeval - number of function evaluations 5500c52818fSSatish Balay . ngeval - number of gradient evaluations 5510c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 5520c52818fSSatish Balay 5530c52818fSSatish Balay Level: intermediate 5540c52818fSSatish Balay 5550c52818fSSatish Balay Note: 55620f4b53cSBarry Smith If the line search is using the `Tao` objective and gradient 55720f4b53cSBarry Smith routines directly (see `TaoLineSearchUseTaoRoutines()`), then the `Tao` 5580c52818fSSatish Balay is already counting the number of evaluations. 5590c52818fSSatish Balay 56020f4b53cSBarry Smith .seealso: `TaoLineSearch` 5610c52818fSSatish Balay @*/ 562d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 563d71ae5a4SJacob Faibussowitsch { 5640c52818fSSatish Balay PetscFunctionBegin; 5650c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 5660c52818fSSatish Balay *nfeval = ls->nfeval; 5670c52818fSSatish Balay *ngeval = ls->ngeval; 5680c52818fSSatish Balay *nfgeval = ls->nfgeval; 5693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5700c52818fSSatish Balay } 5710c52818fSSatish Balay 5720c52818fSSatish Balay /*@ 573441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 57420f4b53cSBarry Smith the standard `Tao` evaluation routines. 5750c52818fSSatish Balay 5760c52818fSSatish Balay Not Collective 5770c52818fSSatish Balay 5780c52818fSSatish Balay Input Parameter: 57920f4b53cSBarry Smith . ls - the `TaoLineSearch` context 5800c52818fSSatish Balay 5810c52818fSSatish Balay Output Parameter: 58220f4b53cSBarry Smith . flg - `PETSC_TRUE` if the line search is using `Tao` evaluation routines, 58320f4b53cSBarry Smith otherwise `PETSC_FALSE` 5840c52818fSSatish Balay 5850c52818fSSatish Balay Level: developer 58620f4b53cSBarry Smith 58720f4b53cSBarry Smith .seealso: `TaoLineSearch` 5880c52818fSSatish Balay @*/ 589d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 590d71ae5a4SJacob Faibussowitsch { 5910c52818fSSatish Balay PetscFunctionBegin; 5920c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 593f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 5943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5950c52818fSSatish Balay } 5960c52818fSSatish Balay 5970c52818fSSatish Balay /*@C 5980c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 5990c52818fSSatish Balay 600c3339decSBarry Smith Logically Collective 6010c52818fSSatish Balay 602d8d19677SJose E. Roman Input Parameters: 60320f4b53cSBarry Smith + ls - the `TaoLineSearch` context 6040c52818fSSatish Balay . func - the objective function evaluation routine 6050c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6060c52818fSSatish Balay 60720f4b53cSBarry Smith Calling sequence of `func`: 60820f4b53cSBarry Smith + ls - the line search context 60920f4b53cSBarry Smith . x - input vector 6100c52818fSSatish Balay . f - function value 6113b242c63SJacob Faibussowitsch - ctx - (optional) user-defined context 6120c52818fSSatish Balay 61320f4b53cSBarry Smith Level: advanced 6140c52818fSSatish Balay 61520f4b53cSBarry Smith Notes: 6160c52818fSSatish Balay Use this routine only if you want the line search objective 61720f4b53cSBarry Smith evaluation routine to be different from the `Tao`'s objective 6180c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6190c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6200c52818fSSatish Balay 6210c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6220c52818fSSatish Balay line search, application programmers should be wary of overriding the 6230c52818fSSatish Balay default objective routine. 6240c52818fSSatish Balay 6251cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 6260c52818fSSatish Balay @*/ 6273b242c63SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *f, void *ctx), void *ctx) 628d71ae5a4SJacob Faibussowitsch { 6290c52818fSSatish Balay PetscFunctionBegin; 6300c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 6310c52818fSSatish Balay 6320c52818fSSatish Balay ls->ops->computeobjective = func; 6330c52818fSSatish Balay if (ctx) ls->userctx_func = ctx; 6340c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 6353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6360c52818fSSatish Balay } 6370c52818fSSatish Balay 6380c52818fSSatish Balay /*@C 6390c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 6400c52818fSSatish Balay 641c3339decSBarry Smith Logically Collective 6420c52818fSSatish Balay 643d8d19677SJose E. Roman Input Parameters: 64420f4b53cSBarry Smith + ls - the `TaoLineSearch` context 6450c52818fSSatish Balay . func - the gradient evaluation routine 6460c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6470c52818fSSatish Balay 64820f4b53cSBarry Smith Calling sequence of `func`: 64920f4b53cSBarry Smith + ls - the linesearch object 65020f4b53cSBarry Smith . x - input vector 6510c52818fSSatish Balay . g - gradient vector 6523b242c63SJacob Faibussowitsch - ctx - (optional) user-defined context 6530c52818fSSatish Balay 6540c52818fSSatish Balay Level: beginner 6550c52818fSSatish Balay 6560c52818fSSatish Balay Note: 6570c52818fSSatish Balay Use this routine only if you want the line search gradient 65820f4b53cSBarry Smith evaluation routine to be different from the `Tao`'s gradient 6590c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6600c52818fSSatish Balay the line search function and/or function/gradient routine. 6610c52818fSSatish Balay 6620c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 6630c52818fSSatish Balay line search, application programmers should be wary of overriding the 6640c52818fSSatish Balay default gradient routine. 6650c52818fSSatish Balay 6661cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 6670c52818fSSatish Balay @*/ 6683b242c63SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec g, void *ctx), void *ctx) 669d71ae5a4SJacob Faibussowitsch { 6700c52818fSSatish Balay PetscFunctionBegin; 6710c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 6720c52818fSSatish Balay ls->ops->computegradient = func; 6730c52818fSSatish Balay if (ctx) ls->userctx_grad = ctx; 6740c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 6753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6760c52818fSSatish Balay } 6770c52818fSSatish Balay 6780c52818fSSatish Balay /*@C 6790c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 6800c52818fSSatish Balay 681c3339decSBarry Smith Logically Collective 6820c52818fSSatish Balay 683d8d19677SJose E. Roman Input Parameters: 68420f4b53cSBarry Smith + ls - the `TaoLineSearch` context 6850c52818fSSatish Balay . func - the objective and gradient evaluation routine 6860c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6870c52818fSSatish Balay 68820f4b53cSBarry Smith Calling sequence of `func`: 68920f4b53cSBarry Smith + ls - the linesearch object 69020f4b53cSBarry Smith . x - input vector 6910c52818fSSatish Balay . f - function value 6920c52818fSSatish Balay . g - gradient vector 6933b242c63SJacob Faibussowitsch - ctx - (optional) user-defined context 6940c52818fSSatish Balay 6950c52818fSSatish Balay Level: beginner 6960c52818fSSatish Balay 6970c52818fSSatish Balay Note: 6980c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 69920f4b53cSBarry Smith evaluation routines to be different from the `Tao`'s objective 7000c52818fSSatish Balay and gradient evaluation routines. 7010c52818fSSatish Balay 7020c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7030c52818fSSatish Balay line search, application programmers should be wary of overriding the 7040c52818fSSatish Balay default objective routine. 7050c52818fSSatish Balay 7061cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 7070c52818fSSatish Balay @*/ 7083b242c63SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, void *ctx), void *ctx) 709d71ae5a4SJacob Faibussowitsch { 7100c52818fSSatish Balay PetscFunctionBegin; 7110c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 7120c52818fSSatish Balay ls->ops->computeobjectiveandgradient = func; 7130c52818fSSatish Balay if (ctx) ls->userctx_funcgrad = ctx; 7140c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7160c52818fSSatish Balay } 7170c52818fSSatish Balay 7180c52818fSSatish Balay /*@C 7190c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 7200c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 7210c52818fSSatish Balay 722c3339decSBarry Smith Logically Collective 7230c52818fSSatish Balay 724d8d19677SJose E. Roman Input Parameters: 72520f4b53cSBarry Smith + ls - the `TaoLineSearch` context 7260c52818fSSatish Balay . func - the objective and gradient evaluation routine 7270c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7280c52818fSSatish Balay 72920f4b53cSBarry Smith Calling sequence of `func`: 73020f4b53cSBarry Smith + ls - the linesearch context 73120f4b53cSBarry Smith . x - input vector 7320c52818fSSatish Balay . s - step direction 7330c52818fSSatish Balay . f - function value 7340c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 7353b242c63SJacob Faibussowitsch - ctx - (optional) user-defined context 7360c52818fSSatish Balay 73720f4b53cSBarry Smith Level: advanced 73820f4b53cSBarry Smith 73920f4b53cSBarry Smith Notes: 7403b242c63SJacob Faibussowitsch Sometimes it is more efficient to compute the inner product of the gradient and the step 7413b242c63SJacob Faibussowitsch direction than it is to compute the gradient, and this is all the line search typically needs 7423b242c63SJacob Faibussowitsch of the gradient. 7433b242c63SJacob Faibussowitsch 74420f4b53cSBarry Smith The gradient will still need to be computed at the end of the line 7450c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 7460c52818fSSatish Balay routine 7470c52818fSSatish Balay 74820f4b53cSBarry Smith Bounded line searches (those used in bounded optimization algorithms) 7490c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 75020f4b53cSBarry Smith x0 and steplength with `TaoLineSearchGetStartingVector()` and `TaoLineSearchGetStepLength()` 7510c52818fSSatish Balay 7520c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7530c52818fSSatish Balay line search, application programmers should be wary of overriding the 7540c52818fSSatish Balay default objective routine. 7550c52818fSSatish Balay 7561cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()` 7570c52818fSSatish Balay @*/ 7583b242c63SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *f, PetscReal *gts, void *ctx), void *ctx) 759d71ae5a4SJacob Faibussowitsch { 7600c52818fSSatish Balay PetscFunctionBegin; 7610c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 7620c52818fSSatish Balay ls->ops->computeobjectiveandgts = func; 7630c52818fSSatish Balay if (ctx) ls->userctx_funcgts = ctx; 7640c52818fSSatish Balay ls->usegts = PETSC_TRUE; 7650c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7670c52818fSSatish Balay } 7680c52818fSSatish Balay 769cc4c1da9SBarry Smith /*@ 77020f4b53cSBarry Smith TaoLineSearchUseTaoRoutines - Informs the `TaoLineSearch` to use the 77120f4b53cSBarry Smith objective and gradient evaluation routines from the given `Tao` object. The default. 7720c52818fSSatish Balay 773c3339decSBarry Smith Logically Collective 7740c52818fSSatish Balay 775d8d19677SJose E. Roman Input Parameters: 77620f4b53cSBarry Smith + ls - the `TaoLineSearch` context 77720f4b53cSBarry Smith - ts - the `Tao` context with defined objective/gradient evaluation routines 7780c52818fSSatish Balay 7790c52818fSSatish Balay Level: developer 7800c52818fSSatish Balay 7811cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()` 7820c52818fSSatish Balay @*/ 783d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts) 784d71ae5a4SJacob Faibussowitsch { 7850c52818fSSatish Balay PetscFunctionBegin; 7860c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 787064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(ts, TAO_CLASSID, 2); 788441846f8SBarry Smith ls->tao = ts; 7890c52818fSSatish Balay ls->usetaoroutines = PETSC_TRUE; 7903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7910c52818fSSatish Balay } 7920c52818fSSatish Balay 7930c52818fSSatish Balay /*@ 7940c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 7950c52818fSSatish Balay 796c3339decSBarry Smith Collective 7970c52818fSSatish Balay 7980c52818fSSatish Balay Input Parameters: 79920f4b53cSBarry Smith + ls - the `TaoLineSearch` context 8000c52818fSSatish Balay - x - input vector 8010c52818fSSatish Balay 8020c52818fSSatish Balay Output Parameter: 80320f4b53cSBarry Smith . f - Objective value at `x` 8040c52818fSSatish Balay 8050c52818fSSatish Balay Level: developer 8060c52818fSSatish Balay 80720f4b53cSBarry Smith Note: 80820f4b53cSBarry Smith `TaoLineSearchComputeObjective()` is typically used within line searches 80920f4b53cSBarry Smith so most users would not generally call this routine themselves. 81020f4b53cSBarry Smith 8111cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()` 8120c52818fSSatish Balay @*/ 813d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 814d71ae5a4SJacob Faibussowitsch { 8150c52818fSSatish Balay Vec gdummy; 8160c52818fSSatish Balay PetscReal gts; 817f06e3bfaSBarry Smith 8180c52818fSSatish Balay PetscFunctionBegin; 8190c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 8200c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 8214f572ea9SToby Isaac PetscAssertPointer(f, 3); 8220c52818fSSatish Balay PetscCheckSameComm(ls, 1, x, 2); 8230c52818fSSatish Balay if (ls->usetaoroutines) { 8249566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(ls->tao, x, f)); 8250c52818fSSatish Balay } else { 8263c859ba3SBarry 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"); 8279566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 828792fecdfSBarry Smith if (ls->ops->computeobjective) PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func)); 829ef1023bdSBarry Smith else if (ls->ops->computeobjectiveandgradient) { 8309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &gdummy)); 831792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgradient)(ls, x, f, gdummy, ls->userctx_funcgrad)); 8329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gdummy)); 833792fecdfSBarry Smith } else PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, >s, ls->userctx_funcgts)); 8349566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 8350c52818fSSatish Balay } 8360c52818fSSatish Balay ls->nfeval++; 8373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8380c52818fSSatish Balay } 8390c52818fSSatish Balay 8400c52818fSSatish Balay /*@ 8410c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 8420c52818fSSatish Balay 843c3339decSBarry Smith Collective 8440c52818fSSatish Balay 8450c52818fSSatish Balay Input Parameters: 84620f4b53cSBarry Smith + ls - the `TaoLineSearch` context 8470c52818fSSatish Balay - x - input vector 8480c52818fSSatish Balay 849d8d19677SJose E. Roman Output Parameters: 85020f4b53cSBarry Smith + f - Objective value at `x` 85120f4b53cSBarry Smith - g - Gradient vector at `x` 8520c52818fSSatish Balay 8530c52818fSSatish Balay Level: developer 8540c52818fSSatish Balay 85520f4b53cSBarry Smith Note: 85620f4b53cSBarry Smith `TaoLineSearchComputeObjectiveAndGradient()` is typically used within line searches 85720f4b53cSBarry Smith so most users would not generally call this routine themselves. 85820f4b53cSBarry Smith 85942747ad1SJacob Faibussowitsch .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchSetObjectiveRoutine()` 8600c52818fSSatish Balay @*/ 861d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 862d71ae5a4SJacob Faibussowitsch { 8630c52818fSSatish Balay PetscFunctionBegin; 8640c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 8650c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 8664f572ea9SToby Isaac PetscAssertPointer(f, 3); 8670c52818fSSatish Balay PetscValidHeaderSpecific(g, VEC_CLASSID, 4); 8680c52818fSSatish Balay PetscCheckSameComm(ls, 1, x, 2); 8690c52818fSSatish Balay PetscCheckSameComm(ls, 1, g, 4); 8700c52818fSSatish Balay if (ls->usetaoroutines) { 8719566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(ls->tao, x, f, g)); 8720c52818fSSatish Balay } else { 8739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 874792fecdfSBarry Smith if (ls->ops->computeobjectiveandgradient) PetscCallBack("TaoLineSearch callback objective/gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, f, g, ls->userctx_funcgrad)); 875ef1023bdSBarry Smith else { 876792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func)); 877792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad)); 878362febeeSStefano Zampini } 8799566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 8809566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f))); 881f06e3bfaSBarry Smith } 882fbe0838dSJason Sarich ls->nfgeval++; 8833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8840c52818fSSatish Balay } 8850c52818fSSatish Balay 8860c52818fSSatish Balay /*@ 8870c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 8880c52818fSSatish Balay 889c3339decSBarry Smith Collective 8900c52818fSSatish Balay 8910c52818fSSatish Balay Input Parameters: 89220f4b53cSBarry Smith + ls - the `TaoLineSearch` context 8930c52818fSSatish Balay - x - input vector 8940c52818fSSatish Balay 8950c52818fSSatish Balay Output Parameter: 8960c52818fSSatish Balay . g - gradient vector 8970c52818fSSatish Balay 8980c52818fSSatish Balay Level: developer 8990c52818fSSatish Balay 90020f4b53cSBarry Smith Note: 90120f4b53cSBarry Smith `TaoComputeGradient()` is typically used within line searches 90220f4b53cSBarry Smith so most users would not generally call this routine themselves. 90320f4b53cSBarry Smith 9041cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()` 9050c52818fSSatish Balay @*/ 906d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 907d71ae5a4SJacob Faibussowitsch { 9080c52818fSSatish Balay PetscReal fdummy; 909f06e3bfaSBarry Smith 9100c52818fSSatish Balay PetscFunctionBegin; 9110c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 9120c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 9130c52818fSSatish Balay PetscValidHeaderSpecific(g, VEC_CLASSID, 3); 9140c52818fSSatish Balay PetscCheckSameComm(ls, 1, x, 2); 9150c52818fSSatish Balay PetscCheckSameComm(ls, 1, g, 3); 9160c52818fSSatish Balay if (ls->usetaoroutines) { 9179566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(ls->tao, x, g)); 9180c52818fSSatish Balay } else { 9199566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 920792fecdfSBarry Smith if (ls->ops->computegradient) PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad)); 92148a46eb9SPierre Jolivet else PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, &fdummy, g, ls->userctx_funcgrad)); 9229566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 9230c52818fSSatish Balay } 9240c52818fSSatish Balay ls->ngeval++; 9253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9260c52818fSSatish Balay } 9270c52818fSSatish Balay 9280c52818fSSatish Balay /*@ 92920f4b53cSBarry Smith TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and 93020f4b53cSBarry Smith step direction at a given point 9310c52818fSSatish Balay 932c3339decSBarry Smith Collective 9330c52818fSSatish Balay 9340c52818fSSatish Balay Input Parameters: 93520f4b53cSBarry Smith + ls - the `TaoLineSearch` context 9360c52818fSSatish Balay - x - input vector 9370c52818fSSatish Balay 938d8d19677SJose E. Roman Output Parameters: 93920f4b53cSBarry Smith + f - Objective value at `x` 94020f4b53cSBarry Smith - gts - inner product of gradient and step direction at `x` 9410c52818fSSatish Balay 9420c52818fSSatish Balay Level: developer 9430c52818fSSatish Balay 94420f4b53cSBarry Smith Note: 94520f4b53cSBarry Smith `TaoLineSearchComputeObjectiveAndGTS()` is typically used within line searches 94620f4b53cSBarry Smith so most users would not generally call this routine themselves. 94720f4b53cSBarry Smith 9481cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()` 9490c52818fSSatish Balay @*/ 950d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 951d71ae5a4SJacob Faibussowitsch { 9520c52818fSSatish Balay PetscFunctionBegin; 9530c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 9540c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 9554f572ea9SToby Isaac PetscAssertPointer(f, 3); 9564f572ea9SToby Isaac PetscAssertPointer(gts, 4); 9570c52818fSSatish Balay PetscCheckSameComm(ls, 1, x, 2); 9589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 959792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback objective/gts", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, gts, ls->userctx_funcgts)); 9609566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0)); 9619566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f))); 9620c52818fSSatish Balay ls->nfeval++; 9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9640c52818fSSatish Balay } 9650c52818fSSatish Balay 9660c52818fSSatish Balay /*@ 9670c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 9680c52818fSSatish Balay 969c3339decSBarry Smith Collective 9700c52818fSSatish Balay 9710c52818fSSatish Balay Input Parameter: 97220f4b53cSBarry Smith . ls - the `TaoLineSearch` context 9730c52818fSSatish Balay 974d8d19677SJose E. Roman Output Parameters: 9750c52818fSSatish Balay + x - the new solution 97620f4b53cSBarry Smith . f - the objective function value at `x` 97720f4b53cSBarry Smith . g - the gradient at `x` 9780c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 9790c52818fSSatish Balay - reason - the reason why the line search terminated 9800c52818fSSatish Balay 9810c52818fSSatish Balay Level: developer 9820c52818fSSatish Balay 98320f4b53cSBarry Smith .seealso: `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()` 9840c52818fSSatish Balay @*/ 985d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 986d71ae5a4SJacob Faibussowitsch { 9870c52818fSSatish Balay PetscFunctionBegin; 9880c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 9890c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 9904f572ea9SToby Isaac PetscAssertPointer(f, 3); 9910c52818fSSatish Balay PetscValidHeaderSpecific(g, VEC_CLASSID, 4); 9924f572ea9SToby Isaac PetscAssertPointer(reason, 6); 9931baa6e33SBarry Smith if (ls->new_x) PetscCall(VecCopy(ls->new_x, x)); 9940c52818fSSatish Balay *f = ls->new_f; 9951baa6e33SBarry Smith if (ls->new_g) PetscCall(VecCopy(ls->new_g, g)); 99697ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 9970c52818fSSatish Balay *reason = ls->reason; 9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9990c52818fSSatish Balay } 10000c52818fSSatish Balay 10010c52818fSSatish Balay /*@ 10020c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 10030c52818fSSatish Balay search. 10040c52818fSSatish Balay 10050c52818fSSatish Balay Not Collective 10060c52818fSSatish Balay 10070c52818fSSatish Balay Input Parameter: 100820f4b53cSBarry Smith . ls - the `TaoLineSearch` context 10090c52818fSSatish Balay 10100c52818fSSatish Balay Output Parameter: 10110c52818fSSatish Balay . x - The initial point of the line search 10120c52818fSSatish Balay 101320f4b53cSBarry Smith Level: advanced 101420f4b53cSBarry Smith 101520f4b53cSBarry Smith .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStepDirection()` 10160c52818fSSatish Balay @*/ 1017d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 1018d71ae5a4SJacob Faibussowitsch { 10190c52818fSSatish Balay PetscFunctionBegin; 10200c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 102197ab8e72SStefano Zampini if (x) *x = ls->start_x; 10223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10230c52818fSSatish Balay } 10240c52818fSSatish Balay 10250c52818fSSatish Balay /*@ 10260c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 10270c52818fSSatish Balay search. 10280c52818fSSatish Balay 10290c52818fSSatish Balay Not Collective 10300c52818fSSatish Balay 10310c52818fSSatish Balay Input Parameter: 103220f4b53cSBarry Smith . ls - the `TaoLineSearch` context 10330c52818fSSatish Balay 10340c52818fSSatish Balay Output Parameter: 10350c52818fSSatish Balay . s - the step direction of the line search 10360c52818fSSatish Balay 10370c52818fSSatish Balay Level: advanced 103820f4b53cSBarry Smith 103920f4b53cSBarry Smith .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()` 10400c52818fSSatish Balay @*/ 1041d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 1042d71ae5a4SJacob Faibussowitsch { 10430c52818fSSatish Balay PetscFunctionBegin; 10440c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 104597ab8e72SStefano Zampini if (s) *s = ls->stepdirection; 10463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10470c52818fSSatish Balay } 10480c52818fSSatish Balay 10490c52818fSSatish Balay /*@ 10500c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 10510c52818fSSatish Balay 10520c52818fSSatish Balay Not Collective 10530c52818fSSatish Balay 10540c52818fSSatish Balay Input Parameter: 105520f4b53cSBarry Smith . ls - the `TaoLineSearch` context 10560c52818fSSatish Balay 10570c52818fSSatish Balay Output Parameter: 1058e056e8ceSJacob Faibussowitsch . f_fullstep - the objective value at the full step length 10590c52818fSSatish Balay 10600c52818fSSatish Balay Level: developer 106120f4b53cSBarry Smith 106220f4b53cSBarry Smith .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()` 10630c52818fSSatish Balay @*/ 1064d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 1065d71ae5a4SJacob Faibussowitsch { 10660c52818fSSatish Balay PetscFunctionBegin; 10670c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 10680c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 10693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10700c52818fSSatish Balay } 10710c52818fSSatish Balay 10720c52818fSSatish Balay /*@ 107320f4b53cSBarry Smith TaoLineSearchSetVariableBounds - Sets the upper and lower bounds for a bounded line search 10740c52818fSSatish Balay 1075c3339decSBarry Smith Logically Collective 10760c52818fSSatish Balay 10770c52818fSSatish Balay Input Parameters: 107820f4b53cSBarry Smith + ls - the `TaoLineSearch` context 10790c52818fSSatish Balay . xl - vector of lower bounds 10800c52818fSSatish Balay - xu - vector of upper bounds 10810c52818fSSatish Balay 10820c52818fSSatish Balay Level: beginner 10830c52818fSSatish Balay 108420f4b53cSBarry Smith Note: 108520f4b53cSBarry Smith If the variable bounds are not set with this routine, then 108620f4b53cSBarry Smith `PETSC_NINFINITY` and `PETSC_INFINITY` are assumed 108720f4b53cSBarry Smith 10881cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoSetVariableBounds()`, `TaoLineSearchCreate()` 10890c52818fSSatish Balay @*/ 1090d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls, Vec xl, Vec xu) 1091d71ae5a4SJacob Faibussowitsch { 10920c52818fSSatish Balay PetscFunctionBegin; 10930c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 109476be6f4fSStefano Zampini if (xl) PetscValidHeaderSpecific(xl, VEC_CLASSID, 2); 109576be6f4fSStefano Zampini if (xu) PetscValidHeaderSpecific(xu, VEC_CLASSID, 3); 10969566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xl)); 10979566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xu)); 10989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->lower)); 10999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->upper)); 11000c52818fSSatish Balay ls->lower = xl; 11010c52818fSSatish Balay ls->upper = xu; 110276be6f4fSStefano Zampini ls->bounded = (PetscBool)(xl || xu); 11033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11040c52818fSSatish Balay } 11050c52818fSSatish Balay 11060c52818fSSatish Balay /*@ 11070c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 11080c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 11090c52818fSSatish Balay 1110c3339decSBarry Smith Logically Collective 11110c52818fSSatish Balay 11120c52818fSSatish Balay Input Parameters: 111320f4b53cSBarry Smith + ls - the `TaoLineSearch` context 11140c52818fSSatish Balay - s - the initial step size 11150c52818fSSatish Balay 11160c52818fSSatish Balay Level: intermediate 11170c52818fSSatish Balay 11181cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()` 11190c52818fSSatish Balay @*/ 1120d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls, PetscReal s) 1121d71ae5a4SJacob Faibussowitsch { 11220c52818fSSatish Balay PetscFunctionBegin; 11230c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 1124a39c8e28SStefano Zampini PetscValidLogicalCollectiveReal(ls, s, 2); 11250c52818fSSatish Balay ls->initstep = s; 11263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11270c52818fSSatish Balay } 11280c52818fSSatish Balay 11290c52818fSSatish Balay /*@ 11300c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 11310c52818fSSatish Balay 11320c52818fSSatish Balay Not Collective 11330c52818fSSatish Balay 113420f4b53cSBarry Smith Input Parameter: 113520f4b53cSBarry Smith . ls - the `TaoLineSearch` context 11360c52818fSSatish Balay 113720f4b53cSBarry Smith Output Parameter: 11380c52818fSSatish Balay . s - the current step length 11390c52818fSSatish Balay 114020f4b53cSBarry Smith Level: intermediate 11410c52818fSSatish Balay 11421cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()` 11430c52818fSSatish Balay @*/ 1144d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls, PetscReal *s) 1145d71ae5a4SJacob Faibussowitsch { 11460c52818fSSatish Balay PetscFunctionBegin; 11470c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 11480c52818fSSatish Balay *s = ls->step; 11493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11500c52818fSSatish Balay } 11510c52818fSSatish Balay 11520ebee16dSLisandro Dalcin /*@C 11530c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 11540c52818fSSatish Balay 1155cc4c1da9SBarry Smith Not Collective, No Fortran Support 11560c52818fSSatish Balay 11570c52818fSSatish Balay Input Parameters: 11580c52818fSSatish Balay + sname - name of a new user-defined solver 11590c52818fSSatish Balay - func - routine to Create method context 11600c52818fSSatish Balay 1161e056e8ceSJacob Faibussowitsch Example Usage: 11620c52818fSSatish Balay .vb 11630c52818fSSatish Balay TaoLineSearchRegister("my_linesearch", MyLinesearchCreate); 11640c52818fSSatish Balay .ve 11650c52818fSSatish Balay 11660c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 1167b44f4de4SBarry Smith .vb 1168b44f4de4SBarry Smith TaoLineSearchSetType(ls, "my_linesearch") 1169b44f4de4SBarry Smith .ve 11700c52818fSSatish Balay or at runtime via the option 1171b44f4de4SBarry Smith .vb 1172b44f4de4SBarry Smith -tao_ls_type my_linesearch 1173b44f4de4SBarry Smith .ve 11740c52818fSSatish Balay 11750c52818fSSatish Balay Level: developer 11760c52818fSSatish Balay 117720f4b53cSBarry Smith Note: 117820f4b53cSBarry Smith `TaoLineSearchRegister()` may be called multiple times to add several user-defined solvers. 117920f4b53cSBarry Smith 11801cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch` 11810ebee16dSLisandro Dalcin @*/ 1182d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 1183d71ae5a4SJacob Faibussowitsch { 11840c52818fSSatish Balay PetscFunctionBegin; 11859566063dSJacob Faibussowitsch PetscCall(TaoLineSearchInitializePackage()); 1186835f2295SStefano Zampini PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, func)); 11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11880c52818fSSatish Balay } 11890c52818fSSatish Balay 1190cc4c1da9SBarry Smith /*@ 11910c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 119220f4b53cSBarry Smith for all `TaoLineSearch` options in the database. 11930c52818fSSatish Balay 1194c3339decSBarry Smith Collective 11950c52818fSSatish Balay 11960c52818fSSatish Balay Input Parameters: 119720f4b53cSBarry Smith + ls - the `TaoLineSearch` solver context 1198e056e8ceSJacob Faibussowitsch - p - the prefix string to prepend to all line search requests 11990c52818fSSatish Balay 120020f4b53cSBarry Smith Level: advanced 120120f4b53cSBarry Smith 12020c52818fSSatish Balay Notes: 12030c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12040c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12050c52818fSSatish Balay 120620f4b53cSBarry Smith This is inherited from the `Tao` object so rarely needs to be set 12070c52818fSSatish Balay 12081cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()` 12090c52818fSSatish Balay @*/ 1210d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 1211d71ae5a4SJacob Faibussowitsch { 1212f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls, p); 12130c52818fSSatish Balay } 12140c52818fSSatish Balay 1215cc4c1da9SBarry Smith /*@ 12160c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 121720f4b53cSBarry Smith `TaoLineSearch` options in the database 12180c52818fSSatish Balay 12190c52818fSSatish Balay Not Collective 12200c52818fSSatish Balay 122120f4b53cSBarry Smith Input Parameter: 122220f4b53cSBarry Smith . ls - the `TaoLineSearch` context 12230c52818fSSatish Balay 122420f4b53cSBarry Smith Output Parameter: 1225e056e8ceSJacob Faibussowitsch . p - pointer to the prefix string used is returned 12260c52818fSSatish Balay 12270c52818fSSatish Balay Level: advanced 12280c52818fSSatish Balay 12291cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()` 12300c52818fSSatish Balay @*/ 1231d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 1232d71ae5a4SJacob Faibussowitsch { 1233f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls, p); 12340c52818fSSatish Balay } 12350c52818fSSatish Balay 1236cc4c1da9SBarry Smith /*@ 12370c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 123820f4b53cSBarry Smith `TaoLineSearch` options in the database. 12390c52818fSSatish Balay 1240c3339decSBarry Smith Logically Collective 12410c52818fSSatish Balay 12420c52818fSSatish Balay Input Parameters: 124320f4b53cSBarry Smith + ls - the `TaoLineSearch` context 1244e056e8ceSJacob Faibussowitsch - p - the prefix string to prepend to all `ls` option requests 124520f4b53cSBarry Smith 124620f4b53cSBarry Smith Level: advanced 12470c52818fSSatish Balay 12480c52818fSSatish Balay Notes: 12490c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12500c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12510c52818fSSatish Balay 125220f4b53cSBarry Smith This is inherited from the `Tao` object so rarely needs to be set 125320f4b53cSBarry Smith 12540c52818fSSatish Balay For example, to distinguish between the runtime options for two 12550c52818fSSatish Balay different line searches, one could call 12560c52818fSSatish Balay .vb 12570c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 12580c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 12590c52818fSSatish Balay .ve 12600c52818fSSatish Balay 12610c52818fSSatish Balay This would enable use of different options for each system, such as 12620c52818fSSatish Balay .vb 12630c52818fSSatish Balay -sys1_tao_ls_type mt 12640c52818fSSatish Balay -sys2_tao_ls_type armijo 12650c52818fSSatish Balay .ve 12660c52818fSSatish Balay 12671cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()` 12680c52818fSSatish Balay @*/ 1269d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 1270d71ae5a4SJacob Faibussowitsch { 1271f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls, p); 12720c52818fSSatish Balay } 1273