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