xref: /libCEED/python/ceed.py (revision 39b2de37682296be8460181179eb4e44de5cc3de)
1*39b2de37Sjeremylt# Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2*39b2de37Sjeremylt# the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3*39b2de37Sjeremylt# reserved. See files LICENSE and NOTICE for details.
4*39b2de37Sjeremylt#
5*39b2de37Sjeremylt# This file is part of CEED, a collection of benchmarks, miniapps, software
6*39b2de37Sjeremylt# libraries and APIs for efficient high-order finite element and spectral
7*39b2de37Sjeremylt# element discretizations for exascale applications. For more information and
8*39b2de37Sjeremylt# source code availability see http://github.com/ceed.
9*39b2de37Sjeremylt#
10*39b2de37Sjeremylt# The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11*39b2de37Sjeremylt# a collaborative effort of two U.S. Department of Energy organizations (Office
12*39b2de37Sjeremylt# of Science and the National Nuclear Security Administration) responsible for
13*39b2de37Sjeremylt# the planning and preparation of a capable exascale ecosystem, including
14*39b2de37Sjeremylt# software, applications, hardware, advanced system engineering and early
15*39b2de37Sjeremylt# testbed platforms, in support of the nation's exascale computing imperative.
16*39b2de37Sjeremylt
17*39b2de37Sjeremyltfrom _ceed_cffi import ffi, lib
18*39b2de37Sjeremyltimport sys
19*39b2de37Sjeremyltimport io
20*39b2de37Sjeremyltfrom abc import ABC
21*39b2de37Sjeremyltfrom .ceed_vector import Vector
22*39b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1
23*39b2de37Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, IdentityElemRestriction, BlockedElemRestriction
24*39b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction
25*39b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator
26*39b2de37Sjeremyltfrom .ceed_constants import *
27*39b2de37Sjeremylt
28*39b2de37Sjeremylt# ------------------------------------------------------------------------------
29*39b2de37Sjeremyltclass Ceed():
30*39b2de37Sjeremylt  """Ceed: core components."""
31*39b2de37Sjeremylt  # Attributes
32*39b2de37Sjeremylt  _pointer = ffi.NULL
33*39b2de37Sjeremylt
34*39b2de37Sjeremylt  # Constructor
35*39b2de37Sjeremylt  def __init__(self, resource = "/cpu/self"):
36*39b2de37Sjeremylt    # libCEED object
37*39b2de37Sjeremylt    self._pointer = ffi.new("Ceed *")
38*39b2de37Sjeremylt
39*39b2de37Sjeremylt    # libCEED call
40*39b2de37Sjeremylt    resourceAscii = ffi.new("char[]", resource.encode("ascii"))
41*39b2de37Sjeremylt    lib.CeedInit(resourceAscii, self._pointer)
42*39b2de37Sjeremylt
43*39b2de37Sjeremylt  # Representation
44*39b2de37Sjeremylt  def __repr__(self):
45*39b2de37Sjeremylt    return "<Ceed instance at " + hex(id(self)) + ">"
46*39b2de37Sjeremylt
47*39b2de37Sjeremylt  # Get Resource
48*39b2de37Sjeremylt  def get_resource(self):
49*39b2de37Sjeremylt    """Get the full resource name for a Ceed context.
50*39b2de37Sjeremylt
51*39b2de37Sjeremylt       Returns:
52*39b2de37Sjeremylt         resource: resource name"""
53*39b2de37Sjeremylt
54*39b2de37Sjeremylt    # libCEED call
55*39b2de37Sjeremylt    resource = ffi.new("char **")
56*39b2de37Sjeremylt    lib.CeedGetResource(self._pointer[0], resource)
57*39b2de37Sjeremylt
58*39b2de37Sjeremylt    return ffi.string(resource[0]).decode("UTF-8")
59*39b2de37Sjeremylt
60*39b2de37Sjeremylt  # Preferred MemType
61*39b2de37Sjeremylt  def get_preferred_memtype(self):
62*39b2de37Sjeremylt    """Return Ceed preferred memory type.
63*39b2de37Sjeremylt
64*39b2de37Sjeremylt       Returns:
65*39b2de37Sjeremylt         memtype: Ceed preferred memory type"""
66*39b2de37Sjeremylt
67*39b2de37Sjeremylt    # libCEED call
68*39b2de37Sjeremylt    memtype = ffi.new("CeedMemType *", MEM_HOST)
69*39b2de37Sjeremylt    lib.CeedGetPreferredMemType(self._pointer[0], memtype)
70*39b2de37Sjeremylt
71*39b2de37Sjeremylt    return memtype[0]
72*39b2de37Sjeremylt
73*39b2de37Sjeremylt  # CeedVector
74*39b2de37Sjeremylt  def Vector(self, size):
75*39b2de37Sjeremylt    """Ceed Vector: storing and manipulating vectors.
76*39b2de37Sjeremylt
77*39b2de37Sjeremylt       Args:
78*39b2de37Sjeremylt         size: length of vector
79*39b2de37Sjeremylt
80*39b2de37Sjeremylt       Returns:
81*39b2de37Sjeremylt         vector: Ceed Vector"""
82*39b2de37Sjeremylt
83*39b2de37Sjeremylt    return Vector(self, size)
84*39b2de37Sjeremylt
85*39b2de37Sjeremylt  # CeedElemRestriction
86*39b2de37Sjeremylt  def ElemRestriction(self, nelem, elemsize, nnodes, ncomp, indices,
87*39b2de37Sjeremylt                      memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES):
88*39b2de37Sjeremylt    """Ceed ElemRestriction: restriction from local vectors to elements.
89*39b2de37Sjeremylt
90*39b2de37Sjeremylt       Args:
91*39b2de37Sjeremylt         nelem: number of elements described in the indices array
92*39b2de37Sjeremylt         elemsize: size (number of nodes) per element
93*39b2de37Sjeremylt         nnodes: the number of nodes in the local vector. The input Ceed Vector
94*39b2de37Sjeremylt                   to which the restriction will be applied is of size
95*39b2de37Sjeremylt                   (nnodes * ncomp). This size may include data
96*39b2de37Sjeremylt                   used by other Ceed ElemRestriction objects describing
97*39b2de37Sjeremylt                   different types of elements.
98*39b2de37Sjeremylt         ncomp: number of field components per interpolation node (1 for scalar fields)
99*39b2de37Sjeremylt         *indices: Numpy array of shape [nelem, elemsize]. Row i holds the
100*39b2de37Sjeremylt                      ordered list of the indices (into the input Ceed Vector)
101*39b2de37Sjeremylt                      for the unknowns corresponding to element i, where
102*39b2de37Sjeremylt                      0 <= i < nelem. All indices must be in the range
103*39b2de37Sjeremylt                      [0, nnodes - 1].
104*39b2de37Sjeremylt         **memtype: memory type of the indices array, default CEED_MEM_HOST
105*39b2de37Sjeremylt         **cmode: copy mode for the indices array, default CEED_COPY_VALUES
106*39b2de37Sjeremylt
107*39b2de37Sjeremylt       Returns:
108*39b2de37Sjeremylt         elemrestriction: Ceed ElemRestiction"""
109*39b2de37Sjeremylt
110*39b2de37Sjeremylt    return ElemRestriction(self, nelem, elemsize, nnodes, ncomp, indices,
111*39b2de37Sjeremylt                           memtype=memtype, cmode=cmode)
112*39b2de37Sjeremylt
113*39b2de37Sjeremylt  def IdentityElemRestriction(self, nelem, elemsize, nnodes, ncomp):
114*39b2de37Sjeremylt    """Ceed Identity ElemRestriction: identity restriction from local vectors to elements.
115*39b2de37Sjeremylt
116*39b2de37Sjeremylt       Args:
117*39b2de37Sjeremylt         nelem: number of elements described in the indices array
118*39b2de37Sjeremylt         elemsize: size (number of nodes) per element
119*39b2de37Sjeremylt         nnodes: the number of nodes in the local vector. The input Ceed Vector
120*39b2de37Sjeremylt                   to which the restriction will be applied is of size
121*39b2de37Sjeremylt                   (nnodes * ncomp). This size may include data
122*39b2de37Sjeremylt                   used by other Ceed ElemRestriction objects describing
123*39b2de37Sjeremylt                   different types of elements.
124*39b2de37Sjeremylt         ncomp: number of field components per interpolation node (1 for scalar fields)
125*39b2de37Sjeremylt
126*39b2de37Sjeremylt       Returns:
127*39b2de37Sjeremylt         elemrestriction: Ceed Identity ElemRestiction"""
128*39b2de37Sjeremylt
129*39b2de37Sjeremylt    return IdentityElemRestriction(self, nelem, elemsize, nnodes, ncomp)
130*39b2de37Sjeremylt
131*39b2de37Sjeremylt  def BlockedElemRestriction(self, nelem, elemsize, blksize, nnodes, ncomp,
132*39b2de37Sjeremylt                             indices, memtype=lib.CEED_MEM_HOST,
133*39b2de37Sjeremylt                             cmode=lib.CEED_COPY_VALUES):
134*39b2de37Sjeremylt    """Ceed Blocked ElemRestriction: blocked restriction from local vectors to elements.
135*39b2de37Sjeremylt
136*39b2de37Sjeremylt       Args:
137*39b2de37Sjeremylt         nelem: number of elements described in the indices array
138*39b2de37Sjeremylt         elemsize: size (number of nodes) per element
139*39b2de37Sjeremylt         blksize: number of elements in a block
140*39b2de37Sjeremylt         nnodes: the number of nodes in the local vector. The input Ceed Vector
141*39b2de37Sjeremylt                   to which the restriction will be applied is of size
142*39b2de37Sjeremylt                   (nnodes * ncomp). This size may include data
143*39b2de37Sjeremylt                   used by other Ceed ElemRestriction objects describing
144*39b2de37Sjeremylt                   different types of elements.
145*39b2de37Sjeremylt         ncomp: number of field components per interpolation node (1 for scalar fields)
146*39b2de37Sjeremylt         *indices: Numpy array of shape [nelem, elemsize]. Row i holds the
147*39b2de37Sjeremylt                      ordered list of the indices (into the input Ceed Vector)
148*39b2de37Sjeremylt                      for the unknowns corresponding to element i, where
149*39b2de37Sjeremylt                      0 <= i < nelem. All indices must be in the range
150*39b2de37Sjeremylt                      [0, nnodes - 1]. The backend will permute and pad this
151*39b2de37Sjeremylt                      array to the desired ordering for the blocksize, which is
152*39b2de37Sjeremylt                      typically given by the backend. The default reordering is
153*39b2de37Sjeremylt                      to interlace elements.
154*39b2de37Sjeremylt         **memtype: memory type of the indices array, default CEED_MEM_HOST
155*39b2de37Sjeremylt         **cmode: copy mode for the indices array, default CEED_COPY_VALUES
156*39b2de37Sjeremylt
157*39b2de37Sjeremylt       Returns:
158*39b2de37Sjeremylt         elemrestriction: Ceed Blocked ElemRestiction"""
159*39b2de37Sjeremylt
160*39b2de37Sjeremylt    return BlockedElemRestriction(self, nelem, elemsize, blksize, nnodes,
161*39b2de37Sjeremylt                                  ncomp, indices, memtype=memtype,
162*39b2de37Sjeremylt                                  cmode=cmode)
163*39b2de37Sjeremylt
164*39b2de37Sjeremylt  # CeedBasis
165*39b2de37Sjeremylt  def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
166*39b2de37Sjeremylt                    qref1d, qweight1d):
167*39b2de37Sjeremylt    """Ceed Tensor H1 Basis: finite element tensor-product basis objects for
168*39b2de37Sjeremylt         H^1 discretizations.
169*39b2de37Sjeremylt
170*39b2de37Sjeremylt       Args:
171*39b2de37Sjeremylt         dim: topological dimension
172*39b2de37Sjeremylt         ncomp: number of field components (1 for scalar fields)
173*39b2de37Sjeremylt         P1d: number of nodes in one dimension
174*39b2de37Sjeremylt         Q1d: number of quadrature points in one dimension
175*39b2de37Sjeremylt         *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix expressing the
176*39b2de37Sjeremylt                     values of nodal basis functions at quadrature points
177*39b2de37Sjeremylt         *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix expressing the
178*39b2de37Sjeremylt                   derivatives of nodal basis functions at quadrature points
179*39b2de37Sjeremylt         *qref1d: Numpy array of length Q1d holding the locations of quadrature points
180*39b2de37Sjeremylt                   on the 1D reference element [-1, 1]
181*39b2de37Sjeremylt         *qweight1d: Numpy array of length Q1d holding the quadrature weights on the
182*39b2de37Sjeremylt                      reference element
183*39b2de37Sjeremylt
184*39b2de37Sjeremylt       Returns:
185*39b2de37Sjeremylt         basis: Ceed Basis"""
186*39b2de37Sjeremylt
187*39b2de37Sjeremylt    return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d,
188*39b2de37Sjeremylt                         qref1d, qweight1d)
189*39b2de37Sjeremylt
190*39b2de37Sjeremylt  def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode):
191*39b2de37Sjeremylt    """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange basis
192*39b2de37Sjeremylt         objects for H^1 discretizations.
193*39b2de37Sjeremylt
194*39b2de37Sjeremylt       Args:
195*39b2de37Sjeremylt         dim: topological dimension
196*39b2de37Sjeremylt         ncomp: number of field components (1 for scalar fields)
197*39b2de37Sjeremylt         P: number of Gauss-Lobatto nodes in one dimension.  The
198*39b2de37Sjeremylt              polynomial degree of the resulting Q_k element is k=P-1.
199*39b2de37Sjeremylt         Q: number of quadrature points in one dimension
200*39b2de37Sjeremylt         qmode: distribution of the Q quadrature points (affects order of
201*39b2de37Sjeremylt                  accuracy for the quadrature)
202*39b2de37Sjeremylt
203*39b2de37Sjeremylt       Returns:
204*39b2de37Sjeremylt         basis: Ceed Basis"""
205*39b2de37Sjeremylt
206*39b2de37Sjeremylt    return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode)
207*39b2de37Sjeremylt
208*39b2de37Sjeremylt  def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight):
209*39b2de37Sjeremylt    """Ceed H1 Basis: finite element non tensor-product basis for H^1 discretizations.
210*39b2de37Sjeremylt
211*39b2de37Sjeremylt       Args:
212*39b2de37Sjeremylt         topo: topology of the element, e.g. hypercube, simplex, etc
213*39b2de37Sjeremylt         ncomp: number of field components (1 for scalar fields)
214*39b2de37Sjeremylt         nnodes: total number of nodes
215*39b2de37Sjeremylt         nqpts: total number of quadrature points
216*39b2de37Sjeremylt         *interp: Numpy array holding the row-major (nqpts * nnodes) matrix expressing
217*39b2de37Sjeremylt                   the values of nodal basis functions at quadrature points
218*39b2de37Sjeremylt         *grad: Numpy array holding the row-major (nqpts * dim * nnodes) matrix
219*39b2de37Sjeremylt                 expressing the derivatives of nodal basis functions at
220*39b2de37Sjeremylt                 quadrature points
221*39b2de37Sjeremylt         *qref: Numpy array of length (nqpts * dim) holding the locations of quadrature
222*39b2de37Sjeremylt                 points on the reference element [-1, 1]
223*39b2de37Sjeremylt         *qweight: Numpy array of length nnodes holding the quadrature weights on the
224*39b2de37Sjeremylt                    reference element
225*39b2de37Sjeremylt
226*39b2de37Sjeremylt       Returns:
227*39b2de37Sjeremylt         basis: Ceed Basis"""
228*39b2de37Sjeremylt
229*39b2de37Sjeremylt    return BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight)
230*39b2de37Sjeremylt
231*39b2de37Sjeremylt  # CeedQFunction
232*39b2de37Sjeremylt  def QFunction(self, vlength, f, source):
233*39b2de37Sjeremylt    """Ceed QFunction: point-wise operation at quadrature points for evaluating
234*39b2de37Sjeremylt         volumetric terms.
235*39b2de37Sjeremylt
236*39b2de37Sjeremylt       Args:
237*39b2de37Sjeremylt         vlength: vector length. Caller must ensure that number of quadrature
238*39b2de37Sjeremylt                    points is a multiple of vlength
239*39b2de37Sjeremylt         f: ctypes function pointer to evaluate action at quadrature points
240*39b2de37Sjeremylt         source: absolute path to source of QFunction,
241*39b2de37Sjeremylt           "\\abs_path\\file.h:function_name
242*39b2de37Sjeremylt
243*39b2de37Sjeremylt       Returns:
244*39b2de37Sjeremylt         qfunction: Ceed QFunction"""
245*39b2de37Sjeremylt
246*39b2de37Sjeremylt    return QFunction(self, vlength, f, source)
247*39b2de37Sjeremylt
248*39b2de37Sjeremylt  def QFunctionByName(self, name):
249*39b2de37Sjeremylt    """Ceed QFunction By Name: point-wise operation at quadrature points
250*39b2de37Sjeremylt         from a given gallery, for evaluating volumetric terms.
251*39b2de37Sjeremylt
252*39b2de37Sjeremylt       Args:
253*39b2de37Sjeremylt         name: name of QFunction to use from gallery
254*39b2de37Sjeremylt
255*39b2de37Sjeremylt       Returns:
256*39b2de37Sjeremylt         qfunction: Ceed QFunction By Name"""
257*39b2de37Sjeremylt
258*39b2de37Sjeremylt    return QFunctionByName(self, name)
259*39b2de37Sjeremylt
260*39b2de37Sjeremylt  def IdentityQFunction(self, size, inmode, outmode):
261*39b2de37Sjeremylt    """Ceed Idenity QFunction: identity qfunction operation.
262*39b2de37Sjeremylt
263*39b2de37Sjeremylt       Args:
264*39b2de37Sjeremylt         size: size of the qfunction fields
265*39b2de37Sjeremylt         **inmode: CeedEvalMode for input to Ceed QFunction
266*39b2de37Sjeremylt         **outmode: CeedEvalMode for output to Ceed QFunction
267*39b2de37Sjeremylt
268*39b2de37Sjeremylt       Returns:
269*39b2de37Sjeremylt         qfunction: Ceed Identity QFunction"""
270*39b2de37Sjeremylt
271*39b2de37Sjeremylt    return IdentityQFunction(self, size, inmode, outmode)
272*39b2de37Sjeremylt
273*39b2de37Sjeremylt  # CeedOperator
274*39b2de37Sjeremylt  def Operator(self, qf, dqf=None, qdfT=None):
275*39b2de37Sjeremylt    """Ceed Operator: composed FE-type operations on vectors.
276*39b2de37Sjeremylt
277*39b2de37Sjeremylt       Args:
278*39b2de37Sjeremylt         qf: QFunction defining the action of the operator at quadrature points
279*39b2de37Sjeremylt         **dqf: QFunction defining the action of the Jacobian, default None
280*39b2de37Sjeremylt         **dqfT: QFunction defining the action of the transpose of the Jacobian,
281*39b2de37Sjeremylt                   default None
282*39b2de37Sjeremylt
283*39b2de37Sjeremylt       Returns:
284*39b2de37Sjeremylt         operator: Ceed Operator"""
285*39b2de37Sjeremylt
286*39b2de37Sjeremylt    return Operator(self, qf, dqf, qdfT)
287*39b2de37Sjeremylt
288*39b2de37Sjeremylt  def CompositeOperator(self):
289*39b2de37Sjeremylt    """Composite Ceed Operator: composition of multiple CeedOperators.
290*39b2de37Sjeremylt
291*39b2de37Sjeremylt       Returns:
292*39b2de37Sjeremylt         operator: Ceed Composite Operator"""
293*39b2de37Sjeremylt
294*39b2de37Sjeremylt    return CompositeOperator(self)
295*39b2de37Sjeremylt
296*39b2de37Sjeremylt  # Destructor
297*39b2de37Sjeremylt  def __del__(self):
298*39b2de37Sjeremylt    # libCEED call
299*39b2de37Sjeremylt    lib.CeedDestroy(self._pointer)
300*39b2de37Sjeremylt
301*39b2de37Sjeremylt# ------------------------------------------------------------------------------
302