1*0c52818fSSatish Balay #include "taosolver.h" 2*0c52818fSSatish Balay #include "taolinesearch.h" /*I "taolinesearch.h" I*/ 3*0c52818fSSatish Balay #include "tao-private/taolinesearch_impl.h" 4*0c52818fSSatish Balay 5*0c52818fSSatish Balay PetscBool TaoLineSearchInitialized = PETSC_FALSE; 6*0c52818fSSatish Balay PetscFunctionList TaoLineSearchList = PETSC_NULL; 7*0c52818fSSatish Balay 8*0c52818fSSatish Balay PetscClassId TAOLINESEARCH_CLASSID=0; 9*0c52818fSSatish Balay PetscLogEvent TaoLineSearch_ApplyEvent = 0, TaoLineSearch_EvalEvent=0; 10*0c52818fSSatish Balay 11*0c52818fSSatish Balay #undef __FUNCT__ 12*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchView" 13*0c52818fSSatish Balay /*@C 14*0c52818fSSatish Balay TaoLineSearchView - Prints information about the TaoLineSearch 15*0c52818fSSatish Balay 16*0c52818fSSatish Balay Collective on TaoLineSearch 17*0c52818fSSatish Balay 18*0c52818fSSatish Balay InputParameters: 19*0c52818fSSatish Balay + ls - the TaoSolver context 20*0c52818fSSatish Balay - viewer - visualization context 21*0c52818fSSatish Balay 22*0c52818fSSatish Balay Options Database Key: 23*0c52818fSSatish Balay . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search 24*0c52818fSSatish Balay 25*0c52818fSSatish Balay Notes: 26*0c52818fSSatish Balay The available visualization contexts include 27*0c52818fSSatish Balay + PETSC_VIEWER_STDOUT_SELF - standard output (default) 28*0c52818fSSatish Balay - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 29*0c52818fSSatish Balay output where only the first processor opens 30*0c52818fSSatish Balay the file. All other processors send their 31*0c52818fSSatish Balay data to the first processor to print. 32*0c52818fSSatish Balay 33*0c52818fSSatish Balay Level: beginner 34*0c52818fSSatish Balay 35*0c52818fSSatish Balay .seealso: PetscViewerASCIIOpen() 36*0c52818fSSatish Balay @*/ 37*0c52818fSSatish Balay 38*0c52818fSSatish Balay PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer) 39*0c52818fSSatish Balay { 40*0c52818fSSatish Balay PetscErrorCode ierr; 41*0c52818fSSatish Balay PetscBool isascii, isstring; 42*0c52818fSSatish Balay const TaoLineSearchType type; 43*0c52818fSSatish Balay PetscBool destroyviewer = PETSC_FALSE; 44*0c52818fSSatish Balay PetscFunctionBegin; 45*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 46*0c52818fSSatish Balay if (!viewer) { 47*0c52818fSSatish Balay destroyviewer = PETSC_TRUE; 48*0c52818fSSatish Balay ierr = PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer); CHKERRQ(ierr); 49*0c52818fSSatish Balay } 50*0c52818fSSatish Balay PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 51*0c52818fSSatish Balay PetscCheckSameComm(ls,1,viewer,2); 52*0c52818fSSatish Balay 53*0c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 54*0c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 55*0c52818fSSatish Balay if (isascii) { 56*0c52818fSSatish Balay if (((PetscObject)ls)->prefix) { 57*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:(%s)\n",((PetscObject)ls)->prefix);CHKERRQ(ierr); 58*0c52818fSSatish Balay } else { 59*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:\n");CHKERRQ(ierr); 60*0c52818fSSatish Balay } 61*0c52818fSSatish Balay ierr = PetscViewerASCIIPushTab(viewer); CHKERRQ(ierr); 62*0c52818fSSatish Balay ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr); 63*0c52818fSSatish Balay if (type) { 64*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"type: %s\n",type);CHKERRQ(ierr); 65*0c52818fSSatish Balay } else { 66*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"type: not set yet\n");CHKERRQ(ierr); 67*0c52818fSSatish Balay } 68*0c52818fSSatish Balay if (ls->ops->view) { 69*0c52818fSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 70*0c52818fSSatish Balay ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr); 71*0c52818fSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 72*0c52818fSSatish Balay } 73*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr); 74*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%G, rtol=%G, gtol=%G\n", ls->ftol, ls->rtol,ls->gtol);CHKERRQ(ierr); 75*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr); 76*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr); 77*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr); 78*0c52818fSSatish Balay 79*0c52818fSSatish Balay if (ls->bounded) { 80*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr); 81*0c52818fSSatish Balay } 82*0c52818fSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %D\n",(int)ls->reason); CHKERRQ(ierr); 83*0c52818fSSatish Balay ierr = PetscViewerASCIIPopTab(viewer); CHKERRQ(ierr); 84*0c52818fSSatish Balay 85*0c52818fSSatish Balay } else if (isstring) { 86*0c52818fSSatish Balay ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr); 87*0c52818fSSatish Balay ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 88*0c52818fSSatish Balay } 89*0c52818fSSatish Balay if (destroyviewer) { 90*0c52818fSSatish Balay ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 91*0c52818fSSatish Balay } 92*0c52818fSSatish Balay PetscFunctionReturn(0); 93*0c52818fSSatish Balay 94*0c52818fSSatish Balay } 95*0c52818fSSatish Balay 96*0c52818fSSatish Balay 97*0c52818fSSatish Balay #undef __FUNCT__ 98*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchCreate" 99*0c52818fSSatish Balay /*@C 100*0c52818fSSatish Balay TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 101*0c52818fSSatish Balay line-searches will automatically create one. 102*0c52818fSSatish Balay 103*0c52818fSSatish Balay Collective on MPI_Comm 104*0c52818fSSatish Balay 105*0c52818fSSatish Balay Input Parameter: 106*0c52818fSSatish Balay . comm - MPI communicator 107*0c52818fSSatish Balay 108*0c52818fSSatish Balay Output Parameter: 109*0c52818fSSatish Balay . newls - the new TaoLineSearch context 110*0c52818fSSatish Balay 111*0c52818fSSatish Balay Available methods include: 112*0c52818fSSatish Balay + more-thuente 113*0c52818fSSatish Balay . gpcg 114*0c52818fSSatish Balay - unit - Do not perform any line search 115*0c52818fSSatish Balay 116*0c52818fSSatish Balay 117*0c52818fSSatish Balay Options Database Keys: 118*0c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 119*0c52818fSSatish Balay 120*0c52818fSSatish Balay Level: beginner 121*0c52818fSSatish Balay 122*0c52818fSSatish Balay .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy() 123*0c52818fSSatish Balay @*/ 124*0c52818fSSatish Balay 125*0c52818fSSatish Balay PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 126*0c52818fSSatish Balay { 127*0c52818fSSatish Balay PetscErrorCode ierr; 128*0c52818fSSatish Balay TaoLineSearch ls; 129*0c52818fSSatish Balay 130*0c52818fSSatish Balay PetscFunctionBegin; 131*0c52818fSSatish Balay PetscValidPointer(newls,2); 132*0c52818fSSatish Balay *newls = PETSC_NULL; 133*0c52818fSSatish Balay 134*0c52818fSSatish Balay #ifndef PETSC_USE_DYNAMIC_LIBRARIES 135*0c52818fSSatish Balay ierr = TaoLineSearchInitializePackage(); CHKERRQ(ierr); 136*0c52818fSSatish Balay #endif 137*0c52818fSSatish Balay 138*0c52818fSSatish Balay ierr = PetscHeaderCreate(ls,_p_TaoLineSearch,struct _TaoLineSearchOps, 139*0c52818fSSatish Balay TAOLINESEARCH_CLASSID, "TaoLineSearch",0,0, 140*0c52818fSSatish Balay comm,TaoLineSearchDestroy, TaoLineSearchView); 141*0c52818fSSatish Balay CHKERRQ(ierr); 142*0c52818fSSatish Balay ls->bounded = 0; 143*0c52818fSSatish Balay ls->max_funcs=30; 144*0c52818fSSatish Balay ls->ftol = 0.0001; 145*0c52818fSSatish Balay ls->gtol = 0.9; 146*0c52818fSSatish Balay ls->rtol = 1.0e-10; 147*0c52818fSSatish Balay ls->stepmin=1.0e-20; 148*0c52818fSSatish Balay ls->stepmax=1.0e+20; 149*0c52818fSSatish Balay ls->step=1.0; 150*0c52818fSSatish Balay ls->nfeval=0; 151*0c52818fSSatish Balay ls->ngeval=0; 152*0c52818fSSatish Balay ls->nfgeval=0; 153*0c52818fSSatish Balay 154*0c52818fSSatish Balay ls->ops->computeobjective=0; 155*0c52818fSSatish Balay ls->ops->computegradient=0; 156*0c52818fSSatish Balay ls->ops->computeobjectiveandgradient=0; 157*0c52818fSSatish Balay ls->ops->computeobjectiveandgts=0; 158*0c52818fSSatish Balay ls->ops->setup=0; 159*0c52818fSSatish Balay ls->ops->apply=0; 160*0c52818fSSatish Balay ls->ops->view=0; 161*0c52818fSSatish Balay ls->ops->setfromoptions=0; 162*0c52818fSSatish Balay ls->ops->reset=0; 163*0c52818fSSatish Balay ls->ops->destroy=0; 164*0c52818fSSatish Balay ls->setupcalled=PETSC_FALSE; 165*0c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 166*0c52818fSSatish Balay *newls = ls; 167*0c52818fSSatish Balay PetscFunctionReturn(0); 168*0c52818fSSatish Balay } 169*0c52818fSSatish Balay 170*0c52818fSSatish Balay 171*0c52818fSSatish Balay #undef __FUNCT__ 172*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetUp" 173*0c52818fSSatish Balay /*@ 174*0c52818fSSatish Balay TaoLineSearchSetUp - Sets up the internal data structures for the later use 175*0c52818fSSatish Balay of a Tao solver 176*0c52818fSSatish Balay 177*0c52818fSSatish Balay Collective on ls 178*0c52818fSSatish Balay 179*0c52818fSSatish Balay Input Parameters: 180*0c52818fSSatish Balay . ls - the TaoLineSearch context 181*0c52818fSSatish Balay 182*0c52818fSSatish Balay Notes: 183*0c52818fSSatish Balay The user will not need to explicitly call TaoLineSearchSetUp(), as it will 184*0c52818fSSatish Balay automatically be called in TaoLineSearchSolve(). However, if the user 185*0c52818fSSatish Balay desires to call it explicitly, it should come after TaoLineSearchCreate() 186*0c52818fSSatish Balay but before TaoLineSearchApply(). 187*0c52818fSSatish Balay 188*0c52818fSSatish Balay Level: developer 189*0c52818fSSatish Balay 190*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 191*0c52818fSSatish Balay @*/ 192*0c52818fSSatish Balay 193*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 194*0c52818fSSatish Balay { 195*0c52818fSSatish Balay PetscErrorCode ierr; 196*0c52818fSSatish Balay const char *default_type=TAOLINESEARCH_MT; 197*0c52818fSSatish Balay PetscBool flg; 198*0c52818fSSatish Balay PetscFunctionBegin; 199*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 200*0c52818fSSatish Balay if (ls->setupcalled) PetscFunctionReturn(0); 201*0c52818fSSatish Balay if (!((PetscObject)ls)->type_name) { 202*0c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,default_type); CHKERRQ(ierr); 203*0c52818fSSatish Balay } 204*0c52818fSSatish Balay if (ls->ops->setup) { 205*0c52818fSSatish Balay ierr = (*ls->ops->setup)(ls); CHKERRQ(ierr); 206*0c52818fSSatish Balay } 207*0c52818fSSatish Balay if (ls->usetaoroutines) { 208*0c52818fSSatish Balay ierr = TaoIsObjectiveDefined(ls->taosolver,&flg); CHKERRQ(ierr); 209*0c52818fSSatish Balay ls->hasobjective = flg; 210*0c52818fSSatish Balay ierr = TaoIsGradientDefined(ls->taosolver,&flg); CHKERRQ(ierr); 211*0c52818fSSatish Balay ls->hasgradient = flg; 212*0c52818fSSatish Balay ierr = TaoIsObjectiveAndGradientDefined(ls->taosolver,&flg); CHKERRQ(ierr); 213*0c52818fSSatish Balay ls->hasobjectiveandgradient = flg; 214*0c52818fSSatish Balay } else { 215*0c52818fSSatish Balay if (ls->ops->computeobjective) { 216*0c52818fSSatish Balay ls->hasobjective = PETSC_TRUE; 217*0c52818fSSatish Balay } else { 218*0c52818fSSatish Balay ls->hasobjective = PETSC_FALSE; 219*0c52818fSSatish Balay } 220*0c52818fSSatish Balay if (ls->ops->computegradient) { 221*0c52818fSSatish Balay ls->hasgradient = PETSC_TRUE; 222*0c52818fSSatish Balay } else { 223*0c52818fSSatish Balay ls->hasgradient = PETSC_FALSE; 224*0c52818fSSatish Balay } 225*0c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 226*0c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_TRUE; 227*0c52818fSSatish Balay } else { 228*0c52818fSSatish Balay ls->hasobjectiveandgradient = PETSC_FALSE; 229*0c52818fSSatish Balay } 230*0c52818fSSatish Balay } 231*0c52818fSSatish Balay ls->setupcalled = PETSC_TRUE; 232*0c52818fSSatish Balay PetscFunctionReturn(0); 233*0c52818fSSatish Balay } 234*0c52818fSSatish Balay 235*0c52818fSSatish Balay #undef __FUNCT__ 236*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchReset" 237*0c52818fSSatish Balay /*@ 238*0c52818fSSatish Balay TaoLineSearchReset - Some line searches may carry state information 239*0c52818fSSatish Balay from one TaoLineSearchApply() to the next. This function resets this 240*0c52818fSSatish Balay state information. 241*0c52818fSSatish Balay 242*0c52818fSSatish Balay Collective on TaoLineSearch 243*0c52818fSSatish Balay 244*0c52818fSSatish Balay Input Parameter: 245*0c52818fSSatish Balay . ls - the TaoLineSearch context 246*0c52818fSSatish Balay 247*0c52818fSSatish Balay Level: developer 248*0c52818fSSatish Balay 249*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 250*0c52818fSSatish Balay @*/ 251*0c52818fSSatish Balay PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 252*0c52818fSSatish Balay { 253*0c52818fSSatish Balay PetscErrorCode ierr; 254*0c52818fSSatish Balay PetscFunctionBegin; 255*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 256*0c52818fSSatish Balay if (ls->ops->reset) { 257*0c52818fSSatish Balay ierr= (*ls->ops->reset)(ls); CHKERRQ(ierr); 258*0c52818fSSatish Balay } 259*0c52818fSSatish Balay PetscFunctionReturn(0); 260*0c52818fSSatish Balay } 261*0c52818fSSatish Balay #undef __FUNCT__ 262*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchDestroy" 263*0c52818fSSatish Balay /*@ 264*0c52818fSSatish Balay TaoLineSearchDestroy - Destroys the TAO context that was created with 265*0c52818fSSatish Balay TaoLineSearchCreate() 266*0c52818fSSatish Balay 267*0c52818fSSatish Balay Collective on TaoLineSearch 268*0c52818fSSatish Balay 269*0c52818fSSatish Balay Input Parameter 270*0c52818fSSatish Balay . ls - the TaoLineSearch context 271*0c52818fSSatish Balay 272*0c52818fSSatish Balay Level: beginner 273*0c52818fSSatish Balay 274*0c52818fSSatish Balay .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 275*0c52818fSSatish Balay @*/ 276*0c52818fSSatish Balay PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 277*0c52818fSSatish Balay { 278*0c52818fSSatish Balay PetscErrorCode ierr; 279*0c52818fSSatish Balay PetscFunctionBegin; 280*0c52818fSSatish Balay if (!*ls) PetscFunctionReturn(0); 281*0c52818fSSatish Balay PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 282*0c52818fSSatish Balay if (--((PetscObject)*ls)->refct > 0) {*ls=0; PetscFunctionReturn(0);} 283*0c52818fSSatish Balay if ((*ls)->stepdirection!=PETSC_NULL) { 284*0c52818fSSatish Balay ierr = PetscObjectDereference((PetscObject)(*ls)->stepdirection);CHKERRQ(ierr); 285*0c52818fSSatish Balay (*ls)->stepdirection = PETSC_NULL; 286*0c52818fSSatish Balay } 287*0c52818fSSatish Balay if ((*ls)->start_x != PETSC_NULL) { 288*0c52818fSSatish Balay ierr = PetscObjectDereference((PetscObject)(*ls)->start_x); CHKERRQ(ierr); 289*0c52818fSSatish Balay (*ls)->start_x = PETSC_NULL; 290*0c52818fSSatish Balay } 291*0c52818fSSatish Balay if ((*ls)->ops->destroy) { 292*0c52818fSSatish Balay ierr = (*(*ls)->ops->destroy)(*ls); CHKERRQ(ierr); 293*0c52818fSSatish Balay } 294*0c52818fSSatish Balay ierr = PetscHeaderDestroy(ls); CHKERRQ(ierr); 295*0c52818fSSatish Balay PetscFunctionReturn(0); 296*0c52818fSSatish Balay } 297*0c52818fSSatish Balay 298*0c52818fSSatish Balay 299*0c52818fSSatish Balay #undef __FUNCT__ 300*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchApply" 301*0c52818fSSatish Balay /*@ 302*0c52818fSSatish Balay TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 303*0c52818fSSatish Balay 304*0c52818fSSatish Balay Collective on TaoLineSearch 305*0c52818fSSatish Balay 306*0c52818fSSatish Balay Input Parameters: 307*0c52818fSSatish Balay + ls - the TaoSolver context 308*0c52818fSSatish Balay . x - The current solution (on output x contains the new solution determined by the line search) 309*0c52818fSSatish Balay . f - objective function value at current solution (on output contains the objective function value at new solution) 310*0c52818fSSatish Balay . g - gradient evaluated at x (on output contains the gradient at new solution) 311*0c52818fSSatish Balay - s - search direction 312*0c52818fSSatish Balay 313*0c52818fSSatish Balay Output Parameters: 314*0c52818fSSatish Balay + x - new solution 315*0c52818fSSatish Balay . f - objective function value at x 316*0c52818fSSatish Balay . g - gradient vector at x 317*0c52818fSSatish Balay . steplength - scalar multiplier of s used ( x = x0 + steplength * x ) 318*0c52818fSSatish Balay - reason - reason why the line-search stopped 319*0c52818fSSatish Balay 320*0c52818fSSatish Balay reason will be set to one of: 321*0c52818fSSatish Balay 322*0c52818fSSatish Balay + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 323*0c52818fSSatish Balay . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 324*0c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 325*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 326*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 327*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 328*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 329*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 330*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 331*0c52818fSSatish Balay - TAOLINESEARCH_SUCCESS - successful line search 332*0c52818fSSatish Balay 333*0c52818fSSatish Balay 334*0c52818fSSatish Balay Note: 335*0c52818fSSatish Balay The algorithm developer must set up the TaoLineSearch with calls to 336*0c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoSolverRoutines() 337*0c52818fSSatish Balay 338*0c52818fSSatish Balay Note: 339*0c52818fSSatish Balay You may or may not need to follow this with a call to 340*0c52818fSSatish Balay TaoAddLineSearchCounts(), depending on whether you want these 341*0c52818fSSatish Balay evaluations to count toward the total function/gradient evaluations. 342*0c52818fSSatish Balay 343*0c52818fSSatish Balay Level: beginner 344*0c52818fSSatish Balay 345*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts() 346*0c52818fSSatish Balay @*/ 347*0c52818fSSatish Balay 348*0c52818fSSatish Balay PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchTerminationReason *reason) 349*0c52818fSSatish Balay { 350*0c52818fSSatish Balay PetscErrorCode ierr; 351*0c52818fSSatish Balay PetscViewer viewer; 352*0c52818fSSatish Balay PetscInt low1,low2,low3,high1,high2,high3; 353*0c52818fSSatish Balay PetscBool flg; 354*0c52818fSSatish Balay char filename[PETSC_MAX_PATH_LEN]; 355*0c52818fSSatish Balay 356*0c52818fSSatish Balay PetscFunctionBegin; 357*0c52818fSSatish Balay *reason = TAOLINESEARCH_CONTINUE_ITERATING; 358*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 359*0c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 360*0c52818fSSatish Balay PetscValidScalarPointer(f,3); 361*0c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 362*0c52818fSSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,5); 363*0c52818fSSatish Balay PetscValidPointer(reason,7); 364*0c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 365*0c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,g,4); 366*0c52818fSSatish Balay PetscCheckSameTypeAndComm(x,2,s,5); 367*0c52818fSSatish Balay ierr = VecGetOwnershipRange(x, &low1, &high1); CHKERRQ(ierr); 368*0c52818fSSatish Balay ierr = VecGetOwnershipRange(g, &low2, &high2); CHKERRQ(ierr); 369*0c52818fSSatish Balay ierr = VecGetOwnershipRange(s, &low3, &high3); CHKERRQ(ierr); 370*0c52818fSSatish Balay if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) { 371*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths"); 372*0c52818fSSatish Balay } 373*0c52818fSSatish Balay 374*0c52818fSSatish Balay 375*0c52818fSSatish Balay if (ls->stepdirection) { 376*0c52818fSSatish Balay ierr = PetscObjectDereference((PetscObject)(s)); CHKERRQ(ierr); 377*0c52818fSSatish Balay } 378*0c52818fSSatish Balay ls->stepdirection = s; 379*0c52818fSSatish Balay ierr = PetscObjectReference((PetscObject)s); CHKERRQ(ierr); 380*0c52818fSSatish Balay 381*0c52818fSSatish Balay ierr = TaoLineSearchSetUp(ls); CHKERRQ(ierr); 382*0c52818fSSatish Balay if (!ls->ops->apply) { 383*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine"); 384*0c52818fSSatish Balay } 385*0c52818fSSatish Balay ls->nfeval=0; 386*0c52818fSSatish Balay ls->ngeval=0; 387*0c52818fSSatish Balay ls->nfgeval=0; 388*0c52818fSSatish Balay /* Check parameter values */ 389*0c52818fSSatish Balay if (ls->ftol < 0.0) { 390*0c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: ftol (%G) < 0\n",ls->ftol); CHKERRQ(ierr); 391*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 392*0c52818fSSatish Balay } 393*0c52818fSSatish Balay if (ls->rtol < 0.0) { 394*0c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: rtol (%G) < 0\n",ls->rtol); CHKERRQ(ierr); 395*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 396*0c52818fSSatish Balay } 397*0c52818fSSatish Balay 398*0c52818fSSatish Balay if (ls->gtol < 0.0) { 399*0c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: gtol (%G) < 0\n",ls->gtol); CHKERRQ(ierr); 400*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 401*0c52818fSSatish Balay } 402*0c52818fSSatish Balay if (ls->stepmin < 0.0) { 403*0c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%G) < 0\n",ls->stepmin); CHKERRQ(ierr); 404*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 405*0c52818fSSatish Balay } 406*0c52818fSSatish Balay if (ls->stepmax < ls->stepmin) { 407*0c52818fSSatish Balay ierr = PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%G) > stepmax (%G)\n",ls->stepmin,ls->stepmax); CHKERRQ(ierr); 408*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 409*0c52818fSSatish Balay } 410*0c52818fSSatish Balay if (ls->max_funcs < 0) { 411*0c52818fSSatish Balay ierr = PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs); CHKERRQ(ierr); 412*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 413*0c52818fSSatish Balay } 414*0c52818fSSatish Balay if (PetscIsInfOrNanReal(*f)) { 415*0c52818fSSatish Balay ierr = PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%G)\n",*f); CHKERRQ(ierr); 416*0c52818fSSatish Balay *reason=TAOLINESEARCH_FAILED_INFORNAN; 417*0c52818fSSatish Balay } 418*0c52818fSSatish Balay 419*0c52818fSSatish Balay 420*0c52818fSSatish Balay if (x != ls->start_x) { 421*0c52818fSSatish Balay ierr = PetscObjectReference((PetscObject)x); 422*0c52818fSSatish Balay if (ls->start_x) { 423*0c52818fSSatish Balay ierr = VecDestroy(&ls->start_x); CHKERRQ(ierr); 424*0c52818fSSatish Balay } 425*0c52818fSSatish Balay ls->start_x = x; 426*0c52818fSSatish Balay } 427*0c52818fSSatish Balay 428*0c52818fSSatish Balay 429*0c52818fSSatish Balay 430*0c52818fSSatish Balay 431*0c52818fSSatish Balay ierr = PetscLogEventBegin(TaoLineSearch_ApplyEvent,ls,0,0,0); CHKERRQ(ierr); 432*0c52818fSSatish Balay ierr = (*ls->ops->apply)(ls,x,f,g,s); CHKERRQ(ierr); 433*0c52818fSSatish Balay ierr = PetscLogEventEnd(TaoLineSearch_ApplyEvent, ls, 0,0,0); CHKERRQ(ierr); 434*0c52818fSSatish Balay *reason=ls->reason; 435*0c52818fSSatish Balay ls->new_f = *f; 436*0c52818fSSatish Balay 437*0c52818fSSatish Balay if (steplength) { 438*0c52818fSSatish Balay *steplength=ls->step; 439*0c52818fSSatish Balay } 440*0c52818fSSatish Balay 441*0c52818fSSatish Balay ierr = PetscOptionsGetString(((PetscObject)ls)->prefix,"-tao_ls_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 442*0c52818fSSatish Balay if (ls->viewls && !PetscPreLoadingOn) { 443*0c52818fSSatish Balay ierr = PetscViewerASCIIOpen(((PetscObject)ls)->comm,filename,&viewer); CHKERRQ(ierr); 444*0c52818fSSatish Balay ierr = TaoLineSearchView(ls,viewer); CHKERRQ(ierr); 445*0c52818fSSatish Balay ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 446*0c52818fSSatish Balay } 447*0c52818fSSatish Balay PetscFunctionReturn(0); 448*0c52818fSSatish Balay } 449*0c52818fSSatish Balay 450*0c52818fSSatish Balay #undef __FUNCT__ 451*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetType" 452*0c52818fSSatish Balay /*@C 453*0c52818fSSatish Balay TaoLineSearchSetType - Sets the algorithm used in a line search 454*0c52818fSSatish Balay 455*0c52818fSSatish Balay Collective on TaoLineSearch 456*0c52818fSSatish Balay 457*0c52818fSSatish Balay Input Parameters: 458*0c52818fSSatish Balay + ls - the TaoLineSearch context 459*0c52818fSSatish Balay - type - a known method 460*0c52818fSSatish Balay 461*0c52818fSSatish Balay 462*0c52818fSSatish Balay Available methods include: 463*0c52818fSSatish Balay + more-thuente 464*0c52818fSSatish Balay . gpcg 465*0c52818fSSatish Balay - unit - Do not perform any line search 466*0c52818fSSatish Balay 467*0c52818fSSatish Balay 468*0c52818fSSatish Balay Options Database Keys: 469*0c52818fSSatish Balay . -tao_ls_type - select which method TAO should use 470*0c52818fSSatish Balay 471*0c52818fSSatish Balay Level: beginner 472*0c52818fSSatish Balay 473*0c52818fSSatish Balay 474*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply() 475*0c52818fSSatish Balay 476*0c52818fSSatish Balay @*/ 477*0c52818fSSatish Balay 478*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, const TaoLineSearchType type) 479*0c52818fSSatish Balay { 480*0c52818fSSatish Balay PetscErrorCode ierr; 481*0c52818fSSatish Balay PetscErrorCode (*r)(TaoLineSearch); 482*0c52818fSSatish Balay PetscBool flg; 483*0c52818fSSatish Balay 484*0c52818fSSatish Balay PetscFunctionBegin; 485*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 486*0c52818fSSatish Balay PetscValidCharPointer(type,2); 487*0c52818fSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg); CHKERRQ(ierr); 488*0c52818fSSatish Balay if (flg) PetscFunctionReturn(0); 489*0c52818fSSatish Balay 490*0c52818fSSatish Balay ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r); CHKERRQ(ierr); 491*0c52818fSSatish Balay if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 492*0c52818fSSatish Balay if (ls->ops->destroy) { 493*0c52818fSSatish Balay ierr = (*(ls)->ops->destroy)(ls); CHKERRQ(ierr); 494*0c52818fSSatish Balay } 495*0c52818fSSatish Balay ls->max_funcs=30; 496*0c52818fSSatish Balay ls->ftol = 0.0001; 497*0c52818fSSatish Balay ls->gtol = 0.9; 498*0c52818fSSatish Balay ls->rtol = 1.0e-10; 499*0c52818fSSatish Balay ls->stepmin=1.0e-20; 500*0c52818fSSatish Balay ls->stepmax=1.0e+20; 501*0c52818fSSatish Balay 502*0c52818fSSatish Balay 503*0c52818fSSatish Balay ls->nfeval=0; 504*0c52818fSSatish Balay ls->ngeval=0; 505*0c52818fSSatish Balay ls->nfgeval=0; 506*0c52818fSSatish Balay ls->ops->setup=0; 507*0c52818fSSatish Balay ls->ops->apply=0; 508*0c52818fSSatish Balay ls->ops->view=0; 509*0c52818fSSatish Balay ls->ops->setfromoptions=0; 510*0c52818fSSatish Balay ls->ops->destroy=0; 511*0c52818fSSatish Balay ls->setupcalled = PETSC_FALSE; 512*0c52818fSSatish Balay ierr = (*r)(ls); CHKERRQ(ierr); 513*0c52818fSSatish Balay ierr = PetscObjectChangeTypeName((PetscObject)ls, type); CHKERRQ(ierr); 514*0c52818fSSatish Balay PetscFunctionReturn(0); 515*0c52818fSSatish Balay } 516*0c52818fSSatish Balay 517*0c52818fSSatish Balay #undef __FUNCT__ 518*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetFromOptions" 519*0c52818fSSatish Balay /*@ 520*0c52818fSSatish Balay TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 521*0c52818fSSatish Balay options. 522*0c52818fSSatish Balay 523*0c52818fSSatish Balay Collective on TaoLineSearch 524*0c52818fSSatish Balay 525*0c52818fSSatish Balay Input Paremeter: 526*0c52818fSSatish Balay . ls - the TaoLineSearch context 527*0c52818fSSatish Balay 528*0c52818fSSatish Balay Options Database Keys: 529*0c52818fSSatish Balay + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 530*0c52818fSSatish Balay . -tao_ls_ftol <tol> - tolerance for sufficient decrease 531*0c52818fSSatish Balay . -tao_ls_gtol <tol> - tolerance for curvature condition 532*0c52818fSSatish Balay . -tao_ls_rtol <tol> - relative tolerance for acceptable step 533*0c52818fSSatish Balay . -tao_ls_stepmin <step> - minimum steplength allowed 534*0c52818fSSatish Balay . -tao_ls_stepmax <step> - maximum steplength allowed 535*0c52818fSSatish Balay . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 536*0c52818fSSatish Balay - -tao_ls_view - display line-search results to standard output 537*0c52818fSSatish Balay 538*0c52818fSSatish Balay Level: beginner 539*0c52818fSSatish Balay @*/ 540*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 541*0c52818fSSatish Balay { 542*0c52818fSSatish Balay PetscErrorCode ierr; 543*0c52818fSSatish Balay 544*0c52818fSSatish Balay const char *default_type=TAOLINESEARCH_MT; 545*0c52818fSSatish Balay char type[256]; 546*0c52818fSSatish Balay PetscBool flg; 547*0c52818fSSatish Balay PetscFunctionBegin; 548*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 549*0c52818fSSatish Balay 550*0c52818fSSatish Balay ierr = PetscObjectOptionsBegin((PetscObject)ls); CHKERRQ(ierr); 551*0c52818fSSatish Balay { 552*0c52818fSSatish Balay if (!TaoLineSearchInitialized) { 553*0c52818fSSatish Balay ierr = TaoLineSearchInitializePackage(); CHKERRQ(ierr); 554*0c52818fSSatish Balay } 555*0c52818fSSatish Balay if (((PetscObject)ls)->type_name) { 556*0c52818fSSatish Balay default_type = ((PetscObject)ls)->type_name; 557*0c52818fSSatish Balay } 558*0c52818fSSatish Balay /* Check for type from options */ 559*0c52818fSSatish Balay ierr = PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg); CHKERRQ(ierr); 560*0c52818fSSatish Balay if (flg) { 561*0c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,type); CHKERRQ(ierr); 562*0c52818fSSatish Balay } else if (!((PetscObject)ls)->type_name) { 563*0c52818fSSatish Balay ierr = TaoLineSearchSetType(ls,default_type); 564*0c52818fSSatish Balay } 565*0c52818fSSatish Balay 566*0c52818fSSatish Balay ierr = PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search", 567*0c52818fSSatish Balay "",ls->max_funcs,&ls->max_funcs,0);CHKERRQ(ierr); 568*0c52818fSSatish Balay ierr = PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","", 569*0c52818fSSatish Balay ls->ftol,&ls->ftol,0);CHKERRQ(ierr); 570*0c52818fSSatish Balay ierr = PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","", 571*0c52818fSSatish Balay ls->gtol,&ls->gtol,0);CHKERRQ(ierr); 572*0c52818fSSatish Balay ierr = PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","", 573*0c52818fSSatish Balay ls->rtol,&ls->rtol,0);CHKERRQ(ierr); 574*0c52818fSSatish Balay ierr = PetscOptionsReal("-tao_ls_stepmin","lower bound for step","", 575*0c52818fSSatish Balay ls->stepmin,&ls->stepmin,0);CHKERRQ(ierr); 576*0c52818fSSatish Balay ierr = PetscOptionsReal("-tao_ls_stepmax","upper bound for step","", 577*0c52818fSSatish Balay ls->stepmax,&ls->stepmax,0);CHKERRQ(ierr); 578*0c52818fSSatish Balay ierr = PetscOptionsBool("-tao_ls_view","view TaoLineSearch info after each line search has completed","TaoLineSearchView",PETSC_FALSE,&ls->viewls,PETSC_NULL);CHKERRQ(ierr); 579*0c52818fSSatish Balay 580*0c52818fSSatish Balay if (ls->ops->setfromoptions) { 581*0c52818fSSatish Balay ierr = (*ls->ops->setfromoptions)(ls); CHKERRQ(ierr); 582*0c52818fSSatish Balay } 583*0c52818fSSatish Balay } 584*0c52818fSSatish Balay ierr = PetscOptionsEnd(); CHKERRQ(ierr); 585*0c52818fSSatish Balay 586*0c52818fSSatish Balay 587*0c52818fSSatish Balay PetscFunctionReturn(0); 588*0c52818fSSatish Balay } 589*0c52818fSSatish Balay 590*0c52818fSSatish Balay #undef __FUNCT__ 591*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetType" 592*0c52818fSSatish Balay /*@C 593*0c52818fSSatish Balay TaoLineSearchGetType - Gets the current line search algorithm 594*0c52818fSSatish Balay 595*0c52818fSSatish Balay Not Collective 596*0c52818fSSatish Balay 597*0c52818fSSatish Balay Input Parameter: 598*0c52818fSSatish Balay . ls - the TaoLineSearch context 599*0c52818fSSatish Balay 600*0c52818fSSatish Balay Output Paramter: 601*0c52818fSSatish Balay . type - the line search algorithm in effect 602*0c52818fSSatish Balay 603*0c52818fSSatish Balay Level: developer 604*0c52818fSSatish Balay 605*0c52818fSSatish Balay @*/ 606*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, const TaoLineSearchType *type) 607*0c52818fSSatish Balay { 608*0c52818fSSatish Balay PetscFunctionBegin; 609*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 610*0c52818fSSatish Balay PetscValidPointer(type,2); 611*0c52818fSSatish Balay *type = ((PetscObject)ls)->type_name; 612*0c52818fSSatish Balay PetscFunctionReturn(0); 613*0c52818fSSatish Balay } 614*0c52818fSSatish Balay 615*0c52818fSSatish Balay #undef __FUNCT__ 616*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetNumberFunctionEvaluations" 617*0c52818fSSatish Balay /*@ 618*0c52818fSSatish Balay TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation 619*0c52818fSSatish Balay routines used by the line search in last application (not cumulative). 620*0c52818fSSatish Balay 621*0c52818fSSatish Balay Not Collective 622*0c52818fSSatish Balay 623*0c52818fSSatish Balay Input Parameter: 624*0c52818fSSatish Balay . ls - the TaoLineSearch context 625*0c52818fSSatish Balay 626*0c52818fSSatish Balay Output Parameters: 627*0c52818fSSatish Balay + nfeval - number of function evaluations 628*0c52818fSSatish Balay . ngeval - number of gradient evaluations 629*0c52818fSSatish Balay - nfgeval - number of function/gradient evaluations 630*0c52818fSSatish Balay 631*0c52818fSSatish Balay Level: intermediate 632*0c52818fSSatish Balay 633*0c52818fSSatish Balay Note: 634*0c52818fSSatish Balay If the line search is using the TaoSolver objective and gradient 635*0c52818fSSatish Balay routines directly (see TaoLineSearchUseTaoSolverRoutines()), then TAO 636*0c52818fSSatish Balay is already counting the number of evaluations. 637*0c52818fSSatish Balay 638*0c52818fSSatish Balay @*/ 639*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval) 640*0c52818fSSatish Balay { 641*0c52818fSSatish Balay PetscFunctionBegin; 642*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 643*0c52818fSSatish Balay *nfeval = ls->nfeval; 644*0c52818fSSatish Balay *ngeval = ls->ngeval; 645*0c52818fSSatish Balay *nfgeval = ls->nfgeval; 646*0c52818fSSatish Balay PetscFunctionReturn(0); 647*0c52818fSSatish Balay } 648*0c52818fSSatish Balay 649*0c52818fSSatish Balay #undef __FUNCT__ 650*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchIsUsingTaoSolverRoutines" 651*0c52818fSSatish Balay /*@ 652*0c52818fSSatish Balay TaoLineSearchIsUsingTaoSolverRoutines - Checks whether the line search is using 653*0c52818fSSatish Balay TaoSolver evaluation routines. 654*0c52818fSSatish Balay 655*0c52818fSSatish Balay Not Collective 656*0c52818fSSatish Balay 657*0c52818fSSatish Balay Input Parameter: 658*0c52818fSSatish Balay . ls - the TaoLineSearch context 659*0c52818fSSatish Balay 660*0c52818fSSatish Balay Output Parameter: 661*0c52818fSSatish Balay . flg - PETSC_TRUE if the line search is using TaoSolver evaluation routines, 662*0c52818fSSatish Balay otherwise PETSC_FALSE 663*0c52818fSSatish Balay 664*0c52818fSSatish Balay Level: developer 665*0c52818fSSatish Balay @*/ 666*0c52818fSSatish Balay PetscErrorCode TaoLineSearchIsUsingTaoSolverRoutines(TaoLineSearch ls, PetscBool *flg) 667*0c52818fSSatish Balay { 668*0c52818fSSatish Balay PetscFunctionBegin; 669*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 670*0c52818fSSatish Balay *flg = (ls->usetaoroutines); 671*0c52818fSSatish Balay PetscFunctionReturn(0); 672*0c52818fSSatish Balay } 673*0c52818fSSatish Balay 674*0c52818fSSatish Balay #undef __FUNCT__ 675*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetObjectiveRoutine" 676*0c52818fSSatish Balay /*@C 677*0c52818fSSatish Balay TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search 678*0c52818fSSatish Balay 679*0c52818fSSatish Balay Logically Collective on TaoLineSearch 680*0c52818fSSatish Balay 681*0c52818fSSatish Balay Input Parameter: 682*0c52818fSSatish Balay + ls - the TaoLineSearch context 683*0c52818fSSatish Balay . func - the objective function evaluation routine 684*0c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 685*0c52818fSSatish Balay 686*0c52818fSSatish Balay Calling sequence of func: 687*0c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx); 688*0c52818fSSatish Balay 689*0c52818fSSatish Balay + x - input vector 690*0c52818fSSatish Balay . f - function value 691*0c52818fSSatish Balay - ctx (optional) user-defined context 692*0c52818fSSatish Balay 693*0c52818fSSatish Balay Level: beginner 694*0c52818fSSatish Balay 695*0c52818fSSatish Balay Note: 696*0c52818fSSatish Balay Use this routine only if you want the line search objective 697*0c52818fSSatish Balay evaluation routine to be different from the TaoSolver's objective 698*0c52818fSSatish Balay evaluation routine. If you use this routine you must also set 699*0c52818fSSatish Balay the line search gradient and/or function/gradient routine. 700*0c52818fSSatish Balay 701*0c52818fSSatish Balay Note: 702*0c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 703*0c52818fSSatish Balay line search, application programmers should be wary of overriding the 704*0c52818fSSatish Balay default objective routine. 705*0c52818fSSatish Balay 706*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoSolverRoutines() 707*0c52818fSSatish Balay @*/ 708*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx) 709*0c52818fSSatish Balay { 710*0c52818fSSatish Balay PetscFunctionBegin; 711*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 712*0c52818fSSatish Balay 713*0c52818fSSatish Balay ls->ops->computeobjective=func; 714*0c52818fSSatish Balay if (ctx) ls->userctx_func=ctx; 715*0c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 716*0c52818fSSatish Balay PetscFunctionReturn(0); 717*0c52818fSSatish Balay 718*0c52818fSSatish Balay 719*0c52818fSSatish Balay } 720*0c52818fSSatish Balay 721*0c52818fSSatish Balay 722*0c52818fSSatish Balay #undef __FUNCT__ 723*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetGradientRoutine" 724*0c52818fSSatish Balay /*@C 725*0c52818fSSatish Balay TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search 726*0c52818fSSatish Balay 727*0c52818fSSatish Balay Logically Collective on TaoLineSearch 728*0c52818fSSatish Balay 729*0c52818fSSatish Balay Input Parameter: 730*0c52818fSSatish Balay + ls - the TaoLineSearch context 731*0c52818fSSatish Balay . func - the gradient evaluation routine 732*0c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 733*0c52818fSSatish Balay 734*0c52818fSSatish Balay Calling sequence of func: 735*0c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, Vec g, void *ctx); 736*0c52818fSSatish Balay 737*0c52818fSSatish Balay + x - input vector 738*0c52818fSSatish Balay . g - gradient vector 739*0c52818fSSatish Balay - ctx (optional) user-defined context 740*0c52818fSSatish Balay 741*0c52818fSSatish Balay Level: beginner 742*0c52818fSSatish Balay 743*0c52818fSSatish Balay Note: 744*0c52818fSSatish Balay Use this routine only if you want the line search gradient 745*0c52818fSSatish Balay evaluation routine to be different from the TaoSolver's gradient 746*0c52818fSSatish Balay evaluation routine. If you use this routine you must also set 747*0c52818fSSatish Balay the line search function and/or function/gradient routine. 748*0c52818fSSatish Balay 749*0c52818fSSatish Balay Note: 750*0c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own gradient routine for the 751*0c52818fSSatish Balay line search, application programmers should be wary of overriding the 752*0c52818fSSatish Balay default gradient routine. 753*0c52818fSSatish Balay 754*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoSolverRoutines() 755*0c52818fSSatish Balay @*/ 756*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx) 757*0c52818fSSatish Balay { 758*0c52818fSSatish Balay PetscFunctionBegin; 759*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 760*0c52818fSSatish Balay 761*0c52818fSSatish Balay ls->ops->computegradient=func; 762*0c52818fSSatish Balay if (ctx) ls->userctx_grad=ctx; 763*0c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 764*0c52818fSSatish Balay PetscFunctionReturn(0); 765*0c52818fSSatish Balay } 766*0c52818fSSatish Balay 767*0c52818fSSatish Balay #undef __FUNCT__ 768*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetObjectiveAndGradientRoutine" 769*0c52818fSSatish Balay /*@C 770*0c52818fSSatish Balay TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search 771*0c52818fSSatish Balay 772*0c52818fSSatish Balay Logically Collective on TaoLineSearch 773*0c52818fSSatish Balay 774*0c52818fSSatish Balay Input Parameter: 775*0c52818fSSatish Balay + ls - the TaoLineSearch context 776*0c52818fSSatish Balay . func - the objective and gradient evaluation routine 777*0c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 778*0c52818fSSatish Balay 779*0c52818fSSatish Balay Calling sequence of func: 780*0c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx); 781*0c52818fSSatish Balay 782*0c52818fSSatish Balay + x - input vector 783*0c52818fSSatish Balay . f - function value 784*0c52818fSSatish Balay . g - gradient vector 785*0c52818fSSatish Balay - ctx (optional) user-defined context 786*0c52818fSSatish Balay 787*0c52818fSSatish Balay Level: beginner 788*0c52818fSSatish Balay 789*0c52818fSSatish Balay Note: 790*0c52818fSSatish Balay Use this routine only if you want the line search objective and gradient 791*0c52818fSSatish Balay evaluation routines to be different from the TaoSolver's objective 792*0c52818fSSatish Balay and gradient evaluation routines. 793*0c52818fSSatish Balay 794*0c52818fSSatish Balay Note: 795*0c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 796*0c52818fSSatish Balay line search, application programmers should be wary of overriding the 797*0c52818fSSatish Balay default objective routine. 798*0c52818fSSatish Balay 799*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoSolverRoutines() 800*0c52818fSSatish Balay @*/ 801*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx) 802*0c52818fSSatish Balay { 803*0c52818fSSatish Balay PetscFunctionBegin; 804*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 805*0c52818fSSatish Balay 806*0c52818fSSatish Balay ls->ops->computeobjectiveandgradient=func; 807*0c52818fSSatish Balay if (ctx) ls->userctx_funcgrad=ctx; 808*0c52818fSSatish Balay ls->usetaoroutines = PETSC_FALSE; 809*0c52818fSSatish Balay PetscFunctionReturn(0); 810*0c52818fSSatish Balay } 811*0c52818fSSatish Balay 812*0c52818fSSatish Balay #undef __FUNCT__ 813*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetObjectiveAndGTSRoutine" 814*0c52818fSSatish Balay /*@C 815*0c52818fSSatish Balay TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and 816*0c52818fSSatish Balay (gradient'*stepdirection) evaluation routine for the line search. 817*0c52818fSSatish Balay Sometimes it is more efficient to compute the inner product of the gradient 818*0c52818fSSatish Balay and the step direction than it is to compute the gradient, and this is all 819*0c52818fSSatish Balay the line search typically needs of the gradient. 820*0c52818fSSatish Balay 821*0c52818fSSatish Balay Logically Collective on TaoLineSearch 822*0c52818fSSatish Balay 823*0c52818fSSatish Balay Input Parameter: 824*0c52818fSSatish Balay + ls - the TaoLineSearch context 825*0c52818fSSatish Balay . func - the objective and gradient evaluation routine 826*0c52818fSSatish Balay - ctx - the (optional) user-defined context for private data 827*0c52818fSSatish Balay 828*0c52818fSSatish Balay Calling sequence of func: 829*0c52818fSSatish Balay $ func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx); 830*0c52818fSSatish Balay 831*0c52818fSSatish Balay + x - input vector 832*0c52818fSSatish Balay . s - step direction 833*0c52818fSSatish Balay . f - function value 834*0c52818fSSatish Balay . gts - inner product of gradient and step direction vectors 835*0c52818fSSatish Balay - ctx (optional) user-defined context 836*0c52818fSSatish Balay 837*0c52818fSSatish Balay Note: The gradient will still need to be computed at the end of the line 838*0c52818fSSatish Balay search, so you will still need to set a line search gradient evaluation 839*0c52818fSSatish Balay routine 840*0c52818fSSatish Balay 841*0c52818fSSatish Balay Note: Bounded line searches (those used in bounded optimization algorithms) 842*0c52818fSSatish Balay don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the 843*0c52818fSSatish Balay x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength() 844*0c52818fSSatish Balay 845*0c52818fSSatish Balay Level: advanced 846*0c52818fSSatish Balay 847*0c52818fSSatish Balay Note: 848*0c52818fSSatish Balay Some algorithms (lcl, gpcg) set their own objective routine for the 849*0c52818fSSatish Balay line search, application programmers should be wary of overriding the 850*0c52818fSSatish Balay default objective routine. 851*0c52818fSSatish Balay 852*0c52818fSSatish Balay .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoSolverRoutines() 853*0c52818fSSatish Balay @*/ 854*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx) 855*0c52818fSSatish Balay { 856*0c52818fSSatish Balay PetscFunctionBegin; 857*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 858*0c52818fSSatish Balay 859*0c52818fSSatish Balay ls->ops->computeobjectiveandgts=func; 860*0c52818fSSatish Balay if (ctx) ls->userctx_funcgts=ctx; 861*0c52818fSSatish Balay ls->usegts = PETSC_TRUE; 862*0c52818fSSatish Balay ls->usetaoroutines=PETSC_FALSE; 863*0c52818fSSatish Balay PetscFunctionReturn(0); 864*0c52818fSSatish Balay } 865*0c52818fSSatish Balay 866*0c52818fSSatish Balay #undef __FUNCT__ 867*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchUseTaoSolverRoutines" 868*0c52818fSSatish Balay /*@C 869*0c52818fSSatish Balay TaoLineSearchUseTaoSolverRoutines - Informs the TaoLineSearch to use the 870*0c52818fSSatish Balay objective and gradient evaluation routines from the given TaoSolver object. 871*0c52818fSSatish Balay 872*0c52818fSSatish Balay Logically Collective on TaoLineSearch 873*0c52818fSSatish Balay 874*0c52818fSSatish Balay Input Parameter: 875*0c52818fSSatish Balay + ls - the TaoLineSearch context 876*0c52818fSSatish Balay - ts - the TaoSolver context with defined objective/gradient evaluation routines 877*0c52818fSSatish Balay 878*0c52818fSSatish Balay Level: developer 879*0c52818fSSatish Balay 880*0c52818fSSatish Balay .seealso: TaoLineSearchCreate() 881*0c52818fSSatish Balay @*/ 882*0c52818fSSatish Balay PetscErrorCode TaoLineSearchUseTaoSolverRoutines(TaoLineSearch ls, TaoSolver ts) 883*0c52818fSSatish Balay { 884*0c52818fSSatish Balay PetscFunctionBegin; 885*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 886*0c52818fSSatish Balay PetscValidHeaderSpecific(ts,TAOSOLVER_CLASSID,1); 887*0c52818fSSatish Balay ls->taosolver = ts; 888*0c52818fSSatish Balay ls->usetaoroutines=PETSC_TRUE; 889*0c52818fSSatish Balay PetscFunctionReturn(0); 890*0c52818fSSatish Balay } 891*0c52818fSSatish Balay 892*0c52818fSSatish Balay 893*0c52818fSSatish Balay #undef __FUNCT__ 894*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeObjective" 895*0c52818fSSatish Balay /*@ 896*0c52818fSSatish Balay TaoLineSearchComputeObjective - Computes the objective function value at a given point 897*0c52818fSSatish Balay 898*0c52818fSSatish Balay Collective on TaoLineSearch 899*0c52818fSSatish Balay 900*0c52818fSSatish Balay Input Parameters: 901*0c52818fSSatish Balay + ls - the TaoLineSearch context 902*0c52818fSSatish Balay - x - input vector 903*0c52818fSSatish Balay 904*0c52818fSSatish Balay Output Parameter: 905*0c52818fSSatish Balay . f - Objective value at X 906*0c52818fSSatish Balay 907*0c52818fSSatish Balay Notes: TaoLineSearchComputeObjective() is typically used within line searches 908*0c52818fSSatish Balay so most users would not generally call this routine themselves. 909*0c52818fSSatish Balay 910*0c52818fSSatish Balay Level: developer 911*0c52818fSSatish Balay 912*0c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 913*0c52818fSSatish Balay @*/ 914*0c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 915*0c52818fSSatish Balay { 916*0c52818fSSatish Balay PetscErrorCode ierr; 917*0c52818fSSatish Balay Vec gdummy; 918*0c52818fSSatish Balay PetscReal gts; 919*0c52818fSSatish Balay PetscFunctionBegin; 920*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 921*0c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 922*0c52818fSSatish Balay PetscValidPointer(f,3); 923*0c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 924*0c52818fSSatish Balay if (ls->usetaoroutines) { 925*0c52818fSSatish Balay ierr = TaoComputeObjective(ls->taosolver,x,f); CHKERRQ(ierr); 926*0c52818fSSatish Balay } else { 927*0c52818fSSatish Balay ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 928*0c52818fSSatish Balay if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient 929*0c52818fSSatish Balay && !ls->ops->computeobjectiveandgts) { 930*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 931*0c52818fSSatish Balay } 932*0c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective routine"); 933*0c52818fSSatish Balay CHKMEMQ; 934*0c52818fSSatish Balay if (ls->ops->computeobjective) { 935*0c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func); CHKERRQ(ierr); 936*0c52818fSSatish Balay } else if (ls->ops->computeobjectiveandgradient) { 937*0c52818fSSatish Balay ierr = VecDuplicate(x,&gdummy); CHKERRQ(ierr); 938*0c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad); CHKERRQ(ierr); 939*0c52818fSSatish Balay ierr = VecDestroy(&gdummy); CHKERRQ(ierr); 940*0c52818fSSatish Balay } else { 941*0c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts); CHKERRQ(ierr); 942*0c52818fSSatish Balay } 943*0c52818fSSatish Balay CHKMEMQ; 944*0c52818fSSatish Balay PetscStackPop; 945*0c52818fSSatish Balay ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 946*0c52818fSSatish Balay } 947*0c52818fSSatish Balay ls->nfeval++; 948*0c52818fSSatish Balay PetscFunctionReturn(0); 949*0c52818fSSatish Balay 950*0c52818fSSatish Balay } 951*0c52818fSSatish Balay 952*0c52818fSSatish Balay 953*0c52818fSSatish Balay #undef __FUNCT__ 954*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeObjectiveAndGradient" 955*0c52818fSSatish Balay /*@ 956*0c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 957*0c52818fSSatish Balay 958*0c52818fSSatish Balay Collective on TaoSOlver 959*0c52818fSSatish Balay 960*0c52818fSSatish Balay Input Parameters: 961*0c52818fSSatish Balay + ls - the TaoLineSearch context 962*0c52818fSSatish Balay - x - input vector 963*0c52818fSSatish Balay 964*0c52818fSSatish Balay Output Parameter: 965*0c52818fSSatish Balay + f - Objective value at X 966*0c52818fSSatish Balay - g - Gradient vector at X 967*0c52818fSSatish Balay 968*0c52818fSSatish Balay Notes: TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 969*0c52818fSSatish Balay so most users would not generally call this routine themselves. 970*0c52818fSSatish Balay 971*0c52818fSSatish Balay Level: developer 972*0c52818fSSatish Balay 973*0c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 974*0c52818fSSatish Balay @*/ 975*0c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 976*0c52818fSSatish Balay { 977*0c52818fSSatish Balay PetscErrorCode ierr; 978*0c52818fSSatish Balay PetscFunctionBegin; 979*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 980*0c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 981*0c52818fSSatish Balay PetscValidPointer(f,3); 982*0c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 983*0c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 984*0c52818fSSatish Balay PetscCheckSameComm(ls,1,g,4); 985*0c52818fSSatish Balay if (ls->usetaoroutines) { 986*0c52818fSSatish Balay ierr = TaoComputeObjectiveAndGradient(ls->taosolver,x,f,g); 987*0c52818fSSatish Balay CHKERRQ(ierr); 988*0c52818fSSatish Balay } else { 989*0c52818fSSatish Balay ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 990*0c52818fSSatish Balay if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) { 991*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 992*0c52818fSSatish Balay } 993*0c52818fSSatish Balay if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) { 994*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 995*0c52818fSSatish Balay } 996*0c52818fSSatish Balay 997*0c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gradient routine"); 998*0c52818fSSatish Balay CHKMEMQ; 999*0c52818fSSatish Balay if (ls->ops->computeobjectiveandgradient) { 1000*0c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad); CHKERRQ(ierr); 1001*0c52818fSSatish Balay } else { 1002*0c52818fSSatish Balay ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func); CHKERRQ(ierr); 1003*0c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad); CHKERRQ(ierr); 1004*0c52818fSSatish Balay } 1005*0c52818fSSatish Balay CHKMEMQ; 1006*0c52818fSSatish Balay PetscStackPop; 1007*0c52818fSSatish Balay ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 1008*0c52818fSSatish Balay } 1009*0c52818fSSatish Balay ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",*f);CHKERRQ(ierr); 1010*0c52818fSSatish Balay ls->nfgeval++; 1011*0c52818fSSatish Balay PetscFunctionReturn(0); 1012*0c52818fSSatish Balay } 1013*0c52818fSSatish Balay 1014*0c52818fSSatish Balay #undef __FUNCT__ 1015*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeGradient" 1016*0c52818fSSatish Balay /*@ 1017*0c52818fSSatish Balay TaoLineSearchComputeGradient - Computes the gradient of the objective function 1018*0c52818fSSatish Balay 1019*0c52818fSSatish Balay Collective on TaoLineSearch 1020*0c52818fSSatish Balay 1021*0c52818fSSatish Balay Input Parameters: 1022*0c52818fSSatish Balay + ls - the TaoLineSearch context 1023*0c52818fSSatish Balay - x - input vector 1024*0c52818fSSatish Balay 1025*0c52818fSSatish Balay Output Parameter: 1026*0c52818fSSatish Balay . g - gradient vector 1027*0c52818fSSatish Balay 1028*0c52818fSSatish Balay Notes: TaoComputeGradient() is typically used within line searches 1029*0c52818fSSatish Balay so most users would not generally call this routine themselves. 1030*0c52818fSSatish Balay 1031*0c52818fSSatish Balay Level: developer 1032*0c52818fSSatish Balay 1033*0c52818fSSatish Balay .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 1034*0c52818fSSatish Balay @*/ 1035*0c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 1036*0c52818fSSatish Balay { 1037*0c52818fSSatish Balay PetscErrorCode ierr; 1038*0c52818fSSatish Balay PetscReal fdummy; 1039*0c52818fSSatish Balay PetscFunctionBegin; 1040*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1041*0c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1042*0c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,3); 1043*0c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 1044*0c52818fSSatish Balay PetscCheckSameComm(ls,1,g,3); 1045*0c52818fSSatish Balay if (ls->usetaoroutines) { 1046*0c52818fSSatish Balay ierr = TaoComputeGradient(ls->taosolver,x,g); CHKERRQ(ierr); 1047*0c52818fSSatish Balay } else { 1048*0c52818fSSatish Balay ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 1049*0c52818fSSatish Balay if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) { 1050*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 1051*0c52818fSSatish Balay } 1052*0c52818fSSatish Balay PetscStackPush("TaoLineSearch user gradient routine"); 1053*0c52818fSSatish Balay CHKMEMQ; 1054*0c52818fSSatish Balay if (ls->ops->computegradient) { 1055*0c52818fSSatish Balay ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad); CHKERRQ(ierr); 1056*0c52818fSSatish Balay } else { 1057*0c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad); CHKERRQ(ierr); 1058*0c52818fSSatish Balay } 1059*0c52818fSSatish Balay CHKMEMQ; 1060*0c52818fSSatish Balay PetscStackPop; 1061*0c52818fSSatish Balay ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 1062*0c52818fSSatish Balay } 1063*0c52818fSSatish Balay ls->ngeval++; 1064*0c52818fSSatish Balay PetscFunctionReturn(0); 1065*0c52818fSSatish Balay } 1066*0c52818fSSatish Balay 1067*0c52818fSSatish Balay #undef __FUNCT__ 1068*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchComputeObjectiveAndGTS" 1069*0c52818fSSatish Balay /*@ 1070*0c52818fSSatish Balay TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 1071*0c52818fSSatish Balay 1072*0c52818fSSatish Balay Collective on TaoSolver 1073*0c52818fSSatish Balay 1074*0c52818fSSatish Balay Input Parameters: 1075*0c52818fSSatish Balay + ls - the TaoLineSearch context 1076*0c52818fSSatish Balay - x - input vector 1077*0c52818fSSatish Balay 1078*0c52818fSSatish Balay Output Parameter: 1079*0c52818fSSatish Balay + f - Objective value at X 1080*0c52818fSSatish Balay - gts - inner product of gradient and step direction at X 1081*0c52818fSSatish Balay 1082*0c52818fSSatish Balay Notes: TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 1083*0c52818fSSatish Balay so most users would not generally call this routine themselves. 1084*0c52818fSSatish Balay 1085*0c52818fSSatish Balay Level: developer 1086*0c52818fSSatish Balay 1087*0c52818fSSatish Balay .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 1088*0c52818fSSatish Balay @*/ 1089*0c52818fSSatish Balay PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 1090*0c52818fSSatish Balay { 1091*0c52818fSSatish Balay PetscErrorCode ierr; 1092*0c52818fSSatish Balay PetscFunctionBegin; 1093*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1094*0c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1095*0c52818fSSatish Balay PetscValidPointer(f,3); 1096*0c52818fSSatish Balay PetscValidPointer(gts,4); 1097*0c52818fSSatish Balay PetscCheckSameComm(ls,1,x,2); 1098*0c52818fSSatish Balay ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 1099*0c52818fSSatish Balay if (!ls->ops->computeobjectiveandgts) { 1100*0c52818fSSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 1101*0c52818fSSatish Balay } 1102*0c52818fSSatish Balay 1103*0c52818fSSatish Balay PetscStackPush("TaoLineSearch user objective/gts routine"); 1104*0c52818fSSatish Balay CHKMEMQ; 1105*0c52818fSSatish Balay ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts); CHKERRQ(ierr); 1106*0c52818fSSatish Balay CHKMEMQ; 1107*0c52818fSSatish Balay PetscStackPop; 1108*0c52818fSSatish Balay ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0); CHKERRQ(ierr); 1109*0c52818fSSatish Balay ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",*f);CHKERRQ(ierr); 1110*0c52818fSSatish Balay ls->nfeval++; 1111*0c52818fSSatish Balay PetscFunctionReturn(0); 1112*0c52818fSSatish Balay } 1113*0c52818fSSatish Balay 1114*0c52818fSSatish Balay 1115*0c52818fSSatish Balay #undef __FUNCT__ 1116*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetSolution" 1117*0c52818fSSatish Balay /*@ 1118*0c52818fSSatish Balay TaoLineSearchGetSolution - Returns the solution to the line search 1119*0c52818fSSatish Balay 1120*0c52818fSSatish Balay Collective on TaoLineSearch 1121*0c52818fSSatish Balay 1122*0c52818fSSatish Balay Input Parameter: 1123*0c52818fSSatish Balay . ls - the TaoLineSearch context 1124*0c52818fSSatish Balay 1125*0c52818fSSatish Balay Output Parameter: 1126*0c52818fSSatish Balay + x - the new solution 1127*0c52818fSSatish Balay . f - the objective function value at x 1128*0c52818fSSatish Balay . g - the gradient at x 1129*0c52818fSSatish Balay . steplength - the multiple of the step direction taken by the line search 1130*0c52818fSSatish Balay - reason - the reason why the line search terminated 1131*0c52818fSSatish Balay 1132*0c52818fSSatish Balay reason will be set to one of: 1133*0c52818fSSatish Balay 1134*0c52818fSSatish Balay + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 1135*0c52818fSSatish Balay . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 1136*0c52818fSSatish Balay . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 1137*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 1138*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 1139*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 1140*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 1141*0c52818fSSatish Balay 1142*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 1143*0c52818fSSatish Balay . TAOLINESEARCH_HALTED_OTHER - any other reason 1144*0c52818fSSatish Balay 1145*0c52818fSSatish Balay + TAOLINESEARCH_SUCCESS - successful line search 1146*0c52818fSSatish Balay 1147*0c52818fSSatish Balay Level: developer 1148*0c52818fSSatish Balay 1149*0c52818fSSatish Balay @*/ 1150*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchTerminationReason *reason) 1151*0c52818fSSatish Balay { 1152*0c52818fSSatish Balay PetscErrorCode ierr; 1153*0c52818fSSatish Balay PetscFunctionBegin; 1154*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1155*0c52818fSSatish Balay PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1156*0c52818fSSatish Balay PetscValidPointer(f,3); 1157*0c52818fSSatish Balay PetscValidHeaderSpecific(g,VEC_CLASSID,4); 1158*0c52818fSSatish Balay PetscValidIntPointer(reason,6); 1159*0c52818fSSatish Balay 1160*0c52818fSSatish Balay if (ls->new_x) { 1161*0c52818fSSatish Balay ierr = VecCopy(ls->new_x,x); CHKERRQ(ierr); 1162*0c52818fSSatish Balay } 1163*0c52818fSSatish Balay *f = ls->new_f; 1164*0c52818fSSatish Balay if (ls->new_g) { 1165*0c52818fSSatish Balay ierr = VecCopy(ls->new_g,g); CHKERRQ(ierr); 1166*0c52818fSSatish Balay } 1167*0c52818fSSatish Balay if (steplength) { 1168*0c52818fSSatish Balay *steplength=ls->step; 1169*0c52818fSSatish Balay } 1170*0c52818fSSatish Balay *reason = ls->reason; 1171*0c52818fSSatish Balay PetscFunctionReturn(0); 1172*0c52818fSSatish Balay } 1173*0c52818fSSatish Balay 1174*0c52818fSSatish Balay #undef __FUNCT__ 1175*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetStartingVector" 1176*0c52818fSSatish Balay /*@ 1177*0c52818fSSatish Balay TaoLineSearchGetStartingVector - Gets a the initial point of the line 1178*0c52818fSSatish Balay search. 1179*0c52818fSSatish Balay 1180*0c52818fSSatish Balay Not Collective 1181*0c52818fSSatish Balay 1182*0c52818fSSatish Balay Input Parameter: 1183*0c52818fSSatish Balay . ls - the TaoLineSearch context 1184*0c52818fSSatish Balay 1185*0c52818fSSatish Balay Output Parameter: 1186*0c52818fSSatish Balay . x - The initial point of the line search 1187*0c52818fSSatish Balay 1188*0c52818fSSatish Balay Level: intermediate 1189*0c52818fSSatish Balay @*/ 1190*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 1191*0c52818fSSatish Balay { 1192*0c52818fSSatish Balay PetscFunctionBegin; 1193*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1194*0c52818fSSatish Balay if (x) { 1195*0c52818fSSatish Balay *x = ls->start_x; 1196*0c52818fSSatish Balay } 1197*0c52818fSSatish Balay PetscFunctionReturn(0); 1198*0c52818fSSatish Balay 1199*0c52818fSSatish Balay } 1200*0c52818fSSatish Balay 1201*0c52818fSSatish Balay #undef __FUNCT__ 1202*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetStepDirection" 1203*0c52818fSSatish Balay /*@ 1204*0c52818fSSatish Balay TaoLineSearchGetStepDirection - Gets the step direction of the line 1205*0c52818fSSatish Balay search. 1206*0c52818fSSatish Balay 1207*0c52818fSSatish Balay Not Collective 1208*0c52818fSSatish Balay 1209*0c52818fSSatish Balay Input Parameter: 1210*0c52818fSSatish Balay . ls - the TaoLineSearch context 1211*0c52818fSSatish Balay 1212*0c52818fSSatish Balay Output Parameter: 1213*0c52818fSSatish Balay . s - the step direction of the line search 1214*0c52818fSSatish Balay 1215*0c52818fSSatish Balay Level: advanced 1216*0c52818fSSatish Balay @*/ 1217*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 1218*0c52818fSSatish Balay { 1219*0c52818fSSatish Balay PetscFunctionBegin; 1220*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1221*0c52818fSSatish Balay if (s) { 1222*0c52818fSSatish Balay *s = ls->stepdirection; 1223*0c52818fSSatish Balay } 1224*0c52818fSSatish Balay PetscFunctionReturn(0); 1225*0c52818fSSatish Balay 1226*0c52818fSSatish Balay } 1227*0c52818fSSatish Balay 1228*0c52818fSSatish Balay #undef __FUNCT__ 1229*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetFullStepObjective" 1230*0c52818fSSatish Balay /*@ 1231*0c52818fSSatish Balay TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 1232*0c52818fSSatish Balay 1233*0c52818fSSatish Balay Not Collective 1234*0c52818fSSatish Balay 1235*0c52818fSSatish Balay Input Parameter: 1236*0c52818fSSatish Balay . ls - the TaoLineSearch context 1237*0c52818fSSatish Balay 1238*0c52818fSSatish Balay Output Parameter: 1239*0c52818fSSatish Balay . f - the objective value at the full step length 1240*0c52818fSSatish Balay 1241*0c52818fSSatish Balay Level: developer 1242*0c52818fSSatish Balay @*/ 1243*0c52818fSSatish Balay 1244*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 1245*0c52818fSSatish Balay { 1246*0c52818fSSatish Balay PetscFunctionBegin; 1247*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1248*0c52818fSSatish Balay *f_fullstep = ls->f_fullstep; 1249*0c52818fSSatish Balay PetscFunctionReturn(0); 1250*0c52818fSSatish Balay } 1251*0c52818fSSatish Balay 1252*0c52818fSSatish Balay #undef __FUNCT__ 1253*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetVariableBounds" 1254*0c52818fSSatish Balay /*@ 1255*0c52818fSSatish Balay TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 1256*0c52818fSSatish Balay 1257*0c52818fSSatish Balay Logically Collective on TaoSolver 1258*0c52818fSSatish Balay 1259*0c52818fSSatish Balay Input Parameters: 1260*0c52818fSSatish Balay + ls - the TaoLineSearch context 1261*0c52818fSSatish Balay . xl - vector of lower bounds 1262*0c52818fSSatish Balay - xu - vector of upper bounds 1263*0c52818fSSatish Balay 1264*0c52818fSSatish Balay Note: If the variable bounds are not set with this routine, then 1265*0c52818fSSatish Balay TAO_NINFINITY and TAO_INFINITY are assumed 1266*0c52818fSSatish Balay 1267*0c52818fSSatish Balay Level: beginner 1268*0c52818fSSatish Balay 1269*0c52818fSSatish Balay .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 1270*0c52818fSSatish Balay @*/ 1271*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 1272*0c52818fSSatish Balay { 1273*0c52818fSSatish Balay PetscFunctionBegin; 1274*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1275*0c52818fSSatish Balay PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1276*0c52818fSSatish Balay PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1277*0c52818fSSatish Balay 1278*0c52818fSSatish Balay ls->lower = xl; 1279*0c52818fSSatish Balay ls->upper = xu; 1280*0c52818fSSatish Balay ls->bounded = 1; 1281*0c52818fSSatish Balay 1282*0c52818fSSatish Balay PetscFunctionReturn(0); 1283*0c52818fSSatish Balay 1284*0c52818fSSatish Balay } 1285*0c52818fSSatish Balay 1286*0c52818fSSatish Balay 1287*0c52818fSSatish Balay #undef __FUNCT__ 1288*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetInitialStepLength" 1289*0c52818fSSatish Balay /*@ 1290*0c52818fSSatish Balay TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 1291*0c52818fSSatish Balay search. If this value is not set then 1.0 is assumed. 1292*0c52818fSSatish Balay 1293*0c52818fSSatish Balay Logically Collective on TaoLineSearch 1294*0c52818fSSatish Balay 1295*0c52818fSSatish Balay Input Parameters: 1296*0c52818fSSatish Balay + ls - the TaoLineSearch context 1297*0c52818fSSatish Balay - s - the initial step size 1298*0c52818fSSatish Balay 1299*0c52818fSSatish Balay Level: intermediate 1300*0c52818fSSatish Balay 1301*0c52818fSSatish Balay .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply() 1302*0c52818fSSatish Balay @*/ 1303*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 1304*0c52818fSSatish Balay { 1305*0c52818fSSatish Balay PetscFunctionBegin; 1306*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1307*0c52818fSSatish Balay 1308*0c52818fSSatish Balay ls->initstep = s; 1309*0c52818fSSatish Balay PetscFunctionReturn(0); 1310*0c52818fSSatish Balay } 1311*0c52818fSSatish Balay 1312*0c52818fSSatish Balay #undef __FUNCT__ 1313*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetStepLength" 1314*0c52818fSSatish Balay /*@ 1315*0c52818fSSatish Balay TaoLineSearchGetStepLength - Get the current step length 1316*0c52818fSSatish Balay 1317*0c52818fSSatish Balay Not Collective 1318*0c52818fSSatish Balay 1319*0c52818fSSatish Balay Input Parameters: 1320*0c52818fSSatish Balay . ls - the TaoLineSearch context 1321*0c52818fSSatish Balay 1322*0c52818fSSatish Balay Output Parameters 1323*0c52818fSSatish Balay . s - the current step length 1324*0c52818fSSatish Balay 1325*0c52818fSSatish Balay Level: beginner 1326*0c52818fSSatish Balay 1327*0c52818fSSatish Balay .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply() 1328*0c52818fSSatish Balay @*/ 1329*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 1330*0c52818fSSatish Balay { 1331*0c52818fSSatish Balay PetscFunctionBegin; 1332*0c52818fSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1333*0c52818fSSatish Balay *s = ls->step; 1334*0c52818fSSatish Balay PetscFunctionReturn(0); 1335*0c52818fSSatish Balay } 1336*0c52818fSSatish Balay 1337*0c52818fSSatish Balay #undef __FUNCT__ 1338*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchRegister" 1339*0c52818fSSatish Balay /*MC 1340*0c52818fSSatish Balay TaoLineSearchRegister - Adds a line-search algorithm to the registry 1341*0c52818fSSatish Balay 1342*0c52818fSSatish Balay Not collective 1343*0c52818fSSatish Balay 1344*0c52818fSSatish Balay Input Parameters: 1345*0c52818fSSatish Balay + sname - name of a new user-defined solver 1346*0c52818fSSatish Balay - func - routine to Create method context 1347*0c52818fSSatish Balay 1348*0c52818fSSatish Balay Notes: 1349*0c52818fSSatish Balay TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 1350*0c52818fSSatish Balay 1351*0c52818fSSatish Balay Sample usage: 1352*0c52818fSSatish Balay .vb 1353*0c52818fSSatish Balay TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 1354*0c52818fSSatish Balay .ve 1355*0c52818fSSatish Balay 1356*0c52818fSSatish Balay Then, your solver can be chosen with the procedural interface via 1357*0c52818fSSatish Balay $ TaoLineSearchSetType(ls,"my_linesearch") 1358*0c52818fSSatish Balay or at runtime via the option 1359*0c52818fSSatish Balay $ -tao_ls_type my_linesearch 1360*0c52818fSSatish Balay 1361*0c52818fSSatish Balay Level: developer 1362*0c52818fSSatish Balay 1363*0c52818fSSatish Balay .seealso: TaoLineSearchRegisterDestroy() 1364*0c52818fSSatish Balay M*/ 1365*0c52818fSSatish Balay PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 1366*0c52818fSSatish Balay { 1367*0c52818fSSatish Balay PetscErrorCode ierr; 1368*0c52818fSSatish Balay PetscFunctionBegin; 1369*0c52818fSSatish Balay ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func); CHKERRQ(ierr); 1370*0c52818fSSatish Balay PetscFunctionReturn(0); 1371*0c52818fSSatish Balay } 1372*0c52818fSSatish Balay 1373*0c52818fSSatish Balay #undef __FUNCT__ 1374*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchRegisterDestroy" 1375*0c52818fSSatish Balay /*@C 1376*0c52818fSSatish Balay TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were 1377*0c52818fSSatish Balay registered by TaoLineSearchRegister(). 1378*0c52818fSSatish Balay 1379*0c52818fSSatish Balay Not Collective 1380*0c52818fSSatish Balay 1381*0c52818fSSatish Balay Level: developer 1382*0c52818fSSatish Balay 1383*0c52818fSSatish Balay .seealso: TaoLineSearchRegister() 1384*0c52818fSSatish Balay @*/ 1385*0c52818fSSatish Balay PetscErrorCode TaoLineSearchRegisterDestroy(void) 1386*0c52818fSSatish Balay { 1387*0c52818fSSatish Balay PetscErrorCode ierr; 1388*0c52818fSSatish Balay PetscFunctionBegin; 1389*0c52818fSSatish Balay ierr = PetscFunctionListDestroy(&TaoLineSearchList); CHKERRQ(ierr); 1390*0c52818fSSatish Balay TaoLineSearchInitialized = PETSC_FALSE; 1391*0c52818fSSatish Balay PetscFunctionReturn(0); 1392*0c52818fSSatish Balay } 1393*0c52818fSSatish Balay 1394*0c52818fSSatish Balay 1395*0c52818fSSatish Balay #undef __FUNCT__ 1396*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchAppendOptionsPrefix" 1397*0c52818fSSatish Balay /*@C 1398*0c52818fSSatish Balay TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 1399*0c52818fSSatish Balay for all TaoLineSearch options in the database. 1400*0c52818fSSatish Balay 1401*0c52818fSSatish Balay 1402*0c52818fSSatish Balay Collective on TaoLineSearch 1403*0c52818fSSatish Balay 1404*0c52818fSSatish Balay Input Parameters: 1405*0c52818fSSatish Balay + ls - the TaoLineSearch solver context 1406*0c52818fSSatish Balay - prefix - the prefix string to prepend to all line search requests 1407*0c52818fSSatish Balay 1408*0c52818fSSatish Balay Notes: 1409*0c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1410*0c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 1411*0c52818fSSatish Balay 1412*0c52818fSSatish Balay 1413*0c52818fSSatish Balay Level: advanced 1414*0c52818fSSatish Balay 1415*0c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 1416*0c52818fSSatish Balay @*/ 1417*0c52818fSSatish Balay PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 1418*0c52818fSSatish Balay { 1419*0c52818fSSatish Balay PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 1420*0c52818fSSatish Balay return(0); 1421*0c52818fSSatish Balay } 1422*0c52818fSSatish Balay 1423*0c52818fSSatish Balay #undef __FUNCT__ 1424*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchGetOptionsPrefix" 1425*0c52818fSSatish Balay /*@C 1426*0c52818fSSatish Balay TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1427*0c52818fSSatish Balay TaoLineSearch options in the database 1428*0c52818fSSatish Balay 1429*0c52818fSSatish Balay Not Collective 1430*0c52818fSSatish Balay 1431*0c52818fSSatish Balay Input Parameters: 1432*0c52818fSSatish Balay . ls - the TaoLineSearch context 1433*0c52818fSSatish Balay 1434*0c52818fSSatish Balay Output Parameters: 1435*0c52818fSSatish Balay . prefix - pointer to the prefix string used is returned 1436*0c52818fSSatish Balay 1437*0c52818fSSatish Balay Notes: On the fortran side, the user should pass in a string 'prefix' of 1438*0c52818fSSatish Balay sufficient length to hold the prefix. 1439*0c52818fSSatish Balay 1440*0c52818fSSatish Balay Level: advanced 1441*0c52818fSSatish Balay 1442*0c52818fSSatish Balay .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 1443*0c52818fSSatish Balay @*/ 1444*0c52818fSSatish Balay PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 1445*0c52818fSSatish Balay { 1446*0c52818fSSatish Balay PetscObjectGetOptionsPrefix((PetscObject)ls,p); 1447*0c52818fSSatish Balay return 0; 1448*0c52818fSSatish Balay } 1449*0c52818fSSatish Balay 1450*0c52818fSSatish Balay #undef __FUNCT__ 1451*0c52818fSSatish Balay #define __FUNCT__ "TaoLineSearchSetOptionsPrefix" 1452*0c52818fSSatish Balay /*@C 1453*0c52818fSSatish Balay TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 1454*0c52818fSSatish Balay TaoLineSearch options in the database. 1455*0c52818fSSatish Balay 1456*0c52818fSSatish Balay 1457*0c52818fSSatish Balay Logically Collective on TaoLineSearch 1458*0c52818fSSatish Balay 1459*0c52818fSSatish Balay Input Parameters: 1460*0c52818fSSatish Balay + ls - the TaoLineSearch context 1461*0c52818fSSatish Balay - prefix - the prefix string to prepend to all TAO option requests 1462*0c52818fSSatish Balay 1463*0c52818fSSatish Balay Notes: 1464*0c52818fSSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1465*0c52818fSSatish Balay The first character of all runtime options is AUTOMATICALLY the hyphen. 1466*0c52818fSSatish Balay 1467*0c52818fSSatish Balay For example, to distinguish between the runtime options for two 1468*0c52818fSSatish Balay different line searches, one could call 1469*0c52818fSSatish Balay .vb 1470*0c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 1471*0c52818fSSatish Balay TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 1472*0c52818fSSatish Balay .ve 1473*0c52818fSSatish Balay 1474*0c52818fSSatish Balay This would enable use of different options for each system, such as 1475*0c52818fSSatish Balay .vb 1476*0c52818fSSatish Balay -sys1_tao_ls_type mt 1477*0c52818fSSatish Balay -sys2_tao_ls_type armijo 1478*0c52818fSSatish Balay .ve 1479*0c52818fSSatish Balay 1480*0c52818fSSatish Balay 1481*0c52818fSSatish Balay Level: advanced 1482*0c52818fSSatish Balay 1483*0c52818fSSatish Balay .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 1484*0c52818fSSatish Balay @*/ 1485*0c52818fSSatish Balay 1486*0c52818fSSatish Balay PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 1487*0c52818fSSatish Balay { 1488*0c52818fSSatish Balay PetscObjectSetOptionsPrefix((PetscObject)ls,p); 1489*0c52818fSSatish Balay return(0); 1490*0c52818fSSatish Balay } 1491