1 /++ 2 Base numeric algorithms. 3 4 Reworked part of `std.numeric`. 5 6 License: $(HTTP www.apache.org/licenses/LICENSE-2.0, Apache-2.0) 7 Authors: Ilya Yaroshenko (API, findLocalMin, findRoot extension), Don Clugston (findRoot), Lars Tandle Kyllingstad (diff) 8 +/ 9 module mir.numeric; 10 11 import mir.internal.utility: isFloatingPoint; 12 import mir.math.common; 13 import mir.math.ieee; 14 15 version(D_Exceptions) 16 { 17 private static immutable findRoot_badBounds = new Exception("findRoot/findLocalMin: f(ax) and f(bx) must have opposite signs to bracket the root."); 18 private static immutable findRoot_nanX = new Exception("findRoot/findLocalMin: ax or bx is NaN."); 19 private static immutable findRoot_nanY = new Exception("findRoot/findLocalMin: f(x) returned NaN."); 20 } 21 22 /++ 23 +/ 24 enum mir_find_root_status 25 { 26 /// Success 27 success, 28 /// 29 badBounds, 30 /// 31 nanX, 32 /// 33 nanY, 34 } 35 36 /// ditto 37 alias FindRootStatus = mir_find_root_status; 38 39 /++ 40 +/ 41 struct mir_find_root_result(T) 42 { 43 /// Left bound 44 T ax = 0; 45 /// Rifht bound 46 T bx = 0; 47 /// `f(ax)` or `f(ax).fabs.fmin(T.max / 2).copysign(f(ax))`. 48 T ay = 0; 49 /// `f(bx)` or `f(bx).fabs.fmin(T.max / 2).copysign(f(bx))`. 50 T by = 0; 51 /// Amount of target function calls. 52 uint iterations; 53 54 @safe pure @nogc scope const @property: 55 56 /++ 57 Returns: self 58 Required_versions:`D_Exceptions` 59 Throws: `Exception` if $(LREF FindRootResult.status) isn't $(LREF mir_find_root_status.success). 60 +/ 61 version(D_Exceptions) 62 ref validate() inout return 63 { 64 with(FindRootStatus) final switch(status) 65 { 66 case success: return this; 67 case badBounds: throw findRoot_badBounds; 68 case nanX: throw findRoot_nanX; 69 case nanY: throw findRoot_nanY; 70 } 71 } 72 73 extern(C++) nothrow: 74 75 /++ 76 Returns: $(LREF mir_find_root_status) 77 +/ 78 FindRootStatus status() 79 { 80 with(FindRootStatus) return 81 ax != ax || bx != bx ? nanX : 82 ay != ay || by != by ? nanY : 83 ay.signbit == by.signbit && ay != 0 && by != 0 ? badBounds : 84 success; 85 } 86 87 /++ 88 A bound that corresponds to the minimal absolute function value. 89 90 Returns: `!(ay.fabs > by.fabs) ? ax : bx` 91 +/ 92 T x() 93 { 94 return !(ay.fabs > by.fabs) ? ax : bx; 95 } 96 97 /++ 98 The minimal of absolute function values. 99 100 Returns: `!(ay.fabs > by.fabs) ? ay : by` 101 +/ 102 T y() 103 { 104 return !(ay.fabs > by.fabs) ? ay : by; 105 } 106 } 107 108 /// ditto 109 alias FindRootResult = mir_find_root_result; 110 111 /++ 112 Find root of a real function f(x) by bracketing, allowing the 113 termination condition to be specified. 114 115 Given a function `f` and a range `[a .. b]` such that `f(a)` 116 and `f(b)` have opposite signs or at least one of them equals ±0, 117 returns the value of `x` in 118 the range which is closest to a root of `f(x)`. If `f(x)` 119 has more than one root in the range, one will be chosen 120 arbitrarily. If `f(x)` returns NaN, NaN will be returned; 121 otherwise, this algorithm is guaranteed to succeed. 122 123 Uses an algorithm based on TOMS748, which uses inverse cubic 124 interpolation whenever possible, otherwise reverting to parabolic 125 or secant interpolation. Compared to TOMS748, this implementation 126 improves worst-case performance by a factor of more than 100, and 127 typical performance by a factor of 2. For 80-bit reals, most 128 problems require 8 to 15 calls to `f(x)` to achieve full machine 129 precision. The worst-case performance (pathological cases) is 130 approximately twice the number of bits. 131 132 References: "On Enclosing Simple Roots of Nonlinear Equations", 133 G. Alefeld, F.A. Potra, Yixun Shi, Mathematics of Computation 61, 134 pp733-744 (1993). Fortran code available from $(HTTP 135 www.netlib.org,www.netlib.org) as algorithm TOMS478. 136 137 Params: 138 f = Function to be analyzed. `f(ax)` and `f(bx)` should have opposite signs, or `lowerBound` and/or `upperBound` 139 should be defined to perform initial interval extension. 140 tolerance = Defines an early termination condition. Receives the 141 current upper and lower bounds on the root. The 142 delegate must return `true` when these bounds are 143 acceptable. If this function always returns `false` or 144 it is null, full machine precision will be achieved. 145 ax = Left inner bound of initial range of `f` known to contain the root. 146 bx = Right inner bound of initial range of `f` known to contain the root. Can be equal to `ax`. 147 fax = Value of `f(ax)` (optional). 148 fbx = Value of `f(bx)` (optional). 149 lowerBound = Lower outer bound for interval extension (optional) 150 upperBound = Upper outer bound for interval extension (optional) 151 maxIterations = Appr. maximum allowed number of function calls (inluciding function calls for grid). 152 steps = Number of nodes in the left side and right side regular grids (optional). If it equals to `0` (default), 153 then the algorithm uses `ieeeMean` for searching. The algorithm performs searching if 154 `sgn(fax)` equals to `sgn(fbx)` and at least one outer bound allows to extend the inner initial range. 155 156 Returns: $(LREF FindRootResult) 157 +/ 158 @fmamath 159 FindRootResult!T findRoot(alias f, alias tolerance = null, T)( 160 const T ax, 161 const T bx, 162 const T fax = T.nan, 163 const T fbx = T.nan, 164 const T lowerBound = T.nan, 165 const T upperBound = T.nan, 166 uint maxIterations = T.sizeof * 16, 167 uint steps = 0, 168 ) 169 if ( 170 isFloatingPoint!T && __traits(compiles, T(f(T.init))) && 171 ( 172 is(typeof(tolerance) == typeof(null)) || 173 __traits(compiles, {auto _ = bool(tolerance(T.init, T.init));} 174 ) 175 )) 176 { 177 if (false) // break attributes 178 T y = f(T(1)); 179 scope funInst = delegate(T x) { 180 return T(f(x)); 181 }; 182 scope fun = funInst.trustedAllAttr; 183 184 static if (is(typeof(tolerance) == typeof(null))) 185 { 186 alias tol = tolerance; 187 } 188 else 189 { 190 if (false) // break attributes 191 bool b = tolerance(T(1), T(1)); 192 scope tolInst = delegate(T a, T b) { return bool(tolerance(a, b)); }; 193 scope tol = tolInst.trustedAllAttr; 194 } 195 return findRootImpl(ax, bx, fax, fbx, lowerBound, upperBound, maxIterations, steps, fun, tol); 196 } 197 198 /// 199 // @nogc 200 version(mir_test) @safe unittest 201 { 202 import mir.math.common: log, exp, fabs; 203 204 auto logRoot = findRoot!log(0, double.infinity).validate.x; 205 assert(logRoot == 1); 206 207 auto shift = 1; 208 auto expm1Root = findRoot!(x => exp(x) - shift) 209 (-double.infinity, double.infinity).validate.x; 210 assert(expm1Root == 0); 211 212 auto approxLogRoot = findRoot!(log, (a, b) => fabs(a - b) < 1e-5)(0, double.infinity).validate.x; 213 assert(fabs(approxLogRoot - 1) < 1e-5); 214 } 215 216 /// With adaptive bounds 217 version(mir_test) @safe unittest 218 { 219 import mir.math.common: log, exp, fabs; 220 221 auto logRoot = findRoot!log( 222 10, 10, // assume we have one initial point 223 double.nan, double.nan, // fa, fb aren't provided by user 224 -double.infinity, double.infinity, // all space is available for the bounds extension. 225 ).validate.x; 226 assert(logRoot == 1); 227 228 auto shift = 1; 229 alias expm1Fun = (double x) => exp(x) - shift; 230 auto expm1RootRet = findRoot!expm1Fun 231 ( 232 11, 10, // reversed order for interval is always OK 233 expm1Fun(11), expm1Fun(10), // the order must be the same as bounds 234 0, double.infinity, // space for the bounds extension. 235 ).validate; 236 assert(expm1Fun(expm1RootRet.x) == 0); 237 238 auto approxLogRoot = findRoot!(log, (a, b) => fabs(a - b) < 1e-5)( 239 -1e10, +1e10, 240 double.nan, double.nan, 241 0, double.infinity, 242 ).validate.x; 243 assert(fabs(approxLogRoot - 1) < 1e-5); 244 } 245 246 /// ditto 247 unittest 248 { 249 import core.stdc.tgmath: atan; 250 import mir.math; 251 import std.meta: AliasSeq; 252 253 const double[2][3] boundaries = [ 254 [0.4, 0.6], 255 [1.4, 1.6], 256 [0.1, 2.1]]; 257 258 enum root = 1.0; 259 260 foreach(fun; AliasSeq!( 261 (double x) => x ^^ 2 - root, 262 (double x) => root - x ^^ 2, 263 (double x) => atan(x - root), 264 )) 265 { 266 foreach(ref bounds; boundaries) 267 { 268 auto result = findRoot!fun( 269 bounds[0], bounds[1], 270 double.nan, double.nan, // f(a) and f(b) not provided 271 -double.max, double.max, // user provided outer bounds 272 ); 273 assert(result.validate.x == root); 274 } 275 } 276 277 foreach(ref bounds; boundaries) 278 { 279 auto result = findRoot!(x => sin(x - root))( 280 bounds[0], bounds[1], 281 double.nan, double.nan, // f(a) and f(b) not provided 282 -10, 10, // user provided outer bounds 283 100, // max iterations, 284 10, // steps for grid 285 ); 286 assert(result.validate.x == root); 287 } 288 // single initial point, grid, positive direction 289 auto result = findRoot!((double x ) => sin(x - root))( 290 0.1, 0.1, // initial point, a = b 291 double.nan, double.nan, // f(a) = f(b) not provided 292 -100.0, 100.0, // user provided outer bounds 293 150, // max iterations, 294 100, // number of steps for grid 295 ); 296 assert(result.validate.x == root); 297 298 // single initial point, grid, negative direction 299 result = findRoot!((double x ) => sin(x - root))( 300 0.1, 0.1, // initial point a = b 301 double.nan, double.nan, // f(a) = f(b) not provided 302 100.0, -100.0, // user provided outer bounds, reversed order 303 150, // max iterations, 304 100, // number of steps for grid 305 ); 306 assert(result.validate.x == double(root - PI)); // Left side root! 307 } 308 309 /++ 310 With adaptive bounds and single initial point. 311 Reverse outer bound order controls first step direction 312 in case of `f(a) == f(b)`. 313 +/ 314 unittest 315 { 316 enum root = 1.0; 317 318 // roots are +/- `root` 319 alias fun = (double x) => x * x - root; 320 321 double lowerBound = -10.0; 322 double upperBound = 10.0; 323 324 assert( 325 findRoot!fun( 326 0, 0, // initial interval 327 double.nan, double.nan, 328 lowerBound, upperBound, 329 // positive direction has higher priority 330 ).validate.x == root 331 ); 332 333 assert( 334 findRoot!fun( 335 0, 0, // initial interval 336 double.nan, double.nan, 337 upperBound, lowerBound, 338 // reversed order 339 ).validate.x == -root // other root 340 ); 341 } 342 343 /// $(LREF findRoot) implementations. 344 export @fmamath FindRootResult!float findRootImpl( 345 float ax, 346 float bx, 347 float fax, 348 float fbx, 349 float lowerBound, 350 float upperBound, 351 uint maxIterations, 352 uint steps, 353 scope float delegate(float) @safe pure nothrow @nogc f, 354 scope bool delegate(float, float) @safe pure nothrow @nogc tolerance, //can be null 355 ) @safe pure nothrow @nogc 356 { 357 pragma(inline, false); 358 return findRootImplGen!float(ax, bx, fax, fbx, lowerBound, upperBound, maxIterations, steps, f, tolerance); 359 } 360 361 /// ditto 362 export @fmamath FindRootResult!double findRootImpl( 363 double ax, 364 double bx, 365 double fax, 366 double fbx, 367 double lowerBound, 368 double upperBound, 369 uint maxIterations, 370 uint steps, 371 scope double delegate(double) @safe pure nothrow @nogc f, 372 scope bool delegate(double, double) @safe pure nothrow @nogc tolerance, //can be null 373 ) @safe pure nothrow @nogc 374 { 375 pragma(inline, false); 376 return findRootImplGen!double(ax, bx, fax, fbx, lowerBound, upperBound, maxIterations, steps, f, tolerance); 377 } 378 379 /// ditto 380 export @fmamath FindRootResult!real findRootImpl( 381 real ax, 382 real bx, 383 real fax, 384 real fbx, 385 real lowerBound, 386 real upperBound, 387 uint maxIterations, 388 uint steps, 389 scope real delegate(real) @safe pure nothrow @nogc f, 390 scope bool delegate(real, real) @safe pure nothrow @nogc tolerance, //can be null 391 ) @safe pure nothrow @nogc 392 { 393 pragma(inline, false); 394 return findRootImplGen!real(ax, bx, fax, fbx, lowerBound, upperBound, maxIterations, steps, f, tolerance); 395 } 396 397 private @fmamath FindRootResult!T findRootImplGen(T)( 398 T a, 399 T b, 400 T fa, 401 T fb, 402 T lb, 403 T ub, 404 uint maxIterations, 405 uint steps, 406 scope const T delegate(T) @safe pure nothrow @nogc f, 407 scope const bool delegate(T, T) @safe pure nothrow @nogc tolerance, //can be null 408 ) @safe pure nothrow @nogc 409 if (isFloatingPoint!T) 410 { 411 version(LDC) pragma(inline, true); 412 // Author: Don Clugston. This code is (heavily) modified from TOMS748 413 // (www.netlib.org). The changes to improve the worst-cast performance are 414 // entirely original. 415 416 // Author: Ilya Yaroshenko (Bounds extension logic, 417 // API improvements, infinity and huge numbers handing, compiled code size reduction) 418 419 T d; // [a .. b] is our current bracket. d is the third best guess. 420 T fd; // Value of f at d. 421 bool done = false; // Has a root been found? 422 uint iterations; 423 424 static void swap(ref T a, ref T b) 425 { 426 T t = a; a = b; b = t; 427 } 428 429 bool exit() 430 { 431 pragma(inline, false); 432 return done 433 || iterations >= maxIterations 434 || b == nextUp(a) 435 || tolerance !is null && tolerance(a, b); 436 } 437 438 // Test the function at point c; update brackets accordingly 439 void bracket(T c) 440 { 441 pragma(inline, false); 442 T fc = f(c); 443 iterations++; 444 if (fc == 0 || fc != fc) // Exact solution, or NaN 445 { 446 a = c; 447 fa = fc; 448 d = c; 449 fd = fc; 450 done = true; 451 return; 452 } 453 454 fc = fc.fabs.fmin(T.max / 2).copysign(fc); 455 456 // Determine new enclosing interval 457 if (signbit(fa) != signbit(fc)) 458 { 459 d = b; 460 fd = fb; 461 b = c; 462 fb = fc; 463 } 464 else 465 { 466 d = a; 467 fd = fa; 468 a = c; 469 fa = fc; 470 } 471 } 472 473 /* Perform a secant interpolation. If the result would lie on a or b, or if 474 a and b differ so wildly in magnitude that the result would be meaningless, 475 perform a bisection instead. 476 */ 477 static T secantInterpolate(T a, T b, T fa, T fb) 478 { 479 pragma(inline, false); 480 if (a - b == a && b != 0 481 || b - a == b && a != 0) 482 { 483 // Catastrophic cancellation 484 return ieeeMean(a, b); 485 } 486 // avoid overflow 487 T m = fa - fb; 488 T wa = fa / m; 489 T wb = fb / m; 490 T c = b * wa - a * wb; 491 if (c == a || c == b || c != c || c.fabs == T.infinity) 492 c = a.half + b.half; 493 return c; 494 } 495 496 /* Uses 'numsteps' newton steps to approximate the zero in [a .. b] of the 497 quadratic polynomial interpolating f(x) at a, b, and d. 498 Returns: 499 The approximate zero in [a .. b] of the quadratic polynomial. 500 */ 501 T newtonQuadratic(int numsteps) 502 { 503 // Find the coefficients of the quadratic polynomial. 504 const T a0 = fa; 505 const T a1 = (fb - fa) / (b - a); 506 const T a2 = ((fd - fb) / (d - b) - a1) / (d - a); 507 508 // Determine the starting point of newton steps. 509 T c = a2.signbit != fa.signbit ? a : b; 510 511 // start the safeguarded newton steps. 512 foreach (int i; 0 .. numsteps) 513 { 514 const T pc = a0 + (a1 + a2 * (c - b))*(c - a); 515 const T pdc = a1 + a2*((2 * c) - (a + b)); 516 if (pdc == 0) 517 return a - a0 / a1; 518 else 519 c = c - pc / pdc; 520 } 521 return c; 522 } 523 524 // Starting with the second iteration, higher-order interpolation can 525 // be used. 526 int itnum = 1; // Iteration number 527 int baditer = 1; // Num bisections to take if an iteration is bad. 528 T c, e; // e is our fourth best guess 529 T fe; 530 bool left; 531 532 // Allow a and b to be provided in reverse order 533 if (a > b) 534 { 535 swap(a, b); 536 swap(fa, fb); 537 } 538 539 if (a != a || b != b) 540 { 541 done = true; 542 goto whileloop; 543 } 544 545 if (lb != lb) 546 { 547 lb = a; 548 } 549 550 if (ub != ub) 551 { 552 ub = b; 553 } 554 555 if (lb > ub) 556 { 557 swap(lb, ub); 558 left = true; 559 } 560 561 if (lb == -T.infinity) 562 { 563 lb = -T.max; 564 } 565 566 if (ub == +T.infinity) 567 { 568 ub = +T.max; 569 } 570 571 if (!(lb <= a)) 572 { 573 a = lb; 574 fa = T.nan; 575 } 576 577 if (!(b <= ub)) 578 { 579 a = lb; 580 fa = T.nan; 581 } 582 583 if (fa != fa) 584 { 585 fa = f(a); 586 iterations++; 587 } 588 589 // On the first iteration we take a secant step: 590 if (fa == 0 || fa != fa) 591 { 592 done = true; 593 b = a; 594 fb = fa; 595 goto whileloop; 596 } 597 598 if (fb != fb) 599 { 600 if (a == b) 601 { 602 fb = fa; 603 } 604 else 605 { 606 fb = f(b); 607 iterations++; 608 } 609 } 610 611 if (fb == 0 || fb != fb) 612 { 613 done = true; 614 a = b; 615 fa = fb; 616 goto whileloop; 617 } 618 619 if (fa.fabs < fb.fabs) 620 { 621 left = true; 622 } 623 else 624 if (fa.fabs > fb.fabs) 625 { 626 left = false; 627 } 628 629 // extend inner boundaries 630 if (fa.signbit == fb.signbit) 631 { 632 T lx = a; 633 T ux = b; 634 T ly = fa; 635 T uy = fb; 636 const sb = fa.signbit; 637 638 import mir.ndslice.topology: linspace; 639 640 typeof(linspace!T([2], [lx, lb])) lgrid, ugrid; 641 if (steps) 642 { 643 lgrid = linspace!T([steps + 1], [lx, lb]); 644 lgrid.popFront; 645 ugrid = linspace!T([steps + 1], [ux, ub]); 646 ugrid.popFront; 647 } 648 649 for(T mx;;) 650 { 651 bool lw = lb < lx; 652 bool uw = ux < ub; 653 654 if (!lw && !uw || iterations >= maxIterations) 655 { 656 done = true; 657 goto whileloop; 658 } 659 660 if (lw && (!uw || left)) 661 { 662 if (lgrid.empty) 663 { 664 mx = ieeeMean(lb, lx); 665 if (mx == lx) 666 mx = lb; 667 } 668 else 669 { 670 mx = lgrid.front; 671 lgrid.popFront; 672 if (mx == lx) 673 continue; 674 } 675 T my = f(mx); 676 iterations++; 677 if (my == 0) 678 { 679 a = b = mx; 680 fa = fb = my; 681 done = true; 682 goto whileloop; 683 } 684 if (mx != mx) 685 { 686 lb = mx; 687 continue; 688 } 689 if (my.signbit == sb) 690 { 691 lx = mx; 692 ly = my; 693 if (lgrid.empty) 694 { 695 left = ly.fabs < uy.fabs; 696 } 697 continue; 698 } 699 a = mx; 700 fa = my; 701 b = lx; 702 fb = ly; 703 break; 704 } 705 else 706 { 707 if (ugrid.empty) 708 { 709 mx = ieeeMean(ub, ux); 710 if (mx == ux) 711 mx = ub; 712 } 713 else 714 { 715 mx = ugrid.front; 716 ugrid.popFront; 717 if (mx == ux) 718 continue; 719 } 720 T my = f(mx); 721 iterations++; 722 if (my == 0) 723 { 724 a = b = mx; 725 fa = fb = my; 726 done = true; 727 goto whileloop; 728 } 729 if (mx != mx) 730 { 731 ub = mx; 732 continue; 733 } 734 if (my.signbit == sb) 735 { 736 ux = mx; 737 uy = my; 738 if (ugrid.empty) 739 { 740 left = !(ly.fabs > uy.fabs); 741 } 742 continue; 743 } 744 b = mx; 745 fb = my; 746 a = ux; 747 fa = uy; 748 break; 749 } 750 } 751 } 752 753 fa = fa.fabs.fmin(T.max / 2).copysign(fa); 754 fb = fb.fabs.fmin(T.max / 2).copysign(fb); 755 bracket(secantInterpolate(a, b, fa, fb)); 756 757 whileloop: 758 while (!exit) 759 { 760 T a0 = a, b0 = b; // record the brackets 761 762 // Do two higher-order (cubic or parabolic) interpolation steps. 763 foreach (int QQ; 0 .. 2) 764 { 765 // Cubic inverse interpolation requires that 766 // all four function values fa, fb, fd, and fe are distinct; 767 // otherwise use quadratic interpolation. 768 bool distinct = (fa != fb) && (fa != fd) && (fa != fe) 769 && (fb != fd) && (fb != fe) && (fd != fe); 770 // The first time, cubic interpolation is impossible. 771 if (itnum<2) distinct = false; 772 bool ok = distinct; 773 if (distinct) 774 { 775 // Cubic inverse interpolation of f(x) at a, b, d, and e 776 const q11 = (d - e) * fd / (fe - fd); 777 const q21 = (b - d) * fb / (fd - fb); 778 const q31 = (a - b) * fa / (fb - fa); 779 const d21 = (b - d) * fd / (fd - fb); 780 const d31 = (a - b) * fb / (fb - fa); 781 782 const q22 = (d21 - q11) * fb / (fe - fb); 783 const q32 = (d31 - q21) * fa / (fd - fa); 784 const d32 = (d31 - q21) * fd / (fd - fa); 785 const q33 = (d32 - q22) * fa / (fe - fa); 786 c = a + (q31 + q32 + q33); 787 if (c != c || (c <= a) || (c >= b)) 788 { 789 // DAC: If the interpolation predicts a or b, it's 790 // probable that it's the actual root. Only allow this if 791 // we're already close to the root. 792 if (c == a && (a - b != a || a - b != -b)) 793 { 794 auto down = !(a - b != a); 795 if (down) 796 c = -c; 797 c = c.nextUp; 798 if (down) 799 c = -c; 800 } 801 else 802 { 803 ok = false; 804 } 805 806 } 807 } 808 if (!ok) 809 { 810 // DAC: Alefeld doesn't explain why the number of newton steps 811 // should vary. 812 c = newtonQuadratic(2 + distinct); 813 if (c != c || (c <= a) || (c >= b)) 814 { 815 // Failure, try a secant step: 816 c = secantInterpolate(a, b, fa, fb); 817 } 818 } 819 ++itnum; 820 e = d; 821 fe = fd; 822 bracket(c); 823 if (exit) 824 break whileloop; 825 if (itnum == 2) 826 continue whileloop; 827 } 828 829 // Now we take a double-length secant step: 830 T u; 831 T fu; 832 if (fabs(fa) < fabs(fb)) 833 { 834 u = a; 835 fu = fa; 836 } 837 else 838 { 839 u = b; 840 fu = fb; 841 } 842 c = u - 2 * (fu / (fb - fa)) * (b - a); 843 844 // DAC: If the secant predicts a value equal to an endpoint, it's 845 // probably false. 846 if (c == a || c == b || c != c || fabs(c - u) > (b - a) * 0.5f) 847 { 848 if ((a - b) == a || (b - a) == b) 849 { 850 c = ieeeMean(a, b); 851 } 852 else 853 { 854 c = a.half + b.half; 855 } 856 } 857 e = d; 858 fe = fd; 859 bracket(c); 860 if (exit) 861 break; 862 863 // IMPROVE THE WORST-CASE PERFORMANCE 864 // We must ensure that the bounds reduce by a factor of 2 865 // in binary space! every iteration. If we haven't achieved this 866 // yet, or if we don't yet know what the exponent is, 867 // perform a binary chop. 868 869 if ((a == 0 870 || b == 0 871 || fabs(a) >= 0.5f * fabs(b) 872 && fabs(b) >= 0.5f * fabs(a)) 873 && b - a < 0.25f * (b0 - a0)) 874 { 875 baditer = 1; 876 continue; 877 } 878 879 // DAC: If this happens on consecutive iterations, we probably have a 880 // pathological function. Perform a number of bisections equal to the 881 // total number of consecutive bad iterations. 882 883 if (b - a < 0.25f * (b0 - a0)) 884 baditer = 1; 885 foreach (int QQ; 0 .. baditer) 886 { 887 e = d; 888 fe = fd; 889 bracket(ieeeMean(a, b)); 890 } 891 ++baditer; 892 } 893 return typeof(return)(a, b, fa, fb, iterations); 894 } 895 896 version(mir_test) @safe unittest 897 { 898 import mir.math.constant; 899 900 int numProblems = 0; 901 int numCalls; 902 903 void testFindRoot(real delegate(real) @nogc @safe nothrow pure f , real x1, real x2, int line = __LINE__) //@nogc @safe nothrow pure 904 { 905 //numCalls=0; 906 //++numProblems; 907 assert(x1 == x1 && x2 == x2); 908 auto result = findRoot!f(x1, x2).validate; 909 910 auto flo = f(result.ax); 911 auto fhi = f(result.bx); 912 if (flo != 0) 913 { 914 assert(flo.signbit != fhi.signbit); 915 } 916 } 917 918 // Test functions 919 real cubicfn(real x) @nogc @safe nothrow pure 920 { 921 //++numCalls; 922 if (x>float.max) 923 x = float.max; 924 if (x<-float.max) 925 x = -float.max; 926 // This has a single real root at -59.286543284815 927 return 0.386*x*x*x + 23*x*x + 15.7*x + 525.2; 928 } 929 // Test a function with more than one root. 930 real multisine(real x) { ++numCalls; return sin(x); } 931 testFindRoot( &multisine, 6, 90); 932 testFindRoot(&cubicfn, -100, 100); 933 testFindRoot( &cubicfn, -double.max, real.max); 934 935 936 /* Tests from the paper: 937 * "On Enclosing Simple Roots of Nonlinear Equations", G. Alefeld, F.A. Potra, 938 * Yixun Shi, Mathematics of Computation 61, pp733-744 (1993). 939 */ 940 // Parameters common to many alefeld tests. 941 int n; 942 real ale_a, ale_b; 943 944 int powercalls = 0; 945 946 real power(real x) 947 { 948 ++powercalls; 949 ++numCalls; 950 return pow(x, n) + double.min_normal; 951 } 952 int [] power_nvals = [3, 5, 7, 9, 19, 25]; 953 // Alefeld paper states that pow(x, n) is a very poor case, where bisection 954 // outperforms his method, and gives total numcalls = 955 // 921 for bisection (2.4 calls per bit), 1830 for Alefeld (4.76/bit), 956 // 0.5f624 for brent (6.8/bit) 957 // ... but that is for double, not real80. 958 // This poor performance seems mainly due to catastrophic cancellation, 959 // which is avoided here by the use of ieeeMean(). 960 // I get: 231 (0.48/bit). 961 // IE this is 10X faster in Alefeld's worst case 962 numProblems=0; 963 foreach (k; power_nvals) 964 { 965 n = k; 966 //testFindRoot(&power, -1, 10); 967 } 968 969 int powerProblems = numProblems; 970 971 // Tests from Alefeld paper 972 973 int [9] alefeldSums; 974 real alefeld0(real x) 975 { 976 ++alefeldSums[0]; 977 ++numCalls; 978 real q = sin(x) - x/2; 979 for (int i=1; i<20; ++i) 980 q+=(2*i-5.0)*(2*i-5.0) / ((x-i*i)*(x-i*i)*(x-i*i)); 981 return q; 982 } 983 real alefeld1(real x) 984 { 985 ++numCalls; 986 ++alefeldSums[1]; 987 return ale_a*x + exp(ale_b * x); 988 } 989 real alefeld2(real x) 990 { 991 ++numCalls; 992 ++alefeldSums[2]; 993 return pow(x, n) - ale_a; 994 } 995 real alefeld3(real x) 996 { 997 ++numCalls; 998 ++alefeldSums[3]; 999 return (1.0 +pow(1.0L-n, 2))*x - pow(1.0L-n*x, 2); 1000 } 1001 real alefeld4(real x) 1002 { 1003 ++numCalls; 1004 ++alefeldSums[4]; 1005 return x*x - pow(1-x, n); 1006 } 1007 real alefeld5(real x) 1008 { 1009 ++numCalls; 1010 ++alefeldSums[5]; 1011 return (1+pow(1.0L-n, 4))*x - pow(1.0L-n*x, 4); 1012 } 1013 real alefeld6(real x) 1014 { 1015 ++numCalls; 1016 ++alefeldSums[6]; 1017 return exp(-n*x)*(x-1.01L) + pow(x, n); 1018 } 1019 real alefeld7(real x) 1020 { 1021 ++numCalls; 1022 ++alefeldSums[7]; 1023 return (n*x-1) / ((n-1)*x); 1024 } 1025 1026 numProblems=0; 1027 testFindRoot(&alefeld0, PI_2, PI); 1028 for (n=1; n <= 10; ++n) 1029 { 1030 testFindRoot(&alefeld0, n*n+1e-9L, (n+1)*(n+1)-1e-9L); 1031 } 1032 ale_a = -40; ale_b = -1; 1033 testFindRoot(&alefeld1, -9, 31); 1034 ale_a = -100; ale_b = -2; 1035 testFindRoot(&alefeld1, -9, 31); 1036 ale_a = -200; ale_b = -3; 1037 testFindRoot(&alefeld1, -9, 31); 1038 int [] nvals_3 = [1, 2, 5, 10, 15, 20]; 1039 int [] nvals_5 = [1, 2, 4, 5, 8, 15, 20]; 1040 int [] nvals_6 = [1, 5, 10, 15, 20]; 1041 int [] nvals_7 = [2, 5, 15, 20]; 1042 1043 for (int i=4; i<12; i+=2) 1044 { 1045 n = i; 1046 ale_a = 0.2; 1047 testFindRoot(&alefeld2, 0, 5); 1048 ale_a=1; 1049 testFindRoot(&alefeld2, 0.95, 4.05); 1050 testFindRoot(&alefeld2, 0, 1.5); 1051 } 1052 foreach (i; nvals_3) 1053 { 1054 n=i; 1055 testFindRoot(&alefeld3, 0, 1); 1056 } 1057 foreach (i; nvals_3) 1058 { 1059 n=i; 1060 testFindRoot(&alefeld4, 0, 1); 1061 } 1062 foreach (i; nvals_5) 1063 { 1064 n=i; 1065 testFindRoot(&alefeld5, 0, 1); 1066 } 1067 foreach (i; nvals_6) 1068 { 1069 n=i; 1070 testFindRoot(&alefeld6, 0, 1); 1071 } 1072 foreach (i; nvals_7) 1073 { 1074 n=i; 1075 testFindRoot(&alefeld7, 0.01L, 1); 1076 } 1077 real worstcase(real x) 1078 { 1079 ++numCalls; 1080 return x<0.3*real.max? -0.999e-3: 1.0; 1081 } 1082 testFindRoot(&worstcase, -real.max, real.max); 1083 1084 // just check that the double + float cases compile 1085 findRoot!(x => 0)(-double.max, double.max); 1086 findRoot!(x => -0.0)(-float.max, float.max); 1087 /* 1088 int grandtotal=0; 1089 foreach (calls; alefeldSums) 1090 { 1091 grandtotal+=calls; 1092 } 1093 grandtotal-=2*numProblems; 1094 printf("\nALEFELD TOTAL = %d avg = %f (alefeld avg=19.3 for double)\n", 1095 grandtotal, (1.0*grandtotal)/numProblems); 1096 powercalls -= 2*powerProblems; 1097 printf("POWER TOTAL = %d avg = %f ", powercalls, 1098 (1.0*powercalls)/powerProblems); 1099 */ 1100 //Issue 14231 1101 auto xp = findRoot!(x => x)(0f, 1f); 1102 auto xn = findRoot!(x => x)(-1f, -0f); 1103 } 1104 1105 /++ 1106 +/ 1107 struct FindLocalMinResult(T) 1108 { 1109 /// 1110 T x = 0; 1111 /// 1112 T y = 0; 1113 /// 1114 T error = 0; 1115 1116 @safe pure @nogc scope const @property: 1117 1118 /++ 1119 Returns: self 1120 Required_versions:`D_Exceptions` 1121 Throws: `Exception` if $(LREF FindRootResult.status) isn't $(LREF mir_find_root_status.success). 1122 +/ 1123 version(D_Exceptions) 1124 ref validate() return 1125 { 1126 with(FindRootStatus) final switch(status) 1127 { 1128 case success: return this; 1129 case badBounds: throw findRoot_badBounds; 1130 case nanX: throw findRoot_nanX; 1131 case nanY: throw findRoot_nanY; 1132 } 1133 } 1134 1135 extern(C++) nothrow: 1136 1137 /++ 1138 Returns: $(LREF mir_find_root_status) 1139 +/ 1140 FindRootStatus status() 1141 { 1142 with(FindRootStatus) return 1143 x != x ? nanX : 1144 y != y ? nanY : 1145 success; 1146 } 1147 } 1148 1149 /++ 1150 Find a real minimum of a real function `f(x)` via bracketing. 1151 Given a function `f` and a range `(ax .. bx)`, 1152 returns the value of `x` in the range which is closest to a minimum of `f(x)`. 1153 `f` is never evaluted at the endpoints of `ax` and `bx`. 1154 If `f(x)` has more than one minimum in the range, one will be chosen arbitrarily. 1155 If `f(x)` returns NaN or -Infinity, `(x, f(x), NaN)` will be returned; 1156 otherwise, this algorithm is guaranteed to succeed. 1157 Params: 1158 f = Function to be analyzed 1159 ax = Left bound of initial range of f known to contain the minimum. 1160 bx = Right bound of initial range of f known to contain the minimum. 1161 relTolerance = Relative tolerance. 1162 absTolerance = Absolute tolerance. 1163 Preconditions: 1164 `ax` and `bx` shall be finite reals. $(BR) 1165 `relTolerance` shall be normal positive real. $(BR) 1166 `absTolerance` shall be normal positive real no less then `T.epsilon*2`. 1167 Returns: 1168 A $(LREF FindLocalMinResult) consisting of `x`, `y = f(x)` and `error = 3 * (absTolerance * fabs(x) + relTolerance)`. 1169 The method used is a combination of golden section search and 1170 successive parabolic interpolation. Convergence is never much slower 1171 than that for a Fibonacci search. 1172 References: 1173 "Algorithms for Minimization without Derivatives", Richard Brent, Prentice-Hall, Inc. (1973) 1174 See_Also: $(LREF findRoot), $(REF isNormal, std,math) 1175 +/ 1176 FindLocalMinResult!T findLocalMin(alias f, T)( 1177 const T ax, 1178 const T bx, 1179 const T relTolerance = sqrt(T.epsilon), 1180 const T absTolerance = sqrt(T.epsilon)) 1181 if (isFloatingPoint!T && __traits(compiles, T(f(T.init)))) 1182 { 1183 if (false) // break attributes 1184 { 1185 T y = f(T(123)); 1186 } 1187 scope funInst = delegate(T x) { 1188 return T(f(x)); 1189 }; 1190 scope fun = funInst.trustedAllAttr; 1191 1192 return findLocalMinImpl(ax, bx, relTolerance, absTolerance, fun); 1193 } 1194 1195 @fmamath 1196 private FindLocalMinResult!float findLocalMinImpl( 1197 const float ax, 1198 const float bx, 1199 const float relTolerance, 1200 const float absTolerance, 1201 scope const float delegate(float) @safe pure nothrow @nogc f, 1202 ) @safe pure nothrow @nogc 1203 { 1204 pragma(inline, false); 1205 return findLocalMinImplGen!float(ax, bx, relTolerance, absTolerance, f); 1206 } 1207 1208 @fmamath 1209 private FindLocalMinResult!double findLocalMinImpl( 1210 const double ax, 1211 const double bx, 1212 const double relTolerance, 1213 const double absTolerance, 1214 scope const double delegate(double) @safe pure nothrow @nogc f, 1215 ) @safe pure nothrow @nogc 1216 { 1217 pragma(inline, false); 1218 return findLocalMinImplGen!double(ax, bx, relTolerance, absTolerance, f); 1219 } 1220 1221 @fmamath 1222 private FindLocalMinResult!real findLocalMinImpl( 1223 const real ax, 1224 const real bx, 1225 const real relTolerance, 1226 const real absTolerance, 1227 scope const real delegate(real) @safe pure nothrow @nogc f, 1228 ) @safe pure nothrow @nogc 1229 { 1230 pragma(inline, false); 1231 return findLocalMinImplGen!real(ax, bx, relTolerance, absTolerance, f); 1232 } 1233 1234 @fmamath 1235 private FindLocalMinResult!T findLocalMinImplGen(T)( 1236 const T ax, 1237 const T bx, 1238 const T relTolerance, 1239 const T absTolerance, 1240 scope const T delegate(T) @safe pure nothrow @nogc f, 1241 ) 1242 if (isFloatingPoint!T && __traits(compiles, {T _ = f(T.init);})) 1243 in 1244 { 1245 assert(relTolerance.fabs >= T.min_normal && relTolerance.fabs < T.infinity, "relTolerance is not normal floating point number"); 1246 assert(absTolerance.fabs >= T.min_normal && absTolerance.fabs < T.infinity, "absTolerance is not normal floating point number"); 1247 assert(relTolerance >= 0, "absTolerance is not positive"); 1248 assert(absTolerance >= T.epsilon*2, "absTolerance is not greater then `2*T.epsilon`"); 1249 } 1250 do 1251 { 1252 version(LDC) pragma(inline, true); 1253 // c is the squared inverse of the golden ratio 1254 // (3 - sqrt(5))/2 1255 // Value obtained from Wolfram Alpha. 1256 enum T c = 0x0.61c8864680b583ea0c633f9fa31237p+0L; 1257 enum T cm1 = 0x0.9e3779b97f4a7c15f39cc0605cedc8p+0L; 1258 T tolerance; 1259 T a = ax > bx ? bx : ax; 1260 T b = ax > bx ? ax : bx; 1261 if (a < -T.max) 1262 a = -T.max; 1263 if (b > +T.max) 1264 b = +T.max; 1265 // sequence of declarations suitable for SIMD instructions 1266 T v = a * cm1 + b * c; 1267 assert(v.fabs < T.infinity); 1268 T fv = v == v ? f(v) : v; 1269 if (fv != fv || fv == -T.infinity) 1270 { 1271 return typeof(return)(v, fv, T.init); 1272 } 1273 T w = v; 1274 T fw = fv; 1275 T x = v; 1276 T fx = fv; 1277 size_t i; 1278 for (T d = 0, e = 0;;) 1279 { 1280 i++; 1281 T m = (a + b) * 0.5f; 1282 // This fix is not part of the original algorithm 1283 if (!((m.fabs < T.infinity))) // fix infinity loop. Issue can be reproduced in R. 1284 { 1285 m = a.half + b.half; 1286 } 1287 tolerance = absTolerance * fabs(x) + relTolerance; 1288 const t2 = tolerance * 2; 1289 // check stopping criterion 1290 if (!(fabs(x - m) > t2 - (b - a) * 0.5f)) 1291 { 1292 break; 1293 } 1294 T p = 0; 1295 T q = 0; 1296 T r = 0; 1297 // fit parabola 1298 if (fabs(e) > tolerance) 1299 { 1300 const xw = x - w; 1301 const fxw = fx - fw; 1302 const xv = x - v; 1303 const fxv = fx - fv; 1304 const xwfxv = xw * fxv; 1305 const xvfxw = xv * fxw; 1306 p = xv * xvfxw - xw * xwfxv; 1307 q = (xvfxw - xwfxv) * 2; 1308 if (q > 0) 1309 p = -p; 1310 else 1311 q = -q; 1312 r = e; 1313 e = d; 1314 } 1315 T u; 1316 // a parabolic-interpolation step 1317 if (fabs(p) < fabs(q * r * 0.5f) && p > q * (a - x) && p < q * (b - x)) 1318 { 1319 d = p / q; 1320 u = x + d; 1321 // f must not be evaluated too close to a or b 1322 if (u - a < t2 || b - u < t2) 1323 d = x < m ? tolerance: -tolerance; 1324 } 1325 // a golden-section step 1326 else 1327 { 1328 e = (x < m ? b : a) - x; 1329 d = c * e; 1330 } 1331 // f must not be evaluated too close to x 1332 u = x + (fabs(d) >= tolerance ? d: d > 0 ? tolerance: -tolerance); 1333 const fu = f(u); 1334 if (fu != fu || fu == -T.infinity) 1335 { 1336 return typeof(return)(u, fu, T.init); 1337 } 1338 // update a, b, v, w, and x 1339 if (fu <= fx) 1340 { 1341 (u < x ? b : a) = x; 1342 v = w; fv = fw; 1343 w = x; fw = fx; 1344 x = u; fx = fu; 1345 } 1346 else 1347 { 1348 (u < x ? a : b) = u; 1349 if (fu <= fw || w == x) 1350 { 1351 v = w; fv = fw; 1352 w = u; fw = fu; 1353 } 1354 else if (fu <= fv || v == x || v == w) 1355 { // do not remove this braces 1356 v = u; fv = fu; 1357 } 1358 } 1359 } 1360 return typeof(return)(x, fx, tolerance * 3); 1361 } 1362 1363 /// 1364 //@nogc 1365 version(mir_test) @safe unittest 1366 { 1367 import mir.math.common: approxEqual; 1368 1369 double shift = 4; 1370 1371 auto ret = findLocalMin!(x => (x-shift)^^2)(-1e7, 1e7).validate; 1372 assert(ret.x.approxEqual(shift)); 1373 assert(ret.y.approxEqual(0.0)); 1374 } 1375 1376 /// 1377 version(mir_test) @safe unittest 1378 { 1379 import mir.math.common: approxEqual, log, fabs; 1380 alias AliasSeq(T...) = T; 1381 static foreach (T; AliasSeq!(double, float, real)) 1382 { 1383 { 1384 auto ret = findLocalMin!(x => (x-4)^^2)(T.min_normal, T(1e7)).validate; 1385 assert(ret.x.approxEqual(T(4))); 1386 assert(ret.y.approxEqual(T(0))); 1387 } 1388 { 1389 auto ret = findLocalMin!(x => fabs(x-1))(-T.max/4, T.max/4, T.min_normal, 2*T.epsilon).validate; 1390 assert(approxEqual(ret.x, T(1))); 1391 assert(approxEqual(ret.y, T(0))); 1392 assert(ret.error <= 10 * T.epsilon); 1393 } 1394 { 1395 auto ret = findLocalMin!(x => T.init)(0, 1, T.min_normal, 2*T.epsilon); 1396 assert(ret.status == FindRootStatus.nanY); 1397 } 1398 { 1399 auto ret = findLocalMin!log( 0, 1, T.min_normal, 2*T.epsilon).validate; 1400 assert(ret.error < 3.00001 * ((2*T.epsilon)*fabs(ret.x)+ T.min_normal)); 1401 assert(ret.x >= 0 && ret.x <= ret.error); 1402 } 1403 { 1404 auto ret = findLocalMin!log(0, T.max, T.min_normal, 2*T.epsilon).validate; 1405 assert(ret.y < -18); 1406 assert(ret.error < 5e-08); 1407 assert(ret.x >= 0 && ret.x <= ret.error); 1408 } 1409 { 1410 auto ret = findLocalMin!(x => -fabs(x))(-1, 1, T.min_normal, 2*T.epsilon).validate; 1411 assert(ret.x.fabs.approxEqual(T(1))); 1412 assert(ret.y.fabs.approxEqual(T(1))); 1413 assert(ret.error.approxEqual(T(0))); 1414 } 1415 } 1416 } 1417 1418 // force disabled FMA math 1419 private static auto half(T)(const T x) 1420 { 1421 pragma(inline, false); 1422 return x * 0.5f; 1423 } 1424 1425 private auto trustedAllAttr(T)(T t) @trusted 1426 { 1427 import std.traits; 1428 enum attrs = (functionAttributes!T & ~FunctionAttribute.system) 1429 | FunctionAttribute.pure_ 1430 | FunctionAttribute.safe 1431 | FunctionAttribute.nogc 1432 | FunctionAttribute.nothrow_; 1433 return cast(SetFunctionAttributes!(T, functionLinkage!T, attrs)) t; 1434 } 1435 1436 /++ 1437 Calculate the derivative of a function. 1438 This function uses Ridders' method of extrapolating the results 1439 of finite difference formulas for consecutively smaller step sizes, 1440 with an improved stopping criterion described in the Numerical Recipes 1441 books by Press et al. 1442 1443 This method gives a much higher degree of accuracy in the answer 1444 compared to a single finite difference calculation, but requires 1445 more function evaluations; typically 6 to 12. The maximum number 1446 of function evaluations is $(D 24). 1447 1448 Params: 1449 f = The function of which to take the derivative. 1450 x = The point at which to take the derivative. 1451 h = A "characteristic scale" over which the function changes. 1452 factor = Stepsize is decreased by factor at each iteration. 1453 safe = Return when error is SAFE worse than the best so far. 1454 1455 References: 1456 $(UL 1457 $(LI 1458 C. J. F. Ridders, 1459 $(I Accurate computation of F'(x) and F'(x)F''(x)). 1460 Advances in Engineering Software, vol. 4 (1982), issue 2, p. 75.) 1461 $(LI 1462 W. H. Press, S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1463 $(I Numerical Recipes in C++) (2nd ed.). 1464 Cambridge University Press, 2003.) 1465 ) 1466 +/ 1467 @fmamath 1468 DiffResult!T diff(alias f, T)(const T x, const T h, const T factor = T(2).sqrt, const T safe = 2) 1469 { 1470 if (false) // break attributes 1471 T y = f(T(1)); 1472 scope funInst = delegate(T x) { 1473 return T(f(x)); 1474 }; 1475 scope fun = funInst.trustedAllAttr; 1476 return diffImpl(fun, x, h, factor, safe); 1477 } 1478 1479 ///ditto 1480 DiffResult!T diffImpl(T) 1481 (scope const T delegate(T) @safe pure nothrow @nogc f, const T x, const T h, const T factor = T(2).sqrt, const T safe = 2) 1482 @safe pure nothrow @nogc 1483 in { 1484 assert(h < T.max); 1485 assert(h > T.min_normal); 1486 } 1487 do { 1488 // Set up Romberg tableau. 1489 import mir.ndslice.slice: sliced; 1490 pragma(inline, false); 1491 1492 enum tableauSize = 16; 1493 T[tableauSize ^^ 2] workspace = void; 1494 auto tab = workspace[].sliced(tableauSize, tableauSize); 1495 1496 // From the NR book: Stop when the difference between consecutive 1497 // approximations is bigger than SAFE*error, where error is an 1498 // estimate of the absolute error in the current (best) approximation. 1499 1500 // First approximation: A_0 1501 T result = void; 1502 T error = T.max; 1503 T hh = h; 1504 1505 tab[0,0] = (f(x + h) - f(x - h)) / (2 * h); 1506 foreach (n; 1 .. tableauSize) 1507 { 1508 // Decrease h. 1509 hh /= factor; 1510 1511 // Compute A_n 1512 tab[0, n] = (f(x + hh) - f(x - hh)) / (2 * hh); 1513 1514 T facm = 1; 1515 foreach (m; 1 .. n + 1) 1516 { 1517 facm *= factor ^^ 2; 1518 1519 // Compute B_(n-1), C_(n-2), ... 1520 T upLeft = tab[m - 1, n - 1]; 1521 T up = tab[m - 1, n]; 1522 T current = (facm * up - upLeft) / (facm - 1); 1523 tab[m, n] = current; 1524 1525 // Calculate and check error. 1526 T currentError = fmax(fabs(current - upLeft), fabs(current - up)); 1527 if (currentError <= error) 1528 { 1529 result = current; 1530 error = currentError; 1531 } 1532 } 1533 1534 if (fabs(tab[n, n] - tab[n - 1, n - 1]) >= safe * error) 1535 break; 1536 } 1537 1538 return typeof(return)(result, error); 1539 } 1540 1541 /// 1542 unittest 1543 { 1544 import mir.math.common; 1545 1546 auto f(double x) { return exp(x) / (sin(x) - x ^^ 2); } 1547 auto d(double x) { return -(exp(x) * ((x - 2) * x - sin(x) + cos(x)))/(x^^2 - sin(x))^^2; } 1548 auto r = diff!f(1.0, 0.01); 1549 assert (approxEqual(r.value, d(1))); 1550 } 1551 1552 /++ 1553 +/ 1554 struct DiffResult(T) 1555 if (__traits(isFloating, T)) 1556 { 1557 /// 1558 T value = 0; 1559 /// 1560 T error = 0; 1561 }