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 PetscInt tabs; 54 ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr); 55 ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)ls)->tablevel);CHKERRQ(ierr); 56 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);CHKERRQ(ierr); 57 if (ls->ops->view) { 58 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 59 ierr = (*ls->ops->view)(ls,viewer);CHKERRQ(ierr); 60 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 61 } 62 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 63 ierr = PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);CHKERRQ(ierr); 64 ierr = PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);CHKERRQ(ierr); 65 ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);CHKERRQ(ierr); 66 ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);CHKERRQ(ierr); 67 ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);CHKERRQ(ierr); 68 69 if (ls->bounded) { 70 ierr = PetscViewerASCIIPrintf(viewer,"using variable bounds\n");CHKERRQ(ierr); 71 } 72 ierr = PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);CHKERRQ(ierr); 73 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 74 ierr = PetscViewerASCIISetTab(viewer, tabs);CHKERRQ(ierr); 75 } else if (isstring) { 76 ierr = TaoLineSearchGetType(ls,&type);CHKERRQ(ierr); 77 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 78 } 79 PetscFunctionReturn(0); 80 } 81 82 /*@C 83 TaoLineSearchCreate - Creates a TAO Line Search object. Algorithms in TAO that use 84 line-searches will automatically create one. 85 86 Collective on MPI_Comm 87 88 Input Parameter: 89 . comm - MPI communicator 90 91 Output Parameter: 92 . newls - the new TaoLineSearch context 93 94 Available methods include: 95 + more-thuente 96 . gpcg 97 - unit - Do not perform any line search 98 99 100 Options Database Keys: 101 . -tao_ls_type - select which method TAO should use 102 103 Level: beginner 104 105 .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy() 106 @*/ 107 108 PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls) 109 { 110 PetscErrorCode ierr; 111 TaoLineSearch ls; 112 113 PetscFunctionBegin; 114 PetscValidPointer(newls,2); 115 *newls = NULL; 116 117 ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr); 118 119 ierr = PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);CHKERRQ(ierr); 120 ls->bounded = 0; 121 ls->max_funcs=30; 122 ls->ftol = 0.0001; 123 ls->gtol = 0.9; 124 #if defined(PETSC_USE_REAL_SINGLE) 125 ls->rtol = 1.0e-5; 126 #else 127 ls->rtol = 1.0e-10; 128 #endif 129 ls->stepmin=1.0e-20; 130 ls->stepmax=1.0e+20; 131 ls->step=1.0; 132 ls->nfeval=0; 133 ls->ngeval=0; 134 ls->nfgeval=0; 135 136 ls->ops->computeobjective=0; 137 ls->ops->computegradient=0; 138 ls->ops->computeobjectiveandgradient=0; 139 ls->ops->computeobjectiveandgts=0; 140 ls->ops->setup=0; 141 ls->ops->apply=0; 142 ls->ops->view=0; 143 ls->ops->setfromoptions=0; 144 ls->ops->reset=0; 145 ls->ops->destroy=0; 146 ls->setupcalled=PETSC_FALSE; 147 ls->usetaoroutines=PETSC_FALSE; 148 *newls = ls; 149 PetscFunctionReturn(0); 150 } 151 152 /*@ 153 TaoLineSearchSetUp - Sets up the internal data structures for the later use 154 of a Tao solver 155 156 Collective on ls 157 158 Input Parameters: 159 . ls - the TaoLineSearch context 160 161 Notes: 162 The user will not need to explicitly call TaoLineSearchSetUp(), as it will 163 automatically be called in TaoLineSearchSolve(). However, if the user 164 desires to call it explicitly, it should come after TaoLineSearchCreate() 165 but before TaoLineSearchApply(). 166 167 Level: developer 168 169 .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 170 @*/ 171 172 PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls) 173 { 174 PetscErrorCode ierr; 175 const char *default_type=TAOLINESEARCHMT; 176 PetscBool flg; 177 178 PetscFunctionBegin; 179 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 180 if (ls->setupcalled) PetscFunctionReturn(0); 181 if (!((PetscObject)ls)->type_name) { 182 ierr = TaoLineSearchSetType(ls,default_type);CHKERRQ(ierr); 183 } 184 if (ls->ops->setup) { 185 ierr = (*ls->ops->setup)(ls);CHKERRQ(ierr); 186 } 187 if (ls->usetaoroutines) { 188 ierr = TaoIsObjectiveDefined(ls->tao,&flg);CHKERRQ(ierr); 189 ls->hasobjective = flg; 190 ierr = TaoIsGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 191 ls->hasgradient = flg; 192 ierr = TaoIsObjectiveAndGradientDefined(ls->tao,&flg);CHKERRQ(ierr); 193 ls->hasobjectiveandgradient = flg; 194 } else { 195 if (ls->ops->computeobjective) { 196 ls->hasobjective = PETSC_TRUE; 197 } else { 198 ls->hasobjective = PETSC_FALSE; 199 } 200 if (ls->ops->computegradient) { 201 ls->hasgradient = PETSC_TRUE; 202 } else { 203 ls->hasgradient = PETSC_FALSE; 204 } 205 if (ls->ops->computeobjectiveandgradient) { 206 ls->hasobjectiveandgradient = PETSC_TRUE; 207 } else { 208 ls->hasobjectiveandgradient = PETSC_FALSE; 209 } 210 } 211 ls->setupcalled = PETSC_TRUE; 212 PetscFunctionReturn(0); 213 } 214 215 /*@ 216 TaoLineSearchReset - Some line searches may carry state information 217 from one TaoLineSearchApply() to the next. This function resets this 218 state information. 219 220 Collective on TaoLineSearch 221 222 Input Parameter: 223 . ls - the TaoLineSearch context 224 225 Level: developer 226 227 .seealso: TaoLineSearchCreate(), TaoLineSearchApply() 228 @*/ 229 PetscErrorCode TaoLineSearchReset(TaoLineSearch ls) 230 { 231 PetscErrorCode ierr; 232 233 PetscFunctionBegin; 234 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 235 if (ls->ops->reset) { 236 ierr = (*ls->ops->reset)(ls);CHKERRQ(ierr); 237 } 238 PetscFunctionReturn(0); 239 } 240 241 /*@ 242 TaoLineSearchDestroy - Destroys the TAO context that was created with 243 TaoLineSearchCreate() 244 245 Collective on TaoLineSearch 246 247 Input Parameter 248 . ls - the TaoLineSearch context 249 250 Level: beginner 251 252 .seealse: TaoLineSearchCreate(), TaoLineSearchSolve() 253 @*/ 254 PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls) 255 { 256 PetscErrorCode ierr; 257 258 PetscFunctionBegin; 259 if (!*ls) PetscFunctionReturn(0); 260 PetscValidHeaderSpecific(*ls,TAOLINESEARCH_CLASSID,1); 261 if (--((PetscObject)*ls)->refct > 0) {*ls=0; PetscFunctionReturn(0);} 262 ierr = VecDestroy(&(*ls)->stepdirection);CHKERRQ(ierr); 263 ierr = VecDestroy(&(*ls)->start_x);CHKERRQ(ierr); 264 if ((*ls)->ops->destroy) { 265 ierr = (*(*ls)->ops->destroy)(*ls);CHKERRQ(ierr); 266 } 267 ierr = PetscHeaderDestroy(ls);CHKERRQ(ierr); 268 PetscFunctionReturn(0); 269 } 270 271 /*@ 272 TaoLineSearchApply - Performs a line-search in a given step direction. Criteria for acceptable step length depends on the line-search algorithm chosen 273 274 Collective on TaoLineSearch 275 276 Input Parameters: 277 + ls - the Tao context 278 . x - The current solution (on output x contains the new solution determined by the line search) 279 . f - objective function value at current solution (on output contains the objective function value at new solution) 280 . g - gradient evaluated at x (on output contains the gradient at new solution) 281 - s - search direction 282 283 Output Parameters: 284 + x - new solution 285 . f - objective function value at x 286 . g - gradient vector at x 287 . steplength - scalar multiplier of s used ( x = x0 + steplength * x ) 288 - reason - reason why the line-search stopped 289 290 reason will be set to one of: 291 292 + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 293 . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 294 . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 295 . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 296 . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 297 . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 298 . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 299 . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 300 . TAOLINESEARCH_HALTED_OTHER - any other reason 301 - TAOLINESEARCH_SUCCESS - successful line search 302 303 Note: 304 The algorithm developer must set up the TaoLineSearch with calls to 305 TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines() 306 307 Note: 308 You may or may not need to follow this with a call to 309 TaoAddLineSearchCounts(), depending on whether you want these 310 evaluations to count toward the total function/gradient evaluations. 311 312 Level: beginner 313 314 .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts() 315 @*/ 316 317 PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 318 { 319 PetscErrorCode ierr; 320 PetscInt low1,low2,low3,high1,high2,high3; 321 322 PetscFunctionBegin; 323 *reason = TAOLINESEARCH_CONTINUE_ITERATING; 324 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 325 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 326 PetscValidRealPointer(f,3); 327 PetscValidHeaderSpecific(g,VEC_CLASSID,4); 328 PetscValidHeaderSpecific(s,VEC_CLASSID,5); 329 PetscValidPointer(reason,7); 330 PetscCheckSameComm(ls,1,x,2); 331 PetscCheckSameTypeAndComm(x,2,g,4); 332 PetscCheckSameTypeAndComm(x,2,s,5); 333 ierr = VecGetOwnershipRange(x, &low1, &high1);CHKERRQ(ierr); 334 ierr = VecGetOwnershipRange(g, &low2, &high2);CHKERRQ(ierr); 335 ierr = VecGetOwnershipRange(s, &low3, &high3);CHKERRQ(ierr); 336 if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths"); 337 338 ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr); 339 ierr = VecDestroy(&ls->stepdirection);CHKERRQ(ierr); 340 ls->stepdirection = s; 341 342 ierr = TaoLineSearchSetUp(ls);CHKERRQ(ierr); 343 if (!ls->ops->apply) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine"); 344 ls->nfeval=0; 345 ls->ngeval=0; 346 ls->nfgeval=0; 347 /* Check parameter values */ 348 if (ls->ftol < 0.0) { 349 ierr = PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);CHKERRQ(ierr); 350 *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 351 } 352 if (ls->rtol < 0.0) { 353 ierr = PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);CHKERRQ(ierr); 354 *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 355 } 356 if (ls->gtol < 0.0) { 357 ierr = PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);CHKERRQ(ierr); 358 *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 359 } 360 if (ls->stepmin < 0.0) { 361 ierr = PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);CHKERRQ(ierr); 362 *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 363 } 364 if (ls->stepmax < ls->stepmin) { 365 ierr = PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);CHKERRQ(ierr); 366 *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 367 } 368 if (ls->max_funcs < 0) { 369 ierr = PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);CHKERRQ(ierr); 370 *reason=TAOLINESEARCH_FAILED_BADPARAMETER; 371 } 372 if (PetscIsInfOrNanReal(*f)) { 373 ierr = PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);CHKERRQ(ierr); 374 *reason=TAOLINESEARCH_FAILED_INFORNAN; 375 } 376 377 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 378 ierr = VecDestroy(&ls->start_x);CHKERRQ(ierr); 379 ls->start_x = x; 380 381 ierr = PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0);CHKERRQ(ierr); 382 ierr = (*ls->ops->apply)(ls,x,f,g,s);CHKERRQ(ierr); 383 ierr = PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0);CHKERRQ(ierr); 384 *reason=ls->reason; 385 ls->new_f = *f; 386 387 if (steplength) { 388 *steplength=ls->step; 389 } 390 391 ierr = TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");CHKERRQ(ierr); 392 PetscFunctionReturn(0); 393 } 394 395 /*@C 396 TaoLineSearchSetType - Sets the algorithm used in a line search 397 398 Collective on TaoLineSearch 399 400 Input Parameters: 401 + ls - the TaoLineSearch context 402 - type - a known method 403 404 Available methods include: 405 + more-thuente 406 . gpcg 407 - unit - Do not perform any line search 408 409 410 Options Database Keys: 411 . -tao_ls_type - select which method TAO should use 412 413 Level: beginner 414 415 416 .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply() 417 418 @*/ 419 420 PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type) 421 { 422 PetscErrorCode ierr; 423 PetscErrorCode (*r)(TaoLineSearch); 424 PetscBool flg; 425 426 PetscFunctionBegin; 427 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 428 PetscValidCharPointer(type,2); 429 ierr = PetscObjectTypeCompare((PetscObject)ls, type, &flg);CHKERRQ(ierr); 430 if (flg) PetscFunctionReturn(0); 431 432 ierr = PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);CHKERRQ(ierr); 433 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type); 434 if (ls->ops->destroy) { 435 ierr = (*(ls)->ops->destroy)(ls);CHKERRQ(ierr); 436 } 437 ls->max_funcs=30; 438 ls->ftol = 0.0001; 439 ls->gtol = 0.9; 440 #if defined(PETSC_USE_REAL_SINGLE) 441 ls->rtol = 1.0e-5; 442 #else 443 ls->rtol = 1.0e-10; 444 #endif 445 ls->stepmin=1.0e-20; 446 ls->stepmax=1.0e+20; 447 448 ls->nfeval=0; 449 ls->ngeval=0; 450 ls->nfgeval=0; 451 ls->ops->setup=0; 452 ls->ops->apply=0; 453 ls->ops->view=0; 454 ls->ops->setfromoptions=0; 455 ls->ops->destroy=0; 456 ls->setupcalled = PETSC_FALSE; 457 ierr = (*r)(ls);CHKERRQ(ierr); 458 ierr = PetscObjectChangeTypeName((PetscObject)ls, type);CHKERRQ(ierr); 459 PetscFunctionReturn(0); 460 } 461 462 /*@ 463 TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user 464 options. 465 466 Collective on TaoLineSearch 467 468 Input Paremeter: 469 . ls - the TaoLineSearch context 470 471 Options Database Keys: 472 + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit) 473 . -tao_ls_ftol <tol> - tolerance for sufficient decrease 474 . -tao_ls_gtol <tol> - tolerance for curvature condition 475 . -tao_ls_rtol <tol> - relative tolerance for acceptable step 476 . -tao_ls_stepmin <step> - minimum steplength allowed 477 . -tao_ls_stepmax <step> - maximum steplength allowed 478 . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed 479 - -tao_ls_view - display line-search results to standard output 480 481 Level: beginner 482 @*/ 483 PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls) 484 { 485 PetscErrorCode ierr; 486 const char *default_type=TAOLINESEARCHMT; 487 char type[256]; 488 PetscBool flg; 489 490 PetscFunctionBegin; 491 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 492 ierr = PetscObjectOptionsBegin((PetscObject)ls);CHKERRQ(ierr); 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, 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: 810 TaoLineSearchComputeObjective() is typically used within line searches 811 so most users would not generally call this routine themselves. 812 813 Level: developer 814 815 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 816 @*/ 817 PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 818 { 819 PetscErrorCode ierr; 820 Vec gdummy; 821 PetscReal gts; 822 823 PetscFunctionBegin; 824 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 825 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 826 PetscValidPointer(f,3); 827 PetscCheckSameComm(ls,1,x,2); 828 if (ls->usetaoroutines) { 829 ierr = TaoComputeObjective(ls->tao,x,f);CHKERRQ(ierr); 830 } else { 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 ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 833 PetscStackPush("TaoLineSearch user objective routine"); 834 if (ls->ops->computeobjective) { 835 ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 836 } else if (ls->ops->computeobjectiveandgradient) { 837 ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr); 838 ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr); 839 ierr = VecDestroy(&gdummy);CHKERRQ(ierr); 840 } else { 841 ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts);CHKERRQ(ierr); 842 } 843 PetscStackPop; 844 ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 845 } 846 ls->nfeval++; 847 PetscFunctionReturn(0); 848 } 849 850 /*@ 851 TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 852 853 Collective on Tao 854 855 Input Parameters: 856 + ls - the TaoLineSearch context 857 - x - input vector 858 859 Output Parameter: 860 + f - Objective value at X 861 - g - Gradient vector at X 862 863 Notes: 864 TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 865 so most users would not generally call this routine themselves. 866 867 Level: developer 868 869 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 870 @*/ 871 PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 872 { 873 PetscErrorCode ierr; 874 875 PetscFunctionBegin; 876 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 877 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 878 PetscValidPointer(f,3); 879 PetscValidHeaderSpecific(g,VEC_CLASSID,4); 880 PetscCheckSameComm(ls,1,x,2); 881 PetscCheckSameComm(ls,1,g,4); 882 if (ls->usetaoroutines) { 883 ierr = TaoComputeObjectiveAndGradient(ls->tao,x,f,g);CHKERRQ(ierr); 884 } else { 885 if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 886 if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 887 ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 888 PetscStackPush("TaoLineSearch user objective/gradient routine"); 889 if (ls->ops->computeobjectiveandgradient) { 890 ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr); 891 } else { 892 ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 893 ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 894 } 895 PetscStackPop; 896 ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 897 ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 898 } 899 ls->nfgeval++; 900 PetscFunctionReturn(0); 901 } 902 903 /*@ 904 TaoLineSearchComputeGradient - Computes the gradient of the objective function 905 906 Collective on TaoLineSearch 907 908 Input Parameters: 909 + ls - the TaoLineSearch context 910 - x - input vector 911 912 Output Parameter: 913 . g - gradient vector 914 915 Notes: 916 TaoComputeGradient() is typically used within line searches 917 so most users would not generally call this routine themselves. 918 919 Level: developer 920 921 .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 922 @*/ 923 PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 924 { 925 PetscErrorCode ierr; 926 PetscReal fdummy; 927 928 PetscFunctionBegin; 929 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 930 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 931 PetscValidHeaderSpecific(g,VEC_CLASSID,3); 932 PetscCheckSameComm(ls,1,x,2); 933 PetscCheckSameComm(ls,1,g,3); 934 if (ls->usetaoroutines) { 935 ierr = TaoComputeGradient(ls->tao,x,g);CHKERRQ(ierr); 936 } else { 937 if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 938 ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 939 PetscStackPush("TaoLineSearch user gradient routine"); 940 if (ls->ops->computegradient) { 941 ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 942 } else { 943 ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr); 944 } 945 PetscStackPop; 946 ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 947 } 948 ls->ngeval++; 949 PetscFunctionReturn(0); 950 } 951 952 /*@ 953 TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 954 955 Collective on Tao 956 957 Input Parameters: 958 + ls - the TaoLineSearch context 959 - x - input vector 960 961 Output Parameter: 962 + f - Objective value at X 963 - gts - inner product of gradient and step direction at X 964 965 Notes: 966 TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 967 so most users would not generally call this routine themselves. 968 969 Level: developer 970 971 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 972 @*/ 973 PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 974 { 975 PetscErrorCode ierr; 976 PetscFunctionBegin; 977 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 978 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 979 PetscValidPointer(f,3); 980 PetscValidPointer(gts,4); 981 PetscCheckSameComm(ls,1,x,2); 982 if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 983 ierr = PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 984 PetscStackPush("TaoLineSearch user objective/gts routine"); 985 ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr); 986 PetscStackPop; 987 ierr = PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);CHKERRQ(ierr); 988 ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 989 ls->nfeval++; 990 PetscFunctionReturn(0); 991 } 992 993 /*@ 994 TaoLineSearchGetSolution - Returns the solution to the line search 995 996 Collective on TaoLineSearch 997 998 Input Parameter: 999 . ls - the TaoLineSearch context 1000 1001 Output Parameter: 1002 + x - the new solution 1003 . f - the objective function value at x 1004 . g - the gradient at x 1005 . steplength - the multiple of the step direction taken by the line search 1006 - reason - the reason why the line search terminated 1007 1008 reason will be set to one of: 1009 1010 + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 1011 . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 1012 . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 1013 . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 1014 . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 1015 . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 1016 . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 1017 1018 . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 1019 . TAOLINESEARCH_HALTED_OTHER - any other reason 1020 1021 + TAOLINESEARCH_SUCCESS - successful line search 1022 1023 Level: developer 1024 1025 @*/ 1026 PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 1027 { 1028 PetscErrorCode ierr; 1029 1030 PetscFunctionBegin; 1031 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1032 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1033 PetscValidPointer(f,3); 1034 PetscValidHeaderSpecific(g,VEC_CLASSID,4); 1035 PetscValidIntPointer(reason,6); 1036 if (ls->new_x) { 1037 ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr); 1038 } 1039 *f = ls->new_f; 1040 if (ls->new_g) { 1041 ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr); 1042 } 1043 if (steplength) { 1044 *steplength=ls->step; 1045 } 1046 *reason = ls->reason; 1047 PetscFunctionReturn(0); 1048 } 1049 1050 /*@ 1051 TaoLineSearchGetStartingVector - Gets a the initial point of the line 1052 search. 1053 1054 Not Collective 1055 1056 Input Parameter: 1057 . ls - the TaoLineSearch context 1058 1059 Output Parameter: 1060 . x - The initial point of the line search 1061 1062 Level: intermediate 1063 @*/ 1064 PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 1065 { 1066 PetscFunctionBegin; 1067 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1068 if (x) { 1069 *x = ls->start_x; 1070 } 1071 PetscFunctionReturn(0); 1072 } 1073 1074 /*@ 1075 TaoLineSearchGetStepDirection - Gets the step direction of the line 1076 search. 1077 1078 Not Collective 1079 1080 Input Parameter: 1081 . ls - the TaoLineSearch context 1082 1083 Output Parameter: 1084 . s - the step direction of the line search 1085 1086 Level: advanced 1087 @*/ 1088 PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 1089 { 1090 PetscFunctionBegin; 1091 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1092 if (s) { 1093 *s = ls->stepdirection; 1094 } 1095 PetscFunctionReturn(0); 1096 1097 } 1098 1099 /*@ 1100 TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 1101 1102 Not Collective 1103 1104 Input Parameter: 1105 . ls - the TaoLineSearch context 1106 1107 Output Parameter: 1108 . f - the objective value at the full step length 1109 1110 Level: developer 1111 @*/ 1112 1113 PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 1114 { 1115 PetscFunctionBegin; 1116 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1117 *f_fullstep = ls->f_fullstep; 1118 PetscFunctionReturn(0); 1119 } 1120 1121 /*@ 1122 TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 1123 1124 Logically Collective on Tao 1125 1126 Input Parameters: 1127 + ls - the TaoLineSearch context 1128 . xl - vector of lower bounds 1129 - xu - vector of upper bounds 1130 1131 Note: If the variable bounds are not set with this routine, then 1132 PETSC_NINFINITY and PETSC_INFINITY are assumed 1133 1134 Level: beginner 1135 1136 .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 1137 @*/ 1138 PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 1139 { 1140 PetscFunctionBegin; 1141 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1142 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1143 PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1144 ls->lower = xl; 1145 ls->upper = xu; 1146 ls->bounded = 1; 1147 PetscFunctionReturn(0); 1148 } 1149 1150 /*@ 1151 TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 1152 search. If this value is not set then 1.0 is assumed. 1153 1154 Logically Collective on TaoLineSearch 1155 1156 Input Parameters: 1157 + ls - the TaoLineSearch context 1158 - s - the initial step size 1159 1160 Level: intermediate 1161 1162 .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply() 1163 @*/ 1164 PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 1165 { 1166 PetscFunctionBegin; 1167 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1168 ls->initstep = s; 1169 PetscFunctionReturn(0); 1170 } 1171 1172 /*@ 1173 TaoLineSearchGetStepLength - Get the current step length 1174 1175 Not Collective 1176 1177 Input Parameters: 1178 . ls - the TaoLineSearch context 1179 1180 Output Parameters 1181 . s - the current step length 1182 1183 Level: beginner 1184 1185 .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply() 1186 @*/ 1187 PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 1188 { 1189 PetscFunctionBegin; 1190 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1191 *s = ls->step; 1192 PetscFunctionReturn(0); 1193 } 1194 1195 /*@C 1196 TaoLineSearchRegister - Adds a line-search algorithm to the registry 1197 1198 Not collective 1199 1200 Input Parameters: 1201 + sname - name of a new user-defined solver 1202 - func - routine to Create method context 1203 1204 Notes: 1205 TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 1206 1207 Sample usage: 1208 .vb 1209 TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 1210 .ve 1211 1212 Then, your solver can be chosen with the procedural interface via 1213 $ TaoLineSearchSetType(ls,"my_linesearch") 1214 or at runtime via the option 1215 $ -tao_ls_type my_linesearch 1216 1217 Level: developer 1218 1219 @*/ 1220 PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 1221 { 1222 PetscErrorCode ierr; 1223 PetscFunctionBegin; 1224 ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr); 1225 PetscFunctionReturn(0); 1226 } 1227 1228 /*@C 1229 TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 1230 for all TaoLineSearch options in the database. 1231 1232 1233 Collective on TaoLineSearch 1234 1235 Input Parameters: 1236 + ls - the TaoLineSearch solver context 1237 - prefix - the prefix string to prepend to all line search requests 1238 1239 Notes: 1240 A hyphen (-) must NOT be given at the beginning of the prefix name. 1241 The first character of all runtime options is AUTOMATICALLY the hyphen. 1242 1243 1244 Level: advanced 1245 1246 .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 1247 @*/ 1248 PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 1249 { 1250 return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 1251 } 1252 1253 /*@C 1254 TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1255 TaoLineSearch options in the database 1256 1257 Not Collective 1258 1259 Input Parameters: 1260 . ls - the TaoLineSearch context 1261 1262 Output Parameters: 1263 . prefix - pointer to the prefix string used is returned 1264 1265 Notes: 1266 On the fortran side, the user should pass in a string 'prefix' of 1267 sufficient length to hold the prefix. 1268 1269 Level: advanced 1270 1271 .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 1272 @*/ 1273 PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 1274 { 1275 return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 1276 } 1277 1278 /*@C 1279 TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 1280 TaoLineSearch options in the database. 1281 1282 1283 Logically Collective on TaoLineSearch 1284 1285 Input Parameters: 1286 + ls - the TaoLineSearch context 1287 - prefix - the prefix string to prepend to all TAO option requests 1288 1289 Notes: 1290 A hyphen (-) must NOT be given at the beginning of the prefix name. 1291 The first character of all runtime options is AUTOMATICALLY the hyphen. 1292 1293 For example, to distinguish between the runtime options for two 1294 different line searches, one could call 1295 .vb 1296 TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 1297 TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 1298 .ve 1299 1300 This would enable use of different options for each system, such as 1301 .vb 1302 -sys1_tao_ls_type mt 1303 -sys2_tao_ls_type armijo 1304 .ve 1305 1306 Level: advanced 1307 1308 .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 1309 @*/ 1310 1311 PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 1312 { 1313 return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 1314 } 1315