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