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