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 1939b2de37Sjeremyltimport io 2039b2de37Sjeremyltfrom abc import ABC 2139b2de37Sjeremyltfrom .ceed_vector import Vector 2239b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1 23*69a53589Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, StridedElemRestriction, BlockedElemRestriction, BlockedStridedElemRestriction 2439b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction 2539b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator 2639b2de37Sjeremyltfrom .ceed_constants import * 2739b2de37Sjeremylt 2839b2de37Sjeremylt# ------------------------------------------------------------------------------ 2939b2de37Sjeremyltclass Ceed(): 3039b2de37Sjeremylt """Ceed: core components.""" 3139b2de37Sjeremylt # Attributes 3239b2de37Sjeremylt _pointer = ffi.NULL 3339b2de37Sjeremylt 3439b2de37Sjeremylt # Constructor 3539b2de37Sjeremylt def __init__(self, resource = "/cpu/self"): 3639b2de37Sjeremylt # libCEED object 3739b2de37Sjeremylt self._pointer = ffi.new("Ceed *") 3839b2de37Sjeremylt 3939b2de37Sjeremylt # libCEED call 4039b2de37Sjeremylt resourceAscii = ffi.new("char[]", resource.encode("ascii")) 4139b2de37Sjeremylt lib.CeedInit(resourceAscii, self._pointer) 4239b2de37Sjeremylt 4339b2de37Sjeremylt # Representation 4439b2de37Sjeremylt def __repr__(self): 4539b2de37Sjeremylt return "<Ceed instance at " + hex(id(self)) + ">" 4639b2de37Sjeremylt 4739b2de37Sjeremylt # Get Resource 4839b2de37Sjeremylt def get_resource(self): 4939b2de37Sjeremylt """Get the full resource name for a Ceed context. 5039b2de37Sjeremylt 5139b2de37Sjeremylt Returns: 5239b2de37Sjeremylt resource: resource name""" 5339b2de37Sjeremylt 5439b2de37Sjeremylt # libCEED call 5539b2de37Sjeremylt resource = ffi.new("char **") 5639b2de37Sjeremylt lib.CeedGetResource(self._pointer[0], resource) 5739b2de37Sjeremylt 5839b2de37Sjeremylt return ffi.string(resource[0]).decode("UTF-8") 5939b2de37Sjeremylt 6039b2de37Sjeremylt # Preferred MemType 6139b2de37Sjeremylt def get_preferred_memtype(self): 6239b2de37Sjeremylt """Return Ceed preferred memory type. 6339b2de37Sjeremylt 6439b2de37Sjeremylt Returns: 6539b2de37Sjeremylt memtype: Ceed preferred memory type""" 6639b2de37Sjeremylt 6739b2de37Sjeremylt # libCEED call 6839b2de37Sjeremylt memtype = ffi.new("CeedMemType *", MEM_HOST) 6939b2de37Sjeremylt lib.CeedGetPreferredMemType(self._pointer[0], memtype) 7039b2de37Sjeremylt 7139b2de37Sjeremylt return memtype[0] 7239b2de37Sjeremylt 7339b2de37Sjeremylt # CeedVector 7439b2de37Sjeremylt def Vector(self, size): 7539b2de37Sjeremylt """Ceed Vector: storing and manipulating vectors. 7639b2de37Sjeremylt 7739b2de37Sjeremylt Args: 7839b2de37Sjeremylt size: length of vector 7939b2de37Sjeremylt 8039b2de37Sjeremylt Returns: 8139b2de37Sjeremylt vector: Ceed Vector""" 8239b2de37Sjeremylt 8339b2de37Sjeremylt return Vector(self, size) 8439b2de37Sjeremylt 8539b2de37Sjeremylt # CeedElemRestriction 8639b2de37Sjeremylt def ElemRestriction(self, nelem, elemsize, nnodes, ncomp, indices, 87a8d32208Sjeremylt memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES, 8861dbc9d2Sjeremylt imode=lib.CEED_NONINTERLACED): 8939b2de37Sjeremylt """Ceed ElemRestriction: restriction from local vectors to elements. 9039b2de37Sjeremylt 9139b2de37Sjeremylt Args: 9239b2de37Sjeremylt nelem: number of elements described in the indices array 9339b2de37Sjeremylt elemsize: size (number of nodes) per element 9439b2de37Sjeremylt nnodes: the number of nodes in the local vector. The input Ceed Vector 9539b2de37Sjeremylt to which the restriction will be applied is of size 9639b2de37Sjeremylt (nnodes * ncomp). This size may include data 9739b2de37Sjeremylt used by other Ceed ElemRestriction objects describing 9839b2de37Sjeremylt different types of elements. 9939b2de37Sjeremylt ncomp: number of field components per interpolation node (1 for scalar fields) 10039b2de37Sjeremylt *indices: Numpy array of shape [nelem, elemsize]. Row i holds the 10139b2de37Sjeremylt ordered list of the indices (into the input Ceed Vector) 10239b2de37Sjeremylt for the unknowns corresponding to element i, where 10339b2de37Sjeremylt 0 <= i < nelem. All indices must be in the range 10439b2de37Sjeremylt [0, nnodes - 1]. 10539b2de37Sjeremylt **memtype: memory type of the indices array, default CEED_MEM_HOST 10639b2de37Sjeremylt **cmode: copy mode for the indices array, default CEED_COPY_VALUES 10761dbc9d2Sjeremylt **imode: ordering of the ncomp components, i.e. it specifies 108a8d32208Sjeremylt the ordering of the components of the local vector used 10961dbc9d2Sjeremylt by this CeedElemRestriction. CEED_NONINTERLACED indicates 11061dbc9d2Sjeremylt the component is the outermost index and CEED_INTERLACED 111a8d32208Sjeremylt indicates the component is the innermost index in 11261dbc9d2Sjeremylt ordering of the local vector, default CEED_NONINTERLACED 11339b2de37Sjeremylt 11439b2de37Sjeremylt Returns: 11539b2de37Sjeremylt elemrestriction: Ceed ElemRestiction""" 11639b2de37Sjeremylt 11739b2de37Sjeremylt return ElemRestriction(self, nelem, elemsize, nnodes, ncomp, indices, 11861dbc9d2Sjeremylt memtype=memtype, cmode=cmode, imode=imode) 11939b2de37Sjeremylt 120*69a53589Sjeremylt def StridedElemRestriction(self, nelem, elemsize, nnodes, ncomp, strides): 121*69a53589Sjeremylt """Ceed Identity ElemRestriction: strided restriction from local vectors to elements. 12239b2de37Sjeremylt 12339b2de37Sjeremylt Args: 12439b2de37Sjeremylt nelem: number of elements described in the indices array 12539b2de37Sjeremylt elemsize: size (number of nodes) per element 12639b2de37Sjeremylt nnodes: the number of nodes in the local vector. The input Ceed Vector 12739b2de37Sjeremylt to which the restriction will be applied is of size 12839b2de37Sjeremylt (nnodes * ncomp). This size may include data 12939b2de37Sjeremylt used by other Ceed ElemRestriction objects describing 13039b2de37Sjeremylt different types of elements. 131a8d32208Sjeremylt ncomp: number of field components per interpolation node 132a8d32208Sjeremylt (1 for scalar fields) 133*69a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 134*69a53589Sjeremylt The data for node i, component j, element k in the 135*69a53589Sjeremylt L-vector is given by 136*69a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 13739b2de37Sjeremylt 13839b2de37Sjeremylt Returns: 139*69a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 14039b2de37Sjeremylt 141*69a53589Sjeremylt return StridedElemRestriction(self, nelem, elemsize, nnodes, ncomp, strides) 14239b2de37Sjeremylt 14339b2de37Sjeremylt def BlockedElemRestriction(self, nelem, elemsize, blksize, nnodes, ncomp, 14439b2de37Sjeremylt indices, memtype=lib.CEED_MEM_HOST, 145a8d32208Sjeremylt cmode=lib.CEED_COPY_VALUES, 14661dbc9d2Sjeremylt imode=lib.CEED_NONINTERLACED): 14739b2de37Sjeremylt """Ceed Blocked ElemRestriction: blocked restriction from local vectors to elements. 14839b2de37Sjeremylt 14939b2de37Sjeremylt Args: 15039b2de37Sjeremylt nelem: number of elements described in the indices array 15139b2de37Sjeremylt elemsize: size (number of nodes) per element 15239b2de37Sjeremylt blksize: number of elements in a block 15339b2de37Sjeremylt nnodes: the number of nodes in the local vector. The input Ceed Vector 15439b2de37Sjeremylt to which the restriction will be applied is of size 15539b2de37Sjeremylt (nnodes * ncomp). This size may include data 15639b2de37Sjeremylt used by other Ceed ElemRestriction objects describing 15739b2de37Sjeremylt different types of elements. 15839b2de37Sjeremylt ncomp: number of field components per interpolation node (1 for scalar fields) 15939b2de37Sjeremylt *indices: Numpy array of shape [nelem, elemsize]. Row i holds the 16039b2de37Sjeremylt ordered list of the indices (into the input Ceed Vector) 16139b2de37Sjeremylt for the unknowns corresponding to element i, where 16239b2de37Sjeremylt 0 <= i < nelem. All indices must be in the range 16339b2de37Sjeremylt [0, nnodes - 1]. The backend will permute and pad this 16439b2de37Sjeremylt array to the desired ordering for the blocksize, which is 16539b2de37Sjeremylt typically given by the backend. The default reordering is 16639b2de37Sjeremylt to interlace elements. 16739b2de37Sjeremylt **memtype: memory type of the indices array, default CEED_MEM_HOST 16839b2de37Sjeremylt **cmode: copy mode for the indices array, default CEED_COPY_VALUES 16961dbc9d2Sjeremylt **imode: ordering of the ncomp components, i.e. it specifies 170a8d32208Sjeremylt the ordering of the components of the local vector used 17161dbc9d2Sjeremylt by this CeedElemRestriction. CEED_NONINTERLACED indicates 17261dbc9d2Sjeremylt the component is the outermost index and CEED_INTERLACED 173a8d32208Sjeremylt indicates the component is the innermost index in 17461dbc9d2Sjeremylt ordering of the local vector, default CEED_NONINTERLACED 17539b2de37Sjeremylt 17639b2de37Sjeremylt Returns: 17739b2de37Sjeremylt elemrestriction: Ceed Blocked ElemRestiction""" 17839b2de37Sjeremylt 17939b2de37Sjeremylt return BlockedElemRestriction(self, nelem, elemsize, blksize, nnodes, 18039b2de37Sjeremylt ncomp, indices, memtype=memtype, 18161dbc9d2Sjeremylt cmode=cmode, imode=imode) 18239b2de37Sjeremylt 183*69a53589Sjeremylt def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, nnodes, 184*69a53589Sjeremylt ncomp, strides): 185*69a53589Sjeremylt """Ceed Blocked Strided ElemRestriction: blocked and strided restriction from local vectors to elements. 186*69a53589Sjeremylt 187*69a53589Sjeremylt Args: 188*69a53589Sjeremylt nelem: number of elements described in the indices array 189*69a53589Sjeremylt elemsize: size (number of nodes) per element 190*69a53589Sjeremylt blksize: number of elements in a block 191*69a53589Sjeremylt nnodes: the number of nodes in the local vector. The input Ceed Vector 192*69a53589Sjeremylt to which the restriction will be applied is of size 193*69a53589Sjeremylt (nnodes * ncomp). This size may include data 194*69a53589Sjeremylt used by other Ceed ElemRestriction objects describing 195*69a53589Sjeremylt different types of elements. 196*69a53589Sjeremylt ncomp: number of field components per interpolation node 197*69a53589Sjeremylt (1 for scalar fields) 198*69a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 199*69a53589Sjeremylt The data for node i, component j, element k in the 200*69a53589Sjeremylt L-vector is given by 201*69a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 202*69a53589Sjeremylt 203*69a53589Sjeremylt Returns: 204*69a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 205*69a53589Sjeremylt 206*69a53589Sjeremylt return BlockedStridedElemRestriction(self, nelem, elemsize, blksize, 207*69a53589Sjeremylt nnodes, ncomp, strides) 208*69a53589Sjeremylt 20939b2de37Sjeremylt # CeedBasis 21039b2de37Sjeremylt def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 21139b2de37Sjeremylt qref1d, qweight1d): 21239b2de37Sjeremylt """Ceed Tensor H1 Basis: finite element tensor-product basis objects for 21339b2de37Sjeremylt H^1 discretizations. 21439b2de37Sjeremylt 21539b2de37Sjeremylt Args: 21639b2de37Sjeremylt dim: topological dimension 21739b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 21839b2de37Sjeremylt P1d: number of nodes in one dimension 21939b2de37Sjeremylt Q1d: number of quadrature points in one dimension 22039b2de37Sjeremylt *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix expressing the 22139b2de37Sjeremylt values of nodal basis functions at quadrature points 22239b2de37Sjeremylt *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix expressing the 22339b2de37Sjeremylt derivatives of nodal basis functions at quadrature points 22439b2de37Sjeremylt *qref1d: Numpy array of length Q1d holding the locations of quadrature points 22539b2de37Sjeremylt on the 1D reference element [-1, 1] 22639b2de37Sjeremylt *qweight1d: Numpy array of length Q1d holding the quadrature weights on the 22739b2de37Sjeremylt reference element 22839b2de37Sjeremylt 22939b2de37Sjeremylt Returns: 23039b2de37Sjeremylt basis: Ceed Basis""" 23139b2de37Sjeremylt 23239b2de37Sjeremylt return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 23339b2de37Sjeremylt qref1d, qweight1d) 23439b2de37Sjeremylt 23539b2de37Sjeremylt def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode): 23639b2de37Sjeremylt """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange basis 23739b2de37Sjeremylt objects for H^1 discretizations. 23839b2de37Sjeremylt 23939b2de37Sjeremylt Args: 24039b2de37Sjeremylt dim: topological dimension 24139b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 24239b2de37Sjeremylt P: number of Gauss-Lobatto nodes in one dimension. The 24339b2de37Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 24439b2de37Sjeremylt Q: number of quadrature points in one dimension 24539b2de37Sjeremylt qmode: distribution of the Q quadrature points (affects order of 24639b2de37Sjeremylt accuracy for the quadrature) 24739b2de37Sjeremylt 24839b2de37Sjeremylt Returns: 24939b2de37Sjeremylt basis: Ceed Basis""" 25039b2de37Sjeremylt 25139b2de37Sjeremylt return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode) 25239b2de37Sjeremylt 25339b2de37Sjeremylt def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight): 25439b2de37Sjeremylt """Ceed H1 Basis: finite element non tensor-product basis for H^1 discretizations. 25539b2de37Sjeremylt 25639b2de37Sjeremylt Args: 25739b2de37Sjeremylt topo: topology of the element, e.g. hypercube, simplex, etc 25839b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 25939b2de37Sjeremylt nnodes: total number of nodes 26039b2de37Sjeremylt nqpts: total number of quadrature points 26139b2de37Sjeremylt *interp: Numpy array holding the row-major (nqpts * nnodes) matrix expressing 26239b2de37Sjeremylt the values of nodal basis functions at quadrature points 26339b2de37Sjeremylt *grad: Numpy array holding the row-major (nqpts * dim * nnodes) matrix 26439b2de37Sjeremylt expressing the derivatives of nodal basis functions at 26539b2de37Sjeremylt quadrature points 26639b2de37Sjeremylt *qref: Numpy array of length (nqpts * dim) holding the locations of quadrature 26739b2de37Sjeremylt points on the reference element [-1, 1] 26839b2de37Sjeremylt *qweight: Numpy array of length nnodes holding the quadrature weights on the 26939b2de37Sjeremylt reference element 27039b2de37Sjeremylt 27139b2de37Sjeremylt Returns: 27239b2de37Sjeremylt basis: Ceed Basis""" 27339b2de37Sjeremylt 27439b2de37Sjeremylt return BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight) 27539b2de37Sjeremylt 27639b2de37Sjeremylt # CeedQFunction 27739b2de37Sjeremylt def QFunction(self, vlength, f, source): 27839b2de37Sjeremylt """Ceed QFunction: point-wise operation at quadrature points for evaluating 27939b2de37Sjeremylt volumetric terms. 28039b2de37Sjeremylt 28139b2de37Sjeremylt Args: 28239b2de37Sjeremylt vlength: vector length. Caller must ensure that number of quadrature 28339b2de37Sjeremylt points is a multiple of vlength 28439b2de37Sjeremylt f: ctypes function pointer to evaluate action at quadrature points 28539b2de37Sjeremylt source: absolute path to source of QFunction, 28639b2de37Sjeremylt "\\abs_path\\file.h:function_name 28739b2de37Sjeremylt 28839b2de37Sjeremylt Returns: 28939b2de37Sjeremylt qfunction: Ceed QFunction""" 29039b2de37Sjeremylt 29139b2de37Sjeremylt return QFunction(self, vlength, f, source) 29239b2de37Sjeremylt 29339b2de37Sjeremylt def QFunctionByName(self, name): 29439b2de37Sjeremylt """Ceed QFunction By Name: point-wise operation at quadrature points 29539b2de37Sjeremylt from a given gallery, for evaluating volumetric terms. 29639b2de37Sjeremylt 29739b2de37Sjeremylt Args: 29839b2de37Sjeremylt name: name of QFunction to use from gallery 29939b2de37Sjeremylt 30039b2de37Sjeremylt Returns: 30139b2de37Sjeremylt qfunction: Ceed QFunction By Name""" 30239b2de37Sjeremylt 30339b2de37Sjeremylt return QFunctionByName(self, name) 30439b2de37Sjeremylt 30539b2de37Sjeremylt def IdentityQFunction(self, size, inmode, outmode): 30639b2de37Sjeremylt """Ceed Idenity QFunction: identity qfunction operation. 30739b2de37Sjeremylt 30839b2de37Sjeremylt Args: 30939b2de37Sjeremylt size: size of the qfunction fields 31039b2de37Sjeremylt **inmode: CeedEvalMode for input to Ceed QFunction 31139b2de37Sjeremylt **outmode: CeedEvalMode for output to Ceed QFunction 31239b2de37Sjeremylt 31339b2de37Sjeremylt Returns: 31439b2de37Sjeremylt qfunction: Ceed Identity QFunction""" 31539b2de37Sjeremylt 31639b2de37Sjeremylt return IdentityQFunction(self, size, inmode, outmode) 31739b2de37Sjeremylt 31839b2de37Sjeremylt # CeedOperator 31939b2de37Sjeremylt def Operator(self, qf, dqf=None, qdfT=None): 32039b2de37Sjeremylt """Ceed Operator: composed FE-type operations on vectors. 32139b2de37Sjeremylt 32239b2de37Sjeremylt Args: 32339b2de37Sjeremylt qf: QFunction defining the action of the operator at quadrature points 32439b2de37Sjeremylt **dqf: QFunction defining the action of the Jacobian, default None 32539b2de37Sjeremylt **dqfT: QFunction defining the action of the transpose of the Jacobian, 32639b2de37Sjeremylt default None 32739b2de37Sjeremylt 32839b2de37Sjeremylt Returns: 32939b2de37Sjeremylt operator: Ceed Operator""" 33039b2de37Sjeremylt 33139b2de37Sjeremylt return Operator(self, qf, dqf, qdfT) 33239b2de37Sjeremylt 33339b2de37Sjeremylt def CompositeOperator(self): 33439b2de37Sjeremylt """Composite Ceed Operator: composition of multiple CeedOperators. 33539b2de37Sjeremylt 33639b2de37Sjeremylt Returns: 33739b2de37Sjeremylt operator: Ceed Composite Operator""" 33839b2de37Sjeremylt 33939b2de37Sjeremylt return CompositeOperator(self) 34039b2de37Sjeremylt 34139b2de37Sjeremylt # Destructor 34239b2de37Sjeremylt def __del__(self): 34339b2de37Sjeremylt # libCEED call 34439b2de37Sjeremylt lib.CeedDestroy(self._pointer) 34539b2de37Sjeremylt 34639b2de37Sjeremylt# ------------------------------------------------------------------------------ 347