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