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