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