diff --git a/FIAT/__init__.py b/FIAT/__init__.py index f5c1677ce..2c87b6426 100644 --- a/FIAT/__init__.py +++ b/FIAT/__init__.py @@ -21,6 +21,7 @@ from FIAT.lagrange import Lagrange from FIAT.gauss_lobatto_legendre import GaussLobattoLegendre from FIAT.gauss_legendre import GaussLegendre +from FIAT.gauss_radau import GaussRadau from FIAT.morley import Morley from FIAT.nedelec import Nedelec from FIAT.nedelec_second_kind import NedelecSecondKind @@ -65,6 +66,7 @@ "Lagrange": Lagrange, "Gauss-Lobatto-Legendre": GaussLobattoLegendre, "Gauss-Legendre": GaussLegendre, + "Gauss-Radau": GaussRadau, "Morley": Morley, "Nedelec 1st kind H(curl)": Nedelec, "Nedelec 2nd kind H(curl)": NedelecSecondKind, diff --git a/FIAT/gauss_radau.py b/FIAT/gauss_radau.py new file mode 100644 index 000000000..89191756d --- /dev/null +++ b/FIAT/gauss_radau.py @@ -0,0 +1,36 @@ +# Copyright (C) 2015 Imperial College London and others. +# +# This file is part of FIAT (https://www.fenicsproject.org) +# +# SPDX-License-Identifier: LGPL-3.0-or-later +# +# Written by Robert C. Kirby (robert_kirby@baylor.edu), 2020 + + +from FIAT import finite_element, polynomial_set, dual_set, functional, quadrature +from FIAT.reference_element import LINE + + +class GaussRadauDualSet(dual_set.DualSet): + """The dual basis for 1D discontinuous elements with nodes at the + Gauss-Radau points.""" + def __init__(self, ref_el, degree, right=True): + # Do DG connectivity because it's bonkers to do one-sided assembly even + # though we have an endpoint in the point set! + entity_ids = {0: {0: [], 1: []}, + 1: {0: list(range(0, degree+1))}} + lr = quadrature.RadauQuadratureLineRule(ref_el, degree+1, right) + nodes = [functional.PointEvaluation(ref_el, x) for x in lr.pts] + + super(GaussRadauDualSet, self).__init__(nodes, ref_el, entity_ids) + + +class GaussRadau(finite_element.CiarletElement): + """1D discontinuous element with nodes at the Gauss-Radau points.""" + def __init__(self, ref_el, degree): + if ref_el.shape != LINE: + raise ValueError("Gauss-Radau elements are only defined in one dimension.") + poly_set = polynomial_set.ONPolynomialSet(ref_el, degree) + dual = GaussRadauDualSet(ref_el, degree) + formdegree = ref_el.get_spatial_dimension() # n-form + super(GaussRadau, self).__init__(poly_set, dual, degree, formdegree) diff --git a/FIAT/quadrature.py b/FIAT/quadrature.py index 46ea01cb5..a93b73c1a 100644 --- a/FIAT/quadrature.py +++ b/FIAT/quadrature.py @@ -123,6 +123,63 @@ def __init__(self, ref_el, m): QuadratureRule.__init__(self, ref_el, xs, ws) +class RadauQuadratureLineRule(QuadratureRule): + """Produce the Gauss--Radau quadrature rules on the interval using + an adaptation of Winkel's Matlab code. + + The quadrature rule uses m points for a degree of precision of 2m-1. + """ + def __init__(self, ref_el, m, right=True): + assert m >= 1 + N = m - 1 + # Use Chebyshev-Gauss-Radau nodes as initial guess for LGR nodes + x = -numpy.cos(2 * numpy.pi * numpy.linspace(0, N, m) / (2 * N + 1)) + + P = numpy.zeros((N + 1, N + 2)) + + xold = 2 + + free = numpy.arange(1, N + 1, dtype='int') + + while numpy.max(numpy.abs(x - xold)) > 5e-16: + xold = x.copy() + + P[0, :] = (-1) ** numpy.arange(0, N + 2) + P[free, 0] = 1 + P[free, 1] = x[free] + + for k in range(2, N + 2): + P[free, k] = ((2 * k - 1) * x[free] * P[free, k - 1] - (k - 1) * P[free, k - 2]) / k + + x[free] = xold[free] - ((1 - xold[free]) / (N + 1)) * (P[free, N] + P[free, N + 1]) / (P[free, N] - P[free, N + 1]) + + # The Legendre-Gauss-Radau Vandermonde + P = P[:, :-1] + # Compute the weights + w = numpy.zeros(N + 1) + w[0] = 2 / (N + 1) ** 2 + w[free] = (1 - x[free])/((N + 1) * P[free, -1])**2 + + if right: + x = numpy.flip(-x) + w = numpy.flip(w) + + xs_ref = x + ws_ref = w + + A, b = reference_element.make_affine_mapping(((-1.,), (1.)), + ref_el.get_vertices()) + + mapping = lambda x: numpy.dot(A, x) + b + + scale = numpy.linalg.det(A) + + xs = tuple([tuple(mapping(x_ref)[0]) for x_ref in xs_ref]) + ws = tuple([scale * w for w in ws_ref]) + + QuadratureRule.__init__(self, ref_el, xs, ws) + + class CollapsedQuadratureTriangleRule(QuadratureRule): """Implements the collapsed quadrature rules defined in Karniadakis & Sherwin by mapping products of Gauss-Jacobi rules diff --git a/test/unit/test_gauss_radau.py b/test/unit/test_gauss_radau.py new file mode 100644 index 000000000..78d0ea8ed --- /dev/null +++ b/test/unit/test_gauss_radau.py @@ -0,0 +1,48 @@ +# Copyright (C) 2016 Imperial College London and others +# +# This file is part of FIAT. +# +# FIAT is free software: you can redistribute it and/or modify +# it under the terms of the GNU Lesser General Public License as published by +# the Free Software Foundation, either version 3 of the License, or +# (at your option) any later version. +# +# FIAT is distributed in the hope that it will be useful, +# but WITHOUT ANY WARRANTY; without even the implied warranty of +# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +# GNU Lesser General Public License for more details. +# +# You should have received a copy of the GNU Lesser General Public License +# along with FIAT. If not, see . +# +# Authors: +# +# Robert Kirby, based on work of David A. Ham +# + +import pytest +import numpy as np + + +@pytest.mark.parametrize("degree", range(1, 7)) +def test_gll_basis_values(degree): + """Ensure that integrating a simple monomial produces the expected results.""" + from FIAT import ufc_simplex, GaussRadau, make_quadrature + + s = ufc_simplex(1) + q = make_quadrature(s, degree + 1) + + fe = GaussRadau(s, degree) + tab = fe.tabulate(0, q.pts)[(0,)] + + for test_degree in range(degree + 1): + coefs = [n(lambda x: x[0]**test_degree) for n in fe.dual.nodes] + integral = np.dot(coefs, np.dot(tab, q.wts)) + reference = np.dot([x[0]**test_degree + for x in q.pts], q.wts) + assert np.allclose(integral, reference, rtol=1e-14) + + +if __name__ == '__main__': + import os + pytest.main(os.path.abspath(__file__)) diff --git a/test/unit/test_quadrature.py b/test/unit/test_quadrature.py index e77d8fa51..8b90d19ba 100644 --- a/test/unit/test_quadrature.py +++ b/test/unit/test_quadrature.py @@ -193,6 +193,18 @@ def test_gauss_lobatto_legendre_quadrature(interval, points, degree): assert numpy.round(q.integrate(lambda x: x[0]**degree) - 1./(degree+1), 14) == 0. +@pytest.mark.parametrize(("points, degree"), tuple((p, d) + for p in range(2, 10) + for d in range(2*p - 1))) +def test_radau_legendre_quadrature(interval, points, degree): + """Check that the quadrature rules correctly integrate all the right + polynomial degrees.""" + + q = FIAT.quadrature.RadauQuadratureLineRule(interval, points) + + assert numpy.round(q.integrate(lambda x: x[0]**degree) - 1./(degree+1), 14) == 0. + + @pytest.mark.parametrize(("points, degree"), tuple((p, d) for p in range(2, 10) for d in range(2*p)))