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