1ba92ff59SBarry Smith #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h> 30c52818fSSatish Balay 46c23d075SBarry Smith PetscFunctionList TaoLineSearchList = NULL; 50c52818fSSatish Balay 60c52818fSSatish Balay PetscClassId TAOLINESEARCH_CLASSID=0; 70ebee16dSLisandro Dalcin 80ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Apply; 90ebee16dSLisandro Dalcin PetscLogEvent TAOLINESEARCH_Eval; 100c52818fSSatish Balay 110c52818fSSatish Balay /*@C 12fe2efc57SMark TaoLineSearchViewFromOptions - View from Options 13fe2efc57SMark 14fe2efc57SMark Collective on TaoLineSearch 15fe2efc57SMark 16fe2efc57SMark Input Parameters: 17fe2efc57SMark + A - the Tao context 18736c3998SJose E. Roman . obj - Optional object 19736c3998SJose E. Roman - name - command line option 20fe2efc57SMark 21fe2efc57SMark Level: intermediate 22fe2efc57SMark .seealso: TaoLineSearch, TaoLineSearchView, PetscObjectViewFromOptions(), TaoLineSearchCreate() 23fe2efc57SMark @*/ 24fe2efc57SMark PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A,PetscObject obj,const char name[]) 25fe2efc57SMark { 26fe2efc57SMark PetscFunctionBegin; 27fe2efc57SMark PetscValidHeaderSpecific(A,TAOLINESEARCH_CLASSID,1); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 29fe2efc57SMark PetscFunctionReturn(0); 30fe2efc57SMark } 31fe2efc57SMark 32fe2efc57SMark /*@C 330c52818fSSatish Balay TaoLineSearchView - Prints information about the TaoLineSearch 340c52818fSSatish Balay 350c52818fSSatish Balay Collective on TaoLineSearch 360c52818fSSatish Balay 370c52818fSSatish Balay InputParameters: 38441846f8SBarry Smith + ls - the Tao context 390c52818fSSatish Balay - viewer - visualization context 400c52818fSSatish Balay 410c52818fSSatish Balay Options Database Key: 420c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search 430c52818fSSatish Balay 440c52818fSSatish Balay Notes: 450c52818fSSatish Balay The available visualization contexts include 460c52818fSSatish Balay + PETSC_VIEWER_STDOUT_SELF - standard output (default) 470c52818fSSatish Balay - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 480c52818fSSatish Balay output where only the first processor opens 490c52818fSSatish Balay the file. All other processors send their 500c52818fSSatish Balay data to the first processor to print. 510c52818fSSatish Balay 520c52818fSSatish Balay Level: beginner 530c52818fSSatish Balay 540c52818fSSatish Balay .seealso: PetscViewerASCIIOpen() 550c52818fSSatish Balay @*/ 560c52818fSSatish Balay 570c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 580c52818fSSatish Balay { 590c52818fSSatish Balay PetscBool isascii, isstring; 60dedfbcbeSJed Brown TaoLineSearchType type; 61f06e3bfaSBarry Smith 620c52818fSSatish Balay PetscFunctionBegin; 630c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 640c52818fSSatish Balay if (!viewer) { 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer)); 660c52818fSSatish Balay } 670c52818fSSatish Balay PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 680c52818fSSatish Balay PetscCheckSameComm(ls,1,viewer,2); 690c52818fSSatish Balay 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 720c52818fSSatish Balay if (isascii) { 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer)); 7463b15415SAlp Dener if (ls->ops->view) { 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 765f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->view)(ls,viewer)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 780c52818fSSatish Balay } 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval)); 850c52818fSSatish Balay 860c52818fSSatish Balay if (ls->bounded) { 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"using variable bounds\n")); 880c52818fSSatish Balay } 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 910c52818fSSatish Balay } else if (isstring) { 925f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchGetType(ls,&type)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringSPrintf(viewer," %-3.3s",type)); 940c52818fSSatish Balay } 950c52818fSSatish Balay PetscFunctionReturn(0); 960c52818fSSatish Balay } 970c52818fSSatish Balay 980c52818fSSatish Balay /*@C 990c52818fSSatish Balay TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 1000c52818fSSatish Balay line-searches will automatically create one. 1010c52818fSSatish Balay 102d083f849SBarry Smith Collective 1030c52818fSSatish Balay 1040c52818fSSatish Balay Input Parameter: 1050c52818fSSatish Balay . comm - MPI communicator 1060c52818fSSatish Balay 1070c52818fSSatish Balay Output Parameter: 1080c52818fSSatish Balay . newls - the new TaoLineSearch context 1090c52818fSSatish Balay 1100c52818fSSatish Balay Available methods include: 111147403d9SBarry Smith + more-thuente - the More-Thuente method 112147403d9SBarry Smith . gpcg - the GPCG method 1130c52818fSSatish Balay - unit - Do not perform any line search 1140c52818fSSatish Balay 1150c52818fSSatish Balay Options Database Keys: 1160c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 1170c52818fSSatish Balay 1180c52818fSSatish Balay Level: beginner 1190c52818fSSatish Balay 1200c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy() 1210c52818fSSatish Balay @*/ 1220c52818fSSatish Balay 1230c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 1240c52818fSSatish Balay { 1250c52818fSSatish Balay TaoLineSearch ls; 1260c52818fSSatish Balay 1270c52818fSSatish Balay PetscFunctionBegin; 1280c52818fSSatish Balay PetscValidPointer(newls,2); 1296c23d075SBarry Smith *newls = NULL; 1300c52818fSSatish Balay 1315f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchInitializePackage()); 1320c52818fSSatish Balay 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView)); 1340c52818fSSatish Balay ls->bounded = 0; 1350c52818fSSatish Balay ls->max_funcs=30; 1360c52818fSSatish Balay ls->ftol = 0.0001; 1370c52818fSSatish Balay ls->gtol = 0.9; 1386f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1396f4723b1SBarry Smith ls->rtol = 1.0e-5; 1406f4723b1SBarry Smith #else 1410c52818fSSatish Balay ls->rtol = 1.0e-10; 1426f4723b1SBarry Smith #endif 1430c52818fSSatish Balay ls->stepmin=1.0e-20; 1440c52818fSSatish Balay ls->stepmax=1.0e+20; 1450c52818fSSatish Balay ls->step=1.0; 1460c52818fSSatish Balay ls->nfeval=0; 1470c52818fSSatish Balay ls->ngeval=0; 1480c52818fSSatish Balay ls->nfgeval=0; 1490c52818fSSatish Balay 15083c8fe1dSLisandro Dalcin ls->ops->computeobjective = NULL; 15183c8fe1dSLisandro Dalcin ls->ops->computegradient = NULL; 15283c8fe1dSLisandro Dalcin ls->ops->computeobjectiveandgradient = NULL; 15383c8fe1dSLisandro Dalcin ls->ops->computeobjectiveandgts = NULL; 15483c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 15583c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 15683c8fe1dSLisandro Dalcin ls->ops->view = NULL; 15783c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 15883c8fe1dSLisandro Dalcin ls->ops->reset = NULL; 15983c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 16083c8fe1dSLisandro Dalcin ls->ops->monitor = NULL; 1612a0dac07SAlp Dener ls->usemonitor=PETSC_FALSE; 1620c52818fSSatish Balay ls->setupcalled=PETSC_FALSE; 1630c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 1640c52818fSSatish Balay *newls = ls; 1650c52818fSSatish Balay PetscFunctionReturn(0); 1660c52818fSSatish Balay } 1670c52818fSSatish Balay 1680c52818fSSatish Balay /*@ 1690c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 1700c52818fSSatish Balay of a Tao solver 1710c52818fSSatish Balay 1720c52818fSSatish Balay Collective on ls 1730c52818fSSatish Balay 1740c52818fSSatish Balay Input Parameters: 1750c52818fSSatish Balay . ls - the TaoLineSearch context 1760c52818fSSatish Balay 1770c52818fSSatish Balay Notes: 1780c52818fSSatish Balay The user will not need to explicitly call TaoLineSearchSetUp(), as it will 1790c52818fSSatish Balay automatically be called in TaoLineSearchSolve(). However, if the user 1800c52818fSSatish Balay desires to call it explicitly, it should come after TaoLineSearchCreate() 1810c52818fSSatish Balay but before TaoLineSearchApply(). 1820c52818fSSatish Balay 1830c52818fSSatish Balay Level: developer 1840c52818fSSatish Balay 1850c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 1860c52818fSSatish Balay @*/ 1870c52818fSSatish Balay 1880c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 1890c52818fSSatish Balay { 1908caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 1910c52818fSSatish Balay PetscBool flg; 192f06e3bfaSBarry Smith 1930c52818fSSatish Balay PetscFunctionBegin; 1940c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1950c52818fSSatish Balay if (ls->setupcalled) PetscFunctionReturn(0); 1960c52818fSSatish Balay if (!((PetscObject)ls)->type_name) { 1975f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchSetType(ls,default_type)); 1980c52818fSSatish Balay } 1990c52818fSSatish Balay if (ls->ops->setup) { 2005f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->setup)(ls)); 2010c52818fSSatish Balay } 2020c52818fSSatish Balay if (ls->usetaoroutines) { 2035f80ce2aSJacob Faibussowitsch CHKERRQ(TaoIsObjectiveDefined(ls->tao,&flg)); 2040c52818fSSatish Balay ls->hasobjective = flg; 2055f80ce2aSJacob Faibussowitsch CHKERRQ(TaoIsGradientDefined(ls->tao,&flg)); 2060c52818fSSatish Balay ls->hasgradient = flg; 2075f80ce2aSJacob Faibussowitsch CHKERRQ(TaoIsObjectiveAndGradientDefined(ls->tao,&flg)); 2080c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 2090c52818fSSatish Balay } else { 2100c52818fSSatish Balay if (ls->ops->computeobjective) { 2110c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 2120c52818fSSatish Balay } else { 2130c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 2140c52818fSSatish Balay } 2150c52818fSSatish Balay if (ls->ops->computegradient) { 2160c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 2170c52818fSSatish Balay } else { 2180c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 2190c52818fSSatish Balay } 2200c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 2210c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 2220c52818fSSatish Balay } else { 2230c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 2240c52818fSSatish Balay } 2250c52818fSSatish Balay } 2260c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 2270c52818fSSatish Balay PetscFunctionReturn(0); 2280c52818fSSatish Balay } 2290c52818fSSatish Balay 2300c52818fSSatish Balay /*@ 2310c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 2320c52818fSSatish Balay from one TaoLineSearchApply() to the next. This function resets this 2330c52818fSSatish Balay state information. 2340c52818fSSatish Balay 2350c52818fSSatish Balay Collective on TaoLineSearch 2360c52818fSSatish Balay 2370c52818fSSatish Balay Input Parameter: 2380c52818fSSatish Balay . ls - the TaoLineSearch context 2390c52818fSSatish Balay 2400c52818fSSatish Balay Level: developer 2410c52818fSSatish Balay 2420c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 2430c52818fSSatish Balay @*/ 2440c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 2450c52818fSSatish Balay { 2460c52818fSSatish Balay PetscFunctionBegin; 2470c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 2480c52818fSSatish Balay if (ls->ops->reset) { 2495f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->reset)(ls)); 2500c52818fSSatish Balay } 2510c52818fSSatish Balay PetscFunctionReturn(0); 2520c52818fSSatish Balay } 253f06e3bfaSBarry Smith 2540c52818fSSatish Balay /*@ 2550c52818fSSatish Balay TaoLineSearchDestroy - Destroys the TAO context that was created with 2560c52818fSSatish Balay TaoLineSearchCreate() 2570c52818fSSatish Balay 2580c52818fSSatish Balay Collective on TaoLineSearch 2590c52818fSSatish Balay 2607a7aea1fSJed Brown Input Parameter: 2610c52818fSSatish Balay . ls - the TaoLineSearch context 2620c52818fSSatish Balay 2630c52818fSSatish Balay Level: beginner 2640c52818fSSatish Balay 2650c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 2660c52818fSSatish Balay @*/ 2670c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 2680c52818fSSatish Balay { 2690c52818fSSatish Balay PetscFunctionBegin; 2700c52818fSSatish Balay if (!*ls) PetscFunctionReturn(0); 2710c52818fSSatish Balay PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 27283c8fe1dSLisandro Dalcin if (--((PetscObject)*ls)->refct > 0) {*ls = NULL; PetscFunctionReturn(0);} 2735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(*ls)->stepdirection)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(*ls)->start_x)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(*ls)->upper)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(*ls)->lower)); 2770c52818fSSatish Balay if ((*ls)->ops->destroy) { 2785f80ce2aSJacob Faibussowitsch CHKERRQ((*(*ls)->ops->destroy)(*ls)); 2790c52818fSSatish Balay } 2802a0dac07SAlp Dener if ((*ls)->usemonitor) { 2815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&(*ls)->viewer)); 2822a0dac07SAlp Dener } 2835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHeaderDestroy(ls)); 2840c52818fSSatish Balay PetscFunctionReturn(0); 2850c52818fSSatish Balay } 2860c52818fSSatish Balay 2870c52818fSSatish Balay /*@ 2880c52818fSSatish Balay TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 2890c52818fSSatish Balay 2900c52818fSSatish Balay Collective on TaoLineSearch 2910c52818fSSatish Balay 2920c52818fSSatish Balay Input Parameters: 293441846f8SBarry Smith + ls - the Tao context 2940c52818fSSatish Balay - s - search direction 2950c52818fSSatish Balay 2966b867d5aSJose E. Roman Input/Output Parameters: 2976b867d5aSJose E. Roman 2980c52818fSSatish Balay Output Parameters: 299f1a722f8SMatthew G. Knepley + x - On input the current solution, on output x contains the new solution determined by the line search 300f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution 301f1a722f8SMatthew G. Knepley . g - On input the gradient evaluated at x, on output contains the gradient at new solution 302f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x) 3030c52818fSSatish Balay - reason - reason why the line-search stopped 3040c52818fSSatish Balay 30597bb3fdcSJose E. Roman Notes: 3060c52818fSSatish Balay reason will be set to one of: 3070c52818fSSatish Balay 3080c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 3090c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 3100c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 3110c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 3120c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 3130c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 3140c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 3150c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 3160c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 3170c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search 3180c52818fSSatish Balay 3190c52818fSSatish Balay The algorithm developer must set up the TaoLineSearch with calls to 320441846f8SBarry Smith TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines() 3210c52818fSSatish Balay 3220c52818fSSatish Balay You may or may not need to follow this with a call to 3230c52818fSSatish Balay TaoAddLineSearchCounts(), depending on whether you want these 3240c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 3250c52818fSSatish Balay 3260c52818fSSatish Balay Level: beginner 3270c52818fSSatish Balay 3280c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts() 3290c52818fSSatish Balay @*/ 3300c52818fSSatish Balay 331e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 3320c52818fSSatish Balay { 3330c52818fSSatish Balay PetscInt low1,low2,low3,high1,high2,high3; 3340c52818fSSatish Balay 3350c52818fSSatish Balay PetscFunctionBegin; 3360c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 3370c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3383f6ba705SLisandro Dalcin PetscValidRealPointer(f,3); 3390c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 3400c52818fSSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,5); 3410c52818fSSatish Balay PetscValidPointer(reason,7); 3420c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 3430c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,g,4); 3440c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,s,5); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x, &low1, &high1)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(g, &low2, &high2)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(s, &low3, &high3)); 3483c859ba3SBarry Smith PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible vector local lengths"); 3490c52818fSSatish Balay 35097ab8e72SStefano Zampini *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)s)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ls->stepdirection)); 353050fc7a3SBarry Smith ls->stepdirection = s; 3540c52818fSSatish Balay 3555f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchSetUp(ls)); 3563c859ba3SBarry Smith PetscCheck(ls->ops->apply,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine"); 3570c52818fSSatish Balay ls->nfeval = 0; 3580c52818fSSatish Balay ls->ngeval = 0; 3590c52818fSSatish Balay ls->nfgeval = 0; 3600c52818fSSatish Balay /* Check parameter values */ 3610c52818fSSatish Balay if (ls->ftol < 0.0) { 3625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol)); 3630c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3640c52818fSSatish Balay } 3650c52818fSSatish Balay if (ls->rtol < 0.0) { 3665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol)); 3670c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3680c52818fSSatish Balay } 3690c52818fSSatish Balay if (ls->gtol < 0.0) { 3705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol)); 3710c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3720c52818fSSatish Balay } 3730c52818fSSatish Balay if (ls->stepmin < 0.0) { 3745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin)); 3750c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3760c52818fSSatish Balay } 3770c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 3785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax)); 3790c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3800c52818fSSatish Balay } 3810c52818fSSatish Balay if (ls->max_funcs < 0) { 3825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n",ls->max_funcs)); 3830c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3840c52818fSSatish Balay } 3850c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 3865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f)); 3870c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_INFORNAN; 3880c52818fSSatish Balay } 3890c52818fSSatish Balay 3905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)x)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ls->start_x)); 3920c52818fSSatish Balay ls->start_x = x; 393050fc7a3SBarry Smith 3945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->apply)(ls,x,f,g,s)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0)); 3970c52818fSSatish Balay *reason = ls->reason; 3980c52818fSSatish Balay ls->new_f = *f; 3990c52818fSSatish Balay 40097ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 4010c52818fSSatish Balay 4025f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view")); 4030c52818fSSatish Balay PetscFunctionReturn(0); 4040c52818fSSatish Balay } 4050c52818fSSatish Balay 4060c52818fSSatish Balay /*@C 4070c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 4080c52818fSSatish Balay 4090c52818fSSatish Balay Collective on TaoLineSearch 4100c52818fSSatish Balay 4110c52818fSSatish Balay Input Parameters: 4120c52818fSSatish Balay + ls - the TaoLineSearch context 413820824deSAlp Dener - type - the TaoLineSearchType selection 4140c52818fSSatish Balay 4150c52818fSSatish Balay Available methods include: 416820824deSAlp Dener + more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition 417820824deSAlp Dener . armijo - simple backtracking line search enforcing only the sufficient decrease condition 418820824deSAlp Dener - unit - do not perform a line search and always accept unit step length 4190c52818fSSatish Balay 4200c52818fSSatish Balay Options Database Keys: 421820824deSAlp Dener . -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime 4220c52818fSSatish Balay 4230c52818fSSatish Balay Level: beginner 4240c52818fSSatish Balay 4250c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply() 4260c52818fSSatish Balay 4270c52818fSSatish Balay @*/ 4280c52818fSSatish Balay 429dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 4300c52818fSSatish Balay { 4310c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 4320c52818fSSatish Balay PetscBool flg; 4330c52818fSSatish Balay 4340c52818fSSatish Balay PetscFunctionBegin; 4350c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4360c52818fSSatish Balay PetscValidCharPointer(type,2); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)ls, type, &flg)); 4380c52818fSSatish Balay if (flg) PetscFunctionReturn(0); 4390c52818fSSatish Balay 4405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r)); 4413c859ba3SBarry Smith PetscCheck(r,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 4420c52818fSSatish Balay if (ls->ops->destroy) { 4435f80ce2aSJacob Faibussowitsch CHKERRQ((*(ls)->ops->destroy)(ls)); 4440c52818fSSatish Balay } 4450c52818fSSatish Balay ls->max_funcs = 30; 4460c52818fSSatish Balay ls->ftol = 0.0001; 4470c52818fSSatish Balay ls->gtol = 0.9; 4486f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 4496f4723b1SBarry Smith ls->rtol = 1.0e-5; 4506f4723b1SBarry Smith #else 4510c52818fSSatish Balay ls->rtol = 1.0e-10; 4526f4723b1SBarry Smith #endif 4530c52818fSSatish Balay ls->stepmin = 1.0e-20; 4540c52818fSSatish Balay ls->stepmax = 1.0e+20; 4550c52818fSSatish Balay 4560c52818fSSatish Balay ls->nfeval = 0; 4570c52818fSSatish Balay ls->ngeval = 0; 4580c52818fSSatish Balay ls->nfgeval = 0; 45983c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 46083c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 46183c8fe1dSLisandro Dalcin ls->ops->view = NULL; 46283c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 46383c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 4640c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4655f80ce2aSJacob Faibussowitsch CHKERRQ((*r)(ls)); 4665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)ls, type)); 4670c52818fSSatish Balay PetscFunctionReturn(0); 4680c52818fSSatish Balay } 4690c52818fSSatish Balay 4702a0dac07SAlp Dener /*@C 4712a0dac07SAlp Dener TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the 4722a0dac07SAlp Dener iteration number, step length, and function value before calling the implementation 4732a0dac07SAlp Dener specific monitor. 4742a0dac07SAlp Dener 4752a0dac07SAlp Dener Input Parameters: 4762a0dac07SAlp Dener + ls - the TaoLineSearch context 4772a0dac07SAlp Dener . its - the current iterate number (>=0) 4782a0dac07SAlp Dener . f - the current objective function value 4792a0dac07SAlp Dener - step - the step length 4802a0dac07SAlp Dener 4812a0dac07SAlp Dener Options Database Key: 4822a0dac07SAlp Dener . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output 4832a0dac07SAlp Dener 4842a0dac07SAlp Dener Level: developer 4852a0dac07SAlp Dener 4862a0dac07SAlp Dener @*/ 4872a0dac07SAlp Dener PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step) 4882a0dac07SAlp Dener { 4892a0dac07SAlp Dener PetscInt tabs; 4902a0dac07SAlp Dener 4912a0dac07SAlp Dener PetscFunctionBegin; 4922a0dac07SAlp Dener PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4932a0dac07SAlp Dener if (ls->usemonitor) { 4945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIGetTab(ls->viewer, &tabs)); 4955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel)); 4965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(ls->viewer, "%3D LS", its)); 4975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f)); 4985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step)); 4992a0dac07SAlp Dener if (ls->ops->monitor && its > 0) { 5005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3)); 5015f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->monitor)(ls)); 5022a0dac07SAlp Dener } 5035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISetTab(ls->viewer, tabs)); 5042a0dac07SAlp Dener } 5052a0dac07SAlp Dener PetscFunctionReturn(0); 5062a0dac07SAlp Dener } 5072a0dac07SAlp Dener 5080c52818fSSatish Balay /*@ 5090c52818fSSatish Balay TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 5100c52818fSSatish Balay options. 5110c52818fSSatish Balay 5120c52818fSSatish Balay Collective on TaoLineSearch 5130c52818fSSatish Balay 51401d2d390SJose E. Roman Input Parameter: 5150c52818fSSatish Balay . ls - the TaoLineSearch context 5160c52818fSSatish Balay 5170c52818fSSatish Balay Options Database Keys: 5180c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 5190c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 5200c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 5210c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 5220c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 5230c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 5240c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 5250c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 5260c52818fSSatish Balay 5270c52818fSSatish Balay Level: beginner 5280c52818fSSatish Balay @*/ 5290c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 5300c52818fSSatish Balay { 5310c52818fSSatish Balay PetscErrorCode ierr; 5328caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 5332a0dac07SAlp Dener char type[256],monfilename[PETSC_MAX_PATH_LEN]; 5342a0dac07SAlp Dener PetscViewer monviewer; 5350c52818fSSatish Balay PetscBool flg; 536f06e3bfaSBarry Smith 5370c52818fSSatish Balay PetscFunctionBegin; 5380c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5390c52818fSSatish Balay ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr); 5400c52818fSSatish Balay if (((PetscObject)ls)->type_name) { 5410c52818fSSatish Balay default_type = ((PetscObject)ls)->type_name; 5420c52818fSSatish Balay } 5430c52818fSSatish Balay /* Check for type from options */ 5445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg)); 5450c52818fSSatish Balay if (flg) { 5465f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchSetType(ls,type)); 5470c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 5485f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchSetType(ls,default_type)); 5490c52818fSSatish Balay } 5500c52818fSSatish Balay 5515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL)); 5525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL)); 5545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL)); 5555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL)); 5565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL)); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,sizeof(monfilename),&flg)); 5582a0dac07SAlp Dener if (flg) { 5595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer)); 5602a0dac07SAlp Dener ls->viewer = monviewer; 5612a0dac07SAlp Dener ls->usemonitor = PETSC_TRUE; 5622a0dac07SAlp Dener } 5630c52818fSSatish Balay if (ls->ops->setfromoptions) { 5645f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->setfromoptions)(PetscOptionsObject,ls)); 5650c52818fSSatish Balay } 5660c52818fSSatish Balay ierr = PetscOptionsEnd();CHKERRQ(ierr); 5670c52818fSSatish Balay PetscFunctionReturn(0); 5680c52818fSSatish Balay } 5690c52818fSSatish Balay 5700c52818fSSatish Balay /*@C 5710c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5720c52818fSSatish Balay 5730c52818fSSatish Balay Not Collective 5740c52818fSSatish Balay 5750c52818fSSatish Balay Input Parameter: 5760c52818fSSatish Balay . ls - the TaoLineSearch context 5770c52818fSSatish Balay 578fd292e60Sprj- Output Parameter: 5790c52818fSSatish Balay . type - the line search algorithm in effect 5800c52818fSSatish Balay 5810c52818fSSatish Balay Level: developer 5820c52818fSSatish Balay 5830c52818fSSatish Balay @*/ 584dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 5850c52818fSSatish Balay { 5860c52818fSSatish Balay PetscFunctionBegin; 5870c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5880c52818fSSatish Balay PetscValidPointer(type,2); 5890c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 5900c52818fSSatish Balay PetscFunctionReturn(0); 5910c52818fSSatish Balay } 5920c52818fSSatish Balay 5930c52818fSSatish Balay /*@ 5940c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 5950c52818fSSatish Balay routines used by the line search in last application (not cumulative). 5960c52818fSSatish Balay 5970c52818fSSatish Balay Not Collective 5980c52818fSSatish Balay 5990c52818fSSatish Balay Input Parameter: 6000c52818fSSatish Balay . ls - the TaoLineSearch context 6010c52818fSSatish Balay 6020c52818fSSatish Balay Output Parameters: 6030c52818fSSatish Balay + nfeval - number of function evaluations 6040c52818fSSatish Balay . ngeval - number of gradient evaluations 6050c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 6060c52818fSSatish Balay 6070c52818fSSatish Balay Level: intermediate 6080c52818fSSatish Balay 6090c52818fSSatish Balay Note: 610441846f8SBarry Smith If the line search is using the Tao objective and gradient 611441846f8SBarry Smith routines directly (see TaoLineSearchUseTaoRoutines()), then TAO 6120c52818fSSatish Balay is already counting the number of evaluations. 6130c52818fSSatish Balay 6140c52818fSSatish Balay @*/ 6150c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 6160c52818fSSatish Balay { 6170c52818fSSatish Balay PetscFunctionBegin; 6180c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6190c52818fSSatish Balay *nfeval = ls->nfeval; 6200c52818fSSatish Balay *ngeval = ls->ngeval; 6210c52818fSSatish Balay *nfgeval = ls->nfgeval; 6220c52818fSSatish Balay PetscFunctionReturn(0); 6230c52818fSSatish Balay } 6240c52818fSSatish Balay 6250c52818fSSatish Balay /*@ 626441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 627441846f8SBarry Smith Tao evaluation routines. 6280c52818fSSatish Balay 6290c52818fSSatish Balay Not Collective 6300c52818fSSatish Balay 6310c52818fSSatish Balay Input Parameter: 6320c52818fSSatish Balay . ls - the TaoLineSearch context 6330c52818fSSatish Balay 6340c52818fSSatish Balay Output Parameter: 635441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines, 6360c52818fSSatish Balay otherwise PETSC_FALSE 6370c52818fSSatish Balay 6380c52818fSSatish Balay Level: developer 6390c52818fSSatish Balay @*/ 640441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 6410c52818fSSatish Balay { 6420c52818fSSatish Balay PetscFunctionBegin; 6430c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 644f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 6450c52818fSSatish Balay PetscFunctionReturn(0); 6460c52818fSSatish Balay } 6470c52818fSSatish Balay 6480c52818fSSatish Balay /*@C 6490c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 6500c52818fSSatish Balay 6510c52818fSSatish Balay Logically Collective on TaoLineSearch 6520c52818fSSatish Balay 653d8d19677SJose E. Roman Input Parameters: 6540c52818fSSatish Balay + ls - the TaoLineSearch context 6550c52818fSSatish Balay . func - the objective function evaluation routine 6560c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6570c52818fSSatish Balay 6580c52818fSSatish Balay Calling sequence of func: 6590c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 6600c52818fSSatish Balay 6610c52818fSSatish Balay + x - input vector 6620c52818fSSatish Balay . f - function value 6630c52818fSSatish Balay - ctx (optional) user-defined context 6640c52818fSSatish Balay 6650c52818fSSatish Balay Level: beginner 6660c52818fSSatish Balay 6670c52818fSSatish Balay Note: 6680c52818fSSatish Balay Use this routine only if you want the line search objective 669441846f8SBarry Smith evaluation routine to be different from the Tao's objective 6700c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6710c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6720c52818fSSatish Balay 6730c52818fSSatish Balay Note: 6740c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6750c52818fSSatish Balay line search, application programmers should be wary of overriding the 6760c52818fSSatish Balay default objective routine. 6770c52818fSSatish Balay 678441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 6790c52818fSSatish Balay @*/ 6800c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx) 6810c52818fSSatish Balay { 6820c52818fSSatish Balay PetscFunctionBegin; 6830c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6840c52818fSSatish Balay 6850c52818fSSatish Balay ls->ops->computeobjective=func; 6860c52818fSSatish Balay if (ctx) ls->userctx_func=ctx; 6870c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 6880c52818fSSatish Balay PetscFunctionReturn(0); 6890c52818fSSatish Balay } 6900c52818fSSatish Balay 6910c52818fSSatish Balay /*@C 6920c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 6930c52818fSSatish Balay 6940c52818fSSatish Balay Logically Collective on TaoLineSearch 6950c52818fSSatish Balay 696d8d19677SJose E. Roman Input Parameters: 6970c52818fSSatish Balay + ls - the TaoLineSearch context 6980c52818fSSatish Balay . func - the gradient evaluation routine 6990c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7000c52818fSSatish Balay 7010c52818fSSatish Balay Calling sequence of func: 7020c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, Vec g, void *ctx); 7030c52818fSSatish Balay 7040c52818fSSatish Balay + x - input vector 7050c52818fSSatish Balay . g - gradient vector 7060c52818fSSatish Balay - ctx (optional) user-defined context 7070c52818fSSatish Balay 7080c52818fSSatish Balay Level: beginner 7090c52818fSSatish Balay 7100c52818fSSatish Balay Note: 7110c52818fSSatish Balay Use this routine only if you want the line search gradient 712441846f8SBarry Smith evaluation routine to be different from the Tao's gradient 7130c52818fSSatish Balay evaluation routine. If you use this routine you must also set 7140c52818fSSatish Balay the line search function and/or function/gradient routine. 7150c52818fSSatish Balay 7160c52818fSSatish Balay Note: 7170c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 7180c52818fSSatish Balay line search, application programmers should be wary of overriding the 7190c52818fSSatish Balay default gradient routine. 7200c52818fSSatish Balay 721441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 7220c52818fSSatish Balay @*/ 7230c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx) 7240c52818fSSatish Balay { 7250c52818fSSatish Balay PetscFunctionBegin; 7260c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7270c52818fSSatish Balay ls->ops->computegradient=func; 7280c52818fSSatish Balay if (ctx) ls->userctx_grad=ctx; 7290c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 7300c52818fSSatish Balay PetscFunctionReturn(0); 7310c52818fSSatish Balay } 7320c52818fSSatish Balay 7330c52818fSSatish Balay /*@C 7340c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 7350c52818fSSatish Balay 7360c52818fSSatish Balay Logically Collective on TaoLineSearch 7370c52818fSSatish Balay 738d8d19677SJose E. Roman Input Parameters: 7390c52818fSSatish Balay + ls - the TaoLineSearch context 7400c52818fSSatish Balay . func - the objective and gradient evaluation routine 7410c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7420c52818fSSatish Balay 7430c52818fSSatish Balay Calling sequence of func: 7440c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 7450c52818fSSatish Balay 7460c52818fSSatish Balay + x - input vector 7470c52818fSSatish Balay . f - function value 7480c52818fSSatish Balay . g - gradient vector 7490c52818fSSatish Balay - ctx (optional) user-defined context 7500c52818fSSatish Balay 7510c52818fSSatish Balay Level: beginner 7520c52818fSSatish Balay 7530c52818fSSatish Balay Note: 7540c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 755441846f8SBarry Smith evaluation routines to be different from the Tao's objective 7560c52818fSSatish Balay and gradient evaluation routines. 7570c52818fSSatish Balay 7580c52818fSSatish Balay Note: 7590c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7600c52818fSSatish Balay line search, application programmers should be wary of overriding the 7610c52818fSSatish Balay default objective routine. 7620c52818fSSatish Balay 763441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines() 7640c52818fSSatish Balay @*/ 7650c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx) 7660c52818fSSatish Balay { 7670c52818fSSatish Balay PetscFunctionBegin; 7680c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7690c52818fSSatish Balay ls->ops->computeobjectiveandgradient=func; 7700c52818fSSatish Balay if (ctx) ls->userctx_funcgrad=ctx; 7710c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7720c52818fSSatish Balay PetscFunctionReturn(0); 7730c52818fSSatish Balay } 7740c52818fSSatish Balay 7750c52818fSSatish Balay /*@C 7760c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 7770c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 7780c52818fSSatish Balay Sometimes it is more efficient to compute the inner product of the gradient 7790c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 7800c52818fSSatish Balay the line search typically needs of the gradient. 7810c52818fSSatish Balay 7820c52818fSSatish Balay Logically Collective on TaoLineSearch 7830c52818fSSatish Balay 784d8d19677SJose E. Roman Input Parameters: 7850c52818fSSatish Balay + ls - the TaoLineSearch context 7860c52818fSSatish Balay . func - the objective and gradient evaluation routine 7870c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7880c52818fSSatish Balay 7890c52818fSSatish Balay Calling sequence of func: 7900c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 7910c52818fSSatish Balay 7920c52818fSSatish Balay + x - input vector 7930c52818fSSatish Balay . s - step direction 7940c52818fSSatish Balay . f - function value 7950c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 7960c52818fSSatish Balay - ctx (optional) user-defined context 7970c52818fSSatish Balay 7980c52818fSSatish Balay Note: The gradient will still need to be computed at the end of the line 7990c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 8000c52818fSSatish Balay routine 8010c52818fSSatish Balay 8020c52818fSSatish Balay Note: Bounded line searches (those used in bounded optimization algorithms) 8030c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 8040c52818fSSatish Balay x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength() 8050c52818fSSatish Balay 8060c52818fSSatish Balay Level: advanced 8070c52818fSSatish Balay 8080c52818fSSatish Balay Note: 8090c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 8100c52818fSSatish Balay line search, application programmers should be wary of overriding the 8110c52818fSSatish Balay default objective routine. 8120c52818fSSatish Balay 813441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines() 8140c52818fSSatish Balay @*/ 8150c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx) 8160c52818fSSatish Balay { 8170c52818fSSatish Balay PetscFunctionBegin; 8180c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8190c52818fSSatish Balay ls->ops->computeobjectiveandgts=func; 8200c52818fSSatish Balay if (ctx) ls->userctx_funcgts=ctx; 8210c52818fSSatish Balay ls->usegts = PETSC_TRUE; 8220c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 8230c52818fSSatish Balay PetscFunctionReturn(0); 8240c52818fSSatish Balay } 8250c52818fSSatish Balay 8260c52818fSSatish Balay /*@C 827441846f8SBarry Smith TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the 828441846f8SBarry Smith objective and gradient evaluation routines from the given Tao object. 8290c52818fSSatish Balay 8300c52818fSSatish Balay Logically Collective on TaoLineSearch 8310c52818fSSatish Balay 832d8d19677SJose E. Roman Input Parameters: 8330c52818fSSatish Balay + ls - the TaoLineSearch context 834441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines 8350c52818fSSatish Balay 8360c52818fSSatish Balay Level: developer 8370c52818fSSatish Balay 8380c52818fSSatish Balay .seealso: TaoLineSearchCreate() 8390c52818fSSatish Balay @*/ 840441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts) 8410c52818fSSatish Balay { 8420c52818fSSatish Balay PetscFunctionBegin; 8430c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 844064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(ts,TAO_CLASSID,2); 845441846f8SBarry Smith ls->tao = ts; 8460c52818fSSatish Balay ls->usetaoroutines = PETSC_TRUE; 8470c52818fSSatish Balay PetscFunctionReturn(0); 8480c52818fSSatish Balay } 8490c52818fSSatish Balay 8500c52818fSSatish Balay /*@ 8510c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 8520c52818fSSatish Balay 8530c52818fSSatish Balay Collective on TaoLineSearch 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 86295452b02SPatrick Sanan Notes: 86395452b02SPatrick Sanan TaoLineSearchComputeObjective() is typically used within line searches 8640c52818fSSatish Balay so most users would not generally call this routine themselves. 8650c52818fSSatish Balay 8660c52818fSSatish Balay Level: developer 8670c52818fSSatish Balay 8680c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 8690c52818fSSatish Balay @*/ 8700c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 8710c52818fSSatish Balay { 8720c52818fSSatish Balay Vec gdummy; 8730c52818fSSatish Balay PetscReal gts; 874f06e3bfaSBarry Smith 8750c52818fSSatish Balay PetscFunctionBegin; 8760c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8770c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 878*dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 8790c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8800c52818fSSatish Balay if (ls->usetaoroutines) { 8815f80ce2aSJacob Faibussowitsch CHKERRQ(TaoComputeObjective(ls->tao,x,f)); 8820c52818fSSatish Balay } else { 8833c859ba3SBarry Smith PetscCheck(ls->ops->computeobjective || ls->ops->computeobjectiveandgradient || ls->ops->computeobjectiveandgts,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 8845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 8850c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective routine"); 8860c52818fSSatish Balay if (ls->ops->computeobjective) { 8875f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjective)(ls,x,f,ls->userctx_func)); 8880c52818fSSatish Balay } else if (ls->ops->computeobjectiveandgradient) { 8895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&gdummy)); 8905f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad)); 8915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&gdummy)); 8920c52818fSSatish Balay } else { 8935f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts)); 8940c52818fSSatish Balay } 8950c52818fSSatish Balay PetscStackPop; 8965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 8970c52818fSSatish Balay } 8980c52818fSSatish Balay ls->nfeval++; 8990c52818fSSatish Balay PetscFunctionReturn(0); 9000c52818fSSatish Balay } 9010c52818fSSatish Balay 9020c52818fSSatish Balay /*@ 9030c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 9040c52818fSSatish Balay 905441846f8SBarry Smith Collective on Tao 9060c52818fSSatish Balay 9070c52818fSSatish Balay Input Parameters: 9080c52818fSSatish Balay + ls - the TaoLineSearch context 9090c52818fSSatish Balay - x - input vector 9100c52818fSSatish Balay 911d8d19677SJose E. Roman Output Parameters: 9120c52818fSSatish Balay + f - Objective value at X 9130c52818fSSatish Balay - g - Gradient vector at X 9140c52818fSSatish Balay 91595452b02SPatrick Sanan Notes: 91695452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGradient() 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: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 9220c52818fSSatish Balay @*/ 9230c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 9240c52818fSSatish Balay { 9250c52818fSSatish Balay PetscFunctionBegin; 9260c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9270c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 928*dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 9290c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 9300c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9310c52818fSSatish Balay PetscCheckSameComm(ls,1,g,4); 9320c52818fSSatish Balay if (ls->usetaoroutines) { 9335f80ce2aSJacob Faibussowitsch CHKERRQ(TaoComputeObjectiveAndGradient(ls->tao,x,f,g)); 9340c52818fSSatish Balay } else { 9355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 9360c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 937362febeeSStefano Zampini PetscStackPush("TaoLineSearch user objective/gradient routine"); 9385f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad)); 9390c52818fSSatish Balay PetscStackPop; 940362febeeSStefano Zampini } else { 9413c859ba3SBarry Smith PetscCheck(ls->ops->computeobjective,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 942362febeeSStefano Zampini PetscStackPush("TaoLineSearch user objective routine"); 9435f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjective)(ls,x,f,ls->userctx_func)); 944362febeeSStefano Zampini PetscStackPop; 9453c859ba3SBarry Smith PetscCheck(ls->ops->computegradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 946362febeeSStefano Zampini PetscStackPush("TaoLineSearch user gradient routine"); 9475f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computegradient)(ls,x,g,ls->userctx_grad)); 948362febeeSStefano Zampini PetscStackPop; 949362febeeSStefano Zampini } 9505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 9515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f))); 952f06e3bfaSBarry Smith } 953fbe0838dSJason Sarich ls->nfgeval++; 9540c52818fSSatish Balay PetscFunctionReturn(0); 9550c52818fSSatish Balay } 9560c52818fSSatish Balay 9570c52818fSSatish Balay /*@ 9580c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 9590c52818fSSatish Balay 9600c52818fSSatish Balay Collective on TaoLineSearch 9610c52818fSSatish Balay 9620c52818fSSatish Balay Input Parameters: 9630c52818fSSatish Balay + ls - the TaoLineSearch context 9640c52818fSSatish Balay - x - input vector 9650c52818fSSatish Balay 9660c52818fSSatish Balay Output Parameter: 9670c52818fSSatish Balay . g - gradient vector 9680c52818fSSatish Balay 96995452b02SPatrick Sanan Notes: 97095452b02SPatrick Sanan TaoComputeGradient() is typically used within line searches 9710c52818fSSatish Balay so most users would not generally call this routine themselves. 9720c52818fSSatish Balay 9730c52818fSSatish Balay Level: developer 9740c52818fSSatish Balay 9750c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 9760c52818fSSatish Balay @*/ 9770c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 9780c52818fSSatish Balay { 9790c52818fSSatish Balay PetscReal fdummy; 980f06e3bfaSBarry Smith 9810c52818fSSatish Balay PetscFunctionBegin; 9820c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9830c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9840c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,3); 9850c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9860c52818fSSatish Balay PetscCheckSameComm(ls,1,g,3); 9870c52818fSSatish Balay if (ls->usetaoroutines) { 9885f80ce2aSJacob Faibussowitsch CHKERRQ(TaoComputeGradient(ls->tao,x,g)); 9890c52818fSSatish Balay } else { 9905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 9910c52818fSSatish Balay PetscStackPush("TaoLineSearch user gradient routine"); 9920c52818fSSatish Balay if (ls->ops->computegradient) { 9935f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computegradient)(ls,x,g,ls->userctx_grad)); 9940c52818fSSatish Balay } else { 9953c859ba3SBarry Smith PetscCheck(ls->ops->computeobjectiveandgradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 9965f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad)); 9970c52818fSSatish Balay } 9980c52818fSSatish Balay PetscStackPop; 9995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 10000c52818fSSatish Balay } 10010c52818fSSatish Balay ls->ngeval++; 10020c52818fSSatish Balay PetscFunctionReturn(0); 10030c52818fSSatish Balay } 10040c52818fSSatish Balay 10050c52818fSSatish Balay /*@ 10060c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 10070c52818fSSatish Balay 1008441846f8SBarry Smith Collective on Tao 10090c52818fSSatish Balay 10100c52818fSSatish Balay Input Parameters: 10110c52818fSSatish Balay + ls - the TaoLineSearch context 10120c52818fSSatish Balay - x - input vector 10130c52818fSSatish Balay 1014d8d19677SJose E. Roman Output Parameters: 10150c52818fSSatish Balay + f - Objective value at X 10160c52818fSSatish Balay - gts - inner product of gradient and step direction at X 10170c52818fSSatish Balay 101895452b02SPatrick Sanan Notes: 101995452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 10200c52818fSSatish Balay so most users would not generally call this routine themselves. 10210c52818fSSatish Balay 10220c52818fSSatish Balay Level: developer 10230c52818fSSatish Balay 10240c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 10250c52818fSSatish Balay @*/ 10260c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 10270c52818fSSatish Balay { 10280c52818fSSatish Balay PetscFunctionBegin; 10290c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10300c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1031*dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 1032*dadcf809SJacob Faibussowitsch PetscValidRealPointer(gts,4); 10330c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 10343c859ba3SBarry Smith PetscCheck(ls->ops->computeobjectiveandgts,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 10355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0)); 10360c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gts routine"); 10375f80ce2aSJacob Faibussowitsch CHKERRQ((*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts)); 10380c52818fSSatish Balay PetscStackPop; 10395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0)); 10405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f))); 10410c52818fSSatish Balay ls->nfeval++; 10420c52818fSSatish Balay PetscFunctionReturn(0); 10430c52818fSSatish Balay } 10440c52818fSSatish Balay 10450c52818fSSatish Balay /*@ 10460c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 10470c52818fSSatish Balay 10480c52818fSSatish Balay Collective on TaoLineSearch 10490c52818fSSatish Balay 10500c52818fSSatish Balay Input Parameter: 10510c52818fSSatish Balay . ls - the TaoLineSearch context 10520c52818fSSatish Balay 1053d8d19677SJose E. Roman Output Parameters: 10540c52818fSSatish Balay + x - the new solution 10550c52818fSSatish Balay . f - the objective function value at x 10560c52818fSSatish Balay . g - the gradient at x 10570c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 10580c52818fSSatish Balay - reason - the reason why the line search terminated 10590c52818fSSatish Balay 10600c52818fSSatish Balay reason will be set to one of: 10610c52818fSSatish Balay 10620c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 10630c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 10640c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 10650c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 10660c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 10670c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 10680c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 10690c52818fSSatish Balay 10700c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 10710c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 10720c52818fSSatish Balay 1073a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search 10740c52818fSSatish Balay 10750c52818fSSatish Balay Level: developer 10760c52818fSSatish Balay 10770c52818fSSatish Balay @*/ 1078e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 10790c52818fSSatish Balay { 10800c52818fSSatish Balay PetscFunctionBegin; 10810c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10820c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1083*dadcf809SJacob Faibussowitsch PetscValidRealPointer(f,3); 10840c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 10850c52818fSSatish Balay PetscValidIntPointer(reason,6); 10860c52818fSSatish Balay if (ls->new_x) { 10875f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(ls->new_x,x)); 10880c52818fSSatish Balay } 10890c52818fSSatish Balay *f = ls->new_f; 10900c52818fSSatish Balay if (ls->new_g) { 10915f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(ls->new_g,g)); 10920c52818fSSatish Balay } 109397ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 10940c52818fSSatish Balay *reason = ls->reason; 10950c52818fSSatish Balay PetscFunctionReturn(0); 10960c52818fSSatish Balay } 10970c52818fSSatish Balay 10980c52818fSSatish Balay /*@ 10990c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 11000c52818fSSatish Balay search. 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 . x - The initial point of the line search 11090c52818fSSatish Balay 11100c52818fSSatish Balay Level: intermediate 11110c52818fSSatish Balay @*/ 11120c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 11130c52818fSSatish Balay { 11140c52818fSSatish Balay PetscFunctionBegin; 11150c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 111697ab8e72SStefano Zampini if (x) *x = ls->start_x; 11170c52818fSSatish Balay PetscFunctionReturn(0); 11180c52818fSSatish Balay } 11190c52818fSSatish Balay 11200c52818fSSatish Balay /*@ 11210c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 11220c52818fSSatish Balay search. 11230c52818fSSatish Balay 11240c52818fSSatish Balay Not Collective 11250c52818fSSatish Balay 11260c52818fSSatish Balay Input Parameter: 11270c52818fSSatish Balay . ls - the TaoLineSearch context 11280c52818fSSatish Balay 11290c52818fSSatish Balay Output Parameter: 11300c52818fSSatish Balay . s - the step direction of the line search 11310c52818fSSatish Balay 11320c52818fSSatish Balay Level: advanced 11330c52818fSSatish Balay @*/ 11340c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 11350c52818fSSatish Balay { 11360c52818fSSatish Balay PetscFunctionBegin; 11370c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 113897ab8e72SStefano Zampini if (s) *s = ls->stepdirection; 11390c52818fSSatish Balay PetscFunctionReturn(0); 11400c52818fSSatish Balay } 11410c52818fSSatish Balay 11420c52818fSSatish Balay /*@ 11430c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 11440c52818fSSatish Balay 11450c52818fSSatish Balay Not Collective 11460c52818fSSatish Balay 11470c52818fSSatish Balay Input Parameter: 11480c52818fSSatish Balay . ls - the TaoLineSearch context 11490c52818fSSatish Balay 11500c52818fSSatish Balay Output Parameter: 11510c52818fSSatish Balay . f - the objective value at the full step length 11520c52818fSSatish Balay 11530c52818fSSatish Balay Level: developer 11540c52818fSSatish Balay @*/ 11550c52818fSSatish Balay 11560c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 11570c52818fSSatish Balay { 11580c52818fSSatish Balay PetscFunctionBegin; 11590c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11600c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 11610c52818fSSatish Balay PetscFunctionReturn(0); 11620c52818fSSatish Balay } 11630c52818fSSatish Balay 11640c52818fSSatish Balay /*@ 11650c52818fSSatish Balay TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 11660c52818fSSatish Balay 1167441846f8SBarry Smith Logically Collective on Tao 11680c52818fSSatish Balay 11690c52818fSSatish Balay Input Parameters: 11700c52818fSSatish Balay + ls - the TaoLineSearch context 11710c52818fSSatish Balay . xl - vector of lower bounds 11720c52818fSSatish Balay - xu - vector of upper bounds 11730c52818fSSatish Balay 11740c52818fSSatish Balay Note: If the variable bounds are not set with this routine, then 1175e270355aSBarry Smith PETSC_NINFINITY and PETSC_INFINITY are assumed 11760c52818fSSatish Balay 11770c52818fSSatish Balay Level: beginner 11780c52818fSSatish Balay 11790c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 11800c52818fSSatish Balay @*/ 11810c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 11820c52818fSSatish Balay { 11830c52818fSSatish Balay PetscFunctionBegin; 11840c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11850c52818fSSatish Balay PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 11860c52818fSSatish Balay PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 11875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)xl)); 11885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)xu)); 11895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ls->lower)); 11905f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ls->upper)); 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 12277a7aea1fSJed Brown 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 PetscFunctionBegin; 12705f80ce2aSJacob Faibussowitsch CHKERRQ(TaoLineSearchInitializePackage()); 12715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func)); 12720c52818fSSatish Balay PetscFunctionReturn(0); 12730c52818fSSatish Balay } 12740c52818fSSatish Balay 12750c52818fSSatish Balay /*@C 12760c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 12770c52818fSSatish Balay for all TaoLineSearch options in the database. 12780c52818fSSatish Balay 12790c52818fSSatish Balay Collective on TaoLineSearch 12800c52818fSSatish Balay 12810c52818fSSatish Balay Input Parameters: 12820c52818fSSatish Balay + ls - the TaoLineSearch solver context 12830c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 12840c52818fSSatish Balay 12850c52818fSSatish Balay Notes: 12860c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 12870c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 12880c52818fSSatish Balay 12890c52818fSSatish Balay Level: advanced 12900c52818fSSatish Balay 12910c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 12920c52818fSSatish Balay @*/ 12930c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 12940c52818fSSatish Balay { 1295f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 12960c52818fSSatish Balay } 12970c52818fSSatish Balay 12980c52818fSSatish Balay /*@C 12990c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 13000c52818fSSatish Balay TaoLineSearch options in the database 13010c52818fSSatish Balay 13020c52818fSSatish Balay Not Collective 13030c52818fSSatish Balay 13040c52818fSSatish Balay Input Parameters: 13050c52818fSSatish Balay . ls - the TaoLineSearch context 13060c52818fSSatish Balay 13070c52818fSSatish Balay Output Parameters: 13080c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 13090c52818fSSatish Balay 131095452b02SPatrick Sanan Notes: 131195452b02SPatrick Sanan On the fortran side, the user should pass in a string 'prefix' of 13120c52818fSSatish Balay sufficient length to hold the prefix. 13130c52818fSSatish Balay 13140c52818fSSatish Balay Level: advanced 13150c52818fSSatish Balay 13160c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 13170c52818fSSatish Balay @*/ 13180c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 13190c52818fSSatish Balay { 1320f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 13210c52818fSSatish Balay } 13220c52818fSSatish Balay 13230c52818fSSatish Balay /*@C 13240c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 13250c52818fSSatish Balay TaoLineSearch options in the database. 13260c52818fSSatish Balay 13270c52818fSSatish Balay Logically Collective on TaoLineSearch 13280c52818fSSatish Balay 13290c52818fSSatish Balay Input Parameters: 13300c52818fSSatish Balay + ls - the TaoLineSearch context 13310c52818fSSatish Balay - prefix - the prefix string to prepend to all TAO option requests 13320c52818fSSatish Balay 13330c52818fSSatish Balay Notes: 13340c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 13350c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 13360c52818fSSatish Balay 13370c52818fSSatish Balay For example, to distinguish between the runtime options for two 13380c52818fSSatish Balay different line searches, one could call 13390c52818fSSatish Balay .vb 13400c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 13410c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 13420c52818fSSatish Balay .ve 13430c52818fSSatish Balay 13440c52818fSSatish Balay This would enable use of different options for each system, such as 13450c52818fSSatish Balay .vb 13460c52818fSSatish Balay -sys1_tao_ls_type mt 13470c52818fSSatish Balay -sys2_tao_ls_type armijo 13480c52818fSSatish Balay .ve 13490c52818fSSatish Balay 13500c52818fSSatish Balay Level: advanced 13510c52818fSSatish Balay 13520c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 13530c52818fSSatish Balay @*/ 13540c52818fSSatish Balay 13550c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 13560c52818fSSatish Balay { 1357f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 13580c52818fSSatish Balay } 1359