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 1220f4b53cSBarry Smith TaoLineSearchViewFromOptions - View a `TaoLineSearch` object based on values in the options database 13fe2efc57SMark 14c3339decSBarry Smith Collective 15fe2efc57SMark 16fe2efc57SMark Input Parameters: 1720f4b53cSBarry Smith + A - the `Tao` context 18736c3998SJose E. Roman . obj - Optional object 19736c3998SJose E. Roman - name - command line option 20fe2efc57SMark 21fe2efc57SMark Level: intermediate 2220f4b53cSBarry Smith 2320f4b53cSBarry Smith Note: 2420f4b53cSBarry Smith See `PetscObjectViewFromOptions()` for available viewer options 2520f4b53cSBarry Smith 26*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchView()`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()` 27fe2efc57SMark @*/ 28d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A, PetscObject obj, const char name[]) 29d71ae5a4SJacob Faibussowitsch { 30fe2efc57SMark PetscFunctionBegin; 31fe2efc57SMark PetscValidHeaderSpecific(A, TAOLINESEARCH_CLASSID, 1); 329566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34fe2efc57SMark } 35fe2efc57SMark 36fe2efc57SMark /*@C 3720f4b53cSBarry Smith TaoLineSearchView - Prints information about the `TaoLineSearch` 380c52818fSSatish Balay 39c3339decSBarry Smith Collective 400c52818fSSatish Balay 410c52818fSSatish Balay InputParameters: 4220f4b53cSBarry Smith + ls - the `TaoLineSearch` context 430c52818fSSatish Balay - viewer - visualization context 440c52818fSSatish Balay 450c52818fSSatish Balay Options Database Key: 4620f4b53cSBarry Smith . -tao_ls_view - Calls `TaoLineSearchView()` at the end of each line search 4720f4b53cSBarry Smith 4820f4b53cSBarry Smith Level: beginner 490c52818fSSatish Balay 500c52818fSSatish Balay Notes: 510c52818fSSatish Balay The available visualization contexts include 5220f4b53cSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 5320f4b53cSBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 540c52818fSSatish Balay output where only the first processor opens 550c52818fSSatish Balay the file. All other processors send their 560c52818fSSatish Balay data to the first processor to print. 570c52818fSSatish Balay 58*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `PetscViewerASCIIOpen()`, `TaoLineSearchViewFromOptions()` 590c52818fSSatish Balay @*/ 60d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 61d71ae5a4SJacob Faibussowitsch { 620c52818fSSatish Balay PetscBool isascii, isstring; 63dedfbcbeSJed Brown TaoLineSearchType type; 64f06e3bfaSBarry Smith 650c52818fSSatish Balay PetscFunctionBegin; 660c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 6748a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer)); 680c52818fSSatish Balay PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 690c52818fSSatish Balay PetscCheckSameComm(ls, 1, viewer, 2); 700c52818fSSatish Balay 719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 730c52818fSSatish Balay if (isascii) { 749566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer)); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 76dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, view, viewer); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "maximum function evaluations=%" PetscInt_FMT "\n", ls->max_funcs)); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "tolerances: ftol=%g, rtol=%g, gtol=%g\n", (double)ls->ftol, (double)ls->rtol, (double)ls->gtol)); 8163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT "\n", ls->nfeval)); 8263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT "\n", ls->ngeval)); 8363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT "\n", ls->nfgeval)); 840c52818fSSatish Balay 8548a46eb9SPierre Jolivet if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(viewer, "using variable bounds\n")); 869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Termination reason: %d\n", (int)ls->reason)); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 880c52818fSSatish Balay } else if (isstring) { 899566063dSJacob Faibussowitsch PetscCall(TaoLineSearchGetType(ls, &type)); 909566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type)); 910c52818fSSatish Balay } 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 930c52818fSSatish Balay } 940c52818fSSatish Balay 950c52818fSSatish Balay /*@C 9620f4b53cSBarry Smith TaoLineSearchCreate - Creates a `TaoLineSearch` object. Algorithms in `Tao` that use 9720f4b53cSBarry Smith line-searches will automatically create one so this all is rarely needed 980c52818fSSatish Balay 99d083f849SBarry Smith Collective 1000c52818fSSatish Balay 1010c52818fSSatish Balay Input Parameter: 1020c52818fSSatish Balay . comm - MPI communicator 1030c52818fSSatish Balay 1040c52818fSSatish Balay Output Parameter: 10520f4b53cSBarry Smith . newls - the new `TaoLineSearch` context 1060c52818fSSatish Balay 10720f4b53cSBarry Smith Options Database Key: 10820f4b53cSBarry Smith . -tao_ls_type - select which method `Tao` should use 1090c52818fSSatish Balay 11020f4b53cSBarry Smith Level: developer 1110c52818fSSatish Balay 112*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()` 1130c52818fSSatish Balay @*/ 114d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 115d71ae5a4SJacob Faibussowitsch { 1160c52818fSSatish Balay TaoLineSearch ls; 1170c52818fSSatish Balay 1180c52818fSSatish Balay PetscFunctionBegin; 1190c52818fSSatish Balay PetscValidPointer(newls, 2); 1209566063dSJacob Faibussowitsch PetscCall(TaoLineSearchInitializePackage()); 1210c52818fSSatish Balay 1229566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(ls, TAOLINESEARCH_CLASSID, "TaoLineSearch", "Linesearch", "Tao", comm, TaoLineSearchDestroy, TaoLineSearchView)); 1230c52818fSSatish Balay ls->max_funcs = 30; 1240c52818fSSatish Balay ls->ftol = 0.0001; 1250c52818fSSatish Balay ls->gtol = 0.9; 1266f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1276f4723b1SBarry Smith ls->rtol = 1.0e-5; 1286f4723b1SBarry Smith #else 1290c52818fSSatish Balay ls->rtol = 1.0e-10; 1306f4723b1SBarry Smith #endif 1310c52818fSSatish Balay ls->stepmin = 1.0e-20; 1320c52818fSSatish Balay ls->stepmax = 1.0e+20; 1330c52818fSSatish Balay ls->step = 1.0; 134a39c8e28SStefano Zampini ls->initstep = 1.0; 1350c52818fSSatish Balay *newls = ls; 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1370c52818fSSatish Balay } 1380c52818fSSatish Balay 1390c52818fSSatish Balay /*@ 1400c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 14120f4b53cSBarry Smith of a `TaoLineSearch` 1420c52818fSSatish Balay 143c3339decSBarry Smith Collective 1440c52818fSSatish Balay 14520f4b53cSBarry Smith Input Parameter: 14620f4b53cSBarry Smith . ls - the `TaoLineSearch` context 1470c52818fSSatish Balay 1480c52818fSSatish Balay Level: developer 1490c52818fSSatish Balay 15020f4b53cSBarry Smith Note: 15120f4b53cSBarry Smith The user will not need to explicitly call `TaoLineSearchSetUp()`, as it will 15220f4b53cSBarry Smith automatically be called in `TaoLineSearchSolve()`. However, if the user 15320f4b53cSBarry Smith desires to call it explicitly, it should come after `TaoLineSearchCreate()` 15420f4b53cSBarry Smith but before `TaoLineSearchApply()`. 1550c52818fSSatish Balay 156*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()` 15720f4b53cSBarry Smith @*/ 158d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 159d71ae5a4SJacob Faibussowitsch { 1608caf6e8cSBarry Smith const char *default_type = TAOLINESEARCHMT; 1610c52818fSSatish Balay PetscBool flg; 162f06e3bfaSBarry Smith 1630c52818fSSatish Balay PetscFunctionBegin; 1640c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 1653ba16761SJacob Faibussowitsch if (ls->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 16648a46eb9SPierre Jolivet if (!((PetscObject)ls)->type_name) PetscCall(TaoLineSearchSetType(ls, default_type)); 167dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, setup); 1680c52818fSSatish Balay if (ls->usetaoroutines) { 1699566063dSJacob Faibussowitsch PetscCall(TaoIsObjectiveDefined(ls->tao, &flg)); 1700c52818fSSatish Balay ls->hasobjective = flg; 1719566063dSJacob Faibussowitsch PetscCall(TaoIsGradientDefined(ls->tao, &flg)); 1720c52818fSSatish Balay ls->hasgradient = flg; 1739566063dSJacob Faibussowitsch PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao, &flg)); 1740c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 1750c52818fSSatish Balay } else { 1760c52818fSSatish Balay if (ls->ops->computeobjective) { 1770c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 1780c52818fSSatish Balay } else { 1790c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 1800c52818fSSatish Balay } 1810c52818fSSatish Balay if (ls->ops->computegradient) { 1820c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 1830c52818fSSatish Balay } else { 1840c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 1850c52818fSSatish Balay } 1860c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 1870c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 1880c52818fSSatish Balay } else { 1890c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 1900c52818fSSatish Balay } 1910c52818fSSatish Balay } 1920c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1940c52818fSSatish Balay } 1950c52818fSSatish Balay 1960c52818fSSatish Balay /*@ 1970c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 19820f4b53cSBarry Smith from one `TaoLineSearchApply()` to the next. This function resets this 1990c52818fSSatish Balay state information. 2000c52818fSSatish Balay 201c3339decSBarry Smith Collective 2020c52818fSSatish Balay 2030c52818fSSatish Balay Input Parameter: 20420f4b53cSBarry Smith . ls - the `TaoLineSearch` context 2050c52818fSSatish Balay 2060c52818fSSatish Balay Level: developer 2070c52818fSSatish Balay 208*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()` 2090c52818fSSatish Balay @*/ 210d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 211d71ae5a4SJacob Faibussowitsch { 2120c52818fSSatish Balay PetscFunctionBegin; 2130c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 214dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, reset); 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2160c52818fSSatish Balay } 217f06e3bfaSBarry Smith 2180c52818fSSatish Balay /*@ 21920f4b53cSBarry Smith TaoLineSearchDestroy - Destroys the `TaoLineSearch` context that was created with 22020f4b53cSBarry Smith `TaoLineSearchCreate()` 2210c52818fSSatish Balay 222c3339decSBarry Smith Collective 2230c52818fSSatish Balay 2247a7aea1fSJed Brown Input Parameter: 22520f4b53cSBarry Smith . ls - the `TaoLineSearch` context 2260c52818fSSatish Balay 22720f4b53cSBarry Smith Level: developer 2280c52818fSSatish Balay 22920f4b53cSBarry Smith .seealse: `TaoLineSearchCreate()`, `TaoLineSearchSolve()` 2300c52818fSSatish Balay @*/ 231d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 232d71ae5a4SJacob Faibussowitsch { 2330c52818fSSatish Balay PetscFunctionBegin; 2343ba16761SJacob Faibussowitsch if (!*ls) PetscFunctionReturn(PETSC_SUCCESS); 2350c52818fSSatish Balay PetscValidHeaderSpecific(*ls, TAOLINESEARCH_CLASSID, 1); 2369371c9d4SSatish Balay if (--((PetscObject)*ls)->refct > 0) { 2379371c9d4SSatish Balay *ls = NULL; 2383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2399371c9d4SSatish Balay } 2409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->stepdirection)); 2419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->start_x)); 2429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->upper)); 2439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->lower)); 24448a46eb9SPierre Jolivet if ((*ls)->ops->destroy) PetscCall((*(*ls)->ops->destroy)(*ls)); 24548a46eb9SPierre Jolivet if ((*ls)->usemonitor) PetscCall(PetscViewerDestroy(&(*ls)->viewer)); 2469566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(ls)); 2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2480c52818fSSatish Balay } 2490c52818fSSatish Balay 2500c52818fSSatish Balay /*@ 25120f4b53cSBarry Smith TaoLineSearchApply - Performs a line-search in a given step direction. 25220f4b53cSBarry Smith Criteria for acceptable step length depends on the line-search algorithm chosen 2530c52818fSSatish Balay 254c3339decSBarry Smith Collective 2550c52818fSSatish Balay 2560c52818fSSatish Balay Input Parameters: 25720f4b53cSBarry Smith + ls - the `TaoLineSearch` context 2580c52818fSSatish Balay - s - search direction 2590c52818fSSatish Balay 2600c52818fSSatish Balay Output Parameters: 26120f4b53cSBarry Smith + x - On input the current solution, on output `x` contains the new solution determined by the line search 262f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution 26320f4b53cSBarry Smith . g - On input the gradient evaluated at `x`, on output contains the gradient at new solution 264f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x) 26520f4b53cSBarry Smith - reason - `TaoLineSearchConvergedReason` reason why the line-search stopped 26620f4b53cSBarry Smith 26720f4b53cSBarry Smith Level: advanced 2680c52818fSSatish Balay 26997bb3fdcSJose E. Roman Notes: 27020f4b53cSBarry Smith The algorithm developer must set up the `TaoLineSearch` with calls to 27120f4b53cSBarry Smith `TaoLineSearchSetObjectiveRoutine()` and `TaoLineSearchSetGradientRoutine()`, 27220f4b53cSBarry Smith `TaoLineSearchSetObjectiveAndGradientRoutine()`, or `TaoLineSearchUseTaoRoutines()`. 27320f4b53cSBarry Smith The latter is done automatically by default and thus requires no user input. 2740c52818fSSatish Balay 2750c52818fSSatish Balay You may or may not need to follow this with a call to 27620f4b53cSBarry Smith `TaoAddLineSearchCounts()`, depending on whether you want these 2770c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 2780c52818fSSatish Balay 279*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearchConvergedReason`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, 28020f4b53cSBarry Smith `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()` 2810c52818fSSatish Balay @*/ 282d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 283d71ae5a4SJacob Faibussowitsch { 2840c52818fSSatish Balay PetscInt low1, low2, low3, high1, high2, high3; 2850c52818fSSatish Balay 2860c52818fSSatish Balay PetscFunctionBegin; 2870c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 2880c52818fSSatish Balay PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 2893f6ba705SLisandro Dalcin PetscValidRealPointer(f, 3); 2900c52818fSSatish Balay PetscValidHeaderSpecific(g, VEC_CLASSID, 4); 2910c52818fSSatish Balay PetscValidHeaderSpecific(s, VEC_CLASSID, 5); 2920c52818fSSatish Balay PetscValidPointer(reason, 7); 2930c52818fSSatish Balay PetscCheckSameComm(ls, 1, x, 2); 2940c52818fSSatish Balay PetscCheckSameTypeAndComm(x, 2, g, 4); 2950c52818fSSatish Balay PetscCheckSameTypeAndComm(x, 2, s, 5); 2969566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low1, &high1)); 2979566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(g, &low2, &high2)); 2989566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(s, &low3, &high3)); 2993c859ba3SBarry Smith PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible vector local lengths"); 3000c52818fSSatish Balay 30197ab8e72SStefano Zampini *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3029566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s)); 3039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->stepdirection)); 304050fc7a3SBarry Smith ls->stepdirection = s; 3050c52818fSSatish Balay 3069566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetUp(ls)); 3070c52818fSSatish Balay ls->nfeval = 0; 3080c52818fSSatish Balay ls->ngeval = 0; 3090c52818fSSatish Balay ls->nfgeval = 0; 3100c52818fSSatish Balay /* Check parameter values */ 3110c52818fSSatish Balay if (ls->ftol < 0.0) { 3129566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: ftol (%g) < 0\n", (double)ls->ftol)); 3130c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3140c52818fSSatish Balay } 3150c52818fSSatish Balay if (ls->rtol < 0.0) { 3169566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: rtol (%g) < 0\n", (double)ls->rtol)); 3170c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3180c52818fSSatish Balay } 3190c52818fSSatish Balay if (ls->gtol < 0.0) { 3209566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: gtol (%g) < 0\n", (double)ls->gtol)); 3210c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3220c52818fSSatish Balay } 3230c52818fSSatish Balay if (ls->stepmin < 0.0) { 3249566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) < 0\n", (double)ls->stepmin)); 3250c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3260c52818fSSatish Balay } 3270c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 3289566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n", (double)ls->stepmin, (double)ls->stepmax)); 3290c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3300c52818fSSatish Balay } 3310c52818fSSatish Balay if (ls->max_funcs < 0) { 3329566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n", ls->max_funcs)); 3330c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3340c52818fSSatish Balay } 3350c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 3369566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Initial Line Search Function Value is Inf or Nan (%g)\n", (double)*f)); 3370c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_INFORNAN; 3380c52818fSSatish Balay } 3390c52818fSSatish Balay 3409566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)x)); 3419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->start_x)); 3420c52818fSSatish Balay ls->start_x = x; 343050fc7a3SBarry Smith 3449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply, ls, 0, 0, 0)); 345dbbe0bcdSBarry Smith PetscUseTypeMethod(ls, apply, x, f, g, s); 3469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0, 0, 0)); 3470c52818fSSatish Balay *reason = ls->reason; 3480c52818fSSatish Balay ls->new_f = *f; 3490c52818fSSatish Balay 35097ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 3510c52818fSSatish Balay 3529566063dSJacob Faibussowitsch PetscCall(TaoLineSearchViewFromOptions(ls, NULL, "-tao_ls_view")); 3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3540c52818fSSatish Balay } 3550c52818fSSatish Balay 3560c52818fSSatish Balay /*@C 3570c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 3580c52818fSSatish Balay 359c3339decSBarry Smith Collective 3600c52818fSSatish Balay 3610c52818fSSatish Balay Input Parameters: 36220f4b53cSBarry Smith + ls - the `TaoLineSearch` context 36320f4b53cSBarry Smith - type - the `TaoLineSearchType` selection 3640c52818fSSatish Balay 36520f4b53cSBarry Smith Options Database Key: 36620f4b53cSBarry Smith . -tao_ls_type <type> - select which method Tao should use at runtime 3670c52818fSSatish Balay 3680c52818fSSatish Balay Level: beginner 3690c52818fSSatish Balay 370*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchCreate()`, `TaoLineSearchGetType()`, 37120f4b53cSBarry Smith `TaoLineSearchApply()` 3720c52818fSSatish Balay @*/ 373d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 374d71ae5a4SJacob Faibussowitsch { 3750c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 3760c52818fSSatish Balay PetscBool flg; 3770c52818fSSatish Balay 3780c52818fSSatish Balay PetscFunctionBegin; 3790c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 3800c52818fSSatish Balay PetscValidCharPointer(type, 2); 3819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg)); 3823ba16761SJacob Faibussowitsch if (flg) PetscFunctionReturn(PETSC_SUCCESS); 3830c52818fSSatish Balay 3849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TaoLineSearchList, type, (void (**)(void)) & r)); 3853c859ba3SBarry Smith PetscCheck(r, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested TaoLineSearch type %s", type); 386dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, destroy); 3870c52818fSSatish Balay ls->max_funcs = 30; 3880c52818fSSatish Balay ls->ftol = 0.0001; 3890c52818fSSatish Balay ls->gtol = 0.9; 3906f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 3916f4723b1SBarry Smith ls->rtol = 1.0e-5; 3926f4723b1SBarry Smith #else 3930c52818fSSatish Balay ls->rtol = 1.0e-10; 3946f4723b1SBarry Smith #endif 3950c52818fSSatish Balay ls->stepmin = 1.0e-20; 3960c52818fSSatish Balay ls->stepmax = 1.0e+20; 3970c52818fSSatish Balay 3980c52818fSSatish Balay ls->nfeval = 0; 3990c52818fSSatish Balay ls->ngeval = 0; 4000c52818fSSatish Balay ls->nfgeval = 0; 40183c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 40283c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 40383c8fe1dSLisandro Dalcin ls->ops->view = NULL; 40483c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 40583c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 4060c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4079566063dSJacob Faibussowitsch PetscCall((*r)(ls)); 4089566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type)); 4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4100c52818fSSatish Balay } 4110c52818fSSatish Balay 4122a0dac07SAlp Dener /*@C 4132a0dac07SAlp Dener TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the 4142a0dac07SAlp Dener iteration number, step length, and function value before calling the implementation 4152a0dac07SAlp Dener specific monitor. 4162a0dac07SAlp Dener 4172a0dac07SAlp Dener Input Parameters: 41820f4b53cSBarry Smith + ls - the `TaoLineSearch` context 4192a0dac07SAlp Dener . its - the current iterate number (>=0) 4202a0dac07SAlp Dener . f - the current objective function value 4212a0dac07SAlp Dener - step - the step length 4222a0dac07SAlp Dener 4232a0dac07SAlp Dener Options Database Key: 4242a0dac07SAlp Dener . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output 4252a0dac07SAlp Dener 4262a0dac07SAlp Dener Level: developer 4272a0dac07SAlp Dener 42820f4b53cSBarry Smith .seealso: `TaoLineSearch` 4292a0dac07SAlp Dener @*/ 430d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step) 431d71ae5a4SJacob Faibussowitsch { 4322a0dac07SAlp Dener PetscInt tabs; 4332a0dac07SAlp Dener 4342a0dac07SAlp Dener PetscFunctionBegin; 4352a0dac07SAlp Dener PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 4362a0dac07SAlp Dener if (ls->usemonitor) { 4379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs)); 4389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel)); 43963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its)); 4409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f)); 4419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step)); 4422a0dac07SAlp Dener if (ls->ops->monitor && its > 0) { 4439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3)); 444dbbe0bcdSBarry Smith PetscUseTypeMethod(ls, monitor); 4452a0dac07SAlp Dener } 4469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs)); 4472a0dac07SAlp Dener } 4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4492a0dac07SAlp Dener } 4502a0dac07SAlp Dener 4510c52818fSSatish Balay /*@ 45220f4b53cSBarry Smith TaoLineSearchSetFromOptions - Sets various `TaoLineSearch` parameters from user 4530c52818fSSatish Balay options. 4540c52818fSSatish Balay 455c3339decSBarry Smith Collective 4560c52818fSSatish Balay 45701d2d390SJose E. Roman Input Parameter: 45820f4b53cSBarry Smith . ls - the `TaoLineSearch` context 4590c52818fSSatish Balay 4600c52818fSSatish Balay Options Database Keys: 46120f4b53cSBarry Smith + -tao_ls_type <type> - The algorithm that `TaoLineSearch` uses (more-thuente, gpcg, unit) 4620c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 4630c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 4640c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 465a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - initial steplength allowed 4660c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 4670c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 4680c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 4690c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 4700c52818fSSatish Balay 4710c52818fSSatish Balay Level: beginner 47220f4b53cSBarry Smith 47320f4b53cSBarry Smith .seealso: `TaoLineSearch` 4740c52818fSSatish Balay @*/ 475d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 476d71ae5a4SJacob Faibussowitsch { 4778caf6e8cSBarry Smith const char *default_type = TAOLINESEARCHMT; 4782a0dac07SAlp Dener char type[256], monfilename[PETSC_MAX_PATH_LEN]; 4792a0dac07SAlp Dener PetscViewer monviewer; 4800c52818fSSatish Balay PetscBool flg; 481f06e3bfaSBarry Smith 4820c52818fSSatish Balay PetscFunctionBegin; 4830c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 484d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)ls); 485ad540459SPierre Jolivet if (((PetscObject)ls)->type_name) default_type = ((PetscObject)ls)->type_name; 4860c52818fSSatish Balay /* Check for type from options */ 4879566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-tao_ls_type", "Tao Line Search type", "TaoLineSearchSetType", TaoLineSearchList, default_type, type, 256, &flg)); 4880c52818fSSatish Balay if (flg) { 4899566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls, type)); 4900c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 4919566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls, default_type)); 4920c52818fSSatish Balay } 4930c52818fSSatish Balay 4949566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_max_funcs", "max function evals in line search", "", ls->max_funcs, &ls->max_funcs, NULL)); 4959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_ftol", "tol for sufficient decrease", "", ls->ftol, &ls->ftol, NULL)); 4969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_gtol", "tol for curvature condition", "", ls->gtol, &ls->gtol, NULL)); 4979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_rtol", "relative tol for acceptable step", "", ls->rtol, &ls->rtol, NULL)); 4989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_stepmin", "lower bound for step", "", ls->stepmin, &ls->stepmin, NULL)); 4999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_stepmax", "upper bound for step", "", ls->stepmax, &ls->stepmax, NULL)); 500a39c8e28SStefano Zampini PetscCall(PetscOptionsReal("-tao_ls_stepinit", "initial step", "", ls->initstep, &ls->initstep, NULL)); 5019566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-tao_ls_monitor", "enable the basic monitor", "TaoLineSearchSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 5022a0dac07SAlp Dener if (flg) { 5039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls), monfilename, &monviewer)); 5042a0dac07SAlp Dener ls->viewer = monviewer; 5052a0dac07SAlp Dener ls->usemonitor = PETSC_TRUE; 5062a0dac07SAlp Dener } 507dbbe0bcdSBarry Smith PetscTryTypeMethod(ls, setfromoptions, PetscOptionsObject); 508d0609cedSBarry Smith PetscOptionsEnd(); 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5100c52818fSSatish Balay } 5110c52818fSSatish Balay 5120c52818fSSatish Balay /*@C 5130c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5140c52818fSSatish Balay 5150c52818fSSatish Balay Not Collective 5160c52818fSSatish Balay 5170c52818fSSatish Balay Input Parameter: 51820f4b53cSBarry Smith . ls - the `TaoLineSearch` context 5190c52818fSSatish Balay 520fd292e60Sprj- Output Parameter: 5210c52818fSSatish Balay . type - the line search algorithm in effect 5220c52818fSSatish Balay 5230c52818fSSatish Balay Level: developer 5240c52818fSSatish Balay 52520f4b53cSBarry Smith .seealso: `TaoLineSearch` 5260c52818fSSatish Balay @*/ 527d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 528d71ae5a4SJacob Faibussowitsch { 5290c52818fSSatish Balay PetscFunctionBegin; 5300c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 5310c52818fSSatish Balay PetscValidPointer(type, 2); 5320c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 5333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5340c52818fSSatish Balay } 5350c52818fSSatish Balay 5360c52818fSSatish Balay /*@ 5370c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 5380c52818fSSatish Balay routines used by the line search in last application (not cumulative). 5390c52818fSSatish Balay 5400c52818fSSatish Balay Not Collective 5410c52818fSSatish Balay 5420c52818fSSatish Balay Input Parameter: 54320f4b53cSBarry Smith . ls - the `TaoLineSearch` context 5440c52818fSSatish Balay 5450c52818fSSatish Balay Output Parameters: 5460c52818fSSatish Balay + nfeval - number of function evaluations 5470c52818fSSatish Balay . ngeval - number of gradient evaluations 5480c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 5490c52818fSSatish Balay 5500c52818fSSatish Balay Level: intermediate 5510c52818fSSatish Balay 5520c52818fSSatish Balay Note: 55320f4b53cSBarry Smith If the line search is using the `Tao` objective and gradient 55420f4b53cSBarry Smith routines directly (see `TaoLineSearchUseTaoRoutines()`), then the `Tao` 5550c52818fSSatish Balay is already counting the number of evaluations. 5560c52818fSSatish Balay 55720f4b53cSBarry Smith .seealso: `TaoLineSearch` 5580c52818fSSatish Balay @*/ 559d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 560d71ae5a4SJacob Faibussowitsch { 5610c52818fSSatish Balay PetscFunctionBegin; 5620c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 5630c52818fSSatish Balay *nfeval = ls->nfeval; 5640c52818fSSatish Balay *ngeval = ls->ngeval; 5650c52818fSSatish Balay *nfgeval = ls->nfgeval; 5663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5670c52818fSSatish Balay } 5680c52818fSSatish Balay 5690c52818fSSatish Balay /*@ 570441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 57120f4b53cSBarry Smith the standard `Tao` evaluation routines. 5720c52818fSSatish Balay 5730c52818fSSatish Balay Not Collective 5740c52818fSSatish Balay 5750c52818fSSatish Balay Input Parameter: 57620f4b53cSBarry Smith . ls - the `TaoLineSearch` context 5770c52818fSSatish Balay 5780c52818fSSatish Balay Output Parameter: 57920f4b53cSBarry Smith . flg - `PETSC_TRUE` if the line search is using `Tao` evaluation routines, 58020f4b53cSBarry Smith otherwise `PETSC_FALSE` 5810c52818fSSatish Balay 5820c52818fSSatish Balay Level: developer 58320f4b53cSBarry Smith 58420f4b53cSBarry Smith .seealso: `TaoLineSearch` 5850c52818fSSatish Balay @*/ 586d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 587d71ae5a4SJacob Faibussowitsch { 5880c52818fSSatish Balay PetscFunctionBegin; 5890c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 590f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 5913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5920c52818fSSatish Balay } 5930c52818fSSatish Balay 5940c52818fSSatish Balay /*@C 5950c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 5960c52818fSSatish Balay 597c3339decSBarry Smith Logically Collective 5980c52818fSSatish Balay 599d8d19677SJose E. Roman Input Parameters: 60020f4b53cSBarry Smith + ls - the `TaoLineSearch` context 6010c52818fSSatish Balay . func - the objective function evaluation routine 6020c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6030c52818fSSatish Balay 60420f4b53cSBarry Smith Calling sequence of `func`: 60520f4b53cSBarry Smith $ PetscErrorCode func(TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 60620f4b53cSBarry Smith + ls - the line search context 60720f4b53cSBarry Smith . x - input vector 6080c52818fSSatish Balay . f - function value 6090c52818fSSatish Balay - ctx (optional) user-defined context 6100c52818fSSatish Balay 61120f4b53cSBarry Smith Level: advanced 6120c52818fSSatish Balay 61320f4b53cSBarry Smith Notes: 6140c52818fSSatish Balay Use this routine only if you want the line search objective 61520f4b53cSBarry Smith evaluation routine to be different from the `Tao`'s objective 6160c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6170c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6180c52818fSSatish Balay 6190c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6200c52818fSSatish Balay line search, application programmers should be wary of overriding the 6210c52818fSSatish Balay default objective routine. 6220c52818fSSatish Balay 623*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 6240c52818fSSatish Balay @*/ 625d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *, void *), void *ctx) 626d71ae5a4SJacob Faibussowitsch { 6270c52818fSSatish Balay PetscFunctionBegin; 6280c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 6290c52818fSSatish Balay 6300c52818fSSatish Balay ls->ops->computeobjective = func; 6310c52818fSSatish Balay if (ctx) ls->userctx_func = ctx; 6320c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 6333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6340c52818fSSatish Balay } 6350c52818fSSatish Balay 6360c52818fSSatish Balay /*@C 6370c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 6380c52818fSSatish Balay 639c3339decSBarry Smith Logically Collective 6400c52818fSSatish Balay 641d8d19677SJose E. Roman Input Parameters: 64220f4b53cSBarry Smith + ls - the `TaoLineSearch` context 6430c52818fSSatish Balay . func - the gradient evaluation routine 6440c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6450c52818fSSatish Balay 64620f4b53cSBarry Smith Calling sequence of `func`: 64720f4b53cSBarry Smith $ PetscErrorCode func(TaoLinesearch ls, Vec x, Vec g, void *ctx); 64820f4b53cSBarry Smith + ls - the linesearch object 64920f4b53cSBarry Smith . x - input vector 6500c52818fSSatish Balay . g - gradient vector 6510c52818fSSatish Balay - ctx (optional) user-defined context 6520c52818fSSatish Balay 6530c52818fSSatish Balay Level: beginner 6540c52818fSSatish Balay 6550c52818fSSatish Balay Note: 6560c52818fSSatish Balay Use this routine only if you want the line search gradient 65720f4b53cSBarry Smith evaluation routine to be different from the `Tao`'s gradient 6580c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6590c52818fSSatish Balay the line search function and/or function/gradient routine. 6600c52818fSSatish Balay 6610c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 6620c52818fSSatish Balay line search, application programmers should be wary of overriding the 6630c52818fSSatish Balay default gradient routine. 6640c52818fSSatish Balay 665*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 6660c52818fSSatish Balay @*/ 667d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec g, void *), void *ctx) 668d71ae5a4SJacob Faibussowitsch { 6690c52818fSSatish Balay PetscFunctionBegin; 6700c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 6710c52818fSSatish Balay ls->ops->computegradient = func; 6720c52818fSSatish Balay if (ctx) ls->userctx_grad = ctx; 6730c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 6743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6750c52818fSSatish Balay } 6760c52818fSSatish Balay 6770c52818fSSatish Balay /*@C 6780c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 6790c52818fSSatish Balay 680c3339decSBarry Smith Logically Collective 6810c52818fSSatish Balay 682d8d19677SJose E. Roman Input Parameters: 68320f4b53cSBarry Smith + ls - the `TaoLineSearch` context 6840c52818fSSatish Balay . func - the objective and gradient evaluation routine 6850c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6860c52818fSSatish Balay 68720f4b53cSBarry Smith Calling sequence of `func`: 68820f4b53cSBarry Smith $ PetscErrorCode func(TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 68920f4b53cSBarry Smith + ls - the linesearch object 69020f4b53cSBarry Smith . x - input vector 6910c52818fSSatish Balay . f - function value 6920c52818fSSatish Balay . g - gradient vector 6930c52818fSSatish Balay - 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 706*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 7070c52818fSSatish Balay @*/ 708d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void *), 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 Sometimes it is more efficient to compute the inner product of the gradient 7220c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 7230c52818fSSatish Balay the line search typically needs of the gradient. 7240c52818fSSatish Balay 725c3339decSBarry Smith Logically Collective 7260c52818fSSatish Balay 727d8d19677SJose E. Roman Input Parameters: 72820f4b53cSBarry Smith + ls - the `TaoLineSearch` context 7290c52818fSSatish Balay . func - the objective and gradient evaluation routine 7300c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7310c52818fSSatish Balay 73220f4b53cSBarry Smith Calling sequence of `func`: 73320f4b53cSBarry Smith $ PetscErrorCode func(TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 73420f4b53cSBarry Smith + ls - the linesearch context 73520f4b53cSBarry Smith . x - input vector 7360c52818fSSatish Balay . s - step direction 7370c52818fSSatish Balay . f - function value 7380c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 7390c52818fSSatish Balay - ctx (optional) user-defined context 7400c52818fSSatish Balay 74120f4b53cSBarry Smith Level: advanced 74220f4b53cSBarry Smith 74320f4b53cSBarry Smith Notes: 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 756*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()` 7570c52818fSSatish Balay @*/ 758d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void *), 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 7690c52818fSSatish Balay /*@C 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 781*1cc06b55SBarry 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 811*1cc06b55SBarry 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); 821dadcf809SJacob Faibussowitsch PetscValidRealPointer(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 859*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `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); 866dadcf809SJacob Faibussowitsch PetscValidRealPointer(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 904*1cc06b55SBarry 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 948*1cc06b55SBarry 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); 955dadcf809SJacob Faibussowitsch PetscValidRealPointer(f, 3); 956dadcf809SJacob Faibussowitsch PetscValidRealPointer(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); 990dadcf809SJacob Faibussowitsch PetscValidRealPointer(f, 3); 9910c52818fSSatish Balay PetscValidHeaderSpecific(g, VEC_CLASSID, 4); 9920c52818fSSatish Balay PetscValidIntPointer(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: 10580c52818fSSatish Balay . f - the objective value at the full step length 10590c52818fSSatish Balay 10600c52818fSSatish Balay Level: developer 106120f4b53cSBarry Smith 106220f4b53cSBarry Smith .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()` 10630c52818fSSatish Balay @*/ 10640c52818fSSatish Balay 1065d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 1066d71ae5a4SJacob Faibussowitsch { 10670c52818fSSatish Balay PetscFunctionBegin; 10680c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 10690c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 10703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10710c52818fSSatish Balay } 10720c52818fSSatish Balay 10730c52818fSSatish Balay /*@ 107420f4b53cSBarry Smith TaoLineSearchSetVariableBounds - Sets the upper and lower bounds for a bounded line search 10750c52818fSSatish Balay 1076c3339decSBarry Smith Logically Collective 10770c52818fSSatish Balay 10780c52818fSSatish Balay Input Parameters: 107920f4b53cSBarry Smith + ls - the `TaoLineSearch` context 10800c52818fSSatish Balay . xl - vector of lower bounds 10810c52818fSSatish Balay - xu - vector of upper bounds 10820c52818fSSatish Balay 10830c52818fSSatish Balay Level: beginner 10840c52818fSSatish Balay 108520f4b53cSBarry Smith Note: 108620f4b53cSBarry Smith If the variable bounds are not set with this routine, then 108720f4b53cSBarry Smith `PETSC_NINFINITY` and `PETSC_INFINITY` are assumed 108820f4b53cSBarry Smith 1089*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoSetVariableBounds()`, `TaoLineSearchCreate()` 10900c52818fSSatish Balay @*/ 1091d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls, Vec xl, Vec xu) 1092d71ae5a4SJacob Faibussowitsch { 10930c52818fSSatish Balay PetscFunctionBegin; 10940c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 109576be6f4fSStefano Zampini if (xl) PetscValidHeaderSpecific(xl, VEC_CLASSID, 2); 109676be6f4fSStefano Zampini if (xu) PetscValidHeaderSpecific(xu, VEC_CLASSID, 3); 10979566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xl)); 10989566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xu)); 10999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->lower)); 11009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->upper)); 11010c52818fSSatish Balay ls->lower = xl; 11020c52818fSSatish Balay ls->upper = xu; 110376be6f4fSStefano Zampini ls->bounded = (PetscBool)(xl || xu); 11043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11050c52818fSSatish Balay } 11060c52818fSSatish Balay 11070c52818fSSatish Balay /*@ 11080c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 11090c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 11100c52818fSSatish Balay 1111c3339decSBarry Smith Logically Collective 11120c52818fSSatish Balay 11130c52818fSSatish Balay Input Parameters: 111420f4b53cSBarry Smith + ls - the `TaoLineSearch` context 11150c52818fSSatish Balay - s - the initial step size 11160c52818fSSatish Balay 11170c52818fSSatish Balay Level: intermediate 11180c52818fSSatish Balay 1119*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()` 11200c52818fSSatish Balay @*/ 1121d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls, PetscReal s) 1122d71ae5a4SJacob Faibussowitsch { 11230c52818fSSatish Balay PetscFunctionBegin; 11240c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 1125a39c8e28SStefano Zampini PetscValidLogicalCollectiveReal(ls, s, 2); 11260c52818fSSatish Balay ls->initstep = s; 11273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11280c52818fSSatish Balay } 11290c52818fSSatish Balay 11300c52818fSSatish Balay /*@ 11310c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 11320c52818fSSatish Balay 11330c52818fSSatish Balay Not Collective 11340c52818fSSatish Balay 113520f4b53cSBarry Smith Input Parameter: 113620f4b53cSBarry Smith . ls - the `TaoLineSearch` context 11370c52818fSSatish Balay 113820f4b53cSBarry Smith Output Parameter: 11390c52818fSSatish Balay . s - the current step length 11400c52818fSSatish Balay 114120f4b53cSBarry Smith Level: intermediate 11420c52818fSSatish Balay 1143*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()` 11440c52818fSSatish Balay @*/ 1145d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls, PetscReal *s) 1146d71ae5a4SJacob Faibussowitsch { 11470c52818fSSatish Balay PetscFunctionBegin; 11480c52818fSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 11490c52818fSSatish Balay *s = ls->step; 11503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11510c52818fSSatish Balay } 11520c52818fSSatish Balay 11530ebee16dSLisandro Dalcin /*@C 11540c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 11550c52818fSSatish Balay 115620f4b53cSBarry Smith Not Collective 11570c52818fSSatish Balay 11580c52818fSSatish Balay Input Parameters: 11590c52818fSSatish Balay + sname - name of a new user-defined solver 11600c52818fSSatish Balay - func - routine to Create method context 11610c52818fSSatish Balay 11620c52818fSSatish Balay Sample usage: 11630c52818fSSatish Balay .vb 11640c52818fSSatish Balay TaoLineSearchRegister("my_linesearch", MyLinesearchCreate); 11650c52818fSSatish Balay .ve 11660c52818fSSatish Balay 11670c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 11680c52818fSSatish Balay $ TaoLineSearchSetType(ls, "my_linesearch") 11690c52818fSSatish Balay or at runtime via the option 11700c52818fSSatish Balay $ -tao_ls_type my_linesearch 11710c52818fSSatish Balay 11720c52818fSSatish Balay Level: developer 11730c52818fSSatish Balay 117420f4b53cSBarry Smith Note: 117520f4b53cSBarry Smith `TaoLineSearchRegister()` may be called multiple times to add several user-defined solvers. 117620f4b53cSBarry Smith 1177*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch` 11780ebee16dSLisandro Dalcin @*/ 1179d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 1180d71ae5a4SJacob Faibussowitsch { 11810c52818fSSatish Balay PetscFunctionBegin; 11829566063dSJacob Faibussowitsch PetscCall(TaoLineSearchInitializePackage()); 11839566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func)); 11843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11850c52818fSSatish Balay } 11860c52818fSSatish Balay 11870c52818fSSatish Balay /*@C 11880c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 118920f4b53cSBarry Smith for all `TaoLineSearch` options in the database. 11900c52818fSSatish Balay 1191c3339decSBarry Smith Collective 11920c52818fSSatish Balay 11930c52818fSSatish Balay Input Parameters: 119420f4b53cSBarry Smith + ls - the `TaoLineSearch` solver context 11950c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 11960c52818fSSatish Balay 119720f4b53cSBarry Smith Level: advanced 119820f4b53cSBarry Smith 11990c52818fSSatish Balay Notes: 12000c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12010c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12020c52818fSSatish Balay 120320f4b53cSBarry Smith This is inherited from the `Tao` object so rarely needs to be set 12040c52818fSSatish Balay 1205*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()` 12060c52818fSSatish Balay @*/ 1207d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 1208d71ae5a4SJacob Faibussowitsch { 1209f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls, p); 12100c52818fSSatish Balay } 12110c52818fSSatish Balay 12120c52818fSSatish Balay /*@C 12130c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 121420f4b53cSBarry Smith `TaoLineSearch` options in the database 12150c52818fSSatish Balay 12160c52818fSSatish Balay Not Collective 12170c52818fSSatish Balay 121820f4b53cSBarry Smith Input Parameter: 121920f4b53cSBarry Smith . ls - the `TaoLineSearch` context 12200c52818fSSatish Balay 122120f4b53cSBarry Smith Output Parameter: 12220c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 12230c52818fSSatish Balay 12240c52818fSSatish Balay Level: advanced 12250c52818fSSatish Balay 122620f4b53cSBarry Smith Fortran Note: 122720f4b53cSBarry Smith The user should pass in a string 'prefix' of 122820f4b53cSBarry Smith sufficient length to hold the prefix. 122920f4b53cSBarry Smith 1230*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()` 12310c52818fSSatish Balay @*/ 1232d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 1233d71ae5a4SJacob Faibussowitsch { 1234f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls, p); 12350c52818fSSatish Balay } 12360c52818fSSatish Balay 12370c52818fSSatish Balay /*@C 12380c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 123920f4b53cSBarry Smith `TaoLineSearch` options in the database. 12400c52818fSSatish Balay 1241c3339decSBarry Smith Logically Collective 12420c52818fSSatish Balay 12430c52818fSSatish Balay Input Parameters: 124420f4b53cSBarry Smith + ls - the `TaoLineSearch` context 124520f4b53cSBarry Smith - prefix - the prefix string to prepend to all `ls` option requests 124620f4b53cSBarry Smith 124720f4b53cSBarry Smith Level: advanced 12480c52818fSSatish Balay 12490c52818fSSatish Balay Notes: 12500c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12510c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12520c52818fSSatish Balay 125320f4b53cSBarry Smith This is inherited from the `Tao` object so rarely needs to be set 125420f4b53cSBarry Smith 12550c52818fSSatish Balay For example, to distinguish between the runtime options for two 12560c52818fSSatish Balay different line searches, one could call 12570c52818fSSatish Balay .vb 12580c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 12590c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 12600c52818fSSatish Balay .ve 12610c52818fSSatish Balay 12620c52818fSSatish Balay This would enable use of different options for each system, such as 12630c52818fSSatish Balay .vb 12640c52818fSSatish Balay -sys1_tao_ls_type mt 12650c52818fSSatish Balay -sys2_tao_ls_type armijo 12660c52818fSSatish Balay .ve 12670c52818fSSatish Balay 1268*1cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()` 12690c52818fSSatish Balay @*/ 1270d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 1271d71ae5a4SJacob Faibussowitsch { 1272f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls, p); 12730c52818fSSatish Balay } 1274