xref: /libCEED/python/ceed_basis.py (revision 5a21df3975c43d9e44497a4a1e70e8ac7dca4bbf)
1*5a21df39Svaleriabarra# Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2*5a21df39Svaleriabarra# the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3*5a21df39Svaleriabarra# reserved. See files LICENSE and NOTICE for details.
4*5a21df39Svaleriabarra#
5*5a21df39Svaleriabarra# This file is part of CEED, a collection of benchmarks, miniapps, software
6*5a21df39Svaleriabarra# libraries and APIs for efficient high-order finite element and spectral
7*5a21df39Svaleriabarra# element discretizations for exascale applications. For more information and
8*5a21df39Svaleriabarra# source code availability see http://github.com/ceed.
9*5a21df39Svaleriabarra#
10*5a21df39Svaleriabarra# The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11*5a21df39Svaleriabarra# a collaborative effort of two U.S. Department of Energy organizations (Office
12*5a21df39Svaleriabarra# of Science and the National Nuclear Security Administration) responsible for
13*5a21df39Svaleriabarra# the planning and preparation of a capable exascale ecosystem, including
14*5a21df39Svaleriabarra# software, applications, hardware, advanced system engineering and early
15*5a21df39Svaleriabarra# testbed platforms, in support of the nation's exascale computing imperative.
16*5a21df39Svaleriabarra
17*5a21df39Svaleriabarrafrom _ceed_cffi import ffi, lib
18*5a21df39Svaleriabarraimport tempfile
19*5a21df39Svaleriabarraimport numpy as np
20*5a21df39Svaleriabarrafrom abc import ABC
21*5a21df39Svaleriabarrafrom .ceed_constants import TRANSPOSE, NOTRANSPOSE
22*5a21df39Svaleriabarra
23*5a21df39Svaleriabarra# ------------------------------------------------------------------------------
24*5a21df39Svaleriabarraclass Basis(ABC):
25*5a21df39Svaleriabarra  """Ceed Basis: finite element basis objects."""
26*5a21df39Svaleriabarra
27*5a21df39Svaleriabarra  # Attributes
28*5a21df39Svaleriabarra  _ceed = ffi.NULL
29*5a21df39Svaleriabarra  _pointer = ffi.NULL
30*5a21df39Svaleriabarra
31*5a21df39Svaleriabarra  # Representation
32*5a21df39Svaleriabarra  def __repr__(self):
33*5a21df39Svaleriabarra    return "<CeedBasis instance at " + hex(id(self)) + ">"
34*5a21df39Svaleriabarra
35*5a21df39Svaleriabarra  # String conversion for print() to stdout
36*5a21df39Svaleriabarra  def __str__(self):
37*5a21df39Svaleriabarra    """View a Basis via print()."""
38*5a21df39Svaleriabarra
39*5a21df39Svaleriabarra    # libCEED call
40*5a21df39Svaleriabarra    with tempfile.NamedTemporaryFile() as key_file:
41*5a21df39Svaleriabarra      with open(key_file.name, 'r+') as stream_file:
42*5a21df39Svaleriabarra        stream = ffi.cast("FILE *", stream_file)
43*5a21df39Svaleriabarra
44*5a21df39Svaleriabarra        lib.CeedBasisView(self._pointer[0], stream)
45*5a21df39Svaleriabarra
46*5a21df39Svaleriabarra        stream_file.seek(0)
47*5a21df39Svaleriabarra        out_string = stream_file.read()
48*5a21df39Svaleriabarra
49*5a21df39Svaleriabarra    return out_string
50*5a21df39Svaleriabarra
51*5a21df39Svaleriabarra  # Apply Basis
52*5a21df39Svaleriabarra  def apply(self, nelem, emode, u, v, tmode=NOTRANSPOSE):
53*5a21df39Svaleriabarra    """Apply basis evaluation from nodes to quadrature points or vice versa.
54*5a21df39Svaleriabarra
55*5a21df39Svaleriabarra       Args:
56*5a21df39Svaleriabarra         nelem: the number of elements to apply the basis evaluation to;
57*5a21df39Svaleriabarra                  the backend will specify the ordering in a
58*5a21df39Svaleriabarra                  BlockedElemRestriction
59*5a21df39Svaleriabarra         emode: basis evaluation mode
60*5a21df39Svaleriabarra         u: input vector
61*5a21df39Svaleriabarra         v: output vector
62*5a21df39Svaleriabarra         **tmode: CEED_NOTRANSPOSE to evaluate from nodes to quadrature
63*5a21df39Svaleriabarra                    points, CEED_TRANSPOSE to apply the transpose, mapping
64*5a21df39Svaleriabarra                    from quadrature points to nodes; default CEED_NOTRANSPOSE"""
65*5a21df39Svaleriabarra
66*5a21df39Svaleriabarra    # libCEED call
67*5a21df39Svaleriabarra    lib.CeedBasisApply(self._pointer[0], nelem, tmode, emode,
68*5a21df39Svaleriabarra                       u._pointer[0], v._pointer[0])
69*5a21df39Svaleriabarra
70*5a21df39Svaleriabarra  # Transpose a Basis
71*5a21df39Svaleriabarra  @property
72*5a21df39Svaleriabarra  def T(self):
73*5a21df39Svaleriabarra    """Transpose a Basis."""
74*5a21df39Svaleriabarra
75*5a21df39Svaleriabarra    return TransposeBasis(self)
76*5a21df39Svaleriabarra
77*5a21df39Svaleriabarra  # Transpose a Basis
78*5a21df39Svaleriabarra  @property
79*5a21df39Svaleriabarra  def transpose(self):
80*5a21df39Svaleriabarra    """Transpose a Basis."""
81*5a21df39Svaleriabarra
82*5a21df39Svaleriabarra    return TransposeBasis(self)
83*5a21df39Svaleriabarra
84*5a21df39Svaleriabarra  # Get number of nodes
85*5a21df39Svaleriabarra  def get_num_nodes(self):
86*5a21df39Svaleriabarra    """Get total number of nodes (in dim dimensions) of a Basis.
87*5a21df39Svaleriabarra
88*5a21df39Svaleriabarra       Returns:
89*5a21df39Svaleriabarra         num_nodes: total number of nodes"""
90*5a21df39Svaleriabarra
91*5a21df39Svaleriabarra    # Setup argument
92*5a21df39Svaleriabarra    p_pointer = ffi.new("CeedInt *")
93*5a21df39Svaleriabarra
94*5a21df39Svaleriabarra    # libCEED call
95*5a21df39Svaleriabarra    lib.CeedBasisGetNumNodes(self._pointer[0], p_pointer)
96*5a21df39Svaleriabarra
97*5a21df39Svaleriabarra    return p_pointer[0]
98*5a21df39Svaleriabarra
99*5a21df39Svaleriabarra  # Get number of quadrature points
100*5a21df39Svaleriabarra  def get_num_quadrature_points(self):
101*5a21df39Svaleriabarra    """Get total number of quadrature points (in dim dimensions) of a Basis.
102*5a21df39Svaleriabarra
103*5a21df39Svaleriabarra       Returns:
104*5a21df39Svaleriabarra         num_qpts: total number of quadrature points"""
105*5a21df39Svaleriabarra
106*5a21df39Svaleriabarra    # Setup argument
107*5a21df39Svaleriabarra    q_pointer = ffi.new("CeedInt *")
108*5a21df39Svaleriabarra
109*5a21df39Svaleriabarra    # libCEED call
110*5a21df39Svaleriabarra    lib.CeedBasisGetNumQuadraturePoints(self._pointer[0], q_pointer)
111*5a21df39Svaleriabarra
112*5a21df39Svaleriabarra    return q_pointer[0]
113*5a21df39Svaleriabarra
114*5a21df39Svaleriabarra  # Gauss quadrature
115*5a21df39Svaleriabarra  @staticmethod
116*5a21df39Svaleriabarra  def gauss_quadrature(q):
117*5a21df39Svaleriabarra    """Construct a Gauss-Legendre quadrature.
118*5a21df39Svaleriabarra
119*5a21df39Svaleriabarra       Args:
120*5a21df39Svaleriabarra         Q: number of quadrature points (integrates polynomials of
121*5a21df39Svaleriabarra              degree 2*Q-1 exactly)
122*5a21df39Svaleriabarra
123*5a21df39Svaleriabarra       Returns:
124*5a21df39Svaleriabarra         (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
125*5a21df39Svaleriabarra                                and array of length Q to hold the weights"""
126*5a21df39Svaleriabarra
127*5a21df39Svaleriabarra    # Setup arguments
128*5a21df39Svaleriabarra    qref1d = np.empty(q, dtype="float64")
129*5a21df39Svaleriabarra    qweight1d = np.empy(q, dtype="float64")
130*5a21df39Svaleriabarra
131*5a21df39Svaleriabarra    qref1d_pointer = ffi.new("CeedScalar *")
132*5a21df39Svaleriabarra    qref1d_pointer = ffi.cast("CeedScalar *", qref1d.__array_interface__['data'][0])
133*5a21df39Svaleriabarra
134*5a21df39Svaleriabarra    qweight1d_pointer = ffi.new("CeedScalar *")
135*5a21df39Svaleriabarra    qweight1d_pointer = ffi.cast("CeedScalar *", qweight1d.__array_interface__['data'][0])
136*5a21df39Svaleriabarra
137*5a21df39Svaleriabarra    # libCEED call
138*5a21df39Svaleriabarra    lib.CeedGaussQuadrature(q, qref1d_pointer, qweight1d_pointer)
139*5a21df39Svaleriabarra
140*5a21df39Svaleriabarra    return qref1d, qweight1d
141*5a21df39Svaleriabarra
142*5a21df39Svaleriabarra  # Lobatto quadrature
143*5a21df39Svaleriabarra  @staticmethod
144*5a21df39Svaleriabarra  def lobatto_quadrature(q):
145*5a21df39Svaleriabarra    """Construct a Gauss-Legendre-Lobatto quadrature.
146*5a21df39Svaleriabarra
147*5a21df39Svaleriabarra       Args:
148*5a21df39Svaleriabarra         q: number of quadrature points (integrates polynomials of
149*5a21df39Svaleriabarra              degree 2*Q-3 exactly)
150*5a21df39Svaleriabarra
151*5a21df39Svaleriabarra       Returns:
152*5a21df39Svaleriabarra         (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
153*5a21df39Svaleriabarra                                and array of length Q to hold the weights"""
154*5a21df39Svaleriabarra
155*5a21df39Svaleriabarra    # Setup arguments
156*5a21df39Svaleriabarra    qref1d = np.empty(q, dtype="float64")
157*5a21df39Svaleriabarra    qref1d_pointer = ffi.new("CeedScalar *")
158*5a21df39Svaleriabarra    qref1d_pointer = ffi.cast("CeedScalar *", qref1d.__array_interface__['data'][0])
159*5a21df39Svaleriabarra
160*5a21df39Svaleriabarra    qweight1d = np.empy(q, dtype="float64")
161*5a21df39Svaleriabarra    qweight1d_pointer = ffi.new("CeedScalar *")
162*5a21df39Svaleriabarra    qweight1d_pointer = ffi.cast("CeedScalar *", qweight1d.__array_interface__['data'][0])
163*5a21df39Svaleriabarra
164*5a21df39Svaleriabarra    # libCEED call
165*5a21df39Svaleriabarra    lib.CeedLobattoQuadrature(q, qref1d_pointer, qweight1d_pointer)
166*5a21df39Svaleriabarra
167*5a21df39Svaleriabarra    return qref1d, qweight1d
168*5a21df39Svaleriabarra
169*5a21df39Svaleriabarra  # QR factorization
170*5a21df39Svaleriabarra  @staticmethod
171*5a21df39Svaleriabarra  def qr_factorization(ceed, mat, tau, m, n):
172*5a21df39Svaleriabarra    """Return QR Factorization of a matrix.
173*5a21df39Svaleriabarra
174*5a21df39Svaleriabarra       Args:
175*5a21df39Svaleriabarra         ceed: Ceed context currently in use
176*5a21df39Svaleriabarra         *mat: Numpy array holding the row-major matrix to be factorized in place
177*5a21df39Svaleriabarra         *tau: Numpy array to hold the vector of lengt m of scaling factors
178*5a21df39Svaleriabarra         m: number of rows
179*5a21df39Svaleriabarra         n: numbef of columns"""
180*5a21df39Svaleriabarra
181*5a21df39Svaleriabarra    # Setup arguments
182*5a21df39Svaleriabarra    mat_pointer = ffi.new("CeedScalar *")
183*5a21df39Svaleriabarra    mat_pointer = ffi.cast("CeedScalar *", mat.__array_interface__['data'][0])
184*5a21df39Svaleriabarra
185*5a21df39Svaleriabarra    tau_pointer = ffi.new("CeedScalar *")
186*5a21df39Svaleriabarra    tau_pointer = ffi.cast("CeedScalar *", tau.__array_interface__['data'][0])
187*5a21df39Svaleriabarra
188*5a21df39Svaleriabarra    # libCEED call
189*5a21df39Svaleriabarra    lib.CeedQRFactorization(ceed._pointer[0], mat_pointer, tau_pointer, m, n)
190*5a21df39Svaleriabarra
191*5a21df39Svaleriabarra    return mat, tau
192*5a21df39Svaleriabarra
193*5a21df39Svaleriabarra  # Symmetric Schur decomposition
194*5a21df39Svaleriabarra  @staticmethod
195*5a21df39Svaleriabarra  def symmetric_schur_decomposition(ceed, mat, n):
196*5a21df39Svaleriabarra    """Return symmetric Schur decomposition of a symmetric matrix
197*5a21df39Svaleriabarra         via symmetric QR factorization.
198*5a21df39Svaleriabarra
199*5a21df39Svaleriabarra       Args:
200*5a21df39Svaleriabarra         ceed: Ceed context currently in use
201*5a21df39Svaleriabarra         *mat: Numpy array holding the row-major matrix to be factorized in place
202*5a21df39Svaleriabarra         n: number of rows/columns
203*5a21df39Svaleriabarra
204*5a21df39Svaleriabarra       Returns:
205*5a21df39Svaleriabarra         lbda: Numpy array of length n holding eigenvalues"""
206*5a21df39Svaleriabarra
207*5a21df39Svaleriabarra    # Setup arguments
208*5a21df39Svaleriabarra    mat_pointer = ffi.new("CeedScalar *")
209*5a21df39Svaleriabarra    mat_pointer = ffi.cast("CeedScalar *", mat.__array_interface__['data'][0])
210*5a21df39Svaleriabarra
211*5a21df39Svaleriabarra    lbda = np.empty(n, dtype="float64")
212*5a21df39Svaleriabarra    l_pointer = ffi.new("CeedScalar *")
213*5a21df39Svaleriabarra    l_pointer = ffi.cast("CeedScalar *", lbda.__array_interface__['data'][0])
214*5a21df39Svaleriabarra
215*5a21df39Svaleriabarra    # libCEED call
216*5a21df39Svaleriabarra    lib.CeedSymmetricSchurDecomposition(ceed._pointer[0], mat_pointer, l_pointer, n)
217*5a21df39Svaleriabarra
218*5a21df39Svaleriabarra    return lbda
219*5a21df39Svaleriabarra
220*5a21df39Svaleriabarra  # Simultaneous Diagonalization
221*5a21df39Svaleriabarra  @staticmethod
222*5a21df39Svaleriabarra  def simultaneous_diagonalization(ceed, matA, matB, n):
223*5a21df39Svaleriabarra    """Return Simultaneous Diagonalization of two matrices.
224*5a21df39Svaleriabarra
225*5a21df39Svaleriabarra       Args:
226*5a21df39Svaleriabarra         ceed: Ceed context currently in use
227*5a21df39Svaleriabarra         *matA: Numpy array holding the row-major matrix to be factorized with
228*5a21df39Svaleriabarra                  eigenvalues
229*5a21df39Svaleriabarra         *matB: Numpy array holding the row-major matrix to be factorized to identity
230*5a21df39Svaleriabarra         n: number of rows/columns
231*5a21df39Svaleriabarra
232*5a21df39Svaleriabarra       Returns:
233*5a21df39Svaleriabarra         (x, lbda): Numpy array holding the row-major orthogonal matrix and
234*5a21df39Svaleriabarra                      Numpy array holding the vector of length n of generalized
235*5a21df39Svaleriabarra                      eigenvalues"""
236*5a21df39Svaleriabarra
237*5a21df39Svaleriabarra    # Setup arguments
238*5a21df39Svaleriabarra    matA_pointer = ffi.new("CeedScalar *")
239*5a21df39Svaleriabarra    matA_pointer = ffi.cast("CeedScalar *", matA.__array_interface__['data'][0])
240*5a21df39Svaleriabarra
241*5a21df39Svaleriabarra    matB_pointer = ffi.new("CeedScalar *")
242*5a21df39Svaleriabarra    matB_pointer = ffi.cast("CeedScalar *", matB.__array_interface__['data'][0])
243*5a21df39Svaleriabarra
244*5a21df39Svaleriabarra    lbda = np.empty(n, dtype="float64")
245*5a21df39Svaleriabarra    l_pointer = ffi.new("CeedScalar *")
246*5a21df39Svaleriabarra    l_pointer = ffi.cast("CeedScalar *", lbda.__array_interface__['data'][0])
247*5a21df39Svaleriabarra
248*5a21df39Svaleriabarra    x = np.empty(n*n, dtype="float64")
249*5a21df39Svaleriabarra    x_pointer = ffi.new("CeedScalar *")
250*5a21df39Svaleriabarra    x_pointer = ffi.cast("CeedScalar *", x.__array_interface__['data'][0])
251*5a21df39Svaleriabarra
252*5a21df39Svaleriabarra    # libCEED call
253*5a21df39Svaleriabarra    lib.CeedSimultaneousDiagonalization(ceed._pointer[0], matA_pointer, matB_pointer,
254*5a21df39Svaleriabarra                                        x_pointer, l_pointer, n)
255*5a21df39Svaleriabarra
256*5a21df39Svaleriabarra    return x, lbda
257*5a21df39Svaleriabarra
258*5a21df39Svaleriabarra  # Destructor
259*5a21df39Svaleriabarra  def __del__(self):
260*5a21df39Svaleriabarra    # libCEED call
261*5a21df39Svaleriabarra    lib.CeedBasisDestroy(self._pointer)
262*5a21df39Svaleriabarra
263*5a21df39Svaleriabarra# ------------------------------------------------------------------------------
264*5a21df39Svaleriabarraclass BasisTensorH1(Basis):
265*5a21df39Svaleriabarra  """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
266*5a21df39Svaleriabarra       H^1 discretizations."""
267*5a21df39Svaleriabarra
268*5a21df39Svaleriabarra  # Constructor
269*5a21df39Svaleriabarra  def __init__(self, ceed, dim, ncomp, P1d, Q1d, interp1d, grad1d,
270*5a21df39Svaleriabarra               qref1d, qweight1d):
271*5a21df39Svaleriabarra
272*5a21df39Svaleriabarra    # Setup arguments
273*5a21df39Svaleriabarra    self._pointer = ffi.new("CeedBasis *")
274*5a21df39Svaleriabarra
275*5a21df39Svaleriabarra    self._ceed = ceed
276*5a21df39Svaleriabarra
277*5a21df39Svaleriabarra    interp1d_pointer = ffi.new("CeedScalar *")
278*5a21df39Svaleriabarra    interp1d_pointer = ffi.cast("CeedScalar *", interp1d.__array_interface__['data'][0])
279*5a21df39Svaleriabarra
280*5a21df39Svaleriabarra    grad1d_pointer = ffi.new("CeedScalar *")
281*5a21df39Svaleriabarra    grad1d_pointer = ffi.cast("CeedScalar *", grad1d.__array_interface__['data'][0])
282*5a21df39Svaleriabarra
283*5a21df39Svaleriabarra    qref1d_pointer = ffi.new("CeedScalar *")
284*5a21df39Svaleriabarra    qref1d_pointer = ffi.cast("CeedScalar *", qref1d.__array_interface__['data'][0])
285*5a21df39Svaleriabarra
286*5a21df39Svaleriabarra    qweight1d_pointer = ffi.new("CeedScalar *")
287*5a21df39Svaleriabarra    qweight1d_pointer = ffi.cast("CeedScalar *", qweight1d.__array_interface__['data'][0])
288*5a21df39Svaleriabarra
289*5a21df39Svaleriabarra    # libCEED call
290*5a21df39Svaleriabarra    lib.CeedBasisCreateTensorH1(self._ceed._pointer[0], dim, ncomp, P1d, Q1d,
291*5a21df39Svaleriabarra                                interp1d_pointer, grad1d_pointer, qref1d_pointer,
292*5a21df39Svaleriabarra                                qweight1d_pointer, self._pointer)
293*5a21df39Svaleriabarra
294*5a21df39Svaleriabarra# ------------------------------------------------------------------------------
295*5a21df39Svaleriabarraclass BasisTensorH1Lagrange(Basis):
296*5a21df39Svaleriabarra  """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange basis
297*5a21df39Svaleriabarra       objects for H^1 discretizations."""
298*5a21df39Svaleriabarra
299*5a21df39Svaleriabarra  # Constructor
300*5a21df39Svaleriabarra  def __init__(self, ceed, dim, ncomp, P, Q, qmode):
301*5a21df39Svaleriabarra
302*5a21df39Svaleriabarra    # Setup arguments
303*5a21df39Svaleriabarra    self._pointer = ffi.new("CeedBasis *")
304*5a21df39Svaleriabarra
305*5a21df39Svaleriabarra    self._ceed = ceed
306*5a21df39Svaleriabarra
307*5a21df39Svaleriabarra    # libCEED call
308*5a21df39Svaleriabarra    lib.CeedBasisCreateTensorH1Lagrange(self._ceed._pointer[0], dim, ncomp, P,
309*5a21df39Svaleriabarra                                        Q, qmode, self._pointer)
310*5a21df39Svaleriabarra
311*5a21df39Svaleriabarra# ------------------------------------------------------------------------------
312*5a21df39Svaleriabarraclass BasisH1(Basis):
313*5a21df39Svaleriabarra  """Ceed H1 Basis: finite element non tensor-product basis for H^1 discretizations."""
314*5a21df39Svaleriabarra
315*5a21df39Svaleriabarra  # Constructor
316*5a21df39Svaleriabarra  def __init__(self, ceed, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight):
317*5a21df39Svaleriabarra
318*5a21df39Svaleriabarra    # Setup arguments
319*5a21df39Svaleriabarra    self._pointer = ffi.new("CeedBasis *")
320*5a21df39Svaleriabarra
321*5a21df39Svaleriabarra    self._ceed = ceed
322*5a21df39Svaleriabarra
323*5a21df39Svaleriabarra    interp_pointer = ffi.new("CeedScalar *")
324*5a21df39Svaleriabarra    interp_pointer = ffi.cast("CeedScalar *", interp.__array_interface__['data'][0])
325*5a21df39Svaleriabarra
326*5a21df39Svaleriabarra    grad_pointer = ffi.new("CeedScalar *")
327*5a21df39Svaleriabarra    grad_pointer = ffi.cast("CeedScalar *", grad.__array_interface__['data'][0])
328*5a21df39Svaleriabarra
329*5a21df39Svaleriabarra    qref_pointer = ffi.new("CeedScalar *")
330*5a21df39Svaleriabarra    qref_pointer = ffi.cast("CeedScalar *", qref.__array_interface__['data'][0])
331*5a21df39Svaleriabarra
332*5a21df39Svaleriabarra    qweight_pointer = ffi.new("CeedScalar *")
333*5a21df39Svaleriabarra    qweight_pointer = ffi.cast("CeedScalar *", qweight.__array_interface__['data'][0])
334*5a21df39Svaleriabarra
335*5a21df39Svaleriabarra    # libCEED call
336*5a21df39Svaleriabarra    lib.CeedBasisCreateH1(self._ceed._pointer[0], topo, ncomp, nnodes, nqpts,
337*5a21df39Svaleriabarra                          interp_pointer, grad_pointer, qref_pointer,
338*5a21df39Svaleriabarra                          qweight_pointer, self._pointer)
339*5a21df39Svaleriabarra
340*5a21df39Svaleriabarra# ------------------------------------------------------------------------------
341*5a21df39Svaleriabarraclass TransposeBasis():
342*5a21df39Svaleriabarra  """Transpose Ceed Basis: transpose of finite element tensor-product basis objects."""
343*5a21df39Svaleriabarra
344*5a21df39Svaleriabarra  # Attributes
345*5a21df39Svaleriabarra  _basis = None
346*5a21df39Svaleriabarra
347*5a21df39Svaleriabarra  # Constructor
348*5a21df39Svaleriabarra  def __init__(self, basis):
349*5a21df39Svaleriabarra
350*5a21df39Svaleriabarra    # Reference basis
351*5a21df39Svaleriabarra    self._basis = basis
352*5a21df39Svaleriabarra
353*5a21df39Svaleriabarra  # Representation
354*5a21df39Svaleriabarra  def __repr__(self):
355*5a21df39Svaleriabarra    return "<Transpose CeedBasis instance at " + hex(id(self)) + ">"
356*5a21df39Svaleriabarra
357*5a21df39Svaleriabarra  # Apply Transpose Basis
358*5a21df39Svaleriabarra  def apply(self, nelem, emode, u, v):
359*5a21df39Svaleriabarra    """Apply basis evaluation from quadrature points to nodes.
360*5a21df39Svaleriabarra
361*5a21df39Svaleriabarra       Args:
362*5a21df39Svaleriabarra         nelem: the number of elements to apply the basis evaluation to;
363*5a21df39Svaleriabarra                  the backend will specify the ordering in a
364*5a21df39Svaleriabarra                  Blocked ElemRestriction
365*5a21df39Svaleriabarra         **emode: basis evaluation mode
366*5a21df39Svaleriabarra         u: input vector
367*5a21df39Svaleriabarra         v: output vector"""
368*5a21df39Svaleriabarra
369*5a21df39Svaleriabarra    # libCEED call
370*5a21df39Svaleriabarra    self._basis.apply(nelem, emode, u, v, tmode=TRANSPOSE)
371*5a21df39Svaleriabarra
372*5a21df39Svaleriabarra# ------------------------------------------------------------------------------
373