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: 814 TaoLineSearchComputeObjective() is typically used within line searches 815 so most users would not generally call this routine themselves. 816 817 Level: developer 818 819 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 820 @*/ 821 PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f) 822 { 823 PetscErrorCode ierr; 824 Vec gdummy; 825 PetscReal gts; 826 827 PetscFunctionBegin; 828 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 829 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 830 PetscValidPointer(f,3); 831 PetscCheckSameComm(ls,1,x,2); 832 if (ls->usetaoroutines) { 833 ierr = TaoComputeObjective(ls->tao,x,f);CHKERRQ(ierr); 834 } else { 835 ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 836 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"); 837 PetscStackPush("TaoLineSearch user objective routine"); 838 if (ls->ops->computeobjective) { 839 ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 840 } else if (ls->ops->computeobjectiveandgradient) { 841 ierr = VecDuplicate(x,&gdummy);CHKERRQ(ierr); 842 ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);CHKERRQ(ierr); 843 ierr = VecDestroy(&gdummy);CHKERRQ(ierr); 844 } else { 845 ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,>s,ls->userctx_funcgts);CHKERRQ(ierr); 846 } 847 PetscStackPop; 848 ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 849 } 850 ls->nfeval++; 851 PetscFunctionReturn(0); 852 } 853 854 /*@ 855 TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point 856 857 Collective on Tao 858 859 Input Parameters: 860 + ls - the TaoLineSearch context 861 - x - input vector 862 863 Output Parameter: 864 + f - Objective value at X 865 - g - Gradient vector at X 866 867 Notes: 868 TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches 869 so most users would not generally call this routine themselves. 870 871 Level: developer 872 873 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 874 @*/ 875 PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g) 876 { 877 PetscErrorCode ierr; 878 879 PetscFunctionBegin; 880 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 881 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 882 PetscValidPointer(f,3); 883 PetscValidHeaderSpecific(g,VEC_CLASSID,4); 884 PetscCheckSameComm(ls,1,x,2); 885 PetscCheckSameComm(ls,1,g,4); 886 if (ls->usetaoroutines) { 887 ierr = TaoComputeObjectiveAndGradient(ls->tao,x,f,g);CHKERRQ(ierr); 888 } else { 889 ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 890 if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set"); 891 if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set"); 892 893 PetscStackPush("TaoLineSearch user objective/gradient routine"); 894 if (ls->ops->computeobjectiveandgradient) { 895 ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);CHKERRQ(ierr); 896 } else { 897 ierr = (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);CHKERRQ(ierr); 898 ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 899 } 900 PetscStackPop; 901 ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 902 ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 903 } 904 ls->nfgeval++; 905 PetscFunctionReturn(0); 906 } 907 908 /*@ 909 TaoLineSearchComputeGradient - Computes the gradient of the objective function 910 911 Collective on TaoLineSearch 912 913 Input Parameters: 914 + ls - the TaoLineSearch context 915 - x - input vector 916 917 Output Parameter: 918 . g - gradient vector 919 920 Notes: 921 TaoComputeGradient() is typically used within line searches 922 so most users would not generally call this routine themselves. 923 924 Level: developer 925 926 .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient() 927 @*/ 928 PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g) 929 { 930 PetscErrorCode ierr; 931 PetscReal fdummy; 932 933 PetscFunctionBegin; 934 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 935 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 936 PetscValidHeaderSpecific(g,VEC_CLASSID,3); 937 PetscCheckSameComm(ls,1,x,2); 938 PetscCheckSameComm(ls,1,g,3); 939 if (ls->usetaoroutines) { 940 ierr = TaoComputeGradient(ls->tao,x,g);CHKERRQ(ierr); 941 } else { 942 ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 943 if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set"); 944 PetscStackPush("TaoLineSearch user gradient routine"); 945 if (ls->ops->computegradient) { 946 ierr = (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);CHKERRQ(ierr); 947 } else { 948 ierr = (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);CHKERRQ(ierr); 949 } 950 PetscStackPop; 951 ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 952 } 953 ls->ngeval++; 954 PetscFunctionReturn(0); 955 } 956 957 /*@ 958 TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point 959 960 Collective on Tao 961 962 Input Parameters: 963 + ls - the TaoLineSearch context 964 - x - input vector 965 966 Output Parameter: 967 + f - Objective value at X 968 - gts - inner product of gradient and step direction at X 969 970 Notes: 971 TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches 972 so most users would not generally call this routine themselves. 973 974 Level: developer 975 976 .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine() 977 @*/ 978 PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts) 979 { 980 PetscErrorCode ierr; 981 PetscFunctionBegin; 982 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 983 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 984 PetscValidPointer(f,3); 985 PetscValidPointer(gts,4); 986 PetscCheckSameComm(ls,1,x,2); 987 ierr = PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 988 if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set"); 989 PetscStackPush("TaoLineSearch user objective/gts routine"); 990 ierr = (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);CHKERRQ(ierr); 991 PetscStackPop; 992 ierr = PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);CHKERRQ(ierr); 993 ierr = PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));CHKERRQ(ierr); 994 ls->nfeval++; 995 PetscFunctionReturn(0); 996 } 997 998 /*@ 999 TaoLineSearchGetSolution - Returns the solution to the line search 1000 1001 Collective on TaoLineSearch 1002 1003 Input Parameter: 1004 . ls - the TaoLineSearch context 1005 1006 Output Parameter: 1007 + x - the new solution 1008 . f - the objective function value at x 1009 . g - the gradient at x 1010 . steplength - the multiple of the step direction taken by the line search 1011 - reason - the reason why the line search terminated 1012 1013 reason will be set to one of: 1014 1015 + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value 1016 . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter 1017 . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction 1018 . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached 1019 . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound 1020 . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound 1021 . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance 1022 1023 . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search 1024 . TAOLINESEARCH_HALTED_OTHER - any other reason 1025 1026 + TAOLINESEARCH_SUCCESS - successful line search 1027 1028 Level: developer 1029 1030 @*/ 1031 PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason) 1032 { 1033 PetscErrorCode ierr; 1034 1035 PetscFunctionBegin; 1036 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1037 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1038 PetscValidPointer(f,3); 1039 PetscValidHeaderSpecific(g,VEC_CLASSID,4); 1040 PetscValidIntPointer(reason,6); 1041 if (ls->new_x) { 1042 ierr = VecCopy(ls->new_x,x);CHKERRQ(ierr); 1043 } 1044 *f = ls->new_f; 1045 if (ls->new_g) { 1046 ierr = VecCopy(ls->new_g,g);CHKERRQ(ierr); 1047 } 1048 if (steplength) { 1049 *steplength=ls->step; 1050 } 1051 *reason = ls->reason; 1052 PetscFunctionReturn(0); 1053 } 1054 1055 /*@ 1056 TaoLineSearchGetStartingVector - Gets a the initial point of the line 1057 search. 1058 1059 Not Collective 1060 1061 Input Parameter: 1062 . ls - the TaoLineSearch context 1063 1064 Output Parameter: 1065 . x - The initial point of the line search 1066 1067 Level: intermediate 1068 @*/ 1069 PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x) 1070 { 1071 PetscFunctionBegin; 1072 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1073 if (x) { 1074 *x = ls->start_x; 1075 } 1076 PetscFunctionReturn(0); 1077 } 1078 1079 /*@ 1080 TaoLineSearchGetStepDirection - Gets the step direction of the line 1081 search. 1082 1083 Not Collective 1084 1085 Input Parameter: 1086 . ls - the TaoLineSearch context 1087 1088 Output Parameter: 1089 . s - the step direction of the line search 1090 1091 Level: advanced 1092 @*/ 1093 PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s) 1094 { 1095 PetscFunctionBegin; 1096 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1097 if (s) { 1098 *s = ls->stepdirection; 1099 } 1100 PetscFunctionReturn(0); 1101 1102 } 1103 1104 /*@ 1105 TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms. 1106 1107 Not Collective 1108 1109 Input Parameter: 1110 . ls - the TaoLineSearch context 1111 1112 Output Parameter: 1113 . f - the objective value at the full step length 1114 1115 Level: developer 1116 @*/ 1117 1118 PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep) 1119 { 1120 PetscFunctionBegin; 1121 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1122 *f_fullstep = ls->f_fullstep; 1123 PetscFunctionReturn(0); 1124 } 1125 1126 /*@ 1127 TaoLineSearchSetVariableBounds - Sets the upper and lower bounds. 1128 1129 Logically Collective on Tao 1130 1131 Input Parameters: 1132 + ls - the TaoLineSearch context 1133 . xl - vector of lower bounds 1134 - xu - vector of upper bounds 1135 1136 Note: If the variable bounds are not set with this routine, then 1137 PETSC_NINFINITY and PETSC_INFINITY are assumed 1138 1139 Level: beginner 1140 1141 .seealso: TaoSetVariableBounds(), TaoLineSearchCreate() 1142 @*/ 1143 PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu) 1144 { 1145 PetscFunctionBegin; 1146 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1147 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1148 PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1149 ls->lower = xl; 1150 ls->upper = xu; 1151 ls->bounded = 1; 1152 PetscFunctionReturn(0); 1153 } 1154 1155 /*@ 1156 TaoLineSearchSetInitialStepLength - Sets the initial step length of a line 1157 search. If this value is not set then 1.0 is assumed. 1158 1159 Logically Collective on TaoLineSearch 1160 1161 Input Parameters: 1162 + ls - the TaoLineSearch context 1163 - s - the initial step size 1164 1165 Level: intermediate 1166 1167 .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply() 1168 @*/ 1169 PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s) 1170 { 1171 PetscFunctionBegin; 1172 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1173 ls->initstep = s; 1174 PetscFunctionReturn(0); 1175 } 1176 1177 /*@ 1178 TaoLineSearchGetStepLength - Get the current step length 1179 1180 Not Collective 1181 1182 Input Parameters: 1183 . ls - the TaoLineSearch context 1184 1185 Output Parameters 1186 . s - the current step length 1187 1188 Level: beginner 1189 1190 .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply() 1191 @*/ 1192 PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s) 1193 { 1194 PetscFunctionBegin; 1195 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 1196 *s = ls->step; 1197 PetscFunctionReturn(0); 1198 } 1199 1200 /*MC 1201 TaoLineSearchRegister - Adds a line-search algorithm to the registry 1202 1203 Not collective 1204 1205 Input Parameters: 1206 + sname - name of a new user-defined solver 1207 - func - routine to Create method context 1208 1209 Notes: 1210 TaoLineSearchRegister() may be called multiple times to add several user-defined solvers. 1211 1212 Sample usage: 1213 .vb 1214 TaoLineSearchRegister("my_linesearch",MyLinesearchCreate); 1215 .ve 1216 1217 Then, your solver can be chosen with the procedural interface via 1218 $ TaoLineSearchSetType(ls,"my_linesearch") 1219 or at runtime via the option 1220 $ -tao_ls_type my_linesearch 1221 1222 Level: developer 1223 1224 .seealso: TaoLineSearchRegisterDestroy() 1225 M*/ 1226 PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch)) 1227 { 1228 PetscErrorCode ierr; 1229 PetscFunctionBegin; 1230 ierr = PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);CHKERRQ(ierr); 1231 PetscFunctionReturn(0); 1232 } 1233 1234 /*@C 1235 TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were 1236 registered by TaoLineSearchRegister(). 1237 1238 Not Collective 1239 1240 Level: developer 1241 1242 .seealso: TaoLineSearchRegister() 1243 @*/ 1244 PetscErrorCode TaoLineSearchRegisterDestroy(void) 1245 { 1246 PetscErrorCode ierr; 1247 PetscFunctionBegin; 1248 ierr = PetscFunctionListDestroy(&TaoLineSearchList);CHKERRQ(ierr); 1249 TaoLineSearchInitialized = PETSC_FALSE; 1250 PetscFunctionReturn(0); 1251 } 1252 1253 /*@C 1254 TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching 1255 for all TaoLineSearch options in the database. 1256 1257 1258 Collective on TaoLineSearch 1259 1260 Input Parameters: 1261 + ls - the TaoLineSearch solver context 1262 - prefix - the prefix string to prepend to all line search requests 1263 1264 Notes: 1265 A hyphen (-) must NOT be given at the beginning of the prefix name. 1266 The first character of all runtime options is AUTOMATICALLY the hyphen. 1267 1268 1269 Level: advanced 1270 1271 .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 1272 @*/ 1273 PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[]) 1274 { 1275 return PetscObjectAppendOptionsPrefix((PetscObject)ls,p); 1276 } 1277 1278 /*@C 1279 TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1280 TaoLineSearch options in the database 1281 1282 Not Collective 1283 1284 Input Parameters: 1285 . ls - the TaoLineSearch context 1286 1287 Output Parameters: 1288 . prefix - pointer to the prefix string used is returned 1289 1290 Notes: 1291 On the fortran side, the user should pass in a string 'prefix' of 1292 sufficient length to hold the prefix. 1293 1294 Level: advanced 1295 1296 .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix() 1297 @*/ 1298 PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[]) 1299 { 1300 return PetscObjectGetOptionsPrefix((PetscObject)ls,p); 1301 } 1302 1303 /*@C 1304 TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all 1305 TaoLineSearch options in the database. 1306 1307 1308 Logically Collective on TaoLineSearch 1309 1310 Input Parameters: 1311 + ls - the TaoLineSearch context 1312 - prefix - the prefix string to prepend to all TAO option requests 1313 1314 Notes: 1315 A hyphen (-) must NOT be given at the beginning of the prefix name. 1316 The first character of all runtime options is AUTOMATICALLY the hyphen. 1317 1318 For example, to distinguish between the runtime options for two 1319 different line searches, one could call 1320 .vb 1321 TaoLineSearchSetOptionsPrefix(ls1,"sys1_") 1322 TaoLineSearchSetOptionsPrefix(ls2,"sys2_") 1323 .ve 1324 1325 This would enable use of different options for each system, such as 1326 .vb 1327 -sys1_tao_ls_type mt 1328 -sys2_tao_ls_type armijo 1329 .ve 1330 1331 Level: advanced 1332 1333 .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix() 1334 @*/ 1335 1336 PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[]) 1337 { 1338 return PetscObjectSetOptionsPrefix((PetscObject)ls,p); 1339 } 1340