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