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