1ba92ff59SBarry Smith #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h> 30c52818fSSatish Balay 46c23d075SBarry Smith PetscFunctionList TaoLineSearchList = NULL; 50c52818fSSatish Balay 60c52818fSSatish Balay PetscClassId TAOLINESEARCH_CLASSID=0; 70ebee16dSLisandro Dalcin 80ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Apply; 90ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Eval; 100c52818fSSatish Balay 110c52818fSSatish Balay /*@C 12fe2efc57SMark TaoLineSearchViewFromOptions - View from Options 13fe2efc57SMark 14fe2efc57SMark Collective on TaoLineSearch 15fe2efc57SMark 16fe2efc57SMark Input Parameters: 17fe2efc57SMark + A - the Tao context 18736c3998SJose E. Roman . obj - Optional object 19736c3998SJose E. Roman - name - command line option 20fe2efc57SMark 21fe2efc57SMark Level: intermediate 22db781477SPatrick Sanan .seealso: `TaoLineSearch`, `TaoLineSearchView`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()` 23fe2efc57SMark @*/ 24fe2efc57SMark PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A,PetscObject obj,const char name[]) 25fe2efc57SMark { 26fe2efc57SMark PetscFunctionBegin; 27fe2efc57SMark PetscValidHeaderSpecific(A,TAOLINESEARCH_CLASSID,1); 289566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 29fe2efc57SMark PetscFunctionReturn(0); 30fe2efc57SMark } 31fe2efc57SMark 32fe2efc57SMark /*@C 330c52818fSSatish Balay TaoLineSearchView - Prints information about the TaoLineSearch 340c52818fSSatish Balay 350c52818fSSatish Balay Collective on TaoLineSearch 360c52818fSSatish Balay 370c52818fSSatish Balay InputParameters: 38441846f8SBarry Smith + ls - the Tao context 390c52818fSSatish Balay - viewer - visualization context 400c52818fSSatish Balay 410c52818fSSatish Balay Options Database Key: 420c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search 430c52818fSSatish Balay 440c52818fSSatish Balay Notes: 450c52818fSSatish Balay The available visualization contexts include 460c52818fSSatish Balay + PETSC_VIEWER_STDOUT_SELF - standard output (default) 470c52818fSSatish Balay - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 480c52818fSSatish Balay output where only the first processor opens 490c52818fSSatish Balay the file. All other processors send their 500c52818fSSatish Balay data to the first processor to print. 510c52818fSSatish Balay 520c52818fSSatish Balay Level: beginner 530c52818fSSatish Balay 54db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()` 550c52818fSSatish Balay @*/ 560c52818fSSatish Balay 570c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 580c52818fSSatish Balay { 590c52818fSSatish Balay PetscBool isascii, isstring; 60dedfbcbeSJed Brown TaoLineSearchType type; 61f06e3bfaSBarry Smith 620c52818fSSatish Balay PetscFunctionBegin; 630c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 640c52818fSSatish Balay if (!viewer) { 659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer)); 660c52818fSSatish Balay } 670c52818fSSatish Balay PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 680c52818fSSatish Balay PetscCheckSameComm(ls,1,viewer,2); 690c52818fSSatish Balay 709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 720c52818fSSatish Balay if (isascii) { 739566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 75*dbbe0bcdSBarry Smith PetscTryTypeMethod(ls,view ,viewer); 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%" PetscInt_FMT "\n",ls->max_funcs)); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol)); 8063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%" PetscInt_FMT "\n",ls->nfeval)); 8163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%" PetscInt_FMT "\n",ls->ngeval)); 8263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%" PetscInt_FMT "\n",ls->nfgeval)); 830c52818fSSatish Balay 840c52818fSSatish Balay if (ls->bounded) { 859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"using variable bounds\n")); 860c52818fSSatish Balay } 879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason)); 889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 890c52818fSSatish Balay } else if (isstring) { 909566063dSJacob Faibussowitsch PetscCall(TaoLineSearchGetType(ls,&type)); 919566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer," %-3.3s",type)); 920c52818fSSatish Balay } 930c52818fSSatish Balay PetscFunctionReturn(0); 940c52818fSSatish Balay } 950c52818fSSatish Balay 960c52818fSSatish Balay /*@C 970c52818fSSatish Balay TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 980c52818fSSatish Balay line-searches will automatically create one. 990c52818fSSatish Balay 100d083f849SBarry Smith Collective 1010c52818fSSatish Balay 1020c52818fSSatish Balay Input Parameter: 1030c52818fSSatish Balay . comm - MPI communicator 1040c52818fSSatish Balay 1050c52818fSSatish Balay Output Parameter: 1060c52818fSSatish Balay . newls - the new TaoLineSearch context 1070c52818fSSatish Balay 1080c52818fSSatish Balay Available methods include: 109147403d9SBarry Smith + more-thuente - the More-Thuente method 110147403d9SBarry Smith . gpcg - the GPCG method 1110c52818fSSatish Balay - unit - Do not perform any line search 1120c52818fSSatish Balay 1130c52818fSSatish Balay Options Database Keys: 1140c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 1150c52818fSSatish Balay 1160c52818fSSatish Balay Level: beginner 1170c52818fSSatish Balay 118db781477SPatrick Sanan .seealso: `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()` 1190c52818fSSatish Balay @*/ 1200c52818fSSatish Balay 1210c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 1220c52818fSSatish Balay { 1230c52818fSSatish Balay TaoLineSearch ls; 1240c52818fSSatish Balay 1250c52818fSSatish Balay PetscFunctionBegin; 1260c52818fSSatish Balay PetscValidPointer(newls,2); 1279566063dSJacob Faibussowitsch PetscCall(TaoLineSearchInitializePackage()); 1280c52818fSSatish Balay 1299566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView)); 1300c52818fSSatish Balay ls->max_funcs = 30; 1310c52818fSSatish Balay ls->ftol = 0.0001; 1320c52818fSSatish Balay ls->gtol = 0.9; 1336f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1346f4723b1SBarry Smith ls->rtol = 1.0e-5; 1356f4723b1SBarry Smith #else 1360c52818fSSatish Balay ls->rtol = 1.0e-10; 1376f4723b1SBarry Smith #endif 1380c52818fSSatish Balay ls->stepmin = 1.0e-20; 1390c52818fSSatish Balay ls->stepmax = 1.0e+20; 1400c52818fSSatish Balay ls->step = 1.0; 141a39c8e28SStefano Zampini ls->initstep = 1.0; 1420c52818fSSatish Balay *newls = ls; 1430c52818fSSatish Balay PetscFunctionReturn(0); 1440c52818fSSatish Balay } 1450c52818fSSatish Balay 1460c52818fSSatish Balay /*@ 1470c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 1480c52818fSSatish Balay of a Tao solver 1490c52818fSSatish Balay 1500c52818fSSatish Balay Collective on ls 1510c52818fSSatish Balay 1520c52818fSSatish Balay Input Parameters: 1530c52818fSSatish Balay . ls - the TaoLineSearch context 1540c52818fSSatish Balay 1550c52818fSSatish Balay Notes: 1560c52818fSSatish Balay The user will not need to explicitly call TaoLineSearchSetUp(), as it will 1570c52818fSSatish Balay automatically be called in TaoLineSearchSolve(). However, if the user 1580c52818fSSatish Balay desires to call it explicitly, it should come after TaoLineSearchCreate() 1590c52818fSSatish Balay but before TaoLineSearchApply(). 1600c52818fSSatish Balay 1610c52818fSSatish Balay Level: developer 1620c52818fSSatish Balay 163db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()` 1640c52818fSSatish Balay @*/ 1650c52818fSSatish Balay 1660c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 1670c52818fSSatish Balay { 1688caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 1690c52818fSSatish Balay PetscBool flg; 170f06e3bfaSBarry Smith 1710c52818fSSatish Balay PetscFunctionBegin; 1720c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1730c52818fSSatish Balay if (ls->setupcalled) PetscFunctionReturn(0); 1740c52818fSSatish Balay if (!((PetscObject)ls)->type_name) { 1759566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls,default_type)); 1760c52818fSSatish Balay } 177*dbbe0bcdSBarry Smith PetscTryTypeMethod(ls,setup); 1780c52818fSSatish Balay if (ls->usetaoroutines) { 1799566063dSJacob Faibussowitsch PetscCall(TaoIsObjectiveDefined(ls->tao,&flg)); 1800c52818fSSatish Balay ls->hasobjective = flg; 1819566063dSJacob Faibussowitsch PetscCall(TaoIsGradientDefined(ls->tao,&flg)); 1820c52818fSSatish Balay ls->hasgradient = flg; 1839566063dSJacob Faibussowitsch PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao,&flg)); 1840c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 1850c52818fSSatish Balay } else { 1860c52818fSSatish Balay if (ls->ops->computeobjective) { 1870c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 1880c52818fSSatish Balay } else { 1890c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 1900c52818fSSatish Balay } 1910c52818fSSatish Balay if (ls->ops->computegradient) { 1920c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 1930c52818fSSatish Balay } else { 1940c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 1950c52818fSSatish Balay } 1960c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 1970c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 1980c52818fSSatish Balay } else { 1990c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 2000c52818fSSatish Balay } 2010c52818fSSatish Balay } 2020c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 2030c52818fSSatish Balay PetscFunctionReturn(0); 2040c52818fSSatish Balay } 2050c52818fSSatish Balay 2060c52818fSSatish Balay /*@ 2070c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 2080c52818fSSatish Balay from one TaoLineSearchApply() to the next. This function resets this 2090c52818fSSatish Balay state information. 2100c52818fSSatish Balay 2110c52818fSSatish Balay Collective on TaoLineSearch 2120c52818fSSatish Balay 2130c52818fSSatish Balay Input Parameter: 2140c52818fSSatish Balay . ls - the TaoLineSearch context 2150c52818fSSatish Balay 2160c52818fSSatish Balay Level: developer 2170c52818fSSatish Balay 218db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchApply()` 2190c52818fSSatish Balay @*/ 2200c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 2210c52818fSSatish Balay { 2220c52818fSSatish Balay PetscFunctionBegin; 2230c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 224*dbbe0bcdSBarry Smith PetscTryTypeMethod(ls,reset); 2250c52818fSSatish Balay PetscFunctionReturn(0); 2260c52818fSSatish Balay } 227f06e3bfaSBarry Smith 2280c52818fSSatish Balay /*@ 2290c52818fSSatish Balay TaoLineSearchDestroy - Destroys the TAO context that was created with 2300c52818fSSatish Balay TaoLineSearchCreate() 2310c52818fSSatish Balay 2320c52818fSSatish Balay Collective on TaoLineSearch 2330c52818fSSatish Balay 2347a7aea1fSJed Brown Input Parameter: 2350c52818fSSatish Balay . ls - the TaoLineSearch context 2360c52818fSSatish Balay 2370c52818fSSatish Balay Level: beginner 2380c52818fSSatish Balay 2390c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 2400c52818fSSatish Balay @*/ 2410c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 2420c52818fSSatish Balay { 2430c52818fSSatish Balay PetscFunctionBegin; 2440c52818fSSatish Balay if (!*ls) PetscFunctionReturn(0); 2450c52818fSSatish Balay PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 24683c8fe1dSLisandro Dalcin if (--((PetscObject)*ls)->refct > 0) {*ls = NULL; PetscFunctionReturn(0);} 2479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->stepdirection)); 2489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->start_x)); 2499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->upper)); 2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ls)->lower)); 2510c52818fSSatish Balay if ((*ls)->ops->destroy) { 2529566063dSJacob Faibussowitsch PetscCall((*(*ls)->ops->destroy)(*ls)); 2530c52818fSSatish Balay } 2542a0dac07SAlp Dener if ((*ls)->usemonitor) { 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*ls)->viewer)); 2562a0dac07SAlp Dener } 2579566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(ls)); 2580c52818fSSatish Balay PetscFunctionReturn(0); 2590c52818fSSatish Balay } 2600c52818fSSatish Balay 2610c52818fSSatish Balay /*@ 2620c52818fSSatish Balay TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 2630c52818fSSatish Balay 2640c52818fSSatish Balay Collective on TaoLineSearch 2650c52818fSSatish Balay 2660c52818fSSatish Balay Input Parameters: 267441846f8SBarry Smith + ls - the Tao context 2680c52818fSSatish Balay - s - search direction 2690c52818fSSatish Balay 2706b867d5aSJose E. Roman Input/Output Parameters: 2716b867d5aSJose E. Roman 2720c52818fSSatish Balay Output Parameters: 273f1a722f8SMatthew G. Knepley + x - On input the current solution, on output x contains the new solution determined by the line search 274f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution 275f1a722f8SMatthew G. Knepley . g - On input the gradient evaluated at x, on output contains the gradient at new solution 276f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x) 2770c52818fSSatish Balay - reason - reason why the line-search stopped 2780c52818fSSatish Balay 27997bb3fdcSJose E. Roman Notes: 2800c52818fSSatish Balay reason will be set to one of: 2810c52818fSSatish Balay 2820c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 2830c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 2840c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 2850c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 2860c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 2870c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 2880c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 2890c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 2900c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 2910c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search 2920c52818fSSatish Balay 2930c52818fSSatish Balay The algorithm developer must set up the TaoLineSearch with calls to 294441846f8SBarry Smith TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines() 2950c52818fSSatish Balay 2960c52818fSSatish Balay You may or may not need to follow this with a call to 2970c52818fSSatish Balay TaoAddLineSearchCounts(), depending on whether you want these 2980c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 2990c52818fSSatish Balay 3000c52818fSSatish Balay Level: beginner 3010c52818fSSatish Balay 302db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()` 3030c52818fSSatish Balay @*/ 3040c52818fSSatish Balay 305e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 3060c52818fSSatish Balay { 3070c52818fSSatish Balay PetscInt low1,low2,low3,high1,high2,high3; 3080c52818fSSatish Balay 3090c52818fSSatish Balay PetscFunctionBegin; 3100c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 3110c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3123f6ba705SLisandro Dalcin PetscValidRealPointer(f,3); 3130c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 3140c52818fSSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,5); 3150c52818fSSatish Balay PetscValidPointer(reason,7); 3160c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 3170c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,g,4); 3180c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,s,5); 3199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low1, &high1)); 3209566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(g, &low2, &high2)); 3219566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(s, &low3, &high3)); 3223c859ba3SBarry Smith PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible vector local lengths"); 3230c52818fSSatish Balay 32497ab8e72SStefano Zampini *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3259566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s)); 3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->stepdirection)); 327050fc7a3SBarry Smith ls->stepdirection = s; 3280c52818fSSatish Balay 3299566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetUp(ls)); 3300c52818fSSatish Balay ls->nfeval = 0; 3310c52818fSSatish Balay ls->ngeval = 0; 3320c52818fSSatish Balay ls->nfgeval = 0; 3330c52818fSSatish Balay /* Check parameter values */ 3340c52818fSSatish Balay if (ls->ftol < 0.0) { 3359566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol)); 3360c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3370c52818fSSatish Balay } 3380c52818fSSatish Balay if (ls->rtol < 0.0) { 3399566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol)); 3400c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3410c52818fSSatish Balay } 3420c52818fSSatish Balay if (ls->gtol < 0.0) { 3439566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol)); 3440c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3450c52818fSSatish Balay } 3460c52818fSSatish Balay if (ls->stepmin < 0.0) { 3479566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin)); 3480c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3490c52818fSSatish Balay } 3500c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 3519566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax)); 3520c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3530c52818fSSatish Balay } 3540c52818fSSatish Balay if (ls->max_funcs < 0) { 3559566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n",ls->max_funcs)); 3560c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3570c52818fSSatish Balay } 3580c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 3599566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f)); 3600c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_INFORNAN; 3610c52818fSSatish Balay } 3620c52818fSSatish Balay 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)x)); 3649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->start_x)); 3650c52818fSSatish Balay ls->start_x = x; 366050fc7a3SBarry Smith 3679566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0)); 368*dbbe0bcdSBarry Smith PetscUseTypeMethod(ls,apply ,x,f,g,s); 3699566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0)); 3700c52818fSSatish Balay *reason = ls->reason; 3710c52818fSSatish Balay ls->new_f = *f; 3720c52818fSSatish Balay 37397ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 3740c52818fSSatish Balay 3759566063dSJacob Faibussowitsch PetscCall(TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view")); 3760c52818fSSatish Balay PetscFunctionReturn(0); 3770c52818fSSatish Balay } 3780c52818fSSatish Balay 3790c52818fSSatish Balay /*@C 3800c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 3810c52818fSSatish Balay 3820c52818fSSatish Balay Collective on TaoLineSearch 3830c52818fSSatish Balay 3840c52818fSSatish Balay Input Parameters: 3850c52818fSSatish Balay + ls - the TaoLineSearch context 386820824deSAlp Dener - type - the TaoLineSearchType selection 3870c52818fSSatish Balay 3880c52818fSSatish Balay Available methods include: 389820824deSAlp Dener + more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition 390820824deSAlp Dener . armijo - simple backtracking line search enforcing only the sufficient decrease condition 391820824deSAlp Dener - unit - do not perform a line search and always accept unit step length 3920c52818fSSatish Balay 3930c52818fSSatish Balay Options Database Keys: 394820824deSAlp Dener . -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime 3950c52818fSSatish Balay 3960c52818fSSatish Balay Level: beginner 3970c52818fSSatish Balay 398db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchGetType()`, `TaoLineSearchApply()` 3990c52818fSSatish Balay 4000c52818fSSatish Balay @*/ 4010c52818fSSatish Balay 402dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 4030c52818fSSatish Balay { 4040c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 4050c52818fSSatish Balay PetscBool flg; 4060c52818fSSatish Balay 4070c52818fSSatish Balay PetscFunctionBegin; 4080c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4090c52818fSSatish Balay PetscValidCharPointer(type,2); 4109566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg)); 4110c52818fSSatish Balay if (flg) PetscFunctionReturn(0); 4120c52818fSSatish Balay 4139566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r)); 4143c859ba3SBarry Smith PetscCheck(r,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 415*dbbe0bcdSBarry Smith PetscTryTypeMethod(ls,destroy); 4160c52818fSSatish Balay ls->max_funcs = 30; 4170c52818fSSatish Balay ls->ftol = 0.0001; 4180c52818fSSatish Balay ls->gtol = 0.9; 4196f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 4206f4723b1SBarry Smith ls->rtol = 1.0e-5; 4216f4723b1SBarry Smith #else 4220c52818fSSatish Balay ls->rtol = 1.0e-10; 4236f4723b1SBarry Smith #endif 4240c52818fSSatish Balay ls->stepmin = 1.0e-20; 4250c52818fSSatish Balay ls->stepmax = 1.0e+20; 4260c52818fSSatish Balay 4270c52818fSSatish Balay ls->nfeval = 0; 4280c52818fSSatish Balay ls->ngeval = 0; 4290c52818fSSatish Balay ls->nfgeval = 0; 43083c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 43183c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 43283c8fe1dSLisandro Dalcin ls->ops->view = NULL; 43383c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 43483c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 4350c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4369566063dSJacob Faibussowitsch PetscCall((*r)(ls)); 4379566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type)); 4380c52818fSSatish Balay PetscFunctionReturn(0); 4390c52818fSSatish Balay } 4400c52818fSSatish Balay 4412a0dac07SAlp Dener /*@C 4422a0dac07SAlp Dener TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the 4432a0dac07SAlp Dener iteration number, step length, and function value before calling the implementation 4442a0dac07SAlp Dener specific monitor. 4452a0dac07SAlp Dener 4462a0dac07SAlp Dener Input Parameters: 4472a0dac07SAlp Dener + ls - the TaoLineSearch context 4482a0dac07SAlp Dener . its - the current iterate number (>=0) 4492a0dac07SAlp Dener . f - the current objective function value 4502a0dac07SAlp Dener - step - the step length 4512a0dac07SAlp Dener 4522a0dac07SAlp Dener Options Database Key: 4532a0dac07SAlp Dener . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output 4542a0dac07SAlp Dener 4552a0dac07SAlp Dener Level: developer 4562a0dac07SAlp Dener 4572a0dac07SAlp Dener @*/ 4582a0dac07SAlp Dener PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step) 4592a0dac07SAlp Dener { 4602a0dac07SAlp Dener PetscInt tabs; 4612a0dac07SAlp Dener 4622a0dac07SAlp Dener PetscFunctionBegin; 4632a0dac07SAlp Dener PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4642a0dac07SAlp Dener if (ls->usemonitor) { 4659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs)); 4669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel)); 46763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its)); 4689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f)); 4699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step)); 4702a0dac07SAlp Dener if (ls->ops->monitor && its > 0) { 4719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3)); 472*dbbe0bcdSBarry Smith PetscUseTypeMethod(ls,monitor); 4732a0dac07SAlp Dener } 4749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs)); 4752a0dac07SAlp Dener } 4762a0dac07SAlp Dener PetscFunctionReturn(0); 4772a0dac07SAlp Dener } 4782a0dac07SAlp Dener 4790c52818fSSatish Balay /*@ 4800c52818fSSatish Balay TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 4810c52818fSSatish Balay options. 4820c52818fSSatish Balay 4830c52818fSSatish Balay Collective on TaoLineSearch 4840c52818fSSatish Balay 48501d2d390SJose E. Roman Input Parameter: 4860c52818fSSatish Balay . ls - the TaoLineSearch context 4870c52818fSSatish Balay 4880c52818fSSatish Balay Options Database Keys: 4890c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 4900c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 4910c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 4920c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 493a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - initial steplength allowed 4940c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 4950c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 4960c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 4970c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 4980c52818fSSatish Balay 4990c52818fSSatish Balay Level: beginner 5000c52818fSSatish Balay @*/ 5010c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 5020c52818fSSatish Balay { 5038caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 5042a0dac07SAlp Dener char type[256],monfilename[PETSC_MAX_PATH_LEN]; 5052a0dac07SAlp Dener PetscViewer monviewer; 5060c52818fSSatish Balay PetscBool flg; 507f06e3bfaSBarry Smith 5080c52818fSSatish Balay PetscFunctionBegin; 5090c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 510d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)ls); 5110c52818fSSatish Balay if (((PetscObject)ls)->type_name) { 5120c52818fSSatish Balay default_type = ((PetscObject)ls)->type_name; 5130c52818fSSatish Balay } 5140c52818fSSatish Balay /* Check for type from options */ 5159566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg)); 5160c52818fSSatish Balay if (flg) { 5179566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls,type)); 5180c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 5199566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(ls,default_type)); 5200c52818fSSatish Balay } 5210c52818fSSatish Balay 5229566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL)); 5239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL)); 5249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL)); 5259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL)); 5269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL)); 5279566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL)); 528a39c8e28SStefano Zampini PetscCall(PetscOptionsReal("-tao_ls_stepinit","initial step","",ls->initstep,&ls->initstep,NULL)); 5299566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,sizeof(monfilename),&flg)); 5302a0dac07SAlp Dener if (flg) { 5319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer)); 5322a0dac07SAlp Dener ls->viewer = monviewer; 5332a0dac07SAlp Dener ls->usemonitor = PETSC_TRUE; 5342a0dac07SAlp Dener } 535*dbbe0bcdSBarry Smith PetscTryTypeMethod(ls,setfromoptions,PetscOptionsObject); 536d0609cedSBarry Smith PetscOptionsEnd(); 5370c52818fSSatish Balay PetscFunctionReturn(0); 5380c52818fSSatish Balay } 5390c52818fSSatish Balay 5400c52818fSSatish Balay /*@C 5410c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5420c52818fSSatish Balay 5430c52818fSSatish Balay Not Collective 5440c52818fSSatish Balay 5450c52818fSSatish Balay Input Parameter: 5460c52818fSSatish Balay . ls - the TaoLineSearch context 5470c52818fSSatish Balay 548fd292e60Sprj- Output Parameter: 5490c52818fSSatish Balay . type - the line search algorithm in effect 5500c52818fSSatish Balay 5510c52818fSSatish Balay Level: developer 5520c52818fSSatish Balay 5530c52818fSSatish Balay @*/ 554dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 5550c52818fSSatish Balay { 5560c52818fSSatish Balay PetscFunctionBegin; 5570c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5580c52818fSSatish Balay PetscValidPointer(type,2); 5590c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 5600c52818fSSatish Balay PetscFunctionReturn(0); 5610c52818fSSatish Balay } 5620c52818fSSatish Balay 5630c52818fSSatish Balay /*@ 5640c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 5650c52818fSSatish Balay routines used by the line search in last application (not cumulative). 5660c52818fSSatish Balay 5670c52818fSSatish Balay Not Collective 5680c52818fSSatish Balay 5690c52818fSSatish Balay Input Parameter: 5700c52818fSSatish Balay . ls - the TaoLineSearch context 5710c52818fSSatish Balay 5720c52818fSSatish Balay Output Parameters: 5730c52818fSSatish Balay + nfeval - number of function evaluations 5740c52818fSSatish Balay . ngeval - number of gradient evaluations 5750c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 5760c52818fSSatish Balay 5770c52818fSSatish Balay Level: intermediate 5780c52818fSSatish Balay 5790c52818fSSatish Balay Note: 580441846f8SBarry Smith If the line search is using the Tao objective and gradient 581441846f8SBarry Smith routines directly (see TaoLineSearchUseTaoRoutines()), then TAO 5820c52818fSSatish Balay is already counting the number of evaluations. 5830c52818fSSatish Balay 5840c52818fSSatish Balay @*/ 5850c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 5860c52818fSSatish Balay { 5870c52818fSSatish Balay PetscFunctionBegin; 5880c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5890c52818fSSatish Balay *nfeval = ls->nfeval; 5900c52818fSSatish Balay *ngeval = ls->ngeval; 5910c52818fSSatish Balay *nfgeval = ls->nfgeval; 5920c52818fSSatish Balay PetscFunctionReturn(0); 5930c52818fSSatish Balay } 5940c52818fSSatish Balay 5950c52818fSSatish Balay /*@ 596441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 597441846f8SBarry Smith Tao evaluation routines. 5980c52818fSSatish Balay 5990c52818fSSatish Balay Not Collective 6000c52818fSSatish Balay 6010c52818fSSatish Balay Input Parameter: 6020c52818fSSatish Balay . ls - the TaoLineSearch context 6030c52818fSSatish Balay 6040c52818fSSatish Balay Output Parameter: 605441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines, 6060c52818fSSatish Balay otherwise PETSC_FALSE 6070c52818fSSatish Balay 6080c52818fSSatish Balay Level: developer 6090c52818fSSatish Balay @*/ 610441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 6110c52818fSSatish Balay { 6120c52818fSSatish Balay PetscFunctionBegin; 6130c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 614f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 6150c52818fSSatish Balay PetscFunctionReturn(0); 6160c52818fSSatish Balay } 6170c52818fSSatish Balay 6180c52818fSSatish Balay /*@C 6190c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 6200c52818fSSatish Balay 6210c52818fSSatish Balay Logically Collective on TaoLineSearch 6220c52818fSSatish Balay 623d8d19677SJose E. Roman Input Parameters: 6240c52818fSSatish Balay + ls - the TaoLineSearch context 6250c52818fSSatish Balay . func - the objective function evaluation routine 6260c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6270c52818fSSatish Balay 6280c52818fSSatish Balay Calling sequence of func: 6290c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 6300c52818fSSatish Balay 6310c52818fSSatish Balay + x - input vector 6320c52818fSSatish Balay . f - function value 6330c52818fSSatish Balay - ctx (optional) user-defined context 6340c52818fSSatish Balay 6350c52818fSSatish Balay Level: beginner 6360c52818fSSatish Balay 6370c52818fSSatish Balay Note: 6380c52818fSSatish Balay Use this routine only if you want the line search objective 639441846f8SBarry Smith evaluation routine to be different from the Tao's objective 6400c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6410c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6420c52818fSSatish Balay 6430c52818fSSatish Balay Note: 6440c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6450c52818fSSatish Balay line search, application programmers should be wary of overriding the 6460c52818fSSatish Balay default objective routine. 6470c52818fSSatish Balay 648db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 6490c52818fSSatish Balay @*/ 6500c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx) 6510c52818fSSatish Balay { 6520c52818fSSatish Balay PetscFunctionBegin; 6530c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6540c52818fSSatish Balay 6550c52818fSSatish Balay ls->ops->computeobjective=func; 6560c52818fSSatish Balay if (ctx) ls->userctx_func=ctx; 6570c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 6580c52818fSSatish Balay PetscFunctionReturn(0); 6590c52818fSSatish Balay } 6600c52818fSSatish Balay 6610c52818fSSatish Balay /*@C 6620c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 6630c52818fSSatish Balay 6640c52818fSSatish Balay Logically Collective on TaoLineSearch 6650c52818fSSatish Balay 666d8d19677SJose E. Roman Input Parameters: 6670c52818fSSatish Balay + ls - the TaoLineSearch context 6680c52818fSSatish Balay . func - the gradient evaluation routine 6690c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6700c52818fSSatish Balay 6710c52818fSSatish Balay Calling sequence of func: 6720c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, Vec g, void *ctx); 6730c52818fSSatish Balay 6740c52818fSSatish Balay + x - input vector 6750c52818fSSatish Balay . g - gradient vector 6760c52818fSSatish Balay - ctx (optional) user-defined context 6770c52818fSSatish Balay 6780c52818fSSatish Balay Level: beginner 6790c52818fSSatish Balay 6800c52818fSSatish Balay Note: 6810c52818fSSatish Balay Use this routine only if you want the line search gradient 682441846f8SBarry Smith evaluation routine to be different from the Tao's gradient 6830c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6840c52818fSSatish Balay the line search function and/or function/gradient routine. 6850c52818fSSatish Balay 6860c52818fSSatish Balay Note: 6870c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 6880c52818fSSatish Balay line search, application programmers should be wary of overriding the 6890c52818fSSatish Balay default gradient routine. 6900c52818fSSatish Balay 691db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 6920c52818fSSatish Balay @*/ 6930c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx) 6940c52818fSSatish Balay { 6950c52818fSSatish Balay PetscFunctionBegin; 6960c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6970c52818fSSatish Balay ls->ops->computegradient=func; 6980c52818fSSatish Balay if (ctx) ls->userctx_grad=ctx; 6990c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 7000c52818fSSatish Balay PetscFunctionReturn(0); 7010c52818fSSatish Balay } 7020c52818fSSatish Balay 7030c52818fSSatish Balay /*@C 7040c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 7050c52818fSSatish Balay 7060c52818fSSatish Balay Logically Collective on TaoLineSearch 7070c52818fSSatish Balay 708d8d19677SJose E. Roman Input Parameters: 7090c52818fSSatish Balay + ls - the TaoLineSearch context 7100c52818fSSatish Balay . func - the objective and gradient evaluation routine 7110c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7120c52818fSSatish Balay 7130c52818fSSatish Balay Calling sequence of func: 7140c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 7150c52818fSSatish Balay 7160c52818fSSatish Balay + x - input vector 7170c52818fSSatish Balay . f - function value 7180c52818fSSatish Balay . g - gradient vector 7190c52818fSSatish Balay - ctx (optional) user-defined context 7200c52818fSSatish Balay 7210c52818fSSatish Balay Level: beginner 7220c52818fSSatish Balay 7230c52818fSSatish Balay Note: 7240c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 725441846f8SBarry Smith evaluation routines to be different from the Tao's objective 7260c52818fSSatish Balay and gradient evaluation routines. 7270c52818fSSatish Balay 7280c52818fSSatish Balay Note: 7290c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7300c52818fSSatish Balay line search, application programmers should be wary of overriding the 7310c52818fSSatish Balay default objective routine. 7320c52818fSSatish Balay 733db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()` 7340c52818fSSatish Balay @*/ 7350c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx) 7360c52818fSSatish Balay { 7370c52818fSSatish Balay PetscFunctionBegin; 7380c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7390c52818fSSatish Balay ls->ops->computeobjectiveandgradient=func; 7400c52818fSSatish Balay if (ctx) ls->userctx_funcgrad=ctx; 7410c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7420c52818fSSatish Balay PetscFunctionReturn(0); 7430c52818fSSatish Balay } 7440c52818fSSatish Balay 7450c52818fSSatish Balay /*@C 7460c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 7470c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 7480c52818fSSatish Balay Sometimes it is more efficient to compute the inner product of the gradient 7490c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 7500c52818fSSatish Balay the line search typically needs of the gradient. 7510c52818fSSatish Balay 7520c52818fSSatish Balay Logically Collective on TaoLineSearch 7530c52818fSSatish Balay 754d8d19677SJose E. Roman Input Parameters: 7550c52818fSSatish Balay + ls - the TaoLineSearch context 7560c52818fSSatish Balay . func - the objective and gradient evaluation routine 7570c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7580c52818fSSatish Balay 7590c52818fSSatish Balay Calling sequence of func: 7600c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 7610c52818fSSatish Balay 7620c52818fSSatish Balay + x - input vector 7630c52818fSSatish Balay . s - step direction 7640c52818fSSatish Balay . f - function value 7650c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 7660c52818fSSatish Balay - ctx (optional) user-defined context 7670c52818fSSatish Balay 7680c52818fSSatish Balay Note: The gradient will still need to be computed at the end of the line 7690c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 7700c52818fSSatish Balay routine 7710c52818fSSatish Balay 7720c52818fSSatish Balay Note: Bounded line searches (those used in bounded optimization algorithms) 7730c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 7740c52818fSSatish Balay x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength() 7750c52818fSSatish Balay 7760c52818fSSatish Balay Level: advanced 7770c52818fSSatish Balay 7780c52818fSSatish Balay Note: 7790c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7800c52818fSSatish Balay line search, application programmers should be wary of overriding the 7810c52818fSSatish Balay default objective routine. 7820c52818fSSatish Balay 783db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()` 7840c52818fSSatish Balay @*/ 7850c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx) 7860c52818fSSatish Balay { 7870c52818fSSatish Balay PetscFunctionBegin; 7880c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7890c52818fSSatish Balay ls->ops->computeobjectiveandgts=func; 7900c52818fSSatish Balay if (ctx) ls->userctx_funcgts=ctx; 7910c52818fSSatish Balay ls->usegts = PETSC_TRUE; 7920c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7930c52818fSSatish Balay PetscFunctionReturn(0); 7940c52818fSSatish Balay } 7950c52818fSSatish Balay 7960c52818fSSatish Balay /*@C 797441846f8SBarry Smith TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the 798441846f8SBarry Smith objective and gradient evaluation routines from the given Tao object. 7990c52818fSSatish Balay 8000c52818fSSatish Balay Logically Collective on TaoLineSearch 8010c52818fSSatish Balay 802d8d19677SJose E. Roman Input Parameters: 8030c52818fSSatish Balay + ls - the TaoLineSearch context 804441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines 8050c52818fSSatish Balay 8060c52818fSSatish Balay Level: developer 8070c52818fSSatish Balay 808db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()` 8090c52818fSSatish Balay @*/ 810441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts) 8110c52818fSSatish Balay { 8120c52818fSSatish Balay PetscFunctionBegin; 8130c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 814064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(ts,TAO_CLASSID,2); 815441846f8SBarry Smith ls->tao = ts; 8160c52818fSSatish Balay ls->usetaoroutines = PETSC_TRUE; 8170c52818fSSatish Balay PetscFunctionReturn(0); 8180c52818fSSatish Balay } 8190c52818fSSatish Balay 8200c52818fSSatish Balay /*@ 8210c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 8220c52818fSSatish Balay 8230c52818fSSatish Balay Collective on TaoLineSearch 8240c52818fSSatish Balay 8250c52818fSSatish Balay Input Parameters: 8260c52818fSSatish Balay + ls - the TaoLineSearch context 8270c52818fSSatish Balay - x - input vector 8280c52818fSSatish Balay 8290c52818fSSatish Balay Output Parameter: 8300c52818fSSatish Balay . f - Objective value at X 8310c52818fSSatish Balay 83295452b02SPatrick Sanan Notes: 83395452b02SPatrick Sanan TaoLineSearchComputeObjective() is typically used within line searches 8340c52818fSSatish Balay so most users would not generally call this routine themselves. 8350c52818fSSatish Balay 8360c52818fSSatish Balay Level: developer 8370c52818fSSatish Balay 838db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()` 8390c52818fSSatish Balay @*/ 8400c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 8410c52818fSSatish Balay { 8420c52818fSSatish Balay Vec gdummy; 8430c52818fSSatish Balay PetscReal gts; 844f06e3bfaSBarry Smith 8450c52818fSSatish Balay PetscFunctionBegin; 8460c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8470c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 848dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 8490c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8500c52818fSSatish Balay if (ls->usetaoroutines) { 8519566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(ls->tao,x,f)); 8520c52818fSSatish Balay } else { 8533c859ba3SBarry 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"); 8549566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 855792fecdfSBarry Smith if (ls->ops->computeobjective) PetscCallBack("TaoLineSearch callback objective",(*ls->ops->computeobjective)(ls,x,f,ls->userctx_func)); 856ef1023bdSBarry Smith else if (ls->ops->computeobjectiveandgradient) { 8579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&gdummy)); 858792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback objective",(*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad)); 8599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gdummy)); 860792fecdfSBarry Smith } else PetscCallBack("TaoLineSearch callback objective",(*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts)); 8619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 8620c52818fSSatish Balay } 8630c52818fSSatish Balay ls->nfeval++; 8640c52818fSSatish Balay PetscFunctionReturn(0); 8650c52818fSSatish Balay } 8660c52818fSSatish Balay 8670c52818fSSatish Balay /*@ 8680c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 8690c52818fSSatish Balay 870441846f8SBarry Smith Collective on Tao 8710c52818fSSatish Balay 8720c52818fSSatish Balay Input Parameters: 8730c52818fSSatish Balay + ls - the TaoLineSearch context 8740c52818fSSatish Balay - x - input vector 8750c52818fSSatish Balay 876d8d19677SJose E. Roman Output Parameters: 8770c52818fSSatish Balay + f - Objective value at X 8780c52818fSSatish Balay - g - Gradient vector at X 8790c52818fSSatish Balay 88095452b02SPatrick Sanan Notes: 88195452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 8820c52818fSSatish Balay so most users would not generally call this routine themselves. 8830c52818fSSatish Balay 8840c52818fSSatish Balay Level: developer 8850c52818fSSatish Balay 886db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()` 8870c52818fSSatish Balay @*/ 8880c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 8890c52818fSSatish Balay { 8900c52818fSSatish Balay PetscFunctionBegin; 8910c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8920c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 893dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 8940c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 8950c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8960c52818fSSatish Balay PetscCheckSameComm(ls,1,g,4); 8970c52818fSSatish Balay if (ls->usetaoroutines) { 8989566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(ls->tao,x,f,g)); 8990c52818fSSatish Balay } else { 9009566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 901792fecdfSBarry Smith if (ls->ops->computeobjectiveandgradient) PetscCallBack("TaoLineSearch callback objective/gradient",(*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad)); 902ef1023bdSBarry Smith else { 903792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback objective",(*ls->ops->computeobjective)(ls,x,f,ls->userctx_func)); 904792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback gradient",(*ls->ops->computegradient)(ls,x,g,ls->userctx_grad)); 905362febeeSStefano Zampini } 9069566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 9079566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f))); 908f06e3bfaSBarry Smith } 909fbe0838dSJason Sarich ls->nfgeval++; 9100c52818fSSatish Balay PetscFunctionReturn(0); 9110c52818fSSatish Balay } 9120c52818fSSatish Balay 9130c52818fSSatish Balay /*@ 9140c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 9150c52818fSSatish Balay 9160c52818fSSatish Balay Collective on TaoLineSearch 9170c52818fSSatish Balay 9180c52818fSSatish Balay Input Parameters: 9190c52818fSSatish Balay + ls - the TaoLineSearch context 9200c52818fSSatish Balay - x - input vector 9210c52818fSSatish Balay 9220c52818fSSatish Balay Output Parameter: 9230c52818fSSatish Balay . g - gradient vector 9240c52818fSSatish Balay 92595452b02SPatrick Sanan Notes: 92695452b02SPatrick Sanan TaoComputeGradient() is typically used within line searches 9270c52818fSSatish Balay so most users would not generally call this routine themselves. 9280c52818fSSatish Balay 9290c52818fSSatish Balay Level: developer 9300c52818fSSatish Balay 931db781477SPatrick Sanan .seealso: `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()` 9320c52818fSSatish Balay @*/ 9330c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 9340c52818fSSatish Balay { 9350c52818fSSatish Balay PetscReal fdummy; 936f06e3bfaSBarry Smith 9370c52818fSSatish Balay PetscFunctionBegin; 9380c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9390c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9400c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,3); 9410c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9420c52818fSSatish Balay PetscCheckSameComm(ls,1,g,3); 9430c52818fSSatish Balay if (ls->usetaoroutines) { 9449566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(ls->tao,x,g)); 9450c52818fSSatish Balay } else { 9469566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 947792fecdfSBarry Smith if (ls->ops->computegradient) PetscCallBack("TaoLineSearch callback gradient",(*ls->ops->computegradient)(ls,x,g,ls->userctx_grad)); 948ef1023bdSBarry Smith else { 949792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback gradient",(*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad)); 9500c52818fSSatish Balay } 9519566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 9520c52818fSSatish Balay } 9530c52818fSSatish Balay ls->ngeval++; 9540c52818fSSatish Balay PetscFunctionReturn(0); 9550c52818fSSatish Balay } 9560c52818fSSatish Balay 9570c52818fSSatish Balay /*@ 9580c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 9590c52818fSSatish Balay 960441846f8SBarry Smith Collective on Tao 9610c52818fSSatish Balay 9620c52818fSSatish Balay Input Parameters: 9630c52818fSSatish Balay + ls - the TaoLineSearch context 9640c52818fSSatish Balay - x - input vector 9650c52818fSSatish Balay 966d8d19677SJose E. Roman Output Parameters: 9670c52818fSSatish Balay + f - Objective value at X 9680c52818fSSatish Balay - gts - inner product of gradient and step direction at X 9690c52818fSSatish Balay 97095452b02SPatrick Sanan Notes: 97195452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 9720c52818fSSatish Balay so most users would not generally call this routine themselves. 9730c52818fSSatish Balay 9740c52818fSSatish Balay Level: developer 9750c52818fSSatish Balay 976db781477SPatrick Sanan .seealso: `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()` 9770c52818fSSatish Balay @*/ 9780c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 9790c52818fSSatish Balay { 9800c52818fSSatish Balay PetscFunctionBegin; 9810c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9820c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 983dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 984dadcf809SJacob Faibussowitsch PetscValidRealPointer(gts,4); 9850c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 987792fecdfSBarry Smith PetscCallBack("TaoLineSearch callback objective/gts",(*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts)); 9889566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 9899566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f))); 9900c52818fSSatish Balay ls->nfeval++; 9910c52818fSSatish Balay PetscFunctionReturn(0); 9920c52818fSSatish Balay } 9930c52818fSSatish Balay 9940c52818fSSatish Balay /*@ 9950c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 9960c52818fSSatish Balay 9970c52818fSSatish Balay Collective on TaoLineSearch 9980c52818fSSatish Balay 9990c52818fSSatish Balay Input Parameter: 10000c52818fSSatish Balay . ls - the TaoLineSearch context 10010c52818fSSatish Balay 1002d8d19677SJose E. Roman Output Parameters: 10030c52818fSSatish Balay + x - the new solution 10040c52818fSSatish Balay . f - the objective function value at x 10050c52818fSSatish Balay . g - the gradient at x 10060c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 10070c52818fSSatish Balay - reason - the reason why the line search terminated 10080c52818fSSatish Balay 10090c52818fSSatish Balay reason will be set to one of: 10100c52818fSSatish Balay 10110c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 10120c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 10130c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 10140c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 10150c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 10160c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 10170c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 10180c52818fSSatish Balay 10190c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 10200c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 10210c52818fSSatish Balay 1022a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search 10230c52818fSSatish Balay 10240c52818fSSatish Balay Level: developer 10250c52818fSSatish Balay 10260c52818fSSatish Balay @*/ 1027e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 10280c52818fSSatish Balay { 10290c52818fSSatish Balay PetscFunctionBegin; 10300c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10310c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1032dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 10330c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 10340c52818fSSatish Balay PetscValidIntPointer(reason,6); 10351baa6e33SBarry Smith if (ls->new_x) PetscCall(VecCopy(ls->new_x,x)); 10360c52818fSSatish Balay *f = ls->new_f; 10371baa6e33SBarry Smith if (ls->new_g) PetscCall(VecCopy(ls->new_g,g)); 103897ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 10390c52818fSSatish Balay *reason = ls->reason; 10400c52818fSSatish Balay PetscFunctionReturn(0); 10410c52818fSSatish Balay } 10420c52818fSSatish Balay 10430c52818fSSatish Balay /*@ 10440c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 10450c52818fSSatish Balay search. 10460c52818fSSatish Balay 10470c52818fSSatish Balay Not Collective 10480c52818fSSatish Balay 10490c52818fSSatish Balay Input Parameter: 10500c52818fSSatish Balay . ls - the TaoLineSearch context 10510c52818fSSatish Balay 10520c52818fSSatish Balay Output Parameter: 10530c52818fSSatish Balay . x - The initial point of the line search 10540c52818fSSatish Balay 10550c52818fSSatish Balay Level: intermediate 10560c52818fSSatish Balay @*/ 10570c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 10580c52818fSSatish Balay { 10590c52818fSSatish Balay PetscFunctionBegin; 10600c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 106197ab8e72SStefano Zampini if (x) *x = ls->start_x; 10620c52818fSSatish Balay PetscFunctionReturn(0); 10630c52818fSSatish Balay } 10640c52818fSSatish Balay 10650c52818fSSatish Balay /*@ 10660c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 10670c52818fSSatish Balay search. 10680c52818fSSatish Balay 10690c52818fSSatish Balay Not Collective 10700c52818fSSatish Balay 10710c52818fSSatish Balay Input Parameter: 10720c52818fSSatish Balay . ls - the TaoLineSearch context 10730c52818fSSatish Balay 10740c52818fSSatish Balay Output Parameter: 10750c52818fSSatish Balay . s - the step direction of the line search 10760c52818fSSatish Balay 10770c52818fSSatish Balay Level: advanced 10780c52818fSSatish Balay @*/ 10790c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 10800c52818fSSatish Balay { 10810c52818fSSatish Balay PetscFunctionBegin; 10820c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 108397ab8e72SStefano Zampini if (s) *s = ls->stepdirection; 10840c52818fSSatish Balay PetscFunctionReturn(0); 10850c52818fSSatish Balay } 10860c52818fSSatish Balay 10870c52818fSSatish Balay /*@ 10880c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 10890c52818fSSatish Balay 10900c52818fSSatish Balay Not Collective 10910c52818fSSatish Balay 10920c52818fSSatish Balay Input Parameter: 10930c52818fSSatish Balay . ls - the TaoLineSearch context 10940c52818fSSatish Balay 10950c52818fSSatish Balay Output Parameter: 10960c52818fSSatish Balay . f - the objective value at the full step length 10970c52818fSSatish Balay 10980c52818fSSatish Balay Level: developer 10990c52818fSSatish Balay @*/ 11000c52818fSSatish Balay 11010c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 11020c52818fSSatish Balay { 11030c52818fSSatish Balay PetscFunctionBegin; 11040c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11050c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 11060c52818fSSatish Balay PetscFunctionReturn(0); 11070c52818fSSatish Balay } 11080c52818fSSatish Balay 11090c52818fSSatish Balay /*@ 11100c52818fSSatish Balay TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 11110c52818fSSatish Balay 1112441846f8SBarry Smith Logically Collective on Tao 11130c52818fSSatish Balay 11140c52818fSSatish Balay Input Parameters: 11150c52818fSSatish Balay + ls - the TaoLineSearch context 11160c52818fSSatish Balay . xl - vector of lower bounds 11170c52818fSSatish Balay - xu - vector of upper bounds 11180c52818fSSatish Balay 11190c52818fSSatish Balay Note: If the variable bounds are not set with this routine, then 1120e270355aSBarry Smith PETSC_NINFINITY and PETSC_INFINITY are assumed 11210c52818fSSatish Balay 11220c52818fSSatish Balay Level: beginner 11230c52818fSSatish Balay 1124db781477SPatrick Sanan .seealso: `TaoSetVariableBounds()`, `TaoLineSearchCreate()` 11250c52818fSSatish Balay @*/ 11260c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 11270c52818fSSatish Balay { 11280c52818fSSatish Balay PetscFunctionBegin; 11290c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 113076be6f4fSStefano Zampini if (xl) PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 113176be6f4fSStefano Zampini if (xu) PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 11329566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xl)); 11339566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xu)); 11349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->lower)); 11359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ls->upper)); 11360c52818fSSatish Balay ls->lower = xl; 11370c52818fSSatish Balay ls->upper = xu; 113876be6f4fSStefano Zampini ls->bounded = (PetscBool)(xl || xu); 11390c52818fSSatish Balay PetscFunctionReturn(0); 11400c52818fSSatish Balay } 11410c52818fSSatish Balay 11420c52818fSSatish Balay /*@ 11430c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 11440c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 11450c52818fSSatish Balay 11460c52818fSSatish Balay Logically Collective on TaoLineSearch 11470c52818fSSatish Balay 11480c52818fSSatish Balay Input Parameters: 11490c52818fSSatish Balay + ls - the TaoLineSearch context 11500c52818fSSatish Balay - s - the initial step size 11510c52818fSSatish Balay 11520c52818fSSatish Balay Level: intermediate 11530c52818fSSatish Balay 1154db781477SPatrick Sanan .seealso: `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()` 11550c52818fSSatish Balay @*/ 11560c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 11570c52818fSSatish Balay { 11580c52818fSSatish Balay PetscFunctionBegin; 11590c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1160a39c8e28SStefano Zampini PetscValidLogicalCollectiveReal(ls,s,2); 11610c52818fSSatish Balay ls->initstep = s; 11620c52818fSSatish Balay PetscFunctionReturn(0); 11630c52818fSSatish Balay } 11640c52818fSSatish Balay 11650c52818fSSatish Balay /*@ 11660c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 11670c52818fSSatish Balay 11680c52818fSSatish Balay Not Collective 11690c52818fSSatish Balay 11700c52818fSSatish Balay Input Parameters: 11710c52818fSSatish Balay . ls - the TaoLineSearch context 11720c52818fSSatish Balay 11737a7aea1fSJed Brown Output Parameters: 11740c52818fSSatish Balay . s - the current step length 11750c52818fSSatish Balay 11760c52818fSSatish Balay Level: beginner 11770c52818fSSatish Balay 1178db781477SPatrick Sanan .seealso: `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()` 11790c52818fSSatish Balay @*/ 11800c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 11810c52818fSSatish Balay { 11820c52818fSSatish Balay PetscFunctionBegin; 11830c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11840c52818fSSatish Balay *s = ls->step; 11850c52818fSSatish Balay PetscFunctionReturn(0); 11860c52818fSSatish Balay } 11870c52818fSSatish Balay 11880ebee16dSLisandro Dalcin /*@C 11890c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 11900c52818fSSatish Balay 11910c52818fSSatish Balay Not collective 11920c52818fSSatish Balay 11930c52818fSSatish Balay Input Parameters: 11940c52818fSSatish Balay + sname - name of a new user-defined solver 11950c52818fSSatish Balay - func - routine to Create method context 11960c52818fSSatish Balay 11970c52818fSSatish Balay Notes: 11980c52818fSSatish Balay TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 11990c52818fSSatish Balay 12000c52818fSSatish Balay Sample usage: 12010c52818fSSatish Balay .vb 12020c52818fSSatish Balay TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 12030c52818fSSatish Balay .ve 12040c52818fSSatish Balay 12050c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 12060c52818fSSatish Balay $ TaoLineSearchSetType(ls,"my_linesearch") 12070c52818fSSatish Balay or at runtime via the option 12080c52818fSSatish Balay $ -tao_ls_type my_linesearch 12090c52818fSSatish Balay 12100c52818fSSatish Balay Level: developer 12110c52818fSSatish Balay 12120ebee16dSLisandro Dalcin @*/ 12130c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 12140c52818fSSatish Balay { 12150c52818fSSatish Balay PetscFunctionBegin; 12169566063dSJacob Faibussowitsch PetscCall(TaoLineSearchInitializePackage()); 12179566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func)); 12180c52818fSSatish Balay PetscFunctionReturn(0); 12190c52818fSSatish Balay } 12200c52818fSSatish Balay 12210c52818fSSatish Balay /*@C 12220c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 12230c52818fSSatish Balay for all TaoLineSearch options in the database. 12240c52818fSSatish Balay 12250c52818fSSatish Balay Collective on TaoLineSearch 12260c52818fSSatish Balay 12270c52818fSSatish Balay Input Parameters: 12280c52818fSSatish Balay + ls - the TaoLineSearch solver context 12290c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 12300c52818fSSatish Balay 12310c52818fSSatish Balay Notes: 12320c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12330c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12340c52818fSSatish Balay 12350c52818fSSatish Balay Level: advanced 12360c52818fSSatish Balay 1237db781477SPatrick Sanan .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()` 12380c52818fSSatish Balay @*/ 12390c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 12400c52818fSSatish Balay { 1241f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 12420c52818fSSatish Balay } 12430c52818fSSatish Balay 12440c52818fSSatish Balay /*@C 12450c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 12460c52818fSSatish Balay TaoLineSearch options in the database 12470c52818fSSatish Balay 12480c52818fSSatish Balay Not Collective 12490c52818fSSatish Balay 12500c52818fSSatish Balay Input Parameters: 12510c52818fSSatish Balay . ls - the TaoLineSearch context 12520c52818fSSatish Balay 12530c52818fSSatish Balay Output Parameters: 12540c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 12550c52818fSSatish Balay 125695452b02SPatrick Sanan Notes: 125795452b02SPatrick Sanan On the fortran side, the user should pass in a string 'prefix' of 12580c52818fSSatish Balay sufficient length to hold the prefix. 12590c52818fSSatish Balay 12600c52818fSSatish Balay Level: advanced 12610c52818fSSatish Balay 1262db781477SPatrick Sanan .seealso: `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()` 12630c52818fSSatish Balay @*/ 12640c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 12650c52818fSSatish Balay { 1266f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 12670c52818fSSatish Balay } 12680c52818fSSatish Balay 12690c52818fSSatish Balay /*@C 12700c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 12710c52818fSSatish Balay TaoLineSearch options in the database. 12720c52818fSSatish Balay 12730c52818fSSatish Balay Logically Collective on TaoLineSearch 12740c52818fSSatish Balay 12750c52818fSSatish Balay Input Parameters: 12760c52818fSSatish Balay + ls - the TaoLineSearch context 12770c52818fSSatish Balay - prefix - the prefix string to prepend to all TAO option requests 12780c52818fSSatish Balay 12790c52818fSSatish Balay Notes: 12800c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12810c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12820c52818fSSatish Balay 12830c52818fSSatish Balay For example, to distinguish between the runtime options for two 12840c52818fSSatish Balay different line searches, one could call 12850c52818fSSatish Balay .vb 12860c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 12870c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 12880c52818fSSatish Balay .ve 12890c52818fSSatish Balay 12900c52818fSSatish Balay This would enable use of different options for each system, such as 12910c52818fSSatish Balay .vb 12920c52818fSSatish Balay -sys1_tao_ls_type mt 12930c52818fSSatish Balay -sys2_tao_ls_type armijo 12940c52818fSSatish Balay .ve 12950c52818fSSatish Balay 12960c52818fSSatish Balay Level: advanced 12970c52818fSSatish Balay 1298db781477SPatrick Sanan .seealso: `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()` 12990c52818fSSatish Balay @*/ 13000c52818fSSatish Balay 13010c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 13020c52818fSSatish Balay { 1303f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 13040c52818fSSatish Balay } 1305