xref: /libCEED/python/ceed_basis.py (revision 9ba83ac0e4b1fca39d6fa6737a318a9f0cbc172d)
1*9ba83ac0SJeremy L Thompson# Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors
23d8e8822SJeremy L Thompson# All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
35a21df39Svaleriabarra#
43d8e8822SJeremy L Thompson# SPDX-License-Identifier: BSD-2-Clause
55a21df39Svaleriabarra#
63d8e8822SJeremy L Thompson# This file is part of CEED:  http://github.com/ceed
75a21df39Svaleriabarra
85a21df39Svaleriabarrafrom _ceed_cffi import ffi, lib
95a21df39Svaleriabarraimport tempfile
105a21df39Svaleriabarraimport numpy as np
115a21df39Svaleriabarrafrom abc import ABC
1280a9ef05SNatalie Beamsfrom .ceed_constants import TRANSPOSE, NOTRANSPOSE, scalar_types
135a21df39Svaleriabarra
145a21df39Svaleriabarra# ------------------------------------------------------------------------------
157a7b0fa3SJed Brown
167a7b0fa3SJed Brown
175a21df39Svaleriabarraclass Basis(ABC):
185a21df39Svaleriabarra    """Ceed Basis: finite element basis objects."""
195a21df39Svaleriabarra
205a21df39Svaleriabarra    # Representation
215a21df39Svaleriabarra    def __repr__(self):
225a21df39Svaleriabarra        return "<CeedBasis instance at " + hex(id(self)) + ">"
235a21df39Svaleriabarra
245a21df39Svaleriabarra    # String conversion for print() to stdout
255a21df39Svaleriabarra    def __str__(self):
265a21df39Svaleriabarra        """View a Basis via print()."""
275a21df39Svaleriabarra
285a21df39Svaleriabarra        # libCEED call
295a21df39Svaleriabarra        with tempfile.NamedTemporaryFile() as key_file:
305a21df39Svaleriabarra            with open(key_file.name, 'r+') as stream_file:
315a21df39Svaleriabarra                stream = ffi.cast("FILE *", stream_file)
325a21df39Svaleriabarra
33477729cfSJeremy L Thompson                err_code = lib.CeedBasisView(self._pointer[0], stream)
34477729cfSJeremy L Thompson                self._ceed._check_error(err_code)
355a21df39Svaleriabarra
365a21df39Svaleriabarra                stream_file.seek(0)
375a21df39Svaleriabarra                out_string = stream_file.read()
385a21df39Svaleriabarra
395a21df39Svaleriabarra        return out_string
405a21df39Svaleriabarra
415a21df39Svaleriabarra    # Apply Basis
425a21df39Svaleriabarra    def apply(self, nelem, emode, u, v, tmode=NOTRANSPOSE):
435a21df39Svaleriabarra        """Apply basis evaluation from nodes to quadrature points or vice versa.
445a21df39Svaleriabarra
455a21df39Svaleriabarra           Args:
465a21df39Svaleriabarra             nelem: the number of elements to apply the basis evaluation to;
475a21df39Svaleriabarra                      the backend will specify the ordering in a
485a21df39Svaleriabarra                      Blocked ElemRestriction
495a21df39Svaleriabarra             emode: basis evaluation mode
505a21df39Svaleriabarra             u: input vector
515a21df39Svaleriabarra             v: output vector
525a21df39Svaleriabarra             **tmode: CEED_NOTRANSPOSE to evaluate from nodes to quadrature
535a21df39Svaleriabarra                        points, CEED_TRANSPOSE to apply the transpose, mapping
545a21df39Svaleriabarra                        from quadrature points to nodes; default CEED_NOTRANSPOSE"""
555a21df39Svaleriabarra
565a21df39Svaleriabarra        # libCEED call
57477729cfSJeremy L Thompson        err_code = lib.CeedBasisApply(self._pointer[0], nelem, tmode, emode,
585a21df39Svaleriabarra                                      u._pointer[0], v._pointer[0])
59477729cfSJeremy L Thompson        self._ceed._check_error(err_code)
605a21df39Svaleriabarra
615a21df39Svaleriabarra    # Transpose a Basis
625a21df39Svaleriabarra    @property
635a21df39Svaleriabarra    def T(self):
645a21df39Svaleriabarra        """Transpose a Basis."""
655a21df39Svaleriabarra
665a21df39Svaleriabarra        return TransposeBasis(self)
675a21df39Svaleriabarra
685a21df39Svaleriabarra    # Transpose a Basis
695a21df39Svaleriabarra    @property
705a21df39Svaleriabarra    def transpose(self):
715a21df39Svaleriabarra        """Transpose a Basis."""
725a21df39Svaleriabarra
735a21df39Svaleriabarra        return TransposeBasis(self)
745a21df39Svaleriabarra
755a21df39Svaleriabarra    # Get number of nodes
765a21df39Svaleriabarra    def get_num_nodes(self):
775a21df39Svaleriabarra        """Get total number of nodes (in dim dimensions) of a Basis.
785a21df39Svaleriabarra
795a21df39Svaleriabarra           Returns:
805a21df39Svaleriabarra             num_nodes: total number of nodes"""
815a21df39Svaleriabarra
825a21df39Svaleriabarra        # Setup argument
835a21df39Svaleriabarra        p_pointer = ffi.new("CeedInt *")
845a21df39Svaleriabarra
855a21df39Svaleriabarra        # libCEED call
86477729cfSJeremy L Thompson        err_code = lib.CeedBasisGetNumNodes(self._pointer[0], p_pointer)
87477729cfSJeremy L Thompson        self._ceed._check_error(err_code)
885a21df39Svaleriabarra
895a21df39Svaleriabarra        return p_pointer[0]
905a21df39Svaleriabarra
915a21df39Svaleriabarra    # Get number of quadrature points
925a21df39Svaleriabarra    def get_num_quadrature_points(self):
935a21df39Svaleriabarra        """Get total number of quadrature points (in dim dimensions) of a Basis.
945a21df39Svaleriabarra
955a21df39Svaleriabarra           Returns:
965a21df39Svaleriabarra             num_qpts: total number of quadrature points"""
975a21df39Svaleriabarra
985a21df39Svaleriabarra        # Setup argument
995a21df39Svaleriabarra        q_pointer = ffi.new("CeedInt *")
1005a21df39Svaleriabarra
1015a21df39Svaleriabarra        # libCEED call
102477729cfSJeremy L Thompson        err_code = lib.CeedBasisGetNumQuadraturePoints(
103477729cfSJeremy L Thompson            self._pointer[0], q_pointer)
104477729cfSJeremy L Thompson        self._ceed._check_error(err_code)
1055a21df39Svaleriabarra
1065a21df39Svaleriabarra        return q_pointer[0]
1075a21df39Svaleriabarra
1085a21df39Svaleriabarra    # Destructor
1095a21df39Svaleriabarra    def __del__(self):
1105a21df39Svaleriabarra        # libCEED call
111477729cfSJeremy L Thompson        err_code = lib.CeedBasisDestroy(self._pointer)
112477729cfSJeremy L Thompson        self._ceed._check_error(err_code)
1135a21df39Svaleriabarra
1145a21df39Svaleriabarra# ------------------------------------------------------------------------------
1157a7b0fa3SJed Brown
1167a7b0fa3SJed Brown
1175a21df39Svaleriabarraclass BasisTensorH1(Basis):
1185a21df39Svaleriabarra    """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
1195a21df39Svaleriabarra         H^1 discretizations."""
1205a21df39Svaleriabarra
1215a21df39Svaleriabarra    # Constructor
1225a21df39Svaleriabarra    def __init__(self, ceed, dim, ncomp, P1d, Q1d, interp1d, grad1d,
1235a21df39Svaleriabarra                 qref1d, qweight1d):
1245a21df39Svaleriabarra
1255a21df39Svaleriabarra        # Setup arguments
1265a21df39Svaleriabarra        self._pointer = ffi.new("CeedBasis *")
1275a21df39Svaleriabarra
1285a21df39Svaleriabarra        self._ceed = ceed
1295a21df39Svaleriabarra
1305a21df39Svaleriabarra        interp1d_pointer = ffi.new("CeedScalar *")
1317a7b0fa3SJed Brown        interp1d_pointer = ffi.cast(
1327a7b0fa3SJed Brown            "CeedScalar *",
1337a7b0fa3SJed Brown            interp1d.__array_interface__['data'][0])
1345a21df39Svaleriabarra
1355a21df39Svaleriabarra        grad1d_pointer = ffi.new("CeedScalar *")
1367a7b0fa3SJed Brown        grad1d_pointer = ffi.cast(
1377a7b0fa3SJed Brown            "CeedScalar *",
1387a7b0fa3SJed Brown            grad1d.__array_interface__['data'][0])
1395a21df39Svaleriabarra
1405a21df39Svaleriabarra        qref1d_pointer = ffi.new("CeedScalar *")
1417a7b0fa3SJed Brown        qref1d_pointer = ffi.cast(
1427a7b0fa3SJed Brown            "CeedScalar *",
1437a7b0fa3SJed Brown            qref1d.__array_interface__['data'][0])
1445a21df39Svaleriabarra
1455a21df39Svaleriabarra        qweight1d_pointer = ffi.new("CeedScalar *")
1467a7b0fa3SJed Brown        qweight1d_pointer = ffi.cast(
1477a7b0fa3SJed Brown            "CeedScalar *",
1487a7b0fa3SJed Brown            qweight1d.__array_interface__['data'][0])
1495a21df39Svaleriabarra
1505a21df39Svaleriabarra        # libCEED call
151477729cfSJeremy L Thompson        err_code = lib.CeedBasisCreateTensorH1(self._ceed._pointer[0], dim, ncomp,
152477729cfSJeremy L Thompson                                               P1d, Q1d, interp1d_pointer,
153477729cfSJeremy L Thompson                                               grad1d_pointer, qref1d_pointer,
1545a21df39Svaleriabarra                                               qweight1d_pointer, self._pointer)
155477729cfSJeremy L Thompson        self._ceed._check_error(err_code)
1565a21df39Svaleriabarra
15766ce82e0Sjeremylt    # Get 1D interpolation matrix
15866ce82e0Sjeremylt    def get_interp_1d(self):
159d99fa3c5SJeremy L Thompson        """Return 1D interpolation matrix of a tensor product Basis.
160d99fa3c5SJeremy L Thompson
161d99fa3c5SJeremy L Thompson           Returns:
162d99fa3c5SJeremy L Thompson             *array: Numpy array"""
163d99fa3c5SJeremy L Thompson
164d99fa3c5SJeremy L Thompson        # Compute the length of the array
165d99fa3c5SJeremy L Thompson        nnodes_pointer = ffi.new("CeedInt *")
166d99fa3c5SJeremy L Thompson        lib.CeedBasisGetNumNodes1D(self._pointer[0], nnodes_pointer)
167d99fa3c5SJeremy L Thompson        nqpts_pointer = ffi.new("CeedInt *")
168d99fa3c5SJeremy L Thompson        lib.CeedBasisGetNumQuadraturePoints1D(self._pointer[0], nqpts_pointer)
169d99fa3c5SJeremy L Thompson        length = nnodes_pointer[0] * nqpts_pointer[0]
170d99fa3c5SJeremy L Thompson
171d99fa3c5SJeremy L Thompson        # Setup the pointer's pointer
172d99fa3c5SJeremy L Thompson        array_pointer = ffi.new("CeedScalar **")
173d99fa3c5SJeremy L Thompson
174d99fa3c5SJeremy L Thompson        # libCEED call
175d99fa3c5SJeremy L Thompson        lib.CeedBasisGetInterp1D(self._pointer[0], array_pointer)
176d99fa3c5SJeremy L Thompson
177d99fa3c5SJeremy L Thompson        # Return array created from buffer
178d99fa3c5SJeremy L Thompson        # Create buffer object from returned pointer
179d99fa3c5SJeremy L Thompson        buff = ffi.buffer(
180d99fa3c5SJeremy L Thompson            array_pointer[0],
181d99fa3c5SJeremy L Thompson            ffi.sizeof("CeedScalar") *
182d99fa3c5SJeremy L Thompson            length)
183d99fa3c5SJeremy L Thompson        # return read only Numpy array
18480a9ef05SNatalie Beams        ret = np.frombuffer(buff, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
185d99fa3c5SJeremy L Thompson        ret.flags['WRITEABLE'] = False
186d99fa3c5SJeremy L Thompson        return ret
187d99fa3c5SJeremy L Thompson
18866ce82e0Sjeremylt    # Get 1D gradient matrix
18966ce82e0Sjeremylt    def get_grad_1d(self):
19066ce82e0Sjeremylt        """Return 1D gradient matrix of a tensor product Basis.
19166ce82e0Sjeremylt
19266ce82e0Sjeremylt           Returns:
19366ce82e0Sjeremylt             *array: Numpy array"""
19466ce82e0Sjeremylt
19566ce82e0Sjeremylt        # Compute the length of the array
19666ce82e0Sjeremylt        nnodes_pointer = ffi.new("CeedInt *")
19766ce82e0Sjeremylt        lib.CeedBasisGetNumNodes1D(self._pointer[0], nnodes_pointer)
19866ce82e0Sjeremylt        nqpts_pointer = ffi.new("CeedInt *")
19966ce82e0Sjeremylt        lib.CeedBasisGetNumQuadraturePoints1D(self._pointer[0], nqpts_pointer)
20066ce82e0Sjeremylt        length = nnodes_pointer[0] * nqpts_pointer[0]
20166ce82e0Sjeremylt
20266ce82e0Sjeremylt        # Setup the pointer's pointer
20366ce82e0Sjeremylt        array_pointer = ffi.new("CeedScalar **")
20466ce82e0Sjeremylt
20566ce82e0Sjeremylt        # libCEED call
20666ce82e0Sjeremylt        lib.CeedBasisGetGrad1D(self._pointer[0], array_pointer)
20766ce82e0Sjeremylt
20866ce82e0Sjeremylt        # Return array created from buffer
20966ce82e0Sjeremylt        # Create buffer object from returned pointer
21066ce82e0Sjeremylt        buff = ffi.buffer(
21166ce82e0Sjeremylt            array_pointer[0],
21266ce82e0Sjeremylt            ffi.sizeof("CeedScalar") *
21366ce82e0Sjeremylt            length)
21466ce82e0Sjeremylt        # return read only Numpy array
21566ce82e0Sjeremylt        ret = np.frombuffer(buff, dtype="float64")
21666ce82e0Sjeremylt        ret.flags['WRITEABLE'] = False
21766ce82e0Sjeremylt        return ret
21866ce82e0Sjeremylt
21966ce82e0Sjeremylt    # Get 1D quadrature weights matrix
22066ce82e0Sjeremylt    def get_q_weight_1d(self):
22166ce82e0Sjeremylt        """Return 1D quadrature weights matrix of a tensor product Basis.
22266ce82e0Sjeremylt
22366ce82e0Sjeremylt           Returns:
22466ce82e0Sjeremylt             *array: Numpy array"""
22566ce82e0Sjeremylt
22666ce82e0Sjeremylt        # Compute the length of the array
22766ce82e0Sjeremylt        nqpts_pointer = ffi.new("CeedInt *")
22866ce82e0Sjeremylt        lib.CeedBasisGetNumQuadraturePoints1D(self._pointer[0], nqpts_pointer)
22966ce82e0Sjeremylt        length = nqpts_pointer[0]
23066ce82e0Sjeremylt
23166ce82e0Sjeremylt        # Setup the pointer's pointer
23266ce82e0Sjeremylt        array_pointer = ffi.new("CeedScalar **")
23366ce82e0Sjeremylt
23466ce82e0Sjeremylt        # libCEED call
23566ce82e0Sjeremylt        lib.CeedBasisGetQWeights(self._pointer[0], array_pointer)
23666ce82e0Sjeremylt
23766ce82e0Sjeremylt        # Return array created from buffer
23866ce82e0Sjeremylt        # Create buffer object from returned pointer
23966ce82e0Sjeremylt        buff = ffi.buffer(
24066ce82e0Sjeremylt            array_pointer[0],
24166ce82e0Sjeremylt            ffi.sizeof("CeedScalar") *
24266ce82e0Sjeremylt            length)
24366ce82e0Sjeremylt        # return read only Numpy array
24466ce82e0Sjeremylt        ret = np.frombuffer(buff, dtype="float64")
24566ce82e0Sjeremylt        ret.flags['WRITEABLE'] = False
24666ce82e0Sjeremylt        return ret
24766ce82e0Sjeremylt
24866ce82e0Sjeremylt    # Get 1D quadrature points matrix
24966ce82e0Sjeremylt    def get_q_ref_1d(self):
25066ce82e0Sjeremylt        """Return 1D quadrature points matrix of a tensor product Basis.
25166ce82e0Sjeremylt
25266ce82e0Sjeremylt           Returns:
25366ce82e0Sjeremylt             *array: Numpy array"""
25466ce82e0Sjeremylt
25566ce82e0Sjeremylt        # Compute the length of the array
25666ce82e0Sjeremylt        nqpts_pointer = ffi.new("CeedInt *")
25766ce82e0Sjeremylt        lib.CeedBasisGetNumQuadraturePoints1D(self._pointer[0], nqpts_pointer)
25866ce82e0Sjeremylt        length = nqpts_pointer[0]
25966ce82e0Sjeremylt
26066ce82e0Sjeremylt        # Setup the pointer's pointer
26166ce82e0Sjeremylt        array_pointer = ffi.new("CeedScalar **")
26266ce82e0Sjeremylt
26366ce82e0Sjeremylt        # libCEED call
26466ce82e0Sjeremylt        lib.CeedBasisGetQRef(self._pointer[0], array_pointer)
26566ce82e0Sjeremylt
26666ce82e0Sjeremylt        # Return array created from buffer
26766ce82e0Sjeremylt        # Create buffer object from returned pointer
26866ce82e0Sjeremylt        buff = ffi.buffer(
26966ce82e0Sjeremylt            array_pointer[0],
27066ce82e0Sjeremylt            ffi.sizeof("CeedScalar") *
27166ce82e0Sjeremylt            length)
27266ce82e0Sjeremylt        # return read only Numpy array
27366ce82e0Sjeremylt        ret = np.frombuffer(buff, dtype="float64")
27466ce82e0Sjeremylt        ret.flags['WRITEABLE'] = False
27566ce82e0Sjeremylt        return ret
27666ce82e0Sjeremylt
277d99fa3c5SJeremy L Thompson
2785a21df39Svaleriabarra# ------------------------------------------------------------------------------
2797a7b0fa3SJed Brown
2807a7b0fa3SJed Brown
281d99fa3c5SJeremy L Thompsonclass BasisTensorH1Lagrange(BasisTensorH1):
2825a21df39Svaleriabarra    """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange basis
2835a21df39Svaleriabarra         objects for H^1 discretizations."""
2845a21df39Svaleriabarra
2855a21df39Svaleriabarra    # Constructor
2865a21df39Svaleriabarra    def __init__(self, ceed, dim, ncomp, P, Q, qmode):
2875a21df39Svaleriabarra
2885a21df39Svaleriabarra        # Setup arguments
2895a21df39Svaleriabarra        self._pointer = ffi.new("CeedBasis *")
2905a21df39Svaleriabarra
2915a21df39Svaleriabarra        self._ceed = ceed
2925a21df39Svaleriabarra
2935a21df39Svaleriabarra        # libCEED call
294477729cfSJeremy L Thompson        err_code = lib.CeedBasisCreateTensorH1Lagrange(self._ceed._pointer[0], dim,
295477729cfSJeremy L Thompson                                                       ncomp, P, Q, qmode, self._pointer)
296477729cfSJeremy L Thompson        self._ceed._check_error(err_code)
2975a21df39Svaleriabarra
2985a21df39Svaleriabarra# ------------------------------------------------------------------------------
2997a7b0fa3SJed Brown
3007a7b0fa3SJed Brown
3015a21df39Svaleriabarraclass BasisH1(Basis):
3025a21df39Svaleriabarra    """Ceed H1 Basis: finite element non tensor-product basis for H^1 discretizations."""
3035a21df39Svaleriabarra
3045a21df39Svaleriabarra    # Constructor
3057a7b0fa3SJed Brown    def __init__(self, ceed, topo, ncomp, nnodes,
3067a7b0fa3SJed Brown                 nqpts, interp, grad, qref, qweight):
3075a21df39Svaleriabarra
3085a21df39Svaleriabarra        # Setup arguments
3095a21df39Svaleriabarra        self._pointer = ffi.new("CeedBasis *")
3105a21df39Svaleriabarra
3115a21df39Svaleriabarra        self._ceed = ceed
3125a21df39Svaleriabarra
3135a21df39Svaleriabarra        interp_pointer = ffi.new("CeedScalar *")
3147a7b0fa3SJed Brown        interp_pointer = ffi.cast(
3157a7b0fa3SJed Brown            "CeedScalar *",
3167a7b0fa3SJed Brown            interp.__array_interface__['data'][0])
3175a21df39Svaleriabarra
3185a21df39Svaleriabarra        grad_pointer = ffi.new("CeedScalar *")
3197a7b0fa3SJed Brown        grad_pointer = ffi.cast(
3207a7b0fa3SJed Brown            "CeedScalar *",
3217a7b0fa3SJed Brown            grad.__array_interface__['data'][0])
3225a21df39Svaleriabarra
3235a21df39Svaleriabarra        qref_pointer = ffi.new("CeedScalar *")
3247a7b0fa3SJed Brown        qref_pointer = ffi.cast(
3257a7b0fa3SJed Brown            "CeedScalar *",
3267a7b0fa3SJed Brown            qref.__array_interface__['data'][0])
3275a21df39Svaleriabarra
3285a21df39Svaleriabarra        qweight_pointer = ffi.new("CeedScalar *")
3297a7b0fa3SJed Brown        qweight_pointer = ffi.cast(
3307a7b0fa3SJed Brown            "CeedScalar *",
3317a7b0fa3SJed Brown            qweight.__array_interface__['data'][0])
3325a21df39Svaleriabarra
3335a21df39Svaleriabarra        # libCEED call
334477729cfSJeremy L Thompson        err_code = lib.CeedBasisCreateH1(self._ceed._pointer[0], topo, ncomp,
335477729cfSJeremy L Thompson                                         nnodes, nqpts, interp_pointer,
336477729cfSJeremy L Thompson                                         grad_pointer, qref_pointer,
3375a21df39Svaleriabarra                                         qweight_pointer, self._pointer)
3385a21df39Svaleriabarra
3395a21df39Svaleriabarra# ------------------------------------------------------------------------------
3407a7b0fa3SJed Brown
3417a7b0fa3SJed Brown
34211b88ddaSSebastian Grimbergclass BasisHdiv(Basis):
34311b88ddaSSebastian Grimberg    """Ceed Hdiv Basis: finite element non tensor-product basis for H(div) discretizations."""
34411b88ddaSSebastian Grimberg
34511b88ddaSSebastian Grimberg    # Constructor
34611b88ddaSSebastian Grimberg    def __init__(self, ceed, topo, ncomp, nnodes,
34711b88ddaSSebastian Grimberg                 nqpts, interp, grad, qref, qweight):
34811b88ddaSSebastian Grimberg
34911b88ddaSSebastian Grimberg        # Setup arguments
35011b88ddaSSebastian Grimberg        self._pointer = ffi.new("CeedBasis *")
35111b88ddaSSebastian Grimberg
35211b88ddaSSebastian Grimberg        self._ceed = ceed
35311b88ddaSSebastian Grimberg
35411b88ddaSSebastian Grimberg        interp_pointer = ffi.new("CeedScalar *")
35511b88ddaSSebastian Grimberg        interp_pointer = ffi.cast(
35611b88ddaSSebastian Grimberg            "CeedScalar *",
35711b88ddaSSebastian Grimberg            interp.__array_interface__['data'][0])
35811b88ddaSSebastian Grimberg
35911b88ddaSSebastian Grimberg        div_pointer = ffi.new("CeedScalar *")
36011b88ddaSSebastian Grimberg        div_pointer = ffi.cast(
36111b88ddaSSebastian Grimberg            "CeedScalar *",
36211b88ddaSSebastian Grimberg            grad.__array_interface__['data'][0])
36311b88ddaSSebastian Grimberg
36411b88ddaSSebastian Grimberg        qref_pointer = ffi.new("CeedScalar *")
36511b88ddaSSebastian Grimberg        qref_pointer = ffi.cast(
36611b88ddaSSebastian Grimberg            "CeedScalar *",
36711b88ddaSSebastian Grimberg            qref.__array_interface__['data'][0])
36811b88ddaSSebastian Grimberg
36911b88ddaSSebastian Grimberg        qweight_pointer = ffi.new("CeedScalar *")
37011b88ddaSSebastian Grimberg        qweight_pointer = ffi.cast(
37111b88ddaSSebastian Grimberg            "CeedScalar *",
37211b88ddaSSebastian Grimberg            qweight.__array_interface__['data'][0])
37311b88ddaSSebastian Grimberg
37411b88ddaSSebastian Grimberg        # libCEED call
37511b88ddaSSebastian Grimberg        err_code = lib.CeedBasisCreateHdiv(self._ceed._pointer[0], topo, ncomp,
37611b88ddaSSebastian Grimberg                                           nnodes, nqpts, interp_pointer,
37711b88ddaSSebastian Grimberg                                           div_pointer, qref_pointer,
37811b88ddaSSebastian Grimberg                                           qweight_pointer, self._pointer)
37911b88ddaSSebastian Grimberg
38011b88ddaSSebastian Grimberg# ------------------------------------------------------------------------------
38111b88ddaSSebastian Grimberg
38211b88ddaSSebastian Grimberg
38311b88ddaSSebastian Grimbergclass BasisHcurl(Basis):
38411b88ddaSSebastian Grimberg    """Ceed Hcurl Basis: finite element non tensor-product basis for H(curl) discretizations."""
38511b88ddaSSebastian Grimberg
38611b88ddaSSebastian Grimberg    # Constructor
38711b88ddaSSebastian Grimberg    def __init__(self, ceed, topo, ncomp, nnodes,
38811b88ddaSSebastian Grimberg                 nqpts, interp, grad, qref, qweight):
38911b88ddaSSebastian Grimberg
39011b88ddaSSebastian Grimberg        # Setup arguments
39111b88ddaSSebastian Grimberg        self._pointer = ffi.new("CeedBasis *")
39211b88ddaSSebastian Grimberg
39311b88ddaSSebastian Grimberg        self._ceed = ceed
39411b88ddaSSebastian Grimberg
39511b88ddaSSebastian Grimberg        interp_pointer = ffi.new("CeedScalar *")
39611b88ddaSSebastian Grimberg        interp_pointer = ffi.cast(
39711b88ddaSSebastian Grimberg            "CeedScalar *",
39811b88ddaSSebastian Grimberg            interp.__array_interface__['data'][0])
39911b88ddaSSebastian Grimberg
40011b88ddaSSebastian Grimberg        curl_pointer = ffi.new("CeedScalar *")
40111b88ddaSSebastian Grimberg        curl_pointer = ffi.cast(
40211b88ddaSSebastian Grimberg            "CeedScalar *",
40311b88ddaSSebastian Grimberg            grad.__array_interface__['data'][0])
40411b88ddaSSebastian Grimberg
40511b88ddaSSebastian Grimberg        qref_pointer = ffi.new("CeedScalar *")
40611b88ddaSSebastian Grimberg        qref_pointer = ffi.cast(
40711b88ddaSSebastian Grimberg            "CeedScalar *",
40811b88ddaSSebastian Grimberg            qref.__array_interface__['data'][0])
40911b88ddaSSebastian Grimberg
41011b88ddaSSebastian Grimberg        qweight_pointer = ffi.new("CeedScalar *")
41111b88ddaSSebastian Grimberg        qweight_pointer = ffi.cast(
41211b88ddaSSebastian Grimberg            "CeedScalar *",
41311b88ddaSSebastian Grimberg            qweight.__array_interface__['data'][0])
41411b88ddaSSebastian Grimberg
41511b88ddaSSebastian Grimberg        # libCEED call
41611b88ddaSSebastian Grimberg        err_code = lib.CeedBasisCreateHcurl(self._ceed._pointer[0], topo, ncomp,
41711b88ddaSSebastian Grimberg                                            nnodes, nqpts, interp_pointer,
41811b88ddaSSebastian Grimberg                                            curl_pointer, qref_pointer,
41911b88ddaSSebastian Grimberg                                            qweight_pointer, self._pointer)
42011b88ddaSSebastian Grimberg
42111b88ddaSSebastian Grimberg# ------------------------------------------------------------------------------
42211b88ddaSSebastian Grimberg
42311b88ddaSSebastian Grimberg
4245a21df39Svaleriabarraclass TransposeBasis():
4255a21df39Svaleriabarra    """Transpose Ceed Basis: transpose of finite element tensor-product basis objects."""
4265a21df39Svaleriabarra
4275a21df39Svaleriabarra    # Attributes
4285a21df39Svaleriabarra    _basis = None
4295a21df39Svaleriabarra
4305a21df39Svaleriabarra    # Constructor
4315a21df39Svaleriabarra    def __init__(self, basis):
4325a21df39Svaleriabarra
4335a21df39Svaleriabarra        # Reference basis
4345a21df39Svaleriabarra        self._basis = basis
4355a21df39Svaleriabarra
4365a21df39Svaleriabarra    # Representation
4375a21df39Svaleriabarra    def __repr__(self):
4385a21df39Svaleriabarra        return "<Transpose CeedBasis instance at " + hex(id(self)) + ">"
4395a21df39Svaleriabarra
4405a21df39Svaleriabarra    # Apply Transpose Basis
4415a21df39Svaleriabarra    def apply(self, nelem, emode, u, v):
4425a21df39Svaleriabarra        """Apply basis evaluation from quadrature points to nodes.
4435a21df39Svaleriabarra
4445a21df39Svaleriabarra           Args:
4455a21df39Svaleriabarra             nelem: the number of elements to apply the basis evaluation to;
4465a21df39Svaleriabarra                      the backend will specify the ordering in a
4475a21df39Svaleriabarra                      Blocked ElemRestriction
4485a21df39Svaleriabarra             **emode: basis evaluation mode
4495a21df39Svaleriabarra             u: input vector
4505a21df39Svaleriabarra             v: output vector"""
4515a21df39Svaleriabarra
4525a21df39Svaleriabarra        # libCEED call
4535a21df39Svaleriabarra        self._basis.apply(nelem, emode, u, v, tmode=TRANSPOSE)
4545a21df39Svaleriabarra
4555a21df39Svaleriabarra# ------------------------------------------------------------------------------
456