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