Report a bug
If you spot a problem with this page, click here to create a GitHub issue.
Improve this page
Quickly fork, edit online, and submit a pull request for this page. Requires a signed-in GitHub account. This works well for small changes. If you'd like to make larger changes you may want to consider using a local clone.

mir.quadrature

This module contains betterC compatible quadrature computation routines.
@nogc size_t gaussHermiteQuadrature(T)(scope T[] x, scope T[] w, scope T[] work);
Gauss-Hermite Quadrature
Parameters:
T[] x (out) user-allocated quadrature nodes in ascending order length of N
T[] w (out) user-allocated corresponding quadrature weights length of N
T[] work (temp) user-allocated workspace length of greate or equal to (N + 1) ^^ 2
Returns:
0 on success, xSTEQR LAPACK code on numerical error.
Examples:
import mir.math.common;

auto n = 5;
auto x = new double[n];
auto w = new double[n];
auto work = new double[(n + 1) ^^ 2];

gaussHermiteQuadrature(x, w, work);

static immutable xc =
   [-2.02018287,
    -0.95857246,
     0.        ,
     0.95857246,
     2.02018287];

static immutable wc =
   [0.01995324,
    0.39361932,
    0.94530872,
    0.39361932,
    0.01995324];

foreach (i; 0 .. n)
{
    assert(x[i].approxEqual(xc[i]));
    assert(w[i].approxEqual(wc[i]));
}
@nogc size_t gaussJacobiQuadrature(T)(scope T[] x, scope T[] w, scope T[] work, const T alpha, const T beta);
Gauss-Jacobi Quadrature
Parameters:
T[] x (out) user-allocated quadrature nodes in ascending order length of N
T[] w (out) user-allocated corresponding quadrature weights length of N
T[] work (temp) user-allocated workspace length of greate or equal to (N + 1) ^^ 2
T alpha parameter '> -1'
T beta parameter '> -1'
Returns:
0 on success, xSTEQR LAPACK code on numerical error.
Examples:
import mir.math.common;

auto n = 5;
auto x = new double[n];
auto w = new double[n];
auto work = new double[(n + 1) ^^ 2];

gaussJacobiQuadrature(x, w, work, 2.3, 3.6);

static immutable xc =
   [-0.6553677 ,
    -0.29480426,
     0.09956621,
     0.47584565,
     0.78356514];

static immutable wc =
   [0.02262392,
    0.19871672,
    0.43585107,
    0.32146619,
    0.0615342 ];

foreach (i; 0 .. n)
{
    assert(x[i].approxEqual(xc[i]));
    assert(w[i].approxEqual(wc[i]));
}
@nogc size_t gaussLaguerreQuadrature(T)(scope T[] x, scope T[] w, scope T[] work, T alpha = 0);
Gauss-Laguerre Quadrature
Parameters:
T[] x (out) user-allocated quadrature nodes in ascending order length of N
T[] w (out) user-allocated corresponding quadrature weights length of N
T[] work (temp) user-allocated workspace length of greate or equal to (N + 1) ^^ 2
T alpha (optional) parameter '> -1'
Returns:
0 on success, xSTEQR LAPACK code on numerical error.
Examples:
import mir.math.common;

auto n = 5;
auto x = new double[n];
auto w = new double[n];
auto work = new double[(n + 1) ^^ 2];

gaussLaguerreQuadrature(x, w, work);

static immutable xc =
   [ 0.26356032,
     1.41340306,
     3.59642577,
     7.08581001,
    12.64080084];

static immutable wc =
   [5.21755611e-01,
    3.98666811e-01,
    7.59424497e-02,
    3.61175868e-03,
    2.33699724e-05];

foreach (i; 0 .. n)
{
    assert(x[i].approxEqual(xc[i]));
    assert(w[i].approxEqual(wc[i]));
}
@nogc size_t gaussLegendreQuadrature(T)(scope T[] x, scope T[] w, scope T[] work);
Gauss-Legendre Quadrature
Parameters:
T[] x (out) user-allocated quadrature nodes in ascending order length of N
T[] w (out) user-allocated corresponding quadrature weights length of N
T[] work (temp) user-allocated workspace length of greate or equal to (N + 1) ^^ 2
Returns:
0 on success, xSTEQR LAPACK code on numerical error.
Examples:
import mir.math.common;

auto n = 5;
auto x = new double[n];
auto w = new double[n];
auto work = new double[(n + 1) ^^ 2];

gaussLegendreQuadrature(x, w, work);

static immutable xc =
   [-0.90617985,
    -0.53846931,
     0.        ,
     0.53846931,
     0.90617985];

static immutable wc =
   [0.23692689,
    0.47862867,
    0.56888889,
    0.47862867,
    0.23692689];

foreach (i; 0 .. n)
{
    assert(x[i].approxEqual(xc[i]));
    assert(w[i].approxEqual(wc[i]));
}
@nogc size_t gaussLobattoQuadrature(T)(scope T[] x, scope T[] w, scope T[] work, const T alpha = 0, const T beta = 0);
Gauss-Lobatto Quadrature
Parameters:
T[] x (out) user-allocated quadrature nodes in ascending order length of N
T[] w (out) user-allocated corresponding quadrature weights length of N
T[] work (temp) user-allocated workspace length of greate or equal to N ^^ 2
T alpha (optional) parameter '> -1'
T beta (optional) parameter '> -1'
Returns:
0 on success, xSTEQR LAPACK code on numerical error.
Examples:
import mir.math.common;

auto n = 7;
auto x = new double[n];
auto w = new double[n];
auto work = new double[n ^^ 2];

gaussLobattoQuadrature(x, w, work);

static immutable xc =
   [-1,
    -sqrt(5.0 / 11 + 2.0 / 11 * sqrt(5.0 / 3)),
    -sqrt(5.0 / 11 - 2.0 / 11 * sqrt(5.0 / 3)),
     0.        ,
    +sqrt(5.0 / 11 - 2.0 / 11 * sqrt(5.0 / 3)),
    +sqrt(5.0 / 11 + 2.0 / 11 * sqrt(5.0 / 3)),
     1];

static immutable wc =
   [1.0 / 21,
    (124.0 - 7 * sqrt(15.0)) / 350,
    (124.0 + 7 * sqrt(15.0)) / 350,
    256.0 / 525,
    (124.0 + 7 * sqrt(15.0)) / 350,
    (124.0 - 7 * sqrt(15.0)) / 350,
    1.0 / 21];

foreach (i; 0 .. n)
    assert(x[i].approxEqual(xc[i]));
foreach (i; 0 .. n)
    assert(w[i].approxEqual(wc[i]));