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