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