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 120c52818fSSatish Balay TaoLineSearchView - Prints information about the TaoLineSearch 130c52818fSSatish Balay 140c52818fSSatish Balay Collective on TaoLineSearch 150c52818fSSatish Balay 160c52818fSSatish Balay InputParameters: 17441846f8SBarry Smith + ls - the Tao context 180c52818fSSatish Balay - viewer - visualization context 190c52818fSSatish Balay 200c52818fSSatish Balay Options Database Key: 210c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search 220c52818fSSatish Balay 230c52818fSSatish Balay Notes: 240c52818fSSatish Balay The available visualization contexts include 250c52818fSSatish Balay + PETSC_VIEWER_STDOUT_SELF - standard output (default) 260c52818fSSatish Balay - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 270c52818fSSatish Balay output where only the first processor opens 280c52818fSSatish Balay the file. All other processors send their 290c52818fSSatish Balay data to the first processor to print. 300c52818fSSatish Balay 310c52818fSSatish Balay Level: beginner 320c52818fSSatish Balay 330c52818fSSatish Balay .seealso: PetscViewerASCIIOpen() 340c52818fSSatish Balay @*/ 350c52818fSSatish Balay 360c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 370c52818fSSatish Balay { 380c52818fSSatish Balay PetscErrorCode ierr; 390c52818fSSatish Balay PetscBool isascii, isstring; 40dedfbcbeSJed Brown TaoLineSearchType type; 41f06e3bfaSBarry Smith 420c52818fSSatish Balay PetscFunctionBegin; 430c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 440c52818fSSatish Balay if (!viewer) { 450c52818fSSatish Balay ierr = PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);CHKERRQ(ierr); 460c52818fSSatish Balay } 470c52818fSSatish Balay PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 480c52818fSSatish Balay PetscCheckSameComm(ls,1,viewer,2); 490c52818fSSatish Balay 500c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 510c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 520c52818fSSatish Balay if (isascii) { 5363b15415SAlp Dener ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);CHKERRQ(ierr); 5463b15415SAlp Dener if (ls->ops->view) { 5563b15415SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 5663b15415SAlp Dener ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr); 5763b15415SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 580c52818fSSatish Balay } 590c52818fSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 600c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr); 61f06e3bfaSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);CHKERRQ(ierr); 620c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr); 630c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr); 640c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr); 650c52818fSSatish Balay 660c52818fSSatish Balay if (ls->bounded) { 670c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr); 680c52818fSSatish Balay } 69de196d00SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);CHKERRQ(ierr); 700c52818fSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 710c52818fSSatish Balay } else if (isstring) { 720c52818fSSatish Balay ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr); 730c52818fSSatish Balay ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 740c52818fSSatish Balay } 750c52818fSSatish Balay PetscFunctionReturn(0); 760c52818fSSatish Balay } 770c52818fSSatish Balay 780c52818fSSatish Balay /*@C 790c52818fSSatish Balay TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 800c52818fSSatish Balay line-searches will automatically create one. 810c52818fSSatish Balay 820c52818fSSatish Balay Collective on MPI_Comm 830c52818fSSatish Balay 840c52818fSSatish Balay Input Parameter: 850c52818fSSatish Balay . comm - MPI communicator 860c52818fSSatish Balay 870c52818fSSatish Balay Output Parameter: 880c52818fSSatish Balay . newls - the new TaoLineSearch context 890c52818fSSatish Balay 900c52818fSSatish Balay Available methods include: 910c52818fSSatish Balay + more-thuente 920c52818fSSatish Balay . gpcg 930c52818fSSatish Balay - unit - Do not perform any line search 940c52818fSSatish Balay 950c52818fSSatish Balay 960c52818fSSatish Balay Options Database Keys: 970c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 980c52818fSSatish Balay 990c52818fSSatish Balay Level: beginner 1000c52818fSSatish Balay 1010c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy() 1020c52818fSSatish Balay @*/ 1030c52818fSSatish Balay 1040c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 1050c52818fSSatish Balay { 1060c52818fSSatish Balay PetscErrorCode ierr; 1070c52818fSSatish Balay TaoLineSearch ls; 1080c52818fSSatish Balay 1090c52818fSSatish Balay PetscFunctionBegin; 1100c52818fSSatish Balay PetscValidPointer(newls,2); 1116c23d075SBarry Smith *newls = NULL; 1120c52818fSSatish Balay 1130c52818fSSatish Balay ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr); 1140c52818fSSatish Balay 11573107ff1SLisandro Dalcin ierr = PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);CHKERRQ(ierr); 1160c52818fSSatish Balay ls->bounded = 0; 1170c52818fSSatish Balay ls->max_funcs=30; 1180c52818fSSatish Balay ls->ftol = 0.0001; 1190c52818fSSatish Balay ls->gtol = 0.9; 1206f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1216f4723b1SBarry Smith ls->rtol = 1.0e-5; 1226f4723b1SBarry Smith #else 1230c52818fSSatish Balay ls->rtol = 1.0e-10; 1246f4723b1SBarry Smith #endif 1250c52818fSSatish Balay ls->stepmin=1.0e-20; 1260c52818fSSatish Balay ls->stepmax=1.0e+20; 1270c52818fSSatish Balay ls->step=1.0; 1280c52818fSSatish Balay ls->nfeval=0; 1290c52818fSSatish Balay ls->ngeval=0; 1300c52818fSSatish Balay ls->nfgeval=0; 1310c52818fSSatish Balay 1320c52818fSSatish Balay ls->ops->computeobjective=0; 1330c52818fSSatish Balay ls->ops->computegradient=0; 1340c52818fSSatish Balay ls->ops->computeobjectiveandgradient=0; 1350c52818fSSatish Balay ls->ops->computeobjectiveandgts=0; 1360c52818fSSatish Balay ls->ops->setup=0; 1370c52818fSSatish Balay ls->ops->apply=0; 1380c52818fSSatish Balay ls->ops->view=0; 1390c52818fSSatish Balay ls->ops->setfromoptions=0; 1400c52818fSSatish Balay ls->ops->reset=0; 1410c52818fSSatish Balay ls->ops->destroy=0; 1422a0dac07SAlp Dener ls->ops->monitor=0; 1432a0dac07SAlp Dener ls->usemonitor=PETSC_FALSE; 1440c52818fSSatish Balay ls->setupcalled=PETSC_FALSE; 1450c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 1460c52818fSSatish Balay *newls = ls; 1470c52818fSSatish Balay PetscFunctionReturn(0); 1480c52818fSSatish Balay } 1490c52818fSSatish Balay 1500c52818fSSatish Balay /*@ 1510c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 1520c52818fSSatish Balay of a Tao solver 1530c52818fSSatish Balay 1540c52818fSSatish Balay Collective on ls 1550c52818fSSatish Balay 1560c52818fSSatish Balay Input Parameters: 1570c52818fSSatish Balay . ls - the TaoLineSearch context 1580c52818fSSatish Balay 1590c52818fSSatish Balay Notes: 1600c52818fSSatish Balay The user will not need to explicitly call TaoLineSearchSetUp(), as it will 1610c52818fSSatish Balay automatically be called in TaoLineSearchSolve(). However, if the user 1620c52818fSSatish Balay desires to call it explicitly, it should come after TaoLineSearchCreate() 1630c52818fSSatish Balay but before TaoLineSearchApply(). 1640c52818fSSatish Balay 1650c52818fSSatish Balay Level: developer 1660c52818fSSatish Balay 1670c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 1680c52818fSSatish Balay @*/ 1690c52818fSSatish Balay 1700c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 1710c52818fSSatish Balay { 1720c52818fSSatish Balay PetscErrorCode ierr; 1738caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 1740c52818fSSatish Balay PetscBool flg; 175f06e3bfaSBarry Smith 1760c52818fSSatish Balay PetscFunctionBegin; 1770c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1780c52818fSSatish Balay if (ls->setupcalled) PetscFunctionReturn(0); 1790c52818fSSatish Balay if (!((PetscObject)ls)->type_name) { 1800c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 1810c52818fSSatish Balay } 1820c52818fSSatish Balay if (ls->ops->setup) { 1830c52818fSSatish Balay ierr = (*ls->ops->setup)(ls);CHKERRQ(ierr); 1840c52818fSSatish Balay } 1850c52818fSSatish Balay if (ls->usetaoroutines) { 186441846f8SBarry Smith ierr = TaoIsObjectiveDefined(ls->tao,&flg);CHKERRQ(ierr); 1870c52818fSSatish Balay ls->hasobjective = flg; 188441846f8SBarry Smith ierr = TaoIsGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 1890c52818fSSatish Balay ls->hasgradient = flg; 190441846f8SBarry Smith ierr = TaoIsObjectiveAndGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 1910c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 1920c52818fSSatish Balay } else { 1930c52818fSSatish Balay if (ls->ops->computeobjective) { 1940c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 1950c52818fSSatish Balay } else { 1960c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 1970c52818fSSatish Balay } 1980c52818fSSatish Balay if (ls->ops->computegradient) { 1990c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 2000c52818fSSatish Balay } else { 2010c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 2020c52818fSSatish Balay } 2030c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 2040c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 2050c52818fSSatish Balay } else { 2060c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 2070c52818fSSatish Balay } 2080c52818fSSatish Balay } 2090c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 2100c52818fSSatish Balay PetscFunctionReturn(0); 2110c52818fSSatish Balay } 2120c52818fSSatish Balay 2130c52818fSSatish Balay /*@ 2140c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 2150c52818fSSatish Balay from one TaoLineSearchApply() to the next. This function resets this 2160c52818fSSatish Balay state information. 2170c52818fSSatish Balay 2180c52818fSSatish Balay Collective on TaoLineSearch 2190c52818fSSatish Balay 2200c52818fSSatish Balay Input Parameter: 2210c52818fSSatish Balay . ls - the TaoLineSearch context 2220c52818fSSatish Balay 2230c52818fSSatish Balay Level: developer 2240c52818fSSatish Balay 2250c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 2260c52818fSSatish Balay @*/ 2270c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 2280c52818fSSatish Balay { 2290c52818fSSatish Balay PetscErrorCode ierr; 230050fc7a3SBarry Smith 2310c52818fSSatish Balay PetscFunctionBegin; 2320c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 2330c52818fSSatish Balay if (ls->ops->reset) { 2340c52818fSSatish Balay ierr = (*ls->ops->reset)(ls);CHKERRQ(ierr); 2350c52818fSSatish Balay } 2360c52818fSSatish Balay PetscFunctionReturn(0); 2370c52818fSSatish Balay } 238f06e3bfaSBarry Smith 2390c52818fSSatish Balay /*@ 2400c52818fSSatish Balay TaoLineSearchDestroy - Destroys the TAO context that was created with 2410c52818fSSatish Balay TaoLineSearchCreate() 2420c52818fSSatish Balay 2430c52818fSSatish Balay Collective on TaoLineSearch 2440c52818fSSatish Balay 2450c52818fSSatish Balay Input Parameter 2460c52818fSSatish Balay . ls - the TaoLineSearch context 2470c52818fSSatish Balay 2480c52818fSSatish Balay Level: beginner 2490c52818fSSatish Balay 2500c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 2510c52818fSSatish Balay @*/ 2520c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 2530c52818fSSatish Balay { 2540c52818fSSatish Balay PetscErrorCode ierr; 255050fc7a3SBarry Smith 2560c52818fSSatish Balay PetscFunctionBegin; 2570c52818fSSatish Balay if (!*ls) PetscFunctionReturn(0); 2580c52818fSSatish Balay PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 2590c52818fSSatish Balay if (--((PetscObject)*ls)->refct > 0) {*ls=0; PetscFunctionReturn(0);} 260050fc7a3SBarry Smith ierr = VecDestroy(&(*ls)->stepdirection);CHKERRQ(ierr); 261050fc7a3SBarry Smith ierr = VecDestroy(&(*ls)->start_x);CHKERRQ(ierr); 2620c52818fSSatish Balay if ((*ls)->ops->destroy) { 2630c52818fSSatish Balay ierr = (*(*ls)->ops->destroy)(*ls);CHKERRQ(ierr); 2640c52818fSSatish Balay } 2652a0dac07SAlp Dener if ((*ls)->usemonitor) { 2662a0dac07SAlp Dener ierr = PetscViewerDestroy(&(*ls)->viewer);CHKERRQ(ierr); 2672a0dac07SAlp Dener } 2680c52818fSSatish Balay ierr = PetscHeaderDestroy(ls);CHKERRQ(ierr); 2690c52818fSSatish Balay PetscFunctionReturn(0); 2700c52818fSSatish Balay } 2710c52818fSSatish Balay 2720c52818fSSatish Balay /*@ 2730c52818fSSatish Balay TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 2740c52818fSSatish Balay 2750c52818fSSatish Balay Collective on TaoLineSearch 2760c52818fSSatish Balay 2770c52818fSSatish Balay Input Parameters: 278441846f8SBarry Smith + ls - the Tao context 2790c52818fSSatish Balay . x - The current solution (on output x contains the new solution determined by the line search) 2800c52818fSSatish Balay . f - objective function value at current solution (on output contains the objective function value at new solution) 2810c52818fSSatish Balay . g - gradient evaluated at x (on output contains the gradient at new solution) 2820c52818fSSatish Balay - s - search direction 2830c52818fSSatish Balay 2840c52818fSSatish Balay Output Parameters: 2850c52818fSSatish Balay + x - new solution 2860c52818fSSatish Balay . f - objective function value at x 2870c52818fSSatish Balay . g - gradient vector at x 2880c52818fSSatish Balay . steplength - scalar multiplier of s used ( x = x0 + steplength * x ) 2890c52818fSSatish Balay - reason - reason why the line-search stopped 2900c52818fSSatish Balay 2910c52818fSSatish Balay reason will be set to one of: 2920c52818fSSatish Balay 2930c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 2940c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 2950c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 2960c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 2970c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 2980c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 2990c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 3000c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 3010c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 3020c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search 3030c52818fSSatish Balay 3040c52818fSSatish Balay Note: 3050c52818fSSatish Balay The algorithm developer must set up the TaoLineSearch with calls to 306441846f8SBarry Smith TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines() 3070c52818fSSatish Balay 3080c52818fSSatish Balay Note: 3090c52818fSSatish Balay You may or may not need to follow this with a call to 3100c52818fSSatish Balay TaoAddLineSearchCounts(), depending on whether you want these 3110c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 3120c52818fSSatish Balay 3130c52818fSSatish Balay Level: beginner 3140c52818fSSatish Balay 3150c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts() 3160c52818fSSatish Balay @*/ 3170c52818fSSatish Balay 318e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 3190c52818fSSatish Balay { 3200c52818fSSatish Balay PetscErrorCode ierr; 3210c52818fSSatish Balay PetscInt low1,low2,low3,high1,high2,high3; 3220c52818fSSatish Balay 3230c52818fSSatish Balay PetscFunctionBegin; 3240c52818fSSatish Balay *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3250c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 3260c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3273f6ba705SLisandro Dalcin PetscValidRealPointer(f,3); 3280c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 3290c52818fSSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,5); 3300c52818fSSatish Balay PetscValidPointer(reason,7); 3310c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 3320c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,g,4); 3330c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,s,5); 3340c52818fSSatish Balay ierr = VecGetOwnershipRange(x, &low1, &high1);CHKERRQ(ierr); 3350c52818fSSatish Balay ierr = VecGetOwnershipRange(g, &low2, &high2);CHKERRQ(ierr); 3360c52818fSSatish Balay ierr = VecGetOwnershipRange(s, &low3, &high3);CHKERRQ(ierr); 337f06e3bfaSBarry Smith if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths"); 3380c52818fSSatish Balay 3390c52818fSSatish Balay ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr); 340050fc7a3SBarry Smith ierr = VecDestroy(&ls->stepdirection);CHKERRQ(ierr); 341050fc7a3SBarry Smith ls->stepdirection = s; 3420c52818fSSatish Balay 3430c52818fSSatish Balay ierr = TaoLineSearchSetUp(ls);CHKERRQ(ierr); 344f06e3bfaSBarry Smith if (!ls->ops->apply) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine"); 3450c52818fSSatish Balay ls->nfeval=0; 3460c52818fSSatish Balay ls->ngeval=0; 3470c52818fSSatish Balay ls->nfgeval=0; 3480c52818fSSatish Balay /* Check parameter values */ 3490c52818fSSatish Balay if (ls->ftol < 0.0) { 350f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);CHKERRQ(ierr); 3510c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3520c52818fSSatish Balay } 3530c52818fSSatish Balay if (ls->rtol < 0.0) { 354f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);CHKERRQ(ierr); 3550c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3560c52818fSSatish Balay } 3570c52818fSSatish Balay if (ls->gtol < 0.0) { 358f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);CHKERRQ(ierr); 3590c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3600c52818fSSatish Balay } 3610c52818fSSatish Balay if (ls->stepmin < 0.0) { 362f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);CHKERRQ(ierr); 3630c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3640c52818fSSatish Balay } 3650c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 366f06e3bfaSBarry Smith ierr = PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);CHKERRQ(ierr); 3670c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3680c52818fSSatish Balay } 3690c52818fSSatish Balay if (ls->max_funcs < 0) { 3700c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);CHKERRQ(ierr); 3710c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3720c52818fSSatish Balay } 3730c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 374f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);CHKERRQ(ierr); 3750c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_INFORNAN; 3760c52818fSSatish Balay } 3770c52818fSSatish Balay 378302440fdSBarry Smith ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 3790c52818fSSatish Balay ierr = VecDestroy(&ls->start_x);CHKERRQ(ierr); 3800c52818fSSatish Balay ls->start_x = x; 381050fc7a3SBarry Smith 3820ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0);CHKERRQ(ierr); 3830c52818fSSatish Balay ierr = (*ls->ops->apply)(ls,x,f,g,s);CHKERRQ(ierr); 3840ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0);CHKERRQ(ierr); 3850c52818fSSatish Balay *reason=ls->reason; 3860c52818fSSatish Balay ls->new_f = *f; 3870c52818fSSatish Balay 3880c52818fSSatish Balay if (steplength) { 3890c52818fSSatish Balay *steplength=ls->step; 3900c52818fSSatish Balay } 3910c52818fSSatish Balay 392fbe0838dSJason Sarich ierr = TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");CHKERRQ(ierr); 3930c52818fSSatish Balay PetscFunctionReturn(0); 3940c52818fSSatish Balay } 3950c52818fSSatish Balay 3960c52818fSSatish Balay /*@C 3970c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 3980c52818fSSatish Balay 3990c52818fSSatish Balay Collective on TaoLineSearch 4000c52818fSSatish Balay 4010c52818fSSatish Balay Input Parameters: 4020c52818fSSatish Balay + ls - the TaoLineSearch context 4030c52818fSSatish Balay - type - a known method 4040c52818fSSatish Balay 4050c52818fSSatish Balay Available methods include: 4060c52818fSSatish Balay + more-thuente 4070c52818fSSatish Balay . gpcg 4080c52818fSSatish Balay - unit - Do not perform any line search 4090c52818fSSatish Balay 4100c52818fSSatish Balay 4110c52818fSSatish Balay Options Database Keys: 4120c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 4130c52818fSSatish Balay 4140c52818fSSatish Balay Level: beginner 4150c52818fSSatish Balay 4160c52818fSSatish Balay 4170c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply() 4180c52818fSSatish Balay 4190c52818fSSatish Balay @*/ 4200c52818fSSatish Balay 421dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 4220c52818fSSatish Balay { 4230c52818fSSatish Balay PetscErrorCode ierr; 4240c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 4250c52818fSSatish Balay PetscBool flg; 4260c52818fSSatish Balay 4270c52818fSSatish Balay PetscFunctionBegin; 4280c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4290c52818fSSatish Balay PetscValidCharPointer(type,2); 4300c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg);CHKERRQ(ierr); 4310c52818fSSatish Balay if (flg) PetscFunctionReturn(0); 4320c52818fSSatish Balay 4330c52818fSSatish Balay ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);CHKERRQ(ierr); 4340c52818fSSatish Balay if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 4350c52818fSSatish Balay if (ls->ops->destroy) { 4360c52818fSSatish Balay ierr = (*(ls)->ops->destroy)(ls);CHKERRQ(ierr); 4370c52818fSSatish Balay } 4380c52818fSSatish Balay ls->max_funcs=30; 4390c52818fSSatish Balay ls->ftol = 0.0001; 4400c52818fSSatish Balay ls->gtol = 0.9; 4416f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 4426f4723b1SBarry Smith ls->rtol = 1.0e-5; 4436f4723b1SBarry Smith #else 4440c52818fSSatish Balay ls->rtol = 1.0e-10; 4456f4723b1SBarry Smith #endif 4460c52818fSSatish Balay ls->stepmin=1.0e-20; 4470c52818fSSatish Balay ls->stepmax=1.0e+20; 4480c52818fSSatish Balay 4490c52818fSSatish Balay ls->nfeval=0; 4500c52818fSSatish Balay ls->ngeval=0; 4510c52818fSSatish Balay ls->nfgeval=0; 4520c52818fSSatish Balay ls->ops->setup=0; 4530c52818fSSatish Balay ls->ops->apply=0; 4540c52818fSSatish Balay ls->ops->view=0; 4550c52818fSSatish Balay ls->ops->setfromoptions=0; 4560c52818fSSatish Balay ls->ops->destroy=0; 4570c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4580c52818fSSatish Balay ierr = (*r)(ls);CHKERRQ(ierr); 4590c52818fSSatish Balay ierr = PetscObjectChangeTypeName((PetscObject)ls, type);CHKERRQ(ierr); 4600c52818fSSatish Balay PetscFunctionReturn(0); 4610c52818fSSatish Balay } 4620c52818fSSatish Balay 4632a0dac07SAlp Dener /*@C 4642a0dac07SAlp Dener TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the 4652a0dac07SAlp Dener iteration number, step length, and function value before calling the implementation 4662a0dac07SAlp Dener specific monitor. 4672a0dac07SAlp Dener 4682a0dac07SAlp Dener Input Parameters: 4692a0dac07SAlp Dener + ls - the TaoLineSearch context 4702a0dac07SAlp Dener . its - the current iterate number (>=0) 4712a0dac07SAlp Dener . f - the current objective function value 4722a0dac07SAlp Dener - step - the step length 4732a0dac07SAlp Dener 4742a0dac07SAlp Dener Options Database Key: 4752a0dac07SAlp Dener . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output 4762a0dac07SAlp Dener 4772a0dac07SAlp Dener Level: developer 4782a0dac07SAlp Dener 4792a0dac07SAlp Dener @*/ 4802a0dac07SAlp Dener PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step) 4812a0dac07SAlp Dener { 4822a0dac07SAlp Dener PetscErrorCode ierr; 4832a0dac07SAlp Dener PetscInt tabs; 4842a0dac07SAlp Dener 4852a0dac07SAlp Dener PetscFunctionBegin; 4862a0dac07SAlp Dener PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4872a0dac07SAlp Dener if (ls->usemonitor) { 4882a0dac07SAlp Dener ierr = PetscViewerASCIIGetTab(ls->viewer, &tabs);CHKERRQ(ierr); 4892a0dac07SAlp Dener ierr = PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel);CHKERRQ(ierr); 4902a0dac07SAlp Dener ierr = PetscViewerASCIIPrintf(ls->viewer, "%3D LS", its);CHKERRQ(ierr); 4912a0dac07SAlp Dener ierr = PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f);CHKERRQ(ierr); 492ad362606SAlp Dener ierr = PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step);CHKERRQ(ierr); 4932a0dac07SAlp Dener if (ls->ops->monitor && its > 0) { 4942a0dac07SAlp Dener ierr = PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3);CHKERRQ(ierr); 4952a0dac07SAlp Dener ierr = (*ls->ops->monitor)(ls);CHKERRQ(ierr); 4962a0dac07SAlp Dener } 4972a0dac07SAlp Dener ierr = PetscViewerASCIISetTab(ls->viewer, tabs);CHKERRQ(ierr); 4982a0dac07SAlp Dener } 4992a0dac07SAlp Dener PetscFunctionReturn(0); 5002a0dac07SAlp Dener } 5012a0dac07SAlp Dener 5020c52818fSSatish Balay /*@ 5030c52818fSSatish Balay TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 5040c52818fSSatish Balay options. 5050c52818fSSatish Balay 5060c52818fSSatish Balay Collective on TaoLineSearch 5070c52818fSSatish Balay 5080c52818fSSatish Balay Input Paremeter: 5090c52818fSSatish Balay . ls - the TaoLineSearch context 5100c52818fSSatish Balay 5110c52818fSSatish Balay Options Database Keys: 5120c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 5130c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 5140c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 5150c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 5160c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 5170c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 5180c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 5190c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 5200c52818fSSatish Balay 5210c52818fSSatish Balay Level: beginner 5220c52818fSSatish Balay @*/ 5230c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 5240c52818fSSatish Balay { 5250c52818fSSatish Balay PetscErrorCode ierr; 5268caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 5272a0dac07SAlp Dener char type[256],monfilename[PETSC_MAX_PATH_LEN]; 5282a0dac07SAlp Dener PetscViewer monviewer; 5290c52818fSSatish Balay PetscBool flg; 530f06e3bfaSBarry Smith 5310c52818fSSatish Balay PetscFunctionBegin; 5320c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5330c52818fSSatish Balay ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr); 5340c52818fSSatish Balay if (((PetscObject)ls)->type_name) { 5350c52818fSSatish Balay default_type = ((PetscObject)ls)->type_name; 5360c52818fSSatish Balay } 5370c52818fSSatish Balay /* Check for type from options */ 5380c52818fSSatish Balay ierr = PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);CHKERRQ(ierr); 5390c52818fSSatish Balay if (flg) { 5400c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,type);CHKERRQ(ierr); 5410c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 542302440fdSBarry Smith ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 5430c52818fSSatish Balay } 5440c52818fSSatish Balay 54594ae4db5SBarry Smith ierr = PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);CHKERRQ(ierr); 54694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);CHKERRQ(ierr); 54794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);CHKERRQ(ierr); 54894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);CHKERRQ(ierr); 54994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);CHKERRQ(ierr); 55094ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);CHKERRQ(ierr); 5512a0dac07SAlp Dener ierr = PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 5522a0dac07SAlp Dener if (flg) { 5532a0dac07SAlp Dener ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer);CHKERRQ(ierr); 5542a0dac07SAlp Dener ls->viewer = monviewer; 5552a0dac07SAlp Dener ls->usemonitor = PETSC_TRUE; 5562a0dac07SAlp Dener } 5570c52818fSSatish Balay if (ls->ops->setfromoptions) { 5583a004c28SBarry Smith ierr = (*ls->ops->setfromoptions)(PetscOptionsObject,ls);CHKERRQ(ierr); 5590c52818fSSatish Balay } 5600c52818fSSatish Balay ierr = PetscOptionsEnd();CHKERRQ(ierr); 5610c52818fSSatish Balay PetscFunctionReturn(0); 5620c52818fSSatish Balay } 5630c52818fSSatish Balay 5640c52818fSSatish Balay /*@C 5650c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5660c52818fSSatish Balay 5670c52818fSSatish Balay Not Collective 5680c52818fSSatish Balay 5690c52818fSSatish Balay Input Parameter: 5700c52818fSSatish Balay . ls - the TaoLineSearch context 5710c52818fSSatish Balay 5720c52818fSSatish Balay Output Paramter: 5730c52818fSSatish Balay . type - the line search algorithm in effect 5740c52818fSSatish Balay 5750c52818fSSatish Balay Level: developer 5760c52818fSSatish Balay 5770c52818fSSatish Balay @*/ 578dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 5790c52818fSSatish Balay { 5800c52818fSSatish Balay PetscFunctionBegin; 5810c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5820c52818fSSatish Balay PetscValidPointer(type,2); 5830c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 5840c52818fSSatish Balay PetscFunctionReturn(0); 5850c52818fSSatish Balay } 5860c52818fSSatish Balay 5870c52818fSSatish Balay /*@ 5880c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 5890c52818fSSatish Balay routines used by the line search in last application (not cumulative). 5900c52818fSSatish Balay 5910c52818fSSatish Balay Not Collective 5920c52818fSSatish Balay 5930c52818fSSatish Balay Input Parameter: 5940c52818fSSatish Balay . ls - the TaoLineSearch context 5950c52818fSSatish Balay 5960c52818fSSatish Balay Output Parameters: 5970c52818fSSatish Balay + nfeval - number of function evaluations 5980c52818fSSatish Balay . ngeval - number of gradient evaluations 5990c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 6000c52818fSSatish Balay 6010c52818fSSatish Balay Level: intermediate 6020c52818fSSatish Balay 6030c52818fSSatish Balay Note: 604441846f8SBarry Smith If the line search is using the Tao objective and gradient 605441846f8SBarry Smith routines directly (see TaoLineSearchUseTaoRoutines()), then TAO 6060c52818fSSatish Balay is already counting the number of evaluations. 6070c52818fSSatish Balay 6080c52818fSSatish Balay @*/ 6090c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 6100c52818fSSatish Balay { 6110c52818fSSatish Balay PetscFunctionBegin; 6120c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6130c52818fSSatish Balay *nfeval = ls->nfeval; 6140c52818fSSatish Balay *ngeval = ls->ngeval; 6150c52818fSSatish Balay *nfgeval = ls->nfgeval; 6160c52818fSSatish Balay PetscFunctionReturn(0); 6170c52818fSSatish Balay } 6180c52818fSSatish Balay 6190c52818fSSatish Balay /*@ 620441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 621441846f8SBarry Smith Tao evaluation routines. 6220c52818fSSatish Balay 6230c52818fSSatish Balay Not Collective 6240c52818fSSatish Balay 6250c52818fSSatish Balay Input Parameter: 6260c52818fSSatish Balay . ls - the TaoLineSearch context 6270c52818fSSatish Balay 6280c52818fSSatish Balay Output Parameter: 629441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines, 6300c52818fSSatish Balay otherwise PETSC_FALSE 6310c52818fSSatish Balay 6320c52818fSSatish Balay Level: developer 6330c52818fSSatish Balay @*/ 634441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 6350c52818fSSatish Balay { 6360c52818fSSatish Balay PetscFunctionBegin; 6370c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 638f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 6390c52818fSSatish Balay PetscFunctionReturn(0); 6400c52818fSSatish Balay } 6410c52818fSSatish Balay 6420c52818fSSatish Balay /*@C 6430c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 6440c52818fSSatish Balay 6450c52818fSSatish Balay Logically Collective on TaoLineSearch 6460c52818fSSatish Balay 6470c52818fSSatish Balay Input Parameter: 6480c52818fSSatish Balay + ls - the TaoLineSearch context 6490c52818fSSatish Balay . func - the objective function evaluation routine 6500c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6510c52818fSSatish Balay 6520c52818fSSatish Balay Calling sequence of func: 6530c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 6540c52818fSSatish Balay 6550c52818fSSatish Balay + x - input vector 6560c52818fSSatish Balay . f - function value 6570c52818fSSatish Balay - ctx (optional) user-defined context 6580c52818fSSatish Balay 6590c52818fSSatish Balay Level: beginner 6600c52818fSSatish Balay 6610c52818fSSatish Balay Note: 6620c52818fSSatish Balay Use this routine only if you want the line search objective 663441846f8SBarry Smith evaluation routine to be different from the Tao's objective 6640c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6650c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6660c52818fSSatish Balay 6670c52818fSSatish Balay Note: 6680c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6690c52818fSSatish Balay line search, application programmers should be wary of overriding the 6700c52818fSSatish Balay default objective routine. 6710c52818fSSatish Balay 672441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 6730c52818fSSatish Balay @*/ 6740c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx) 6750c52818fSSatish Balay { 6760c52818fSSatish Balay PetscFunctionBegin; 6770c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6780c52818fSSatish Balay 6790c52818fSSatish Balay ls->ops->computeobjective=func; 6800c52818fSSatish Balay if (ctx) ls->userctx_func=ctx; 6810c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 6820c52818fSSatish Balay PetscFunctionReturn(0); 6830c52818fSSatish Balay } 6840c52818fSSatish Balay 6850c52818fSSatish Balay /*@C 6860c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 6870c52818fSSatish Balay 6880c52818fSSatish Balay Logically Collective on TaoLineSearch 6890c52818fSSatish Balay 6900c52818fSSatish Balay Input Parameter: 6910c52818fSSatish Balay + ls - the TaoLineSearch context 6920c52818fSSatish Balay . func - the gradient evaluation routine 6930c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6940c52818fSSatish Balay 6950c52818fSSatish Balay Calling sequence of func: 6960c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, Vec g, void *ctx); 6970c52818fSSatish Balay 6980c52818fSSatish Balay + x - input vector 6990c52818fSSatish Balay . g - gradient vector 7000c52818fSSatish Balay - ctx (optional) user-defined context 7010c52818fSSatish Balay 7020c52818fSSatish Balay Level: beginner 7030c52818fSSatish Balay 7040c52818fSSatish Balay Note: 7050c52818fSSatish Balay Use this routine only if you want the line search gradient 706441846f8SBarry Smith evaluation routine to be different from the Tao's gradient 7070c52818fSSatish Balay evaluation routine. If you use this routine you must also set 7080c52818fSSatish Balay the line search function and/or function/gradient routine. 7090c52818fSSatish Balay 7100c52818fSSatish Balay Note: 7110c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 7120c52818fSSatish Balay line search, application programmers should be wary of overriding the 7130c52818fSSatish Balay default gradient routine. 7140c52818fSSatish Balay 715441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 7160c52818fSSatish Balay @*/ 7170c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx) 7180c52818fSSatish Balay { 7190c52818fSSatish Balay PetscFunctionBegin; 7200c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7210c52818fSSatish Balay ls->ops->computegradient=func; 7220c52818fSSatish Balay if (ctx) ls->userctx_grad=ctx; 7230c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 7240c52818fSSatish Balay PetscFunctionReturn(0); 7250c52818fSSatish Balay } 7260c52818fSSatish Balay 7270c52818fSSatish Balay /*@C 7280c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 7290c52818fSSatish Balay 7300c52818fSSatish Balay Logically Collective on TaoLineSearch 7310c52818fSSatish Balay 7320c52818fSSatish Balay Input Parameter: 7330c52818fSSatish Balay + ls - the TaoLineSearch context 7340c52818fSSatish Balay . func - the objective and gradient evaluation routine 7350c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7360c52818fSSatish Balay 7370c52818fSSatish Balay Calling sequence of func: 7380c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 7390c52818fSSatish Balay 7400c52818fSSatish Balay + x - input vector 7410c52818fSSatish Balay . f - function value 7420c52818fSSatish Balay . g - gradient vector 7430c52818fSSatish Balay - ctx (optional) user-defined context 7440c52818fSSatish Balay 7450c52818fSSatish Balay Level: beginner 7460c52818fSSatish Balay 7470c52818fSSatish Balay Note: 7480c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 749441846f8SBarry Smith evaluation routines to be different from the Tao's objective 7500c52818fSSatish Balay and gradient evaluation routines. 7510c52818fSSatish Balay 7520c52818fSSatish Balay Note: 7530c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7540c52818fSSatish Balay line search, application programmers should be wary of overriding the 7550c52818fSSatish Balay default objective routine. 7560c52818fSSatish Balay 757441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines() 7580c52818fSSatish Balay @*/ 7590c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx) 7600c52818fSSatish Balay { 7610c52818fSSatish Balay PetscFunctionBegin; 7620c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7630c52818fSSatish Balay ls->ops->computeobjectiveandgradient=func; 7640c52818fSSatish Balay if (ctx) ls->userctx_funcgrad=ctx; 7650c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7660c52818fSSatish Balay PetscFunctionReturn(0); 7670c52818fSSatish Balay } 7680c52818fSSatish Balay 7690c52818fSSatish Balay /*@C 7700c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 7710c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 7720c52818fSSatish Balay Sometimes it is more efficient to compute the inner product of the gradient 7730c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 7740c52818fSSatish Balay the line search typically needs of the gradient. 7750c52818fSSatish Balay 7760c52818fSSatish Balay Logically Collective on TaoLineSearch 7770c52818fSSatish Balay 7780c52818fSSatish Balay Input Parameter: 7790c52818fSSatish Balay + ls - the TaoLineSearch context 7800c52818fSSatish Balay . func - the objective and gradient evaluation routine 7810c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7820c52818fSSatish Balay 7830c52818fSSatish Balay Calling sequence of func: 7840c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 7850c52818fSSatish Balay 7860c52818fSSatish Balay + x - input vector 7870c52818fSSatish Balay . s - step direction 7880c52818fSSatish Balay . f - function value 7890c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 7900c52818fSSatish Balay - ctx (optional) user-defined context 7910c52818fSSatish Balay 7920c52818fSSatish Balay Note: The gradient will still need to be computed at the end of the line 7930c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 7940c52818fSSatish Balay routine 7950c52818fSSatish Balay 7960c52818fSSatish Balay Note: Bounded line searches (those used in bounded optimization algorithms) 7970c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 7980c52818fSSatish Balay x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength() 7990c52818fSSatish Balay 8000c52818fSSatish Balay Level: advanced 8010c52818fSSatish Balay 8020c52818fSSatish Balay Note: 8030c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 8040c52818fSSatish Balay line search, application programmers should be wary of overriding the 8050c52818fSSatish Balay default objective routine. 8060c52818fSSatish Balay 807441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines() 8080c52818fSSatish Balay @*/ 8090c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx) 8100c52818fSSatish Balay { 8110c52818fSSatish Balay PetscFunctionBegin; 8120c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8130c52818fSSatish Balay ls->ops->computeobjectiveandgts=func; 8140c52818fSSatish Balay if (ctx) ls->userctx_funcgts=ctx; 8150c52818fSSatish Balay ls->usegts = PETSC_TRUE; 8160c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 8170c52818fSSatish Balay PetscFunctionReturn(0); 8180c52818fSSatish Balay } 8190c52818fSSatish Balay 8200c52818fSSatish Balay /*@C 821441846f8SBarry Smith TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the 822441846f8SBarry Smith objective and gradient evaluation routines from the given Tao object. 8230c52818fSSatish Balay 8240c52818fSSatish Balay Logically Collective on TaoLineSearch 8250c52818fSSatish Balay 8260c52818fSSatish Balay Input Parameter: 8270c52818fSSatish Balay + ls - the TaoLineSearch context 828441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines 8290c52818fSSatish Balay 8300c52818fSSatish Balay Level: developer 8310c52818fSSatish Balay 8320c52818fSSatish Balay .seealso: TaoLineSearchCreate() 8330c52818fSSatish Balay @*/ 834441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts) 8350c52818fSSatish Balay { 8360c52818fSSatish Balay PetscFunctionBegin; 8370c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 838441846f8SBarry Smith PetscValidHeaderSpecific(ts,TAO_CLASSID,1); 839441846f8SBarry Smith ls->tao = ts; 8400c52818fSSatish Balay ls->usetaoroutines=PETSC_TRUE; 8410c52818fSSatish Balay PetscFunctionReturn(0); 8420c52818fSSatish Balay } 8430c52818fSSatish Balay 8440c52818fSSatish Balay /*@ 8450c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 8460c52818fSSatish Balay 8470c52818fSSatish Balay Collective on TaoLineSearch 8480c52818fSSatish Balay 8490c52818fSSatish Balay Input Parameters: 8500c52818fSSatish Balay + ls - the TaoLineSearch context 8510c52818fSSatish Balay - x - input vector 8520c52818fSSatish Balay 8530c52818fSSatish Balay Output Parameter: 8540c52818fSSatish Balay . f - Objective value at X 8550c52818fSSatish Balay 85695452b02SPatrick Sanan Notes: 85795452b02SPatrick Sanan TaoLineSearchComputeObjective() is typically used within line searches 8580c52818fSSatish Balay so most users would not generally call this routine themselves. 8590c52818fSSatish Balay 8600c52818fSSatish Balay Level: developer 8610c52818fSSatish Balay 8620c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 8630c52818fSSatish Balay @*/ 8640c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 8650c52818fSSatish Balay { 8660c52818fSSatish Balay PetscErrorCode ierr; 8670c52818fSSatish Balay Vec gdummy; 8680c52818fSSatish Balay PetscReal gts; 869f06e3bfaSBarry Smith 8700c52818fSSatish Balay PetscFunctionBegin; 8710c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8720c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 8730c52818fSSatish Balay PetscValidPointer(f,3); 8740c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8750c52818fSSatish Balay if (ls->usetaoroutines) { 876441846f8SBarry Smith ierr = TaoComputeObjective(ls->tao,x,f);CHKERRQ(ierr); 8770c52818fSSatish Balay } else { 878f06e3bfaSBarry Smith if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient && !ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 8790ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 8800c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective routine"); 8810c52818fSSatish Balay if (ls->ops->computeobjective) { 8820c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 8830c52818fSSatish Balay } else if (ls->ops->computeobjectiveandgradient) { 8840c52818fSSatish Balay ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr); 8850c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr); 8860c52818fSSatish Balay ierr = VecDestroy(&gdummy);CHKERRQ(ierr); 8870c52818fSSatish Balay } else { 8880c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts);CHKERRQ(ierr); 8890c52818fSSatish Balay } 8900c52818fSSatish Balay PetscStackPop; 8910ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 8920c52818fSSatish Balay } 8930c52818fSSatish Balay ls->nfeval++; 8940c52818fSSatish Balay PetscFunctionReturn(0); 8950c52818fSSatish Balay } 8960c52818fSSatish Balay 8970c52818fSSatish Balay /*@ 8980c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 8990c52818fSSatish Balay 900441846f8SBarry Smith Collective on Tao 9010c52818fSSatish Balay 9020c52818fSSatish Balay Input Parameters: 9030c52818fSSatish Balay + ls - the TaoLineSearch context 9040c52818fSSatish Balay - x - input vector 9050c52818fSSatish Balay 9060c52818fSSatish Balay Output Parameter: 9070c52818fSSatish Balay + f - Objective value at X 9080c52818fSSatish Balay - g - Gradient vector at X 9090c52818fSSatish Balay 91095452b02SPatrick Sanan Notes: 91195452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 9120c52818fSSatish Balay so most users would not generally call this routine themselves. 9130c52818fSSatish Balay 9140c52818fSSatish Balay Level: developer 9150c52818fSSatish Balay 9160c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 9170c52818fSSatish Balay @*/ 9180c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 9190c52818fSSatish Balay { 9200c52818fSSatish Balay PetscErrorCode ierr; 921f06e3bfaSBarry Smith 9220c52818fSSatish Balay PetscFunctionBegin; 9230c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9240c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9250c52818fSSatish Balay PetscValidPointer(f,3); 9260c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 9270c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9280c52818fSSatish Balay PetscCheckSameComm(ls,1,g,4); 9290c52818fSSatish Balay if (ls->usetaoroutines) { 930441846f8SBarry Smith ierr = TaoComputeObjectiveAndGradient(ls->tao,x,f,g);CHKERRQ(ierr); 9310c52818fSSatish Balay } else { 932f06e3bfaSBarry Smith if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 933f06e3bfaSBarry Smith if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 9340ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9350c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gradient routine"); 9360c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 9370c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr); 9380c52818fSSatish Balay } else { 9390c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 9400c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 9410c52818fSSatish Balay } 9420c52818fSSatish Balay PetscStackPop; 9430ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 944335036cbSBarry Smith ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 945f06e3bfaSBarry Smith } 946fbe0838dSJason Sarich ls->nfgeval++; 9470c52818fSSatish Balay PetscFunctionReturn(0); 9480c52818fSSatish Balay } 9490c52818fSSatish Balay 9500c52818fSSatish Balay /*@ 9510c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 9520c52818fSSatish Balay 9530c52818fSSatish Balay Collective on TaoLineSearch 9540c52818fSSatish Balay 9550c52818fSSatish Balay Input Parameters: 9560c52818fSSatish Balay + ls - the TaoLineSearch context 9570c52818fSSatish Balay - x - input vector 9580c52818fSSatish Balay 9590c52818fSSatish Balay Output Parameter: 9600c52818fSSatish Balay . g - gradient vector 9610c52818fSSatish Balay 96295452b02SPatrick Sanan Notes: 96395452b02SPatrick Sanan TaoComputeGradient() is typically used within line searches 9640c52818fSSatish Balay so most users would not generally call this routine themselves. 9650c52818fSSatish Balay 9660c52818fSSatish Balay Level: developer 9670c52818fSSatish Balay 9680c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 9690c52818fSSatish Balay @*/ 9700c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 9710c52818fSSatish Balay { 9720c52818fSSatish Balay PetscErrorCode ierr; 9730c52818fSSatish Balay PetscReal fdummy; 974f06e3bfaSBarry Smith 9750c52818fSSatish Balay PetscFunctionBegin; 9760c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9770c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9780c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,3); 9790c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9800c52818fSSatish Balay PetscCheckSameComm(ls,1,g,3); 9810c52818fSSatish Balay if (ls->usetaoroutines) { 982441846f8SBarry Smith ierr = TaoComputeGradient(ls->tao,x,g);CHKERRQ(ierr); 9830c52818fSSatish Balay } else { 984f06e3bfaSBarry Smith if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 9850ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9860c52818fSSatish Balay PetscStackPush("TaoLineSearch user gradient routine"); 9870c52818fSSatish Balay if (ls->ops->computegradient) { 9880c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 9890c52818fSSatish Balay } else { 9900c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr); 9910c52818fSSatish Balay } 9920c52818fSSatish Balay PetscStackPop; 9930ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9940c52818fSSatish Balay } 9950c52818fSSatish Balay ls->ngeval++; 9960c52818fSSatish Balay PetscFunctionReturn(0); 9970c52818fSSatish Balay } 9980c52818fSSatish Balay 9990c52818fSSatish Balay /*@ 10000c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 10010c52818fSSatish Balay 1002441846f8SBarry Smith Collective on Tao 10030c52818fSSatish Balay 10040c52818fSSatish Balay Input Parameters: 10050c52818fSSatish Balay + ls - the TaoLineSearch context 10060c52818fSSatish Balay - x - input vector 10070c52818fSSatish Balay 10080c52818fSSatish Balay Output Parameter: 10090c52818fSSatish Balay + f - Objective value at X 10100c52818fSSatish Balay - gts - inner product of gradient and step direction at X 10110c52818fSSatish Balay 101295452b02SPatrick Sanan Notes: 101395452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 10140c52818fSSatish Balay so most users would not generally call this routine themselves. 10150c52818fSSatish Balay 10160c52818fSSatish Balay Level: developer 10170c52818fSSatish Balay 10180c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 10190c52818fSSatish Balay @*/ 10200c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 10210c52818fSSatish Balay { 10220c52818fSSatish Balay PetscErrorCode ierr; 10230c52818fSSatish Balay PetscFunctionBegin; 10240c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10250c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 10260c52818fSSatish Balay PetscValidPointer(f,3); 10270c52818fSSatish Balay PetscValidPointer(gts,4); 10280c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 1029f06e3bfaSBarry Smith if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 10300ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 10310c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gts routine"); 10320c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr); 10330c52818fSSatish Balay PetscStackPop; 10340ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 1035335036cbSBarry Smith ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 10360c52818fSSatish Balay ls->nfeval++; 10370c52818fSSatish Balay PetscFunctionReturn(0); 10380c52818fSSatish Balay } 10390c52818fSSatish Balay 10400c52818fSSatish Balay /*@ 10410c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 10420c52818fSSatish Balay 10430c52818fSSatish Balay Collective on TaoLineSearch 10440c52818fSSatish Balay 10450c52818fSSatish Balay Input Parameter: 10460c52818fSSatish Balay . ls - the TaoLineSearch context 10470c52818fSSatish Balay 10480c52818fSSatish Balay Output Parameter: 10490c52818fSSatish Balay + x - the new solution 10500c52818fSSatish Balay . f - the objective function value at x 10510c52818fSSatish Balay . g - the gradient at x 10520c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 10530c52818fSSatish Balay - reason - the reason why the line search terminated 10540c52818fSSatish Balay 10550c52818fSSatish Balay reason will be set to one of: 10560c52818fSSatish Balay 10570c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 10580c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 10590c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 10600c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 10610c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 10620c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 10630c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 10640c52818fSSatish Balay 10650c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 10660c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 10670c52818fSSatish Balay 1068*a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search 10690c52818fSSatish Balay 10700c52818fSSatish Balay Level: developer 10710c52818fSSatish Balay 10720c52818fSSatish Balay @*/ 1073e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 10740c52818fSSatish Balay { 10750c52818fSSatish Balay PetscErrorCode ierr; 1076f06e3bfaSBarry Smith 10770c52818fSSatish Balay PetscFunctionBegin; 10780c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10790c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 10800c52818fSSatish Balay PetscValidPointer(f,3); 10810c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 10820c52818fSSatish Balay PetscValidIntPointer(reason,6); 10830c52818fSSatish Balay if (ls->new_x) { 10840c52818fSSatish Balay ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr); 10850c52818fSSatish Balay } 10860c52818fSSatish Balay *f = ls->new_f; 10870c52818fSSatish Balay if (ls->new_g) { 10880c52818fSSatish Balay ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr); 10890c52818fSSatish Balay } 10900c52818fSSatish Balay if (steplength) { 10910c52818fSSatish Balay *steplength=ls->step; 10920c52818fSSatish Balay } 10930c52818fSSatish Balay *reason = ls->reason; 10940c52818fSSatish Balay PetscFunctionReturn(0); 10950c52818fSSatish Balay } 10960c52818fSSatish Balay 10970c52818fSSatish Balay /*@ 10980c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 10990c52818fSSatish Balay search. 11000c52818fSSatish Balay 11010c52818fSSatish Balay Not Collective 11020c52818fSSatish Balay 11030c52818fSSatish Balay Input Parameter: 11040c52818fSSatish Balay . ls - the TaoLineSearch context 11050c52818fSSatish Balay 11060c52818fSSatish Balay Output Parameter: 11070c52818fSSatish Balay . x - The initial point of the line search 11080c52818fSSatish Balay 11090c52818fSSatish Balay Level: intermediate 11100c52818fSSatish Balay @*/ 11110c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 11120c52818fSSatish Balay { 11130c52818fSSatish Balay PetscFunctionBegin; 11140c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11150c52818fSSatish Balay if (x) { 11160c52818fSSatish Balay *x = ls->start_x; 11170c52818fSSatish Balay } 11180c52818fSSatish Balay PetscFunctionReturn(0); 11190c52818fSSatish Balay } 11200c52818fSSatish Balay 11210c52818fSSatish Balay /*@ 11220c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 11230c52818fSSatish Balay search. 11240c52818fSSatish Balay 11250c52818fSSatish Balay Not Collective 11260c52818fSSatish Balay 11270c52818fSSatish Balay Input Parameter: 11280c52818fSSatish Balay . ls - the TaoLineSearch context 11290c52818fSSatish Balay 11300c52818fSSatish Balay Output Parameter: 11310c52818fSSatish Balay . s - the step direction of the line search 11320c52818fSSatish Balay 11330c52818fSSatish Balay Level: advanced 11340c52818fSSatish Balay @*/ 11350c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 11360c52818fSSatish Balay { 11370c52818fSSatish Balay PetscFunctionBegin; 11380c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11390c52818fSSatish Balay if (s) { 11400c52818fSSatish Balay *s = ls->stepdirection; 11410c52818fSSatish Balay } 11420c52818fSSatish Balay PetscFunctionReturn(0); 11430c52818fSSatish Balay 11440c52818fSSatish Balay } 11450c52818fSSatish Balay 11460c52818fSSatish Balay /*@ 11470c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 11480c52818fSSatish Balay 11490c52818fSSatish Balay Not Collective 11500c52818fSSatish Balay 11510c52818fSSatish Balay Input Parameter: 11520c52818fSSatish Balay . ls - the TaoLineSearch context 11530c52818fSSatish Balay 11540c52818fSSatish Balay Output Parameter: 11550c52818fSSatish Balay . f - the objective value at the full step length 11560c52818fSSatish Balay 11570c52818fSSatish Balay Level: developer 11580c52818fSSatish Balay @*/ 11590c52818fSSatish Balay 11600c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 11610c52818fSSatish Balay { 11620c52818fSSatish Balay PetscFunctionBegin; 11630c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11640c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 11650c52818fSSatish Balay PetscFunctionReturn(0); 11660c52818fSSatish Balay } 11670c52818fSSatish Balay 11680c52818fSSatish Balay /*@ 11690c52818fSSatish Balay TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 11700c52818fSSatish Balay 1171441846f8SBarry Smith Logically Collective on Tao 11720c52818fSSatish Balay 11730c52818fSSatish Balay Input Parameters: 11740c52818fSSatish Balay + ls - the TaoLineSearch context 11750c52818fSSatish Balay . xl - vector of lower bounds 11760c52818fSSatish Balay - xu - vector of upper bounds 11770c52818fSSatish Balay 11780c52818fSSatish Balay Note: If the variable bounds are not set with this routine, then 1179e270355aSBarry Smith PETSC_NINFINITY and PETSC_INFINITY are assumed 11800c52818fSSatish Balay 11810c52818fSSatish Balay Level: beginner 11820c52818fSSatish Balay 11830c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 11840c52818fSSatish Balay @*/ 11850c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 11860c52818fSSatish Balay { 11870c52818fSSatish Balay PetscFunctionBegin; 11880c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11890c52818fSSatish Balay PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 11900c52818fSSatish Balay PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 11910c52818fSSatish Balay ls->lower = xl; 11920c52818fSSatish Balay ls->upper = xu; 11930c52818fSSatish Balay ls->bounded = 1; 11940c52818fSSatish Balay PetscFunctionReturn(0); 11950c52818fSSatish Balay } 11960c52818fSSatish Balay 11970c52818fSSatish Balay /*@ 11980c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 11990c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 12000c52818fSSatish Balay 12010c52818fSSatish Balay Logically Collective on TaoLineSearch 12020c52818fSSatish Balay 12030c52818fSSatish Balay Input Parameters: 12040c52818fSSatish Balay + ls - the TaoLineSearch context 12050c52818fSSatish Balay - s - the initial step size 12060c52818fSSatish Balay 12070c52818fSSatish Balay Level: intermediate 12080c52818fSSatish Balay 12090c52818fSSatish Balay .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply() 12100c52818fSSatish Balay @*/ 12110c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 12120c52818fSSatish Balay { 12130c52818fSSatish Balay PetscFunctionBegin; 12140c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 12150c52818fSSatish Balay ls->initstep = s; 12160c52818fSSatish Balay PetscFunctionReturn(0); 12170c52818fSSatish Balay } 12180c52818fSSatish Balay 12190c52818fSSatish Balay /*@ 12200c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 12210c52818fSSatish Balay 12220c52818fSSatish Balay Not Collective 12230c52818fSSatish Balay 12240c52818fSSatish Balay Input Parameters: 12250c52818fSSatish Balay . ls - the TaoLineSearch context 12260c52818fSSatish Balay 12270c52818fSSatish Balay Output Parameters 12280c52818fSSatish Balay . s - the current step length 12290c52818fSSatish Balay 12300c52818fSSatish Balay Level: beginner 12310c52818fSSatish Balay 12320c52818fSSatish Balay .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply() 12330c52818fSSatish Balay @*/ 12340c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 12350c52818fSSatish Balay { 12360c52818fSSatish Balay PetscFunctionBegin; 12370c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 12380c52818fSSatish Balay *s = ls->step; 12390c52818fSSatish Balay PetscFunctionReturn(0); 12400c52818fSSatish Balay } 12410c52818fSSatish Balay 12420ebee16dSLisandro Dalcin /*@C 12430c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 12440c52818fSSatish Balay 12450c52818fSSatish Balay Not collective 12460c52818fSSatish Balay 12470c52818fSSatish Balay Input Parameters: 12480c52818fSSatish Balay + sname - name of a new user-defined solver 12490c52818fSSatish Balay - func - routine to Create method context 12500c52818fSSatish Balay 12510c52818fSSatish Balay Notes: 12520c52818fSSatish Balay TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 12530c52818fSSatish Balay 12540c52818fSSatish Balay Sample usage: 12550c52818fSSatish Balay .vb 12560c52818fSSatish Balay TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 12570c52818fSSatish Balay .ve 12580c52818fSSatish Balay 12590c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 12600c52818fSSatish Balay $ TaoLineSearchSetType(ls,"my_linesearch") 12610c52818fSSatish Balay or at runtime via the option 12620c52818fSSatish Balay $ -tao_ls_type my_linesearch 12630c52818fSSatish Balay 12640c52818fSSatish Balay Level: developer 12650c52818fSSatish Balay 12660ebee16dSLisandro Dalcin @*/ 12670c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 12680c52818fSSatish Balay { 12690c52818fSSatish Balay PetscErrorCode ierr; 12700c52818fSSatish Balay PetscFunctionBegin; 12711d36bdfdSBarry Smith ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr); 12720c52818fSSatish Balay ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr); 12730c52818fSSatish Balay PetscFunctionReturn(0); 12740c52818fSSatish Balay } 12750c52818fSSatish Balay 12760c52818fSSatish Balay /*@C 12770c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 12780c52818fSSatish Balay for all TaoLineSearch options in the database. 12790c52818fSSatish Balay 12800c52818fSSatish Balay 12810c52818fSSatish Balay Collective on TaoLineSearch 12820c52818fSSatish Balay 12830c52818fSSatish Balay Input Parameters: 12840c52818fSSatish Balay + ls - the TaoLineSearch solver context 12850c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 12860c52818fSSatish Balay 12870c52818fSSatish Balay Notes: 12880c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12890c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12900c52818fSSatish Balay 12910c52818fSSatish Balay 12920c52818fSSatish Balay Level: advanced 12930c52818fSSatish Balay 12940c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 12950c52818fSSatish Balay @*/ 12960c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 12970c52818fSSatish Balay { 1298f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 12990c52818fSSatish Balay } 13000c52818fSSatish Balay 13010c52818fSSatish Balay /*@C 13020c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 13030c52818fSSatish Balay TaoLineSearch options in the database 13040c52818fSSatish Balay 13050c52818fSSatish Balay Not Collective 13060c52818fSSatish Balay 13070c52818fSSatish Balay Input Parameters: 13080c52818fSSatish Balay . ls - the TaoLineSearch context 13090c52818fSSatish Balay 13100c52818fSSatish Balay Output Parameters: 13110c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 13120c52818fSSatish Balay 131395452b02SPatrick Sanan Notes: 131495452b02SPatrick Sanan On the fortran side, the user should pass in a string 'prefix' of 13150c52818fSSatish Balay sufficient length to hold the prefix. 13160c52818fSSatish Balay 13170c52818fSSatish Balay Level: advanced 13180c52818fSSatish Balay 13190c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 13200c52818fSSatish Balay @*/ 13210c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 13220c52818fSSatish Balay { 1323f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 13240c52818fSSatish Balay } 13250c52818fSSatish Balay 13260c52818fSSatish Balay /*@C 13270c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 13280c52818fSSatish Balay TaoLineSearch options in the database. 13290c52818fSSatish Balay 13300c52818fSSatish Balay 13310c52818fSSatish Balay Logically Collective on TaoLineSearch 13320c52818fSSatish Balay 13330c52818fSSatish Balay Input Parameters: 13340c52818fSSatish Balay + ls - the TaoLineSearch context 13350c52818fSSatish Balay - prefix - the prefix string to prepend to all TAO option requests 13360c52818fSSatish Balay 13370c52818fSSatish Balay Notes: 13380c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 13390c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 13400c52818fSSatish Balay 13410c52818fSSatish Balay For example, to distinguish between the runtime options for two 13420c52818fSSatish Balay different line searches, one could call 13430c52818fSSatish Balay .vb 13440c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 13450c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 13460c52818fSSatish Balay .ve 13470c52818fSSatish Balay 13480c52818fSSatish Balay This would enable use of different options for each system, such as 13490c52818fSSatish Balay .vb 13500c52818fSSatish Balay -sys1_tao_ls_type mt 13510c52818fSSatish Balay -sys2_tao_ls_type armijo 13520c52818fSSatish Balay .ve 13530c52818fSSatish Balay 13540c52818fSSatish Balay Level: advanced 13550c52818fSSatish Balay 13560c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 13570c52818fSSatish Balay @*/ 13580c52818fSSatish Balay 13590c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 13600c52818fSSatish Balay { 1361f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 13620c52818fSSatish Balay } 1363