xref: /libCEED/python/ceed.py (revision 80a9ef0545a39c00cdcaab1ca26f8053604f3120)
139b2de37Sjeremylt# Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
239b2de37Sjeremylt# the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
339b2de37Sjeremylt# reserved. See files LICENSE and NOTICE for details.
439b2de37Sjeremylt#
539b2de37Sjeremylt# This file is part of CEED, a collection of benchmarks, miniapps, software
639b2de37Sjeremylt# libraries and APIs for efficient high-order finite element and spectral
739b2de37Sjeremylt# element discretizations for exascale applications. For more information and
839b2de37Sjeremylt# source code availability see http://github.com/ceed.
939b2de37Sjeremylt#
1039b2de37Sjeremylt# The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
1139b2de37Sjeremylt# a collaborative effort of two U.S. Department of Energy organizations (Office
1239b2de37Sjeremylt# of Science and the National Nuclear Security Administration) responsible for
1339b2de37Sjeremylt# the planning and preparation of a capable exascale ecosystem, including
1439b2de37Sjeremylt# software, applications, hardware, advanced system engineering and early
1539b2de37Sjeremylt# testbed platforms, in support of the nation's exascale computing imperative.
1639b2de37Sjeremylt
1739b2de37Sjeremyltfrom _ceed_cffi import ffi, lib
1839b2de37Sjeremyltimport sys
19477729cfSJeremy L Thompsonimport os
2039b2de37Sjeremyltimport io
21e15f9bd0SJeremy L Thompsonimport numpy as np
220a0da059Sjeremyltimport tempfile
2339b2de37Sjeremyltfrom abc import ABC
2439b2de37Sjeremyltfrom .ceed_vector import Vector
2539b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1
2669a53589Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, StridedElemRestriction, BlockedElemRestriction, BlockedStridedElemRestriction
2739b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction
28777ff853SJeremy L Thompsonfrom .ceed_qfunctioncontext import QFunctionContext
2939b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator
3039b2de37Sjeremyltfrom .ceed_constants import *
3139b2de37Sjeremylt
3239b2de37Sjeremylt# ------------------------------------------------------------------------------
337a7b0fa3SJed Brown
347a7b0fa3SJed Brown
3539b2de37Sjeremyltclass Ceed():
3639b2de37Sjeremylt    """Ceed: core components."""
3739b2de37Sjeremylt
3839b2de37Sjeremylt    # Constructor
3974b533adSJed Brown    def __init__(self, resource="/cpu/self", on_error="store"):
4039b2de37Sjeremylt        # libCEED object
4139b2de37Sjeremylt        self._pointer = ffi.new("Ceed *")
4239b2de37Sjeremylt
4339b2de37Sjeremylt        # libCEED call
4439b2de37Sjeremylt        resourceAscii = ffi.new("char[]", resource.encode("ascii"))
45477729cfSJeremy L Thompson        os.environ["CEED_ERROR_HANDLER"] = "return"
46477729cfSJeremy L Thompson        err_code = lib.CeedInit(resourceAscii, self._pointer)
47477729cfSJeremy L Thompson        if err_code:
48477729cfSJeremy L Thompson            raise Exception("Error initializing backend resource: " + resource)
4974b533adSJed Brown        error_handlers = dict(
5074b533adSJed Brown            store="CeedErrorStore",
5174b533adSJed Brown            abort="CeedErrorAbort",
5274b533adSJed Brown        )
53477729cfSJeremy L Thompson        lib.CeedSetErrorHandler(
54477729cfSJeremy L Thompson            self._pointer[0], ffi.addressof(
5574b533adSJed Brown                lib, error_handlers[on_error]))
5639b2de37Sjeremylt
5739b2de37Sjeremylt    # Representation
5839b2de37Sjeremylt    def __repr__(self):
5939b2de37Sjeremylt        return "<Ceed instance at " + hex(id(self)) + ">"
6039b2de37Sjeremylt
610a0da059Sjeremylt    # String conversion for print() to stdout
620a0da059Sjeremylt    def __str__(self):
630a0da059Sjeremylt        """View a Ceed via print()."""
640a0da059Sjeremylt
650a0da059Sjeremylt        # libCEED call
660a0da059Sjeremylt        with tempfile.NamedTemporaryFile() as key_file:
670a0da059Sjeremylt            with open(key_file.name, 'r+') as stream_file:
680a0da059Sjeremylt                stream = ffi.cast("FILE *", stream_file)
690a0da059Sjeremylt
70477729cfSJeremy L Thompson                err_code = lib.CeedView(self._pointer[0], stream)
71477729cfSJeremy L Thompson                self._check_error(err_code)
720a0da059Sjeremylt
730a0da059Sjeremylt                stream_file.seek(0)
740a0da059Sjeremylt                out_string = stream_file.read()
750a0da059Sjeremylt
760a0da059Sjeremylt        return out_string
770a0da059Sjeremylt
78477729cfSJeremy L Thompson    # Error handler
79477729cfSJeremy L Thompson    def _check_error(self, err_code):
80477729cfSJeremy L Thompson        """Check return code and retrieve error message for non-zero code"""
81e15f9bd0SJeremy L Thompson        if (err_code != lib.CEED_ERROR_SUCCESS):
82477729cfSJeremy L Thompson            message = ffi.new("char **")
83477729cfSJeremy L Thompson            lib.CeedGetErrorMessage(self._pointer[0], message)
84477729cfSJeremy L Thompson            raise Exception(ffi.string(message[0]).decode("UTF-8"))
85477729cfSJeremy L Thompson
8639b2de37Sjeremylt    # Get Resource
8739b2de37Sjeremylt    def get_resource(self):
8839b2de37Sjeremylt        """Get the full resource name for a Ceed context.
8939b2de37Sjeremylt
9039b2de37Sjeremylt           Returns:
9139b2de37Sjeremylt             resource: resource name"""
9239b2de37Sjeremylt
9339b2de37Sjeremylt        # libCEED call
9439b2de37Sjeremylt        resource = ffi.new("char **")
95477729cfSJeremy L Thompson        err_code = lib.CeedGetResource(self._pointer[0], resource)
96477729cfSJeremy L Thompson        self._check_error(err_code)
9739b2de37Sjeremylt
9839b2de37Sjeremylt        return ffi.string(resource[0]).decode("UTF-8")
9939b2de37Sjeremylt
10039b2de37Sjeremylt    # Preferred MemType
10139b2de37Sjeremylt    def get_preferred_memtype(self):
10239b2de37Sjeremylt        """Return Ceed preferred memory type.
10339b2de37Sjeremylt
10439b2de37Sjeremylt           Returns:
10539b2de37Sjeremylt             memtype: Ceed preferred memory type"""
10639b2de37Sjeremylt
10739b2de37Sjeremylt        # libCEED call
10839b2de37Sjeremylt        memtype = ffi.new("CeedMemType *", MEM_HOST)
109477729cfSJeremy L Thompson        err_code = lib.CeedGetPreferredMemType(self._pointer[0], memtype)
110477729cfSJeremy L Thompson        self._check_error(err_code)
11139b2de37Sjeremylt
11239b2de37Sjeremylt        return memtype[0]
11339b2de37Sjeremylt
114*80a9ef05SNatalie Beams    # Convenience function to get CeedScalar type
115*80a9ef05SNatalie Beams    def scalar_type(self):
116*80a9ef05SNatalie Beams        """Return dtype string for CeedScalar.
117*80a9ef05SNatalie Beams
118*80a9ef05SNatalie Beams           Returns:
119*80a9ef05SNatalie Beams             np_dtype: String naming numpy data type corresponding to CeedScalar"""
120*80a9ef05SNatalie Beams
121*80a9ef05SNatalie Beams        return scalar_types[lib.CEED_SCALAR_TYPE]
122*80a9ef05SNatalie Beams
123e15f9bd0SJeremy L Thompson    # --- Basis utility functions ---
124e15f9bd0SJeremy L Thompson
125e15f9bd0SJeremy L Thompson    # Gauss quadrature
126e15f9bd0SJeremy L Thompson    def gauss_quadrature(self, q):
127e15f9bd0SJeremy L Thompson        """Construct a Gauss-Legendre quadrature.
128e15f9bd0SJeremy L Thompson
129e15f9bd0SJeremy L Thompson           Args:
130e15f9bd0SJeremy L Thompson             Q: number of quadrature points (integrates polynomials of
131e15f9bd0SJeremy L Thompson                  degree 2*Q-1 exactly)
132e15f9bd0SJeremy L Thompson
133e15f9bd0SJeremy L Thompson           Returns:
134e15f9bd0SJeremy L Thompson             (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
135e15f9bd0SJeremy L Thompson                                    and array of length Q to hold the weights"""
136e15f9bd0SJeremy L Thompson
137e15f9bd0SJeremy L Thompson        # Setup arguments
138*80a9ef05SNatalie Beams        qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
139*80a9ef05SNatalie Beams        qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
140e15f9bd0SJeremy L Thompson
141e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.new("CeedScalar *")
142e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.cast(
143e15f9bd0SJeremy L Thompson            "CeedScalar *",
144e15f9bd0SJeremy L Thompson            qref1d.__array_interface__['data'][0])
145e15f9bd0SJeremy L Thompson
146e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.new("CeedScalar *")
147e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.cast(
148e15f9bd0SJeremy L Thompson            "CeedScalar *",
149e15f9bd0SJeremy L Thompson            qweight1d.__array_interface__['data'][0])
150e15f9bd0SJeremy L Thompson
151e15f9bd0SJeremy L Thompson        # libCEED call
152e15f9bd0SJeremy L Thompson        err_code = lib.CeedGaussQuadrature(q, qref1d_pointer, qweight1d_pointer)
153e15f9bd0SJeremy L Thompson        self._check_error(err_code)
154e15f9bd0SJeremy L Thompson
155e15f9bd0SJeremy L Thompson        return qref1d, qweight1d
156e15f9bd0SJeremy L Thompson
157e15f9bd0SJeremy L Thompson    # Lobatto quadrature
158e15f9bd0SJeremy L Thompson    def lobatto_quadrature(self, q):
159e15f9bd0SJeremy L Thompson        """Construct a Gauss-Legendre-Lobatto quadrature.
160e15f9bd0SJeremy L Thompson
161e15f9bd0SJeremy L Thompson           Args:
162e15f9bd0SJeremy L Thompson             q: number of quadrature points (integrates polynomials of
163e15f9bd0SJeremy L Thompson                  degree 2*Q-3 exactly)
164e15f9bd0SJeremy L Thompson
165e15f9bd0SJeremy L Thompson           Returns:
166e15f9bd0SJeremy L Thompson             (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
167e15f9bd0SJeremy L Thompson                                    and array of length Q to hold the weights"""
168e15f9bd0SJeremy L Thompson
169e15f9bd0SJeremy L Thompson        # Setup arguments
170*80a9ef05SNatalie Beams        qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
171e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.new("CeedScalar *")
172e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.cast(
173e15f9bd0SJeremy L Thompson            "CeedScalar *",
174e15f9bd0SJeremy L Thompson            qref1d.__array_interface__['data'][0])
175e15f9bd0SJeremy L Thompson
176*80a9ef05SNatalie Beams        qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
177e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.new("CeedScalar *")
178e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.cast(
179e15f9bd0SJeremy L Thompson            "CeedScalar *",
180e15f9bd0SJeremy L Thompson            qweight1d.__array_interface__['data'][0])
181e15f9bd0SJeremy L Thompson
182e15f9bd0SJeremy L Thompson        # libCEED call
183e15f9bd0SJeremy L Thompson        err_code = lib.CeedLobattoQuadrature(
184e15f9bd0SJeremy L Thompson            q, qref1d_pointer, qweight1d_pointer)
185e15f9bd0SJeremy L Thompson        self._check_error(err_code)
186e15f9bd0SJeremy L Thompson
187e15f9bd0SJeremy L Thompson        return qref1d, qweight1d
188e15f9bd0SJeremy L Thompson
189e15f9bd0SJeremy L Thompson    # QR factorization
190e15f9bd0SJeremy L Thompson    def qr_factorization(self, mat, tau, m, n):
191e15f9bd0SJeremy L Thompson        """Return QR Factorization of a matrix.
192e15f9bd0SJeremy L Thompson
193e15f9bd0SJeremy L Thompson           Args:
194e15f9bd0SJeremy L Thompson             ceed: Ceed context currently in use
195e15f9bd0SJeremy L Thompson             *mat: Numpy array holding the row-major matrix to be factorized in place
196e15f9bd0SJeremy L Thompson             *tau: Numpy array to hold the vector of lengt m of scaling factors
197e15f9bd0SJeremy L Thompson             m: number of rows
198e15f9bd0SJeremy L Thompson             n: numbef of columns"""
199e15f9bd0SJeremy L Thompson
200e15f9bd0SJeremy L Thompson        # Setup arguments
201e15f9bd0SJeremy L Thompson        mat_pointer = ffi.new("CeedScalar *")
202e15f9bd0SJeremy L Thompson        mat_pointer = ffi.cast(
203e15f9bd0SJeremy L Thompson            "CeedScalar *",
204e15f9bd0SJeremy L Thompson            mat.__array_interface__['data'][0])
205e15f9bd0SJeremy L Thompson
206e15f9bd0SJeremy L Thompson        tau_pointer = ffi.new("CeedScalar *")
207e15f9bd0SJeremy L Thompson        tau_pointer = ffi.cast(
208e15f9bd0SJeremy L Thompson            "CeedScalar *",
209e15f9bd0SJeremy L Thompson            tau.__array_interface__['data'][0])
210e15f9bd0SJeremy L Thompson
211e15f9bd0SJeremy L Thompson        # libCEED call
212e15f9bd0SJeremy L Thompson        err_code = lib.CeedQRFactorization(
213e15f9bd0SJeremy L Thompson            self._pointer[0], mat_pointer, tau_pointer, m, n)
214e15f9bd0SJeremy L Thompson        self._check_error(err_code)
215e15f9bd0SJeremy L Thompson
216e15f9bd0SJeremy L Thompson        return mat, tau
217e15f9bd0SJeremy L Thompson
218e15f9bd0SJeremy L Thompson    # Symmetric Schur decomposition
219e15f9bd0SJeremy L Thompson    def symmetric_schur_decomposition(self, mat, n):
220e15f9bd0SJeremy L Thompson        """Return symmetric Schur decomposition of a symmetric matrix
221e15f9bd0SJeremy L Thompson             via symmetric QR factorization.
222e15f9bd0SJeremy L Thompson
223e15f9bd0SJeremy L Thompson           Args:
224e15f9bd0SJeremy L Thompson             ceed: Ceed context currently in use
225e15f9bd0SJeremy L Thompson             *mat: Numpy array holding the row-major matrix to be factorized in place
226e15f9bd0SJeremy L Thompson             n: number of rows/columns
227e15f9bd0SJeremy L Thompson
228e15f9bd0SJeremy L Thompson           Returns:
229e15f9bd0SJeremy L Thompson             lbda: Numpy array of length n holding eigenvalues"""
230e15f9bd0SJeremy L Thompson
231e15f9bd0SJeremy L Thompson        # Setup arguments
232e15f9bd0SJeremy L Thompson        mat_pointer = ffi.new("CeedScalar *")
233e15f9bd0SJeremy L Thompson        mat_pointer = ffi.cast(
234e15f9bd0SJeremy L Thompson            "CeedScalar *",
235e15f9bd0SJeremy L Thompson            mat.__array_interface__['data'][0])
236e15f9bd0SJeremy L Thompson
237*80a9ef05SNatalie Beams        lbda = np.empty(n, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
238e15f9bd0SJeremy L Thompson        l_pointer = ffi.new("CeedScalar *")
239e15f9bd0SJeremy L Thompson        l_pointer = ffi.cast(
240e15f9bd0SJeremy L Thompson            "CeedScalar *",
241e15f9bd0SJeremy L Thompson            lbda.__array_interface__['data'][0])
242e15f9bd0SJeremy L Thompson
243e15f9bd0SJeremy L Thompson        # libCEED call
244e15f9bd0SJeremy L Thompson        err_code = lib.CeedSymmetricSchurDecomposition(
245e15f9bd0SJeremy L Thompson            self._pointer[0], mat_pointer, l_pointer, n)
246e15f9bd0SJeremy L Thompson        self._check_error(err_code)
247e15f9bd0SJeremy L Thompson
248e15f9bd0SJeremy L Thompson        return lbda
249e15f9bd0SJeremy L Thompson
250e15f9bd0SJeremy L Thompson    # Simultaneous Diagonalization
251e15f9bd0SJeremy L Thompson    def simultaneous_diagonalization(self, matA, matB, n):
252e15f9bd0SJeremy L Thompson        """Return Simultaneous Diagonalization of two matrices.
253e15f9bd0SJeremy L Thompson
254e15f9bd0SJeremy L Thompson           Args:
255e15f9bd0SJeremy L Thompson             ceed: Ceed context currently in use
256e15f9bd0SJeremy L Thompson             *matA: Numpy array holding the row-major matrix to be factorized with
257e15f9bd0SJeremy L Thompson                      eigenvalues
258e15f9bd0SJeremy L Thompson             *matB: Numpy array holding the row-major matrix to be factorized to identity
259e15f9bd0SJeremy L Thompson             n: number of rows/columns
260e15f9bd0SJeremy L Thompson
261e15f9bd0SJeremy L Thompson           Returns:
262e15f9bd0SJeremy L Thompson             (x, lbda): Numpy array holding the row-major orthogonal matrix and
263e15f9bd0SJeremy L Thompson                          Numpy array holding the vector of length n of generalized
264e15f9bd0SJeremy L Thompson                          eigenvalues"""
265e15f9bd0SJeremy L Thompson
266e15f9bd0SJeremy L Thompson        # Setup arguments
267e15f9bd0SJeremy L Thompson        matA_pointer = ffi.new("CeedScalar *")
268e15f9bd0SJeremy L Thompson        matA_pointer = ffi.cast(
269e15f9bd0SJeremy L Thompson            "CeedScalar *",
270e15f9bd0SJeremy L Thompson            matA.__array_interface__['data'][0])
271e15f9bd0SJeremy L Thompson
272e15f9bd0SJeremy L Thompson        matB_pointer = ffi.new("CeedScalar *")
273e15f9bd0SJeremy L Thompson        matB_pointer = ffi.cast(
274e15f9bd0SJeremy L Thompson            "CeedScalar *",
275e15f9bd0SJeremy L Thompson            matB.__array_interface__['data'][0])
276e15f9bd0SJeremy L Thompson
277*80a9ef05SNatalie Beams        lbda = np.empty(n, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
278e15f9bd0SJeremy L Thompson        l_pointer = ffi.new("CeedScalar *")
279e15f9bd0SJeremy L Thompson        l_pointer = ffi.cast(
280e15f9bd0SJeremy L Thompson            "CeedScalar *",
281e15f9bd0SJeremy L Thompson            lbda.__array_interface__['data'][0])
282e15f9bd0SJeremy L Thompson
283*80a9ef05SNatalie Beams        x = np.empty(n * n, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
284e15f9bd0SJeremy L Thompson        x_pointer = ffi.new("CeedScalar *")
285e15f9bd0SJeremy L Thompson        x_pointer = ffi.cast("CeedScalar *", x.__array_interface__['data'][0])
286e15f9bd0SJeremy L Thompson
287e15f9bd0SJeremy L Thompson        # libCEED call
288e15f9bd0SJeremy L Thompson        err_code = lib.CeedSimultaneousDiagonalization(self._pointer[0], matA_pointer, matB_pointer,
289e15f9bd0SJeremy L Thompson                                                       x_pointer, l_pointer, n)
290e15f9bd0SJeremy L Thompson        self._check_error(err_code)
291e15f9bd0SJeremy L Thompson
292e15f9bd0SJeremy L Thompson        return x, lbda
293e15f9bd0SJeremy L Thompson
294e15f9bd0SJeremy L Thompson    # --- libCEED objects ---
295e15f9bd0SJeremy L Thompson
29639b2de37Sjeremylt    # CeedVector
29739b2de37Sjeremylt    def Vector(self, size):
29839b2de37Sjeremylt        """Ceed Vector: storing and manipulating vectors.
29939b2de37Sjeremylt
30039b2de37Sjeremylt           Args:
30139b2de37Sjeremylt             size: length of vector
30239b2de37Sjeremylt
30339b2de37Sjeremylt           Returns:
30439b2de37Sjeremylt             vector: Ceed Vector"""
30539b2de37Sjeremylt
30639b2de37Sjeremylt        return Vector(self, size)
30739b2de37Sjeremylt
30839b2de37Sjeremylt    # CeedElemRestriction
309d979a051Sjeremylt    def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets,
310d979a051Sjeremylt                        memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
31139b2de37Sjeremylt        """Ceed ElemRestriction: restriction from local vectors to elements.
31239b2de37Sjeremylt
31339b2de37Sjeremylt           Args:
314d979a051Sjeremylt             nelem: number of elements described by the restriction
31539b2de37Sjeremylt             elemsize: size (number of nodes) per element
316d979a051Sjeremylt             ncomp: number of field components per interpolation node
317d979a051Sjeremylt                      (1 for scalar fields)
318d979a051Sjeremylt             compstride: Stride between components for the same L-vector "node".
319d979a051Sjeremylt                           Data for node i, component k can be found in the
320d979a051Sjeremylt                           L-vector at index [offsets[i] + k*compstride].
321d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
322d979a051Sjeremylt                       the elements and fields given by this restriction.
323d979a051Sjeremylt             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
324d979a051Sjeremylt                         ordered list of the offsets (into the input Ceed Vector)
32539b2de37Sjeremylt                         for the unknowns corresponding to element i, where
326d979a051Sjeremylt                         0 <= i < nelem. All offsets must be in the range
327d979a051Sjeremylt                         [0, lsize - 1].
328d979a051Sjeremylt             **memtype: memory type of the offsets array, default CEED_MEM_HOST
329d979a051Sjeremylt             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
33039b2de37Sjeremylt
33139b2de37Sjeremylt           Returns:
33239b2de37Sjeremylt             elemrestriction: Ceed ElemRestiction"""
33339b2de37Sjeremylt
334d979a051Sjeremylt        return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize,
335d979a051Sjeremylt                               offsets, memtype=memtype, cmode=cmode)
33639b2de37Sjeremylt
337d979a051Sjeremylt    def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides):
338d979a051Sjeremylt        """Ceed Identity ElemRestriction: strided restriction from local vectors
339d979a051Sjeremylt             to elements.
34039b2de37Sjeremylt
34139b2de37Sjeremylt           Args:
342d979a051Sjeremylt             nelem: number of elements described by the restriction
34339b2de37Sjeremylt             elemsize: size (number of nodes) per element
344a8d32208Sjeremylt             ncomp: number of field components per interpolation node
345a8d32208Sjeremylt                      (1 for scalar fields)
346d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
347d979a051Sjeremylt                      the elements and fields given by this restriction.
34869a53589Sjeremylt             *strides: Array for strides between [nodes, components, elements].
34969a53589Sjeremylt                         The data for node i, component j, element k in the
35069a53589Sjeremylt                         L-vector is given by
35169a53589Sjeremylt                           i*strides[0] + j*strides[1] + k*strides[2]
35239b2de37Sjeremylt
35339b2de37Sjeremylt           Returns:
35469a53589Sjeremylt             elemrestriction: Ceed Strided ElemRestiction"""
35539b2de37Sjeremylt
3567a7b0fa3SJed Brown        return StridedElemRestriction(
357d979a051Sjeremylt            self, nelem, elemsize, ncomp, lsize, strides)
35839b2de37Sjeremylt
359d979a051Sjeremylt    def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride,
360d979a051Sjeremylt                               lsize, offsets, memtype=lib.CEED_MEM_HOST,
361d979a051Sjeremylt                               cmode=lib.CEED_COPY_VALUES):
362d979a051Sjeremylt        """Ceed Blocked ElemRestriction: blocked restriction from local vectors
363d979a051Sjeremylt             to elements.
36439b2de37Sjeremylt
36539b2de37Sjeremylt           Args:
366d979a051Sjeremylt             nelem: number of elements described by the restriction
36739b2de37Sjeremylt             elemsize: size (number of nodes) per element
36839b2de37Sjeremylt             blksize: number of elements in a block
369d979a051Sjeremylt             ncomp: number of field components per interpolation node
370d979a051Sjeremylt                       (1 for scalar fields)
371d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
372d979a051Sjeremylt                      the elements and fields given by this restriction.
373d979a051Sjeremylt             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
374d979a051Sjeremylt                         ordered list of the offsets (into the input Ceed Vector)
37539b2de37Sjeremylt                         for the unknowns corresponding to element i, where
376d979a051Sjeremylt                         0 <= i < nelem. All offsets must be in the range
377d979a051Sjeremylt                         [0, lsize - 1]. The backend will permute and pad this
37839b2de37Sjeremylt                         array to the desired ordering for the blocksize, which is
37939b2de37Sjeremylt                         typically given by the backend. The default reordering is
38039b2de37Sjeremylt                         to interlace elements.
381d979a051Sjeremylt             **memtype: memory type of the offsets array, default CEED_MEM_HOST
382d979a051Sjeremylt             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
38339b2de37Sjeremylt
38439b2de37Sjeremylt           Returns:
38539b2de37Sjeremylt             elemrestriction: Ceed Blocked ElemRestiction"""
38639b2de37Sjeremylt
387d979a051Sjeremylt        return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp,
388d979a051Sjeremylt                                      compstride, lsize, offsets,
389d979a051Sjeremylt                                      memtype=memtype, cmode=cmode)
39039b2de37Sjeremylt
391d979a051Sjeremylt    def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp,
392d979a051Sjeremylt                                      lsize, strides):
393d979a051Sjeremylt        """Ceed Blocked Strided ElemRestriction: blocked and strided restriction
394d979a051Sjeremylt             from local vectors to elements.
39569a53589Sjeremylt
39669a53589Sjeremylt           Args:
397d979a051Sjeremylt             nelem: number of elements described in the restriction
39869a53589Sjeremylt             elemsize: size (number of nodes) per element
39969a53589Sjeremylt             blksize: number of elements in a block
40069a53589Sjeremylt             ncomp: number of field components per interpolation node
40169a53589Sjeremylt                      (1 for scalar fields)
402d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
403d979a051Sjeremylt                      the elements and fields given by this restriction.
40469a53589Sjeremylt             *strides: Array for strides between [nodes, components, elements].
40569a53589Sjeremylt                         The data for node i, component j, element k in the
40669a53589Sjeremylt                         L-vector is given by
40769a53589Sjeremylt                           i*strides[0] + j*strides[1] + k*strides[2]
40869a53589Sjeremylt
40969a53589Sjeremylt           Returns:
41069a53589Sjeremylt             elemrestriction: Ceed Strided ElemRestiction"""
41169a53589Sjeremylt
41269a53589Sjeremylt        return BlockedStridedElemRestriction(self, nelem, elemsize, blksize,
413d979a051Sjeremylt                                             ncomp, lsize, strides)
41469a53589Sjeremylt
41539b2de37Sjeremylt    # CeedBasis
41639b2de37Sjeremylt    def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
41739b2de37Sjeremylt                      qref1d, qweight1d):
41839b2de37Sjeremylt        """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
41939b2de37Sjeremylt             H^1 discretizations.
42039b2de37Sjeremylt
42139b2de37Sjeremylt           Args:
42239b2de37Sjeremylt             dim: topological dimension
42339b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
42439b2de37Sjeremylt             P1d: number of nodes in one dimension
42539b2de37Sjeremylt             Q1d: number of quadrature points in one dimension
426d979a051Sjeremylt             *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix
427d979a051Sjeremylt                          expressing the values of nodal basis functions at
428d979a051Sjeremylt                          quadrature points
429d979a051Sjeremylt             *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix
430d979a051Sjeremylt                        expressing the derivatives of nodal basis functions at
431d979a051Sjeremylt                        quadrature points
432d979a051Sjeremylt             *qref1d: Numpy array of length Q1d holding the locations of
433d979a051Sjeremylt                        quadrature points on the 1D reference element [-1, 1]
434d979a051Sjeremylt             *qweight1d: Numpy array of length Q1d holding the quadrature
435d979a051Sjeremylt                           weights on the reference element
43639b2de37Sjeremylt
43739b2de37Sjeremylt           Returns:
43839b2de37Sjeremylt             basis: Ceed Basis"""
43939b2de37Sjeremylt
44039b2de37Sjeremylt        return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
44139b2de37Sjeremylt                             qref1d, qweight1d)
44239b2de37Sjeremylt
44339b2de37Sjeremylt    def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode):
444d979a051Sjeremylt        """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange
445d979a051Sjeremylt             basis objects for H^1 discretizations.
44639b2de37Sjeremylt
44739b2de37Sjeremylt           Args:
44839b2de37Sjeremylt             dim: topological dimension
44939b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
45039b2de37Sjeremylt             P: number of Gauss-Lobatto nodes in one dimension.  The
45139b2de37Sjeremylt                  polynomial degree of the resulting Q_k element is k=P-1.
45239b2de37Sjeremylt             Q: number of quadrature points in one dimension
45339b2de37Sjeremylt             qmode: distribution of the Q quadrature points (affects order of
45439b2de37Sjeremylt                      accuracy for the quadrature)
45539b2de37Sjeremylt
45639b2de37Sjeremylt           Returns:
45739b2de37Sjeremylt             basis: Ceed Basis"""
45839b2de37Sjeremylt
45939b2de37Sjeremylt        return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode)
46039b2de37Sjeremylt
46139b2de37Sjeremylt    def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight):
462d979a051Sjeremylt        """Ceed H1 Basis: finite element non tensor-product basis for H^1
463d979a051Sjeremylt             discretizations.
46439b2de37Sjeremylt
46539b2de37Sjeremylt           Args:
46639b2de37Sjeremylt             topo: topology of the element, e.g. hypercube, simplex, etc
46739b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
46839b2de37Sjeremylt             nnodes: total number of nodes
46939b2de37Sjeremylt             nqpts: total number of quadrature points
470d979a051Sjeremylt             *interp: Numpy array holding the row-major (nqpts * nnodes) matrix
471d979a051Sjeremylt                       expressing the values of nodal basis functions at
47239b2de37Sjeremylt                       quadrature points
473d979a051Sjeremylt             *grad: Numpy array holding the row-major (nqpts * dim * nnodes)
474d979a051Sjeremylt                     matrix expressing the derivatives of nodal basis functions
475d979a051Sjeremylt                     at quadrature points
476d979a051Sjeremylt             *qref: Numpy array of length (nqpts * dim) holding the locations of
477d979a051Sjeremylt                     quadrature points on the reference element [-1, 1]
478d979a051Sjeremylt             *qweight: Numpy array of length nnodes holding the quadrature
479d979a051Sjeremylt                        weights on the reference element
48039b2de37Sjeremylt
48139b2de37Sjeremylt           Returns:
48239b2de37Sjeremylt             basis: Ceed Basis"""
48339b2de37Sjeremylt
4847a7b0fa3SJed Brown        return BasisH1(self, topo, ncomp, nnodes, nqpts,
4857a7b0fa3SJed Brown                       interp, grad, qref, qweight)
48639b2de37Sjeremylt
48739b2de37Sjeremylt    # CeedQFunction
48839b2de37Sjeremylt    def QFunction(self, vlength, f, source):
489d979a051Sjeremylt        """Ceed QFunction: point-wise operation at quadrature points for
490d979a051Sjeremylt             evaluating volumetric terms.
49139b2de37Sjeremylt
49239b2de37Sjeremylt           Args:
49339b2de37Sjeremylt             vlength: vector length. Caller must ensure that number of quadrature
49439b2de37Sjeremylt                        points is a multiple of vlength
49539b2de37Sjeremylt             f: ctypes function pointer to evaluate action at quadrature points
49639b2de37Sjeremylt             source: absolute path to source of QFunction,
49739b2de37Sjeremylt               "\\abs_path\\file.h:function_name
49839b2de37Sjeremylt
49939b2de37Sjeremylt           Returns:
50039b2de37Sjeremylt             qfunction: Ceed QFunction"""
50139b2de37Sjeremylt
50239b2de37Sjeremylt        return QFunction(self, vlength, f, source)
50339b2de37Sjeremylt
50439b2de37Sjeremylt    def QFunctionByName(self, name):
50539b2de37Sjeremylt        """Ceed QFunction By Name: point-wise operation at quadrature points
50639b2de37Sjeremylt             from a given gallery, for evaluating volumetric terms.
50739b2de37Sjeremylt
50839b2de37Sjeremylt           Args:
50939b2de37Sjeremylt             name: name of QFunction to use from gallery
51039b2de37Sjeremylt
51139b2de37Sjeremylt           Returns:
51239b2de37Sjeremylt             qfunction: Ceed QFunction By Name"""
51339b2de37Sjeremylt
51439b2de37Sjeremylt        return QFunctionByName(self, name)
51539b2de37Sjeremylt
51639b2de37Sjeremylt    def IdentityQFunction(self, size, inmode, outmode):
51739b2de37Sjeremylt        """Ceed Idenity QFunction: identity qfunction operation.
51839b2de37Sjeremylt
51939b2de37Sjeremylt           Args:
52039b2de37Sjeremylt             size: size of the qfunction fields
52139b2de37Sjeremylt             **inmode: CeedEvalMode for input to Ceed QFunction
52239b2de37Sjeremylt             **outmode: CeedEvalMode for output to Ceed QFunction
52339b2de37Sjeremylt
52439b2de37Sjeremylt           Returns:
52539b2de37Sjeremylt             qfunction: Ceed Identity QFunction"""
52639b2de37Sjeremylt
52739b2de37Sjeremylt        return IdentityQFunction(self, size, inmode, outmode)
52839b2de37Sjeremylt
529777ff853SJeremy L Thompson    def QFunctionContext(self):
530777ff853SJeremy L Thompson        """Ceed QFunction Context: stores Ceed QFunction user context data.
531777ff853SJeremy L Thompson
532777ff853SJeremy L Thompson           Returns:
533777ff853SJeremy L Thompson             userContext: Ceed QFunction Context"""
534777ff853SJeremy L Thompson
535777ff853SJeremy L Thompson        return QFunctionContext(self)
536777ff853SJeremy L Thompson
53739b2de37Sjeremylt    # CeedOperator
53839b2de37Sjeremylt    def Operator(self, qf, dqf=None, qdfT=None):
53939b2de37Sjeremylt        """Ceed Operator: composed FE-type operations on vectors.
54039b2de37Sjeremylt
54139b2de37Sjeremylt           Args:
542d979a051Sjeremylt             qf: QFunction defining the action of the operator at quadrature
543d979a051Sjeremylt                   points
54439b2de37Sjeremylt             **dqf: QFunction defining the action of the Jacobian, default None
545d979a051Sjeremylt             **dqfT: QFunction defining the action of the transpose of the
546d979a051Sjeremylt                       Jacobian, default None
54739b2de37Sjeremylt
54839b2de37Sjeremylt           Returns:
54939b2de37Sjeremylt             operator: Ceed Operator"""
55039b2de37Sjeremylt
55139b2de37Sjeremylt        return Operator(self, qf, dqf, qdfT)
55239b2de37Sjeremylt
55339b2de37Sjeremylt    def CompositeOperator(self):
55439b2de37Sjeremylt        """Composite Ceed Operator: composition of multiple CeedOperators.
55539b2de37Sjeremylt
55639b2de37Sjeremylt           Returns:
55739b2de37Sjeremylt             operator: Ceed Composite Operator"""
55839b2de37Sjeremylt
55939b2de37Sjeremylt        return CompositeOperator(self)
56039b2de37Sjeremylt
56139b2de37Sjeremylt    # Destructor
56239b2de37Sjeremylt    def __del__(self):
56339b2de37Sjeremylt        # libCEED call
56439b2de37Sjeremylt        lib.CeedDestroy(self._pointer)
56539b2de37Sjeremylt
56639b2de37Sjeremylt# ------------------------------------------------------------------------------
567