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