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 PetscErrorCode ierr; 15 TaoLineSearch_MT *mt = (TaoLineSearch_MT*)(ls->data); 16 17 PetscFunctionBegin; 18 ierr = PetscObjectDereference((PetscObject)mt->x);CHKERRQ(ierr); 19 ierr = VecDestroy(&mt->work);CHKERRQ(ierr); 20 ierr = PetscFree(ls->data);CHKERRQ(ierr); 21 PetscFunctionReturn(0); 22 } 23 24 static PetscErrorCode TaoLineSearchSetFromOptions_MT(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls) 25 { 26 PetscFunctionBegin; 27 PetscFunctionReturn(0); 28 } 29 30 static PetscErrorCode TaoLineSearchMonitor_MT(TaoLineSearch ls) 31 { 32 TaoLineSearch_MT *mt = (TaoLineSearch_MT*)ls->data; 33 PetscErrorCode ierr; 34 35 PetscFunctionBegin; 36 ierr = PetscViewerASCIIPrintf(ls->viewer, "stx: %g, fx: %g, dgx: %g\n", (double)mt->stx, (double)mt->fx, (double)mt->dgx);CHKERRQ(ierr); 37 ierr = PetscViewerASCIIPrintf(ls->viewer, "sty: %g, fy: %g, dgy: %g\n", (double)mt->sty, (double)mt->fy, (double)mt->dgy);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 static PetscErrorCode TaoLineSearchApply_MT(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 42 { 43 PetscErrorCode ierr; 44 TaoLineSearch_MT *mt = (TaoLineSearch_MT*)(ls->data); 45 PetscReal xtrapf = 4.0; 46 PetscReal finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym; 47 PetscReal dgx, dgy, dg, dg2, fx, fy, stx, sty, dgtest; 48 PetscReal ftest1=0.0, ftest2=0.0; 49 PetscInt i, stage1,n1,n2,nn1,nn2; 50 PetscReal bstepmin1, bstepmin2, bstepmax; 51 PetscBool g_computed = PETSC_FALSE; /* to prevent extra gradient computation */ 52 53 PetscFunctionBegin; 54 ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 55 ierr = TaoLineSearchMonitor(ls, 0, *f, 0.0);CHKERRQ(ierr); 56 /* Check work vector */ 57 if (!mt->work) { 58 ierr = VecDuplicate(x,&mt->work);CHKERRQ(ierr); 59 mt->x = x; 60 ierr = PetscObjectReference((PetscObject)mt->x);CHKERRQ(ierr); 61 } else if (x != mt->x) { 62 ierr = VecDestroy(&mt->work);CHKERRQ(ierr); 63 ierr = VecDuplicate(x,&mt->work);CHKERRQ(ierr); 64 ierr = PetscObjectDereference((PetscObject)mt->x);CHKERRQ(ierr); 65 mt->x = x; 66 ierr = PetscObjectReference((PetscObject)mt->x);CHKERRQ(ierr); 67 } 68 69 if (ls->bounded) { 70 /* Compute step length needed to make all variables equal a bound */ 71 /* Compute the smallest steplength that will make one nonbinding variable 72 equal the bound */ 73 ierr = VecGetLocalSize(ls->upper,&n1);CHKERRQ(ierr); 74 ierr = VecGetLocalSize(mt->x, &n2);CHKERRQ(ierr); 75 ierr = VecGetSize(ls->upper,&nn1);CHKERRQ(ierr); 76 ierr = VecGetSize(mt->x,&nn2);CHKERRQ(ierr); 77 PetscCheck(n1 == n2 && nn1 == nn2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Variable vector not compatible with bounds vector"); 78 ierr = VecScale(s,-1.0);CHKERRQ(ierr); 79 ierr = VecBoundGradientProjection(s,x,ls->lower,ls->upper,s);CHKERRQ(ierr); 80 ierr = VecScale(s,-1.0);CHKERRQ(ierr); 81 ierr = VecStepBoundInfo(x,s,ls->lower,ls->upper,&bstepmin1,&bstepmin2,&bstepmax);CHKERRQ(ierr); 82 ls->stepmax = PetscMin(bstepmax,1.0e15); 83 } 84 85 ierr = VecDot(g,s,&dginit);CHKERRQ(ierr); 86 if (PetscIsInfOrNanReal(dginit)) { 87 ierr = PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)dginit);CHKERRQ(ierr); 88 ls->reason = TAOLINESEARCH_FAILED_INFORNAN; 89 PetscFunctionReturn(0); 90 } 91 if (dginit >= 0.0) { 92 ierr = PetscInfo(ls,"Initial Line Search step * g is not descent direction (%g)\n",(double)dginit);CHKERRQ(ierr); 93 ls->reason = TAOLINESEARCH_FAILED_ASCENT; 94 PetscFunctionReturn(0); 95 } 96 97 /* Initialization */ 98 mt->bracket = 0; 99 stage1 = 1; 100 finit = *f; 101 dgtest = ls->ftol * dginit; 102 width = ls->stepmax - ls->stepmin; 103 width1 = width * 2.0; 104 ierr = VecCopy(x,mt->work);CHKERRQ(ierr); 105 /* Variable dictionary: 106 stx, fx, dgx - the step, function, and derivative at the best step 107 sty, fy, dgy - the step, function, and derivative at the other endpoint 108 of the interval of uncertainty 109 step, f, dg - the step, function, and derivative at the current step */ 110 111 stx = 0.0; 112 fx = finit; 113 dgx = dginit; 114 sty = 0.0; 115 fy = finit; 116 dgy = dginit; 117 118 ls->step = ls->initstep; 119 for (i=0; i< ls->max_funcs; i++) { 120 /* Set min and max steps to correspond to the interval of uncertainty */ 121 if (mt->bracket) { 122 ls->stepmin = PetscMin(stx,sty); 123 ls->stepmax = PetscMax(stx,sty); 124 } else { 125 ls->stepmin = stx; 126 ls->stepmax = ls->step + xtrapf * (ls->step - stx); 127 } 128 129 /* Force the step to be within the bounds */ 130 ls->step = PetscMax(ls->step,ls->stepmin); 131 ls->step = PetscMin(ls->step,ls->stepmax); 132 133 /* If an unusual termination is to occur, then let step be the lowest 134 point obtained thus far */ 135 if ((stx!=0) && (((mt->bracket) && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || ((mt->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) || 136 ((ls->nfeval+ls->nfgeval) >= ls->max_funcs - 1) || (mt->infoc == 0))) { 137 ls->step = stx; 138 } 139 140 ierr = VecCopy(x,mt->work);CHKERRQ(ierr); 141 ierr = VecAXPY(mt->work,ls->step,s);CHKERRQ(ierr); /* W = X + step*S */ 142 143 if (ls->bounded) { 144 ierr = VecMedian(ls->lower, mt->work, ls->upper, mt->work);CHKERRQ(ierr); 145 } 146 if (ls->usegts) { 147 ierr = TaoLineSearchComputeObjectiveAndGTS(ls,mt->work,f,&dg);CHKERRQ(ierr); 148 g_computed = PETSC_FALSE; 149 } else { 150 ierr = TaoLineSearchComputeObjectiveAndGradient(ls,mt->work,f,g);CHKERRQ(ierr); 151 g_computed = PETSC_TRUE; 152 if (ls->bounded) { 153 ierr = VecDot(g,x,&dg);CHKERRQ(ierr); 154 ierr = VecDot(g,mt->work,&dg2);CHKERRQ(ierr); 155 dg = (dg2 - dg)/ls->step; 156 } else { 157 ierr = VecDot(g,s,&dg);CHKERRQ(ierr); 158 } 159 } 160 161 /* update bracketing parameters in the MT context for printouts in monitor */ 162 mt->stx = stx; 163 mt->fx = fx; 164 mt->dgx = dgx; 165 mt->sty = sty; 166 mt->fy = fy; 167 mt->dgy = dgy; 168 ierr = TaoLineSearchMonitor(ls, i+1, *f, ls->step);CHKERRQ(ierr); 169 170 if (i == 0) ls->f_fullstep=*f; 171 172 if (PetscIsInfOrNanReal(*f) || PetscIsInfOrNanReal(dg)) { 173 /* User provided compute function generated Not-a-Number, assume 174 domain violation and set function value and directional 175 derivative to infinity. */ 176 *f = PETSC_INFINITY; 177 dg = PETSC_INFINITY; 178 } 179 180 ftest1 = finit + ls->step * dgtest; 181 if (ls->bounded) ftest2 = finit + ls->step * dgtest * ls->ftol; 182 183 /* Convergence testing */ 184 if (((*f - ftest1 <= 1.0e-10 * PetscAbsReal(finit)) && (PetscAbsReal(dg) + ls->gtol*dginit <= 0.0))) { 185 ierr = PetscInfo(ls, "Line search success: Sufficient decrease and directional deriv conditions hold\n");CHKERRQ(ierr); 186 ls->reason = TAOLINESEARCH_SUCCESS; 187 break; 188 } 189 190 /* Check Armijo if beyond the first breakpoint */ 191 if (ls->bounded && (*f <= ftest2) && (ls->step >= bstepmin2)) { 192 ierr = PetscInfo(ls,"Line search success: Sufficient decrease.\n");CHKERRQ(ierr); 193 ls->reason = TAOLINESEARCH_SUCCESS; 194 break; 195 } 196 197 /* Checks for bad cases */ 198 if (((mt->bracket) && (ls->step <= ls->stepmin||ls->step >= ls->stepmax)) || (!mt->infoc)) { 199 ierr = PetscInfo(ls,"Rounding errors may prevent further progress. May not be a step satisfying\n");CHKERRQ(ierr); 200 ierr = PetscInfo(ls,"sufficient decrease and curvature conditions. Tolerances may be too small.\n");CHKERRQ(ierr); 201 ls->reason = TAOLINESEARCH_HALTED_OTHER; 202 break; 203 } 204 if ((ls->step == ls->stepmax) && (*f <= ftest1) && (dg <= dgtest)) { 205 ierr = PetscInfo(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax);CHKERRQ(ierr); 206 ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND; 207 break; 208 } 209 if ((ls->step == ls->stepmin) && (*f >= ftest1) && (dg >= dgtest)) { 210 ierr = PetscInfo(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin);CHKERRQ(ierr); 211 ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND; 212 break; 213 } 214 if ((mt->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)) { 215 ierr = PetscInfo(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol);CHKERRQ(ierr); 216 ls->reason = TAOLINESEARCH_HALTED_RTOL; 217 break; 218 } 219 220 /* In the first stage, we seek a step for which the modified function 221 has a nonpositive value and nonnegative derivative */ 222 if ((stage1) && (*f <= ftest1) && (dg >= dginit * PetscMin(ls->ftol, ls->gtol))) stage1 = 0; 223 224 /* A modified function is used to predict the step only if we 225 have not obtained a step for which the modified function has a 226 nonpositive function value and nonnegative derivative, and if a 227 lower function value has been obtained but the decrease is not 228 sufficient */ 229 230 if ((stage1) && (*f <= fx) && (*f > ftest1)) { 231 fm = *f - ls->step * dgtest; /* Define modified function */ 232 fxm = fx - stx * dgtest; /* and derivatives */ 233 fym = fy - sty * dgtest; 234 dgm = dg - dgtest; 235 dgxm = dgx - dgtest; 236 dgym = dgy - dgtest; 237 238 /* if (dgxm * (ls->step - stx) >= 0.0) */ 239 /* Update the interval of uncertainty and compute the new step */ 240 ierr = Tao_mcstep(ls,&stx,&fxm,&dgxm,&sty,&fym,&dgym,&ls->step,&fm,&dgm);CHKERRQ(ierr); 241 242 fx = fxm + stx * dgtest; /* Reset the function and */ 243 fy = fym + sty * dgtest; /* gradient values */ 244 dgx = dgxm + dgtest; 245 dgy = dgym + dgtest; 246 } else { 247 /* Update the interval of uncertainty and compute the new step */ 248 ierr = Tao_mcstep(ls,&stx,&fx,&dgx,&sty,&fy,&dgy,&ls->step,f,&dg);CHKERRQ(ierr); 249 } 250 251 /* Force a sufficient decrease in the interval of uncertainty */ 252 if (mt->bracket) { 253 if (PetscAbsReal(sty - stx) >= 0.66 * width1) ls->step = stx + 0.5*(sty - stx); 254 width1 = width; 255 width = PetscAbsReal(sty - stx); 256 } 257 } 258 if ((ls->nfeval+ls->nfgeval) > ls->max_funcs) { 259 ierr = PetscInfo(ls,"Number of line search function evals (%D) > maximum (%D)\n",(ls->nfeval+ls->nfgeval),ls->max_funcs);CHKERRQ(ierr); 260 ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 261 } 262 263 /* Finish computations */ 264 ierr = PetscInfo(ls,"%D function evals in line search, step = %g\n",(ls->nfeval+ls->nfgeval),(double)ls->step);CHKERRQ(ierr); 265 266 /* Set new solution vector and compute gradient if needed */ 267 ierr = VecCopy(mt->work,x);CHKERRQ(ierr); 268 if (!g_computed) { 269 ierr = TaoLineSearchComputeGradient(ls,mt->work,g);CHKERRQ(ierr); 270 } 271 PetscFunctionReturn(0); 272 } 273 274 /*MC 275 TAOLINESEARCHMT - Line-search type with cubic interpolation that satisfies both the sufficient decrease and 276 curvature conditions. This method can take step lengths greater than 1. 277 278 More-Thuente line-search can be selected with "-tao_ls_type more-thuente". 279 280 References: 281 . * - JORGE J. MORE AND DAVID J. THUENTE, LINE SEARCH ALGORITHMS WITH GUARANTEED SUFFICIENT DECREASE. 282 ACM Trans. Math. Software 20, no. 3 (1994): 286-307. 283 284 Level: developer 285 286 .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchApply() 287 288 .keywords: Tao, linesearch 289 M*/ 290 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_MT(TaoLineSearch ls) 291 { 292 PetscErrorCode ierr; 293 TaoLineSearch_MT *ctx; 294 295 PetscFunctionBegin; 296 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 297 ierr = PetscNewLog(ls,&ctx);CHKERRQ(ierr); 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(0); 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)) { 410 stpf = stpc; 411 } else { 412 stpf = stpc + 0.5*(stpq - stpc); 413 } 414 mtP->bracket = 1; 415 } else if (sgnd < 0.0) { 416 /* Case 2: A lower function value and derivatives of 417 opposite sign. The minimum is bracketed. If the cubic 418 step is closer to stx than the quadratic (secant) step, 419 the cubic step is taken, else the quadratic step is taken. */ 420 421 mtP->infoc = 2; 422 bound = 0; 423 theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp; 424 s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx)); 425 s = PetscMax(s,PetscAbsReal(*dp)); 426 gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s)); 427 if (*stp > *stx) gamma1 = -gamma1; 428 p = (gamma1 - *dp) + theta; 429 q = ((gamma1 - *dp) + gamma1) + *dx; 430 r = p/q; 431 stpc = *stp + r*(*stx - *stp); 432 stpq = *stp + (*dp/(*dp-*dx))*(*stx - *stp); 433 434 if (PetscAbsReal(stpc-*stp) > PetscAbsReal(stpq-*stp)) { 435 stpf = stpc; 436 } else { 437 stpf = stpq; 438 } 439 mtP->bracket = 1; 440 } else if (PetscAbsReal(*dp) < PetscAbsReal(*dx)) { 441 /* Case 3: A lower function value, derivatives of the 442 same sign, and the magnitude of the derivative decreases. 443 The cubic step is only used if the cubic tends to infinity 444 in the direction of the step or if the minimum of the cubic 445 is beyond stp. Otherwise the cubic step is defined to be 446 either stepmin or stepmax. The quadratic (secant) step is also 447 computed and if the minimum is bracketed then the step 448 closest to stx is taken, else the step farthest away is taken. */ 449 450 mtP->infoc = 3; 451 bound = 1; 452 theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp; 453 s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx)); 454 s = PetscMax(s,PetscAbsReal(*dp)); 455 456 /* The case gamma1 = 0 only arises if the cubic does not tend 457 to infinity in the direction of the step. */ 458 gamma1 = s*PetscSqrtScalar(PetscMax(0.0,PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s))); 459 if (*stp > *stx) gamma1 = -gamma1; 460 p = (gamma1 - *dp) + theta; 461 q = (gamma1 + (*dx - *dp)) + gamma1; 462 r = p/q; 463 if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r*(*stx - *stp); 464 else if (*stp > *stx) stpc = ls->stepmax; 465 else stpc = ls->stepmin; 466 stpq = *stp + (*dp/(*dp-*dx)) * (*stx - *stp); 467 468 if (mtP->bracket) { 469 if (PetscAbsReal(*stp-stpc) < PetscAbsReal(*stp-stpq)) { 470 stpf = stpc; 471 } else { 472 stpf = stpq; 473 } 474 } else { 475 if (PetscAbsReal(*stp-stpc) > PetscAbsReal(*stp-stpq)) { 476 stpf = stpc; 477 } else { 478 stpf = stpq; 479 } 480 } 481 } else { 482 /* Case 4: A lower function value, derivatives of the 483 same sign, and the magnitude of the derivative does 484 not decrease. If the minimum is not bracketed, the step 485 is either stpmin or stpmax, else the cubic step is taken. */ 486 487 mtP->infoc = 4; 488 bound = 0; 489 if (mtP->bracket) { 490 theta = 3*(*fp - *fy)/(*sty - *stp) + *dy + *dp; 491 s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dy)); 492 s = PetscMax(s,PetscAbsReal(*dp)); 493 gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dy/s)*(*dp/s)); 494 if (*stp > *sty) gamma1 = -gamma1; 495 p = (gamma1 - *dp) + theta; 496 q = ((gamma1 - *dp) + gamma1) + *dy; 497 r = p/q; 498 stpc = *stp + r*(*sty - *stp); 499 stpf = stpc; 500 } else if (*stp > *stx) { 501 stpf = ls->stepmax; 502 } else { 503 stpf = ls->stepmin; 504 } 505 } 506 507 /* Update the interval of uncertainty. This update does not 508 depend on the new step or the case analysis above. */ 509 510 if (*fp > *fx) { 511 *sty = *stp; 512 *fy = *fp; 513 *dy = *dp; 514 } else { 515 if (sgnd < 0.0) { 516 *sty = *stx; 517 *fy = *fx; 518 *dy = *dx; 519 } 520 *stx = *stp; 521 *fx = *fp; 522 *dx = *dp; 523 } 524 525 /* Compute the new step and safeguard it. */ 526 stpf = PetscMin(ls->stepmax,stpf); 527 stpf = PetscMax(ls->stepmin,stpf); 528 *stp = stpf; 529 if (mtP->bracket && bound) { 530 if (*sty > *stx) { 531 *stp = PetscMin(*stx+0.66*(*sty-*stx),*stp); 532 } else { 533 *stp = PetscMax(*stx+0.66*(*sty-*stx),*stp); 534 } 535 } 536 PetscFunctionReturn(0); 537 } 538