1 2 #include <petsc/private/taolinesearchimpl.h> 3 #include <../src/tao/linesearch/impls/morethuente/morethuente.h> 4 5 /* 6 This algorithm is taken from More' and Thuente, "Line search algorithms 7 with guaranteed sufficient decrease", Argonne National Laboratory, 8 Technical Report MCS-P330-1092. 9 */ 10 11 static PetscErrorCode Tao_mcstep(TaoLineSearch ls, PetscReal *stx, PetscReal *fx, PetscReal *dx, PetscReal *sty, PetscReal *fy, PetscReal *dy, PetscReal *stp, PetscReal *fp, PetscReal *dp); 12 13 static PetscErrorCode TaoLineSearchDestroy_MT(TaoLineSearch ls) 14 { 15 TaoLineSearch_MT *mt = (TaoLineSearch_MT *)(ls->data); 16 17 PetscFunctionBegin; 18 PetscCall(PetscObjectDereference((PetscObject)mt->x)); 19 PetscCall(VecDestroy(&mt->work)); 20 PetscCall(PetscFree(ls->data)); 21 PetscFunctionReturn(PETSC_SUCCESS); 22 } 23 24 static PetscErrorCode TaoLineSearchSetFromOptions_MT(TaoLineSearch ls, PetscOptionItems *PetscOptionsObject) 25 { 26 PetscFunctionBegin; 27 PetscFunctionReturn(PETSC_SUCCESS); 28 } 29 30 static PetscErrorCode TaoLineSearchMonitor_MT(TaoLineSearch ls) 31 { 32 TaoLineSearch_MT *mt = (TaoLineSearch_MT *)ls->data; 33 34 PetscFunctionBegin; 35 PetscCall(PetscViewerASCIIPrintf(ls->viewer, "stx: %g, fx: %g, dgx: %g\n", (double)mt->stx, (double)mt->fx, (double)mt->dgx)); 36 PetscCall(PetscViewerASCIIPrintf(ls->viewer, "sty: %g, fy: %g, dgy: %g\n", (double)mt->sty, (double)mt->fy, (double)mt->dgy)); 37 PetscFunctionReturn(PETSC_SUCCESS); 38 } 39 40 static PetscErrorCode TaoLineSearchApply_MT(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 41 { 42 TaoLineSearch_MT *mt = (TaoLineSearch_MT *)(ls->data); 43 PetscReal xtrapf = 4.0; 44 PetscReal finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym; 45 PetscReal dgx, dgy, dg, dg2, fx, fy, stx, sty, dgtest; 46 PetscReal ftest1 = 0.0, ftest2 = 0.0; 47 PetscInt i, stage1, n1, n2, nn1, nn2; 48 PetscReal bstepmin1, bstepmin2, bstepmax, ostepmin, ostepmax; 49 PetscBool g_computed = PETSC_FALSE; /* to prevent extra gradient computation */ 50 51 PetscFunctionBegin; 52 ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 53 PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0)); 54 /* Check work vector */ 55 if (!mt->work) { 56 PetscCall(VecDuplicate(x, &mt->work)); 57 mt->x = x; 58 PetscCall(PetscObjectReference((PetscObject)mt->x)); 59 } else if (x != mt->x) { 60 PetscCall(VecDestroy(&mt->work)); 61 PetscCall(VecDuplicate(x, &mt->work)); 62 PetscCall(PetscObjectDereference((PetscObject)mt->x)); 63 mt->x = x; 64 PetscCall(PetscObjectReference((PetscObject)mt->x)); 65 } 66 67 ostepmax = ls->stepmax; 68 ostepmin = ls->stepmin; 69 70 if (ls->bounded) { 71 /* Compute step length needed to make all variables equal a bound */ 72 /* Compute the smallest steplength that will make one nonbinding variable 73 equal the bound */ 74 PetscCall(VecGetLocalSize(ls->upper, &n1)); 75 PetscCall(VecGetLocalSize(mt->x, &n2)); 76 PetscCall(VecGetSize(ls->upper, &nn1)); 77 PetscCall(VecGetSize(mt->x, &nn2)); 78 PetscCheck(n1 == n2 && nn1 == nn2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Variable vector not compatible with bounds vector"); 79 PetscCall(VecScale(s, -1.0)); 80 PetscCall(VecBoundGradientProjection(s, x, ls->lower, ls->upper, s)); 81 PetscCall(VecScale(s, -1.0)); 82 PetscCall(VecStepBoundInfo(x, s, ls->lower, ls->upper, &bstepmin1, &bstepmin2, &bstepmax)); 83 ls->stepmax = PetscMin(bstepmax, ls->stepmax); 84 } 85 86 PetscCall(VecDot(g, s, &dginit)); 87 if (PetscIsInfOrNanReal(dginit)) { 88 PetscCall(PetscInfo(ls, "Initial Line Search step * g is Inf or Nan (%g)\n", (double)dginit)); 89 ls->reason = TAOLINESEARCH_FAILED_INFORNAN; 90 PetscFunctionReturn(PETSC_SUCCESS); 91 } 92 if (dginit >= 0.0) { 93 PetscCall(PetscInfo(ls, "Initial Line Search step * g is not descent direction (%g)\n", (double)dginit)); 94 ls->reason = TAOLINESEARCH_FAILED_ASCENT; 95 PetscFunctionReturn(PETSC_SUCCESS); 96 } 97 98 /* Initialization */ 99 mt->bracket = 0; 100 stage1 = 1; 101 finit = *f; 102 dgtest = ls->ftol * dginit; 103 width = ls->stepmax - ls->stepmin; 104 width1 = width * 2.0; 105 PetscCall(VecCopy(x, mt->work)); 106 /* Variable dictionary: 107 stx, fx, dgx - the step, function, and derivative at the best step 108 sty, fy, dgy - the step, function, and derivative at the other endpoint 109 of the interval of uncertainty 110 step, f, dg - the step, function, and derivative at the current step */ 111 112 stx = 0.0; 113 fx = finit; 114 dgx = dginit; 115 sty = 0.0; 116 fy = finit; 117 dgy = dginit; 118 119 ls->step = ls->initstep; 120 for (i = 0; i < ls->max_funcs; i++) { 121 /* Set min and max steps to correspond to the interval of uncertainty */ 122 if (mt->bracket) { 123 ls->stepmin = PetscMin(stx, sty); 124 ls->stepmax = PetscMax(stx, sty); 125 } else { 126 ls->stepmin = stx; 127 ls->stepmax = ls->step + xtrapf * (ls->step - stx); 128 } 129 130 /* Force the step to be within the bounds */ 131 ls->step = PetscMax(ls->step, ls->stepmin); 132 ls->step = PetscMin(ls->step, ls->stepmax); 133 134 /* If an unusual termination is to occur, then let step be the lowest 135 point obtained thus far */ 136 if (stx != 0 && ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) || (ls->nfeval + ls->nfgeval >= ls->max_funcs - 1) || mt->infoc == 0)) 137 ls->step = stx; 138 139 PetscCall(VecWAXPY(mt->work, ls->step, s, x)); /* W = X + step*S */ 140 141 if (ls->step == 0.0) { 142 PetscCall(PetscInfo(ls, "Step size is zero.\n")); 143 ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND; 144 break; 145 } 146 147 if (ls->bounded) PetscCall(VecMedian(ls->lower, mt->work, ls->upper, mt->work)); 148 if (ls->usegts) { 149 PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls, mt->work, f, &dg)); 150 g_computed = PETSC_FALSE; 151 } else { 152 PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, mt->work, f, g)); 153 g_computed = PETSC_TRUE; 154 if (ls->bounded) { 155 PetscCall(VecDot(g, x, &dg)); 156 PetscCall(VecDot(g, mt->work, &dg2)); 157 dg = (dg2 - dg) / ls->step; 158 } else { 159 PetscCall(VecDot(g, s, &dg)); 160 } 161 } 162 163 /* update bracketing parameters in the MT context for printouts in monitor */ 164 mt->stx = stx; 165 mt->fx = fx; 166 mt->dgx = dgx; 167 mt->sty = sty; 168 mt->fy = fy; 169 mt->dgy = dgy; 170 PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step)); 171 172 if (i == 0) ls->f_fullstep = *f; 173 174 if (PetscIsInfOrNanReal(*f) || PetscIsInfOrNanReal(dg)) { 175 /* User provided compute function generated Not-a-Number, assume 176 domain violation and set function value and directional 177 derivative to infinity. */ 178 *f = PETSC_INFINITY; 179 dg = PETSC_INFINITY; 180 } 181 182 ftest1 = finit + ls->step * dgtest; 183 if (ls->bounded) ftest2 = finit + ls->step * dgtest * ls->ftol; 184 185 /* Convergence testing */ 186 if ((*f - ftest1 <= PETSC_SMALL * PetscAbsReal(finit)) && (PetscAbsReal(dg) + ls->gtol * dginit <= 0.0)) { 187 PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease and directional deriv conditions hold\n")); 188 ls->reason = TAOLINESEARCH_SUCCESS; 189 break; 190 } 191 192 /* Check Armijo if beyond the first breakpoint */ 193 if (ls->bounded && *f <= ftest2 && ls->step >= bstepmin2) { 194 PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease.\n")); 195 ls->reason = TAOLINESEARCH_SUCCESS; 196 break; 197 } 198 199 /* Checks for bad cases */ 200 if ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || !mt->infoc) { 201 PetscCall(PetscInfo(ls, "Rounding errors may prevent further progress. May not be a step satisfying\nsufficient decrease and curvature conditions. Tolerances may be too small.\n")); 202 ls->reason = TAOLINESEARCH_HALTED_OTHER; 203 break; 204 } 205 if (ls->step == ls->stepmax && *f <= ftest1 && dg <= dgtest) { 206 PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax)); 207 ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND; 208 break; 209 } 210 if (ls->step == ls->stepmin && *f >= ftest1 && dg >= dgtest) { 211 PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin)); 212 ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND; 213 break; 214 } 215 if (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) { 216 PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol)); 217 ls->reason = TAOLINESEARCH_HALTED_RTOL; 218 break; 219 } 220 221 /* In the first stage, we seek a step for which the modified function 222 has a nonpositive value and nonnegative derivative */ 223 if (stage1 && *f <= ftest1 && dg >= dginit * PetscMin(ls->ftol, ls->gtol)) stage1 = 0; 224 225 /* A modified function is used to predict the step only if we 226 have not obtained a step for which the modified function has a 227 nonpositive function value and nonnegative derivative, and if a 228 lower function value has been obtained but the decrease is not 229 sufficient */ 230 231 if (stage1 && *f <= fx && *f > ftest1) { 232 fm = *f - ls->step * dgtest; /* Define modified function */ 233 fxm = fx - stx * dgtest; /* and derivatives */ 234 fym = fy - sty * dgtest; 235 dgm = dg - dgtest; 236 dgxm = dgx - dgtest; 237 dgym = dgy - dgtest; 238 239 /* if (dgxm * (ls->step - stx) >= 0.0) */ 240 /* Update the interval of uncertainty and compute the new step */ 241 PetscCall(Tao_mcstep(ls, &stx, &fxm, &dgxm, &sty, &fym, &dgym, &ls->step, &fm, &dgm)); 242 243 fx = fxm + stx * dgtest; /* Reset the function and */ 244 fy = fym + sty * dgtest; /* gradient values */ 245 dgx = dgxm + dgtest; 246 dgy = dgym + dgtest; 247 } else { 248 /* Update the interval of uncertainty and compute the new step */ 249 PetscCall(Tao_mcstep(ls, &stx, &fx, &dgx, &sty, &fy, &dgy, &ls->step, f, &dg)); 250 } 251 252 /* Force a sufficient decrease in the interval of uncertainty */ 253 if (mt->bracket) { 254 if (PetscAbsReal(sty - stx) >= 0.66 * width1) ls->step = stx + 0.5 * (sty - stx); 255 width1 = width; 256 width = PetscAbsReal(sty - stx); 257 } 258 } 259 if (ls->nfeval + ls->nfgeval > ls->max_funcs) { 260 PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs)); 261 ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 262 } 263 ls->stepmax = ostepmax; 264 ls->stepmin = ostepmin; 265 266 /* Finish computations */ 267 PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step)); 268 269 /* Set new solution vector and compute gradient if needed */ 270 PetscCall(VecCopy(mt->work, x)); 271 if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, mt->work, g)); 272 PetscFunctionReturn(PETSC_SUCCESS); 273 } 274 275 /*MC 276 TAOLINESEARCHMT - Line-search type with cubic interpolation that satisfies both the sufficient decrease and 277 curvature conditions. This method can take step lengths greater than 1. 278 279 More-Thuente line-search can be selected with "-tao_ls_type more-thuente". 280 281 References: 282 . * - JORGE J. MORE AND DAVID J. THUENTE, LINE SEARCH ALGORITHMS WITH GUARANTEED SUFFICIENT DECREASE. 283 ACM Trans. Math. Software 20, no. 3 (1994): 286-307. 284 285 Level: developer 286 287 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()` 288 289 .keywords: Tao, linesearch 290 M*/ 291 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_MT(TaoLineSearch ls) 292 { 293 TaoLineSearch_MT *ctx; 294 295 PetscFunctionBegin; 296 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 297 PetscCall(PetscNew(&ctx)); 298 ctx->bracket = 0; 299 ctx->infoc = 1; 300 ls->data = (void *)ctx; 301 ls->initstep = 1.0; 302 ls->ops->setup = NULL; 303 ls->ops->reset = NULL; 304 ls->ops->apply = TaoLineSearchApply_MT; 305 ls->ops->destroy = TaoLineSearchDestroy_MT; 306 ls->ops->setfromoptions = TaoLineSearchSetFromOptions_MT; 307 ls->ops->monitor = TaoLineSearchMonitor_MT; 308 PetscFunctionReturn(PETSC_SUCCESS); 309 } 310 311 /* 312 The subroutine mcstep is taken from the work of Jorge Nocedal. 313 this is a variant of More' and Thuente's routine. 314 315 subroutine mcstep 316 317 the purpose of mcstep is to compute a safeguarded step for 318 a linesearch and to update an interval of uncertainty for 319 a minimizer of the function. 320 321 the parameter stx contains the step with the least function 322 value. the parameter stp contains the current step. it is 323 assumed that the derivative at stx is negative in the 324 direction of the step. if bracket is set true then a 325 minimizer has been bracketed in an interval of uncertainty 326 with endpoints stx and sty. 327 328 the subroutine statement is 329 330 subroutine mcstep(stx,fx,dx,sty,fy,dy,stp,fp,dp,bracket, 331 stpmin,stpmax,info) 332 333 where 334 335 stx, fx, and dx are variables which specify the step, 336 the function, and the derivative at the best step obtained 337 so far. The derivative must be negative in the direction 338 of the step, that is, dx and stp-stx must have opposite 339 signs. On output these parameters are updated appropriately. 340 341 sty, fy, and dy are variables which specify the step, 342 the function, and the derivative at the other endpoint of 343 the interval of uncertainty. On output these parameters are 344 updated appropriately. 345 346 stp, fp, and dp are variables which specify the step, 347 the function, and the derivative at the current step. 348 If bracket is set true then on input stp must be 349 between stx and sty. On output stp is set to the new step. 350 351 bracket is a logical variable which specifies if a minimizer 352 has been bracketed. If the minimizer has not been bracketed 353 then on input bracket must be set false. If the minimizer 354 is bracketed then on output bracket is set true. 355 356 stpmin and stpmax are input variables which specify lower 357 and upper bounds for the step. 358 359 info is an integer output variable set as follows: 360 if info = 1,2,3,4,5, then the step has been computed 361 according to one of the five cases below. otherwise 362 info = 0, and this indicates improper input parameters. 363 364 subprograms called 365 366 fortran-supplied ... abs,max,min,sqrt 367 368 argonne national laboratory. minpack project. june 1983 369 jorge j. more', david j. thuente 370 371 */ 372 373 static PetscErrorCode Tao_mcstep(TaoLineSearch ls, PetscReal *stx, PetscReal *fx, PetscReal *dx, PetscReal *sty, PetscReal *fy, PetscReal *dy, PetscReal *stp, PetscReal *fp, PetscReal *dp) 374 { 375 TaoLineSearch_MT *mtP = (TaoLineSearch_MT *)ls->data; 376 PetscReal gamma1, p, q, r, s, sgnd, stpc, stpf, stpq, theta; 377 PetscInt bound; 378 379 PetscFunctionBegin; 380 /* Check the input parameters for errors */ 381 mtP->infoc = 0; 382 PetscCheck(!mtP->bracket || (*stp > PetscMin(*stx, *sty) && *stp < PetscMax(*stx, *sty)), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad stp in bracket"); 383 PetscCheck(*dx * (*stp - *stx) < 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "dx * (stp-stx) >= 0.0"); 384 PetscCheck(ls->stepmax >= ls->stepmin, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "stepmax > stepmin"); 385 386 /* Determine if the derivatives have opposite sign */ 387 sgnd = *dp * (*dx / PetscAbsReal(*dx)); 388 389 if (*fp > *fx) { 390 /* Case 1: a higher function value. 391 The minimum is bracketed. If the cubic step is closer 392 to stx than the quadratic step, the cubic step is taken, 393 else the average of the cubic and quadratic steps is taken. */ 394 395 mtP->infoc = 1; 396 bound = 1; 397 theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp; 398 s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx)); 399 s = PetscMax(s, PetscAbsReal(*dp)); 400 gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s)); 401 if (*stp < *stx) gamma1 = -gamma1; 402 /* Can p be 0? Check */ 403 p = (gamma1 - *dx) + theta; 404 q = ((gamma1 - *dx) + gamma1) + *dp; 405 r = p / q; 406 stpc = *stx + r * (*stp - *stx); 407 stpq = *stx + ((*dx / ((*fx - *fp) / (*stp - *stx) + *dx)) * 0.5) * (*stp - *stx); 408 409 if (PetscAbsReal(stpc - *stx) < PetscAbsReal(stpq - *stx)) stpf = stpc; 410 else stpf = stpc + 0.5 * (stpq - stpc); 411 mtP->bracket = 1; 412 } else if (sgnd < 0.0) { 413 /* Case 2: A lower function value and derivatives of 414 opposite sign. The minimum is bracketed. If the cubic 415 step is closer to stx than the quadratic (secant) step, 416 the cubic step is taken, else the quadratic step is taken. */ 417 418 mtP->infoc = 2; 419 bound = 0; 420 theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp; 421 s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx)); 422 s = PetscMax(s, PetscAbsReal(*dp)); 423 gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s)); 424 if (*stp > *stx) gamma1 = -gamma1; 425 p = (gamma1 - *dp) + theta; 426 q = ((gamma1 - *dp) + gamma1) + *dx; 427 r = p / q; 428 stpc = *stp + r * (*stx - *stp); 429 stpq = *stp + (*dp / (*dp - *dx)) * (*stx - *stp); 430 431 if (PetscAbsReal(stpc - *stp) > PetscAbsReal(stpq - *stp)) stpf = stpc; 432 else stpf = stpq; 433 mtP->bracket = 1; 434 } else if (PetscAbsReal(*dp) < PetscAbsReal(*dx)) { 435 /* Case 3: A lower function value, derivatives of the 436 same sign, and the magnitude of the derivative decreases. 437 The cubic step is only used if the cubic tends to infinity 438 in the direction of the step or if the minimum of the cubic 439 is beyond stp. Otherwise the cubic step is defined to be 440 either stepmin or stepmax. The quadratic (secant) step is also 441 computed and if the minimum is bracketed then the step 442 closest to stx is taken, else the step farthest away is taken. */ 443 444 mtP->infoc = 3; 445 bound = 1; 446 theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp; 447 s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx)); 448 s = PetscMax(s, PetscAbsReal(*dp)); 449 450 /* The case gamma1 = 0 only arises if the cubic does not tend 451 to infinity in the direction of the step. */ 452 gamma1 = s * PetscSqrtScalar(PetscMax(0.0, PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s))); 453 if (*stp > *stx) gamma1 = -gamma1; 454 p = (gamma1 - *dp) + theta; 455 q = (gamma1 + (*dx - *dp)) + gamma1; 456 r = p / q; 457 if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r * (*stx - *stp); 458 else if (*stp > *stx) stpc = ls->stepmax; 459 else stpc = ls->stepmin; 460 stpq = *stp + (*dp / (*dp - *dx)) * (*stx - *stp); 461 462 if (mtP->bracket) { 463 if (PetscAbsReal(*stp - stpc) < PetscAbsReal(*stp - stpq)) stpf = stpc; 464 else stpf = stpq; 465 } else { 466 if (PetscAbsReal(*stp - stpc) > PetscAbsReal(*stp - stpq)) stpf = stpc; 467 else stpf = stpq; 468 } 469 } else { 470 /* Case 4: A lower function value, derivatives of the 471 same sign, and the magnitude of the derivative does 472 not decrease. If the minimum is not bracketed, the step 473 is either stpmin or stpmax, else the cubic step is taken. */ 474 475 mtP->infoc = 4; 476 bound = 0; 477 if (mtP->bracket) { 478 theta = 3 * (*fp - *fy) / (*sty - *stp) + *dy + *dp; 479 s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dy)); 480 s = PetscMax(s, PetscAbsReal(*dp)); 481 gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dy / s) * (*dp / s)); 482 if (*stp > *sty) gamma1 = -gamma1; 483 p = (gamma1 - *dp) + theta; 484 q = ((gamma1 - *dp) + gamma1) + *dy; 485 r = p / q; 486 stpc = *stp + r * (*sty - *stp); 487 stpf = stpc; 488 } else if (*stp > *stx) { 489 stpf = ls->stepmax; 490 } else { 491 stpf = ls->stepmin; 492 } 493 } 494 495 /* Update the interval of uncertainty. This update does not 496 depend on the new step or the case analysis above. */ 497 498 if (*fp > *fx) { 499 *sty = *stp; 500 *fy = *fp; 501 *dy = *dp; 502 } else { 503 if (sgnd < 0.0) { 504 *sty = *stx; 505 *fy = *fx; 506 *dy = *dx; 507 } 508 *stx = *stp; 509 *fx = *fp; 510 *dx = *dp; 511 } 512 513 /* Compute the new step and safeguard it. */ 514 stpf = PetscMin(ls->stepmax, stpf); 515 stpf = PetscMax(ls->stepmin, stpf); 516 *stp = stpf; 517 if (mtP->bracket && bound) { 518 if (*sty > *stx) *stp = PetscMin(*stx + 0.66 * (*sty - *stx), *stp); 519 else *stp = PetscMax(*stx + 0.66 * (*sty - *stx), *stp); 520 } 521 PetscFunctionReturn(PETSC_SUCCESS); 522 } 523