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; 7*0ebee16dSLisandro Dalcin 8*0ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Apply; 9*0ebee16dSLisandro 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 PetscInt tabs; 5463b15415SAlp Dener ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr); 5563b15415SAlp Dener ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)ls)->tablevel);CHKERRQ(ierr); 5663b15415SAlp Dener ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);CHKERRQ(ierr); 5763b15415SAlp Dener if (ls->ops->view) { 5863b15415SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 5963b15415SAlp Dener ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr); 6063b15415SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 610c52818fSSatish Balay } 620c52818fSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 630c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr); 64f06e3bfaSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);CHKERRQ(ierr); 650c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr); 660c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr); 670c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr); 680c52818fSSatish Balay 690c52818fSSatish Balay if (ls->bounded) { 700c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr); 710c52818fSSatish Balay } 72de196d00SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);CHKERRQ(ierr); 730c52818fSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7463b15415SAlp Dener ierr = PetscViewerASCIISetTab(viewer, tabs);CHKERRQ(ierr); 750c52818fSSatish Balay } else if (isstring) { 760c52818fSSatish Balay ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr); 770c52818fSSatish Balay ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 780c52818fSSatish Balay } 790c52818fSSatish Balay PetscFunctionReturn(0); 800c52818fSSatish Balay } 810c52818fSSatish Balay 820c52818fSSatish Balay /*@C 830c52818fSSatish Balay TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 840c52818fSSatish Balay line-searches will automatically create one. 850c52818fSSatish Balay 860c52818fSSatish Balay Collective on MPI_Comm 870c52818fSSatish Balay 880c52818fSSatish Balay Input Parameter: 890c52818fSSatish Balay . comm - MPI communicator 900c52818fSSatish Balay 910c52818fSSatish Balay Output Parameter: 920c52818fSSatish Balay . newls - the new TaoLineSearch context 930c52818fSSatish Balay 940c52818fSSatish Balay Available methods include: 950c52818fSSatish Balay + more-thuente 960c52818fSSatish Balay . gpcg 970c52818fSSatish Balay - unit - Do not perform any line search 980c52818fSSatish Balay 990c52818fSSatish Balay 1000c52818fSSatish Balay Options Database Keys: 1010c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 1020c52818fSSatish Balay 1030c52818fSSatish Balay Level: beginner 1040c52818fSSatish Balay 1050c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy() 1060c52818fSSatish Balay @*/ 1070c52818fSSatish Balay 1080c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 1090c52818fSSatish Balay { 1100c52818fSSatish Balay PetscErrorCode ierr; 1110c52818fSSatish Balay TaoLineSearch ls; 1120c52818fSSatish Balay 1130c52818fSSatish Balay PetscFunctionBegin; 1140c52818fSSatish Balay PetscValidPointer(newls,2); 1156c23d075SBarry Smith *newls = NULL; 1160c52818fSSatish Balay 1170c52818fSSatish Balay ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr); 1180c52818fSSatish Balay 11973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);CHKERRQ(ierr); 1200c52818fSSatish Balay ls->bounded = 0; 1210c52818fSSatish Balay ls->max_funcs=30; 1220c52818fSSatish Balay ls->ftol = 0.0001; 1230c52818fSSatish Balay ls->gtol = 0.9; 1246f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1256f4723b1SBarry Smith ls->rtol = 1.0e-5; 1266f4723b1SBarry Smith #else 1270c52818fSSatish Balay ls->rtol = 1.0e-10; 1286f4723b1SBarry Smith #endif 1290c52818fSSatish Balay ls->stepmin=1.0e-20; 1300c52818fSSatish Balay ls->stepmax=1.0e+20; 1310c52818fSSatish Balay ls->step=1.0; 1320c52818fSSatish Balay ls->nfeval=0; 1330c52818fSSatish Balay ls->ngeval=0; 1340c52818fSSatish Balay ls->nfgeval=0; 1350c52818fSSatish Balay 1360c52818fSSatish Balay ls->ops->computeobjective=0; 1370c52818fSSatish Balay ls->ops->computegradient=0; 1380c52818fSSatish Balay ls->ops->computeobjectiveandgradient=0; 1390c52818fSSatish Balay ls->ops->computeobjectiveandgts=0; 1400c52818fSSatish Balay ls->ops->setup=0; 1410c52818fSSatish Balay ls->ops->apply=0; 1420c52818fSSatish Balay ls->ops->view=0; 1430c52818fSSatish Balay ls->ops->setfromoptions=0; 1440c52818fSSatish Balay ls->ops->reset=0; 1450c52818fSSatish Balay ls->ops->destroy=0; 1460c52818fSSatish Balay ls->setupcalled=PETSC_FALSE; 1470c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 1480c52818fSSatish Balay *newls = ls; 1490c52818fSSatish Balay PetscFunctionReturn(0); 1500c52818fSSatish Balay } 1510c52818fSSatish Balay 1520c52818fSSatish Balay /*@ 1530c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 1540c52818fSSatish Balay of a Tao solver 1550c52818fSSatish Balay 1560c52818fSSatish Balay Collective on ls 1570c52818fSSatish Balay 1580c52818fSSatish Balay Input Parameters: 1590c52818fSSatish Balay . ls - the TaoLineSearch context 1600c52818fSSatish Balay 1610c52818fSSatish Balay Notes: 1620c52818fSSatish Balay The user will not need to explicitly call TaoLineSearchSetUp(), as it will 1630c52818fSSatish Balay automatically be called in TaoLineSearchSolve(). However, if the user 1640c52818fSSatish Balay desires to call it explicitly, it should come after TaoLineSearchCreate() 1650c52818fSSatish Balay but before TaoLineSearchApply(). 1660c52818fSSatish Balay 1670c52818fSSatish Balay Level: developer 1680c52818fSSatish Balay 1690c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 1700c52818fSSatish Balay @*/ 1710c52818fSSatish Balay 1720c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 1730c52818fSSatish Balay { 1740c52818fSSatish Balay PetscErrorCode ierr; 1758caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 1760c52818fSSatish Balay PetscBool flg; 177f06e3bfaSBarry Smith 1780c52818fSSatish Balay PetscFunctionBegin; 1790c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1800c52818fSSatish Balay if (ls->setupcalled) PetscFunctionReturn(0); 1810c52818fSSatish Balay if (!((PetscObject)ls)->type_name) { 1820c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 1830c52818fSSatish Balay } 1840c52818fSSatish Balay if (ls->ops->setup) { 1850c52818fSSatish Balay ierr = (*ls->ops->setup)(ls);CHKERRQ(ierr); 1860c52818fSSatish Balay } 1870c52818fSSatish Balay if (ls->usetaoroutines) { 188441846f8SBarry Smith ierr = TaoIsObjectiveDefined(ls->tao,&flg);CHKERRQ(ierr); 1890c52818fSSatish Balay ls->hasobjective = flg; 190441846f8SBarry Smith ierr = TaoIsGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 1910c52818fSSatish Balay ls->hasgradient = flg; 192441846f8SBarry Smith ierr = TaoIsObjectiveAndGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 1930c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 1940c52818fSSatish Balay } else { 1950c52818fSSatish Balay if (ls->ops->computeobjective) { 1960c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 1970c52818fSSatish Balay } else { 1980c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 1990c52818fSSatish Balay } 2000c52818fSSatish Balay if (ls->ops->computegradient) { 2010c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 2020c52818fSSatish Balay } else { 2030c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 2040c52818fSSatish Balay } 2050c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 2060c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 2070c52818fSSatish Balay } else { 2080c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 2090c52818fSSatish Balay } 2100c52818fSSatish Balay } 2110c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 2120c52818fSSatish Balay PetscFunctionReturn(0); 2130c52818fSSatish Balay } 2140c52818fSSatish Balay 2150c52818fSSatish Balay /*@ 2160c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 2170c52818fSSatish Balay from one TaoLineSearchApply() to the next. This function resets this 2180c52818fSSatish Balay state information. 2190c52818fSSatish Balay 2200c52818fSSatish Balay Collective on TaoLineSearch 2210c52818fSSatish Balay 2220c52818fSSatish Balay Input Parameter: 2230c52818fSSatish Balay . ls - the TaoLineSearch context 2240c52818fSSatish Balay 2250c52818fSSatish Balay Level: developer 2260c52818fSSatish Balay 2270c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 2280c52818fSSatish Balay @*/ 2290c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 2300c52818fSSatish Balay { 2310c52818fSSatish Balay PetscErrorCode ierr; 232050fc7a3SBarry Smith 2330c52818fSSatish Balay PetscFunctionBegin; 2340c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 2350c52818fSSatish Balay if (ls->ops->reset) { 2360c52818fSSatish Balay ierr = (*ls->ops->reset)(ls);CHKERRQ(ierr); 2370c52818fSSatish Balay } 2380c52818fSSatish Balay PetscFunctionReturn(0); 2390c52818fSSatish Balay } 240f06e3bfaSBarry Smith 2410c52818fSSatish Balay /*@ 2420c52818fSSatish Balay TaoLineSearchDestroy - Destroys the TAO context that was created with 2430c52818fSSatish Balay TaoLineSearchCreate() 2440c52818fSSatish Balay 2450c52818fSSatish Balay Collective on TaoLineSearch 2460c52818fSSatish Balay 2470c52818fSSatish Balay Input Parameter 2480c52818fSSatish Balay . ls - the TaoLineSearch context 2490c52818fSSatish Balay 2500c52818fSSatish Balay Level: beginner 2510c52818fSSatish Balay 2520c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 2530c52818fSSatish Balay @*/ 2540c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 2550c52818fSSatish Balay { 2560c52818fSSatish Balay PetscErrorCode ierr; 257050fc7a3SBarry Smith 2580c52818fSSatish Balay PetscFunctionBegin; 2590c52818fSSatish Balay if (!*ls) PetscFunctionReturn(0); 2600c52818fSSatish Balay PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 2610c52818fSSatish Balay if (--((PetscObject)*ls)->refct > 0) {*ls=0; PetscFunctionReturn(0);} 262050fc7a3SBarry Smith ierr = VecDestroy(&(*ls)->stepdirection);CHKERRQ(ierr); 263050fc7a3SBarry Smith ierr = VecDestroy(&(*ls)->start_x);CHKERRQ(ierr); 2640c52818fSSatish Balay if ((*ls)->ops->destroy) { 2650c52818fSSatish Balay ierr = (*(*ls)->ops->destroy)(*ls);CHKERRQ(ierr); 2660c52818fSSatish Balay } 2670c52818fSSatish Balay ierr = PetscHeaderDestroy(ls);CHKERRQ(ierr); 2680c52818fSSatish Balay PetscFunctionReturn(0); 2690c52818fSSatish Balay } 2700c52818fSSatish Balay 2710c52818fSSatish Balay /*@ 2720c52818fSSatish Balay TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 2730c52818fSSatish Balay 2740c52818fSSatish Balay Collective on TaoLineSearch 2750c52818fSSatish Balay 2760c52818fSSatish Balay Input Parameters: 277441846f8SBarry Smith + ls - the Tao context 2780c52818fSSatish Balay . x - The current solution (on output x contains the new solution determined by the line search) 2790c52818fSSatish Balay . f - objective function value at current solution (on output contains the objective function value at new solution) 2800c52818fSSatish Balay . g - gradient evaluated at x (on output contains the gradient at new solution) 2810c52818fSSatish Balay - s - search direction 2820c52818fSSatish Balay 2830c52818fSSatish Balay Output Parameters: 2840c52818fSSatish Balay + x - new solution 2850c52818fSSatish Balay . f - objective function value at x 2860c52818fSSatish Balay . g - gradient vector at x 2870c52818fSSatish Balay . steplength - scalar multiplier of s used ( x = x0 + steplength * x ) 2880c52818fSSatish Balay - reason - reason why the line-search stopped 2890c52818fSSatish Balay 2900c52818fSSatish Balay reason will be set to one of: 2910c52818fSSatish Balay 2920c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 2930c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 2940c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 2950c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 2960c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 2970c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 2980c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 2990c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 3000c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 3010c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search 3020c52818fSSatish Balay 3030c52818fSSatish Balay Note: 3040c52818fSSatish Balay The algorithm developer must set up the TaoLineSearch with calls to 305441846f8SBarry Smith TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines() 3060c52818fSSatish Balay 3070c52818fSSatish Balay Note: 3080c52818fSSatish Balay You may or may not need to follow this with a call to 3090c52818fSSatish Balay TaoAddLineSearchCounts(), depending on whether you want these 3100c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 3110c52818fSSatish Balay 3120c52818fSSatish Balay Level: beginner 3130c52818fSSatish Balay 3140c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts() 3150c52818fSSatish Balay @*/ 3160c52818fSSatish Balay 317e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 3180c52818fSSatish Balay { 3190c52818fSSatish Balay PetscErrorCode ierr; 3200c52818fSSatish Balay PetscInt low1,low2,low3,high1,high2,high3; 3210c52818fSSatish Balay 3220c52818fSSatish Balay PetscFunctionBegin; 3230c52818fSSatish Balay *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3240c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 3250c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3263f6ba705SLisandro Dalcin PetscValidRealPointer(f,3); 3270c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 3280c52818fSSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,5); 3290c52818fSSatish Balay PetscValidPointer(reason,7); 3300c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 3310c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,g,4); 3320c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,s,5); 3330c52818fSSatish Balay ierr = VecGetOwnershipRange(x, &low1, &high1);CHKERRQ(ierr); 3340c52818fSSatish Balay ierr = VecGetOwnershipRange(g, &low2, &high2);CHKERRQ(ierr); 3350c52818fSSatish Balay ierr = VecGetOwnershipRange(s, &low3, &high3);CHKERRQ(ierr); 336f06e3bfaSBarry Smith if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths"); 3370c52818fSSatish Balay 3380c52818fSSatish Balay ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr); 339050fc7a3SBarry Smith ierr = VecDestroy(&ls->stepdirection);CHKERRQ(ierr); 340050fc7a3SBarry Smith ls->stepdirection = s; 3410c52818fSSatish Balay 3420c52818fSSatish Balay ierr = TaoLineSearchSetUp(ls);CHKERRQ(ierr); 343f06e3bfaSBarry Smith if (!ls->ops->apply) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine"); 3440c52818fSSatish Balay ls->nfeval=0; 3450c52818fSSatish Balay ls->ngeval=0; 3460c52818fSSatish Balay ls->nfgeval=0; 3470c52818fSSatish Balay /* Check parameter values */ 3480c52818fSSatish Balay if (ls->ftol < 0.0) { 349f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);CHKERRQ(ierr); 3500c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3510c52818fSSatish Balay } 3520c52818fSSatish Balay if (ls->rtol < 0.0) { 353f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);CHKERRQ(ierr); 3540c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3550c52818fSSatish Balay } 3560c52818fSSatish Balay if (ls->gtol < 0.0) { 357f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);CHKERRQ(ierr); 3580c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3590c52818fSSatish Balay } 3600c52818fSSatish Balay if (ls->stepmin < 0.0) { 361f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);CHKERRQ(ierr); 3620c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3630c52818fSSatish Balay } 3640c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 365f06e3bfaSBarry Smith ierr = PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);CHKERRQ(ierr); 3660c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3670c52818fSSatish Balay } 3680c52818fSSatish Balay if (ls->max_funcs < 0) { 3690c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);CHKERRQ(ierr); 3700c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 3710c52818fSSatish Balay } 3720c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 373f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);CHKERRQ(ierr); 3740c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_INFORNAN; 3750c52818fSSatish Balay } 3760c52818fSSatish Balay 377302440fdSBarry Smith ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 3780c52818fSSatish Balay ierr = VecDestroy(&ls->start_x);CHKERRQ(ierr); 3790c52818fSSatish Balay ls->start_x = x; 380050fc7a3SBarry Smith 381*0ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0);CHKERRQ(ierr); 3820c52818fSSatish Balay ierr = (*ls->ops->apply)(ls,x,f,g,s);CHKERRQ(ierr); 383*0ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0);CHKERRQ(ierr); 3840c52818fSSatish Balay *reason=ls->reason; 3850c52818fSSatish Balay ls->new_f = *f; 3860c52818fSSatish Balay 3870c52818fSSatish Balay if (steplength) { 3880c52818fSSatish Balay *steplength=ls->step; 3890c52818fSSatish Balay } 3900c52818fSSatish Balay 391fbe0838dSJason Sarich ierr = TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");CHKERRQ(ierr); 3920c52818fSSatish Balay PetscFunctionReturn(0); 3930c52818fSSatish Balay } 3940c52818fSSatish Balay 3950c52818fSSatish Balay /*@C 3960c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 3970c52818fSSatish Balay 3980c52818fSSatish Balay Collective on TaoLineSearch 3990c52818fSSatish Balay 4000c52818fSSatish Balay Input Parameters: 4010c52818fSSatish Balay + ls - the TaoLineSearch context 4020c52818fSSatish Balay - type - a known method 4030c52818fSSatish Balay 4040c52818fSSatish Balay Available methods include: 4050c52818fSSatish Balay + more-thuente 4060c52818fSSatish Balay . gpcg 4070c52818fSSatish Balay - unit - Do not perform any line search 4080c52818fSSatish Balay 4090c52818fSSatish Balay 4100c52818fSSatish Balay Options Database Keys: 4110c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 4120c52818fSSatish Balay 4130c52818fSSatish Balay Level: beginner 4140c52818fSSatish Balay 4150c52818fSSatish Balay 4160c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply() 4170c52818fSSatish Balay 4180c52818fSSatish Balay @*/ 4190c52818fSSatish Balay 420dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 4210c52818fSSatish Balay { 4220c52818fSSatish Balay PetscErrorCode ierr; 4230c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 4240c52818fSSatish Balay PetscBool flg; 4250c52818fSSatish Balay 4260c52818fSSatish Balay PetscFunctionBegin; 4270c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4280c52818fSSatish Balay PetscValidCharPointer(type,2); 4290c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg);CHKERRQ(ierr); 4300c52818fSSatish Balay if (flg) PetscFunctionReturn(0); 4310c52818fSSatish Balay 4320c52818fSSatish Balay ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);CHKERRQ(ierr); 4330c52818fSSatish Balay if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 4340c52818fSSatish Balay if (ls->ops->destroy) { 4350c52818fSSatish Balay ierr = (*(ls)->ops->destroy)(ls);CHKERRQ(ierr); 4360c52818fSSatish Balay } 4370c52818fSSatish Balay ls->max_funcs=30; 4380c52818fSSatish Balay ls->ftol = 0.0001; 4390c52818fSSatish Balay ls->gtol = 0.9; 4406f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 4416f4723b1SBarry Smith ls->rtol = 1.0e-5; 4426f4723b1SBarry Smith #else 4430c52818fSSatish Balay ls->rtol = 1.0e-10; 4446f4723b1SBarry Smith #endif 4450c52818fSSatish Balay ls->stepmin=1.0e-20; 4460c52818fSSatish Balay ls->stepmax=1.0e+20; 4470c52818fSSatish Balay 4480c52818fSSatish Balay ls->nfeval=0; 4490c52818fSSatish Balay ls->ngeval=0; 4500c52818fSSatish Balay ls->nfgeval=0; 4510c52818fSSatish Balay ls->ops->setup=0; 4520c52818fSSatish Balay ls->ops->apply=0; 4530c52818fSSatish Balay ls->ops->view=0; 4540c52818fSSatish Balay ls->ops->setfromoptions=0; 4550c52818fSSatish Balay ls->ops->destroy=0; 4560c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4570c52818fSSatish Balay ierr = (*r)(ls);CHKERRQ(ierr); 4580c52818fSSatish Balay ierr = PetscObjectChangeTypeName((PetscObject)ls, type);CHKERRQ(ierr); 4590c52818fSSatish Balay PetscFunctionReturn(0); 4600c52818fSSatish Balay } 4610c52818fSSatish Balay 4620c52818fSSatish Balay /*@ 4630c52818fSSatish Balay TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 4640c52818fSSatish Balay options. 4650c52818fSSatish Balay 4660c52818fSSatish Balay Collective on TaoLineSearch 4670c52818fSSatish Balay 4680c52818fSSatish Balay Input Paremeter: 4690c52818fSSatish Balay . ls - the TaoLineSearch context 4700c52818fSSatish Balay 4710c52818fSSatish Balay Options Database Keys: 4720c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 4730c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 4740c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 4750c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 4760c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 4770c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 4780c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 4790c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 4800c52818fSSatish Balay 4810c52818fSSatish Balay Level: beginner 4820c52818fSSatish Balay @*/ 4830c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 4840c52818fSSatish Balay { 4850c52818fSSatish Balay PetscErrorCode ierr; 4868caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 4870c52818fSSatish Balay char type[256]; 4880c52818fSSatish Balay PetscBool flg; 489f06e3bfaSBarry Smith 4900c52818fSSatish Balay PetscFunctionBegin; 4910c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4920c52818fSSatish Balay ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr); 4930c52818fSSatish Balay if (((PetscObject)ls)->type_name) { 4940c52818fSSatish Balay default_type = ((PetscObject)ls)->type_name; 4950c52818fSSatish Balay } 4960c52818fSSatish Balay /* Check for type from options */ 4970c52818fSSatish Balay ierr = PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);CHKERRQ(ierr); 4980c52818fSSatish Balay if (flg) { 4990c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,type);CHKERRQ(ierr); 5000c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 501302440fdSBarry Smith ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 5020c52818fSSatish Balay } 5030c52818fSSatish Balay 50494ae4db5SBarry Smith ierr = PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);CHKERRQ(ierr); 50594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);CHKERRQ(ierr); 50694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);CHKERRQ(ierr); 50794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);CHKERRQ(ierr); 50894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);CHKERRQ(ierr); 50994ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);CHKERRQ(ierr); 5100c52818fSSatish Balay if (ls->ops->setfromoptions) { 5113a004c28SBarry Smith ierr = (*ls->ops->setfromoptions)(PetscOptionsObject,ls);CHKERRQ(ierr); 5120c52818fSSatish Balay } 5130c52818fSSatish Balay ierr = PetscOptionsEnd();CHKERRQ(ierr); 5140c52818fSSatish Balay PetscFunctionReturn(0); 5150c52818fSSatish Balay } 5160c52818fSSatish Balay 5170c52818fSSatish Balay /*@C 5180c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5190c52818fSSatish Balay 5200c52818fSSatish Balay Not Collective 5210c52818fSSatish Balay 5220c52818fSSatish Balay Input Parameter: 5230c52818fSSatish Balay . ls - the TaoLineSearch context 5240c52818fSSatish Balay 5250c52818fSSatish Balay Output Paramter: 5260c52818fSSatish Balay . type - the line search algorithm in effect 5270c52818fSSatish Balay 5280c52818fSSatish Balay Level: developer 5290c52818fSSatish Balay 5300c52818fSSatish Balay @*/ 531dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 5320c52818fSSatish Balay { 5330c52818fSSatish Balay PetscFunctionBegin; 5340c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5350c52818fSSatish Balay PetscValidPointer(type,2); 5360c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 5370c52818fSSatish Balay PetscFunctionReturn(0); 5380c52818fSSatish Balay } 5390c52818fSSatish Balay 5400c52818fSSatish Balay /*@ 5410c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 5420c52818fSSatish Balay routines used by the line search in last application (not cumulative). 5430c52818fSSatish Balay 5440c52818fSSatish Balay Not Collective 5450c52818fSSatish Balay 5460c52818fSSatish Balay Input Parameter: 5470c52818fSSatish Balay . ls - the TaoLineSearch context 5480c52818fSSatish Balay 5490c52818fSSatish Balay Output Parameters: 5500c52818fSSatish Balay + nfeval - number of function evaluations 5510c52818fSSatish Balay . ngeval - number of gradient evaluations 5520c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 5530c52818fSSatish Balay 5540c52818fSSatish Balay Level: intermediate 5550c52818fSSatish Balay 5560c52818fSSatish Balay Note: 557441846f8SBarry Smith If the line search is using the Tao objective and gradient 558441846f8SBarry Smith routines directly (see TaoLineSearchUseTaoRoutines()), then TAO 5590c52818fSSatish Balay is already counting the number of evaluations. 5600c52818fSSatish Balay 5610c52818fSSatish Balay @*/ 5620c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 5630c52818fSSatish Balay { 5640c52818fSSatish Balay PetscFunctionBegin; 5650c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5660c52818fSSatish Balay *nfeval = ls->nfeval; 5670c52818fSSatish Balay *ngeval = ls->ngeval; 5680c52818fSSatish Balay *nfgeval = ls->nfgeval; 5690c52818fSSatish Balay PetscFunctionReturn(0); 5700c52818fSSatish Balay } 5710c52818fSSatish Balay 5720c52818fSSatish Balay /*@ 573441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 574441846f8SBarry Smith Tao evaluation routines. 5750c52818fSSatish Balay 5760c52818fSSatish Balay Not Collective 5770c52818fSSatish Balay 5780c52818fSSatish Balay Input Parameter: 5790c52818fSSatish Balay . ls - the TaoLineSearch context 5800c52818fSSatish Balay 5810c52818fSSatish Balay Output Parameter: 582441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines, 5830c52818fSSatish Balay otherwise PETSC_FALSE 5840c52818fSSatish Balay 5850c52818fSSatish Balay Level: developer 5860c52818fSSatish Balay @*/ 587441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 5880c52818fSSatish Balay { 5890c52818fSSatish Balay PetscFunctionBegin; 5900c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 591f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 5920c52818fSSatish Balay PetscFunctionReturn(0); 5930c52818fSSatish Balay } 5940c52818fSSatish Balay 5950c52818fSSatish Balay /*@C 5960c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 5970c52818fSSatish Balay 5980c52818fSSatish Balay Logically Collective on TaoLineSearch 5990c52818fSSatish Balay 6000c52818fSSatish Balay Input Parameter: 6010c52818fSSatish Balay + ls - the TaoLineSearch context 6020c52818fSSatish Balay . func - the objective function evaluation routine 6030c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6040c52818fSSatish Balay 6050c52818fSSatish Balay Calling sequence of func: 6060c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 6070c52818fSSatish Balay 6080c52818fSSatish Balay + x - input vector 6090c52818fSSatish Balay . f - function value 6100c52818fSSatish Balay - ctx (optional) user-defined context 6110c52818fSSatish Balay 6120c52818fSSatish Balay Level: beginner 6130c52818fSSatish Balay 6140c52818fSSatish Balay Note: 6150c52818fSSatish Balay Use this routine only if you want the line search objective 616441846f8SBarry Smith evaluation routine to be different from the Tao's objective 6170c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6180c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6190c52818fSSatish Balay 6200c52818fSSatish Balay Note: 6210c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6220c52818fSSatish Balay line search, application programmers should be wary of overriding the 6230c52818fSSatish Balay default objective routine. 6240c52818fSSatish Balay 625441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 6260c52818fSSatish Balay @*/ 6270c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx) 6280c52818fSSatish Balay { 6290c52818fSSatish Balay PetscFunctionBegin; 6300c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6310c52818fSSatish Balay 6320c52818fSSatish Balay ls->ops->computeobjective=func; 6330c52818fSSatish Balay if (ctx) ls->userctx_func=ctx; 6340c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 6350c52818fSSatish Balay PetscFunctionReturn(0); 6360c52818fSSatish Balay } 6370c52818fSSatish Balay 6380c52818fSSatish Balay /*@C 6390c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 6400c52818fSSatish Balay 6410c52818fSSatish Balay Logically Collective on TaoLineSearch 6420c52818fSSatish Balay 6430c52818fSSatish Balay Input Parameter: 6440c52818fSSatish Balay + ls - the TaoLineSearch context 6450c52818fSSatish Balay . func - the gradient evaluation routine 6460c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6470c52818fSSatish Balay 6480c52818fSSatish Balay Calling sequence of func: 6490c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, Vec g, void *ctx); 6500c52818fSSatish Balay 6510c52818fSSatish Balay + x - input vector 6520c52818fSSatish Balay . g - gradient vector 6530c52818fSSatish Balay - ctx (optional) user-defined context 6540c52818fSSatish Balay 6550c52818fSSatish Balay Level: beginner 6560c52818fSSatish Balay 6570c52818fSSatish Balay Note: 6580c52818fSSatish Balay Use this routine only if you want the line search gradient 659441846f8SBarry Smith evaluation routine to be different from the Tao's gradient 6600c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6610c52818fSSatish Balay the line search function and/or function/gradient routine. 6620c52818fSSatish Balay 6630c52818fSSatish Balay Note: 6640c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 6650c52818fSSatish Balay line search, application programmers should be wary of overriding the 6660c52818fSSatish Balay default gradient routine. 6670c52818fSSatish Balay 668441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 6690c52818fSSatish Balay @*/ 6700c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx) 6710c52818fSSatish Balay { 6720c52818fSSatish Balay PetscFunctionBegin; 6730c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6740c52818fSSatish Balay ls->ops->computegradient=func; 6750c52818fSSatish Balay if (ctx) ls->userctx_grad=ctx; 6760c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 6770c52818fSSatish Balay PetscFunctionReturn(0); 6780c52818fSSatish Balay } 6790c52818fSSatish Balay 6800c52818fSSatish Balay /*@C 6810c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 6820c52818fSSatish Balay 6830c52818fSSatish Balay Logically Collective on TaoLineSearch 6840c52818fSSatish Balay 6850c52818fSSatish Balay Input Parameter: 6860c52818fSSatish Balay + ls - the TaoLineSearch context 6870c52818fSSatish Balay . func - the objective and gradient evaluation routine 6880c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6890c52818fSSatish Balay 6900c52818fSSatish Balay Calling sequence of func: 6910c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 6920c52818fSSatish Balay 6930c52818fSSatish Balay + x - input vector 6940c52818fSSatish Balay . f - function value 6950c52818fSSatish Balay . g - gradient vector 6960c52818fSSatish Balay - ctx (optional) user-defined context 6970c52818fSSatish Balay 6980c52818fSSatish Balay Level: beginner 6990c52818fSSatish Balay 7000c52818fSSatish Balay Note: 7010c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 702441846f8SBarry Smith evaluation routines to be different from the Tao's objective 7030c52818fSSatish Balay and gradient evaluation routines. 7040c52818fSSatish Balay 7050c52818fSSatish Balay Note: 7060c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7070c52818fSSatish Balay line search, application programmers should be wary of overriding the 7080c52818fSSatish Balay default objective routine. 7090c52818fSSatish Balay 710441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines() 7110c52818fSSatish Balay @*/ 7120c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx) 7130c52818fSSatish Balay { 7140c52818fSSatish Balay PetscFunctionBegin; 7150c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7160c52818fSSatish Balay ls->ops->computeobjectiveandgradient=func; 7170c52818fSSatish Balay if (ctx) ls->userctx_funcgrad=ctx; 7180c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7190c52818fSSatish Balay PetscFunctionReturn(0); 7200c52818fSSatish Balay } 7210c52818fSSatish Balay 7220c52818fSSatish Balay /*@C 7230c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 7240c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 7250c52818fSSatish Balay Sometimes it is more efficient to compute the inner product of the gradient 7260c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 7270c52818fSSatish Balay the line search typically needs of the gradient. 7280c52818fSSatish Balay 7290c52818fSSatish Balay Logically Collective on TaoLineSearch 7300c52818fSSatish Balay 7310c52818fSSatish Balay Input Parameter: 7320c52818fSSatish Balay + ls - the TaoLineSearch context 7330c52818fSSatish Balay . func - the objective and gradient evaluation routine 7340c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7350c52818fSSatish Balay 7360c52818fSSatish Balay Calling sequence of func: 7370c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 7380c52818fSSatish Balay 7390c52818fSSatish Balay + x - input vector 7400c52818fSSatish Balay . s - step direction 7410c52818fSSatish Balay . f - function value 7420c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 7430c52818fSSatish Balay - ctx (optional) user-defined context 7440c52818fSSatish Balay 7450c52818fSSatish Balay Note: The gradient will still need to be computed at the end of the line 7460c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 7470c52818fSSatish Balay routine 7480c52818fSSatish Balay 7490c52818fSSatish Balay Note: Bounded line searches (those used in bounded optimization algorithms) 7500c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 7510c52818fSSatish Balay x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength() 7520c52818fSSatish Balay 7530c52818fSSatish Balay Level: advanced 7540c52818fSSatish Balay 7550c52818fSSatish Balay Note: 7560c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7570c52818fSSatish Balay line search, application programmers should be wary of overriding the 7580c52818fSSatish Balay default objective routine. 7590c52818fSSatish Balay 760441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines() 7610c52818fSSatish Balay @*/ 7620c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx) 7630c52818fSSatish Balay { 7640c52818fSSatish Balay PetscFunctionBegin; 7650c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7660c52818fSSatish Balay ls->ops->computeobjectiveandgts=func; 7670c52818fSSatish Balay if (ctx) ls->userctx_funcgts=ctx; 7680c52818fSSatish Balay ls->usegts = PETSC_TRUE; 7690c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 7700c52818fSSatish Balay PetscFunctionReturn(0); 7710c52818fSSatish Balay } 7720c52818fSSatish Balay 7730c52818fSSatish Balay /*@C 774441846f8SBarry Smith TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the 775441846f8SBarry Smith objective and gradient evaluation routines from the given Tao object. 7760c52818fSSatish Balay 7770c52818fSSatish Balay Logically Collective on TaoLineSearch 7780c52818fSSatish Balay 7790c52818fSSatish Balay Input Parameter: 7800c52818fSSatish Balay + ls - the TaoLineSearch context 781441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines 7820c52818fSSatish Balay 7830c52818fSSatish Balay Level: developer 7840c52818fSSatish Balay 7850c52818fSSatish Balay .seealso: TaoLineSearchCreate() 7860c52818fSSatish Balay @*/ 787441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts) 7880c52818fSSatish Balay { 7890c52818fSSatish Balay PetscFunctionBegin; 7900c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 791441846f8SBarry Smith PetscValidHeaderSpecific(ts,TAO_CLASSID,1); 792441846f8SBarry Smith ls->tao = ts; 7930c52818fSSatish Balay ls->usetaoroutines=PETSC_TRUE; 7940c52818fSSatish Balay PetscFunctionReturn(0); 7950c52818fSSatish Balay } 7960c52818fSSatish Balay 7970c52818fSSatish Balay /*@ 7980c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 7990c52818fSSatish Balay 8000c52818fSSatish Balay Collective on TaoLineSearch 8010c52818fSSatish Balay 8020c52818fSSatish Balay Input Parameters: 8030c52818fSSatish Balay + ls - the TaoLineSearch context 8040c52818fSSatish Balay - x - input vector 8050c52818fSSatish Balay 8060c52818fSSatish Balay Output Parameter: 8070c52818fSSatish Balay . f - Objective value at X 8080c52818fSSatish Balay 80995452b02SPatrick Sanan Notes: 81095452b02SPatrick Sanan TaoLineSearchComputeObjective() is typically used within line searches 8110c52818fSSatish Balay so most users would not generally call this routine themselves. 8120c52818fSSatish Balay 8130c52818fSSatish Balay Level: developer 8140c52818fSSatish Balay 8150c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 8160c52818fSSatish Balay @*/ 8170c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 8180c52818fSSatish Balay { 8190c52818fSSatish Balay PetscErrorCode ierr; 8200c52818fSSatish Balay Vec gdummy; 8210c52818fSSatish Balay PetscReal gts; 822f06e3bfaSBarry Smith 8230c52818fSSatish Balay PetscFunctionBegin; 8240c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8250c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 8260c52818fSSatish Balay PetscValidPointer(f,3); 8270c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8280c52818fSSatish Balay if (ls->usetaoroutines) { 829441846f8SBarry Smith ierr = TaoComputeObjective(ls->tao,x,f);CHKERRQ(ierr); 8300c52818fSSatish Balay } else { 831f06e3bfaSBarry 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"); 832*0ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 8330c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective routine"); 8340c52818fSSatish Balay if (ls->ops->computeobjective) { 8350c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 8360c52818fSSatish Balay } else if (ls->ops->computeobjectiveandgradient) { 8370c52818fSSatish Balay ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr); 8380c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr); 8390c52818fSSatish Balay ierr = VecDestroy(&gdummy);CHKERRQ(ierr); 8400c52818fSSatish Balay } else { 8410c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts);CHKERRQ(ierr); 8420c52818fSSatish Balay } 8430c52818fSSatish Balay PetscStackPop; 844*0ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 8450c52818fSSatish Balay } 8460c52818fSSatish Balay ls->nfeval++; 8470c52818fSSatish Balay PetscFunctionReturn(0); 8480c52818fSSatish Balay } 8490c52818fSSatish Balay 8500c52818fSSatish Balay /*@ 8510c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 8520c52818fSSatish Balay 853441846f8SBarry Smith Collective on Tao 8540c52818fSSatish Balay 8550c52818fSSatish Balay Input Parameters: 8560c52818fSSatish Balay + ls - the TaoLineSearch context 8570c52818fSSatish Balay - x - input vector 8580c52818fSSatish Balay 8590c52818fSSatish Balay Output Parameter: 8600c52818fSSatish Balay + f - Objective value at X 8610c52818fSSatish Balay - g - Gradient vector at X 8620c52818fSSatish Balay 86395452b02SPatrick Sanan Notes: 86495452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 8650c52818fSSatish Balay so most users would not generally call this routine themselves. 8660c52818fSSatish Balay 8670c52818fSSatish Balay Level: developer 8680c52818fSSatish Balay 8690c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 8700c52818fSSatish Balay @*/ 8710c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 8720c52818fSSatish Balay { 8730c52818fSSatish Balay PetscErrorCode ierr; 874f06e3bfaSBarry Smith 8750c52818fSSatish Balay PetscFunctionBegin; 8760c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8770c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 8780c52818fSSatish Balay PetscValidPointer(f,3); 8790c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 8800c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8810c52818fSSatish Balay PetscCheckSameComm(ls,1,g,4); 8820c52818fSSatish Balay if (ls->usetaoroutines) { 883441846f8SBarry Smith ierr = TaoComputeObjectiveAndGradient(ls->tao,x,f,g);CHKERRQ(ierr); 8840c52818fSSatish Balay } else { 885f06e3bfaSBarry Smith if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 886f06e3bfaSBarry Smith if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 887*0ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 8880c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gradient routine"); 8890c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 8900c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr); 8910c52818fSSatish Balay } else { 8920c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 8930c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 8940c52818fSSatish Balay } 8950c52818fSSatish Balay PetscStackPop; 896*0ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 897335036cbSBarry Smith ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 898f06e3bfaSBarry Smith } 899fbe0838dSJason Sarich ls->nfgeval++; 9000c52818fSSatish Balay PetscFunctionReturn(0); 9010c52818fSSatish Balay } 9020c52818fSSatish Balay 9030c52818fSSatish Balay /*@ 9040c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 9050c52818fSSatish Balay 9060c52818fSSatish Balay Collective on TaoLineSearch 9070c52818fSSatish Balay 9080c52818fSSatish Balay Input Parameters: 9090c52818fSSatish Balay + ls - the TaoLineSearch context 9100c52818fSSatish Balay - x - input vector 9110c52818fSSatish Balay 9120c52818fSSatish Balay Output Parameter: 9130c52818fSSatish Balay . g - gradient vector 9140c52818fSSatish Balay 91595452b02SPatrick Sanan Notes: 91695452b02SPatrick Sanan TaoComputeGradient() is typically used within line searches 9170c52818fSSatish Balay so most users would not generally call this routine themselves. 9180c52818fSSatish Balay 9190c52818fSSatish Balay Level: developer 9200c52818fSSatish Balay 9210c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 9220c52818fSSatish Balay @*/ 9230c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 9240c52818fSSatish Balay { 9250c52818fSSatish Balay PetscErrorCode ierr; 9260c52818fSSatish Balay PetscReal fdummy; 927f06e3bfaSBarry Smith 9280c52818fSSatish Balay PetscFunctionBegin; 9290c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9300c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9310c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,3); 9320c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9330c52818fSSatish Balay PetscCheckSameComm(ls,1,g,3); 9340c52818fSSatish Balay if (ls->usetaoroutines) { 935441846f8SBarry Smith ierr = TaoComputeGradient(ls->tao,x,g);CHKERRQ(ierr); 9360c52818fSSatish Balay } else { 937f06e3bfaSBarry Smith if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 938*0ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9390c52818fSSatish Balay PetscStackPush("TaoLineSearch user gradient routine"); 9400c52818fSSatish Balay if (ls->ops->computegradient) { 9410c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 9420c52818fSSatish Balay } else { 9430c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr); 9440c52818fSSatish Balay } 9450c52818fSSatish Balay PetscStackPop; 946*0ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9470c52818fSSatish Balay } 9480c52818fSSatish Balay ls->ngeval++; 9490c52818fSSatish Balay PetscFunctionReturn(0); 9500c52818fSSatish Balay } 9510c52818fSSatish Balay 9520c52818fSSatish Balay /*@ 9530c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 9540c52818fSSatish Balay 955441846f8SBarry Smith Collective on Tao 9560c52818fSSatish Balay 9570c52818fSSatish Balay Input Parameters: 9580c52818fSSatish Balay + ls - the TaoLineSearch context 9590c52818fSSatish Balay - x - input vector 9600c52818fSSatish Balay 9610c52818fSSatish Balay Output Parameter: 9620c52818fSSatish Balay + f - Objective value at X 9630c52818fSSatish Balay - gts - inner product of gradient and step direction at X 9640c52818fSSatish Balay 96595452b02SPatrick Sanan Notes: 96695452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 9670c52818fSSatish Balay so most users would not generally call this routine themselves. 9680c52818fSSatish Balay 9690c52818fSSatish Balay Level: developer 9700c52818fSSatish Balay 9710c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 9720c52818fSSatish Balay @*/ 9730c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 9740c52818fSSatish Balay { 9750c52818fSSatish Balay PetscErrorCode ierr; 9760c52818fSSatish Balay PetscFunctionBegin; 9770c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9780c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9790c52818fSSatish Balay PetscValidPointer(f,3); 9800c52818fSSatish Balay PetscValidPointer(gts,4); 9810c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 982f06e3bfaSBarry Smith if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 983*0ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9840c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gts routine"); 9850c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr); 9860c52818fSSatish Balay PetscStackPop; 987*0ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 988335036cbSBarry Smith ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 9890c52818fSSatish Balay ls->nfeval++; 9900c52818fSSatish Balay PetscFunctionReturn(0); 9910c52818fSSatish Balay } 9920c52818fSSatish Balay 9930c52818fSSatish Balay /*@ 9940c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 9950c52818fSSatish Balay 9960c52818fSSatish Balay Collective on TaoLineSearch 9970c52818fSSatish Balay 9980c52818fSSatish Balay Input Parameter: 9990c52818fSSatish Balay . ls - the TaoLineSearch context 10000c52818fSSatish Balay 10010c52818fSSatish Balay Output Parameter: 10020c52818fSSatish Balay + x - the new solution 10030c52818fSSatish Balay . f - the objective function value at x 10040c52818fSSatish Balay . g - the gradient at x 10050c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 10060c52818fSSatish Balay - reason - the reason why the line search terminated 10070c52818fSSatish Balay 10080c52818fSSatish Balay reason will be set to one of: 10090c52818fSSatish Balay 10100c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 10110c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 10120c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 10130c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 10140c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 10150c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 10160c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 10170c52818fSSatish Balay 10180c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 10190c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 10200c52818fSSatish Balay 10210c52818fSSatish Balay + TAOLINESEARCH_SUCCESS - successful line search 10220c52818fSSatish Balay 10230c52818fSSatish Balay Level: developer 10240c52818fSSatish Balay 10250c52818fSSatish Balay @*/ 1026e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 10270c52818fSSatish Balay { 10280c52818fSSatish Balay PetscErrorCode ierr; 1029f06e3bfaSBarry Smith 10300c52818fSSatish Balay PetscFunctionBegin; 10310c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10320c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 10330c52818fSSatish Balay PetscValidPointer(f,3); 10340c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 10350c52818fSSatish Balay PetscValidIntPointer(reason,6); 10360c52818fSSatish Balay if (ls->new_x) { 10370c52818fSSatish Balay ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr); 10380c52818fSSatish Balay } 10390c52818fSSatish Balay *f = ls->new_f; 10400c52818fSSatish Balay if (ls->new_g) { 10410c52818fSSatish Balay ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr); 10420c52818fSSatish Balay } 10430c52818fSSatish Balay if (steplength) { 10440c52818fSSatish Balay *steplength=ls->step; 10450c52818fSSatish Balay } 10460c52818fSSatish Balay *reason = ls->reason; 10470c52818fSSatish Balay PetscFunctionReturn(0); 10480c52818fSSatish Balay } 10490c52818fSSatish Balay 10500c52818fSSatish Balay /*@ 10510c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 10520c52818fSSatish Balay search. 10530c52818fSSatish Balay 10540c52818fSSatish Balay Not Collective 10550c52818fSSatish Balay 10560c52818fSSatish Balay Input Parameter: 10570c52818fSSatish Balay . ls - the TaoLineSearch context 10580c52818fSSatish Balay 10590c52818fSSatish Balay Output Parameter: 10600c52818fSSatish Balay . x - The initial point of the line search 10610c52818fSSatish Balay 10620c52818fSSatish Balay Level: intermediate 10630c52818fSSatish Balay @*/ 10640c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 10650c52818fSSatish Balay { 10660c52818fSSatish Balay PetscFunctionBegin; 10670c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10680c52818fSSatish Balay if (x) { 10690c52818fSSatish Balay *x = ls->start_x; 10700c52818fSSatish Balay } 10710c52818fSSatish Balay PetscFunctionReturn(0); 10720c52818fSSatish Balay } 10730c52818fSSatish Balay 10740c52818fSSatish Balay /*@ 10750c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 10760c52818fSSatish Balay search. 10770c52818fSSatish Balay 10780c52818fSSatish Balay Not Collective 10790c52818fSSatish Balay 10800c52818fSSatish Balay Input Parameter: 10810c52818fSSatish Balay . ls - the TaoLineSearch context 10820c52818fSSatish Balay 10830c52818fSSatish Balay Output Parameter: 10840c52818fSSatish Balay . s - the step direction of the line search 10850c52818fSSatish Balay 10860c52818fSSatish Balay Level: advanced 10870c52818fSSatish Balay @*/ 10880c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 10890c52818fSSatish Balay { 10900c52818fSSatish Balay PetscFunctionBegin; 10910c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10920c52818fSSatish Balay if (s) { 10930c52818fSSatish Balay *s = ls->stepdirection; 10940c52818fSSatish Balay } 10950c52818fSSatish Balay PetscFunctionReturn(0); 10960c52818fSSatish Balay 10970c52818fSSatish Balay } 10980c52818fSSatish Balay 10990c52818fSSatish Balay /*@ 11000c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 11010c52818fSSatish Balay 11020c52818fSSatish Balay Not Collective 11030c52818fSSatish Balay 11040c52818fSSatish Balay Input Parameter: 11050c52818fSSatish Balay . ls - the TaoLineSearch context 11060c52818fSSatish Balay 11070c52818fSSatish Balay Output Parameter: 11080c52818fSSatish Balay . f - the objective value at the full step length 11090c52818fSSatish Balay 11100c52818fSSatish Balay Level: developer 11110c52818fSSatish Balay @*/ 11120c52818fSSatish Balay 11130c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 11140c52818fSSatish Balay { 11150c52818fSSatish Balay PetscFunctionBegin; 11160c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11170c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 11180c52818fSSatish Balay PetscFunctionReturn(0); 11190c52818fSSatish Balay } 11200c52818fSSatish Balay 11210c52818fSSatish Balay /*@ 11220c52818fSSatish Balay TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 11230c52818fSSatish Balay 1124441846f8SBarry Smith Logically Collective on Tao 11250c52818fSSatish Balay 11260c52818fSSatish Balay Input Parameters: 11270c52818fSSatish Balay + ls - the TaoLineSearch context 11280c52818fSSatish Balay . xl - vector of lower bounds 11290c52818fSSatish Balay - xu - vector of upper bounds 11300c52818fSSatish Balay 11310c52818fSSatish Balay Note: If the variable bounds are not set with this routine, then 1132e270355aSBarry Smith PETSC_NINFINITY and PETSC_INFINITY are assumed 11330c52818fSSatish Balay 11340c52818fSSatish Balay Level: beginner 11350c52818fSSatish Balay 11360c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 11370c52818fSSatish Balay @*/ 11380c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 11390c52818fSSatish Balay { 11400c52818fSSatish Balay PetscFunctionBegin; 11410c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11420c52818fSSatish Balay PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 11430c52818fSSatish Balay PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 11440c52818fSSatish Balay ls->lower = xl; 11450c52818fSSatish Balay ls->upper = xu; 11460c52818fSSatish Balay ls->bounded = 1; 11470c52818fSSatish Balay PetscFunctionReturn(0); 11480c52818fSSatish Balay } 11490c52818fSSatish Balay 11500c52818fSSatish Balay /*@ 11510c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 11520c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 11530c52818fSSatish Balay 11540c52818fSSatish Balay Logically Collective on TaoLineSearch 11550c52818fSSatish Balay 11560c52818fSSatish Balay Input Parameters: 11570c52818fSSatish Balay + ls - the TaoLineSearch context 11580c52818fSSatish Balay - s - the initial step size 11590c52818fSSatish Balay 11600c52818fSSatish Balay Level: intermediate 11610c52818fSSatish Balay 11620c52818fSSatish Balay .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply() 11630c52818fSSatish Balay @*/ 11640c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 11650c52818fSSatish Balay { 11660c52818fSSatish Balay PetscFunctionBegin; 11670c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11680c52818fSSatish Balay ls->initstep = s; 11690c52818fSSatish Balay PetscFunctionReturn(0); 11700c52818fSSatish Balay } 11710c52818fSSatish Balay 11720c52818fSSatish Balay /*@ 11730c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 11740c52818fSSatish Balay 11750c52818fSSatish Balay Not Collective 11760c52818fSSatish Balay 11770c52818fSSatish Balay Input Parameters: 11780c52818fSSatish Balay . ls - the TaoLineSearch context 11790c52818fSSatish Balay 11800c52818fSSatish Balay Output Parameters 11810c52818fSSatish Balay . s - the current step length 11820c52818fSSatish Balay 11830c52818fSSatish Balay Level: beginner 11840c52818fSSatish Balay 11850c52818fSSatish Balay .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply() 11860c52818fSSatish Balay @*/ 11870c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 11880c52818fSSatish Balay { 11890c52818fSSatish Balay PetscFunctionBegin; 11900c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11910c52818fSSatish Balay *s = ls->step; 11920c52818fSSatish Balay PetscFunctionReturn(0); 11930c52818fSSatish Balay } 11940c52818fSSatish Balay 1195*0ebee16dSLisandro Dalcin /*@C 11960c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 11970c52818fSSatish Balay 11980c52818fSSatish Balay Not collective 11990c52818fSSatish Balay 12000c52818fSSatish Balay Input Parameters: 12010c52818fSSatish Balay + sname - name of a new user-defined solver 12020c52818fSSatish Balay - func - routine to Create method context 12030c52818fSSatish Balay 12040c52818fSSatish Balay Notes: 12050c52818fSSatish Balay TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 12060c52818fSSatish Balay 12070c52818fSSatish Balay Sample usage: 12080c52818fSSatish Balay .vb 12090c52818fSSatish Balay TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 12100c52818fSSatish Balay .ve 12110c52818fSSatish Balay 12120c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 12130c52818fSSatish Balay $ TaoLineSearchSetType(ls,"my_linesearch") 12140c52818fSSatish Balay or at runtime via the option 12150c52818fSSatish Balay $ -tao_ls_type my_linesearch 12160c52818fSSatish Balay 12170c52818fSSatish Balay Level: developer 12180c52818fSSatish Balay 1219*0ebee16dSLisandro Dalcin @*/ 12200c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 12210c52818fSSatish Balay { 12220c52818fSSatish Balay PetscErrorCode ierr; 12230c52818fSSatish Balay PetscFunctionBegin; 12240c52818fSSatish Balay ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr); 12250c52818fSSatish Balay PetscFunctionReturn(0); 12260c52818fSSatish Balay } 12270c52818fSSatish Balay 12280c52818fSSatish Balay /*@C 12290c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 12300c52818fSSatish Balay for all TaoLineSearch options in the database. 12310c52818fSSatish Balay 12320c52818fSSatish Balay 12330c52818fSSatish Balay Collective on TaoLineSearch 12340c52818fSSatish Balay 12350c52818fSSatish Balay Input Parameters: 12360c52818fSSatish Balay + ls - the TaoLineSearch solver context 12370c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 12380c52818fSSatish Balay 12390c52818fSSatish Balay Notes: 12400c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12410c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12420c52818fSSatish Balay 12430c52818fSSatish Balay 12440c52818fSSatish Balay Level: advanced 12450c52818fSSatish Balay 12460c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 12470c52818fSSatish Balay @*/ 12480c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 12490c52818fSSatish Balay { 1250f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 12510c52818fSSatish Balay } 12520c52818fSSatish Balay 12530c52818fSSatish Balay /*@C 12540c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 12550c52818fSSatish Balay TaoLineSearch options in the database 12560c52818fSSatish Balay 12570c52818fSSatish Balay Not Collective 12580c52818fSSatish Balay 12590c52818fSSatish Balay Input Parameters: 12600c52818fSSatish Balay . ls - the TaoLineSearch context 12610c52818fSSatish Balay 12620c52818fSSatish Balay Output Parameters: 12630c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 12640c52818fSSatish Balay 126595452b02SPatrick Sanan Notes: 126695452b02SPatrick Sanan On the fortran side, the user should pass in a string 'prefix' of 12670c52818fSSatish Balay sufficient length to hold the prefix. 12680c52818fSSatish Balay 12690c52818fSSatish Balay Level: advanced 12700c52818fSSatish Balay 12710c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 12720c52818fSSatish Balay @*/ 12730c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 12740c52818fSSatish Balay { 1275f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 12760c52818fSSatish Balay } 12770c52818fSSatish Balay 12780c52818fSSatish Balay /*@C 12790c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 12800c52818fSSatish Balay TaoLineSearch options in the database. 12810c52818fSSatish Balay 12820c52818fSSatish Balay 12830c52818fSSatish Balay Logically Collective on TaoLineSearch 12840c52818fSSatish Balay 12850c52818fSSatish Balay Input Parameters: 12860c52818fSSatish Balay + ls - the TaoLineSearch context 12870c52818fSSatish Balay - prefix - the prefix string to prepend to all TAO option requests 12880c52818fSSatish Balay 12890c52818fSSatish Balay Notes: 12900c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12910c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12920c52818fSSatish Balay 12930c52818fSSatish Balay For example, to distinguish between the runtime options for two 12940c52818fSSatish Balay different line searches, one could call 12950c52818fSSatish Balay .vb 12960c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 12970c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 12980c52818fSSatish Balay .ve 12990c52818fSSatish Balay 13000c52818fSSatish Balay This would enable use of different options for each system, such as 13010c52818fSSatish Balay .vb 13020c52818fSSatish Balay -sys1_tao_ls_type mt 13030c52818fSSatish Balay -sys2_tao_ls_type armijo 13040c52818fSSatish Balay .ve 13050c52818fSSatish Balay 13060c52818fSSatish Balay Level: advanced 13070c52818fSSatish Balay 13080c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 13090c52818fSSatish Balay @*/ 13100c52818fSSatish Balay 13110c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 13120c52818fSSatish Balay { 1313f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 13140c52818fSSatish Balay } 1315