From 9fca26188a4533ff86e006a8fe61b8edbb73ed8b Mon Sep 17 00:00:00 2001 From: Pablo Brubeck Date: Wed, 8 Feb 2023 16:50:17 +0000 Subject: [PATCH] Add hierarchical and FDM elements (#34) * Add two FDM discontinuous elements * add hierarchical high-order elements * add FIAT.Legendre and FIAT.IntegratedLegendre * implement barycentric interpolation as a PolynomialSet * use sympy to simplify tabulation on sympy objects --- FIAT/__init__.py | 3 +- FIAT/barycentric_interpolation.py | 113 +++++++++++---- FIAT/fdm_element.py | 232 ++++++++++++++++++++---------- FIAT/gauss_legendre.py | 31 ++-- FIAT/gauss_lobatto_legendre.py | 31 ++-- FIAT/hierarchical.py | 85 +++++++++++ test/unit/test_fdm.py | 63 +++++--- test/unit/test_hierarchical.py | 73 ++++++++++ 8 files changed, 460 insertions(+), 171 deletions(-) create mode 100644 FIAT/hierarchical.py create mode 100644 test/unit/test_hierarchical.py diff --git a/FIAT/__init__.py b/FIAT/__init__.py index 38c15614..534a3216 100644 --- a/FIAT/__init__.py +++ b/FIAT/__init__.py @@ -47,7 +47,8 @@ from FIAT.restricted import RestrictedElement # noqa: F401 from FIAT.quadrature_element import QuadratureElement # noqa: F401 from FIAT.kong_mulder_veldhuizen import KongMulderVeldhuizen # noqa: F401 -from FIAT.fdm_element import FDMLagrange, FDMHermite # noqa: F401 +from FIAT.hierarchical import Legendre, IntegratedLegendre # noqa: F401 +from FIAT.fdm_element import FDMLagrange, FDMDiscontinuousLagrange, FDMQuadrature, FDMBrokenH1, FDMBrokenL2, FDMHermite # noqa: F401 # Important functionality from FIAT.quadrature import make_quadrature # noqa: F401 diff --git a/FIAT/barycentric_interpolation.py b/FIAT/barycentric_interpolation.py index 42d9fc5b..64c6913e 100644 --- a/FIAT/barycentric_interpolation.py +++ b/FIAT/barycentric_interpolation.py @@ -7,39 +7,92 @@ # Written by Pablo D. Brubeck (brubeck@protonmail.com), 2021 import numpy +from FIAT import reference_element, expansions, polynomial_set +from FIAT.functional import index_iterator -def barycentric_interpolation(xsrc, xdst, order=0): - """Return tabulations of a 1D Lagrange nodal basis via the second barycentric interpolation formula +def make_dmat(x): + """Returns Lagrange differentiation matrix and barycentric weights + associated with x[j].""" + dmat = numpy.add.outer(-x, x) + numpy.fill_diagonal(dmat, 1.0) + wts = numpy.prod(dmat, axis=0) + numpy.reciprocal(wts, out=wts) + numpy.divide(numpy.divide.outer(wts, wts), dmat, out=dmat) + numpy.fill_diagonal(dmat, dmat.diagonal() - numpy.sum(dmat, axis=0)) + return dmat, wts - See Berrut and Trefethen (2004) https://doi.org/10.1137/S0036144502417715 Eq. (4.2) & (9.4) - :arg xsrc: a :class:`numpy.array` with the nodes defining the Lagrange polynomial basis - :arg xdst: a :class:`numpy.array` with the interpolation points - :arg order: the integer order of differentiation - :returns: dict of tabulations up to the given order (in the same format as :meth:`~.CiarletElement.tabulate`) +class LagrangeLineExpansionSet(expansions.LineExpansionSet): + """Evaluates a Lagrange basis on a line reference element + via the second barycentric interpolation formula. See Berrut and Trefethen (2004) + https://doi.org/10.1137/S0036144502417715 Eq. (4.2) & (9.4) """ - # w = barycentric weights - # D = spectral differentiation matrix (D.T : u(xsrc) -> u'(xsrc)) - # I = barycentric interpolation matrix (I.T : u(xsrc) -> u(xdst)) - - D = numpy.add.outer(-xsrc, xsrc) - numpy.fill_diagonal(D, 1.0E0) - w = 1.0E0 / numpy.prod(D, axis=0) - D = numpy.divide.outer(w, w) / D - numpy.fill_diagonal(D, numpy.diag(D) - numpy.sum(D, axis=0)) - - I = numpy.add.outer(-xsrc, xdst) - idx = numpy.argwhere(numpy.isclose(I, 0.0E0, 0.0E0)) - I[idx[:, 0], idx[:, 1]] = 1.0E0 - I = 1.0E0 / I - I *= w[:, None] - I[:, idx[:, 1]] = 0.0E0 - I[idx[:, 0], idx[:, 1]] = 1.0E0 - I = (1.0E0 / numpy.sum(I, axis=0)) * I - - tabulation = {(0,): I} - for k in range(order): - tabulation[(k+1,)] = numpy.dot(D, tabulation[(k,)]) - return tabulation + def __init__(self, ref_el, pts): + self.nodes = numpy.array(pts).flatten() + self.dmat, self.weights = make_dmat(self.nodes) + super(LagrangeLineExpansionSet, self).__init__(ref_el) + + def get_num_members(self, n): + return len(self.nodes) + + def tabulate(self, n, pts): + assert n == len(self.nodes)-1 + results = numpy.add.outer(-self.nodes, numpy.array(pts).flatten()) + with numpy.errstate(divide='ignore', invalid='ignore'): + numpy.reciprocal(results, out=results) + numpy.multiply(results, self.weights[:, None], out=results) + numpy.multiply(1.0 / numpy.sum(results, axis=0), results, out=results) + + results[results != results] = 1.0 + if results.dtype == object: + from sympy import simplify + results = numpy.array(list(map(simplify, results))) + return results + + def tabulate_derivatives(self, n, pts): + return numpy.dot(self.dmat, self.tabulate(n, pts)) + + +class LagrangePolynomialSet(polynomial_set.PolynomialSet): + + def __init__(self, ref_el, pts, shape=tuple()): + degree = len(pts) - 1 + if shape == tuple(): + num_components = 1 + else: + flat_shape = numpy.ravel(shape) + num_components = numpy.prod(flat_shape) + num_exp_functions = expansions.polynomial_dimension(ref_el, degree) + num_members = num_components * num_exp_functions + embedded_degree = degree + expansion_set = get_expansion_set(ref_el, pts) + + # set up coefficients + if shape == tuple(): + coeffs = numpy.eye(num_members) + else: + coeffs_shape = tuple([num_members] + list(shape) + [num_exp_functions]) + coeffs = numpy.zeros(coeffs_shape, "d") + # use functional's index_iterator function + cur_bf = 0 + for idx in index_iterator(shape): + n = expansions.polynomial_dimension(ref_el, embedded_degree) + for exp_bf in range(n): + cur_idx = tuple([cur_bf] + list(idx) + [exp_bf]) + coeffs[cur_idx] = 1.0 + cur_bf += 1 + + dmats = [numpy.transpose(expansion_set.dmat)] + super(LagrangePolynomialSet, self).__init__(ref_el, degree, embedded_degree, + expansion_set, coeffs, dmats) + + +def get_expansion_set(ref_el, pts): + """Returns an ExpansionSet instance appopriate for the given + reference element.""" + if ref_el.get_shape() == reference_element.LINE: + return LagrangeLineExpansionSet(ref_el, pts) + else: + raise Exception("Unknown reference element type.") diff --git a/FIAT/fdm_element.py b/FIAT/fdm_element.py index b0485c6d..a045c551 100644 --- a/FIAT/fdm_element.py +++ b/FIAT/fdm_element.py @@ -6,13 +6,14 @@ # # Written by Pablo D. Brubeck (brubeck@protonmail.com), 2021 +import abc import numpy -from FIAT import finite_element, polynomial_set, dual_set, functional, quadrature +from FIAT import finite_element, dual_set, functional, quadrature from FIAT.reference_element import LINE -from FIAT.barycentric_interpolation import barycentric_interpolation from FIAT.lagrange import make_entity_permutations -from FIAT.gauss_lobatto_legendre import GaussLobattoLegendre +from FIAT.hierarchical import IntegratedLegendre +from FIAT.barycentric_interpolation import LagrangePolynomialSet from FIAT.P0 import P0Dual @@ -20,108 +21,189 @@ def sym_eig(A, B): """ A numpy-only implementation of `scipy.linalg.eigh` """ - L = numpy.linalg.cholesky(B) - Linv = numpy.linalg.inv(L) + Linv = numpy.linalg.inv(numpy.linalg.cholesky(B)) C = numpy.dot(Linv, numpy.dot(A, Linv.T)) - Z, W = numpy.linalg.eigh(C) - V = numpy.dot(Linv.T, W) + Z, V = numpy.linalg.eigh(C, "U") + V = numpy.dot(Linv.T, V) + return Z, V + + +def tridiag_eig(A, B): + """ + Same as sym_eig, but assumes that A is already diagonal and B tri-diagonal + """ + a = numpy.reciprocal(A.diagonal()) + numpy.sqrt(a, out=a) + C = numpy.multiply(a, B) + numpy.multiply(C, a[:, None], out=C) + Z, V = numpy.linalg.eigh(C, "U") + numpy.reciprocal(Z, out=Z) + numpy.multiply(numpy.sqrt(Z), V, out=V) + numpy.multiply(V, a[:, None], out=V) return Z, V class FDMDual(dual_set.DualSet): - """The dual basis for 1D CG elements with FDM shape functions.""" - def __init__(self, ref_el, degree, order=1): - bc_nodes = [] - for x in ref_el.get_vertices(): - bc_nodes.append([functional.PointEvaluation(ref_el, x), - *[functional.PointDerivative(ref_el, x, [alpha]) for alpha in range(1, order)]]) - bc_nodes[1].reverse() - k = len(bc_nodes[0]) - idof = slice(k, -k) - bdof = list(range(-k, k)) - bdof = bdof[k:] + bdof[:k] - - # Define the generalized eigenproblem on a GLL element - gll = GaussLobattoLegendre(ref_el, degree) - xhat = numpy.array([list(x.get_point_dict().keys())[0][0] for x in gll.dual_basis()]) + """The dual basis for 1D elements with FDM shape functions.""" + def __init__(self, ref_el, degree, bc_order=1, formdegree=0, orthogonalize=False): + # Define the generalized eigenproblem on a reference element + embedded_degree = degree + formdegree + embedded = IntegratedLegendre(ref_el, embedded_degree) + edim = embedded.space_dimension() + self.embedded = embedded + + vertices = ref_el.get_vertices() + rule = quadrature.GaussLegendreQuadratureLineRule(ref_el, edim) + self.rule = rule + + solve_eig = sym_eig + if bc_order == 1: + solve_eig = tridiag_eig # Tabulate the BC nodes - constraints = gll.tabulate(order-1, ref_el.get_vertices()) - C = numpy.column_stack(list(constraints.values())) - perm = list(range(len(bdof))) - perm = perm[::2] + perm[-1::-2] - C = C[:, perm].T + if bc_order == 0: + C = numpy.empty((0, edim), "d") + else: + constraints = embedded.tabulate(bc_order-1, vertices) + C = numpy.transpose(numpy.column_stack(list(constraints.values()))) + bdof = slice(None, C.shape[0]) + idof = slice(C.shape[0], None) # Tabulate the basis that splits the DOFs into interior and bcs - E = numpy.eye(gll.space_dimension()) + E = numpy.eye(edim) E[bdof, idof] = -C[:, idof] - E[bdof, :] = numpy.dot(numpy.linalg.inv(C[:, bdof]), E[bdof, :]) + E[bdof, :] = numpy.linalg.solve(C[:, bdof], E[bdof, :]) # Assemble the constrained Galerkin matrices on the reference cell - rule = quadrature.GaussLegendreQuadratureLineRule(ref_el, degree+1) - phi = gll.tabulate(order, rule.get_points()) - E0 = numpy.dot(phi[(0, )].T, E) - Ek = numpy.dot(phi[(order, )].T, E) - B = numpy.dot(numpy.multiply(E0.T, rule.get_weights()), E0) - A = numpy.dot(numpy.multiply(Ek.T, rule.get_weights()), Ek) + k = max(1, bc_order) + phi = embedded.tabulate(k, rule.get_points()) + wts = rule.get_weights() + E0 = numpy.dot(E.T, phi[(0, )]) + Ek = numpy.dot(E.T, phi[(k, )]) + B = numpy.dot(numpy.multiply(E0, wts), E0.T) + A = numpy.dot(numpy.multiply(Ek, wts), Ek.T) # Eigenfunctions in the constrained basis S = numpy.eye(A.shape[0]) - if S.shape[0] > len(bdof): - _, Sii = sym_eig(A[idof, idof], B[idof, idof]) + lam = numpy.ones((A.shape[0],)) + if S.shape[0] > C.shape[0]: + lam[idof], Sii = solve_eig(A[idof, idof], B[idof, idof]) S[idof, idof] = Sii S[idof, bdof] = numpy.dot(Sii, numpy.dot(Sii.T, -B[idof, bdof])) - # Eigenfunctions in the Lagrange basis - S = numpy.dot(E, S) - self.gll_points = xhat - self.gll_tabulation = S.T + if orthogonalize: + Abb = numpy.dot(S[:, bdof].T, numpy.dot(A, S[:, bdof])) + Bbb = numpy.dot(S[:, bdof].T, numpy.dot(B, S[:, bdof])) + _, Qbb = sym_eig(Abb, Bbb) + S[:, bdof] = numpy.dot(S[:, bdof], Qbb) - # Interpolate eigenfunctions onto the quadrature points - basis = numpy.dot(S.T, phi[(0, )]) - nodes = bc_nodes[0] + [functional.IntegralMoment(ref_el, rule, f) for f in basis[idof]] + bc_nodes[1] + if formdegree == 0: + # Tabulate eigenbasis + basis = numpy.dot(S.T, E0) + else: + # Tabulate the derivative of the eigenbasis and normalize + if bc_order == 0: + idof = lam > 1.0E-12 + lam[~idof] = 1.0E0 + numpy.reciprocal(lam, out=lam) + numpy.sqrt(lam, out=lam) + numpy.multiply(S, lam, out=S) + basis = numpy.dot(S.T, Ek) - entity_ids = {0: {0: [0], 1: [degree]}, - 1: {0: list(range(1, degree))}} - entity_permutations = {} - entity_permutations[0] = {0: {0: [0]}, 1: {0: [0]}} - entity_permutations[1] = {0: make_entity_permutations(1, degree - 1)} + bc_nodes = [] + if formdegree == 0: + if orthogonalize: + idof = slice(None) + elif bc_order > 0: + bc_nodes += [functional.PointEvaluation(ref_el, x) for x in vertices] + for alpha in range(1, bc_order): + bc_nodes += [functional.PointDerivative(ref_el, x, [alpha]) for x in vertices] + + elif bc_order > 0: + basis[bdof, :] = numpy.sqrt(1.0E0 / ref_el.volume()) + idof = slice(formdegree, None) + + nodes = bc_nodes + [functional.IntegralMoment(ref_el, rule, f) for f in basis[idof]] + + if len(bc_nodes) > 0: + entity_ids = {0: {0: [0], 1: [1]}, + 1: {0: list(range(2, degree+1))}} + entity_permutations = {} + entity_permutations[0] = {0: {0: [0]}, 1: {0: [0]}} + entity_permutations[1] = {0: make_entity_permutations(1, degree - 1)} + else: + entity_ids = {0: {0: [], 1: []}, + 1: {0: list(range(0, degree+1))}} + entity_permutations = {} + entity_permutations[0] = {0: {0: []}, 1: {0: []}} + entity_permutations[1] = {0: make_entity_permutations(1, degree + 1)} super(FDMDual, self).__init__(nodes, ref_el, entity_ids, entity_permutations) -class FDMLagrange(finite_element.CiarletElement): - """1D CG element with interior shape functions that diagonalize the Laplacian.""" - _order = 1 +class FDMFiniteElement(finite_element.CiarletElement): + """1D element that diagonalizes bilinear forms with BCs.""" + + _orthogonalize = False + + @property + @abc.abstractmethod + def _bc_order(self): + pass + + @property + @abc.abstractmethod + def _formdegree(self): + pass def __init__(self, ref_el, degree): if ref_el.shape != LINE: - raise ValueError("FDM elements are only defined in one dimension.") - poly_set = polynomial_set.ONPolynomialSet(ref_el, degree) + raise ValueError("%s is only defined in one dimension." % type(self)) if degree == 0: dual = P0Dual(ref_el) else: - dual = FDMDual(ref_el, degree, order=self._order) - formdegree = 0 - super(FDMLagrange, self).__init__(poly_set, dual, degree, formdegree) - - def tabulate(self, order, points, entity=None): - # This overrides the default with a more numerically stable algorithm - if hasattr(self.dual, "gll_points"): - if entity is None: - entity = (self.ref_el.get_dimension(), 0) - - entity_dim, entity_id = entity - transform = self.ref_el.get_entity_transform(entity_dim, entity_id) - xsrc = self.dual.gll_points - xdst = numpy.array(list(map(transform, points))).flatten() - tabulation = barycentric_interpolation(xsrc, xdst, order=order) - for key in tabulation: - tabulation[key] = numpy.dot(self.dual.gll_tabulation, tabulation[key]) - return tabulation + dual = FDMDual(ref_el, degree, bc_order=self._bc_order, + formdegree=self._formdegree, orthogonalize=self._orthogonalize) + + if self._formdegree == 0: + poly_set = dual.embedded.poly_set else: - return super(FDMLagrange, self).tabulate(order, points, entity) + lr = quadrature.GaussLegendreQuadratureLineRule(ref_el, degree+1) + poly_set = LagrangePolynomialSet(ref_el, lr.get_points()) + super(FDMFiniteElement, self).__init__(poly_set, dual, degree, self._formdegree) + + +class FDMLagrange(FDMFiniteElement): + """1D CG element with interior shape functions that diagonalize the Laplacian.""" + _bc_order = 1 + _formdegree = 0 + + +class FDMDiscontinuousLagrange(FDMFiniteElement): + """1D DG element with derivatives of interior CG FDM shape functions.""" + _bc_order = 1 + _formdegree = 1 + + +class FDMQuadrature(FDMFiniteElement): + """1D DG element with interior CG FDM shape functions and orthogonalized vertex modes.""" + _bc_order = 1 + _formdegree = 0 + _orthogonalize = True + + +class FDMBrokenH1(FDMFiniteElement): + """1D DG element with shape functions that diagonalize the Laplacian.""" + _bc_order = 0 + _formdegree = 0 + + +class FDMBrokenL2(FDMFiniteElement): + """1D DG element with the derivates of DG FDM shape functions.""" + _bc_order = 0 + _formdegree = 1 -class FDMHermite(FDMLagrange): +class FDMHermite(FDMFiniteElement): """1D CG element with interior shape functions that diagonalize the biharmonic operator.""" - _order = 2 + _bc_order = 2 + _formdegree = 0 diff --git a/FIAT/gauss_legendre.py b/FIAT/gauss_legendre.py index 4cf41ce6..a6919c87 100644 --- a/FIAT/gauss_legendre.py +++ b/FIAT/gauss_legendre.py @@ -8,12 +8,10 @@ # # Modified by Pablo D. Brubeck (brubeck@protonmail.com), 2021 -import numpy - -from FIAT import finite_element, polynomial_set, dual_set, functional, quadrature +from FIAT import finite_element, dual_set, functional, quadrature from FIAT.reference_element import LINE -from FIAT.barycentric_interpolation import barycentric_interpolation from FIAT.lagrange import make_entity_permutations +from FIAT.barycentric_interpolation import LagrangePolynomialSet class GaussLegendreDualSet(dual_set.DualSet): @@ -36,25 +34,12 @@ class GaussLegendre(finite_element.CiarletElement): def __init__(self, ref_el, degree): if ref_el.shape != LINE: raise ValueError("Gauss-Legendre elements are only defined in one dimension.") - poly_set = polynomial_set.ONPolynomialSet(ref_el, degree) dual = GaussLegendreDualSet(ref_el, degree) + points = [] + for node in dual.nodes: + # Assert singleton point for each node. + pt, = node.get_point_dict().keys() + points.append(pt) + poly_set = LagrangePolynomialSet(ref_el, points) formdegree = ref_el.get_spatial_dimension() # n-form super(GaussLegendre, self).__init__(poly_set, dual, degree, formdegree) - - def tabulate(self, order, points, entity=None): - # This overrides the default with a more numerically stable algorithm - - if entity is None: - entity = (self.ref_el.get_dimension(), 0) - - entity_dim, entity_id = entity - transform = self.ref_el.get_entity_transform(entity_dim, entity_id) - - xsrc = [] - for node in self.dual.nodes: - # Assert singleton point for each node. - (pt,), = node.get_point_dict().keys() - xsrc.append(pt) - xsrc = numpy.asarray(xsrc) - xdst = numpy.array(list(map(transform, points))).flatten() - return barycentric_interpolation(xsrc, xdst, order=order) diff --git a/FIAT/gauss_lobatto_legendre.py b/FIAT/gauss_lobatto_legendre.py index d6b8b333..1cec6031 100644 --- a/FIAT/gauss_lobatto_legendre.py +++ b/FIAT/gauss_lobatto_legendre.py @@ -8,12 +8,10 @@ # # Modified by Pablo D. Brubeck (brubeck@protonmail.com), 2021 -import numpy - -from FIAT import finite_element, polynomial_set, dual_set, functional, quadrature +from FIAT import finite_element, dual_set, functional, quadrature from FIAT.reference_element import LINE -from FIAT.barycentric_interpolation import barycentric_interpolation from FIAT.lagrange import make_entity_permutations +from FIAT.barycentric_interpolation import LagrangePolynomialSet class GaussLobattoLegendreDualSet(dual_set.DualSet): @@ -36,25 +34,12 @@ class GaussLobattoLegendre(finite_element.CiarletElement): def __init__(self, ref_el, degree): if ref_el.shape != LINE: raise ValueError("Gauss-Lobatto-Legendre elements are only defined in one dimension.") - poly_set = polynomial_set.ONPolynomialSet(ref_el, degree) dual = GaussLobattoLegendreDualSet(ref_el, degree) + points = [] + for node in dual.nodes: + # Assert singleton point for each node. + pt, = node.get_point_dict().keys() + points.append(pt) + poly_set = LagrangePolynomialSet(ref_el, points) formdegree = 0 # 0-form super(GaussLobattoLegendre, self).__init__(poly_set, dual, degree, formdegree) - - def tabulate(self, order, points, entity=None): - # This overrides the default with a more numerically stable algorithm - - if entity is None: - entity = (self.ref_el.get_dimension(), 0) - - entity_dim, entity_id = entity - transform = self.ref_el.get_entity_transform(entity_dim, entity_id) - - xsrc = [] - for node in self.dual.nodes: - # Assert singleton point for each node. - (pt,), = node.get_point_dict().keys() - xsrc.append(pt) - xsrc = numpy.asarray(xsrc) - xdst = numpy.array(list(map(transform, points))).flatten() - return barycentric_interpolation(xsrc, xdst, order=order) diff --git a/FIAT/hierarchical.py b/FIAT/hierarchical.py new file mode 100644 index 00000000..313fd5ae --- /dev/null +++ b/FIAT/hierarchical.py @@ -0,0 +1,85 @@ +# 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 Pablo D. Brubeck (brubeck@protonmail.com), 2022 + +import numpy + +from FIAT import (finite_element, reference_element, + dual_set, functional, quadrature, + jacobi, barycentric_interpolation) +from FIAT.lagrange import make_entity_permutations +from FIAT.barycentric_interpolation import LagrangePolynomialSet + + +class LegendreDual(dual_set.DualSet): + """The dual basis for Legendre elements.""" + def __init__(self, ref_el, degree, rule): + v1 = ref_el.get_vertices() + A, b = reference_element.make_affine_mapping(v1, [(-1.0,), (1.0,)]) + mapping = lambda x: numpy.dot(A, x) + b + xhat = numpy.array([mapping(pt) for pt in rule.pts]) + + basis = jacobi.eval_jacobi_batch(0, 0, degree, xhat) + nodes = [functional.IntegralMoment(ref_el, rule, f) for f in basis] + + entity_ids = {0: {0: [], 1: []}, + 1: {0: list(range(0, degree+1))}} + entity_permutations = {} + entity_permutations[0] = {0: {0: []}, 1: {0: []}} + entity_permutations[1] = {0: make_entity_permutations(1, degree + 1)} + super(LegendreDual, self).__init__(nodes, ref_el, entity_ids, entity_permutations) + + +class Legendre(finite_element.CiarletElement): + """1D discontinuous element with Legendre polynomials.""" + + def __init__(self, ref_el, degree): + if ref_el.shape != reference_element.LINE: + raise ValueError("%s is only defined in one dimension." % type(self)) + rule = quadrature.GaussLegendreQuadratureLineRule(ref_el, degree+1) + poly_set = LagrangePolynomialSet(ref_el, rule.get_points()) + dual = LegendreDual(ref_el, degree, rule) + formdegree = ref_el.get_spatial_dimension() # n-form + super(Legendre, self).__init__(poly_set, dual, degree, formdegree) + + +class IntegratedLegendreDual(dual_set.DualSet): + """The dual basis for integrated Legendre elements.""" + def __init__(self, ref_el, degree, rule): + v1 = ref_el.get_vertices() + A, b = reference_element.make_affine_mapping(v1, [(-1.0,), (1.0,)]) + mapping = lambda x: numpy.dot(A, x) + b + xhat = numpy.array([mapping(pt) for pt in rule.pts]) + + W = rule.get_weights() + D, _ = barycentric_interpolation.make_dmat(numpy.array(rule.pts).flatten()) + P = jacobi.eval_jacobi_batch(0, 0, degree-1, xhat) + basis = numpy.dot(numpy.multiply(P, W), numpy.multiply(D.T, 1.0/W)) + + nodes = [functional.PointEvaluation(ref_el, x) for x in v1] + nodes += [functional.IntegralMoment(ref_el, rule, f) for f in basis[2::2]] + nodes += [functional.IntegralMoment(ref_el, rule, f) for f in basis[1::2]] + + entity_ids = {0: {0: [0], 1: [1]}, + 1: {0: list(range(2, degree+1))}} + entity_permutations = {} + entity_permutations[0] = {0: {0: [0]}, 1: {0: [0]}} + entity_permutations[1] = {0: make_entity_permutations(1, degree - 1)} + super(IntegratedLegendreDual, self).__init__(nodes, ref_el, entity_ids, entity_permutations) + + +class IntegratedLegendre(finite_element.CiarletElement): + """1D continuous element with integrated Legendre polynomials.""" + + def __init__(self, ref_el, degree): + if ref_el.shape != reference_element.LINE: + raise ValueError("%s is only defined in one dimension." % type(self)) + rule = quadrature.GaussLegendreQuadratureLineRule(ref_el, degree+1) + poly_set = LagrangePolynomialSet(ref_el, rule.get_points()) + dual = IntegratedLegendreDual(ref_el, degree, rule) + formdegree = 0 # 0-form + super(IntegratedLegendre, self).__init__(poly_set, dual, degree, formdegree) diff --git a/test/unit/test_fdm.py b/test/unit/test_fdm.py index 591f6ed6..5286124b 100644 --- a/test/unit/test_fdm.py +++ b/test/unit/test_fdm.py @@ -23,15 +23,31 @@ import numpy as np -@pytest.mark.parametrize("degree", range(1, 7)) -def test_fdm_basis_values(degree): +def make_fdm_element(ref_el, family, degree): + from FIAT import FDMLagrange, FDMDiscontinuousLagrange, FDMBrokenH1, FDMBrokenL2, FDMQuadrature + if family == "CG": + return FDMLagrange(ref_el, degree) + elif family == "DG": + return FDMDiscontinuousLagrange(ref_el, degree) + elif family == "BrokenH1": + return FDMBrokenH1(ref_el, degree) + elif family == "BrokenL2": + return FDMBrokenL2(ref_el, degree) + elif family == "Quadrature": + return FDMQuadrature(ref_el, degree) + + +@pytest.mark.parametrize("family, degree", [(f, degree - 1 if f in {"DG", "BrokenL2"} else degree) + for f in ("CG", "DG", "BrokenH1", "BrokenL2", "Quadrature") + for degree in range(1, 7)]) +def test_fdm_basis_values(family, degree): """Ensure that integrating a simple monomial produces the expected results.""" - from FIAT import ufc_simplex, FDMLagrange, make_quadrature + from FIAT import ufc_simplex, make_quadrature s = ufc_simplex(1) q = make_quadrature(s, degree + 1) + fe = make_fdm_element(s, family, degree) - fe = FDMLagrange(s, degree) tab = fe.tabulate(0, q.pts)[(0,)] for test_degree in range(degree + 1): @@ -42,23 +58,32 @@ def test_fdm_basis_values(degree): assert np.allclose(integral, reference, rtol=1e-14) -@pytest.mark.parametrize("degree", range(1, 7)) -def test_sparsity(degree): - from FIAT import ufc_simplex, FDMLagrange, make_quadrature - cell = ufc_simplex(1) - fe = FDMLagrange(cell, degree) +@pytest.mark.parametrize("family, degree", [(f, degree - 1 if f in {"DG", "BrokenL2"} else degree) + for f in ("CG", "DG", "BrokenH1", "BrokenL2", "Quadrature") + for degree in range(1, 7)]) +def test_fdm_sparsity(family, degree): + from FIAT import ufc_simplex, make_quadrature + + s = ufc_simplex(1) + q = make_quadrature(s, degree+1) + fe = make_fdm_element(s, family, degree) + + if family == "CG": + expected = [degree + 3, 5 * degree - 1] + elif family == "DG": + expected = [degree + 1] + elif family == "BrokenH1": + expected = [degree + 1, degree] + elif family == "BrokenL2": + expected = [degree + 1] + elif family == "Quadrature": + expected = [degree + 1, 3 * degree - 1 - (degree == 1)] - rule = make_quadrature(cell, degree+1) - basis = fe.tabulate(1, rule.get_points()) - Jhat = basis[(0,)] - Dhat = basis[(1,)] - what = rule.get_weights() - Ahat = np.dot(np.multiply(Dhat, what), Dhat.T) - Bhat = np.dot(np.multiply(Jhat, what), Jhat.T) nnz = lambda A: A.size - np.sum(np.isclose(A, 0.0E0, rtol=1E-14)) - ndof = fe.space_dimension() - assert nnz(Ahat) == 5*ndof-6 - assert nnz(Bhat) == ndof+2 + moments = lambda v, u: np.dot(np.multiply(v, q.get_weights()), u.T) + tab = fe.tabulate(len(expected)-1, q.get_points()) + for k, ennz in enumerate(expected): + assert nnz(moments(tab[(k, )], tab[(k, )])) == ennz if __name__ == '__main__': diff --git a/test/unit/test_hierarchical.py b/test/unit/test_hierarchical.py new file mode 100644 index 00000000..c802cfd2 --- /dev/null +++ b/test/unit/test_hierarchical.py @@ -0,0 +1,73 @@ +# 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: +# +# Pablo Brubeck + +import pytest +import numpy as np + + +@pytest.mark.parametrize("family, degree", [(f, degree - 1 if f == "DG" else degree) + for f in ("CG", "DG") + for degree in range(1, 7)]) +def test_hierarchical_basis_values(family, degree): + """Ensure that integrating a simple monomial produces the expected results.""" + from FIAT import ufc_simplex, Legendre, IntegratedLegendre, make_quadrature + + s = ufc_simplex(1) + q = make_quadrature(s, degree + 1) + if family == "CG": + fe = IntegratedLegendre(s, degree) + else: + fe = Legendre(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) + + +@pytest.mark.parametrize("family, degree", [(f, degree - 1 if f == "DG" else degree) + for f in ("CG", "DG") + for degree in range(1, 7)]) +def test_hierarchical_sparsity(family, degree): + from FIAT import ufc_simplex, Legendre, IntegratedLegendre, make_quadrature + + s = ufc_simplex(1) + q = make_quadrature(s, degree+1) + if family == "CG": + fe = IntegratedLegendre(s, degree) + expected = [5 * min(degree, 3) + 3 * max(0, degree-3) - 1, degree + 3] + else: + fe = Legendre(s, degree) + expected = [degree + 1] + + nnz = lambda A: A.size - np.sum(np.isclose(A, 0.0E0, rtol=1E-14)) + moments = lambda v, u: np.dot(np.multiply(v, q.get_weights()), u.T) + tab = fe.tabulate(len(expected)-1, q.get_points()) + for k, ennz in enumerate(expected): + assert nnz(moments(tab[(k, )], tab[(k, )])) == ennz + + +if __name__ == '__main__': + import os + pytest.main(os.path.abspath(__file__))