1*3d8e8822SJeremy L Thompson# Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors 2*3d8e8822SJeremy L Thompson# All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 339b2de37Sjeremylt# 4*3d8e8822SJeremy L Thompson# SPDX-License-Identifier: BSD-2-Clause 539b2de37Sjeremylt# 6*3d8e8822SJeremy 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 1639b2de37Sjeremyltfrom .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1 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 # QR factorization 181e15f9bd0SJeremy L Thompson def qr_factorization(self, mat, tau, m, n): 182e15f9bd0SJeremy L Thompson """Return QR Factorization of a matrix. 183e15f9bd0SJeremy L Thompson 184e15f9bd0SJeremy L Thompson Args: 185e15f9bd0SJeremy L Thompson ceed: Ceed context currently in use 186e15f9bd0SJeremy L Thompson *mat: Numpy array holding the row-major matrix to be factorized in place 187e15f9bd0SJeremy L Thompson *tau: Numpy array to hold the vector of lengt m of scaling factors 188e15f9bd0SJeremy L Thompson m: number of rows 189e15f9bd0SJeremy L Thompson n: numbef of columns""" 190e15f9bd0SJeremy L Thompson 191e15f9bd0SJeremy L Thompson # Setup arguments 192e15f9bd0SJeremy L Thompson mat_pointer = ffi.new("CeedScalar *") 193e15f9bd0SJeremy L Thompson mat_pointer = ffi.cast( 194e15f9bd0SJeremy L Thompson "CeedScalar *", 195e15f9bd0SJeremy L Thompson mat.__array_interface__['data'][0]) 196e15f9bd0SJeremy L Thompson 197e15f9bd0SJeremy L Thompson tau_pointer = ffi.new("CeedScalar *") 198e15f9bd0SJeremy L Thompson tau_pointer = ffi.cast( 199e15f9bd0SJeremy L Thompson "CeedScalar *", 200e15f9bd0SJeremy L Thompson tau.__array_interface__['data'][0]) 201e15f9bd0SJeremy L Thompson 202e15f9bd0SJeremy L Thompson # libCEED call 203e15f9bd0SJeremy L Thompson err_code = lib.CeedQRFactorization( 204e15f9bd0SJeremy L Thompson self._pointer[0], mat_pointer, tau_pointer, m, n) 205e15f9bd0SJeremy L Thompson self._check_error(err_code) 206e15f9bd0SJeremy L Thompson 207e15f9bd0SJeremy L Thompson return mat, tau 208e15f9bd0SJeremy L Thompson 209e15f9bd0SJeremy L Thompson # Symmetric Schur decomposition 210e15f9bd0SJeremy L Thompson def symmetric_schur_decomposition(self, mat, n): 211e15f9bd0SJeremy L Thompson """Return symmetric Schur decomposition of a symmetric matrix 212e15f9bd0SJeremy L Thompson via symmetric QR factorization. 213e15f9bd0SJeremy L Thompson 214e15f9bd0SJeremy L Thompson Args: 215e15f9bd0SJeremy L Thompson ceed: Ceed context currently in use 216e15f9bd0SJeremy L Thompson *mat: Numpy array holding the row-major matrix to be factorized in place 217e15f9bd0SJeremy L Thompson n: number of rows/columns 218e15f9bd0SJeremy L Thompson 219e15f9bd0SJeremy L Thompson Returns: 220e15f9bd0SJeremy L Thompson lbda: Numpy array of length n holding eigenvalues""" 221e15f9bd0SJeremy L Thompson 222e15f9bd0SJeremy L Thompson # Setup arguments 223e15f9bd0SJeremy L Thompson mat_pointer = ffi.new("CeedScalar *") 224e15f9bd0SJeremy L Thompson mat_pointer = ffi.cast( 225e15f9bd0SJeremy L Thompson "CeedScalar *", 226e15f9bd0SJeremy L Thompson mat.__array_interface__['data'][0]) 227e15f9bd0SJeremy L Thompson 22880a9ef05SNatalie Beams lbda = np.empty(n, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 229e15f9bd0SJeremy L Thompson l_pointer = ffi.new("CeedScalar *") 230e15f9bd0SJeremy L Thompson l_pointer = ffi.cast( 231e15f9bd0SJeremy L Thompson "CeedScalar *", 232e15f9bd0SJeremy L Thompson lbda.__array_interface__['data'][0]) 233e15f9bd0SJeremy L Thompson 234e15f9bd0SJeremy L Thompson # libCEED call 235e15f9bd0SJeremy L Thompson err_code = lib.CeedSymmetricSchurDecomposition( 236e15f9bd0SJeremy L Thompson self._pointer[0], mat_pointer, l_pointer, n) 237e15f9bd0SJeremy L Thompson self._check_error(err_code) 238e15f9bd0SJeremy L Thompson 239e15f9bd0SJeremy L Thompson return lbda 240e15f9bd0SJeremy L Thompson 241e15f9bd0SJeremy L Thompson # Simultaneous Diagonalization 242e15f9bd0SJeremy L Thompson def simultaneous_diagonalization(self, matA, matB, n): 243e15f9bd0SJeremy L Thompson """Return Simultaneous Diagonalization of two matrices. 244e15f9bd0SJeremy L Thompson 245e15f9bd0SJeremy L Thompson Args: 246e15f9bd0SJeremy L Thompson ceed: Ceed context currently in use 247e15f9bd0SJeremy L Thompson *matA: Numpy array holding the row-major matrix to be factorized with 248e15f9bd0SJeremy L Thompson eigenvalues 249e15f9bd0SJeremy L Thompson *matB: Numpy array holding the row-major matrix to be factorized to identity 250e15f9bd0SJeremy L Thompson n: number of rows/columns 251e15f9bd0SJeremy L Thompson 252e15f9bd0SJeremy L Thompson Returns: 253e15f9bd0SJeremy L Thompson (x, lbda): Numpy array holding the row-major orthogonal matrix and 254e15f9bd0SJeremy L Thompson Numpy array holding the vector of length n of generalized 255e15f9bd0SJeremy L Thompson eigenvalues""" 256e15f9bd0SJeremy L Thompson 257e15f9bd0SJeremy L Thompson # Setup arguments 258e15f9bd0SJeremy L Thompson matA_pointer = ffi.new("CeedScalar *") 259e15f9bd0SJeremy L Thompson matA_pointer = ffi.cast( 260e15f9bd0SJeremy L Thompson "CeedScalar *", 261e15f9bd0SJeremy L Thompson matA.__array_interface__['data'][0]) 262e15f9bd0SJeremy L Thompson 263e15f9bd0SJeremy L Thompson matB_pointer = ffi.new("CeedScalar *") 264e15f9bd0SJeremy L Thompson matB_pointer = ffi.cast( 265e15f9bd0SJeremy L Thompson "CeedScalar *", 266e15f9bd0SJeremy L Thompson matB.__array_interface__['data'][0]) 267e15f9bd0SJeremy L Thompson 26880a9ef05SNatalie Beams lbda = np.empty(n, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 269e15f9bd0SJeremy L Thompson l_pointer = ffi.new("CeedScalar *") 270e15f9bd0SJeremy L Thompson l_pointer = ffi.cast( 271e15f9bd0SJeremy L Thompson "CeedScalar *", 272e15f9bd0SJeremy L Thompson lbda.__array_interface__['data'][0]) 273e15f9bd0SJeremy L Thompson 27480a9ef05SNatalie Beams x = np.empty(n * n, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 275e15f9bd0SJeremy L Thompson x_pointer = ffi.new("CeedScalar *") 276e15f9bd0SJeremy L Thompson x_pointer = ffi.cast("CeedScalar *", x.__array_interface__['data'][0]) 277e15f9bd0SJeremy L Thompson 278e15f9bd0SJeremy L Thompson # libCEED call 279e15f9bd0SJeremy L Thompson err_code = lib.CeedSimultaneousDiagonalization(self._pointer[0], matA_pointer, matB_pointer, 280e15f9bd0SJeremy L Thompson x_pointer, l_pointer, n) 281e15f9bd0SJeremy L Thompson self._check_error(err_code) 282e15f9bd0SJeremy L Thompson 283e15f9bd0SJeremy L Thompson return x, lbda 284e15f9bd0SJeremy L Thompson 285e15f9bd0SJeremy L Thompson # --- libCEED objects --- 286e15f9bd0SJeremy L Thompson 28739b2de37Sjeremylt # CeedVector 28839b2de37Sjeremylt def Vector(self, size): 28939b2de37Sjeremylt """Ceed Vector: storing and manipulating vectors. 29039b2de37Sjeremylt 29139b2de37Sjeremylt Args: 29239b2de37Sjeremylt size: length of vector 29339b2de37Sjeremylt 29439b2de37Sjeremylt Returns: 29539b2de37Sjeremylt vector: Ceed Vector""" 29639b2de37Sjeremylt 29739b2de37Sjeremylt return Vector(self, size) 29839b2de37Sjeremylt 29939b2de37Sjeremylt # CeedElemRestriction 300d979a051Sjeremylt def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets, 301d979a051Sjeremylt memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 30239b2de37Sjeremylt """Ceed ElemRestriction: restriction from local vectors to elements. 30339b2de37Sjeremylt 30439b2de37Sjeremylt Args: 305d979a051Sjeremylt nelem: number of elements described by the restriction 30639b2de37Sjeremylt elemsize: size (number of nodes) per element 307d979a051Sjeremylt ncomp: number of field components per interpolation node 308d979a051Sjeremylt (1 for scalar fields) 309d979a051Sjeremylt compstride: Stride between components for the same L-vector "node". 310d979a051Sjeremylt Data for node i, component k can be found in the 311d979a051Sjeremylt L-vector at index [offsets[i] + k*compstride]. 312d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 313d979a051Sjeremylt the elements and fields given by this restriction. 314d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 315d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 31639b2de37Sjeremylt for the unknowns corresponding to element i, where 317d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 318d979a051Sjeremylt [0, lsize - 1]. 319d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 320d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 32139b2de37Sjeremylt 32239b2de37Sjeremylt Returns: 32339b2de37Sjeremylt elemrestriction: Ceed ElemRestiction""" 32439b2de37Sjeremylt 325d979a051Sjeremylt return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 326d979a051Sjeremylt offsets, memtype=memtype, cmode=cmode) 32739b2de37Sjeremylt 328d979a051Sjeremylt def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides): 329d979a051Sjeremylt """Ceed Identity ElemRestriction: strided restriction from local vectors 330d979a051Sjeremylt to elements. 33139b2de37Sjeremylt 33239b2de37Sjeremylt Args: 333d979a051Sjeremylt nelem: number of elements described by the restriction 33439b2de37Sjeremylt elemsize: size (number of nodes) per element 335a8d32208Sjeremylt ncomp: number of field components per interpolation node 336a8d32208Sjeremylt (1 for scalar fields) 337d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 338d979a051Sjeremylt the elements and fields given by this restriction. 33969a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 34069a53589Sjeremylt The data for node i, component j, element k in the 34169a53589Sjeremylt L-vector is given by 34269a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 34339b2de37Sjeremylt 34439b2de37Sjeremylt Returns: 34569a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 34639b2de37Sjeremylt 3477a7b0fa3SJed Brown return StridedElemRestriction( 348d979a051Sjeremylt self, nelem, elemsize, ncomp, lsize, strides) 34939b2de37Sjeremylt 350d979a051Sjeremylt def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 351d979a051Sjeremylt lsize, offsets, memtype=lib.CEED_MEM_HOST, 352d979a051Sjeremylt cmode=lib.CEED_COPY_VALUES): 353d979a051Sjeremylt """Ceed Blocked ElemRestriction: blocked restriction from local vectors 354d979a051Sjeremylt to elements. 35539b2de37Sjeremylt 35639b2de37Sjeremylt Args: 357d979a051Sjeremylt nelem: number of elements described by the restriction 35839b2de37Sjeremylt elemsize: size (number of nodes) per element 35939b2de37Sjeremylt blksize: number of elements in a block 360d979a051Sjeremylt ncomp: number of field components per interpolation node 361d979a051Sjeremylt (1 for scalar fields) 362d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 363d979a051Sjeremylt the elements and fields given by this restriction. 364d979a051Sjeremylt *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 365d979a051Sjeremylt ordered list of the offsets (into the input Ceed Vector) 36639b2de37Sjeremylt for the unknowns corresponding to element i, where 367d979a051Sjeremylt 0 <= i < nelem. All offsets must be in the range 368d979a051Sjeremylt [0, lsize - 1]. The backend will permute and pad this 36939b2de37Sjeremylt array to the desired ordering for the blocksize, which is 37039b2de37Sjeremylt typically given by the backend. The default reordering is 37139b2de37Sjeremylt to interlace elements. 372d979a051Sjeremylt **memtype: memory type of the offsets array, default CEED_MEM_HOST 373d979a051Sjeremylt **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 37439b2de37Sjeremylt 37539b2de37Sjeremylt Returns: 37639b2de37Sjeremylt elemrestriction: Ceed Blocked ElemRestiction""" 37739b2de37Sjeremylt 378d979a051Sjeremylt return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, 379d979a051Sjeremylt compstride, lsize, offsets, 380d979a051Sjeremylt memtype=memtype, cmode=cmode) 38139b2de37Sjeremylt 382d979a051Sjeremylt def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp, 383d979a051Sjeremylt lsize, strides): 384d979a051Sjeremylt """Ceed Blocked Strided ElemRestriction: blocked and strided restriction 385d979a051Sjeremylt from local vectors to elements. 38669a53589Sjeremylt 38769a53589Sjeremylt Args: 388d979a051Sjeremylt nelem: number of elements described in the restriction 38969a53589Sjeremylt elemsize: size (number of nodes) per element 39069a53589Sjeremylt blksize: number of elements in a block 39169a53589Sjeremylt ncomp: number of field components per interpolation node 39269a53589Sjeremylt (1 for scalar fields) 393d979a051Sjeremylt lsize: The size of the L-vector. This vector may be larger than 394d979a051Sjeremylt the elements and fields given by this restriction. 39569a53589Sjeremylt *strides: Array for strides between [nodes, components, elements]. 39669a53589Sjeremylt The data for node i, component j, element k in the 39769a53589Sjeremylt L-vector is given by 39869a53589Sjeremylt i*strides[0] + j*strides[1] + k*strides[2] 39969a53589Sjeremylt 40069a53589Sjeremylt Returns: 40169a53589Sjeremylt elemrestriction: Ceed Strided ElemRestiction""" 40269a53589Sjeremylt 40369a53589Sjeremylt return BlockedStridedElemRestriction(self, nelem, elemsize, blksize, 404d979a051Sjeremylt ncomp, lsize, strides) 40569a53589Sjeremylt 40639b2de37Sjeremylt # CeedBasis 40739b2de37Sjeremylt def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 40839b2de37Sjeremylt qref1d, qweight1d): 40939b2de37Sjeremylt """Ceed Tensor H1 Basis: finite element tensor-product basis objects for 41039b2de37Sjeremylt H^1 discretizations. 41139b2de37Sjeremylt 41239b2de37Sjeremylt Args: 41339b2de37Sjeremylt dim: topological dimension 41439b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 41539b2de37Sjeremylt P1d: number of nodes in one dimension 41639b2de37Sjeremylt Q1d: number of quadrature points in one dimension 417d979a051Sjeremylt *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix 418d979a051Sjeremylt expressing the values of nodal basis functions at 419d979a051Sjeremylt quadrature points 420d979a051Sjeremylt *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix 421d979a051Sjeremylt expressing the derivatives of nodal basis functions at 422d979a051Sjeremylt quadrature points 423d979a051Sjeremylt *qref1d: Numpy array of length Q1d holding the locations of 424d979a051Sjeremylt quadrature points on the 1D reference element [-1, 1] 425d979a051Sjeremylt *qweight1d: Numpy array of length Q1d holding the quadrature 426d979a051Sjeremylt weights on the reference element 42739b2de37Sjeremylt 42839b2de37Sjeremylt Returns: 42939b2de37Sjeremylt basis: Ceed Basis""" 43039b2de37Sjeremylt 43139b2de37Sjeremylt return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 43239b2de37Sjeremylt qref1d, qweight1d) 43339b2de37Sjeremylt 43439b2de37Sjeremylt def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode): 435d979a051Sjeremylt """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange 436d979a051Sjeremylt basis objects for H^1 discretizations. 43739b2de37Sjeremylt 43839b2de37Sjeremylt Args: 43939b2de37Sjeremylt dim: topological dimension 44039b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 44139b2de37Sjeremylt P: number of Gauss-Lobatto nodes in one dimension. The 44239b2de37Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 44339b2de37Sjeremylt Q: number of quadrature points in one dimension 44439b2de37Sjeremylt qmode: distribution of the Q quadrature points (affects order of 44539b2de37Sjeremylt accuracy for the quadrature) 44639b2de37Sjeremylt 44739b2de37Sjeremylt Returns: 44839b2de37Sjeremylt basis: Ceed Basis""" 44939b2de37Sjeremylt 45039b2de37Sjeremylt return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode) 45139b2de37Sjeremylt 45239b2de37Sjeremylt def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight): 453d979a051Sjeremylt """Ceed H1 Basis: finite element non tensor-product basis for H^1 454d979a051Sjeremylt discretizations. 45539b2de37Sjeremylt 45639b2de37Sjeremylt Args: 45739b2de37Sjeremylt topo: topology of the element, e.g. hypercube, simplex, etc 45839b2de37Sjeremylt ncomp: number of field components (1 for scalar fields) 45939b2de37Sjeremylt nnodes: total number of nodes 46039b2de37Sjeremylt nqpts: total number of quadrature points 461d979a051Sjeremylt *interp: Numpy array holding the row-major (nqpts * nnodes) matrix 462d979a051Sjeremylt expressing the values of nodal basis functions at 46339b2de37Sjeremylt quadrature points 464d979a051Sjeremylt *grad: Numpy array holding the row-major (nqpts * dim * nnodes) 465d979a051Sjeremylt matrix expressing the derivatives of nodal basis functions 466d979a051Sjeremylt at quadrature points 467d979a051Sjeremylt *qref: Numpy array of length (nqpts * dim) holding the locations of 468d979a051Sjeremylt quadrature points on the reference element [-1, 1] 469d979a051Sjeremylt *qweight: Numpy array of length nnodes holding the quadrature 470d979a051Sjeremylt weights on the reference element 47139b2de37Sjeremylt 47239b2de37Sjeremylt Returns: 47339b2de37Sjeremylt basis: Ceed Basis""" 47439b2de37Sjeremylt 4757a7b0fa3SJed Brown return BasisH1(self, topo, ncomp, nnodes, nqpts, 4767a7b0fa3SJed Brown interp, grad, qref, qweight) 47739b2de37Sjeremylt 47839b2de37Sjeremylt # CeedQFunction 47939b2de37Sjeremylt def QFunction(self, vlength, f, source): 480d979a051Sjeremylt """Ceed QFunction: point-wise operation at quadrature points for 481d979a051Sjeremylt evaluating volumetric terms. 48239b2de37Sjeremylt 48339b2de37Sjeremylt Args: 48439b2de37Sjeremylt vlength: vector length. Caller must ensure that number of quadrature 48539b2de37Sjeremylt points is a multiple of vlength 48639b2de37Sjeremylt f: ctypes function pointer to evaluate action at quadrature points 48739b2de37Sjeremylt source: absolute path to source of QFunction, 48839b2de37Sjeremylt "\\abs_path\\file.h:function_name 48939b2de37Sjeremylt 49039b2de37Sjeremylt Returns: 49139b2de37Sjeremylt qfunction: Ceed QFunction""" 49239b2de37Sjeremylt 49339b2de37Sjeremylt return QFunction(self, vlength, f, source) 49439b2de37Sjeremylt 49539b2de37Sjeremylt def QFunctionByName(self, name): 49639b2de37Sjeremylt """Ceed QFunction By Name: point-wise operation at quadrature points 49739b2de37Sjeremylt from a given gallery, for evaluating volumetric terms. 49839b2de37Sjeremylt 49939b2de37Sjeremylt Args: 50039b2de37Sjeremylt name: name of QFunction to use from gallery 50139b2de37Sjeremylt 50239b2de37Sjeremylt Returns: 50339b2de37Sjeremylt qfunction: Ceed QFunction By Name""" 50439b2de37Sjeremylt 50539b2de37Sjeremylt return QFunctionByName(self, name) 50639b2de37Sjeremylt 50739b2de37Sjeremylt def IdentityQFunction(self, size, inmode, outmode): 50839b2de37Sjeremylt """Ceed Idenity QFunction: identity qfunction operation. 50939b2de37Sjeremylt 51039b2de37Sjeremylt Args: 51139b2de37Sjeremylt size: size of the qfunction fields 51239b2de37Sjeremylt **inmode: CeedEvalMode for input to Ceed QFunction 51339b2de37Sjeremylt **outmode: CeedEvalMode for output to Ceed QFunction 51439b2de37Sjeremylt 51539b2de37Sjeremylt Returns: 51639b2de37Sjeremylt qfunction: Ceed Identity QFunction""" 51739b2de37Sjeremylt 51839b2de37Sjeremylt return IdentityQFunction(self, size, inmode, outmode) 51939b2de37Sjeremylt 520777ff853SJeremy L Thompson def QFunctionContext(self): 521777ff853SJeremy L Thompson """Ceed QFunction Context: stores Ceed QFunction user context data. 522777ff853SJeremy L Thompson 523777ff853SJeremy L Thompson Returns: 524777ff853SJeremy L Thompson userContext: Ceed QFunction Context""" 525777ff853SJeremy L Thompson 526777ff853SJeremy L Thompson return QFunctionContext(self) 527777ff853SJeremy L Thompson 52839b2de37Sjeremylt # CeedOperator 52939b2de37Sjeremylt def Operator(self, qf, dqf=None, qdfT=None): 53039b2de37Sjeremylt """Ceed Operator: composed FE-type operations on vectors. 53139b2de37Sjeremylt 53239b2de37Sjeremylt Args: 533d979a051Sjeremylt qf: QFunction defining the action of the operator at quadrature 534d979a051Sjeremylt points 53539b2de37Sjeremylt **dqf: QFunction defining the action of the Jacobian, default None 536d979a051Sjeremylt **dqfT: QFunction defining the action of the transpose of the 537d979a051Sjeremylt Jacobian, default None 53839b2de37Sjeremylt 53939b2de37Sjeremylt Returns: 54039b2de37Sjeremylt operator: Ceed Operator""" 54139b2de37Sjeremylt 54239b2de37Sjeremylt return Operator(self, qf, dqf, qdfT) 54339b2de37Sjeremylt 54439b2de37Sjeremylt def CompositeOperator(self): 54539b2de37Sjeremylt """Composite Ceed Operator: composition of multiple CeedOperators. 54639b2de37Sjeremylt 54739b2de37Sjeremylt Returns: 54839b2de37Sjeremylt operator: Ceed Composite Operator""" 54939b2de37Sjeremylt 55039b2de37Sjeremylt return CompositeOperator(self) 55139b2de37Sjeremylt 55239b2de37Sjeremylt # Destructor 55339b2de37Sjeremylt def __del__(self): 55439b2de37Sjeremylt # libCEED call 55539b2de37Sjeremylt lib.CeedDestroy(self._pointer) 55639b2de37Sjeremylt 55739b2de37Sjeremylt# ------------------------------------------------------------------------------ 558