xref: /libCEED/python/ceed_basis.py (revision 7a7b0fa3029ed576dc3d3fb37afaf50f3689b395)
15a21df39Svaleriabarra# Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
25a21df39Svaleriabarra# the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
35a21df39Svaleriabarra# reserved. See files LICENSE and NOTICE for details.
45a21df39Svaleriabarra#
55a21df39Svaleriabarra# This file is part of CEED, a collection of benchmarks, miniapps, software
65a21df39Svaleriabarra# libraries and APIs for efficient high-order finite element and spectral
75a21df39Svaleriabarra# element discretizations for exascale applications. For more information and
85a21df39Svaleriabarra# source code availability see http://github.com/ceed.
95a21df39Svaleriabarra#
105a21df39Svaleriabarra# The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
115a21df39Svaleriabarra# a collaborative effort of two U.S. Department of Energy organizations (Office
125a21df39Svaleriabarra# of Science and the National Nuclear Security Administration) responsible for
135a21df39Svaleriabarra# the planning and preparation of a capable exascale ecosystem, including
145a21df39Svaleriabarra# software, applications, hardware, advanced system engineering and early
155a21df39Svaleriabarra# testbed platforms, in support of the nation's exascale computing imperative.
165a21df39Svaleriabarra
175a21df39Svaleriabarrafrom _ceed_cffi import ffi, lib
185a21df39Svaleriabarraimport tempfile
195a21df39Svaleriabarraimport numpy as np
205a21df39Svaleriabarrafrom abc import ABC
215a21df39Svaleriabarrafrom .ceed_constants import TRANSPOSE, NOTRANSPOSE
225a21df39Svaleriabarra
235a21df39Svaleriabarra# ------------------------------------------------------------------------------
24*7a7b0fa3SJed Brown
25*7a7b0fa3SJed Brown
265a21df39Svaleriabarraclass Basis(ABC):
275a21df39Svaleriabarra    """Ceed Basis: finite element basis objects."""
285a21df39Svaleriabarra
295a21df39Svaleriabarra    # Attributes
305a21df39Svaleriabarra    _ceed = ffi.NULL
315a21df39Svaleriabarra    _pointer = ffi.NULL
325a21df39Svaleriabarra
335a21df39Svaleriabarra    # Representation
345a21df39Svaleriabarra    def __repr__(self):
355a21df39Svaleriabarra        return "<CeedBasis instance at " + hex(id(self)) + ">"
365a21df39Svaleriabarra
375a21df39Svaleriabarra    # String conversion for print() to stdout
385a21df39Svaleriabarra    def __str__(self):
395a21df39Svaleriabarra        """View a Basis via print()."""
405a21df39Svaleriabarra
415a21df39Svaleriabarra        # libCEED call
425a21df39Svaleriabarra        with tempfile.NamedTemporaryFile() as key_file:
435a21df39Svaleriabarra            with open(key_file.name, 'r+') as stream_file:
445a21df39Svaleriabarra                stream = ffi.cast("FILE *", stream_file)
455a21df39Svaleriabarra
465a21df39Svaleriabarra                lib.CeedBasisView(self._pointer[0], stream)
475a21df39Svaleriabarra
485a21df39Svaleriabarra                stream_file.seek(0)
495a21df39Svaleriabarra                out_string = stream_file.read()
505a21df39Svaleriabarra
515a21df39Svaleriabarra        return out_string
525a21df39Svaleriabarra
535a21df39Svaleriabarra    # Apply Basis
545a21df39Svaleriabarra    def apply(self, nelem, emode, u, v, tmode=NOTRANSPOSE):
555a21df39Svaleriabarra        """Apply basis evaluation from nodes to quadrature points or vice versa.
565a21df39Svaleriabarra
575a21df39Svaleriabarra           Args:
585a21df39Svaleriabarra             nelem: the number of elements to apply the basis evaluation to;
595a21df39Svaleriabarra                      the backend will specify the ordering in a
605a21df39Svaleriabarra                      BlockedElemRestriction
615a21df39Svaleriabarra             emode: basis evaluation mode
625a21df39Svaleriabarra             u: input vector
635a21df39Svaleriabarra             v: output vector
645a21df39Svaleriabarra             **tmode: CEED_NOTRANSPOSE to evaluate from nodes to quadrature
655a21df39Svaleriabarra                        points, CEED_TRANSPOSE to apply the transpose, mapping
665a21df39Svaleriabarra                        from quadrature points to nodes; default CEED_NOTRANSPOSE"""
675a21df39Svaleriabarra
685a21df39Svaleriabarra        # libCEED call
695a21df39Svaleriabarra        lib.CeedBasisApply(self._pointer[0], nelem, tmode, emode,
705a21df39Svaleriabarra                           u._pointer[0], v._pointer[0])
715a21df39Svaleriabarra
725a21df39Svaleriabarra    # Transpose a Basis
735a21df39Svaleriabarra    @property
745a21df39Svaleriabarra    def T(self):
755a21df39Svaleriabarra        """Transpose a Basis."""
765a21df39Svaleriabarra
775a21df39Svaleriabarra        return TransposeBasis(self)
785a21df39Svaleriabarra
795a21df39Svaleriabarra    # Transpose a Basis
805a21df39Svaleriabarra    @property
815a21df39Svaleriabarra    def transpose(self):
825a21df39Svaleriabarra        """Transpose a Basis."""
835a21df39Svaleriabarra
845a21df39Svaleriabarra        return TransposeBasis(self)
855a21df39Svaleriabarra
865a21df39Svaleriabarra    # Get number of nodes
875a21df39Svaleriabarra    def get_num_nodes(self):
885a21df39Svaleriabarra        """Get total number of nodes (in dim dimensions) of a Basis.
895a21df39Svaleriabarra
905a21df39Svaleriabarra           Returns:
915a21df39Svaleriabarra             num_nodes: total number of nodes"""
925a21df39Svaleriabarra
935a21df39Svaleriabarra        # Setup argument
945a21df39Svaleriabarra        p_pointer = ffi.new("CeedInt *")
955a21df39Svaleriabarra
965a21df39Svaleriabarra        # libCEED call
975a21df39Svaleriabarra        lib.CeedBasisGetNumNodes(self._pointer[0], p_pointer)
985a21df39Svaleriabarra
995a21df39Svaleriabarra        return p_pointer[0]
1005a21df39Svaleriabarra
1015a21df39Svaleriabarra    # Get number of quadrature points
1025a21df39Svaleriabarra    def get_num_quadrature_points(self):
1035a21df39Svaleriabarra        """Get total number of quadrature points (in dim dimensions) of a Basis.
1045a21df39Svaleriabarra
1055a21df39Svaleriabarra           Returns:
1065a21df39Svaleriabarra             num_qpts: total number of quadrature points"""
1075a21df39Svaleriabarra
1085a21df39Svaleriabarra        # Setup argument
1095a21df39Svaleriabarra        q_pointer = ffi.new("CeedInt *")
1105a21df39Svaleriabarra
1115a21df39Svaleriabarra        # libCEED call
1125a21df39Svaleriabarra        lib.CeedBasisGetNumQuadraturePoints(self._pointer[0], q_pointer)
1135a21df39Svaleriabarra
1145a21df39Svaleriabarra        return q_pointer[0]
1155a21df39Svaleriabarra
1165a21df39Svaleriabarra    # Gauss quadrature
1175a21df39Svaleriabarra    @staticmethod
1185a21df39Svaleriabarra    def gauss_quadrature(q):
1195a21df39Svaleriabarra        """Construct a Gauss-Legendre quadrature.
1205a21df39Svaleriabarra
1215a21df39Svaleriabarra           Args:
1225a21df39Svaleriabarra             Q: number of quadrature points (integrates polynomials of
1235a21df39Svaleriabarra                  degree 2*Q-1 exactly)
1245a21df39Svaleriabarra
1255a21df39Svaleriabarra           Returns:
1265a21df39Svaleriabarra             (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
1275a21df39Svaleriabarra                                    and array of length Q to hold the weights"""
1285a21df39Svaleriabarra
1295a21df39Svaleriabarra        # Setup arguments
1305a21df39Svaleriabarra        qref1d = np.empty(q, dtype="float64")
131fc2a3161SJed Brown        qweight1d = np.empty(q, dtype="float64")
1325a21df39Svaleriabarra
1335a21df39Svaleriabarra        qref1d_pointer = ffi.new("CeedScalar *")
134*7a7b0fa3SJed Brown        qref1d_pointer = ffi.cast(
135*7a7b0fa3SJed Brown            "CeedScalar *",
136*7a7b0fa3SJed Brown            qref1d.__array_interface__['data'][0])
1375a21df39Svaleriabarra
1385a21df39Svaleriabarra        qweight1d_pointer = ffi.new("CeedScalar *")
139*7a7b0fa3SJed Brown        qweight1d_pointer = ffi.cast(
140*7a7b0fa3SJed Brown            "CeedScalar *",
141*7a7b0fa3SJed Brown            qweight1d.__array_interface__['data'][0])
1425a21df39Svaleriabarra
1435a21df39Svaleriabarra        # libCEED call
1445a21df39Svaleriabarra        lib.CeedGaussQuadrature(q, qref1d_pointer, qweight1d_pointer)
1455a21df39Svaleriabarra
1465a21df39Svaleriabarra        return qref1d, qweight1d
1475a21df39Svaleriabarra
1485a21df39Svaleriabarra    # Lobatto quadrature
1495a21df39Svaleriabarra    @staticmethod
1505a21df39Svaleriabarra    def lobatto_quadrature(q):
1515a21df39Svaleriabarra        """Construct a Gauss-Legendre-Lobatto quadrature.
1525a21df39Svaleriabarra
1535a21df39Svaleriabarra           Args:
1545a21df39Svaleriabarra             q: number of quadrature points (integrates polynomials of
1555a21df39Svaleriabarra                  degree 2*Q-3 exactly)
1565a21df39Svaleriabarra
1575a21df39Svaleriabarra           Returns:
1585a21df39Svaleriabarra             (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
1595a21df39Svaleriabarra                                    and array of length Q to hold the weights"""
1605a21df39Svaleriabarra
1615a21df39Svaleriabarra        # Setup arguments
1625a21df39Svaleriabarra        qref1d = np.empty(q, dtype="float64")
1635a21df39Svaleriabarra        qref1d_pointer = ffi.new("CeedScalar *")
164*7a7b0fa3SJed Brown        qref1d_pointer = ffi.cast(
165*7a7b0fa3SJed Brown            "CeedScalar *",
166*7a7b0fa3SJed Brown            qref1d.__array_interface__['data'][0])
1675a21df39Svaleriabarra
168fc2a3161SJed Brown        qweight1d = np.empty(q, dtype="float64")
1695a21df39Svaleriabarra        qweight1d_pointer = ffi.new("CeedScalar *")
170*7a7b0fa3SJed Brown        qweight1d_pointer = ffi.cast(
171*7a7b0fa3SJed Brown            "CeedScalar *",
172*7a7b0fa3SJed Brown            qweight1d.__array_interface__['data'][0])
1735a21df39Svaleriabarra
1745a21df39Svaleriabarra        # libCEED call
1755a21df39Svaleriabarra        lib.CeedLobattoQuadrature(q, qref1d_pointer, qweight1d_pointer)
1765a21df39Svaleriabarra
1775a21df39Svaleriabarra        return qref1d, qweight1d
1785a21df39Svaleriabarra
1795a21df39Svaleriabarra    # QR factorization
1805a21df39Svaleriabarra    @staticmethod
1815a21df39Svaleriabarra    def qr_factorization(ceed, mat, tau, m, n):
1825a21df39Svaleriabarra        """Return QR Factorization of a matrix.
1835a21df39Svaleriabarra
1845a21df39Svaleriabarra           Args:
1855a21df39Svaleriabarra             ceed: Ceed context currently in use
1865a21df39Svaleriabarra             *mat: Numpy array holding the row-major matrix to be factorized in place
1875a21df39Svaleriabarra             *tau: Numpy array to hold the vector of lengt m of scaling factors
1885a21df39Svaleriabarra             m: number of rows
1895a21df39Svaleriabarra             n: numbef of columns"""
1905a21df39Svaleriabarra
1915a21df39Svaleriabarra        # Setup arguments
1925a21df39Svaleriabarra        mat_pointer = ffi.new("CeedScalar *")
193*7a7b0fa3SJed Brown        mat_pointer = ffi.cast(
194*7a7b0fa3SJed Brown            "CeedScalar *",
195*7a7b0fa3SJed Brown            mat.__array_interface__['data'][0])
1965a21df39Svaleriabarra
1975a21df39Svaleriabarra        tau_pointer = ffi.new("CeedScalar *")
198*7a7b0fa3SJed Brown        tau_pointer = ffi.cast(
199*7a7b0fa3SJed Brown            "CeedScalar *",
200*7a7b0fa3SJed Brown            tau.__array_interface__['data'][0])
2015a21df39Svaleriabarra
2025a21df39Svaleriabarra        # libCEED call
203*7a7b0fa3SJed Brown        lib.CeedQRFactorization(
204*7a7b0fa3SJed Brown            ceed._pointer[0], mat_pointer, tau_pointer, m, n)
2055a21df39Svaleriabarra
2065a21df39Svaleriabarra        return mat, tau
2075a21df39Svaleriabarra
2085a21df39Svaleriabarra    # Symmetric Schur decomposition
2095a21df39Svaleriabarra    @staticmethod
2105a21df39Svaleriabarra    def symmetric_schur_decomposition(ceed, mat, n):
2115a21df39Svaleriabarra        """Return symmetric Schur decomposition of a symmetric matrix
2125a21df39Svaleriabarra             via symmetric QR factorization.
2135a21df39Svaleriabarra
2145a21df39Svaleriabarra           Args:
2155a21df39Svaleriabarra             ceed: Ceed context currently in use
2165a21df39Svaleriabarra             *mat: Numpy array holding the row-major matrix to be factorized in place
2175a21df39Svaleriabarra             n: number of rows/columns
2185a21df39Svaleriabarra
2195a21df39Svaleriabarra           Returns:
2205a21df39Svaleriabarra             lbda: Numpy array of length n holding eigenvalues"""
2215a21df39Svaleriabarra
2225a21df39Svaleriabarra        # Setup arguments
2235a21df39Svaleriabarra        mat_pointer = ffi.new("CeedScalar *")
224*7a7b0fa3SJed Brown        mat_pointer = ffi.cast(
225*7a7b0fa3SJed Brown            "CeedScalar *",
226*7a7b0fa3SJed Brown            mat.__array_interface__['data'][0])
2275a21df39Svaleriabarra
2285a21df39Svaleriabarra        lbda = np.empty(n, dtype="float64")
2295a21df39Svaleriabarra        l_pointer = ffi.new("CeedScalar *")
230*7a7b0fa3SJed Brown        l_pointer = ffi.cast(
231*7a7b0fa3SJed Brown            "CeedScalar *",
232*7a7b0fa3SJed Brown            lbda.__array_interface__['data'][0])
2335a21df39Svaleriabarra
2345a21df39Svaleriabarra        # libCEED call
235*7a7b0fa3SJed Brown        lib.CeedSymmetricSchurDecomposition(
236*7a7b0fa3SJed Brown            ceed._pointer[0], mat_pointer, l_pointer, n)
2375a21df39Svaleriabarra
2385a21df39Svaleriabarra        return lbda
2395a21df39Svaleriabarra
2405a21df39Svaleriabarra    # Simultaneous Diagonalization
2415a21df39Svaleriabarra    @staticmethod
2425a21df39Svaleriabarra    def simultaneous_diagonalization(ceed, matA, matB, n):
2435a21df39Svaleriabarra        """Return Simultaneous Diagonalization of two matrices.
2445a21df39Svaleriabarra
2455a21df39Svaleriabarra           Args:
2465a21df39Svaleriabarra             ceed: Ceed context currently in use
2475a21df39Svaleriabarra             *matA: Numpy array holding the row-major matrix to be factorized with
2485a21df39Svaleriabarra                      eigenvalues
2495a21df39Svaleriabarra             *matB: Numpy array holding the row-major matrix to be factorized to identity
2505a21df39Svaleriabarra             n: number of rows/columns
2515a21df39Svaleriabarra
2525a21df39Svaleriabarra           Returns:
2535a21df39Svaleriabarra             (x, lbda): Numpy array holding the row-major orthogonal matrix and
2545a21df39Svaleriabarra                          Numpy array holding the vector of length n of generalized
2555a21df39Svaleriabarra                          eigenvalues"""
2565a21df39Svaleriabarra
2575a21df39Svaleriabarra        # Setup arguments
2585a21df39Svaleriabarra        matA_pointer = ffi.new("CeedScalar *")
259*7a7b0fa3SJed Brown        matA_pointer = ffi.cast(
260*7a7b0fa3SJed Brown            "CeedScalar *",
261*7a7b0fa3SJed Brown            matA.__array_interface__['data'][0])
2625a21df39Svaleriabarra
2635a21df39Svaleriabarra        matB_pointer = ffi.new("CeedScalar *")
264*7a7b0fa3SJed Brown        matB_pointer = ffi.cast(
265*7a7b0fa3SJed Brown            "CeedScalar *",
266*7a7b0fa3SJed Brown            matB.__array_interface__['data'][0])
2675a21df39Svaleriabarra
2685a21df39Svaleriabarra        lbda = np.empty(n, dtype="float64")
2695a21df39Svaleriabarra        l_pointer = ffi.new("CeedScalar *")
270*7a7b0fa3SJed Brown        l_pointer = ffi.cast(
271*7a7b0fa3SJed Brown            "CeedScalar *",
272*7a7b0fa3SJed Brown            lbda.__array_interface__['data'][0])
2735a21df39Svaleriabarra
2745a21df39Svaleriabarra        x = np.empty(n * n, dtype="float64")
2755a21df39Svaleriabarra        x_pointer = ffi.new("CeedScalar *")
2765a21df39Svaleriabarra        x_pointer = ffi.cast("CeedScalar *", x.__array_interface__['data'][0])
2775a21df39Svaleriabarra
2785a21df39Svaleriabarra        # libCEED call
2795a21df39Svaleriabarra        lib.CeedSimultaneousDiagonalization(ceed._pointer[0], matA_pointer, matB_pointer,
2805a21df39Svaleriabarra                                            x_pointer, l_pointer, n)
2815a21df39Svaleriabarra
2825a21df39Svaleriabarra        return x, lbda
2835a21df39Svaleriabarra
2845a21df39Svaleriabarra    # Destructor
2855a21df39Svaleriabarra    def __del__(self):
2865a21df39Svaleriabarra        # libCEED call
2875a21df39Svaleriabarra        lib.CeedBasisDestroy(self._pointer)
2885a21df39Svaleriabarra
2895a21df39Svaleriabarra# ------------------------------------------------------------------------------
290*7a7b0fa3SJed Brown
291*7a7b0fa3SJed Brown
2925a21df39Svaleriabarraclass BasisTensorH1(Basis):
2935a21df39Svaleriabarra    """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
2945a21df39Svaleriabarra         H^1 discretizations."""
2955a21df39Svaleriabarra
2965a21df39Svaleriabarra    # Constructor
2975a21df39Svaleriabarra    def __init__(self, ceed, dim, ncomp, P1d, Q1d, interp1d, grad1d,
2985a21df39Svaleriabarra                 qref1d, qweight1d):
2995a21df39Svaleriabarra
3005a21df39Svaleriabarra        # Setup arguments
3015a21df39Svaleriabarra        self._pointer = ffi.new("CeedBasis *")
3025a21df39Svaleriabarra
3035a21df39Svaleriabarra        self._ceed = ceed
3045a21df39Svaleriabarra
3055a21df39Svaleriabarra        interp1d_pointer = ffi.new("CeedScalar *")
306*7a7b0fa3SJed Brown        interp1d_pointer = ffi.cast(
307*7a7b0fa3SJed Brown            "CeedScalar *",
308*7a7b0fa3SJed Brown            interp1d.__array_interface__['data'][0])
3095a21df39Svaleriabarra
3105a21df39Svaleriabarra        grad1d_pointer = ffi.new("CeedScalar *")
311*7a7b0fa3SJed Brown        grad1d_pointer = ffi.cast(
312*7a7b0fa3SJed Brown            "CeedScalar *",
313*7a7b0fa3SJed Brown            grad1d.__array_interface__['data'][0])
3145a21df39Svaleriabarra
3155a21df39Svaleriabarra        qref1d_pointer = ffi.new("CeedScalar *")
316*7a7b0fa3SJed Brown        qref1d_pointer = ffi.cast(
317*7a7b0fa3SJed Brown            "CeedScalar *",
318*7a7b0fa3SJed Brown            qref1d.__array_interface__['data'][0])
3195a21df39Svaleriabarra
3205a21df39Svaleriabarra        qweight1d_pointer = ffi.new("CeedScalar *")
321*7a7b0fa3SJed Brown        qweight1d_pointer = ffi.cast(
322*7a7b0fa3SJed Brown            "CeedScalar *",
323*7a7b0fa3SJed Brown            qweight1d.__array_interface__['data'][0])
3245a21df39Svaleriabarra
3255a21df39Svaleriabarra        # libCEED call
3265a21df39Svaleriabarra        lib.CeedBasisCreateTensorH1(self._ceed._pointer[0], dim, ncomp, P1d, Q1d,
3275a21df39Svaleriabarra                                    interp1d_pointer, grad1d_pointer, qref1d_pointer,
3285a21df39Svaleriabarra                                    qweight1d_pointer, self._pointer)
3295a21df39Svaleriabarra
3305a21df39Svaleriabarra# ------------------------------------------------------------------------------
331*7a7b0fa3SJed Brown
332*7a7b0fa3SJed Brown
3335a21df39Svaleriabarraclass BasisTensorH1Lagrange(Basis):
3345a21df39Svaleriabarra    """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange basis
3355a21df39Svaleriabarra         objects for H^1 discretizations."""
3365a21df39Svaleriabarra
3375a21df39Svaleriabarra    # Constructor
3385a21df39Svaleriabarra    def __init__(self, ceed, dim, ncomp, P, Q, qmode):
3395a21df39Svaleriabarra
3405a21df39Svaleriabarra        # Setup arguments
3415a21df39Svaleriabarra        self._pointer = ffi.new("CeedBasis *")
3425a21df39Svaleriabarra
3435a21df39Svaleriabarra        self._ceed = ceed
3445a21df39Svaleriabarra
3455a21df39Svaleriabarra        # libCEED call
3465a21df39Svaleriabarra        lib.CeedBasisCreateTensorH1Lagrange(self._ceed._pointer[0], dim, ncomp, P,
3475a21df39Svaleriabarra                                            Q, qmode, self._pointer)
3485a21df39Svaleriabarra
3495a21df39Svaleriabarra# ------------------------------------------------------------------------------
350*7a7b0fa3SJed Brown
351*7a7b0fa3SJed Brown
3525a21df39Svaleriabarraclass BasisH1(Basis):
3535a21df39Svaleriabarra    """Ceed H1 Basis: finite element non tensor-product basis for H^1 discretizations."""
3545a21df39Svaleriabarra
3555a21df39Svaleriabarra    # Constructor
356*7a7b0fa3SJed Brown    def __init__(self, ceed, topo, ncomp, nnodes,
357*7a7b0fa3SJed Brown                 nqpts, interp, grad, qref, qweight):
3585a21df39Svaleriabarra
3595a21df39Svaleriabarra        # Setup arguments
3605a21df39Svaleriabarra        self._pointer = ffi.new("CeedBasis *")
3615a21df39Svaleriabarra
3625a21df39Svaleriabarra        self._ceed = ceed
3635a21df39Svaleriabarra
3645a21df39Svaleriabarra        interp_pointer = ffi.new("CeedScalar *")
365*7a7b0fa3SJed Brown        interp_pointer = ffi.cast(
366*7a7b0fa3SJed Brown            "CeedScalar *",
367*7a7b0fa3SJed Brown            interp.__array_interface__['data'][0])
3685a21df39Svaleriabarra
3695a21df39Svaleriabarra        grad_pointer = ffi.new("CeedScalar *")
370*7a7b0fa3SJed Brown        grad_pointer = ffi.cast(
371*7a7b0fa3SJed Brown            "CeedScalar *",
372*7a7b0fa3SJed Brown            grad.__array_interface__['data'][0])
3735a21df39Svaleriabarra
3745a21df39Svaleriabarra        qref_pointer = ffi.new("CeedScalar *")
375*7a7b0fa3SJed Brown        qref_pointer = ffi.cast(
376*7a7b0fa3SJed Brown            "CeedScalar *",
377*7a7b0fa3SJed Brown            qref.__array_interface__['data'][0])
3785a21df39Svaleriabarra
3795a21df39Svaleriabarra        qweight_pointer = ffi.new("CeedScalar *")
380*7a7b0fa3SJed Brown        qweight_pointer = ffi.cast(
381*7a7b0fa3SJed Brown            "CeedScalar *",
382*7a7b0fa3SJed Brown            qweight.__array_interface__['data'][0])
3835a21df39Svaleriabarra
3845a21df39Svaleriabarra        # libCEED call
3855a21df39Svaleriabarra        lib.CeedBasisCreateH1(self._ceed._pointer[0], topo, ncomp, nnodes, nqpts,
3865a21df39Svaleriabarra                              interp_pointer, grad_pointer, qref_pointer,
3875a21df39Svaleriabarra                              qweight_pointer, self._pointer)
3885a21df39Svaleriabarra
3895a21df39Svaleriabarra# ------------------------------------------------------------------------------
390*7a7b0fa3SJed Brown
391*7a7b0fa3SJed Brown
3925a21df39Svaleriabarraclass TransposeBasis():
3935a21df39Svaleriabarra    """Transpose Ceed Basis: transpose of finite element tensor-product basis objects."""
3945a21df39Svaleriabarra
3955a21df39Svaleriabarra    # Attributes
3965a21df39Svaleriabarra    _basis = None
3975a21df39Svaleriabarra
3985a21df39Svaleriabarra    # Constructor
3995a21df39Svaleriabarra    def __init__(self, basis):
4005a21df39Svaleriabarra
4015a21df39Svaleriabarra        # Reference basis
4025a21df39Svaleriabarra        self._basis = basis
4035a21df39Svaleriabarra
4045a21df39Svaleriabarra    # Representation
4055a21df39Svaleriabarra    def __repr__(self):
4065a21df39Svaleriabarra        return "<Transpose CeedBasis instance at " + hex(id(self)) + ">"
4075a21df39Svaleriabarra
4085a21df39Svaleriabarra    # Apply Transpose Basis
4095a21df39Svaleriabarra    def apply(self, nelem, emode, u, v):
4105a21df39Svaleriabarra        """Apply basis evaluation from quadrature points to nodes.
4115a21df39Svaleriabarra
4125a21df39Svaleriabarra           Args:
4135a21df39Svaleriabarra             nelem: the number of elements to apply the basis evaluation to;
4145a21df39Svaleriabarra                      the backend will specify the ordering in a
4155a21df39Svaleriabarra                      Blocked ElemRestriction
4165a21df39Svaleriabarra             **emode: basis evaluation mode
4175a21df39Svaleriabarra             u: input vector
4185a21df39Svaleriabarra             v: output vector"""
4195a21df39Svaleriabarra
4205a21df39Svaleriabarra        # libCEED call
4215a21df39Svaleriabarra        self._basis.apply(nelem, emode, u, v, tmode=TRANSPOSE)
4225a21df39Svaleriabarra
4235a21df39Svaleriabarra# ------------------------------------------------------------------------------
424