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