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