xref: /libCEED/python/ceed.py (revision 709403c1d8a3197f71934fa05519eedc9835a633)
13d8e8822SJeremy L Thompson# Copyright (c) 2017-2022, 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.
339b2de37Sjeremylt#
43d8e8822SJeremy L Thompson# SPDX-License-Identifier: BSD-2-Clause
539b2de37Sjeremylt#
63d8e8822SJeremy L Thompson# This file is part of CEED:  http://github.com/ceed
739b2de37Sjeremylt
839b2de37Sjeremyltfrom _ceed_cffi import ffi, lib
939b2de37Sjeremyltimport sys
10477729cfSJeremy L Thompsonimport os
1139b2de37Sjeremyltimport io
12e15f9bd0SJeremy L Thompsonimport numpy as np
130a0da059Sjeremyltimport tempfile
1439b2de37Sjeremyltfrom abc import ABC
1539b2de37Sjeremyltfrom .ceed_vector import Vector
1697c1c57aSSebastian Grimbergfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1, BasisHdiv, BasisHcurl
17*709403c1SSebastian Grimbergfrom .ceed_elemrestriction import ElemRestriction, OrientedElemRestriction, CurlOrientedElemRestriction, StridedElemRestriction
18*709403c1SSebastian Grimbergfrom .ceed_elemrestriction import BlockedElemRestriction, BlockedOrientedElemRestriction, BlockedCurlOrientedElemRestriction, BlockedStridedElemRestriction
1939b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction
20777ff853SJeremy L Thompsonfrom .ceed_qfunctioncontext import QFunctionContext
2139b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator
2239b2de37Sjeremyltfrom .ceed_constants import *
2339b2de37Sjeremylt
2439b2de37Sjeremylt# ------------------------------------------------------------------------------
257a7b0fa3SJed Brown
267a7b0fa3SJed Brown
2739b2de37Sjeremyltclass Ceed():
2839b2de37Sjeremylt    """Ceed: core components."""
2939b2de37Sjeremylt
3039b2de37Sjeremylt    # Constructor
3174b533adSJed Brown    def __init__(self, resource="/cpu/self", on_error="store"):
3239b2de37Sjeremylt        # libCEED object
3339b2de37Sjeremylt        self._pointer = ffi.new("Ceed *")
3439b2de37Sjeremylt
3539b2de37Sjeremylt        # libCEED call
3639b2de37Sjeremylt        resourceAscii = ffi.new("char[]", resource.encode("ascii"))
37477729cfSJeremy L Thompson        os.environ["CEED_ERROR_HANDLER"] = "return"
38477729cfSJeremy L Thompson        err_code = lib.CeedInit(resourceAscii, self._pointer)
39477729cfSJeremy L Thompson        if err_code:
40477729cfSJeremy L Thompson            raise Exception("Error initializing backend resource: " + resource)
4174b533adSJed Brown        error_handlers = dict(
4274b533adSJed Brown            store="CeedErrorStore",
4374b533adSJed Brown            abort="CeedErrorAbort",
4474b533adSJed Brown        )
45477729cfSJeremy L Thompson        lib.CeedSetErrorHandler(
46477729cfSJeremy L Thompson            self._pointer[0], ffi.addressof(
4774b533adSJed Brown                lib, error_handlers[on_error]))
4839b2de37Sjeremylt
4939b2de37Sjeremylt    # Representation
5039b2de37Sjeremylt    def __repr__(self):
5139b2de37Sjeremylt        return "<Ceed instance at " + hex(id(self)) + ">"
5239b2de37Sjeremylt
530a0da059Sjeremylt    # String conversion for print() to stdout
540a0da059Sjeremylt    def __str__(self):
550a0da059Sjeremylt        """View a Ceed via print()."""
560a0da059Sjeremylt
570a0da059Sjeremylt        # libCEED call
580a0da059Sjeremylt        with tempfile.NamedTemporaryFile() as key_file:
590a0da059Sjeremylt            with open(key_file.name, 'r+') as stream_file:
600a0da059Sjeremylt                stream = ffi.cast("FILE *", stream_file)
610a0da059Sjeremylt
62477729cfSJeremy L Thompson                err_code = lib.CeedView(self._pointer[0], stream)
63477729cfSJeremy L Thompson                self._check_error(err_code)
640a0da059Sjeremylt
650a0da059Sjeremylt                stream_file.seek(0)
660a0da059Sjeremylt                out_string = stream_file.read()
670a0da059Sjeremylt
680a0da059Sjeremylt        return out_string
690a0da059Sjeremylt
70477729cfSJeremy L Thompson    # Error handler
71477729cfSJeremy L Thompson    def _check_error(self, err_code):
72477729cfSJeremy L Thompson        """Check return code and retrieve error message for non-zero code"""
73e15f9bd0SJeremy L Thompson        if (err_code != lib.CEED_ERROR_SUCCESS):
74477729cfSJeremy L Thompson            message = ffi.new("char **")
75477729cfSJeremy L Thompson            lib.CeedGetErrorMessage(self._pointer[0], message)
76477729cfSJeremy L Thompson            raise Exception(ffi.string(message[0]).decode("UTF-8"))
77477729cfSJeremy L Thompson
7839b2de37Sjeremylt    # Get Resource
7939b2de37Sjeremylt    def get_resource(self):
8039b2de37Sjeremylt        """Get the full resource name for a Ceed context.
8139b2de37Sjeremylt
8239b2de37Sjeremylt           Returns:
8339b2de37Sjeremylt             resource: resource name"""
8439b2de37Sjeremylt
8539b2de37Sjeremylt        # libCEED call
8639b2de37Sjeremylt        resource = ffi.new("char **")
87477729cfSJeremy L Thompson        err_code = lib.CeedGetResource(self._pointer[0], resource)
88477729cfSJeremy L Thompson        self._check_error(err_code)
8939b2de37Sjeremylt
9039b2de37Sjeremylt        return ffi.string(resource[0]).decode("UTF-8")
9139b2de37Sjeremylt
9239b2de37Sjeremylt    # Preferred MemType
9339b2de37Sjeremylt    def get_preferred_memtype(self):
9439b2de37Sjeremylt        """Return Ceed preferred memory type.
9539b2de37Sjeremylt
9639b2de37Sjeremylt           Returns:
9739b2de37Sjeremylt             memtype: Ceed preferred memory type"""
9839b2de37Sjeremylt
9939b2de37Sjeremylt        # libCEED call
10039b2de37Sjeremylt        memtype = ffi.new("CeedMemType *", MEM_HOST)
101477729cfSJeremy L Thompson        err_code = lib.CeedGetPreferredMemType(self._pointer[0], memtype)
102477729cfSJeremy L Thompson        self._check_error(err_code)
10339b2de37Sjeremylt
10439b2de37Sjeremylt        return memtype[0]
10539b2de37Sjeremylt
10680a9ef05SNatalie Beams    # Convenience function to get CeedScalar type
10780a9ef05SNatalie Beams    def scalar_type(self):
10880a9ef05SNatalie Beams        """Return dtype string for CeedScalar.
10980a9ef05SNatalie Beams
11080a9ef05SNatalie Beams           Returns:
11180a9ef05SNatalie Beams             np_dtype: String naming numpy data type corresponding to CeedScalar"""
11280a9ef05SNatalie Beams
11380a9ef05SNatalie Beams        return scalar_types[lib.CEED_SCALAR_TYPE]
11480a9ef05SNatalie Beams
115e15f9bd0SJeremy L Thompson    # --- Basis utility functions ---
116e15f9bd0SJeremy L Thompson
117e15f9bd0SJeremy L Thompson    # Gauss quadrature
118e15f9bd0SJeremy L Thompson    def gauss_quadrature(self, q):
119e15f9bd0SJeremy L Thompson        """Construct a Gauss-Legendre quadrature.
120e15f9bd0SJeremy L Thompson
121e15f9bd0SJeremy L Thompson           Args:
122e15f9bd0SJeremy L Thompson             Q: number of quadrature points (integrates polynomials of
123e15f9bd0SJeremy L Thompson                  degree 2*Q-1 exactly)
124e15f9bd0SJeremy L Thompson
125e15f9bd0SJeremy L Thompson           Returns:
126e15f9bd0SJeremy L Thompson             (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
127e15f9bd0SJeremy L Thompson                                    and array of length Q to hold the weights"""
128e15f9bd0SJeremy L Thompson
129e15f9bd0SJeremy L Thompson        # Setup arguments
13080a9ef05SNatalie Beams        qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
13180a9ef05SNatalie Beams        qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
132e15f9bd0SJeremy L Thompson
133e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.new("CeedScalar *")
134e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.cast(
135e15f9bd0SJeremy L Thompson            "CeedScalar *",
136e15f9bd0SJeremy L Thompson            qref1d.__array_interface__['data'][0])
137e15f9bd0SJeremy L Thompson
138e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.new("CeedScalar *")
139e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.cast(
140e15f9bd0SJeremy L Thompson            "CeedScalar *",
141e15f9bd0SJeremy L Thompson            qweight1d.__array_interface__['data'][0])
142e15f9bd0SJeremy L Thompson
143e15f9bd0SJeremy L Thompson        # libCEED call
144e15f9bd0SJeremy L Thompson        err_code = lib.CeedGaussQuadrature(q, qref1d_pointer, qweight1d_pointer)
145e15f9bd0SJeremy L Thompson        self._check_error(err_code)
146e15f9bd0SJeremy L Thompson
147e15f9bd0SJeremy L Thompson        return qref1d, qweight1d
148e15f9bd0SJeremy L Thompson
149e15f9bd0SJeremy L Thompson    # Lobatto quadrature
150e15f9bd0SJeremy L Thompson    def lobatto_quadrature(self, q):
151e15f9bd0SJeremy L Thompson        """Construct a Gauss-Legendre-Lobatto quadrature.
152e15f9bd0SJeremy L Thompson
153e15f9bd0SJeremy L Thompson           Args:
154e15f9bd0SJeremy L Thompson             q: number of quadrature points (integrates polynomials of
155e15f9bd0SJeremy L Thompson                  degree 2*Q-3 exactly)
156e15f9bd0SJeremy L Thompson
157e15f9bd0SJeremy L Thompson           Returns:
158e15f9bd0SJeremy L Thompson             (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1]
159e15f9bd0SJeremy L Thompson                                    and array of length Q to hold the weights"""
160e15f9bd0SJeremy L Thompson
161e15f9bd0SJeremy L Thompson        # Setup arguments
16280a9ef05SNatalie Beams        qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
163e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.new("CeedScalar *")
164e15f9bd0SJeremy L Thompson        qref1d_pointer = ffi.cast(
165e15f9bd0SJeremy L Thompson            "CeedScalar *",
166e15f9bd0SJeremy L Thompson            qref1d.__array_interface__['data'][0])
167e15f9bd0SJeremy L Thompson
16880a9ef05SNatalie Beams        qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE])
169e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.new("CeedScalar *")
170e15f9bd0SJeremy L Thompson        qweight1d_pointer = ffi.cast(
171e15f9bd0SJeremy L Thompson            "CeedScalar *",
172e15f9bd0SJeremy L Thompson            qweight1d.__array_interface__['data'][0])
173e15f9bd0SJeremy L Thompson
174e15f9bd0SJeremy L Thompson        # libCEED call
175e15f9bd0SJeremy L Thompson        err_code = lib.CeedLobattoQuadrature(
176e15f9bd0SJeremy L Thompson            q, qref1d_pointer, qweight1d_pointer)
177e15f9bd0SJeremy L Thompson        self._check_error(err_code)
178e15f9bd0SJeremy L Thompson
179e15f9bd0SJeremy L Thompson        return qref1d, qweight1d
180e15f9bd0SJeremy L Thompson
181e15f9bd0SJeremy L Thompson    # --- libCEED objects ---
182e15f9bd0SJeremy L Thompson
18339b2de37Sjeremylt    # CeedVector
18439b2de37Sjeremylt    def Vector(self, size):
18539b2de37Sjeremylt        """Ceed Vector: storing and manipulating vectors.
18639b2de37Sjeremylt
18739b2de37Sjeremylt           Args:
18839b2de37Sjeremylt             size: length of vector
18939b2de37Sjeremylt
19039b2de37Sjeremylt           Returns:
19139b2de37Sjeremylt             vector: Ceed Vector"""
19239b2de37Sjeremylt
19339b2de37Sjeremylt        return Vector(self, size)
19439b2de37Sjeremylt
19539b2de37Sjeremylt    # CeedElemRestriction
196d979a051Sjeremylt    def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets,
197d979a051Sjeremylt                        memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
19839b2de37Sjeremylt        """Ceed ElemRestriction: restriction from local vectors to elements.
19939b2de37Sjeremylt
20039b2de37Sjeremylt           Args:
201d979a051Sjeremylt             nelem: number of elements described by the restriction
20239b2de37Sjeremylt             elemsize: size (number of nodes) per element
203d979a051Sjeremylt             ncomp: number of field components per interpolation node
204d979a051Sjeremylt                      (1 for scalar fields)
205d979a051Sjeremylt             compstride: Stride between components for the same L-vector "node".
206d979a051Sjeremylt                           Data for node i, component k can be found in the
207d979a051Sjeremylt                           L-vector at index [offsets[i] + k*compstride].
208d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
209d979a051Sjeremylt                       the elements and fields given by this restriction.
210d979a051Sjeremylt             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
211d979a051Sjeremylt                         ordered list of the offsets (into the input Ceed Vector)
21239b2de37Sjeremylt                         for the unknowns corresponding to element i, where
213d979a051Sjeremylt                         0 <= i < nelem. All offsets must be in the range
214d979a051Sjeremylt                         [0, lsize - 1].
215d979a051Sjeremylt             **memtype: memory type of the offsets array, default CEED_MEM_HOST
216d979a051Sjeremylt             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
21739b2de37Sjeremylt
21839b2de37Sjeremylt           Returns:
21939b2de37Sjeremylt             elemrestriction: Ceed ElemRestiction"""
22039b2de37Sjeremylt
221d979a051Sjeremylt        return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize,
222d979a051Sjeremylt                               offsets, memtype=memtype, cmode=cmode)
22339b2de37Sjeremylt
224*709403c1SSebastian Grimberg    def OrientedElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize,
225*709403c1SSebastian Grimberg                                offsets, orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
226*709403c1SSebastian Grimberg        """Ceed Oriented ElemRestriction: oriented restriction from local vectors
227*709403c1SSebastian Grimberg             to elements.
228*709403c1SSebastian Grimberg
229*709403c1SSebastian Grimberg           Args:
230*709403c1SSebastian Grimberg             nelem: number of elements described by the restriction
231*709403c1SSebastian Grimberg             elemsize: size (number of nodes) per element
232*709403c1SSebastian Grimberg             ncomp: number of field components per interpolation node
233*709403c1SSebastian Grimberg                      (1 for scalar fields)
234*709403c1SSebastian Grimberg             compstride: Stride between components for the same L-vector "node".
235*709403c1SSebastian Grimberg                           Data for node i, component k can be found in the
236*709403c1SSebastian Grimberg                           L-vector at index [offsets[i] + k*compstride].
237*709403c1SSebastian Grimberg             lsize: The size of the L-vector. This vector may be larger than
238*709403c1SSebastian Grimberg                       the elements and fields given by this restriction.
239*709403c1SSebastian Grimberg             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
240*709403c1SSebastian Grimberg                         ordered list of the offsets (into the input Ceed Vector)
241*709403c1SSebastian Grimberg                         for the unknowns corresponding to element i, where
242*709403c1SSebastian Grimberg                         0 <= i < nelem. All offsets must be in the range
243*709403c1SSebastian Grimberg                         [0, lsize - 1].
244*709403c1SSebastian Grimberg             *orients: Numpy array of shape [nelem, elemsize]. Row i holds the
245*709403c1SSebastian Grimberg                         ordered list of the orientations for the unknowns
246*709403c1SSebastian Grimberg                         corresponding to element i, with bool false used for
247*709403c1SSebastian Grimberg                         positively oriented and true to flip the orientation.
248*709403c1SSebastian Grimberg             **memtype: memory type of the offsets array, default CEED_MEM_HOST
249*709403c1SSebastian Grimberg             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
250*709403c1SSebastian Grimberg
251*709403c1SSebastian Grimberg           Returns:
252*709403c1SSebastian Grimberg             elemrestriction: Ceed Oriented ElemRestiction"""
253*709403c1SSebastian Grimberg
254*709403c1SSebastian Grimberg        return OrientedElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize,
255*709403c1SSebastian Grimberg                                       offsets, orients, memtype=memtype, cmode=cmode)
256*709403c1SSebastian Grimberg
257*709403c1SSebastian Grimberg    def CurlOrientedElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize,
258*709403c1SSebastian Grimberg                                    offsets, curl_orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
259*709403c1SSebastian Grimberg        """Ceed Curl Oriented ElemRestriction: curl-oriented restriction from local
260*709403c1SSebastian Grimberg             vectors to elements.
261*709403c1SSebastian Grimberg
262*709403c1SSebastian Grimberg           Args:
263*709403c1SSebastian Grimberg             nelem: number of elements described by the restriction
264*709403c1SSebastian Grimberg             elemsize: size (number of nodes) per element
265*709403c1SSebastian Grimberg             ncomp: number of field components per interpolation node
266*709403c1SSebastian Grimberg                      (1 for scalar fields)
267*709403c1SSebastian Grimberg             compstride: Stride between components for the same L-vector "node".
268*709403c1SSebastian Grimberg                           Data for node i, component k can be found in the
269*709403c1SSebastian Grimberg                           L-vector at index [offsets[i] + k*compstride].
270*709403c1SSebastian Grimberg             lsize: The size of the L-vector. This vector may be larger than
271*709403c1SSebastian Grimberg                       the elements and fields given by this restriction.
272*709403c1SSebastian Grimberg             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
273*709403c1SSebastian Grimberg                         ordered list of the offsets (into the input Ceed Vector)
274*709403c1SSebastian Grimberg                         for the unknowns corresponding to element i, where
275*709403c1SSebastian Grimberg                         0 <= i < nelem. All offsets must be in the range
276*709403c1SSebastian Grimberg                         [0, lsize - 1].
277*709403c1SSebastian Grimberg             *curl_orients: Numpy array of shape [nelem, 3 * elemsize]. Row i holds
278*709403c1SSebastian Grimberg                              a row-major tridiagonal matrix (curl_orients[i, 0] =
279*709403c1SSebastian Grimberg                              curl_orients[i, 3 * elemsize - 1] = 0, where 0 <= i < nelem)
280*709403c1SSebastian Grimberg                              which is applied to the element unknowns upon restriction.
281*709403c1SSebastian Grimberg             **memtype: memory type of the offsets array, default CEED_MEM_HOST
282*709403c1SSebastian Grimberg             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
283*709403c1SSebastian Grimberg
284*709403c1SSebastian Grimberg           Returns:
285*709403c1SSebastian Grimberg             elemrestriction: Ceed Curl Oriented ElemRestiction"""
286*709403c1SSebastian Grimberg
287*709403c1SSebastian Grimberg        return CurlOrientedElemRestriction(
288*709403c1SSebastian Grimberg            self, nelem, elemsize, ncomp, compstride, lsize, offsets, curl_orients, memtype=memtype, cmode=cmode)
289*709403c1SSebastian Grimberg
290d979a051Sjeremylt    def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides):
291d979a051Sjeremylt        """Ceed Identity ElemRestriction: strided restriction from local vectors
292d979a051Sjeremylt             to elements.
29339b2de37Sjeremylt
29439b2de37Sjeremylt           Args:
295d979a051Sjeremylt             nelem: number of elements described by the restriction
29639b2de37Sjeremylt             elemsize: size (number of nodes) per element
297a8d32208Sjeremylt             ncomp: number of field components per interpolation node
298a8d32208Sjeremylt                      (1 for scalar fields)
299d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
300d979a051Sjeremylt                      the elements and fields given by this restriction.
30169a53589Sjeremylt             *strides: Array for strides between [nodes, components, elements].
30269a53589Sjeremylt                         The data for node i, component j, element k in the
30369a53589Sjeremylt                         L-vector is given by
30469a53589Sjeremylt                           i*strides[0] + j*strides[1] + k*strides[2]
30539b2de37Sjeremylt
30639b2de37Sjeremylt           Returns:
30769a53589Sjeremylt             elemrestriction: Ceed Strided ElemRestiction"""
30839b2de37Sjeremylt
3097a7b0fa3SJed Brown        return StridedElemRestriction(
310d979a051Sjeremylt            self, nelem, elemsize, ncomp, lsize, strides)
31139b2de37Sjeremylt
312d979a051Sjeremylt    def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride,
313d979a051Sjeremylt                               lsize, offsets, memtype=lib.CEED_MEM_HOST,
314d979a051Sjeremylt                               cmode=lib.CEED_COPY_VALUES):
315d979a051Sjeremylt        """Ceed Blocked ElemRestriction: blocked restriction from local vectors
316d979a051Sjeremylt             to elements.
31739b2de37Sjeremylt
31839b2de37Sjeremylt           Args:
319d979a051Sjeremylt             nelem: number of elements described by the restriction
32039b2de37Sjeremylt             elemsize: size (number of nodes) per element
32139b2de37Sjeremylt             blksize: number of elements in a block
322d979a051Sjeremylt             ncomp: number of field components per interpolation node
323d979a051Sjeremylt                       (1 for scalar fields)
324d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
325d979a051Sjeremylt                      the elements and fields given by this restriction.
326d979a051Sjeremylt             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
327d979a051Sjeremylt                         ordered list of the offsets (into the input Ceed Vector)
32839b2de37Sjeremylt                         for the unknowns corresponding to element i, where
329d979a051Sjeremylt                         0 <= i < nelem. All offsets must be in the range
330d979a051Sjeremylt                         [0, lsize - 1]. The backend will permute and pad this
33139b2de37Sjeremylt                         array to the desired ordering for the blocksize, which is
33239b2de37Sjeremylt                         typically given by the backend. The default reordering is
33339b2de37Sjeremylt                         to interlace elements.
334d979a051Sjeremylt             **memtype: memory type of the offsets array, default CEED_MEM_HOST
335d979a051Sjeremylt             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
33639b2de37Sjeremylt
33739b2de37Sjeremylt           Returns:
33839b2de37Sjeremylt             elemrestriction: Ceed Blocked ElemRestiction"""
33939b2de37Sjeremylt
340d979a051Sjeremylt        return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp,
341d979a051Sjeremylt                                      compstride, lsize, offsets,
342d979a051Sjeremylt                                      memtype=memtype, cmode=cmode)
34339b2de37Sjeremylt
344*709403c1SSebastian Grimberg    def BlockedOrientedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride,
345*709403c1SSebastian Grimberg                                       lsize, offsets, orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
346*709403c1SSebastian Grimberg        """Ceed Blocked Oriented ElemRestriction: blocked oriented restriction
347*709403c1SSebastian Grimberg             from local vectors to elements.
348*709403c1SSebastian Grimberg
349*709403c1SSebastian Grimberg           Args:
350*709403c1SSebastian Grimberg             nelem: number of elements described by the restriction
351*709403c1SSebastian Grimberg             elemsize: size (number of nodes) per element
352*709403c1SSebastian Grimberg             blksize: number of elements in a block
353*709403c1SSebastian Grimberg             ncomp: number of field components per interpolation node
354*709403c1SSebastian Grimberg                      (1 for scalar fields)
355*709403c1SSebastian Grimberg             compstride: Stride between components for the same L-vector "node".
356*709403c1SSebastian Grimberg                           Data for node i, component k can be found in the
357*709403c1SSebastian Grimberg                           L-vector at index [offsets[i] + k*compstride].
358*709403c1SSebastian Grimberg             lsize: The size of the L-vector. This vector may be larger than
359*709403c1SSebastian Grimberg                       the elements and fields given by this restriction.
360*709403c1SSebastian Grimberg             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
361*709403c1SSebastian Grimberg                         ordered list of the offsets (into the input Ceed Vector)
362*709403c1SSebastian Grimberg                         for the unknowns corresponding to element i, where
363*709403c1SSebastian Grimberg                         0 <= i < nelem. All offsets must be in the range
364*709403c1SSebastian Grimberg                         [0, lsize - 1]. The backend will permute and pad this
365*709403c1SSebastian Grimberg                         array to the desired ordering for the blocksize, which is
366*709403c1SSebastian Grimberg                         typically given by the backend. The default reordering is
367*709403c1SSebastian Grimberg                         to interlace elements.
368*709403c1SSebastian Grimberg             *orients: Numpy array of shape [nelem, elemsize]. Row i holds the
369*709403c1SSebastian Grimberg                         ordered list of the orientations for the unknowns
370*709403c1SSebastian Grimberg                         corresponding to element i, with ool false is used for
371*709403c1SSebastian Grimberg                         positively oriented and true to flip the orientation.
372*709403c1SSebastian Grimberg             **memtype: memory type of the offsets array, default CEED_MEM_HOST
373*709403c1SSebastian Grimberg             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
374*709403c1SSebastian Grimberg
375*709403c1SSebastian Grimberg           Returns:
376*709403c1SSebastian Grimberg             elemrestriction: Ceed Blocked Oriented ElemRestiction"""
377*709403c1SSebastian Grimberg
378*709403c1SSebastian Grimberg        return BlockedOrientedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, lsize,
379*709403c1SSebastian Grimberg                                              offsets, orients, memtype=memtype, cmode=cmode)
380*709403c1SSebastian Grimberg
381*709403c1SSebastian Grimberg    def BlockedCurlOrientedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride,
382*709403c1SSebastian Grimberg                                           lsize, offsets, curl_orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
383*709403c1SSebastian Grimberg        """Ceed Blocked Curl Oriented ElemRestriction: blocked curl-oriented
384*709403c1SSebastian Grimberg             restriction from local vectors to elements.
385*709403c1SSebastian Grimberg
386*709403c1SSebastian Grimberg           Args:
387*709403c1SSebastian Grimberg             nelem: number of elements described by the restriction
388*709403c1SSebastian Grimberg             elemsize: size (number of nodes) per element
389*709403c1SSebastian Grimberg             blksize: number of elements in a block
390*709403c1SSebastian Grimberg             ncomp: number of field components per interpolation node
391*709403c1SSebastian Grimberg                      (1 for scalar fields)
392*709403c1SSebastian Grimberg             compstride: Stride between components for the same L-vector "node".
393*709403c1SSebastian Grimberg                           Data for node i, component k can be found in the
394*709403c1SSebastian Grimberg                           L-vector at index [offsets[i] + k*compstride].
395*709403c1SSebastian Grimberg             lsize: The size of the L-vector. This vector may be larger than
396*709403c1SSebastian Grimberg                       the elements and fields given by this restriction.
397*709403c1SSebastian Grimberg             *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the
398*709403c1SSebastian Grimberg                         ordered list of the offsets (into the input Ceed Vector)
399*709403c1SSebastian Grimberg                         for the unknowns corresponding to element i, where
400*709403c1SSebastian Grimberg                         0 <= i < nelem. All offsets must be in the range
401*709403c1SSebastian Grimberg                         [0, lsize - 1]. The backend will permute and pad this
402*709403c1SSebastian Grimberg                         array to the desired ordering for the blocksize, which is
403*709403c1SSebastian Grimberg                         typically given by the backend. The default reordering is
404*709403c1SSebastian Grimberg                         to interlace elements.
405*709403c1SSebastian Grimberg             *curl_orients: Numpy array of shape [nelem, 3 * elemsize]. Row i holds
406*709403c1SSebastian Grimberg                              a row-major tridiagonal matrix (curl_orients[i, 0] =
407*709403c1SSebastian Grimberg                              curl_orients[i, 3 * elemsize - 1] = 0, where 0 <= i < nelem)
408*709403c1SSebastian Grimberg                              which is applied to the element unknowns upon restriction.
409*709403c1SSebastian Grimberg             **memtype: memory type of the offsets array, default CEED_MEM_HOST
410*709403c1SSebastian Grimberg             **cmode: copy mode for the offsets array, default CEED_COPY_VALUES
411*709403c1SSebastian Grimberg
412*709403c1SSebastian Grimberg           Returns:
413*709403c1SSebastian Grimberg             elemrestriction: Ceed Blocked Curl Oriented ElemRestiction"""
414*709403c1SSebastian Grimberg
415*709403c1SSebastian Grimberg        return BlockedCurlOrientedElemRestriction(
416*709403c1SSebastian Grimberg            self, nelem, elemsize, blksize, ncomp, compstride, lsize, offsets, curl_orients, memtype=memtype, cmode=cmode)
417*709403c1SSebastian Grimberg
418d979a051Sjeremylt    def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp,
419d979a051Sjeremylt                                      lsize, strides):
420d979a051Sjeremylt        """Ceed Blocked Strided ElemRestriction: blocked and strided restriction
421d979a051Sjeremylt             from local vectors to elements.
42269a53589Sjeremylt
42369a53589Sjeremylt           Args:
424d979a051Sjeremylt             nelem: number of elements described in the restriction
42569a53589Sjeremylt             elemsize: size (number of nodes) per element
42669a53589Sjeremylt             blksize: number of elements in a block
42769a53589Sjeremylt             ncomp: number of field components per interpolation node
42869a53589Sjeremylt                      (1 for scalar fields)
429d979a051Sjeremylt             lsize: The size of the L-vector. This vector may be larger than
430d979a051Sjeremylt                      the elements and fields given by this restriction.
43169a53589Sjeremylt             *strides: Array for strides between [nodes, components, elements].
43269a53589Sjeremylt                         The data for node i, component j, element k in the
43369a53589Sjeremylt                         L-vector is given by
43469a53589Sjeremylt                           i*strides[0] + j*strides[1] + k*strides[2]
43569a53589Sjeremylt
43669a53589Sjeremylt           Returns:
437*709403c1SSebastian Grimberg             elemrestriction: Ceed Blocked Strided ElemRestiction"""
43869a53589Sjeremylt
43969a53589Sjeremylt        return BlockedStridedElemRestriction(self, nelem, elemsize, blksize,
440d979a051Sjeremylt                                             ncomp, lsize, strides)
44169a53589Sjeremylt
44239b2de37Sjeremylt    # CeedBasis
44339b2de37Sjeremylt    def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
44439b2de37Sjeremylt                      qref1d, qweight1d):
44539b2de37Sjeremylt        """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
44639b2de37Sjeremylt             H^1 discretizations.
44739b2de37Sjeremylt
44839b2de37Sjeremylt           Args:
44939b2de37Sjeremylt             dim: topological dimension
45039b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
45139b2de37Sjeremylt             P1d: number of nodes in one dimension
45239b2de37Sjeremylt             Q1d: number of quadrature points in one dimension
453d979a051Sjeremylt             *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix
454d979a051Sjeremylt                          expressing the values of nodal basis functions at
455d979a051Sjeremylt                          quadrature points
456d979a051Sjeremylt             *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix
457d979a051Sjeremylt                        expressing the derivatives of nodal basis functions at
458d979a051Sjeremylt                        quadrature points
459d979a051Sjeremylt             *qref1d: Numpy array of length Q1d holding the locations of
460d979a051Sjeremylt                        quadrature points on the 1D reference element [-1, 1]
461d979a051Sjeremylt             *qweight1d: Numpy array of length Q1d holding the quadrature
462d979a051Sjeremylt                           weights on the reference element
46339b2de37Sjeremylt
46439b2de37Sjeremylt           Returns:
46539b2de37Sjeremylt             basis: Ceed Basis"""
46639b2de37Sjeremylt
46739b2de37Sjeremylt        return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
46839b2de37Sjeremylt                             qref1d, qweight1d)
46939b2de37Sjeremylt
47039b2de37Sjeremylt    def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode):
471d979a051Sjeremylt        """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange
472d979a051Sjeremylt             basis objects for H^1 discretizations.
47339b2de37Sjeremylt
47439b2de37Sjeremylt           Args:
47539b2de37Sjeremylt             dim: topological dimension
47639b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
47739b2de37Sjeremylt             P: number of Gauss-Lobatto nodes in one dimension.  The
47839b2de37Sjeremylt                  polynomial degree of the resulting Q_k element is k=P-1.
47939b2de37Sjeremylt             Q: number of quadrature points in one dimension
48039b2de37Sjeremylt             qmode: distribution of the Q quadrature points (affects order of
48139b2de37Sjeremylt                      accuracy for the quadrature)
48239b2de37Sjeremylt
48339b2de37Sjeremylt           Returns:
48439b2de37Sjeremylt             basis: Ceed Basis"""
48539b2de37Sjeremylt
48639b2de37Sjeremylt        return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode)
48739b2de37Sjeremylt
48839b2de37Sjeremylt    def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight):
489d979a051Sjeremylt        """Ceed H1 Basis: finite element non tensor-product basis for H^1
490d979a051Sjeremylt             discretizations.
49139b2de37Sjeremylt
49239b2de37Sjeremylt           Args:
49339b2de37Sjeremylt             topo: topology of the element, e.g. hypercube, simplex, etc
49439b2de37Sjeremylt             ncomp: number of field components (1 for scalar fields)
49539b2de37Sjeremylt             nnodes: total number of nodes
49639b2de37Sjeremylt             nqpts: total number of quadrature points
497d979a051Sjeremylt             *interp: Numpy array holding the row-major (nqpts * nnodes) matrix
498d979a051Sjeremylt                       expressing the values of nodal basis functions at
49939b2de37Sjeremylt                       quadrature points
50097c1c57aSSebastian Grimberg             *grad: Numpy array holding the row-major (dim * nqpts * nnodes)
501d979a051Sjeremylt                     matrix expressing the derivatives of nodal basis functions
502d979a051Sjeremylt                     at quadrature points
503d979a051Sjeremylt             *qref: Numpy array of length (nqpts * dim) holding the locations of
504d979a051Sjeremylt                     quadrature points on the reference element [-1, 1]
505d979a051Sjeremylt             *qweight: Numpy array of length nnodes holding the quadrature
506d979a051Sjeremylt                        weights on the reference element
50739b2de37Sjeremylt
50839b2de37Sjeremylt           Returns:
50939b2de37Sjeremylt             basis: Ceed Basis"""
51039b2de37Sjeremylt
5117a7b0fa3SJed Brown        return BasisH1(self, topo, ncomp, nnodes, nqpts,
5127a7b0fa3SJed Brown                       interp, grad, qref, qweight)
51339b2de37Sjeremylt
51497c1c57aSSebastian Grimberg    def BasisHdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight):
51597c1c57aSSebastian Grimberg        """Ceed Hdiv Basis: finite element non tensor-product basis for H(div)
51697c1c57aSSebastian Grimberg             discretizations.
51797c1c57aSSebastian Grimberg
51897c1c57aSSebastian Grimberg           Args:
51997c1c57aSSebastian Grimberg             topo: topology of the element, e.g. hypercube, simplex, etc
52097c1c57aSSebastian Grimberg             ncomp: number of field components (1 for scalar fields)
52197c1c57aSSebastian Grimberg             nnodes: total number of nodes
52297c1c57aSSebastian Grimberg             nqpts: total number of quadrature points
52397c1c57aSSebastian Grimberg             *interp: Numpy array holding the row-major (dim * nqpts * nnodes)
52497c1c57aSSebastian Grimberg                       matrix expressing the values of basis functions at
52597c1c57aSSebastian Grimberg                       quadrature points
52697c1c57aSSebastian Grimberg             *div: Numpy array holding the row-major (nqpts * nnodes) matrix
52797c1c57aSSebastian Grimberg                    expressing the divergence of basis functions at
52897c1c57aSSebastian Grimberg                    quadrature points
52997c1c57aSSebastian Grimberg             *qref: Numpy array of length (nqpts * dim) holding the locations of
53097c1c57aSSebastian Grimberg                     quadrature points on the reference element [-1, 1]
53197c1c57aSSebastian Grimberg             *qweight: Numpy array of length nnodes holding the quadrature
53297c1c57aSSebastian Grimberg                        weights on the reference element
53397c1c57aSSebastian Grimberg
53497c1c57aSSebastian Grimberg           Returns:
53597c1c57aSSebastian Grimberg             basis: Ceed Basis"""
53697c1c57aSSebastian Grimberg
53797c1c57aSSebastian Grimberg        return BasisHdiv(self, topo, ncomp, nnodes, nqpts,
53897c1c57aSSebastian Grimberg                         interp, div, qref, qweight)
53997c1c57aSSebastian Grimberg
54097c1c57aSSebastian Grimberg    def BasisHcurl(self, topo, ncomp, nnodes, nqpts,
54197c1c57aSSebastian Grimberg                   interp, curl, qref, qweight):
54297c1c57aSSebastian Grimberg        """Ceed Hcurl Basis: finite element non tensor-product basis for H(curl)
54397c1c57aSSebastian Grimberg             discretizations.
54497c1c57aSSebastian Grimberg
54597c1c57aSSebastian Grimberg           Args:
54697c1c57aSSebastian Grimberg             topo: topology of the element, e.g. hypercube, simplex, etc
54797c1c57aSSebastian Grimberg             ncomp: number of field components (1 for scalar fields)
54897c1c57aSSebastian Grimberg             nnodes: total number of nodes
54997c1c57aSSebastian Grimberg             nqpts: total number of quadrature points
55097c1c57aSSebastian Grimberg             *interp: Numpy array holding the row-major (dim * nqpts * nnodes)
55197c1c57aSSebastian Grimberg                       matrix expressing the values of basis functions at
55297c1c57aSSebastian Grimberg                       quadrature points
55397c1c57aSSebastian Grimberg             *curl: Numpy array holding the row-major (curlcomp * nqpts * nnodes),
55497c1c57aSSebastian Grimberg                     curlcomp = 1 if dim < 3 else dim, matrix expressing the curl
55597c1c57aSSebastian Grimberg                     of basis functions at quadrature points
55697c1c57aSSebastian Grimberg             *qref: Numpy array of length (nqpts * dim) holding the locations of
55797c1c57aSSebastian Grimberg                     quadrature points on the reference element [-1, 1]
55897c1c57aSSebastian Grimberg             *qweight: Numpy array of length nnodes holding the quadrature
55997c1c57aSSebastian Grimberg                        weights on the reference element
56097c1c57aSSebastian Grimberg
56197c1c57aSSebastian Grimberg           Returns:
56297c1c57aSSebastian Grimberg             basis: Ceed Basis"""
56397c1c57aSSebastian Grimberg
56497c1c57aSSebastian Grimberg        return BasisHcurl(self, topo, ncomp, nnodes, nqpts,
56597c1c57aSSebastian Grimberg                          interp, curl, qref, qweight)
56697c1c57aSSebastian Grimberg
56739b2de37Sjeremylt    # CeedQFunction
56839b2de37Sjeremylt    def QFunction(self, vlength, f, source):
569d979a051Sjeremylt        """Ceed QFunction: point-wise operation at quadrature points for
570d979a051Sjeremylt             evaluating volumetric terms.
57139b2de37Sjeremylt
57239b2de37Sjeremylt           Args:
57339b2de37Sjeremylt             vlength: vector length. Caller must ensure that number of quadrature
57439b2de37Sjeremylt                        points is a multiple of vlength
57539b2de37Sjeremylt             f: ctypes function pointer to evaluate action at quadrature points
57639b2de37Sjeremylt             source: absolute path to source of QFunction,
57739b2de37Sjeremylt               "\\abs_path\\file.h:function_name
57839b2de37Sjeremylt
57939b2de37Sjeremylt           Returns:
58039b2de37Sjeremylt             qfunction: Ceed QFunction"""
58139b2de37Sjeremylt
58239b2de37Sjeremylt        return QFunction(self, vlength, f, source)
58339b2de37Sjeremylt
58439b2de37Sjeremylt    def QFunctionByName(self, name):
58539b2de37Sjeremylt        """Ceed QFunction By Name: point-wise operation at quadrature points
58639b2de37Sjeremylt             from a given gallery, for evaluating volumetric terms.
58739b2de37Sjeremylt
58839b2de37Sjeremylt           Args:
58939b2de37Sjeremylt             name: name of QFunction to use from gallery
59039b2de37Sjeremylt
59139b2de37Sjeremylt           Returns:
59239b2de37Sjeremylt             qfunction: Ceed QFunction By Name"""
59339b2de37Sjeremylt
59439b2de37Sjeremylt        return QFunctionByName(self, name)
59539b2de37Sjeremylt
59639b2de37Sjeremylt    def IdentityQFunction(self, size, inmode, outmode):
59739b2de37Sjeremylt        """Ceed Idenity QFunction: identity qfunction operation.
59839b2de37Sjeremylt
59939b2de37Sjeremylt           Args:
60039b2de37Sjeremylt             size: size of the qfunction fields
60139b2de37Sjeremylt             **inmode: CeedEvalMode for input to Ceed QFunction
60239b2de37Sjeremylt             **outmode: CeedEvalMode for output to Ceed QFunction
60339b2de37Sjeremylt
60439b2de37Sjeremylt           Returns:
60539b2de37Sjeremylt             qfunction: Ceed Identity QFunction"""
60639b2de37Sjeremylt
60739b2de37Sjeremylt        return IdentityQFunction(self, size, inmode, outmode)
60839b2de37Sjeremylt
609777ff853SJeremy L Thompson    def QFunctionContext(self):
610777ff853SJeremy L Thompson        """Ceed QFunction Context: stores Ceed QFunction user context data.
611777ff853SJeremy L Thompson
612777ff853SJeremy L Thompson           Returns:
613777ff853SJeremy L Thompson             userContext: Ceed QFunction Context"""
614777ff853SJeremy L Thompson
615777ff853SJeremy L Thompson        return QFunctionContext(self)
616777ff853SJeremy L Thompson
61739b2de37Sjeremylt    # CeedOperator
61839b2de37Sjeremylt    def Operator(self, qf, dqf=None, qdfT=None):
61939b2de37Sjeremylt        """Ceed Operator: composed FE-type operations on vectors.
62039b2de37Sjeremylt
62139b2de37Sjeremylt           Args:
622d979a051Sjeremylt             qf: QFunction defining the action of the operator at quadrature
623d979a051Sjeremylt                   points
62439b2de37Sjeremylt             **dqf: QFunction defining the action of the Jacobian, default None
625d979a051Sjeremylt             **dqfT: QFunction defining the action of the transpose of the
626d979a051Sjeremylt                       Jacobian, default None
62739b2de37Sjeremylt
62839b2de37Sjeremylt           Returns:
62939b2de37Sjeremylt             operator: Ceed Operator"""
63039b2de37Sjeremylt
63139b2de37Sjeremylt        return Operator(self, qf, dqf, qdfT)
63239b2de37Sjeremylt
63339b2de37Sjeremylt    def CompositeOperator(self):
63439b2de37Sjeremylt        """Composite Ceed Operator: composition of multiple CeedOperators.
63539b2de37Sjeremylt
63639b2de37Sjeremylt           Returns:
63739b2de37Sjeremylt             operator: Ceed Composite Operator"""
63839b2de37Sjeremylt
63939b2de37Sjeremylt        return CompositeOperator(self)
64039b2de37Sjeremylt
64139b2de37Sjeremylt    # Destructor
64239b2de37Sjeremylt    def __del__(self):
64339b2de37Sjeremylt        # libCEED call
64439b2de37Sjeremylt        lib.CeedDestroy(self._pointer)
64539b2de37Sjeremylt
64639b2de37Sjeremylt# ------------------------------------------------------------------------------
647