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 PetscErrorCode ierr; 27fe2efc57SMark 28fe2efc57SMark PetscFunctionBegin; 29fe2efc57SMark PetscValidHeaderSpecific(A,TAOLINESEARCH_CLASSID,1); 30fe2efc57SMark ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr); 31fe2efc57SMark PetscFunctionReturn(0); 32fe2efc57SMark } 33fe2efc57SMark 34fe2efc57SMark /*@C 350c52818fSSatish Balay TaoLineSearchView - Prints information about the TaoLineSearch 360c52818fSSatish Balay 370c52818fSSatish Balay Collective on TaoLineSearch 380c52818fSSatish Balay 390c52818fSSatish Balay InputParameters: 40441846f8SBarry Smith + ls - the Tao context 410c52818fSSatish Balay - viewer - visualization context 420c52818fSSatish Balay 430c52818fSSatish Balay Options Database Key: 440c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search 450c52818fSSatish Balay 460c52818fSSatish Balay Notes: 470c52818fSSatish Balay The available visualization contexts include 480c52818fSSatish Balay + PETSC_VIEWER_STDOUT_SELF - standard output (default) 490c52818fSSatish Balay - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 500c52818fSSatish Balay output where only the first processor opens 510c52818fSSatish Balay the file. All other processors send their 520c52818fSSatish Balay data to the first processor to print. 530c52818fSSatish Balay 540c52818fSSatish Balay Level: beginner 550c52818fSSatish Balay 560c52818fSSatish Balay .seealso: PetscViewerASCIIOpen() 570c52818fSSatish Balay @*/ 580c52818fSSatish Balay 590c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 600c52818fSSatish Balay { 610c52818fSSatish Balay PetscErrorCode ierr; 620c52818fSSatish Balay PetscBool isascii, isstring; 63dedfbcbeSJed Brown TaoLineSearchType type; 64f06e3bfaSBarry Smith 650c52818fSSatish Balay PetscFunctionBegin; 660c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 670c52818fSSatish Balay if (!viewer) { 680c52818fSSatish Balay ierr = PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);CHKERRQ(ierr); 690c52818fSSatish Balay } 700c52818fSSatish Balay PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 710c52818fSSatish Balay PetscCheckSameComm(ls,1,viewer,2); 720c52818fSSatish Balay 730c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 740c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 750c52818fSSatish Balay if (isascii) { 7663b15415SAlp Dener ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);CHKERRQ(ierr); 7763b15415SAlp Dener if (ls->ops->view) { 7863b15415SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7963b15415SAlp Dener ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr); 8063b15415SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 810c52818fSSatish Balay } 820c52818fSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 830c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr); 84f06e3bfaSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);CHKERRQ(ierr); 850c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr); 860c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr); 870c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr); 880c52818fSSatish Balay 890c52818fSSatish Balay if (ls->bounded) { 900c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr); 910c52818fSSatish Balay } 92de196d00SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);CHKERRQ(ierr); 930c52818fSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 940c52818fSSatish Balay } else if (isstring) { 950c52818fSSatish Balay ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr); 960c52818fSSatish Balay ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 970c52818fSSatish Balay } 980c52818fSSatish Balay PetscFunctionReturn(0); 990c52818fSSatish Balay } 1000c52818fSSatish Balay 1010c52818fSSatish Balay /*@C 1020c52818fSSatish Balay TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 1030c52818fSSatish Balay line-searches will automatically create one. 1040c52818fSSatish Balay 105d083f849SBarry Smith Collective 1060c52818fSSatish Balay 1070c52818fSSatish Balay Input Parameter: 1080c52818fSSatish Balay . comm - MPI communicator 1090c52818fSSatish Balay 1100c52818fSSatish Balay Output Parameter: 1110c52818fSSatish Balay . newls - the new TaoLineSearch context 1120c52818fSSatish Balay 1130c52818fSSatish Balay Available methods include: 114147403d9SBarry Smith + more-thuente - the More-Thuente method 115147403d9SBarry Smith . gpcg - the GPCG method 1160c52818fSSatish Balay - unit - Do not perform any line search 1170c52818fSSatish Balay 1180c52818fSSatish Balay Options Database Keys: 1190c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 1200c52818fSSatish Balay 1210c52818fSSatish Balay Level: beginner 1220c52818fSSatish Balay 1230c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy() 1240c52818fSSatish Balay @*/ 1250c52818fSSatish Balay 1260c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 1270c52818fSSatish Balay { 1280c52818fSSatish Balay PetscErrorCode ierr; 1290c52818fSSatish Balay TaoLineSearch ls; 1300c52818fSSatish Balay 1310c52818fSSatish Balay PetscFunctionBegin; 1320c52818fSSatish Balay PetscValidPointer(newls,2); 1336c23d075SBarry Smith *newls = NULL; 1340c52818fSSatish Balay 1350c52818fSSatish Balay ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr); 1360c52818fSSatish Balay 13773107ff1SLisandro Dalcin ierr = PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);CHKERRQ(ierr); 1380c52818fSSatish Balay ls->bounded = 0; 1390c52818fSSatish Balay ls->max_funcs=30; 1400c52818fSSatish Balay ls->ftol = 0.0001; 1410c52818fSSatish Balay ls->gtol = 0.9; 1426f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1436f4723b1SBarry Smith ls->rtol = 1.0e-5; 1446f4723b1SBarry Smith #else 1450c52818fSSatish Balay ls->rtol = 1.0e-10; 1466f4723b1SBarry Smith #endif 1470c52818fSSatish Balay ls->stepmin=1.0e-20; 1480c52818fSSatish Balay ls->stepmax=1.0e+20; 1490c52818fSSatish Balay ls->step=1.0; 1500c52818fSSatish Balay ls->nfeval=0; 1510c52818fSSatish Balay ls->ngeval=0; 1520c52818fSSatish Balay ls->nfgeval=0; 1530c52818fSSatish Balay 15483c8fe1dSLisandro Dalcin ls->ops->computeobjective = NULL; 15583c8fe1dSLisandro Dalcin ls->ops->computegradient = NULL; 15683c8fe1dSLisandro Dalcin ls->ops->computeobjectiveandgradient = NULL; 15783c8fe1dSLisandro Dalcin ls->ops->computeobjectiveandgts = NULL; 15883c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 15983c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 16083c8fe1dSLisandro Dalcin ls->ops->view = NULL; 16183c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 16283c8fe1dSLisandro Dalcin ls->ops->reset = NULL; 16383c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 16483c8fe1dSLisandro Dalcin ls->ops->monitor = NULL; 1652a0dac07SAlp Dener ls->usemonitor=PETSC_FALSE; 1660c52818fSSatish Balay ls->setupcalled=PETSC_FALSE; 1670c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 1680c52818fSSatish Balay *newls = ls; 1690c52818fSSatish Balay PetscFunctionReturn(0); 1700c52818fSSatish Balay } 1710c52818fSSatish Balay 1720c52818fSSatish Balay /*@ 1730c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 1740c52818fSSatish Balay of a Tao solver 1750c52818fSSatish Balay 1760c52818fSSatish Balay Collective on ls 1770c52818fSSatish Balay 1780c52818fSSatish Balay Input Parameters: 1790c52818fSSatish Balay . ls - the TaoLineSearch context 1800c52818fSSatish Balay 1810c52818fSSatish Balay Notes: 1820c52818fSSatish Balay The user will not need to explicitly call TaoLineSearchSetUp(), as it will 1830c52818fSSatish Balay automatically be called in TaoLineSearchSolve(). However, if the user 1840c52818fSSatish Balay desires to call it explicitly, it should come after TaoLineSearchCreate() 1850c52818fSSatish Balay but before TaoLineSearchApply(). 1860c52818fSSatish Balay 1870c52818fSSatish Balay Level: developer 1880c52818fSSatish Balay 1890c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 1900c52818fSSatish Balay @*/ 1910c52818fSSatish Balay 1920c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 1930c52818fSSatish Balay { 1940c52818fSSatish Balay PetscErrorCode ierr; 1958caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 1960c52818fSSatish Balay PetscBool flg; 197f06e3bfaSBarry Smith 1980c52818fSSatish Balay PetscFunctionBegin; 1990c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 2000c52818fSSatish Balay if (ls->setupcalled) PetscFunctionReturn(0); 2010c52818fSSatish Balay if (!((PetscObject)ls)->type_name) { 2020c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 2030c52818fSSatish Balay } 2040c52818fSSatish Balay if (ls->ops->setup) { 2050c52818fSSatish Balay ierr = (*ls->ops->setup)(ls);CHKERRQ(ierr); 2060c52818fSSatish Balay } 2070c52818fSSatish Balay if (ls->usetaoroutines) { 208441846f8SBarry Smith ierr = TaoIsObjectiveDefined(ls->tao,&flg);CHKERRQ(ierr); 2090c52818fSSatish Balay ls->hasobjective = flg; 210441846f8SBarry Smith ierr = TaoIsGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 2110c52818fSSatish Balay ls->hasgradient = flg; 212441846f8SBarry Smith ierr = TaoIsObjectiveAndGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 2130c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 2140c52818fSSatish Balay } else { 2150c52818fSSatish Balay if (ls->ops->computeobjective) { 2160c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 2170c52818fSSatish Balay } else { 2180c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 2190c52818fSSatish Balay } 2200c52818fSSatish Balay if (ls->ops->computegradient) { 2210c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 2220c52818fSSatish Balay } else { 2230c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 2240c52818fSSatish Balay } 2250c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 2260c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 2270c52818fSSatish Balay } else { 2280c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 2290c52818fSSatish Balay } 2300c52818fSSatish Balay } 2310c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 2320c52818fSSatish Balay PetscFunctionReturn(0); 2330c52818fSSatish Balay } 2340c52818fSSatish Balay 2350c52818fSSatish Balay /*@ 2360c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 2370c52818fSSatish Balay from one TaoLineSearchApply() to the next. This function resets this 2380c52818fSSatish Balay state information. 2390c52818fSSatish Balay 2400c52818fSSatish Balay Collective on TaoLineSearch 2410c52818fSSatish Balay 2420c52818fSSatish Balay Input Parameter: 2430c52818fSSatish Balay . ls - the TaoLineSearch context 2440c52818fSSatish Balay 2450c52818fSSatish Balay Level: developer 2460c52818fSSatish Balay 2470c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 2480c52818fSSatish Balay @*/ 2490c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 2500c52818fSSatish Balay { 2510c52818fSSatish Balay PetscErrorCode ierr; 252050fc7a3SBarry Smith 2530c52818fSSatish Balay PetscFunctionBegin; 2540c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 2550c52818fSSatish Balay if (ls->ops->reset) { 2560c52818fSSatish Balay ierr = (*ls->ops->reset)(ls);CHKERRQ(ierr); 2570c52818fSSatish Balay } 2580c52818fSSatish Balay PetscFunctionReturn(0); 2590c52818fSSatish Balay } 260f06e3bfaSBarry Smith 2610c52818fSSatish Balay /*@ 2620c52818fSSatish Balay TaoLineSearchDestroy - Destroys the TAO context that was created with 2630c52818fSSatish Balay TaoLineSearchCreate() 2640c52818fSSatish Balay 2650c52818fSSatish Balay Collective on TaoLineSearch 2660c52818fSSatish Balay 2677a7aea1fSJed Brown Input Parameter: 2680c52818fSSatish Balay . ls - the TaoLineSearch context 2690c52818fSSatish Balay 2700c52818fSSatish Balay Level: beginner 2710c52818fSSatish Balay 2720c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 2730c52818fSSatish Balay @*/ 2740c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 2750c52818fSSatish Balay { 2760c52818fSSatish Balay PetscErrorCode ierr; 277050fc7a3SBarry Smith 2780c52818fSSatish Balay PetscFunctionBegin; 2790c52818fSSatish Balay if (!*ls) PetscFunctionReturn(0); 2800c52818fSSatish Balay PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 28183c8fe1dSLisandro Dalcin if (--((PetscObject)*ls)->refct > 0) {*ls = NULL; PetscFunctionReturn(0);} 282050fc7a3SBarry Smith ierr = VecDestroy(&(*ls)->stepdirection);CHKERRQ(ierr); 283050fc7a3SBarry Smith ierr = VecDestroy(&(*ls)->start_x);CHKERRQ(ierr); 28497ab8e72SStefano Zampini ierr = VecDestroy(&(*ls)->upper);CHKERRQ(ierr); 28597ab8e72SStefano Zampini ierr = VecDestroy(&(*ls)->lower);CHKERRQ(ierr); 2860c52818fSSatish Balay if ((*ls)->ops->destroy) { 2870c52818fSSatish Balay ierr = (*(*ls)->ops->destroy)(*ls);CHKERRQ(ierr); 2880c52818fSSatish Balay } 2892a0dac07SAlp Dener if ((*ls)->usemonitor) { 2902a0dac07SAlp Dener ierr = PetscViewerDestroy(&(*ls)->viewer);CHKERRQ(ierr); 2912a0dac07SAlp Dener } 2920c52818fSSatish Balay ierr = PetscHeaderDestroy(ls);CHKERRQ(ierr); 2930c52818fSSatish Balay PetscFunctionReturn(0); 2940c52818fSSatish Balay } 2950c52818fSSatish Balay 2960c52818fSSatish Balay /*@ 2970c52818fSSatish Balay TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 2980c52818fSSatish Balay 2990c52818fSSatish Balay Collective on TaoLineSearch 3000c52818fSSatish Balay 3010c52818fSSatish Balay Input Parameters: 302441846f8SBarry Smith + ls - the Tao context 3030c52818fSSatish Balay - s - search direction 3040c52818fSSatish Balay 3056b867d5aSJose E. Roman Input/Output Parameters: 3066b867d5aSJose E. Roman 3070c52818fSSatish Balay Output Parameters: 308f1a722f8SMatthew G. Knepley + x - On input the current solution, on output x contains the new solution determined by the line search 309f1a722f8SMatthew G. Knepley . f - On input the objective function value at current solution, on output contains the objective function value at new solution 310f1a722f8SMatthew G. Knepley . g - On input the gradient evaluated at x, on output contains the gradient at new solution 311f1a722f8SMatthew G. Knepley . steplength - scalar multiplier of s used ( x = x0 + steplength * x) 3120c52818fSSatish Balay - reason - reason why the line-search stopped 3130c52818fSSatish Balay 31497bb3fdcSJose E. Roman Notes: 3150c52818fSSatish Balay reason will be set to one of: 3160c52818fSSatish Balay 3170c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 3180c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 3190c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 3200c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 3210c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 3220c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 3230c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 3240c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 3250c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 3260c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search 3270c52818fSSatish Balay 3280c52818fSSatish Balay The algorithm developer must set up the TaoLineSearch with calls to 329441846f8SBarry Smith TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines() 3300c52818fSSatish Balay 3310c52818fSSatish Balay You may or may not need to follow this with a call to 3320c52818fSSatish Balay TaoAddLineSearchCounts(), depending on whether you want these 3330c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 3340c52818fSSatish Balay 3350c52818fSSatish Balay Level: beginner 3360c52818fSSatish Balay 3370c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts() 3380c52818fSSatish Balay @*/ 3390c52818fSSatish Balay 340e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 3410c52818fSSatish Balay { 3420c52818fSSatish Balay PetscErrorCode ierr; 3430c52818fSSatish Balay PetscInt low1,low2,low3,high1,high2,high3; 3440c52818fSSatish Balay 3450c52818fSSatish Balay PetscFunctionBegin; 3460c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 3470c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3483f6ba705SLisandro Dalcin PetscValidRealPointer(f,3); 3490c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 3500c52818fSSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,5); 3510c52818fSSatish Balay PetscValidPointer(reason,7); 3520c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 3530c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,g,4); 3540c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,s,5); 3550c52818fSSatish Balay ierr = VecGetOwnershipRange(x, &low1, &high1);CHKERRQ(ierr); 3560c52818fSSatish Balay ierr = VecGetOwnershipRange(g, &low2, &high2);CHKERRQ(ierr); 3570c52818fSSatish Balay ierr = VecGetOwnershipRange(s, &low3, &high3);CHKERRQ(ierr); 358*3c859ba3SBarry Smith PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible vector local lengths"); 3590c52818fSSatish Balay 36097ab8e72SStefano Zampini *reason = TAOLINESEARCH_CONTINUE_ITERATING; 3610c52818fSSatish Balay ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr); 362050fc7a3SBarry Smith ierr = VecDestroy(&ls->stepdirection);CHKERRQ(ierr); 363050fc7a3SBarry Smith ls->stepdirection = s; 3640c52818fSSatish Balay 3650c52818fSSatish Balay ierr = TaoLineSearchSetUp(ls);CHKERRQ(ierr); 366*3c859ba3SBarry Smith PetscCheck(ls->ops->apply,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine"); 3670c52818fSSatish Balay ls->nfeval = 0; 3680c52818fSSatish Balay ls->ngeval = 0; 3690c52818fSSatish Balay ls->nfgeval = 0; 3700c52818fSSatish Balay /* Check parameter values */ 3710c52818fSSatish Balay if (ls->ftol < 0.0) { 3727d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);CHKERRQ(ierr); 3730c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3740c52818fSSatish Balay } 3750c52818fSSatish Balay if (ls->rtol < 0.0) { 3767d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);CHKERRQ(ierr); 3770c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3780c52818fSSatish Balay } 3790c52818fSSatish Balay if (ls->gtol < 0.0) { 3807d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);CHKERRQ(ierr); 3810c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3820c52818fSSatish Balay } 3830c52818fSSatish Balay if (ls->stepmin < 0.0) { 3847d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);CHKERRQ(ierr); 3850c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3860c52818fSSatish Balay } 3870c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 3887d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);CHKERRQ(ierr); 3890c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3900c52818fSSatish Balay } 3910c52818fSSatish Balay if (ls->max_funcs < 0) { 3927d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n",ls->max_funcs);CHKERRQ(ierr); 3930c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_BADPARAMETER; 3940c52818fSSatish Balay } 3950c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 3967d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);CHKERRQ(ierr); 3970c52818fSSatish Balay *reason = TAOLINESEARCH_FAILED_INFORNAN; 3980c52818fSSatish Balay } 3990c52818fSSatish Balay 400302440fdSBarry Smith ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 4010c52818fSSatish Balay ierr = VecDestroy(&ls->start_x);CHKERRQ(ierr); 4020c52818fSSatish Balay ls->start_x = x; 403050fc7a3SBarry Smith 4040ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0);CHKERRQ(ierr); 4050c52818fSSatish Balay ierr = (*ls->ops->apply)(ls,x,f,g,s);CHKERRQ(ierr); 4060ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0);CHKERRQ(ierr); 4070c52818fSSatish Balay *reason = ls->reason; 4080c52818fSSatish Balay ls->new_f = *f; 4090c52818fSSatish Balay 41097ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 4110c52818fSSatish Balay 412fbe0838dSJason Sarich ierr = TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");CHKERRQ(ierr); 4130c52818fSSatish Balay PetscFunctionReturn(0); 4140c52818fSSatish Balay } 4150c52818fSSatish Balay 4160c52818fSSatish Balay /*@C 4170c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 4180c52818fSSatish Balay 4190c52818fSSatish Balay Collective on TaoLineSearch 4200c52818fSSatish Balay 4210c52818fSSatish Balay Input Parameters: 4220c52818fSSatish Balay + ls - the TaoLineSearch context 423820824deSAlp Dener - type - the TaoLineSearchType selection 4240c52818fSSatish Balay 4250c52818fSSatish Balay Available methods include: 426820824deSAlp Dener + more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition 427820824deSAlp Dener . armijo - simple backtracking line search enforcing only the sufficient decrease condition 428820824deSAlp Dener - unit - do not perform a line search and always accept unit step length 4290c52818fSSatish Balay 4300c52818fSSatish Balay Options Database Keys: 431820824deSAlp Dener . -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime 4320c52818fSSatish Balay 4330c52818fSSatish Balay Level: beginner 4340c52818fSSatish Balay 4350c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply() 4360c52818fSSatish Balay 4370c52818fSSatish Balay @*/ 4380c52818fSSatish Balay 439dedfbcbeSJed Brown PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 4400c52818fSSatish Balay { 4410c52818fSSatish Balay PetscErrorCode ierr; 4420c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 4430c52818fSSatish Balay PetscBool flg; 4440c52818fSSatish Balay 4450c52818fSSatish Balay PetscFunctionBegin; 4460c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 4470c52818fSSatish Balay PetscValidCharPointer(type,2); 4480c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg);CHKERRQ(ierr); 4490c52818fSSatish Balay if (flg) PetscFunctionReturn(0); 4500c52818fSSatish Balay 4510c52818fSSatish Balay ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);CHKERRQ(ierr); 452*3c859ba3SBarry Smith PetscCheck(r,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 4530c52818fSSatish Balay if (ls->ops->destroy) { 4540c52818fSSatish Balay ierr = (*(ls)->ops->destroy)(ls);CHKERRQ(ierr); 4550c52818fSSatish Balay } 4560c52818fSSatish Balay ls->max_funcs = 30; 4570c52818fSSatish Balay ls->ftol = 0.0001; 4580c52818fSSatish Balay ls->gtol = 0.9; 4596f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 4606f4723b1SBarry Smith ls->rtol = 1.0e-5; 4616f4723b1SBarry Smith #else 4620c52818fSSatish Balay ls->rtol = 1.0e-10; 4636f4723b1SBarry Smith #endif 4640c52818fSSatish Balay ls->stepmin = 1.0e-20; 4650c52818fSSatish Balay ls->stepmax = 1.0e+20; 4660c52818fSSatish Balay 4670c52818fSSatish Balay ls->nfeval = 0; 4680c52818fSSatish Balay ls->ngeval = 0; 4690c52818fSSatish Balay ls->nfgeval = 0; 47083c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 47183c8fe1dSLisandro Dalcin ls->ops->apply = NULL; 47283c8fe1dSLisandro Dalcin ls->ops->view = NULL; 47383c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 47483c8fe1dSLisandro Dalcin ls->ops->destroy = NULL; 4750c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 4760c52818fSSatish Balay ierr = (*r)(ls);CHKERRQ(ierr); 4770c52818fSSatish Balay ierr = PetscObjectChangeTypeName((PetscObject)ls, type);CHKERRQ(ierr); 4780c52818fSSatish Balay PetscFunctionReturn(0); 4790c52818fSSatish Balay } 4800c52818fSSatish Balay 4812a0dac07SAlp Dener /*@C 4822a0dac07SAlp Dener TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the 4832a0dac07SAlp Dener iteration number, step length, and function value before calling the implementation 4842a0dac07SAlp Dener specific monitor. 4852a0dac07SAlp Dener 4862a0dac07SAlp Dener Input Parameters: 4872a0dac07SAlp Dener + ls - the TaoLineSearch context 4882a0dac07SAlp Dener . its - the current iterate number (>=0) 4892a0dac07SAlp Dener . f - the current objective function value 4902a0dac07SAlp Dener - step - the step length 4912a0dac07SAlp Dener 4922a0dac07SAlp Dener Options Database Key: 4932a0dac07SAlp Dener . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output 4942a0dac07SAlp Dener 4952a0dac07SAlp Dener Level: developer 4962a0dac07SAlp Dener 4972a0dac07SAlp Dener @*/ 4982a0dac07SAlp Dener PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step) 4992a0dac07SAlp Dener { 5002a0dac07SAlp Dener PetscErrorCode ierr; 5012a0dac07SAlp Dener PetscInt tabs; 5022a0dac07SAlp Dener 5032a0dac07SAlp Dener PetscFunctionBegin; 5042a0dac07SAlp Dener PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5052a0dac07SAlp Dener if (ls->usemonitor) { 5062a0dac07SAlp Dener ierr = PetscViewerASCIIGetTab(ls->viewer, &tabs);CHKERRQ(ierr); 5072a0dac07SAlp Dener ierr = PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel);CHKERRQ(ierr); 5082a0dac07SAlp Dener ierr = PetscViewerASCIIPrintf(ls->viewer, "%3D LS", its);CHKERRQ(ierr); 5092a0dac07SAlp Dener ierr = PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f);CHKERRQ(ierr); 510ad362606SAlp Dener ierr = PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step);CHKERRQ(ierr); 5112a0dac07SAlp Dener if (ls->ops->monitor && its > 0) { 5122a0dac07SAlp Dener ierr = PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3);CHKERRQ(ierr); 5132a0dac07SAlp Dener ierr = (*ls->ops->monitor)(ls);CHKERRQ(ierr); 5142a0dac07SAlp Dener } 5152a0dac07SAlp Dener ierr = PetscViewerASCIISetTab(ls->viewer, tabs);CHKERRQ(ierr); 5162a0dac07SAlp Dener } 5172a0dac07SAlp Dener PetscFunctionReturn(0); 5182a0dac07SAlp Dener } 5192a0dac07SAlp Dener 5200c52818fSSatish Balay /*@ 5210c52818fSSatish Balay TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 5220c52818fSSatish Balay options. 5230c52818fSSatish Balay 5240c52818fSSatish Balay Collective on TaoLineSearch 5250c52818fSSatish Balay 52601d2d390SJose E. Roman Input Parameter: 5270c52818fSSatish Balay . ls - the TaoLineSearch context 5280c52818fSSatish Balay 5290c52818fSSatish Balay Options Database Keys: 5300c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 5310c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 5320c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 5330c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 5340c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 5350c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 5360c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 5370c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 5380c52818fSSatish Balay 5390c52818fSSatish Balay Level: beginner 5400c52818fSSatish Balay @*/ 5410c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 5420c52818fSSatish Balay { 5430c52818fSSatish Balay PetscErrorCode ierr; 5448caf6e8cSBarry Smith const char *default_type=TAOLINESEARCHMT; 5452a0dac07SAlp Dener char type[256],monfilename[PETSC_MAX_PATH_LEN]; 5462a0dac07SAlp Dener PetscViewer monviewer; 5470c52818fSSatish Balay PetscBool flg; 548f06e3bfaSBarry Smith 5490c52818fSSatish Balay PetscFunctionBegin; 5500c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 5510c52818fSSatish Balay ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr); 5520c52818fSSatish Balay if (((PetscObject)ls)->type_name) { 5530c52818fSSatish Balay default_type = ((PetscObject)ls)->type_name; 5540c52818fSSatish Balay } 5550c52818fSSatish Balay /* Check for type from options */ 5560c52818fSSatish Balay ierr = PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);CHKERRQ(ierr); 5570c52818fSSatish Balay if (flg) { 5580c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,type);CHKERRQ(ierr); 5590c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 560302440fdSBarry Smith ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 5610c52818fSSatish Balay } 5620c52818fSSatish Balay 56394ae4db5SBarry Smith ierr = PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);CHKERRQ(ierr); 56494ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);CHKERRQ(ierr); 56594ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);CHKERRQ(ierr); 56694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);CHKERRQ(ierr); 56794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);CHKERRQ(ierr); 56894ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);CHKERRQ(ierr); 569589a23caSBarry Smith ierr = PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,sizeof(monfilename),&flg);CHKERRQ(ierr); 5702a0dac07SAlp Dener if (flg) { 5712a0dac07SAlp Dener ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer);CHKERRQ(ierr); 5722a0dac07SAlp Dener ls->viewer = monviewer; 5732a0dac07SAlp Dener ls->usemonitor = PETSC_TRUE; 5742a0dac07SAlp Dener } 5750c52818fSSatish Balay if (ls->ops->setfromoptions) { 5763a004c28SBarry Smith ierr = (*ls->ops->setfromoptions)(PetscOptionsObject,ls);CHKERRQ(ierr); 5770c52818fSSatish Balay } 5780c52818fSSatish Balay ierr = PetscOptionsEnd();CHKERRQ(ierr); 5790c52818fSSatish Balay PetscFunctionReturn(0); 5800c52818fSSatish Balay } 5810c52818fSSatish Balay 5820c52818fSSatish Balay /*@C 5830c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 5840c52818fSSatish Balay 5850c52818fSSatish Balay Not Collective 5860c52818fSSatish Balay 5870c52818fSSatish Balay Input Parameter: 5880c52818fSSatish Balay . ls - the TaoLineSearch context 5890c52818fSSatish Balay 590fd292e60Sprj- Output Parameter: 5910c52818fSSatish Balay . type - the line search algorithm in effect 5920c52818fSSatish Balay 5930c52818fSSatish Balay Level: developer 5940c52818fSSatish Balay 5950c52818fSSatish Balay @*/ 596dedfbcbeSJed Brown PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type) 5970c52818fSSatish Balay { 5980c52818fSSatish Balay PetscFunctionBegin; 5990c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6000c52818fSSatish Balay PetscValidPointer(type,2); 6010c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 6020c52818fSSatish Balay PetscFunctionReturn(0); 6030c52818fSSatish Balay } 6040c52818fSSatish Balay 6050c52818fSSatish Balay /*@ 6060c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 6070c52818fSSatish Balay routines used by the line search in last application (not cumulative). 6080c52818fSSatish Balay 6090c52818fSSatish Balay Not Collective 6100c52818fSSatish Balay 6110c52818fSSatish Balay Input Parameter: 6120c52818fSSatish Balay . ls - the TaoLineSearch context 6130c52818fSSatish Balay 6140c52818fSSatish Balay Output Parameters: 6150c52818fSSatish Balay + nfeval - number of function evaluations 6160c52818fSSatish Balay . ngeval - number of gradient evaluations 6170c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 6180c52818fSSatish Balay 6190c52818fSSatish Balay Level: intermediate 6200c52818fSSatish Balay 6210c52818fSSatish Balay Note: 622441846f8SBarry Smith If the line search is using the Tao objective and gradient 623441846f8SBarry Smith routines directly (see TaoLineSearchUseTaoRoutines()), then TAO 6240c52818fSSatish Balay is already counting the number of evaluations. 6250c52818fSSatish Balay 6260c52818fSSatish Balay @*/ 6270c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 6280c52818fSSatish Balay { 6290c52818fSSatish Balay PetscFunctionBegin; 6300c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6310c52818fSSatish Balay *nfeval = ls->nfeval; 6320c52818fSSatish Balay *ngeval = ls->ngeval; 6330c52818fSSatish Balay *nfgeval = ls->nfgeval; 6340c52818fSSatish Balay PetscFunctionReturn(0); 6350c52818fSSatish Balay } 6360c52818fSSatish Balay 6370c52818fSSatish Balay /*@ 638441846f8SBarry Smith TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using 639441846f8SBarry Smith Tao evaluation routines. 6400c52818fSSatish Balay 6410c52818fSSatish Balay Not Collective 6420c52818fSSatish Balay 6430c52818fSSatish Balay Input Parameter: 6440c52818fSSatish Balay . ls - the TaoLineSearch context 6450c52818fSSatish Balay 6460c52818fSSatish Balay Output Parameter: 647441846f8SBarry Smith . flg - PETSC_TRUE if the line search is using Tao evaluation routines, 6480c52818fSSatish Balay otherwise PETSC_FALSE 6490c52818fSSatish Balay 6500c52818fSSatish Balay Level: developer 6510c52818fSSatish Balay @*/ 652441846f8SBarry Smith PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg) 6530c52818fSSatish Balay { 6540c52818fSSatish Balay PetscFunctionBegin; 6550c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 656f06e3bfaSBarry Smith *flg = ls->usetaoroutines; 6570c52818fSSatish Balay PetscFunctionReturn(0); 6580c52818fSSatish Balay } 6590c52818fSSatish Balay 6600c52818fSSatish Balay /*@C 6610c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 6620c52818fSSatish Balay 6630c52818fSSatish Balay Logically Collective on TaoLineSearch 6640c52818fSSatish Balay 665d8d19677SJose E. Roman Input Parameters: 6660c52818fSSatish Balay + ls - the TaoLineSearch context 6670c52818fSSatish Balay . func - the objective function evaluation routine 6680c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 6690c52818fSSatish Balay 6700c52818fSSatish Balay Calling sequence of func: 6710c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 6720c52818fSSatish Balay 6730c52818fSSatish Balay + x - input vector 6740c52818fSSatish Balay . f - function value 6750c52818fSSatish Balay - ctx (optional) user-defined context 6760c52818fSSatish Balay 6770c52818fSSatish Balay Level: beginner 6780c52818fSSatish Balay 6790c52818fSSatish Balay Note: 6800c52818fSSatish Balay Use this routine only if you want the line search objective 681441846f8SBarry Smith evaluation routine to be different from the Tao's objective 6820c52818fSSatish Balay evaluation routine. If you use this routine you must also set 6830c52818fSSatish Balay the line search gradient and/or function/gradient routine. 6840c52818fSSatish Balay 6850c52818fSSatish Balay Note: 6860c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 6870c52818fSSatish Balay line search, application programmers should be wary of overriding the 6880c52818fSSatish Balay default objective routine. 6890c52818fSSatish Balay 690441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 6910c52818fSSatish Balay @*/ 6920c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx) 6930c52818fSSatish Balay { 6940c52818fSSatish Balay PetscFunctionBegin; 6950c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 6960c52818fSSatish Balay 6970c52818fSSatish Balay ls->ops->computeobjective=func; 6980c52818fSSatish Balay if (ctx) ls->userctx_func=ctx; 6990c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 7000c52818fSSatish Balay PetscFunctionReturn(0); 7010c52818fSSatish Balay } 7020c52818fSSatish Balay 7030c52818fSSatish Balay /*@C 7040c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 7050c52818fSSatish Balay 7060c52818fSSatish Balay Logically Collective on TaoLineSearch 7070c52818fSSatish Balay 708d8d19677SJose E. Roman Input Parameters: 7090c52818fSSatish Balay + ls - the TaoLineSearch context 7100c52818fSSatish Balay . func - the gradient evaluation routine 7110c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7120c52818fSSatish Balay 7130c52818fSSatish Balay Calling sequence of func: 7140c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, Vec g, void *ctx); 7150c52818fSSatish Balay 7160c52818fSSatish Balay + x - input vector 7170c52818fSSatish Balay . g - gradient vector 7180c52818fSSatish Balay - ctx (optional) user-defined context 7190c52818fSSatish Balay 7200c52818fSSatish Balay Level: beginner 7210c52818fSSatish Balay 7220c52818fSSatish Balay Note: 7230c52818fSSatish Balay Use this routine only if you want the line search gradient 724441846f8SBarry Smith evaluation routine to be different from the Tao's gradient 7250c52818fSSatish Balay evaluation routine. If you use this routine you must also set 7260c52818fSSatish Balay the line search function and/or function/gradient routine. 7270c52818fSSatish Balay 7280c52818fSSatish Balay Note: 7290c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 7300c52818fSSatish Balay line search, application programmers should be wary of overriding the 7310c52818fSSatish Balay default gradient routine. 7320c52818fSSatish Balay 733441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines() 7340c52818fSSatish Balay @*/ 7350c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx) 7360c52818fSSatish Balay { 7370c52818fSSatish Balay PetscFunctionBegin; 7380c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7390c52818fSSatish Balay ls->ops->computegradient=func; 7400c52818fSSatish Balay if (ctx) ls->userctx_grad=ctx; 7410c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 7420c52818fSSatish Balay PetscFunctionReturn(0); 7430c52818fSSatish Balay } 7440c52818fSSatish Balay 7450c52818fSSatish Balay /*@C 7460c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 7470c52818fSSatish Balay 7480c52818fSSatish Balay Logically Collective on TaoLineSearch 7490c52818fSSatish Balay 750d8d19677SJose E. Roman Input Parameters: 7510c52818fSSatish Balay + ls - the TaoLineSearch context 7520c52818fSSatish Balay . func - the objective and gradient evaluation routine 7530c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 7540c52818fSSatish Balay 7550c52818fSSatish Balay Calling sequence of func: 7560c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 7570c52818fSSatish Balay 7580c52818fSSatish Balay + x - input vector 7590c52818fSSatish Balay . f - function value 7600c52818fSSatish Balay . g - gradient vector 7610c52818fSSatish Balay - ctx (optional) user-defined context 7620c52818fSSatish Balay 7630c52818fSSatish Balay Level: beginner 7640c52818fSSatish Balay 7650c52818fSSatish Balay Note: 7660c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 767441846f8SBarry Smith evaluation routines to be different from the Tao's objective 7680c52818fSSatish Balay and gradient evaluation routines. 7690c52818fSSatish Balay 7700c52818fSSatish Balay Note: 7710c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 7720c52818fSSatish Balay line search, application programmers should be wary of overriding the 7730c52818fSSatish Balay default objective routine. 7740c52818fSSatish Balay 775441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines() 7760c52818fSSatish Balay @*/ 7770c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx) 7780c52818fSSatish Balay { 7790c52818fSSatish Balay PetscFunctionBegin; 7800c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 7810c52818fSSatish Balay ls->ops->computeobjectiveandgradient=func; 7820c52818fSSatish Balay if (ctx) ls->userctx_funcgrad=ctx; 7830c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 7840c52818fSSatish Balay PetscFunctionReturn(0); 7850c52818fSSatish Balay } 7860c52818fSSatish Balay 7870c52818fSSatish Balay /*@C 7880c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 7890c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 7900c52818fSSatish Balay Sometimes it is more efficient to compute the inner product of the gradient 7910c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 7920c52818fSSatish Balay the line search typically needs of the gradient. 7930c52818fSSatish Balay 7940c52818fSSatish Balay Logically Collective on TaoLineSearch 7950c52818fSSatish Balay 796d8d19677SJose E. Roman Input Parameters: 7970c52818fSSatish Balay + ls - the TaoLineSearch context 7980c52818fSSatish Balay . func - the objective and gradient evaluation routine 7990c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 8000c52818fSSatish Balay 8010c52818fSSatish Balay Calling sequence of func: 8020c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 8030c52818fSSatish Balay 8040c52818fSSatish Balay + x - input vector 8050c52818fSSatish Balay . s - step direction 8060c52818fSSatish Balay . f - function value 8070c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 8080c52818fSSatish Balay - ctx (optional) user-defined context 8090c52818fSSatish Balay 8100c52818fSSatish Balay Note: The gradient will still need to be computed at the end of the line 8110c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 8120c52818fSSatish Balay routine 8130c52818fSSatish Balay 8140c52818fSSatish Balay Note: Bounded line searches (those used in bounded optimization algorithms) 8150c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 8160c52818fSSatish Balay x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength() 8170c52818fSSatish Balay 8180c52818fSSatish Balay Level: advanced 8190c52818fSSatish Balay 8200c52818fSSatish Balay Note: 8210c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 8220c52818fSSatish Balay line search, application programmers should be wary of overriding the 8230c52818fSSatish Balay default objective routine. 8240c52818fSSatish Balay 825441846f8SBarry Smith .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines() 8260c52818fSSatish Balay @*/ 8270c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx) 8280c52818fSSatish Balay { 8290c52818fSSatish Balay PetscFunctionBegin; 8300c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8310c52818fSSatish Balay ls->ops->computeobjectiveandgts=func; 8320c52818fSSatish Balay if (ctx) ls->userctx_funcgts=ctx; 8330c52818fSSatish Balay ls->usegts = PETSC_TRUE; 8340c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 8350c52818fSSatish Balay PetscFunctionReturn(0); 8360c52818fSSatish Balay } 8370c52818fSSatish Balay 8380c52818fSSatish Balay /*@C 839441846f8SBarry Smith TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the 840441846f8SBarry Smith objective and gradient evaluation routines from the given Tao object. 8410c52818fSSatish Balay 8420c52818fSSatish Balay Logically Collective on TaoLineSearch 8430c52818fSSatish Balay 844d8d19677SJose E. Roman Input Parameters: 8450c52818fSSatish Balay + ls - the TaoLineSearch context 846441846f8SBarry Smith - ts - the Tao context with defined objective/gradient evaluation routines 8470c52818fSSatish Balay 8480c52818fSSatish Balay Level: developer 8490c52818fSSatish Balay 8500c52818fSSatish Balay .seealso: TaoLineSearchCreate() 8510c52818fSSatish Balay @*/ 852441846f8SBarry Smith PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts) 8530c52818fSSatish Balay { 8540c52818fSSatish Balay PetscFunctionBegin; 8550c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 856064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(ts,TAO_CLASSID,2); 857441846f8SBarry Smith ls->tao = ts; 8580c52818fSSatish Balay ls->usetaoroutines = PETSC_TRUE; 8590c52818fSSatish Balay PetscFunctionReturn(0); 8600c52818fSSatish Balay } 8610c52818fSSatish Balay 8620c52818fSSatish Balay /*@ 8630c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 8640c52818fSSatish Balay 8650c52818fSSatish Balay Collective on TaoLineSearch 8660c52818fSSatish Balay 8670c52818fSSatish Balay Input Parameters: 8680c52818fSSatish Balay + ls - the TaoLineSearch context 8690c52818fSSatish Balay - x - input vector 8700c52818fSSatish Balay 8710c52818fSSatish Balay Output Parameter: 8720c52818fSSatish Balay . f - Objective value at X 8730c52818fSSatish Balay 87495452b02SPatrick Sanan Notes: 87595452b02SPatrick Sanan TaoLineSearchComputeObjective() is typically used within line searches 8760c52818fSSatish Balay so most users would not generally call this routine themselves. 8770c52818fSSatish Balay 8780c52818fSSatish Balay Level: developer 8790c52818fSSatish Balay 8800c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 8810c52818fSSatish Balay @*/ 8820c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 8830c52818fSSatish Balay { 8840c52818fSSatish Balay PetscErrorCode ierr; 8850c52818fSSatish Balay Vec gdummy; 8860c52818fSSatish Balay PetscReal gts; 887f06e3bfaSBarry Smith 8880c52818fSSatish Balay PetscFunctionBegin; 8890c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 8900c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 8910c52818fSSatish Balay PetscValidPointer(f,3); 8920c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 8930c52818fSSatish Balay if (ls->usetaoroutines) { 894441846f8SBarry Smith ierr = TaoComputeObjective(ls->tao,x,f);CHKERRQ(ierr); 8950c52818fSSatish Balay } else { 896*3c859ba3SBarry 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"); 8970ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 8980c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective routine"); 8990c52818fSSatish Balay if (ls->ops->computeobjective) { 9000c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 9010c52818fSSatish Balay } else if (ls->ops->computeobjectiveandgradient) { 9020c52818fSSatish Balay ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr); 9030c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr); 9040c52818fSSatish Balay ierr = VecDestroy(&gdummy);CHKERRQ(ierr); 9050c52818fSSatish Balay } else { 9060c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts);CHKERRQ(ierr); 9070c52818fSSatish Balay } 9080c52818fSSatish Balay PetscStackPop; 9090ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9100c52818fSSatish Balay } 9110c52818fSSatish Balay ls->nfeval++; 9120c52818fSSatish Balay PetscFunctionReturn(0); 9130c52818fSSatish Balay } 9140c52818fSSatish Balay 9150c52818fSSatish Balay /*@ 9160c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 9170c52818fSSatish Balay 918441846f8SBarry Smith Collective on Tao 9190c52818fSSatish Balay 9200c52818fSSatish Balay Input Parameters: 9210c52818fSSatish Balay + ls - the TaoLineSearch context 9220c52818fSSatish Balay - x - input vector 9230c52818fSSatish Balay 924d8d19677SJose E. Roman Output Parameters: 9250c52818fSSatish Balay + f - Objective value at X 9260c52818fSSatish Balay - g - Gradient vector at X 9270c52818fSSatish Balay 92895452b02SPatrick Sanan Notes: 92995452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 9300c52818fSSatish Balay so most users would not generally call this routine themselves. 9310c52818fSSatish Balay 9320c52818fSSatish Balay Level: developer 9330c52818fSSatish Balay 9340c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 9350c52818fSSatish Balay @*/ 9360c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 9370c52818fSSatish Balay { 9380c52818fSSatish Balay PetscErrorCode ierr; 939f06e3bfaSBarry Smith 9400c52818fSSatish Balay PetscFunctionBegin; 9410c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9420c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 9430c52818fSSatish Balay PetscValidPointer(f,3); 9440c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 9450c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 9460c52818fSSatish Balay PetscCheckSameComm(ls,1,g,4); 9470c52818fSSatish Balay if (ls->usetaoroutines) { 948441846f8SBarry Smith ierr = TaoComputeObjectiveAndGradient(ls->tao,x,f,g);CHKERRQ(ierr); 9490c52818fSSatish Balay } else { 9500ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9510c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 952362febeeSStefano Zampini PetscStackPush("TaoLineSearch user objective/gradient routine"); 9530c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr); 9540c52818fSSatish Balay PetscStackPop; 955362febeeSStefano Zampini } else { 956*3c859ba3SBarry Smith PetscCheck(ls->ops->computeobjective,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 957362febeeSStefano Zampini PetscStackPush("TaoLineSearch user objective routine"); 958362febeeSStefano Zampini ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 959362febeeSStefano Zampini PetscStackPop; 960*3c859ba3SBarry Smith PetscCheck(ls->ops->computegradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 961362febeeSStefano Zampini PetscStackPush("TaoLineSearch user gradient routine"); 962362febeeSStefano Zampini ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 963362febeeSStefano Zampini PetscStackPop; 964362febeeSStefano Zampini } 9650ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 9667d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 967f06e3bfaSBarry Smith } 968fbe0838dSJason Sarich ls->nfgeval++; 9690c52818fSSatish Balay PetscFunctionReturn(0); 9700c52818fSSatish Balay } 9710c52818fSSatish Balay 9720c52818fSSatish Balay /*@ 9730c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 9740c52818fSSatish Balay 9750c52818fSSatish Balay Collective on TaoLineSearch 9760c52818fSSatish Balay 9770c52818fSSatish Balay Input Parameters: 9780c52818fSSatish Balay + ls - the TaoLineSearch context 9790c52818fSSatish Balay - x - input vector 9800c52818fSSatish Balay 9810c52818fSSatish Balay Output Parameter: 9820c52818fSSatish Balay . g - gradient vector 9830c52818fSSatish Balay 98495452b02SPatrick Sanan Notes: 98595452b02SPatrick Sanan TaoComputeGradient() is typically used within line searches 9860c52818fSSatish Balay so most users would not generally call this routine themselves. 9870c52818fSSatish Balay 9880c52818fSSatish Balay Level: developer 9890c52818fSSatish Balay 9900c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 9910c52818fSSatish Balay @*/ 9920c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 9930c52818fSSatish Balay { 9940c52818fSSatish Balay PetscErrorCode ierr; 9950c52818fSSatish Balay PetscReal fdummy; 996f06e3bfaSBarry Smith 9970c52818fSSatish Balay PetscFunctionBegin; 9980c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 9990c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 10000c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,3); 10010c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 10020c52818fSSatish Balay PetscCheckSameComm(ls,1,g,3); 10030c52818fSSatish Balay if (ls->usetaoroutines) { 1004441846f8SBarry Smith ierr = TaoComputeGradient(ls->tao,x,g);CHKERRQ(ierr); 10050c52818fSSatish Balay } else { 10060ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 10070c52818fSSatish Balay PetscStackPush("TaoLineSearch user gradient routine"); 10080c52818fSSatish Balay if (ls->ops->computegradient) { 10090c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 10100c52818fSSatish Balay } else { 1011*3c859ba3SBarry Smith PetscCheck(ls->ops->computeobjectiveandgradient,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 10120c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr); 10130c52818fSSatish Balay } 10140c52818fSSatish Balay PetscStackPop; 10150ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 10160c52818fSSatish Balay } 10170c52818fSSatish Balay ls->ngeval++; 10180c52818fSSatish Balay PetscFunctionReturn(0); 10190c52818fSSatish Balay } 10200c52818fSSatish Balay 10210c52818fSSatish Balay /*@ 10220c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 10230c52818fSSatish Balay 1024441846f8SBarry Smith Collective on Tao 10250c52818fSSatish Balay 10260c52818fSSatish Balay Input Parameters: 10270c52818fSSatish Balay + ls - the TaoLineSearch context 10280c52818fSSatish Balay - x - input vector 10290c52818fSSatish Balay 1030d8d19677SJose E. Roman Output Parameters: 10310c52818fSSatish Balay + f - Objective value at X 10320c52818fSSatish Balay - gts - inner product of gradient and step direction at X 10330c52818fSSatish Balay 103495452b02SPatrick Sanan Notes: 103595452b02SPatrick Sanan TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 10360c52818fSSatish Balay so most users would not generally call this routine themselves. 10370c52818fSSatish Balay 10380c52818fSSatish Balay Level: developer 10390c52818fSSatish Balay 10400c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 10410c52818fSSatish Balay @*/ 10420c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 10430c52818fSSatish Balay { 10440c52818fSSatish Balay PetscErrorCode ierr; 10450c52818fSSatish Balay PetscFunctionBegin; 10460c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 10470c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 10480c52818fSSatish Balay PetscValidPointer(f,3); 10490c52818fSSatish Balay PetscValidPointer(gts,4); 10500c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 1051*3c859ba3SBarry Smith PetscCheck(ls->ops->computeobjectiveandgts,PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 10520ebee16dSLisandro Dalcin ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 10530c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gts routine"); 10540c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr); 10550c52818fSSatish Balay PetscStackPop; 10560ebee16dSLisandro Dalcin ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 10577d3de750SJacob Faibussowitsch ierr = PetscInfo(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 10580c52818fSSatish Balay ls->nfeval++; 10590c52818fSSatish Balay PetscFunctionReturn(0); 10600c52818fSSatish Balay } 10610c52818fSSatish Balay 10620c52818fSSatish Balay /*@ 10630c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 10640c52818fSSatish Balay 10650c52818fSSatish Balay Collective on TaoLineSearch 10660c52818fSSatish Balay 10670c52818fSSatish Balay Input Parameter: 10680c52818fSSatish Balay . ls - the TaoLineSearch context 10690c52818fSSatish Balay 1070d8d19677SJose E. Roman Output Parameters: 10710c52818fSSatish Balay + x - the new solution 10720c52818fSSatish Balay . f - the objective function value at x 10730c52818fSSatish Balay . g - the gradient at x 10740c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 10750c52818fSSatish Balay - reason - the reason why the line search terminated 10760c52818fSSatish Balay 10770c52818fSSatish Balay reason will be set to one of: 10780c52818fSSatish Balay 10790c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 10800c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 10810c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 10820c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 10830c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 10840c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 10850c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 10860c52818fSSatish Balay 10870c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 10880c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 10890c52818fSSatish Balay 1090a2b725a8SWilliam Gropp - TAOLINESEARCH_SUCCESS - successful line search 10910c52818fSSatish Balay 10920c52818fSSatish Balay Level: developer 10930c52818fSSatish Balay 10940c52818fSSatish Balay @*/ 1095e4cb33bbSBarry Smith PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 10960c52818fSSatish Balay { 10970c52818fSSatish Balay PetscErrorCode ierr; 1098f06e3bfaSBarry Smith 10990c52818fSSatish Balay PetscFunctionBegin; 11000c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11010c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 11020c52818fSSatish Balay PetscValidPointer(f,3); 11030c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 11040c52818fSSatish Balay PetscValidIntPointer(reason,6); 11050c52818fSSatish Balay if (ls->new_x) { 11060c52818fSSatish Balay ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr); 11070c52818fSSatish Balay } 11080c52818fSSatish Balay *f = ls->new_f; 11090c52818fSSatish Balay if (ls->new_g) { 11100c52818fSSatish Balay ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr); 11110c52818fSSatish Balay } 111297ab8e72SStefano Zampini if (steplength) *steplength = ls->step; 11130c52818fSSatish Balay *reason = ls->reason; 11140c52818fSSatish Balay PetscFunctionReturn(0); 11150c52818fSSatish Balay } 11160c52818fSSatish Balay 11170c52818fSSatish Balay /*@ 11180c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 11190c52818fSSatish Balay search. 11200c52818fSSatish Balay 11210c52818fSSatish Balay Not Collective 11220c52818fSSatish Balay 11230c52818fSSatish Balay Input Parameter: 11240c52818fSSatish Balay . ls - the TaoLineSearch context 11250c52818fSSatish Balay 11260c52818fSSatish Balay Output Parameter: 11270c52818fSSatish Balay . x - The initial point of the line search 11280c52818fSSatish Balay 11290c52818fSSatish Balay Level: intermediate 11300c52818fSSatish Balay @*/ 11310c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 11320c52818fSSatish Balay { 11330c52818fSSatish Balay PetscFunctionBegin; 11340c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 113597ab8e72SStefano Zampini if (x) *x = ls->start_x; 11360c52818fSSatish Balay PetscFunctionReturn(0); 11370c52818fSSatish Balay } 11380c52818fSSatish Balay 11390c52818fSSatish Balay /*@ 11400c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 11410c52818fSSatish Balay search. 11420c52818fSSatish Balay 11430c52818fSSatish Balay Not Collective 11440c52818fSSatish Balay 11450c52818fSSatish Balay Input Parameter: 11460c52818fSSatish Balay . ls - the TaoLineSearch context 11470c52818fSSatish Balay 11480c52818fSSatish Balay Output Parameter: 11490c52818fSSatish Balay . s - the step direction of the line search 11500c52818fSSatish Balay 11510c52818fSSatish Balay Level: advanced 11520c52818fSSatish Balay @*/ 11530c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 11540c52818fSSatish Balay { 11550c52818fSSatish Balay PetscFunctionBegin; 11560c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 115797ab8e72SStefano Zampini if (s) *s = ls->stepdirection; 11580c52818fSSatish Balay PetscFunctionReturn(0); 11590c52818fSSatish Balay } 11600c52818fSSatish Balay 11610c52818fSSatish Balay /*@ 11620c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 11630c52818fSSatish Balay 11640c52818fSSatish Balay Not Collective 11650c52818fSSatish Balay 11660c52818fSSatish Balay Input Parameter: 11670c52818fSSatish Balay . ls - the TaoLineSearch context 11680c52818fSSatish Balay 11690c52818fSSatish Balay Output Parameter: 11700c52818fSSatish Balay . f - the objective value at the full step length 11710c52818fSSatish Balay 11720c52818fSSatish Balay Level: developer 11730c52818fSSatish Balay @*/ 11740c52818fSSatish Balay 11750c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 11760c52818fSSatish Balay { 11770c52818fSSatish Balay PetscFunctionBegin; 11780c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 11790c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 11800c52818fSSatish Balay PetscFunctionReturn(0); 11810c52818fSSatish Balay } 11820c52818fSSatish Balay 11830c52818fSSatish Balay /*@ 11840c52818fSSatish Balay TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 11850c52818fSSatish Balay 1186441846f8SBarry Smith Logically Collective on Tao 11870c52818fSSatish Balay 11880c52818fSSatish Balay Input Parameters: 11890c52818fSSatish Balay + ls - the TaoLineSearch context 11900c52818fSSatish Balay . xl - vector of lower bounds 11910c52818fSSatish Balay - xu - vector of upper bounds 11920c52818fSSatish Balay 11930c52818fSSatish Balay Note: If the variable bounds are not set with this routine, then 1194e270355aSBarry Smith PETSC_NINFINITY and PETSC_INFINITY are assumed 11950c52818fSSatish Balay 11960c52818fSSatish Balay Level: beginner 11970c52818fSSatish Balay 11980c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 11990c52818fSSatish Balay @*/ 12000c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 12010c52818fSSatish Balay { 120297ab8e72SStefano Zampini PetscErrorCode ierr; 120397ab8e72SStefano Zampini 12040c52818fSSatish Balay PetscFunctionBegin; 12050c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 12060c52818fSSatish Balay PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 12070c52818fSSatish Balay PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 120897ab8e72SStefano Zampini ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr); 120997ab8e72SStefano Zampini ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr); 121097ab8e72SStefano Zampini ierr = VecDestroy(&ls->lower);CHKERRQ(ierr); 121197ab8e72SStefano Zampini ierr = VecDestroy(&ls->upper);CHKERRQ(ierr); 12120c52818fSSatish Balay ls->lower = xl; 12130c52818fSSatish Balay ls->upper = xu; 12140c52818fSSatish Balay ls->bounded = 1; 12150c52818fSSatish Balay PetscFunctionReturn(0); 12160c52818fSSatish Balay } 12170c52818fSSatish Balay 12180c52818fSSatish Balay /*@ 12190c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 12200c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 12210c52818fSSatish Balay 12220c52818fSSatish Balay Logically Collective on TaoLineSearch 12230c52818fSSatish Balay 12240c52818fSSatish Balay Input Parameters: 12250c52818fSSatish Balay + ls - the TaoLineSearch context 12260c52818fSSatish Balay - s - the initial step size 12270c52818fSSatish Balay 12280c52818fSSatish Balay Level: intermediate 12290c52818fSSatish Balay 12300c52818fSSatish Balay .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply() 12310c52818fSSatish Balay @*/ 12320c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 12330c52818fSSatish Balay { 12340c52818fSSatish Balay PetscFunctionBegin; 12350c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 12360c52818fSSatish Balay ls->initstep = s; 12370c52818fSSatish Balay PetscFunctionReturn(0); 12380c52818fSSatish Balay } 12390c52818fSSatish Balay 12400c52818fSSatish Balay /*@ 12410c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 12420c52818fSSatish Balay 12430c52818fSSatish Balay Not Collective 12440c52818fSSatish Balay 12450c52818fSSatish Balay Input Parameters: 12460c52818fSSatish Balay . ls - the TaoLineSearch context 12470c52818fSSatish Balay 12487a7aea1fSJed Brown Output Parameters: 12490c52818fSSatish Balay . s - the current step length 12500c52818fSSatish Balay 12510c52818fSSatish Balay Level: beginner 12520c52818fSSatish Balay 12530c52818fSSatish Balay .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply() 12540c52818fSSatish Balay @*/ 12550c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 12560c52818fSSatish Balay { 12570c52818fSSatish Balay PetscFunctionBegin; 12580c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 12590c52818fSSatish Balay *s = ls->step; 12600c52818fSSatish Balay PetscFunctionReturn(0); 12610c52818fSSatish Balay } 12620c52818fSSatish Balay 12630ebee16dSLisandro Dalcin /*@C 12640c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 12650c52818fSSatish Balay 12660c52818fSSatish Balay Not collective 12670c52818fSSatish Balay 12680c52818fSSatish Balay Input Parameters: 12690c52818fSSatish Balay + sname - name of a new user-defined solver 12700c52818fSSatish Balay - func - routine to Create method context 12710c52818fSSatish Balay 12720c52818fSSatish Balay Notes: 12730c52818fSSatish Balay TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 12740c52818fSSatish Balay 12750c52818fSSatish Balay Sample usage: 12760c52818fSSatish Balay .vb 12770c52818fSSatish Balay TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 12780c52818fSSatish Balay .ve 12790c52818fSSatish Balay 12800c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 12810c52818fSSatish Balay $ TaoLineSearchSetType(ls,"my_linesearch") 12820c52818fSSatish Balay or at runtime via the option 12830c52818fSSatish Balay $ -tao_ls_type my_linesearch 12840c52818fSSatish Balay 12850c52818fSSatish Balay Level: developer 12860c52818fSSatish Balay 12870ebee16dSLisandro Dalcin @*/ 12880c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 12890c52818fSSatish Balay { 12900c52818fSSatish Balay PetscErrorCode ierr; 129197ab8e72SStefano Zampini 12920c52818fSSatish Balay PetscFunctionBegin; 12931d36bdfdSBarry Smith ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr); 12940c52818fSSatish Balay ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr); 12950c52818fSSatish Balay PetscFunctionReturn(0); 12960c52818fSSatish Balay } 12970c52818fSSatish Balay 12980c52818fSSatish Balay /*@C 12990c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 13000c52818fSSatish Balay for all TaoLineSearch options in the database. 13010c52818fSSatish Balay 13020c52818fSSatish Balay Collective on TaoLineSearch 13030c52818fSSatish Balay 13040c52818fSSatish Balay Input Parameters: 13050c52818fSSatish Balay + ls - the TaoLineSearch solver context 13060c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 13070c52818fSSatish Balay 13080c52818fSSatish Balay Notes: 13090c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 13100c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 13110c52818fSSatish Balay 13120c52818fSSatish Balay Level: advanced 13130c52818fSSatish Balay 13140c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 13150c52818fSSatish Balay @*/ 13160c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 13170c52818fSSatish Balay { 1318f06e3bfaSBarry Smith return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 13190c52818fSSatish Balay } 13200c52818fSSatish Balay 13210c52818fSSatish Balay /*@C 13220c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 13230c52818fSSatish Balay TaoLineSearch options in the database 13240c52818fSSatish Balay 13250c52818fSSatish Balay Not Collective 13260c52818fSSatish Balay 13270c52818fSSatish Balay Input Parameters: 13280c52818fSSatish Balay . ls - the TaoLineSearch context 13290c52818fSSatish Balay 13300c52818fSSatish Balay Output Parameters: 13310c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 13320c52818fSSatish Balay 133395452b02SPatrick Sanan Notes: 133495452b02SPatrick Sanan On the fortran side, the user should pass in a string 'prefix' of 13350c52818fSSatish Balay sufficient length to hold the prefix. 13360c52818fSSatish Balay 13370c52818fSSatish Balay Level: advanced 13380c52818fSSatish Balay 13390c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 13400c52818fSSatish Balay @*/ 13410c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 13420c52818fSSatish Balay { 1343f06e3bfaSBarry Smith return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 13440c52818fSSatish Balay } 13450c52818fSSatish Balay 13460c52818fSSatish Balay /*@C 13470c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 13480c52818fSSatish Balay TaoLineSearch options in the database. 13490c52818fSSatish Balay 13500c52818fSSatish Balay Logically Collective on TaoLineSearch 13510c52818fSSatish Balay 13520c52818fSSatish Balay Input Parameters: 13530c52818fSSatish Balay + ls - the TaoLineSearch context 13540c52818fSSatish Balay - prefix - the prefix string to prepend to all TAO option requests 13550c52818fSSatish Balay 13560c52818fSSatish Balay Notes: 13570c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 13580c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 13590c52818fSSatish Balay 13600c52818fSSatish Balay For example, to distinguish between the runtime options for two 13610c52818fSSatish Balay different line searches, one could call 13620c52818fSSatish Balay .vb 13630c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 13640c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 13650c52818fSSatish Balay .ve 13660c52818fSSatish Balay 13670c52818fSSatish Balay This would enable use of different options for each system, such as 13680c52818fSSatish Balay .vb 13690c52818fSSatish Balay -sys1_tao_ls_type mt 13700c52818fSSatish Balay -sys2_tao_ls_type armijo 13710c52818fSSatish Balay .ve 13720c52818fSSatish Balay 13730c52818fSSatish Balay Level: advanced 13740c52818fSSatish Balay 13750c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 13760c52818fSSatish Balay @*/ 13770c52818fSSatish Balay 13780c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 13790c52818fSSatish Balay { 1380f06e3bfaSBarry Smith return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 13810c52818fSSatish Balay } 1382