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