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