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 16*97c1c57aSSebastian Grimbergfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1, BasisHdiv, BasisHcurl 1769a53589Sjeremyltfrom .ceed_elemrestriction import ElemRestriction, StridedElemRestriction, BlockedElemRestriction, BlockedStridedElemRestriction 1839b2de37Sjeremyltfrom .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction 19777ff853SJeremy L Thompsonfrom .ceed_qfunctioncontext import QFunctionContext 2039b2de37Sjeremyltfrom .ceed_operator import Operator, CompositeOperator 2139b2de37Sjeremyltfrom .ceed_constants import * 2239b2de37Sjeremylt 2339b2de37Sjeremylt# ------------------------------------------------------------------------------ 247a7b0fa3SJed Brown 257a7b0fa3SJed Brown 2639b2de37Sjeremyltclass Ceed(): 2739b2de37Sjeremylt """Ceed: core components.""" 2839b2de37Sjeremylt 2939b2de37Sjeremylt # Constructor 3074b533adSJed Brown def __init__(self, resource="/cpu/self", on_error="store"): 3139b2de37Sjeremylt # libCEED object 3239b2de37Sjeremylt self._pointer = ffi.new("Ceed *") 3339b2de37Sjeremylt 3439b2de37Sjeremylt # libCEED call 3539b2de37Sjeremylt resourceAscii = ffi.new("char[]", resource.encode("ascii")) 36477729cfSJeremy L Thompson os.environ["CEED_ERROR_HANDLER"] = "return" 37477729cfSJeremy L Thompson err_code = lib.CeedInit(resourceAscii, self._pointer) 38477729cfSJeremy L Thompson if err_code: 39477729cfSJeremy L Thompson raise Exception("Error initializing backend resource: " + resource) 4074b533adSJed Brown error_handlers = dict( 4174b533adSJed Brown store="CeedErrorStore", 4274b533adSJed Brown abort="CeedErrorAbort", 4374b533adSJed Brown ) 44477729cfSJeremy L Thompson lib.CeedSetErrorHandler( 45477729cfSJeremy L Thompson self._pointer[0], ffi.addressof( 4674b533adSJed Brown lib, error_handlers[on_error])) 4739b2de37Sjeremylt 4839b2de37Sjeremylt # Representation 4939b2de37Sjeremylt def __repr__(self): 5039b2de37Sjeremylt return "<Ceed instance at " + hex(id(self)) + ">" 5139b2de37Sjeremylt 520a0da059Sjeremylt # String conversion for print() to stdout 530a0da059Sjeremylt def __str__(self): 540a0da059Sjeremylt """View a Ceed via print().""" 550a0da059Sjeremylt 560a0da059Sjeremylt # libCEED call 570a0da059Sjeremylt with tempfile.NamedTemporaryFile() as key_file: 580a0da059Sjeremylt with open(key_file.name, 'r+') as stream_file: 590a0da059Sjeremylt stream = ffi.cast("FILE *", stream_file) 600a0da059Sjeremylt 61477729cfSJeremy L Thompson err_code = lib.CeedView(self._pointer[0], stream) 62477729cfSJeremy L Thompson self._check_error(err_code) 630a0da059Sjeremylt 640a0da059Sjeremylt stream_file.seek(0) 650a0da059Sjeremylt out_string = stream_file.read() 660a0da059Sjeremylt 670a0da059Sjeremylt return out_string 680a0da059Sjeremylt 69477729cfSJeremy L Thompson # Error handler 70477729cfSJeremy L Thompson def _check_error(self, err_code): 71477729cfSJeremy L Thompson """Check return code and retrieve error message for non-zero code""" 72e15f9bd0SJeremy L Thompson if (err_code != lib.CEED_ERROR_SUCCESS): 73477729cfSJeremy L Thompson message = ffi.new("char **") 74477729cfSJeremy L Thompson lib.CeedGetErrorMessage(self._pointer[0], message) 75477729cfSJeremy L Thompson raise Exception(ffi.string(message[0]).decode("UTF-8")) 76477729cfSJeremy L Thompson 7739b2de37Sjeremylt # Get Resource 7839b2de37Sjeremylt def get_resource(self): 7939b2de37Sjeremylt """Get the full resource name for a Ceed context. 8039b2de37Sjeremylt 8139b2de37Sjeremylt Returns: 8239b2de37Sjeremylt resource: resource name""" 8339b2de37Sjeremylt 8439b2de37Sjeremylt # libCEED call 8539b2de37Sjeremylt resource = ffi.new("char **") 86477729cfSJeremy L Thompson err_code = lib.CeedGetResource(self._pointer[0], resource) 87477729cfSJeremy L Thompson self._check_error(err_code) 8839b2de37Sjeremylt 8939b2de37Sjeremylt return ffi.string(resource[0]).decode("UTF-8") 9039b2de37Sjeremylt 9139b2de37Sjeremylt # Preferred MemType 9239b2de37Sjeremylt def get_preferred_memtype(self): 9339b2de37Sjeremylt """Return Ceed preferred memory type. 9439b2de37Sjeremylt 9539b2de37Sjeremylt Returns: 9639b2de37Sjeremylt memtype: Ceed preferred memory type""" 9739b2de37Sjeremylt 9839b2de37Sjeremylt # libCEED call 9939b2de37Sjeremylt memtype = ffi.new("CeedMemType *", MEM_HOST) 100477729cfSJeremy L Thompson err_code = lib.CeedGetPreferredMemType(self._pointer[0], memtype) 101477729cfSJeremy L Thompson self._check_error(err_code) 10239b2de37Sjeremylt 10339b2de37Sjeremylt return memtype[0] 10439b2de37Sjeremylt 10580a9ef05SNatalie Beams # Convenience function to get CeedScalar type 10680a9ef05SNatalie Beams def scalar_type(self): 10780a9ef05SNatalie Beams """Return dtype string for CeedScalar. 10880a9ef05SNatalie Beams 10980a9ef05SNatalie Beams Returns: 11080a9ef05SNatalie Beams np_dtype: String naming numpy data type corresponding to CeedScalar""" 11180a9ef05SNatalie Beams 11280a9ef05SNatalie Beams return scalar_types[lib.CEED_SCALAR_TYPE] 11380a9ef05SNatalie Beams 114e15f9bd0SJeremy L Thompson # --- Basis utility functions --- 115e15f9bd0SJeremy L Thompson 116e15f9bd0SJeremy L Thompson # Gauss quadrature 117e15f9bd0SJeremy L Thompson def gauss_quadrature(self, q): 118e15f9bd0SJeremy L Thompson """Construct a Gauss-Legendre quadrature. 119e15f9bd0SJeremy L Thompson 120e15f9bd0SJeremy L Thompson Args: 121e15f9bd0SJeremy L Thompson Q: number of quadrature points (integrates polynomials of 122e15f9bd0SJeremy L Thompson degree 2*Q-1 exactly) 123e15f9bd0SJeremy L Thompson 124e15f9bd0SJeremy L Thompson Returns: 125e15f9bd0SJeremy L Thompson (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1] 126e15f9bd0SJeremy L Thompson and array of length Q to hold the weights""" 127e15f9bd0SJeremy L Thompson 128e15f9bd0SJeremy L Thompson # Setup arguments 12980a9ef05SNatalie Beams qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 13080a9ef05SNatalie Beams qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 131e15f9bd0SJeremy L Thompson 132e15f9bd0SJeremy L Thompson qref1d_pointer = ffi.new("CeedScalar *") 133e15f9bd0SJeremy L Thompson qref1d_pointer = ffi.cast( 134e15f9bd0SJeremy L Thompson "CeedScalar *", 135e15f9bd0SJeremy L Thompson qref1d.__array_interface__['data'][0]) 136e15f9bd0SJeremy L Thompson 137e15f9bd0SJeremy L Thompson qweight1d_pointer = ffi.new("CeedScalar *") 138e15f9bd0SJeremy L Thompson qweight1d_pointer = ffi.cast( 139e15f9bd0SJeremy L Thompson "CeedScalar *", 140e15f9bd0SJeremy L Thompson qweight1d.__array_interface__['data'][0]) 141e15f9bd0SJeremy L Thompson 142e15f9bd0SJeremy L Thompson # libCEED call 143e15f9bd0SJeremy L Thompson err_code = lib.CeedGaussQuadrature(q, qref1d_pointer, qweight1d_pointer) 144e15f9bd0SJeremy L Thompson self._check_error(err_code) 145e15f9bd0SJeremy L Thompson 146e15f9bd0SJeremy L Thompson return qref1d, qweight1d 147e15f9bd0SJeremy L Thompson 148e15f9bd0SJeremy L Thompson # Lobatto quadrature 149e15f9bd0SJeremy L Thompson def lobatto_quadrature(self, q): 150e15f9bd0SJeremy L Thompson """Construct a Gauss-Legendre-Lobatto quadrature. 151e15f9bd0SJeremy L Thompson 152e15f9bd0SJeremy L Thompson Args: 153e15f9bd0SJeremy L Thompson q: number of quadrature points (integrates polynomials of 154e15f9bd0SJeremy L Thompson degree 2*Q-3 exactly) 155e15f9bd0SJeremy L Thompson 156e15f9bd0SJeremy L Thompson Returns: 157e15f9bd0SJeremy L Thompson (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1] 158e15f9bd0SJeremy L Thompson and array of length Q to hold the weights""" 159e15f9bd0SJeremy L Thompson 160e15f9bd0SJeremy L Thompson # Setup arguments 16180a9ef05SNatalie Beams qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 162e15f9bd0SJeremy L Thompson qref1d_pointer = ffi.new("CeedScalar *") 163e15f9bd0SJeremy L Thompson qref1d_pointer = ffi.cast( 164e15f9bd0SJeremy L Thompson "CeedScalar *", 165e15f9bd0SJeremy L Thompson qref1d.__array_interface__['data'][0]) 166e15f9bd0SJeremy L Thompson 16780a9ef05SNatalie Beams qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 168e15f9bd0SJeremy L Thompson qweight1d_pointer = ffi.new("CeedScalar *") 169e15f9bd0SJeremy L Thompson qweight1d_pointer = ffi.cast( 170e15f9bd0SJeremy L Thompson "CeedScalar *", 171e15f9bd0SJeremy L Thompson qweight1d.__array_interface__['data'][0]) 172e15f9bd0SJeremy L Thompson 173e15f9bd0SJeremy L Thompson # libCEED call 174e15f9bd0SJeremy L Thompson err_code = lib.CeedLobattoQuadrature( 175e15f9bd0SJeremy L Thompson q, qref1d_pointer, qweight1d_pointer) 176e15f9bd0SJeremy L Thompson self._check_error(err_code) 177e15f9bd0SJeremy L Thompson 178e15f9bd0SJeremy L Thompson return qref1d, qweight1d 179e15f9bd0SJeremy L Thompson 180e15f9bd0SJeremy L Thompson # --- libCEED objects --- 181e15f9bd0SJeremy L Thompson 18239b2de37Sjeremylt # CeedVector 18339b2de37Sjeremylt def Vector(self, size): 18439b2de37Sjeremylt """Ceed Vector: storing and manipulating vectors. 18539b2de37Sjeremylt 18639b2de37Sjeremylt Args: 18739b2de37Sjeremylt size: length of vector 18839b2de37Sjeremylt 18939b2de37Sjeremylt Returns: 19039b2de37Sjeremylt vector: Ceed Vector""" 19139b2de37Sjeremylt 19239b2de37Sjeremylt return Vector(self, size) 19339b2de37Sjeremylt 19439b2de37Sjeremylt # CeedElemRestriction 195d979a051Sjeremylt def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets, 196d979a051Sjeremylt memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 19739b2de37Sjeremylt """Ceed ElemRestriction: restriction from local vectors to elements. 19839b2de37Sjeremylt 19939b2de37Sjeremylt Args: 200d979a051Sjeremylt nelem: number of elements described by the restriction 20139b2de37Sjeremylt elemsize: size (number of nodes) per element 202d979a051Sjeremylt ncomp: number of field components per interpolation node 203d979a051Sjeremylt (1 for scalar fields) 204d979a051Sjeremylt compstride: Stride between components for the same L-vector "node". 205d979a051Sjeremylt Data for node i, component k can be found in the 206d979a051Sjeremylt L-vector at index [offsets[i] + k*compstride]. 207d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 208d979a051Sjeremylt the elements and fields given by this restriction. 209d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 210d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 21139b2de37Sjeremylt for the unknowns corresponding to element i, where 212d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 213d979a051Sjeremylt [0, lsize - 1]. 214d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 215d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 21639b2de37Sjeremylt 21739b2de37Sjeremylt Returns: 21839b2de37Sjeremylt elemrestriction: Ceed ElemRestiction""" 21939b2de37Sjeremylt 220d979a051Sjeremylt return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 221d979a051Sjeremylt offsets, memtype=memtype, cmode=cmode) 22239b2de37Sjeremylt 223d979a051Sjeremylt def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides): 224d979a051Sjeremylt """Ceed Identity ElemRestriction: strided restriction from local vectors 225d979a051Sjeremylt to elements. 22639b2de37Sjeremylt 22739b2de37Sjeremylt Args: 228d979a051Sjeremylt nelem: number of elements described by the restriction 22939b2de37Sjeremylt elemsize: size (number of nodes) per element 230a8d32208Sjeremylt ncomp: number of field components per interpolation node 231a8d32208Sjeremylt (1 for scalar fields) 232d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 233d979a051Sjeremylt the elements and fields given by this restriction. 23469a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 23569a53589Sjeremylt The data for node i, component j, element k in the 23669a53589Sjeremylt L-vector is given by 23769a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 23839b2de37Sjeremylt 23939b2de37Sjeremylt Returns: 24069a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 24139b2de37Sjeremylt 2427a7b0fa3SJed Brown return StridedElemRestriction( 243d979a051Sjeremylt self, nelem, elemsize, ncomp, lsize, strides) 24439b2de37Sjeremylt 245d979a051Sjeremylt def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 246d979a051Sjeremylt lsize, offsets, memtype=lib.CEED_MEM_HOST, 247d979a051Sjeremylt cmode=lib.CEED_COPY_VALUES): 248d979a051Sjeremylt """Ceed Blocked ElemRestriction: blocked restriction from local vectors 249d979a051Sjeremylt to elements. 25039b2de37Sjeremylt 25139b2de37Sjeremylt Args: 252d979a051Sjeremylt nelem: number of elements described by the restriction 25339b2de37Sjeremylt elemsize: size (number of nodes) per element 25439b2de37Sjeremylt blksize: number of elements in a block 255d979a051Sjeremylt ncomp: number of field components per interpolation node 256d979a051Sjeremylt (1 for scalar fields) 257d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 258d979a051Sjeremylt the elements and fields given by this restriction. 259d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 260d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 26139b2de37Sjeremylt for the unknowns corresponding to element i, where 262d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 263d979a051Sjeremylt [0, lsize - 1]. The backend will permute and pad this 26439b2de37Sjeremylt array to the desired ordering for the blocksize, which is 26539b2de37Sjeremylt typically given by the backend. The default reordering is 26639b2de37Sjeremylt to interlace elements. 267d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 268d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 26939b2de37Sjeremylt 27039b2de37Sjeremylt Returns: 27139b2de37Sjeremylt elemrestriction: Ceed Blocked ElemRestiction""" 27239b2de37Sjeremylt 273d979a051Sjeremylt return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, 274d979a051Sjeremylt compstride, lsize, offsets, 275d979a051Sjeremylt memtype=memtype, cmode=cmode) 27639b2de37Sjeremylt 277d979a051Sjeremylt def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp, 278d979a051Sjeremylt lsize, strides): 279d979a051Sjeremylt """Ceed Blocked Strided ElemRestriction: blocked and strided restriction 280d979a051Sjeremylt from local vectors to elements. 28169a53589Sjeremylt 28269a53589Sjeremylt Args: 283d979a051Sjeremylt nelem: number of elements described in the restriction 28469a53589Sjeremylt elemsize: size (number of nodes) per element 28569a53589Sjeremylt blksize: number of elements in a block 28669a53589Sjeremylt ncomp: number of field components per interpolation node 28769a53589Sjeremylt (1 for scalar fields) 288d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 289d979a051Sjeremylt the elements and fields given by this restriction. 29069a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 29169a53589Sjeremylt The data for node i, component j, element k in the 29269a53589Sjeremylt L-vector is given by 29369a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 29469a53589Sjeremylt 29569a53589Sjeremylt Returns: 29669a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 29769a53589Sjeremylt 29869a53589Sjeremylt return BlockedStridedElemRestriction(self, nelem, elemsize, blksize, 299d979a051Sjeremylt ncomp, lsize, strides) 30069a53589Sjeremylt 30139b2de37Sjeremylt # CeedBasis 30239b2de37Sjeremylt def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 30339b2de37Sjeremylt qref1d, qweight1d): 30439b2de37Sjeremylt """Ceed Tensor H1 Basis: finite element tensor-product basis objects for 30539b2de37Sjeremylt H^1 discretizations. 30639b2de37Sjeremylt 30739b2de37Sjeremylt Args: 30839b2de37Sjeremylt dim: topological dimension 30939b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 31039b2de37Sjeremylt P1d: number of nodes in one dimension 31139b2de37Sjeremylt Q1d: number of quadrature points in one dimension 312d979a051Sjeremylt *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix 313d979a051Sjeremylt expressing the values of nodal basis functions at 314d979a051Sjeremylt quadrature points 315d979a051Sjeremylt *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix 316d979a051Sjeremylt expressing the derivatives of nodal basis functions at 317d979a051Sjeremylt quadrature points 318d979a051Sjeremylt *qref1d: Numpy array of length Q1d holding the locations of 319d979a051Sjeremylt quadrature points on the 1D reference element [-1, 1] 320d979a051Sjeremylt *qweight1d: Numpy array of length Q1d holding the quadrature 321d979a051Sjeremylt weights on the reference element 32239b2de37Sjeremylt 32339b2de37Sjeremylt Returns: 32439b2de37Sjeremylt basis: Ceed Basis""" 32539b2de37Sjeremylt 32639b2de37Sjeremylt return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 32739b2de37Sjeremylt qref1d, qweight1d) 32839b2de37Sjeremylt 32939b2de37Sjeremylt def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode): 330d979a051Sjeremylt """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange 331d979a051Sjeremylt basis objects for H^1 discretizations. 33239b2de37Sjeremylt 33339b2de37Sjeremylt Args: 33439b2de37Sjeremylt dim: topological dimension 33539b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 33639b2de37Sjeremylt P: number of Gauss-Lobatto nodes in one dimension. The 33739b2de37Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 33839b2de37Sjeremylt Q: number of quadrature points in one dimension 33939b2de37Sjeremylt qmode: distribution of the Q quadrature points (affects order of 34039b2de37Sjeremylt accuracy for the quadrature) 34139b2de37Sjeremylt 34239b2de37Sjeremylt Returns: 34339b2de37Sjeremylt basis: Ceed Basis""" 34439b2de37Sjeremylt 34539b2de37Sjeremylt return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode) 34639b2de37Sjeremylt 34739b2de37Sjeremylt def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight): 348d979a051Sjeremylt """Ceed H1 Basis: finite element non tensor-product basis for H^1 349d979a051Sjeremylt discretizations. 35039b2de37Sjeremylt 35139b2de37Sjeremylt Args: 35239b2de37Sjeremylt topo: topology of the element, e.g. hypercube, simplex, etc 35339b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 35439b2de37Sjeremylt nnodes: total number of nodes 35539b2de37Sjeremylt nqpts: total number of quadrature points 356d979a051Sjeremylt *interp: Numpy array holding the row-major (nqpts * nnodes) matrix 357d979a051Sjeremylt expressing the values of nodal basis functions at 35839b2de37Sjeremylt quadrature points 359*97c1c57aSSebastian Grimberg *grad: Numpy array holding the row-major (dim * nqpts * nnodes) 360d979a051Sjeremylt matrix expressing the derivatives of nodal basis functions 361d979a051Sjeremylt at quadrature points 362d979a051Sjeremylt *qref: Numpy array of length (nqpts * dim) holding the locations of 363d979a051Sjeremylt quadrature points on the reference element [-1, 1] 364d979a051Sjeremylt *qweight: Numpy array of length nnodes holding the quadrature 365d979a051Sjeremylt weights on the reference element 36639b2de37Sjeremylt 36739b2de37Sjeremylt Returns: 36839b2de37Sjeremylt basis: Ceed Basis""" 36939b2de37Sjeremylt 3707a7b0fa3SJed Brown return BasisH1(self, topo, ncomp, nnodes, nqpts, 3717a7b0fa3SJed Brown interp, grad, qref, qweight) 37239b2de37Sjeremylt 373*97c1c57aSSebastian Grimberg def BasisHdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight): 374*97c1c57aSSebastian Grimberg """Ceed Hdiv Basis: finite element non tensor-product basis for H(div) 375*97c1c57aSSebastian Grimberg discretizations. 376*97c1c57aSSebastian Grimberg 377*97c1c57aSSebastian Grimberg Args: 378*97c1c57aSSebastian Grimberg topo: topology of the element, e.g. hypercube, simplex, etc 379*97c1c57aSSebastian Grimberg ncomp: number of field components (1 for scalar fields) 380*97c1c57aSSebastian Grimberg nnodes: total number of nodes 381*97c1c57aSSebastian Grimberg nqpts: total number of quadrature points 382*97c1c57aSSebastian Grimberg *interp: Numpy array holding the row-major (dim * nqpts * nnodes) 383*97c1c57aSSebastian Grimberg matrix expressing the values of basis functions at 384*97c1c57aSSebastian Grimberg quadrature points 385*97c1c57aSSebastian Grimberg *div: Numpy array holding the row-major (nqpts * nnodes) matrix 386*97c1c57aSSebastian Grimberg expressing the divergence of basis functions at 387*97c1c57aSSebastian Grimberg quadrature points 388*97c1c57aSSebastian Grimberg *qref: Numpy array of length (nqpts * dim) holding the locations of 389*97c1c57aSSebastian Grimberg quadrature points on the reference element [-1, 1] 390*97c1c57aSSebastian Grimberg *qweight: Numpy array of length nnodes holding the quadrature 391*97c1c57aSSebastian Grimberg weights on the reference element 392*97c1c57aSSebastian Grimberg 393*97c1c57aSSebastian Grimberg Returns: 394*97c1c57aSSebastian Grimberg basis: Ceed Basis""" 395*97c1c57aSSebastian Grimberg 396*97c1c57aSSebastian Grimberg return BasisHdiv(self, topo, ncomp, nnodes, nqpts, 397*97c1c57aSSebastian Grimberg interp, div, qref, qweight) 398*97c1c57aSSebastian Grimberg 399*97c1c57aSSebastian Grimberg def BasisHcurl(self, topo, ncomp, nnodes, nqpts, 400*97c1c57aSSebastian Grimberg interp, curl, qref, qweight): 401*97c1c57aSSebastian Grimberg """Ceed Hcurl Basis: finite element non tensor-product basis for H(curl) 402*97c1c57aSSebastian Grimberg discretizations. 403*97c1c57aSSebastian Grimberg 404*97c1c57aSSebastian Grimberg Args: 405*97c1c57aSSebastian Grimberg topo: topology of the element, e.g. hypercube, simplex, etc 406*97c1c57aSSebastian Grimberg ncomp: number of field components (1 for scalar fields) 407*97c1c57aSSebastian Grimberg nnodes: total number of nodes 408*97c1c57aSSebastian Grimberg nqpts: total number of quadrature points 409*97c1c57aSSebastian Grimberg *interp: Numpy array holding the row-major (dim * nqpts * nnodes) 410*97c1c57aSSebastian Grimberg matrix expressing the values of basis functions at 411*97c1c57aSSebastian Grimberg quadrature points 412*97c1c57aSSebastian Grimberg *curl: Numpy array holding the row-major (curlcomp * nqpts * nnodes), 413*97c1c57aSSebastian Grimberg curlcomp = 1 if dim < 3 else dim, matrix expressing the curl 414*97c1c57aSSebastian Grimberg of basis functions at quadrature points 415*97c1c57aSSebastian Grimberg *qref: Numpy array of length (nqpts * dim) holding the locations of 416*97c1c57aSSebastian Grimberg quadrature points on the reference element [-1, 1] 417*97c1c57aSSebastian Grimberg *qweight: Numpy array of length nnodes holding the quadrature 418*97c1c57aSSebastian Grimberg weights on the reference element 419*97c1c57aSSebastian Grimberg 420*97c1c57aSSebastian Grimberg Returns: 421*97c1c57aSSebastian Grimberg basis: Ceed Basis""" 422*97c1c57aSSebastian Grimberg 423*97c1c57aSSebastian Grimberg return BasisHcurl(self, topo, ncomp, nnodes, nqpts, 424*97c1c57aSSebastian Grimberg interp, curl, qref, qweight) 425*97c1c57aSSebastian Grimberg 42639b2de37Sjeremylt # CeedQFunction 42739b2de37Sjeremylt def QFunction(self, vlength, f, source): 428d979a051Sjeremylt """Ceed QFunction: point-wise operation at quadrature points for 429d979a051Sjeremylt evaluating volumetric terms. 43039b2de37Sjeremylt 43139b2de37Sjeremylt Args: 43239b2de37Sjeremylt vlength: vector length. Caller must ensure that number of quadrature 43339b2de37Sjeremylt points is a multiple of vlength 43439b2de37Sjeremylt f: ctypes function pointer to evaluate action at quadrature points 43539b2de37Sjeremylt source: absolute path to source of QFunction, 43639b2de37Sjeremylt "\\abs_path\\file.h:function_name 43739b2de37Sjeremylt 43839b2de37Sjeremylt Returns: 43939b2de37Sjeremylt qfunction: Ceed QFunction""" 44039b2de37Sjeremylt 44139b2de37Sjeremylt return QFunction(self, vlength, f, source) 44239b2de37Sjeremylt 44339b2de37Sjeremylt def QFunctionByName(self, name): 44439b2de37Sjeremylt """Ceed QFunction By Name: point-wise operation at quadrature points 44539b2de37Sjeremylt from a given gallery, for evaluating volumetric terms. 44639b2de37Sjeremylt 44739b2de37Sjeremylt Args: 44839b2de37Sjeremylt name: name of QFunction to use from gallery 44939b2de37Sjeremylt 45039b2de37Sjeremylt Returns: 45139b2de37Sjeremylt qfunction: Ceed QFunction By Name""" 45239b2de37Sjeremylt 45339b2de37Sjeremylt return QFunctionByName(self, name) 45439b2de37Sjeremylt 45539b2de37Sjeremylt def IdentityQFunction(self, size, inmode, outmode): 45639b2de37Sjeremylt """Ceed Idenity QFunction: identity qfunction operation. 45739b2de37Sjeremylt 45839b2de37Sjeremylt Args: 45939b2de37Sjeremylt size: size of the qfunction fields 46039b2de37Sjeremylt **inmode: CeedEvalMode for input to Ceed QFunction 46139b2de37Sjeremylt **outmode: CeedEvalMode for output to Ceed QFunction 46239b2de37Sjeremylt 46339b2de37Sjeremylt Returns: 46439b2de37Sjeremylt qfunction: Ceed Identity QFunction""" 46539b2de37Sjeremylt 46639b2de37Sjeremylt return IdentityQFunction(self, size, inmode, outmode) 46739b2de37Sjeremylt 468777ff853SJeremy L Thompson def QFunctionContext(self): 469777ff853SJeremy L Thompson """Ceed QFunction Context: stores Ceed QFunction user context data. 470777ff853SJeremy L Thompson 471777ff853SJeremy L Thompson Returns: 472777ff853SJeremy L Thompson userContext: Ceed QFunction Context""" 473777ff853SJeremy L Thompson 474777ff853SJeremy L Thompson return QFunctionContext(self) 475777ff853SJeremy L Thompson 47639b2de37Sjeremylt # CeedOperator 47739b2de37Sjeremylt def Operator(self, qf, dqf=None, qdfT=None): 47839b2de37Sjeremylt """Ceed Operator: composed FE-type operations on vectors. 47939b2de37Sjeremylt 48039b2de37Sjeremylt Args: 481d979a051Sjeremylt qf: QFunction defining the action of the operator at quadrature 482d979a051Sjeremylt points 48339b2de37Sjeremylt **dqf: QFunction defining the action of the Jacobian, default None 484d979a051Sjeremylt **dqfT: QFunction defining the action of the transpose of the 485d979a051Sjeremylt Jacobian, default None 48639b2de37Sjeremylt 48739b2de37Sjeremylt Returns: 48839b2de37Sjeremylt operator: Ceed Operator""" 48939b2de37Sjeremylt 49039b2de37Sjeremylt return Operator(self, qf, dqf, qdfT) 49139b2de37Sjeremylt 49239b2de37Sjeremylt def CompositeOperator(self): 49339b2de37Sjeremylt """Composite Ceed Operator: composition of multiple CeedOperators. 49439b2de37Sjeremylt 49539b2de37Sjeremylt Returns: 49639b2de37Sjeremylt operator: Ceed Composite Operator""" 49739b2de37Sjeremylt 49839b2de37Sjeremylt return CompositeOperator(self) 49939b2de37Sjeremylt 50039b2de37Sjeremylt # Destructor 50139b2de37Sjeremylt def __del__(self): 50239b2de37Sjeremylt # libCEED call 50339b2de37Sjeremylt lib.CeedDestroy(self._pointer) 50439b2de37Sjeremylt 50539b2de37Sjeremylt# ------------------------------------------------------------------------------ 506