139b2de37Sjeremylt# Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 239b2de37Sjeremylt# the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 339b2de37Sjeremylt# reserved. See files LICENSE and NOTICE for details. 439b2de37Sjeremylt# 539b2de37Sjeremylt# This file is part of CEED, a collection of benchmarks, miniapps, software 639b2de37Sjeremylt# libraries and APIs for efficient high-order finite element and spectral 739b2de37Sjeremylt# element discretizations for exascale applications. For more information and 839b2de37Sjeremylt# source code availability see http://github.com/ceed. 939b2de37Sjeremylt# 1039b2de37Sjeremylt# The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 1139b2de37Sjeremylt# a collaborative effort of two U.S. Department of Energy organizations (Office 1239b2de37Sjeremylt# of Science and the National Nuclear Security Administration) responsible for 1339b2de37Sjeremylt# the planning and preparation of a capable exascale ecosystem, including 1439b2de37Sjeremylt# software, applications, hardware, advanced system engineering and early 1539b2de37Sjeremylt# testbed platforms, in support of the nation's exascale computing imperative. 1639b2de37Sjeremylt 1739b2de37Sjeremyltfrom _ceed_cffi import ffi, lib 1839b2de37Sjeremyltimport sys 19477729cfSJeremy L Thompsonimport os 2039b2de37Sjeremyltimport io 210a0da059Sjeremyltimport tempfile 2239b2de37Sjeremyltfrom abc import ABC 2339b2de37Sjeremyltfrom .ceed_vector import Vector 2439b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1 2569a53589Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, StridedElemRestriction, BlockedElemRestriction, BlockedStridedElemRestriction 2639b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction 27777ff853SJeremy L Thompsonfrom .ceed_qfunctioncontext import QFunctionContext 2839b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator 2939b2de37Sjeremyltfrom .ceed_constants import * 3039b2de37Sjeremylt 3139b2de37Sjeremylt# ------------------------------------------------------------------------------ 327a7b0fa3SJed Brown 337a7b0fa3SJed Brown 3439b2de37Sjeremyltclass Ceed(): 3539b2de37Sjeremylt """Ceed: core components.""" 3639b2de37Sjeremylt # Attributes 3739b2de37Sjeremylt _pointer = ffi.NULL 3839b2de37Sjeremylt 3939b2de37Sjeremylt # Constructor 40*74b533adSJed Brown def __init__(self, resource="/cpu/self", on_error="store"): 4139b2de37Sjeremylt # libCEED object 4239b2de37Sjeremylt self._pointer = ffi.new("Ceed *") 4339b2de37Sjeremylt 4439b2de37Sjeremylt # libCEED call 4539b2de37Sjeremylt resourceAscii = ffi.new("char[]", resource.encode("ascii")) 46477729cfSJeremy L Thompson os.environ["CEED_ERROR_HANDLER"] = "return" 47477729cfSJeremy L Thompson err_code = lib.CeedInit(resourceAscii, self._pointer) 48477729cfSJeremy L Thompson if err_code: 49477729cfSJeremy L Thompson raise Exception("Error initializing backend resource: " + resource) 50*74b533adSJed Brown error_handlers = dict( 51*74b533adSJed Brown store="CeedErrorStore", 52*74b533adSJed Brown abort="CeedErrorAbort", 53*74b533adSJed Brown ) 54477729cfSJeremy L Thompson lib.CeedSetErrorHandler( 55477729cfSJeremy L Thompson self._pointer[0], ffi.addressof( 56*74b533adSJed Brown lib, error_handlers[on_error])) 5739b2de37Sjeremylt 5839b2de37Sjeremylt # Representation 5939b2de37Sjeremylt def __repr__(self): 6039b2de37Sjeremylt return "<Ceed instance at " + hex(id(self)) + ">" 6139b2de37Sjeremylt 620a0da059Sjeremylt # String conversion for print() to stdout 630a0da059Sjeremylt def __str__(self): 640a0da059Sjeremylt """View a Ceed via print().""" 650a0da059Sjeremylt 660a0da059Sjeremylt # libCEED call 670a0da059Sjeremylt with tempfile.NamedTemporaryFile() as key_file: 680a0da059Sjeremylt with open(key_file.name, 'r+') as stream_file: 690a0da059Sjeremylt stream = ffi.cast("FILE *", stream_file) 700a0da059Sjeremylt 71477729cfSJeremy L Thompson err_code = lib.CeedView(self._pointer[0], stream) 72477729cfSJeremy L Thompson self._check_error(err_code) 730a0da059Sjeremylt 740a0da059Sjeremylt stream_file.seek(0) 750a0da059Sjeremylt out_string = stream_file.read() 760a0da059Sjeremylt 770a0da059Sjeremylt return out_string 780a0da059Sjeremylt 79477729cfSJeremy L Thompson # Error handler 80477729cfSJeremy L Thompson def _check_error(self, err_code): 81477729cfSJeremy L Thompson """Check return code and retrieve error message for non-zero code""" 82477729cfSJeremy L Thompson if (err_code): 83477729cfSJeremy L Thompson message = ffi.new("char **") 84477729cfSJeremy L Thompson lib.CeedGetErrorMessage(self._pointer[0], message) 85477729cfSJeremy L Thompson raise Exception(ffi.string(message[0]).decode("UTF-8")) 86477729cfSJeremy L Thompson 8739b2de37Sjeremylt # Get Resource 8839b2de37Sjeremylt def get_resource(self): 8939b2de37Sjeremylt """Get the full resource name for a Ceed context. 9039b2de37Sjeremylt 9139b2de37Sjeremylt Returns: 9239b2de37Sjeremylt resource: resource name""" 9339b2de37Sjeremylt 9439b2de37Sjeremylt # libCEED call 9539b2de37Sjeremylt resource = ffi.new("char **") 96477729cfSJeremy L Thompson err_code = lib.CeedGetResource(self._pointer[0], resource) 97477729cfSJeremy L Thompson self._check_error(err_code) 9839b2de37Sjeremylt 9939b2de37Sjeremylt return ffi.string(resource[0]).decode("UTF-8") 10039b2de37Sjeremylt 10139b2de37Sjeremylt # Preferred MemType 10239b2de37Sjeremylt def get_preferred_memtype(self): 10339b2de37Sjeremylt """Return Ceed preferred memory type. 10439b2de37Sjeremylt 10539b2de37Sjeremylt Returns: 10639b2de37Sjeremylt memtype: Ceed preferred memory type""" 10739b2de37Sjeremylt 10839b2de37Sjeremylt # libCEED call 10939b2de37Sjeremylt memtype = ffi.new("CeedMemType *", MEM_HOST) 110477729cfSJeremy L Thompson err_code = lib.CeedGetPreferredMemType(self._pointer[0], memtype) 111477729cfSJeremy L Thompson self._check_error(err_code) 11239b2de37Sjeremylt 11339b2de37Sjeremylt return memtype[0] 11439b2de37Sjeremylt 11539b2de37Sjeremylt # CeedVector 11639b2de37Sjeremylt def Vector(self, size): 11739b2de37Sjeremylt """Ceed Vector: storing and manipulating vectors. 11839b2de37Sjeremylt 11939b2de37Sjeremylt Args: 12039b2de37Sjeremylt size: length of vector 12139b2de37Sjeremylt 12239b2de37Sjeremylt Returns: 12339b2de37Sjeremylt vector: Ceed Vector""" 12439b2de37Sjeremylt 12539b2de37Sjeremylt return Vector(self, size) 12639b2de37Sjeremylt 12739b2de37Sjeremylt # CeedElemRestriction 128d979a051Sjeremylt def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets, 129d979a051Sjeremylt memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 13039b2de37Sjeremylt """Ceed ElemRestriction: restriction from local vectors to elements. 13139b2de37Sjeremylt 13239b2de37Sjeremylt Args: 133d979a051Sjeremylt nelem: number of elements described by the restriction 13439b2de37Sjeremylt elemsize: size (number of nodes) per element 135d979a051Sjeremylt ncomp: number of field components per interpolation node 136d979a051Sjeremylt (1 for scalar fields) 137d979a051Sjeremylt compstride: Stride between components for the same L-vector "node". 138d979a051Sjeremylt Data for node i, component k can be found in the 139d979a051Sjeremylt L-vector at index [offsets[i] + k*compstride]. 140d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 141d979a051Sjeremylt the elements and fields given by this restriction. 142d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 143d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 14439b2de37Sjeremylt for the unknowns corresponding to element i, where 145d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 146d979a051Sjeremylt [0, lsize - 1]. 147d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 148d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 14939b2de37Sjeremylt 15039b2de37Sjeremylt Returns: 15139b2de37Sjeremylt elemrestriction: Ceed ElemRestiction""" 15239b2de37Sjeremylt 153d979a051Sjeremylt return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 154d979a051Sjeremylt offsets, memtype=memtype, cmode=cmode) 15539b2de37Sjeremylt 156d979a051Sjeremylt def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides): 157d979a051Sjeremylt """Ceed Identity ElemRestriction: strided restriction from local vectors 158d979a051Sjeremylt to elements. 15939b2de37Sjeremylt 16039b2de37Sjeremylt Args: 161d979a051Sjeremylt nelem: number of elements described by the restriction 16239b2de37Sjeremylt elemsize: size (number of nodes) per element 163a8d32208Sjeremylt ncomp: number of field components per interpolation node 164a8d32208Sjeremylt (1 for scalar fields) 165d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 166d979a051Sjeremylt the elements and fields given by this restriction. 16769a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 16869a53589Sjeremylt The data for node i, component j, element k in the 16969a53589Sjeremylt L-vector is given by 17069a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 17139b2de37Sjeremylt 17239b2de37Sjeremylt Returns: 17369a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 17439b2de37Sjeremylt 1757a7b0fa3SJed Brown return StridedElemRestriction( 176d979a051Sjeremylt self, nelem, elemsize, ncomp, lsize, strides) 17739b2de37Sjeremylt 178d979a051Sjeremylt def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 179d979a051Sjeremylt lsize, offsets, memtype=lib.CEED_MEM_HOST, 180d979a051Sjeremylt cmode=lib.CEED_COPY_VALUES): 181d979a051Sjeremylt """Ceed Blocked ElemRestriction: blocked restriction from local vectors 182d979a051Sjeremylt to elements. 18339b2de37Sjeremylt 18439b2de37Sjeremylt Args: 185d979a051Sjeremylt nelem: number of elements described by the restriction 18639b2de37Sjeremylt elemsize: size (number of nodes) per element 18739b2de37Sjeremylt blksize: number of elements in a block 188d979a051Sjeremylt ncomp: number of field components per interpolation node 189d979a051Sjeremylt (1 for scalar fields) 190d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 191d979a051Sjeremylt the elements and fields given by this restriction. 192d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 193d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 19439b2de37Sjeremylt for the unknowns corresponding to element i, where 195d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 196d979a051Sjeremylt [0, lsize - 1]. The backend will permute and pad this 19739b2de37Sjeremylt array to the desired ordering for the blocksize, which is 19839b2de37Sjeremylt typically given by the backend. The default reordering is 19939b2de37Sjeremylt to interlace elements. 200d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 201d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 20239b2de37Sjeremylt 20339b2de37Sjeremylt Returns: 20439b2de37Sjeremylt elemrestriction: Ceed Blocked ElemRestiction""" 20539b2de37Sjeremylt 206d979a051Sjeremylt return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, 207d979a051Sjeremylt compstride, lsize, offsets, 208d979a051Sjeremylt memtype=memtype, cmode=cmode) 20939b2de37Sjeremylt 210d979a051Sjeremylt def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp, 211d979a051Sjeremylt lsize, strides): 212d979a051Sjeremylt """Ceed Blocked Strided ElemRestriction: blocked and strided restriction 213d979a051Sjeremylt from local vectors to elements. 21469a53589Sjeremylt 21569a53589Sjeremylt Args: 216d979a051Sjeremylt nelem: number of elements described in the restriction 21769a53589Sjeremylt elemsize: size (number of nodes) per element 21869a53589Sjeremylt blksize: number of elements in a block 21969a53589Sjeremylt ncomp: number of field components per interpolation node 22069a53589Sjeremylt (1 for scalar fields) 221d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 222d979a051Sjeremylt the elements and fields given by this restriction. 22369a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 22469a53589Sjeremylt The data for node i, component j, element k in the 22569a53589Sjeremylt L-vector is given by 22669a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 22769a53589Sjeremylt 22869a53589Sjeremylt Returns: 22969a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 23069a53589Sjeremylt 23169a53589Sjeremylt return BlockedStridedElemRestriction(self, nelem, elemsize, blksize, 232d979a051Sjeremylt ncomp, lsize, strides) 23369a53589Sjeremylt 23439b2de37Sjeremylt # CeedBasis 23539b2de37Sjeremylt def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 23639b2de37Sjeremylt qref1d, qweight1d): 23739b2de37Sjeremylt """Ceed Tensor H1 Basis: finite element tensor-product basis objects for 23839b2de37Sjeremylt H^1 discretizations. 23939b2de37Sjeremylt 24039b2de37Sjeremylt Args: 24139b2de37Sjeremylt dim: topological dimension 24239b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 24339b2de37Sjeremylt P1d: number of nodes in one dimension 24439b2de37Sjeremylt Q1d: number of quadrature points in one dimension 245d979a051Sjeremylt *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix 246d979a051Sjeremylt expressing the values of nodal basis functions at 247d979a051Sjeremylt quadrature points 248d979a051Sjeremylt *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix 249d979a051Sjeremylt expressing the derivatives of nodal basis functions at 250d979a051Sjeremylt quadrature points 251d979a051Sjeremylt *qref1d: Numpy array of length Q1d holding the locations of 252d979a051Sjeremylt quadrature points on the 1D reference element [-1, 1] 253d979a051Sjeremylt *qweight1d: Numpy array of length Q1d holding the quadrature 254d979a051Sjeremylt weights on the reference element 25539b2de37Sjeremylt 25639b2de37Sjeremylt Returns: 25739b2de37Sjeremylt basis: Ceed Basis""" 25839b2de37Sjeremylt 25939b2de37Sjeremylt return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 26039b2de37Sjeremylt qref1d, qweight1d) 26139b2de37Sjeremylt 26239b2de37Sjeremylt def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode): 263d979a051Sjeremylt """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange 264d979a051Sjeremylt basis objects for H^1 discretizations. 26539b2de37Sjeremylt 26639b2de37Sjeremylt Args: 26739b2de37Sjeremylt dim: topological dimension 26839b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 26939b2de37Sjeremylt P: number of Gauss-Lobatto nodes in one dimension. The 27039b2de37Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 27139b2de37Sjeremylt Q: number of quadrature points in one dimension 27239b2de37Sjeremylt qmode: distribution of the Q quadrature points (affects order of 27339b2de37Sjeremylt accuracy for the quadrature) 27439b2de37Sjeremylt 27539b2de37Sjeremylt Returns: 27639b2de37Sjeremylt basis: Ceed Basis""" 27739b2de37Sjeremylt 27839b2de37Sjeremylt return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode) 27939b2de37Sjeremylt 28039b2de37Sjeremylt def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight): 281d979a051Sjeremylt """Ceed H1 Basis: finite element non tensor-product basis for H^1 282d979a051Sjeremylt discretizations. 28339b2de37Sjeremylt 28439b2de37Sjeremylt Args: 28539b2de37Sjeremylt topo: topology of the element, e.g. hypercube, simplex, etc 28639b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 28739b2de37Sjeremylt nnodes: total number of nodes 28839b2de37Sjeremylt nqpts: total number of quadrature points 289d979a051Sjeremylt *interp: Numpy array holding the row-major (nqpts * nnodes) matrix 290d979a051Sjeremylt expressing the values of nodal basis functions at 29139b2de37Sjeremylt quadrature points 292d979a051Sjeremylt *grad: Numpy array holding the row-major (nqpts * dim * nnodes) 293d979a051Sjeremylt matrix expressing the derivatives of nodal basis functions 294d979a051Sjeremylt at quadrature points 295d979a051Sjeremylt *qref: Numpy array of length (nqpts * dim) holding the locations of 296d979a051Sjeremylt quadrature points on the reference element [-1, 1] 297d979a051Sjeremylt *qweight: Numpy array of length nnodes holding the quadrature 298d979a051Sjeremylt weights on the reference element 29939b2de37Sjeremylt 30039b2de37Sjeremylt Returns: 30139b2de37Sjeremylt basis: Ceed Basis""" 30239b2de37Sjeremylt 3037a7b0fa3SJed Brown return BasisH1(self, topo, ncomp, nnodes, nqpts, 3047a7b0fa3SJed Brown interp, grad, qref, qweight) 30539b2de37Sjeremylt 30639b2de37Sjeremylt # CeedQFunction 30739b2de37Sjeremylt def QFunction(self, vlength, f, source): 308d979a051Sjeremylt """Ceed QFunction: point-wise operation at quadrature points for 309d979a051Sjeremylt evaluating volumetric terms. 31039b2de37Sjeremylt 31139b2de37Sjeremylt Args: 31239b2de37Sjeremylt vlength: vector length. Caller must ensure that number of quadrature 31339b2de37Sjeremylt points is a multiple of vlength 31439b2de37Sjeremylt f: ctypes function pointer to evaluate action at quadrature points 31539b2de37Sjeremylt source: absolute path to source of QFunction, 31639b2de37Sjeremylt "\\abs_path\\file.h:function_name 31739b2de37Sjeremylt 31839b2de37Sjeremylt Returns: 31939b2de37Sjeremylt qfunction: Ceed QFunction""" 32039b2de37Sjeremylt 32139b2de37Sjeremylt return QFunction(self, vlength, f, source) 32239b2de37Sjeremylt 32339b2de37Sjeremylt def QFunctionByName(self, name): 32439b2de37Sjeremylt """Ceed QFunction By Name: point-wise operation at quadrature points 32539b2de37Sjeremylt from a given gallery, for evaluating volumetric terms. 32639b2de37Sjeremylt 32739b2de37Sjeremylt Args: 32839b2de37Sjeremylt name: name of QFunction to use from gallery 32939b2de37Sjeremylt 33039b2de37Sjeremylt Returns: 33139b2de37Sjeremylt qfunction: Ceed QFunction By Name""" 33239b2de37Sjeremylt 33339b2de37Sjeremylt return QFunctionByName(self, name) 33439b2de37Sjeremylt 33539b2de37Sjeremylt def IdentityQFunction(self, size, inmode, outmode): 33639b2de37Sjeremylt """Ceed Idenity QFunction: identity qfunction operation. 33739b2de37Sjeremylt 33839b2de37Sjeremylt Args: 33939b2de37Sjeremylt size: size of the qfunction fields 34039b2de37Sjeremylt **inmode: CeedEvalMode for input to Ceed QFunction 34139b2de37Sjeremylt **outmode: CeedEvalMode for output to Ceed QFunction 34239b2de37Sjeremylt 34339b2de37Sjeremylt Returns: 34439b2de37Sjeremylt qfunction: Ceed Identity QFunction""" 34539b2de37Sjeremylt 34639b2de37Sjeremylt return IdentityQFunction(self, size, inmode, outmode) 34739b2de37Sjeremylt 348777ff853SJeremy L Thompson def QFunctionContext(self): 349777ff853SJeremy L Thompson """Ceed QFunction Context: stores Ceed QFunction user context data. 350777ff853SJeremy L Thompson 351777ff853SJeremy L Thompson Returns: 352777ff853SJeremy L Thompson userContext: Ceed QFunction Context""" 353777ff853SJeremy L Thompson 354777ff853SJeremy L Thompson return QFunctionContext(self) 355777ff853SJeremy L Thompson 35639b2de37Sjeremylt # CeedOperator 35739b2de37Sjeremylt def Operator(self, qf, dqf=None, qdfT=None): 35839b2de37Sjeremylt """Ceed Operator: composed FE-type operations on vectors. 35939b2de37Sjeremylt 36039b2de37Sjeremylt Args: 361d979a051Sjeremylt qf: QFunction defining the action of the operator at quadrature 362d979a051Sjeremylt points 36339b2de37Sjeremylt **dqf: QFunction defining the action of the Jacobian, default None 364d979a051Sjeremylt **dqfT: QFunction defining the action of the transpose of the 365d979a051Sjeremylt Jacobian, default None 36639b2de37Sjeremylt 36739b2de37Sjeremylt Returns: 36839b2de37Sjeremylt operator: Ceed Operator""" 36939b2de37Sjeremylt 37039b2de37Sjeremylt return Operator(self, qf, dqf, qdfT) 37139b2de37Sjeremylt 37239b2de37Sjeremylt def CompositeOperator(self): 37339b2de37Sjeremylt """Composite Ceed Operator: composition of multiple CeedOperators. 37439b2de37Sjeremylt 37539b2de37Sjeremylt Returns: 37639b2de37Sjeremylt operator: Ceed Composite Operator""" 37739b2de37Sjeremylt 37839b2de37Sjeremylt return CompositeOperator(self) 37939b2de37Sjeremylt 38039b2de37Sjeremylt # Destructor 38139b2de37Sjeremylt def __del__(self): 38239b2de37Sjeremylt # libCEED call 38339b2de37Sjeremylt lib.CeedDestroy(self._pointer) 38439b2de37Sjeremylt 38539b2de37Sjeremylt# ------------------------------------------------------------------------------ 386