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 200a0da059Sjeremyltimport tempfile 2139b2de37Sjeremyltfrom abc import ABC 2239b2de37Sjeremyltfrom .ceed_vector import Vector 2339b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1 2469a53589Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, StridedElemRestriction, BlockedElemRestriction, BlockedStridedElemRestriction 2539b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction 2639b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator 2739b2de37Sjeremyltfrom .ceed_constants import * 2839b2de37Sjeremylt 2939b2de37Sjeremylt# ------------------------------------------------------------------------------ 307a7b0fa3SJed Brown 317a7b0fa3SJed Brown 3239b2de37Sjeremyltclass Ceed(): 3339b2de37Sjeremylt """Ceed: core components.""" 3439b2de37Sjeremylt # Attributes 3539b2de37Sjeremylt _pointer = ffi.NULL 3639b2de37Sjeremylt 3739b2de37Sjeremylt # Constructor 3839b2de37Sjeremylt def __init__(self, resource="/cpu/self"): 3939b2de37Sjeremylt # libCEED object 4039b2de37Sjeremylt self._pointer = ffi.new("Ceed *") 4139b2de37Sjeremylt 4239b2de37Sjeremylt # libCEED call 4339b2de37Sjeremylt resourceAscii = ffi.new("char[]", resource.encode("ascii")) 4439b2de37Sjeremylt lib.CeedInit(resourceAscii, self._pointer) 4539b2de37Sjeremylt 4639b2de37Sjeremylt # Representation 4739b2de37Sjeremylt def __repr__(self): 4839b2de37Sjeremylt return "<Ceed instance at " + hex(id(self)) + ">" 4939b2de37Sjeremylt 500a0da059Sjeremylt # String conversion for print() to stdout 510a0da059Sjeremylt def __str__(self): 520a0da059Sjeremylt """View a Ceed via print().""" 530a0da059Sjeremylt 540a0da059Sjeremylt # libCEED call 550a0da059Sjeremylt with tempfile.NamedTemporaryFile() as key_file: 560a0da059Sjeremylt with open(key_file.name, 'r+') as stream_file: 570a0da059Sjeremylt stream = ffi.cast("FILE *", stream_file) 580a0da059Sjeremylt 590a0da059Sjeremylt lib.CeedView(self._pointer[0], stream) 600a0da059Sjeremylt 610a0da059Sjeremylt stream_file.seek(0) 620a0da059Sjeremylt out_string = stream_file.read() 630a0da059Sjeremylt 640a0da059Sjeremylt return out_string 650a0da059Sjeremylt 6639b2de37Sjeremylt # Get Resource 6739b2de37Sjeremylt def get_resource(self): 6839b2de37Sjeremylt """Get the full resource name for a Ceed context. 6939b2de37Sjeremylt 7039b2de37Sjeremylt Returns: 7139b2de37Sjeremylt resource: resource name""" 7239b2de37Sjeremylt 7339b2de37Sjeremylt # libCEED call 7439b2de37Sjeremylt resource = ffi.new("char **") 7539b2de37Sjeremylt lib.CeedGetResource(self._pointer[0], resource) 7639b2de37Sjeremylt 7739b2de37Sjeremylt return ffi.string(resource[0]).decode("UTF-8") 7839b2de37Sjeremylt 7939b2de37Sjeremylt # Preferred MemType 8039b2de37Sjeremylt def get_preferred_memtype(self): 8139b2de37Sjeremylt """Return Ceed preferred memory type. 8239b2de37Sjeremylt 8339b2de37Sjeremylt Returns: 8439b2de37Sjeremylt memtype: Ceed preferred memory type""" 8539b2de37Sjeremylt 8639b2de37Sjeremylt # libCEED call 8739b2de37Sjeremylt memtype = ffi.new("CeedMemType *", MEM_HOST) 8839b2de37Sjeremylt lib.CeedGetPreferredMemType(self._pointer[0], memtype) 8939b2de37Sjeremylt 9039b2de37Sjeremylt return memtype[0] 9139b2de37Sjeremylt 9239b2de37Sjeremylt # CeedVector 9339b2de37Sjeremylt def Vector(self, size): 9439b2de37Sjeremylt """Ceed Vector: storing and manipulating vectors. 9539b2de37Sjeremylt 9639b2de37Sjeremylt Args: 9739b2de37Sjeremylt size: length of vector 9839b2de37Sjeremylt 9939b2de37Sjeremylt Returns: 10039b2de37Sjeremylt vector: Ceed Vector""" 10139b2de37Sjeremylt 10239b2de37Sjeremylt return Vector(self, size) 10339b2de37Sjeremylt 10439b2de37Sjeremylt # CeedElemRestriction 105*d979a051Sjeremylt def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets, 106*d979a051Sjeremylt memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 10739b2de37Sjeremylt """Ceed ElemRestriction: restriction from local vectors to elements. 10839b2de37Sjeremylt 10939b2de37Sjeremylt Args: 110*d979a051Sjeremylt nelem: number of elements described by the restriction 11139b2de37Sjeremylt elemsize: size (number of nodes) per element 112*d979a051Sjeremylt ncomp: number of field components per interpolation node 113*d979a051Sjeremylt (1 for scalar fields) 114*d979a051Sjeremylt compstride: Stride between components for the same L-vector "node". 115*d979a051Sjeremylt Data for node i, component k can be found in the 116*d979a051Sjeremylt L-vector at index [offsets[i] + k*compstride]. 117*d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 118*d979a051Sjeremylt the elements and fields given by this restriction. 119*d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 120*d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 12139b2de37Sjeremylt for the unknowns corresponding to element i, where 122*d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 123*d979a051Sjeremylt [0, lsize - 1]. 124*d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 125*d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 12639b2de37Sjeremylt 12739b2de37Sjeremylt Returns: 12839b2de37Sjeremylt elemrestriction: Ceed ElemRestiction""" 12939b2de37Sjeremylt 130*d979a051Sjeremylt return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 131*d979a051Sjeremylt offsets, memtype=memtype, cmode=cmode) 13239b2de37Sjeremylt 133*d979a051Sjeremylt def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides): 134*d979a051Sjeremylt """Ceed Identity ElemRestriction: strided restriction from local vectors 135*d979a051Sjeremylt to elements. 13639b2de37Sjeremylt 13739b2de37Sjeremylt Args: 138*d979a051Sjeremylt nelem: number of elements described by the restriction 13939b2de37Sjeremylt elemsize: size (number of nodes) per element 140a8d32208Sjeremylt ncomp: number of field components per interpolation node 141a8d32208Sjeremylt (1 for scalar fields) 142*d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 143*d979a051Sjeremylt the elements and fields given by this restriction. 14469a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 14569a53589Sjeremylt The data for node i, component j, element k in the 14669a53589Sjeremylt L-vector is given by 14769a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 14839b2de37Sjeremylt 14939b2de37Sjeremylt Returns: 15069a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 15139b2de37Sjeremylt 1527a7b0fa3SJed Brown return StridedElemRestriction( 153*d979a051Sjeremylt self, nelem, elemsize, ncomp, lsize, strides) 15439b2de37Sjeremylt 155*d979a051Sjeremylt def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 156*d979a051Sjeremylt lsize, offsets, memtype=lib.CEED_MEM_HOST, 157*d979a051Sjeremylt cmode=lib.CEED_COPY_VALUES): 158*d979a051Sjeremylt """Ceed Blocked ElemRestriction: blocked restriction from local vectors 159*d979a051Sjeremylt to elements. 16039b2de37Sjeremylt 16139b2de37Sjeremylt Args: 162*d979a051Sjeremylt nelem: number of elements described by the restriction 16339b2de37Sjeremylt elemsize: size (number of nodes) per element 16439b2de37Sjeremylt blksize: number of elements in a block 165*d979a051Sjeremylt ncomp: number of field components per interpolation node 166*d979a051Sjeremylt (1 for scalar fields) 167*d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 168*d979a051Sjeremylt the elements and fields given by this restriction. 169*d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 170*d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 17139b2de37Sjeremylt for the unknowns corresponding to element i, where 172*d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 173*d979a051Sjeremylt [0, lsize - 1]. The backend will permute and pad this 17439b2de37Sjeremylt array to the desired ordering for the blocksize, which is 17539b2de37Sjeremylt typically given by the backend. The default reordering is 17639b2de37Sjeremylt to interlace elements. 177*d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 178*d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 17939b2de37Sjeremylt 18039b2de37Sjeremylt Returns: 18139b2de37Sjeremylt elemrestriction: Ceed Blocked ElemRestiction""" 18239b2de37Sjeremylt 183*d979a051Sjeremylt return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, 184*d979a051Sjeremylt compstride, lsize, offsets, 185*d979a051Sjeremylt memtype=memtype, cmode=cmode) 18639b2de37Sjeremylt 187*d979a051Sjeremylt def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp, 188*d979a051Sjeremylt lsize, strides): 189*d979a051Sjeremylt """Ceed Blocked Strided ElemRestriction: blocked and strided restriction 190*d979a051Sjeremylt from local vectors to elements. 19169a53589Sjeremylt 19269a53589Sjeremylt Args: 193*d979a051Sjeremylt nelem: number of elements described in the restriction 19469a53589Sjeremylt elemsize: size (number of nodes) per element 19569a53589Sjeremylt blksize: number of elements in a block 19669a53589Sjeremylt ncomp: number of field components per interpolation node 19769a53589Sjeremylt (1 for scalar fields) 198*d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 199*d979a051Sjeremylt the elements and fields given by this restriction. 20069a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 20169a53589Sjeremylt The data for node i, component j, element k in the 20269a53589Sjeremylt L-vector is given by 20369a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 20469a53589Sjeremylt 20569a53589Sjeremylt Returns: 20669a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 20769a53589Sjeremylt 20869a53589Sjeremylt return BlockedStridedElemRestriction(self, nelem, elemsize, blksize, 209*d979a051Sjeremylt ncomp, lsize, strides) 21069a53589Sjeremylt 21139b2de37Sjeremylt # CeedBasis 21239b2de37Sjeremylt def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 21339b2de37Sjeremylt qref1d, qweight1d): 21439b2de37Sjeremylt """Ceed Tensor H1 Basis: finite element tensor-product basis objects for 21539b2de37Sjeremylt H^1 discretizations. 21639b2de37Sjeremylt 21739b2de37Sjeremylt Args: 21839b2de37Sjeremylt dim: topological dimension 21939b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 22039b2de37Sjeremylt P1d: number of nodes in one dimension 22139b2de37Sjeremylt Q1d: number of quadrature points in one dimension 222*d979a051Sjeremylt *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix 223*d979a051Sjeremylt expressing the values of nodal basis functions at 224*d979a051Sjeremylt quadrature points 225*d979a051Sjeremylt *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix 226*d979a051Sjeremylt expressing the derivatives of nodal basis functions at 227*d979a051Sjeremylt quadrature points 228*d979a051Sjeremylt *qref1d: Numpy array of length Q1d holding the locations of 229*d979a051Sjeremylt quadrature points on the 1D reference element [-1, 1] 230*d979a051Sjeremylt *qweight1d: Numpy array of length Q1d holding the quadrature 231*d979a051Sjeremylt weights on the reference element 23239b2de37Sjeremylt 23339b2de37Sjeremylt Returns: 23439b2de37Sjeremylt basis: Ceed Basis""" 23539b2de37Sjeremylt 23639b2de37Sjeremylt return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 23739b2de37Sjeremylt qref1d, qweight1d) 23839b2de37Sjeremylt 23939b2de37Sjeremylt def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode): 240*d979a051Sjeremylt """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange 241*d979a051Sjeremylt basis objects for H^1 discretizations. 24239b2de37Sjeremylt 24339b2de37Sjeremylt Args: 24439b2de37Sjeremylt dim: topological dimension 24539b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 24639b2de37Sjeremylt P: number of Gauss-Lobatto nodes in one dimension. The 24739b2de37Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 24839b2de37Sjeremylt Q: number of quadrature points in one dimension 24939b2de37Sjeremylt qmode: distribution of the Q quadrature points (affects order of 25039b2de37Sjeremylt accuracy for the quadrature) 25139b2de37Sjeremylt 25239b2de37Sjeremylt Returns: 25339b2de37Sjeremylt basis: Ceed Basis""" 25439b2de37Sjeremylt 25539b2de37Sjeremylt return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode) 25639b2de37Sjeremylt 25739b2de37Sjeremylt def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight): 258*d979a051Sjeremylt """Ceed H1 Basis: finite element non tensor-product basis for H^1 259*d979a051Sjeremylt discretizations. 26039b2de37Sjeremylt 26139b2de37Sjeremylt Args: 26239b2de37Sjeremylt topo: topology of the element, e.g. hypercube, simplex, etc 26339b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 26439b2de37Sjeremylt nnodes: total number of nodes 26539b2de37Sjeremylt nqpts: total number of quadrature points 266*d979a051Sjeremylt *interp: Numpy array holding the row-major (nqpts * nnodes) matrix 267*d979a051Sjeremylt expressing the values of nodal basis functions at 26839b2de37Sjeremylt quadrature points 269*d979a051Sjeremylt *grad: Numpy array holding the row-major (nqpts * dim * nnodes) 270*d979a051Sjeremylt matrix expressing the derivatives of nodal basis functions 271*d979a051Sjeremylt at quadrature points 272*d979a051Sjeremylt *qref: Numpy array of length (nqpts * dim) holding the locations of 273*d979a051Sjeremylt quadrature points on the reference element [-1, 1] 274*d979a051Sjeremylt *qweight: Numpy array of length nnodes holding the quadrature 275*d979a051Sjeremylt weights on the reference element 27639b2de37Sjeremylt 27739b2de37Sjeremylt Returns: 27839b2de37Sjeremylt basis: Ceed Basis""" 27939b2de37Sjeremylt 2807a7b0fa3SJed Brown return BasisH1(self, topo, ncomp, nnodes, nqpts, 2817a7b0fa3SJed Brown interp, grad, qref, qweight) 28239b2de37Sjeremylt 28339b2de37Sjeremylt # CeedQFunction 28439b2de37Sjeremylt def QFunction(self, vlength, f, source): 285*d979a051Sjeremylt """Ceed QFunction: point-wise operation at quadrature points for 286*d979a051Sjeremylt evaluating volumetric terms. 28739b2de37Sjeremylt 28839b2de37Sjeremylt Args: 28939b2de37Sjeremylt vlength: vector length. Caller must ensure that number of quadrature 29039b2de37Sjeremylt points is a multiple of vlength 29139b2de37Sjeremylt f: ctypes function pointer to evaluate action at quadrature points 29239b2de37Sjeremylt source: absolute path to source of QFunction, 29339b2de37Sjeremylt "\\abs_path\\file.h:function_name 29439b2de37Sjeremylt 29539b2de37Sjeremylt Returns: 29639b2de37Sjeremylt qfunction: Ceed QFunction""" 29739b2de37Sjeremylt 29839b2de37Sjeremylt return QFunction(self, vlength, f, source) 29939b2de37Sjeremylt 30039b2de37Sjeremylt def QFunctionByName(self, name): 30139b2de37Sjeremylt """Ceed QFunction By Name: point-wise operation at quadrature points 30239b2de37Sjeremylt from a given gallery, for evaluating volumetric terms. 30339b2de37Sjeremylt 30439b2de37Sjeremylt Args: 30539b2de37Sjeremylt name: name of QFunction to use from gallery 30639b2de37Sjeremylt 30739b2de37Sjeremylt Returns: 30839b2de37Sjeremylt qfunction: Ceed QFunction By Name""" 30939b2de37Sjeremylt 31039b2de37Sjeremylt return QFunctionByName(self, name) 31139b2de37Sjeremylt 31239b2de37Sjeremylt def IdentityQFunction(self, size, inmode, outmode): 31339b2de37Sjeremylt """Ceed Idenity QFunction: identity qfunction operation. 31439b2de37Sjeremylt 31539b2de37Sjeremylt Args: 31639b2de37Sjeremylt size: size of the qfunction fields 31739b2de37Sjeremylt **inmode: CeedEvalMode for input to Ceed QFunction 31839b2de37Sjeremylt **outmode: CeedEvalMode for output to Ceed QFunction 31939b2de37Sjeremylt 32039b2de37Sjeremylt Returns: 32139b2de37Sjeremylt qfunction: Ceed Identity QFunction""" 32239b2de37Sjeremylt 32339b2de37Sjeremylt return IdentityQFunction(self, size, inmode, outmode) 32439b2de37Sjeremylt 32539b2de37Sjeremylt # CeedOperator 32639b2de37Sjeremylt def Operator(self, qf, dqf=None, qdfT=None): 32739b2de37Sjeremylt """Ceed Operator: composed FE-type operations on vectors. 32839b2de37Sjeremylt 32939b2de37Sjeremylt Args: 330*d979a051Sjeremylt qf: QFunction defining the action of the operator at quadrature 331*d979a051Sjeremylt points 33239b2de37Sjeremylt **dqf: QFunction defining the action of the Jacobian, default None 333*d979a051Sjeremylt **dqfT: QFunction defining the action of the transpose of the 334*d979a051Sjeremylt Jacobian, default None 33539b2de37Sjeremylt 33639b2de37Sjeremylt Returns: 33739b2de37Sjeremylt operator: Ceed Operator""" 33839b2de37Sjeremylt 33939b2de37Sjeremylt return Operator(self, qf, dqf, qdfT) 34039b2de37Sjeremylt 34139b2de37Sjeremylt def CompositeOperator(self): 34239b2de37Sjeremylt """Composite Ceed Operator: composition of multiple CeedOperators. 34339b2de37Sjeremylt 34439b2de37Sjeremylt Returns: 34539b2de37Sjeremylt operator: Ceed Composite Operator""" 34639b2de37Sjeremylt 34739b2de37Sjeremylt return CompositeOperator(self) 34839b2de37Sjeremylt 34939b2de37Sjeremylt # Destructor 35039b2de37Sjeremylt def __del__(self): 35139b2de37Sjeremylt # libCEED call 35239b2de37Sjeremylt lib.CeedDestroy(self._pointer) 35339b2de37Sjeremylt 35439b2de37Sjeremylt# ------------------------------------------------------------------------------ 355