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