예제 #1
0
파일: legendre.py 프로젝트: jonathf/chaospy
def legendre_simple(order):
    """
    Simple Legendre quadrature on the [0, 1] interval.

    Use :func:`chaospy.quadrature.legendre` instead.
    """
    coefficients = chaospy.construct_recurrence_coefficients(
        order=int(order), dist=chaospy.Uniform(-1, 1))
    [abscissas], [weights] = chaospy.coefficients_to_quadrature(coefficients)
    return abscissas * 0.5 + 0.5, weights
예제 #2
0
파일: jacobi.py 프로젝트: jonathf/chaospy
def jacobi(order, alpha, beta, lower=-1, upper=1, physicist=False):
    """
    Gauss-Jacobi quadrature rule.

    Compute the sample points and weights for Gauss-Jacobi quadrature. The
    sample points are the roots of the nth degree Jacobi polynomial. These
    sample points and weights correctly integrate polynomials of degree
    :math:`2N-1` or less.

    Gaussian quadrature come in two variants: physicist and probabilist. For
    Gauss-Jacobi physicist means a weight function
    :math:`(1-x)^\alpha (1+x)^\beta` and
    weights that sum to :math`2^{\alpha+\beta}`, and probabilist means a weight
    function is :math:`B(\alpha, \beta) x^{\alpha-1}(1-x)^{\beta-1}` (where
    :math:`B` is the beta normalizing constant) which sum to 1.

    Args:
        order (int):
            The quadrature order.
        alpha (float):
            First Jakobi shape parameter.
        beta (float):
            Second Jakobi shape parameter.
        lower (float):
            Lower bound for the integration interval.
        upper (float):
            Upper bound for the integration interval.
        physicist (bool):
            Use physicist weights instead of probabilist.

    Returns:
        abscissas (numpy.ndarray):
            The ``order+1`` quadrature points for where to evaluate the model
            function with.
        weights (numpy.ndarray):
            The quadrature weights associated with each abscissas.

    Examples:
        >>> abscissas, weights = chaospy.quadrature.jacobi(3, alpha=2, beta=2)
        >>> abscissas
        array([[-0.69474659, -0.25056281,  0.25056281,  0.69474659]])
        >>> weights
        array([0.09535261, 0.40464739, 0.40464739, 0.09535261])

    See also:
        :func:`chaospy.quadrature.gaussian`

    """
    order = int(order)
    coefficients = chaospy.construct_recurrence_coefficients(
        order=order, dist=chaospy.Beta(alpha + 1, beta + 1, lower, upper))
    [abscissas], [weights] = chaospy.coefficients_to_quadrature(coefficients)
    weights *= 2**(alpha + beta) if physicist else 1
    return abscissas[numpy.newaxis], weights
예제 #3
0
파일: hermite.py 프로젝트: jonathf/chaospy
def hermite(order, mu=0.0, sigma=1.0, physicist=False):
    r"""
    Gauss-Hermite quadrature rule.

    Compute the sample points and weights for Gauss-Hermite quadrature. The
    sample points are the roots of the nth degree Hermite polynomial. These
    sample points and weights correctly integrate polynomials of degree
    :math:`2N-1` or less.

    Gaussian quadrature come in two variants: physicist and probabilist. For
    Gauss-Hermite physicist means a weight function :math:`e^{-x^2}` and
    weights that sum to :math`\sqrt(\pi)`, and probabilist means a weight
    function is :math:`e^{-x^2/2}` and sum to 1.

    Args:
        order (int):
            The quadrature order.
        mu (float):
            Non-centrality parameter.
        sigma (float):
            Scale parameter.
        physicist (bool):
            Use physicist weights instead of probabilist variant.

    Returns:
        abscissas (numpy.ndarray):
            The ``order+1`` quadrature points for where to evaluate the model
            function with.
        weights (numpy.ndarray):
            The quadrature weights associated with each abscissas.

    Examples:
        >>> abscissas, weights = chaospy.quadrature.hermite(3)
        >>> abscissas
        array([[-2.33441422, -0.74196378,  0.74196378,  2.33441422]])
        >>> weights
        array([0.04587585, 0.45412415, 0.45412415, 0.04587585])

    See also:
        :func:`chaospy.quadrature.gaussian`

    """
    order = int(order)
    sigma = float(sigma * 2**-0.5 if physicist else sigma)
    coefficients = chaospy.construct_recurrence_coefficients(
        order=order, dist=chaospy.Normal(0, sigma))
    [abscissas], [weights] = chaospy.coefficients_to_quadrature(coefficients)
    weights = weights * numpy.pi**0.5 if physicist else weights
    if order % 2 == 0:
        abscissas[len(abscissas) // 2] = 0
    abscissas += mu
    return abscissas[numpy.newaxis], weights
예제 #4
0
파일: laguerre.py 프로젝트: jonathf/chaospy
def laguerre(order, alpha=0.0, physicist=False):
    r"""
    Generalized Gauss-Laguerre quadrature rule.

    Compute the sample points and weights for Gauss-Laguerre quadrature. The
    sample points are the roots of the nth degree Laguerre polynomial. These
    sample points and weights correctly integrate polynomials of degree
    :math:`2N-1` or less.

    Gaussian quadrature come in two variants: physicist and probabilist. For
    Gauss-Laguerre physicist means a weight function :math:`x^\alpha e^{-x}`
    and weights that sum to :math`\Gamma(\alpha+1)`, and probabilist means a
    weight function is :math:`x^\alpha e^{-x}` and sum to 1.

    Args:
        order (int):
            The quadrature order.
        alpha (float):
            Shape parameter. Defaults to non-generalized Laguerre if 0.
        physicist (bool):
            Use physicist weights instead of probabilist.

    Returns:
        abscissas (numpy.ndarray):
            The ``order+1`` quadrature points for where to evaluate the model
            function with.
        weights (numpy.ndarray):
            The quadrature weights associated with each abscissas.

    Examples:
        >>> abscissas, weights = chaospy.quadrature.laguerre(2)
        >>> abscissas
        array([[0.41577456, 2.29428036, 6.28994508]])
        >>> weights
        array([0.71109301, 0.27851773, 0.01038926])

    See also:
        :func:`chaospy.quadrature.gaussian`

    """
    order = int(order)
    coefficients = chaospy.construct_recurrence_coefficients(
        order=order, dist=chaospy.Gamma(alpha + 1))
    [abscissas], [weights] = chaospy.coefficients_to_quadrature(coefficients)
    weights *= gamma(alpha + 1) if physicist else 1
    return abscissas[numpy.newaxis], weights
예제 #5
0
def quad_gauss_kronrod(
    order,
    dist,
    recurrence_algorithm="stieltjes",
    rule="clenshaw_curtis",
    tolerance=1e-10,
    scaling=3,
    n_max=5000,
):
    """
    Generate Gauss-Kronrod quadrature abscissas and weights.

    Gauss-Kronrod quadrature is an adaptive method for Gaussian quadrature
    rule. It builds on top of other quadrature rules by extending "pure"
    Gaussian quadrature rules with extra abscissas and new weights such that
    already used abscissas can be reused. For more details, see `Wikipedia
    article
    <https://en.wikipedia.org/wiki/Gauss%E2%80%93Kronrod_quadrature_formula>`_.

    For each order ``N`` taken with ordinary Gaussian quadrature, Gauss-Kronrod
    will create ``2N+1`` abscissas where all of the ``N`` "old" abscissas are
    all interlaced between the "new" ones.

    The algorithm is well suited for any Jacobi scheme, i.e. quadrature
    involving Uniform or Beta distribution, and might work on others as well.
    However, it will not work everywhere. For example `Kahaner and Monegato
    <https://link.springer.com/article/10.1007/BF01590820>`_ showed that higher
    order Gauss-Kronrod quadrature for Gauss-Hermite and Gauss-Laguerre does
    not exist.

    Args:
        order (int):
            The order of the quadrature.
        dist (chaospy.distributions.baseclass.Distribution):
            The distribution which density will be used as weight function.
        recurrence_algorithm (str):
            Name of the algorithm used to generate abscissas and weights. If
            omitted, ``analytical`` will be tried first, and ``stieltjes`` used
            if that fails.
        rule (str):
            In the case of ``lanczos`` or ``stieltjes``, defines the
            proxy-integration scheme.
        tolerance (float):
            The allowed relative error in norm between two quadrature orders
            before method assumes convergence.
        scaling (float):
            A multiplier the adaptive order increases with for each step
            quadrature order is not converged. Use 0 to indicate unit
            increments.
        n_max (int):
            The allowed number of quadrature points to use in approximation.

    Returns:
        (numpy.ndarray, numpy.ndarray):
            abscissas:
                The quadrature points for where to evaluate the model function
                with ``abscissas.shape == (len(dist), N)`` where ``N`` is the
                number of samples.
            weights:
                The quadrature weights with ``weights.shape == (N,)``.

    Raises:
        ValueError:
            Error raised if Kronrod algorithm results in negative recurrence
            coefficients.

    Notes:
        Code is adapted from `quadpy <https://github.com/nschloe/quadpy>`_,
        which adapted his code from `W. Gautschi
        <https://www.cs.purdue.edu/archives/2002/wxg/codes/OPQ.html>`_.
        Algorithm for calculating Kronrod-Jacobi matrices was first published
        Algorithm as proposed by Laurie :cite:`laurie_calculation_1997`.

    Example:
        >>> distribution = chaospy.Uniform(-1, 1)
        >>> abscissas, weights = quad_gauss_kronrod(3, distribution)
        >>> abscissas.round(2)
        array([[-0.98, -0.86, -0.64, -0.34,  0.  ,  0.34,  0.64,  0.86,  0.98]])
        >>> weights.round(3)
        array([0.031, 0.085, 0.133, 0.163, 0.173, 0.163, 0.133, 0.085, 0.031])
    """
    assert not rule.startswith("gauss"), "recursive Gaussian quadrature call"

    length = int(numpy.ceil(3 * (order + 1) / 2.0))
    coefficients = chaospy.construct_recurrence_coefficients(
        order=length,
        dist=dist,
        recurrence_algorithm=recurrence_algorithm,
        rule=rule,
        tolerance=tolerance,
        scaling=scaling,
        n_max=n_max,
    )

    coefficients = [
        kronrod_jacobi(order + 1, coeffs) for coeffs in coefficients
    ]

    abscissas, weights = chaospy.coefficients_to_quadrature(coefficients)
    return combine_quadrature(abscissas, weights)
예제 #6
0
def radau(
    order,
    dist,
    fixed_point=None,
    recurrence_algorithm="stieltjes",
    rule="clenshaw_curtis",
    tolerance=1e-10,
    scaling=3,
    n_max=5000,
):
    """
    Generate the quadrature nodes and weights in Gauss-Radau quadrature.

    Gauss-Radau formula for numerical estimation of integrals. It requires
    :math:`m+1` points and fits all Polynomials to degree :math:`2m`, so it
    effectively fits exactly all Polynomials of degree :math:`2m-1`.

    It allows for a single abscissas to be user defined, while the others are
    built around this point.

    Canonically, Radau is built around Legendre weight function with the fixed
    point at the left end. Not all distributions/fixed point combinations
    allows for the building of a quadrature scheme.

    Args:
        order (int):
            Quadrature order.
        dist (:class:`chaospy.Distribution`):
            The distribution weights to be used to create higher order nodes
            from.
        fixed_point (float):
            Fixed point abscissas assumed to be included in the quadrature. If
            omitted, use distribution lower bound.
        rule (str):
            In the case of ``lanczos`` or ``stieltjes``, defines the
            proxy-integration scheme.
        accuracy (int):
            In the case ``rule`` is used, defines the quadrature order of the
            scheme used. In practice, must be at least as large as ``order``.
        recurrence_algorithm (str):
            Name of the algorithm used to generate abscissas and weights. If
            omitted, ``analytical`` will be tried first, and ``stieltjes`` used
            if that fails.

    Returns:
        abscissas (numpy.ndarray):
            The quadrature points for where to evaluate the model function
            with ``abscissas.shape == (len(dist), N)`` where ``N`` is the
            number of samples.
        weights (numpy.ndarray):
            The quadrature weights with ``weights.shape == (N,)``.

    Example:
        >>> distribution = chaospy.Uniform(-1, 1)
        >>> abscissas, weights = chaospy.quadrature.radau(4, distribution)
        >>> abscissas.round(3)
        array([[-1.   , -0.887, -0.64 , -0.295,  0.094,  0.468,  0.771,  0.955]])
        >>> weights.round(3)
        array([0.016, 0.093, 0.152, 0.188, 0.196, 0.174, 0.125, 0.057])

    """
    if fixed_point is None:
        fixed_point = dist.lower
    else:
        fixed_point = numpy.ones(len(dist)) * fixed_point

    if order == 0:
        return fixed_point.reshape(-1, 1), numpy.ones(1)

    coefficients = chaospy.construct_recurrence_coefficients(
        order=2 * order - 1,
        dist=dist,
        recurrence_algorithm=recurrence_algorithm,
        rule=rule,
        tolerance=tolerance,
        scaling=scaling,
        n_max=n_max,
    )
    coefficients = [
        radau_jakobi(coeffs, point) for point, coeffs in zip(fixed_point, coefficients)
    ]

    abscissas, weights = chaospy.coefficients_to_quadrature(coefficients)

    return combine_quadrature(abscissas, weights)
예제 #7
0
def quad_gauss_radau(
        order,
        dist,
        fixed_point=None,
        recurrence_algorithm="stieltjes",
        rule="clenshaw_curtis",
        tolerance=1e-10,
        scaling=3,
        n_max=5000,
):
    """
    Generate the quadrature nodes and weights in Gauss-Radau quadrature.

    Args:
        order (int):
            Quadrature order.
        dist (chaospy.distributions.baseclass.Distribution):
            The distribution weights to be used to create higher order nodes
            from.
        fixed_point (float):
            Fixed point abscissas assumed to be included in the quadrature. If
            imitted, use distribution lower point ``dist.range()[0]``.
        rule (str):
            In the case of ``lanczos`` or ``stieltjes``, defines the
            proxy-integration scheme.
        accuracy (int):
            In the case ``rule`` is used, defines the quadrature order of the
            scheme used. In practice, must be at least as large as ``order``.
        recurrence_algorithm (str):
            Name of the algorithm used to generate abscissas and weights. If
            omitted, ``analytical`` will be tried first, and ``stieltjes`` used
            if that fails.

    Returns:
        (numpy.ndarray, numpy.ndarray):
            abscissas:
                The quadrature points for where to evaluate the model function
                with ``abscissas.shape == (len(dist), N)`` where ``N`` is the
                number of samples.
            weights:
                The quadrature weights with ``weights.shape == (N,)``.

    Example:
        >>> abscissas, weights = quad_gauss_radau(4, chaospy.Uniform(-1, 1))
        >>> abscissas.round(3)
        array([[-1.   , -0.887, -0.64 , -0.295,  0.094,  0.468,  0.771,  0.955]])
        >>> weights.round(3)
        array([0.016, 0.093, 0.152, 0.188, 0.196, 0.174, 0.125, 0.057])
    """
    assert not rule.startswith("gauss"), "recursive Gaussian quadrature call"
    if fixed_point is None:
        fixed_point = dist.lower
    else:
        fixed_point = numpy.ones(len(dist))*fixed_point

    if order == 0:
        return fixed_point.reshape(-1, 1), numpy.ones(1)

    coefficients = chaospy.construct_recurrence_coefficients(
        order=2*order-1,
        dist=dist,
        recurrence_algorithm=recurrence_algorithm,
        rule=rule,
        tolerance=tolerance,
        scaling=scaling,
        n_max=n_max,
    )
    coefficients = [radau_jakobi(coeffs, point)
                    for point, coeffs in zip(fixed_point, coefficients)]

    abscissas, weights = chaospy.coefficients_to_quadrature(coefficients)

    return combine_quadrature(abscissas, weights)
예제 #8
0
def quad_gauss_lobatto(
    order,
    dist,
    recurrence_algorithm="stieltjes",
    rule="fejer",
    tolerance=1e-10,
    scaling=3,
    n_max=5000,
):
    """
    Generate the abscissas and weights in Gauss-Loboto quadrature.

    Args:
        order (int):
            Quadrature order.
        dist (chaospy.distributions.baseclass.Distribution):
            The distribution weights to be used to create higher order nodes
            from.
        recurrence_algorithm (str):
            Name of the algorithm used to generate abscissas and weights.
        rule (str):
            In the case of ``lanczos`` or ``stieltjes``, defines the
            proxy-integration scheme.
        tolerance (float):
            The allowed relative error in norm between two quadrature orders
            before method assumes convergence.
        scaling (float):
            A multiplier the adaptive order increases with for each step
            quadrature order is not converged. Use 0 to indicate unit
            increments.
        n_max (int):
            The allowed number of quadrature points to use in approximation.

    Returns:
        (numpy.ndarray, numpy.ndarray):
            abscissas:
                The quadrature points for where to evaluate the model function
                with ``abscissas.shape == (len(dist), N)`` where ``N`` is the
                number of samples.
            weights:
                The quadrature weights with ``weights.shape == (N,)``.

    Example:
        >>> abscissas, weights = quad_gauss_lobatto(
        ...     4, chaospy.Uniform(-1, 1))
        >>> abscissas.round(3)
        array([[-1.   , -0.872, -0.592, -0.209,  0.209,  0.592,  0.872,  1.   ]])
        >>> weights.round(3)
        array([0.018, 0.105, 0.171, 0.206, 0.206, 0.171, 0.105, 0.018])
    """
    assert not rule.startswith("gauss"), "recursive Gaussian quadrature call"
    if order == 0:
        return dist.lower.reshape(1, -1), numpy.array([1.])

    coefficients = chaospy.construct_recurrence_coefficients(
        order=2 * order - 1,
        dist=dist,
        recurrence_algorithm=recurrence_algorithm,
        rule=rule,
        tolerance=tolerance,
        scaling=scaling,
        n_max=n_max,
    )
    coefficients = [
        _lobatto(coeffs, (lo, up))
        for coeffs, lo, up in zip(coefficients, dist.lower, dist.upper)
    ]
    abscissas, weights = chaospy.coefficients_to_quadrature(coefficients)

    return combine_quadrature(abscissas, weights)
예제 #9
0
def quad_gauss_legendre(
    order,
    domain=(0, 1),
    recurrence_algorithm="stieltjes",
    rule="clenshaw_curtis",
    tolerance=1e-10,
    scaling=3,
    n_max=5000,
):
    r"""
    Generate the quadrature nodes and weights in Gauss-Legendre quadrature.

    Note that this rule exists to allow for integrating functions with weight
    functions without actually adding the quadrature. Like:

    .. math:
        \int_a^b p(x) f(x) dx \approx \sum_i p(X_i) f(X_i) W_i

    instead of the more traditional:

    .. math:
        \int_a^b p(x) f(x) dx \approx \sum_i f(X_i) W_i

    To get the behavior where the weight function is taken into consideration,
    use :func:`chaospy.quad_gaussian`.

    Args:
        order (int, numpy.ndarray):
            Quadrature order.
        domain (chaospy.distributions.baseclass.Distribution, numpy.ndarray):
            Either distribution or bounding of interval to integrate over.
        recurrence_algorithm (str):
            Name of the algorithm used to generate abscissas and weights.
        rule (str):
            In the case of ``lanczos`` or ``stieltjes``, defines the
            proxy-integration scheme.
        tolerance (float):
            The allowed relative error in norm between two quadrature orders
            before method assumes convergence.
        scaling (float):
            A multiplier the adaptive order increases with for each step
            quadrature order is not converged. Use 0 to indicate unit
            increments.
        n_max (int):
            The allowed number of quadrature points to use in approximation.

    Returns:
        (numpy.ndarray, numpy.ndarray):
            abscissas:
                The quadrature points for where to evaluate the model function
                with ``abscissas.shape == (len(dist), N)`` where ``N`` is the
                number of samples.
            weights:
                The quadrature weights with ``weights.shape == (N,)``.

    Example:
        >>> abscissas, weights = quad_gauss_legendre(3)
        >>> abscissas.round(4)
        array([[0.0694, 0.33  , 0.67  , 0.9306]])
        >>> weights.round(4)
        array([0.1739, 0.3261, 0.3261, 0.1739])
    """
    from ..distributions.baseclass import Distribution
    from ..distributions.collection import Uniform
    if isinstance(domain, Distribution):
        abscissas, weights = quad_gauss_legendre(
            order=order,
            domain=(domain.lower, domain.upper),
            recurrence_algorithm=recurrence_algorithm,
            rule=rule,
            tolerance=tolerance,
            scaling=scaling,
            n_max=n_max,
        )
        eps = 1e-14 * (domain.upper - domain.lower)
        abscissas_ = numpy.clip(abscissas.T, domain.lower + eps,
                                domain.upper - eps).T
        weights *= domain.pdf(abscissas_).flatten()
        weights /= numpy.sum(weights)
        return abscissas, weights

    order = numpy.asarray(order, dtype=int).flatten()
    lower, upper = numpy.array(domain)
    lower = numpy.asarray(lower).flatten()
    upper = numpy.asarray(upper).flatten()

    dim = max(lower.size, upper.size, order.size)
    order = numpy.ones(dim, dtype=int) * order
    lower = numpy.ones(dim) * lower
    upper = numpy.ones(dim) * upper

    coefficients = chaospy.construct_recurrence_coefficients(
        order=numpy.max(order),
        dist=Uniform(0, 1),
        recurrence_algorithm=recurrence_algorithm,
        rule=rule,
        tolerance=tolerance,
        scaling=scaling,
        n_max=n_max,
    )

    abscissas, weights = zip(*[
        chaospy.coefficients_to_quadrature(coefficients[:order_ + 1])
        for order_ in order
    ])
    abscissas = list(numpy.asarray(abscissas).reshape(dim, -1))
    weights = list(numpy.asarray(weights).reshape(dim, -1))

    return combine_quadrature(abscissas, weights, (lower, upper))
예제 #10
0
def quad_gaussian(
    order,
    dist,
    recurrence_algorithm="stieltjes",
    rule="clenshaw_curtis",
    tolerance=1e-10,
    scaling=3,
    n_max=5000,
):
    """
    Create Gaussian quadrature nodes and weights.

    Generating Gaussian quadrature by first generating so called *three terms
    recurrence* coefficients using one various different algorithms. The
    coefficients are them converted to abscissas and weights by constructing
    lower banded Jacobi matrix and calculating eigenvalues and eigenvectors,
    which can be directly translated to abscissas and weights. Construction of
    the coefficients is potentially numerically unstable, there for multiple
    algorithms exists.

    Args:
        dist (chaospy.Distribution):
            The distribution which density will be used as weight function.
        order (int):
            The order of the quadrature.
        recurrence_algorithm (str):
            Name of the algorithm used to generate abscissas and weights.
        rule (str):
            In the case of ``lanczos`` or ``stieltjes``, defines the
            proxy-integration scheme.
        tolerance (float):
            The allowed relative error in norm between two quadrature orders
            before method assumes convergence.
        scaling (float):
            A multiplier the adaptive order increases with for each step
            quadrature order is not converged. Use 0 to indicate unit
            increments.
        n_max (int):
            The allowed number of quadrature points to use in approximation.

    Returns:
        (numpy.ndarray, numpy.ndarray):
            Gaussian quadrature abscissas and weights with shapes
            ``(len(dist), order+1)`` and ``(order+1,)`` respectively.

    Raises:
        NotImplementedError:
            In the case of recurrence algorithm ``analytical``, error is raised
            if the distribution does not implement the three terms recurrence
            algorithm analytically.
        numpy.linalg.LinAlgError:
            For non-canonical random variables, the construction might fail
            because of illegal numerical operations.

    Notes:
        Quadrature constructed as outlined by Walter Gautschi
        :cite:`gautschi_construction_1968`.

    Examples:
        >>> distribution = chaospy.Normal(0, 1)
        >>> abscissas, weights = chaospy.quad_gaussian(
        ...     5, distribution, recurrence_algorithm="stieltjes")
        >>> abscissas.round(4)
        array([[-3.3243, -1.8892, -0.6167,  0.6167,  1.8892,  3.3243]])
        >>> weights.round(4)
        array([0.0026, 0.0886, 0.4088, 0.4088, 0.0886, 0.0026])
        >>> distribution = chaospy.J(chaospy.Uniform(), chaospy.Normal())
        >>> abscissas, weights = chaospy.quad_gaussian(
        ...     2, distribution, recurrence_algorithm="chebyshev")
        >>> abscissas.round(2)
        array([[ 0.11,  0.11,  0.11,  0.5 ,  0.5 ,  0.5 ,  0.89,  0.89,  0.89],
               [-1.73,  0.  ,  1.73, -1.73,  0.  ,  1.73, -1.73,  0.  ,  1.73]])
        >>> weights.round(3)
        array([0.046, 0.185, 0.046, 0.074, 0.296, 0.074, 0.046, 0.185, 0.046])
    """
    coefficients = chaospy.construct_recurrence_coefficients(
        order=order,
        dist=dist,
        recurrence_algorithm=recurrence_algorithm,
        rule=rule,
        tolerance=tolerance,
        scaling=scaling,
        n_max=n_max,
    )
    abscissas, weights = chaospy.coefficients_to_quadrature(coefficients)
    return combine_quadrature(abscissas, weights)