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