Esempio n. 1
0
    def __init__(self, order, dims):
        """
        :arg order: The total degree to which the quadrature rule is exact
            for *interpolation*.
        :arg dims: The number of dimensions for the quadrature rule.
            2 for quadrature on triangles and 3 for tetrahedra.
        """

        if dims == 2:
            from modepy.quadrature.vr_quad_data_tri import triangle_data as table
            ref_volume = 2
        elif dims == 3:
            from modepy.quadrature.vr_quad_data_tet import tetrahedron_data as table
            ref_volume = 4 / 3
        else:
            raise ValueError("invalid dimensionality")

        from modepy.tools import EQUILATERAL_TO_UNIT_MAP
        e2u = EQUILATERAL_TO_UNIT_MAP[dims]

        try:
            order_table = table[order]
        except KeyError:
            raise QuadratureRuleUnavailable

        nodes = e2u(order_table["points"])
        wts = order_table["weights"]
        wts = wts * (ref_volume / np.sum(wts))

        Quadrature.__init__(self, nodes, wts)

        self.exact_to = order_table["quad_degree"]
    def __init__(self, order, dimension):
        """
        :arg order: A parameter correlated with the total degree of polynomials
            that are integrated exactly. (See also :attr:`exact_to`.)
        :arg dimension: The number of dimensions for the quadrature rule.
            Any positive integer.
        """
        s = order
        n = dimension
        d = 2*s + 1

        self.exact_to = d

        if dimension == 0:
            nodes = np.zeros((dimension, 1))
            weights = np.ones(1)

            Quadrature.__init__(self, nodes, weights)
            return

        from pytools import \
                generate_decreasing_nonnegative_tuples_summing_to, \
                generate_unique_permutations, \
                factorial, \
                wandering_element

        points_to_weights = {}

        for i in range(s + 1):
            weight = (-1)**i * 2**(-2*s) \
                    * (d + n - 2*i)**d \
                    / factorial(i) \
                    / factorial(d + n - i)

            for t in generate_decreasing_nonnegative_tuples_summing_to(s - i, n + 1):
                for beta in generate_unique_permutations(t):
                    denominator = d + n - 2*i
                    point = tuple(
                            _simplify_fraction((2*beta_i + 1, denominator))
                            for beta_i in beta)

                    points_to_weights[point] = \
                            points_to_weights.get(point, 0) + weight

        from operator import add

        vertices = ([-1 * np.ones((n,))]
                + [np.array(x)
                    for x in wandering_element(n, landscape=-1, wanderer=1)])

        nodes = []
        weights = []

        dim_factor = 2**n
        for p, w in points_to_weights.items():
            real_p = reduce(add, (a/b * v for (a, b), v in zip(p, vertices)))
            nodes.append(real_p)
            weights.append(dim_factor * w)

        Quadrature.__init__(self, np.array(nodes).T, np.array(weights))
Esempio n. 3
0
    def __init__(self, order, dims):
        """
        :arg order: The total degree to which the quadrature rule is exact.
        :arg dims: The number of dimensions for the quadrature rule.
            2 for quadrature on triangles and 3 for tetrahedra.
        """

        if dims == 2:
            from modepy.quadrature.xg_quad_data import triangle_table as table
        elif dims == 3:
            from modepy.quadrature.xg_quad_data import tetrahedron_table as table
        else:
            raise ValueError("invalid dimensionality")
        try:
            order_table = table[order]
        except KeyError:
            raise QuadratureRuleUnavailable

        from modepy.tools import EQUILATERAL_TO_UNIT_MAP
        e2u = EQUILATERAL_TO_UNIT_MAP[dims]

        nodes = e2u(order_table["points"].T)
        wts = order_table["weights"] * e2u.jacobian

        Quadrature.__init__(self, nodes, wts)

        self.exact_to = order
Esempio n. 4
0
    def __init__(self, order, dims):
        """
        :arg order: The total degree to which the quadrature rule is exact.
        :arg dims: The number of dimensions for the quadrature rule.
            2 for quadrature on triangles and 3 for tetrahedra.
        """

        if dims == 2:
            from modepy.quadrature.xg_quad_data import triangle_table as table
        elif dims == 3:
            from modepy.quadrature.xg_quad_data import tetrahedron_table as table
        else:
            raise ValueError("invalid dimensionality")
        try:
            order_table = table[order]
        except KeyError:
            raise QuadratureRuleUnavailable

        from modepy.tools import EQUILATERAL_TO_UNIT_MAP
        e2u = EQUILATERAL_TO_UNIT_MAP[dims]

        nodes = e2u(order_table["points"].T)
        wts = order_table["weights"]*e2u.jacobian

        Quadrature.__init__(self, nodes, wts)

        self.exact_to = order
Esempio n. 5
0
    def __init__(self, order, dims):
        """
        :arg order: The total degree to which the quadrature rule is exact
            for *interpolation*.
        :arg dims: The number of dimensions for the quadrature rule.
            2 for quadrature on triangles and 3 for tetrahedra.
        """

        if dims == 2:
            from modepy.quadrature.vr_quad_data_tri import triangle_data as table
            ref_volume = 2
        elif dims == 3:
            from modepy.quadrature.vr_quad_data_tet import tetrahedron_data as table
            ref_volume = 4/3
        else:
            raise ValueError("invalid dimensionality")

        from modepy.tools import EQUILATERAL_TO_UNIT_MAP
        e2u = EQUILATERAL_TO_UNIT_MAP[dims]

        try:
            order_table = table[order]
        except KeyError:
            raise QuadratureRuleUnavailable

        nodes = e2u(order_table["points"])
        wts = order_table["weights"]
        wts = wts * (ref_volume/np.sum(wts))

        Quadrature.__init__(self, nodes, wts)

        self.exact_to = order_table["quad_degree"]
Esempio n. 6
0
    def __init__(self, order, dimension):
        """
        :arg order: A parameter correlated with the total degree of polynomials
            that are integrated exactly. (See also :attr:`exact_to`.)
        :arg dimension: The number of dimensions for the quadrature rule.
            Any positive integer.
        """
        s = order
        n = dimension
        d = 2*s+1

        self.exact_to = d

        if dimension == 0:
            nodes = np.zeros((dimension, 1))
            weights = np.ones(1)

            Quadrature.__init__(self, nodes, weights)
            return

        from pytools import \
                generate_decreasing_nonnegative_tuples_summing_to, \
                generate_unique_permutations, \
                factorial, \
                wandering_element

        points_to_weights = {}

        for i in range(s+1):
            weight = (-1)**i * 2**(-2*s) \
                    * (d + n-2*i)**d \
                    / factorial(i) \
                    / factorial(d+n-i)

            for t in generate_decreasing_nonnegative_tuples_summing_to(s-i, n+1):
                for beta in generate_unique_permutations(t):
                    denominator = d+n-2*i
                    point = tuple(
                            _simplify_fraction((2*beta_i+1, denominator))
                            for beta_i in beta)

                    points_to_weights[point] = \
                            points_to_weights.get(point, 0) + weight

        from operator import add

        vertices = [-1 * np.ones((n,))] \
                + [np.array(x)
                        for x in wandering_element(n, landscape=-1, wanderer=1)]

        nodes = []
        weights = []

        dim_factor = 2**n
        for p, w in six.iteritems(points_to_weights):
            real_p = reduce(add, (a/b*v for (a, b), v in zip(p, vertices)))
            nodes.append(real_p)
            weights.append(dim_factor*w)

        Quadrature.__init__(self, np.array(nodes).T, np.array(weights))
Esempio n. 7
0
    def __init__(self, alpha, beta, N, backend=None):  # noqa
        # default backend
        if backend is None:
            backend = 'builtin'
        if backend == 'builtin':
            x, w = self.compute_weights_and_nodes(N, alpha, beta)
        elif backend == 'scipy':
            from scipy.special.orthogonal import roots_jacobi
            x, w = roots_jacobi(N + 1, alpha, beta)
        else:
            raise NotImplementedError("Unsupported backend: %s" % backend)

        self.exact_to = 2*N + 1
        Quadrature.__init__(self, x, w)
Esempio n. 8
0
    def __init__(self, alpha, beta, N, backend=None):  # noqa: N803
        r"""
        :arg backend: Either ``"builtin"`` or ``"scipy"``. When the
            ``"builtin"`` backend is in use, there is an additional
            requirement that :math:`\alpha + \beta \ne -1`, with the exception
            of the Chebyshev quadrature :math:`\alpha = \beta = -1/2`. The
            ``"scipy"`` backend has no such restriction.
        """
        if backend is None:
            backend = "builtin"

        if backend == "builtin":
            x, w = self.compute_weights_and_nodes(N, alpha, beta)
        elif backend == "scipy":
            from scipy.special.orthogonal import roots_jacobi
            x, w = roots_jacobi(N + 1, alpha, beta)
        else:
            raise NotImplementedError("Unsupported backend: %s" % backend)

        self.exact_to = 2 * N + 1
        Quadrature.__init__(self, x, w)
Esempio n. 9
0
 def __init__(self, alpha, beta, N):  # noqa
     x, w = self.compute_weights_and_nodes(N, alpha, beta)
     Quadrature.__init__(self, x, w)
Esempio n. 10
0
 def __init__(self, N):  # noqa
     x, w = _fejer(N, 'cc')
     self.exact_to = N
     Quadrature.__init__(self, x, w)