1# Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors 2# All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3# 4# SPDX-License-Identifier: BSD-2-Clause 5# 6# This file is part of CEED: http://github.com/ceed 7 8from _ceed_cffi import ffi, lib 9import sys 10import os 11import io 12import numpy as np 13import tempfile 14from abc import ABC 15from .ceed_vector import Vector 16from .ceed_basis import BasisTensorH1, BasisTensorH1Lagrange, BasisH1, BasisHdiv, BasisHcurl 17from .ceed_elemrestriction import ElemRestriction, OrientedElemRestriction, CurlOrientedElemRestriction, StridedElemRestriction 18from .ceed_elemrestriction import BlockedElemRestriction, BlockedOrientedElemRestriction, BlockedCurlOrientedElemRestriction, BlockedStridedElemRestriction 19from .ceed_qfunction import QFunction, QFunctionByName, IdentityQFunction 20from .ceed_qfunctioncontext import QFunctionContext 21from .ceed_operator import Operator, CompositeOperator 22from .ceed_constants import * 23 24# ------------------------------------------------------------------------------ 25 26 27class Ceed(): 28 """Ceed: core components.""" 29 30 # Constructor 31 def __init__(self, resource="/cpu/self", on_error="store"): 32 # libCEED object 33 self._pointer = ffi.new("Ceed *") 34 35 # libCEED call 36 resourceAscii = ffi.new("char[]", resource.encode("ascii")) 37 os.environ["CEED_ERROR_HANDLER"] = "return" 38 err_code = lib.CeedInit(resourceAscii, self._pointer) 39 if err_code: 40 raise Exception("Error initializing backend resource: " + resource) 41 error_handlers = dict( 42 store="CeedErrorStore", 43 abort="CeedErrorAbort", 44 ) 45 lib.CeedSetErrorHandler( 46 self._pointer[0], ffi.addressof( 47 lib, error_handlers[on_error])) 48 49 # Representation 50 def __repr__(self): 51 return "<Ceed instance at " + hex(id(self)) + ">" 52 53 # String conversion for print() to stdout 54 def __str__(self): 55 """View a Ceed via print().""" 56 57 # libCEED call 58 with tempfile.NamedTemporaryFile() as key_file: 59 with open(key_file.name, 'r+') as stream_file: 60 stream = ffi.cast("FILE *", stream_file) 61 62 err_code = lib.CeedView(self._pointer[0], stream) 63 self._check_error(err_code) 64 65 stream_file.seek(0) 66 out_string = stream_file.read() 67 68 return out_string 69 70 # Error handler 71 def _check_error(self, err_code): 72 """Check return code and retrieve error message for non-zero code""" 73 if (err_code != lib.CEED_ERROR_SUCCESS): 74 message = ffi.new("char **") 75 lib.CeedGetErrorMessage(self._pointer[0], message) 76 raise Exception(ffi.string(message[0]).decode("UTF-8")) 77 78 # Get Resource 79 def get_resource(self): 80 """Get the full resource name for a Ceed context. 81 82 Returns: 83 resource: resource name""" 84 85 # libCEED call 86 resource = ffi.new("char **") 87 err_code = lib.CeedGetResource(self._pointer[0], resource) 88 self._check_error(err_code) 89 90 return ffi.string(resource[0]).decode("UTF-8") 91 92 # Preferred MemType 93 def get_preferred_memtype(self): 94 """Return Ceed preferred memory type. 95 96 Returns: 97 memtype: Ceed preferred memory type""" 98 99 # libCEED call 100 memtype = ffi.new("CeedMemType *", MEM_HOST) 101 err_code = lib.CeedGetPreferredMemType(self._pointer[0], memtype) 102 self._check_error(err_code) 103 104 return memtype[0] 105 106 # Convenience function to get CeedScalar type 107 def scalar_type(self): 108 """Return dtype string for CeedScalar. 109 110 Returns: 111 np_dtype: String naming numpy data type corresponding to CeedScalar""" 112 113 return scalar_types[lib.CEED_SCALAR_TYPE] 114 115 # --- Basis utility functions --- 116 117 # Gauss quadrature 118 def gauss_quadrature(self, q): 119 """Construct a Gauss-Legendre quadrature. 120 121 Args: 122 Q: number of quadrature points (integrates polynomials of 123 degree 2*Q-1 exactly) 124 125 Returns: 126 (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1] 127 and array of length Q to hold the weights""" 128 129 # Setup arguments 130 qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 131 qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 132 133 qref1d_pointer = ffi.new("CeedScalar *") 134 qref1d_pointer = ffi.cast( 135 "CeedScalar *", 136 qref1d.__array_interface__['data'][0]) 137 138 qweight1d_pointer = ffi.new("CeedScalar *") 139 qweight1d_pointer = ffi.cast( 140 "CeedScalar *", 141 qweight1d.__array_interface__['data'][0]) 142 143 # libCEED call 144 err_code = lib.CeedGaussQuadrature(q, qref1d_pointer, qweight1d_pointer) 145 self._check_error(err_code) 146 147 return qref1d, qweight1d 148 149 # Lobatto quadrature 150 def lobatto_quadrature(self, q): 151 """Construct a Gauss-Legendre-Lobatto quadrature. 152 153 Args: 154 q: number of quadrature points (integrates polynomials of 155 degree 2*Q-3 exactly) 156 157 Returns: 158 (qref1d, qweight1d): array of length Q to hold the abscissa on [-1, 1] 159 and array of length Q to hold the weights""" 160 161 # Setup arguments 162 qref1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 163 qref1d_pointer = ffi.new("CeedScalar *") 164 qref1d_pointer = ffi.cast( 165 "CeedScalar *", 166 qref1d.__array_interface__['data'][0]) 167 168 qweight1d = np.empty(q, dtype=scalar_types[lib.CEED_SCALAR_TYPE]) 169 qweight1d_pointer = ffi.new("CeedScalar *") 170 qweight1d_pointer = ffi.cast( 171 "CeedScalar *", 172 qweight1d.__array_interface__['data'][0]) 173 174 # libCEED call 175 err_code = lib.CeedLobattoQuadrature( 176 q, qref1d_pointer, qweight1d_pointer) 177 self._check_error(err_code) 178 179 return qref1d, qweight1d 180 181 # --- libCEED objects --- 182 183 # CeedVector 184 def Vector(self, size): 185 """Ceed Vector: storing and manipulating vectors. 186 187 Args: 188 size: length of vector 189 190 Returns: 191 vector: Ceed Vector""" 192 193 return Vector(self, size) 194 195 # CeedElemRestriction 196 def ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, offsets, 197 memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 198 """Ceed ElemRestriction: restriction from local vectors to elements. 199 200 Args: 201 nelem: number of elements described by the restriction 202 elemsize: size (number of nodes) per element 203 ncomp: number of field components per interpolation node 204 (1 for scalar fields) 205 compstride: Stride between components for the same L-vector "node". 206 Data for node i, component k can be found in the 207 L-vector at index [offsets[i] + k*compstride]. 208 lsize: The size of the L-vector. This vector may be larger than 209 the elements and fields given by this restriction. 210 *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 211 ordered list of the offsets (into the input Ceed Vector) 212 for the unknowns corresponding to element i, where 213 0 <= i < nelem. All offsets must be in the range 214 [0, lsize - 1]. 215 **memtype: memory type of the offsets array, default CEED_MEM_HOST 216 **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 217 218 Returns: 219 elemrestriction: Ceed ElemRestiction""" 220 221 return ElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 222 offsets, memtype=memtype, cmode=cmode) 223 224 def OrientedElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 225 offsets, orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 226 """Ceed Oriented ElemRestriction: oriented restriction from local vectors 227 to elements. 228 229 Args: 230 nelem: number of elements described by the restriction 231 elemsize: size (number of nodes) per element 232 ncomp: number of field components per interpolation node 233 (1 for scalar fields) 234 compstride: Stride between components for the same L-vector "node". 235 Data for node i, component k can be found in the 236 L-vector at index [offsets[i] + k*compstride]. 237 lsize: The size of the L-vector. This vector may be larger than 238 the elements and fields given by this restriction. 239 *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 240 ordered list of the offsets (into the input Ceed Vector) 241 for the unknowns corresponding to element i, where 242 0 <= i < nelem. All offsets must be in the range 243 [0, lsize - 1]. 244 *orients: Numpy array of shape [nelem, elemsize]. Row i holds the 245 ordered list of the orientations for the unknowns 246 corresponding to element i, with bool false used for 247 positively oriented and true to flip the orientation. 248 **memtype: memory type of the offsets array, default CEED_MEM_HOST 249 **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 250 251 Returns: 252 elemrestriction: Ceed Oriented ElemRestiction""" 253 254 return OrientedElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 255 offsets, orients, memtype=memtype, cmode=cmode) 256 257 def CurlOrientedElemRestriction(self, nelem, elemsize, ncomp, compstride, lsize, 258 offsets, curl_orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 259 """Ceed Curl Oriented ElemRestriction: curl-oriented restriction from local 260 vectors to elements. 261 262 Args: 263 nelem: number of elements described by the restriction 264 elemsize: size (number of nodes) per element 265 ncomp: number of field components per interpolation node 266 (1 for scalar fields) 267 compstride: Stride between components for the same L-vector "node". 268 Data for node i, component k can be found in the 269 L-vector at index [offsets[i] + k*compstride]. 270 lsize: The size of the L-vector. This vector may be larger than 271 the elements and fields given by this restriction. 272 *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 273 ordered list of the offsets (into the input Ceed Vector) 274 for the unknowns corresponding to element i, where 275 0 <= i < nelem. All offsets must be in the range 276 [0, lsize - 1]. 277 *curl_orients: Numpy array of shape [nelem, 3 * elemsize]. Row i holds 278 a row-major tridiagonal matrix (curl_orients[i, 0] = 279 curl_orients[i, 3 * elemsize - 1] = 0, where 0 <= i < nelem) 280 which is applied to the element unknowns upon restriction. 281 **memtype: memory type of the offsets array, default CEED_MEM_HOST 282 **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 283 284 Returns: 285 elemrestriction: Ceed Curl Oriented ElemRestiction""" 286 287 return CurlOrientedElemRestriction( 288 self, nelem, elemsize, ncomp, compstride, lsize, offsets, curl_orients, memtype=memtype, cmode=cmode) 289 290 def StridedElemRestriction(self, nelem, elemsize, ncomp, lsize, strides): 291 """Ceed Identity ElemRestriction: strided restriction from local vectors 292 to elements. 293 294 Args: 295 nelem: number of elements described by the restriction 296 elemsize: size (number of nodes) per element 297 ncomp: number of field components per interpolation node 298 (1 for scalar fields) 299 lsize: The size of the L-vector. This vector may be larger than 300 the elements and fields given by this restriction. 301 *strides: Array for strides between [nodes, components, elements]. 302 The data for node i, component j, element k in the 303 L-vector is given by 304 i*strides[0] + j*strides[1] + k*strides[2] 305 306 Returns: 307 elemrestriction: Ceed Strided ElemRestiction""" 308 309 return StridedElemRestriction( 310 self, nelem, elemsize, ncomp, lsize, strides) 311 312 def BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 313 lsize, offsets, memtype=lib.CEED_MEM_HOST, 314 cmode=lib.CEED_COPY_VALUES): 315 """Ceed Blocked ElemRestriction: blocked restriction from local vectors 316 to elements. 317 318 Args: 319 nelem: number of elements described by the restriction 320 elemsize: size (number of nodes) per element 321 blksize: number of elements in a block 322 ncomp: number of field components per interpolation node 323 (1 for scalar fields) 324 lsize: The size of the L-vector. This vector may be larger than 325 the elements and fields given by this restriction. 326 *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 327 ordered list of the offsets (into the input Ceed Vector) 328 for the unknowns corresponding to element i, where 329 0 <= i < nelem. All offsets must be in the range 330 [0, lsize - 1]. The backend will permute and pad this 331 array to the desired ordering for the blocksize, which is 332 typically given by the backend. The default reordering is 333 to interlace elements. 334 **memtype: memory type of the offsets array, default CEED_MEM_HOST 335 **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 336 337 Returns: 338 elemrestriction: Ceed Blocked ElemRestiction""" 339 340 return BlockedElemRestriction(self, nelem, elemsize, blksize, ncomp, 341 compstride, lsize, offsets, 342 memtype=memtype, cmode=cmode) 343 344 def BlockedOrientedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 345 lsize, offsets, orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 346 """Ceed Blocked Oriented ElemRestriction: blocked oriented restriction 347 from local vectors to elements. 348 349 Args: 350 nelem: number of elements described by the restriction 351 elemsize: size (number of nodes) per element 352 blksize: number of elements in a block 353 ncomp: number of field components per interpolation node 354 (1 for scalar fields) 355 compstride: Stride between components for the same L-vector "node". 356 Data for node i, component k can be found in the 357 L-vector at index [offsets[i] + k*compstride]. 358 lsize: The size of the L-vector. This vector may be larger than 359 the elements and fields given by this restriction. 360 *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 361 ordered list of the offsets (into the input Ceed Vector) 362 for the unknowns corresponding to element i, where 363 0 <= i < nelem. All offsets must be in the range 364 [0, lsize - 1]. The backend will permute and pad this 365 array to the desired ordering for the blocksize, which is 366 typically given by the backend. The default reordering is 367 to interlace elements. 368 *orients: Numpy array of shape [nelem, elemsize]. Row i holds the 369 ordered list of the orientations for the unknowns 370 corresponding to element i, with ool false is used for 371 positively oriented and true to flip the orientation. 372 **memtype: memory type of the offsets array, default CEED_MEM_HOST 373 **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 374 375 Returns: 376 elemrestriction: Ceed Blocked Oriented ElemRestiction""" 377 378 return BlockedOrientedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, lsize, 379 offsets, orients, memtype=memtype, cmode=cmode) 380 381 def BlockedCurlOrientedElemRestriction(self, nelem, elemsize, blksize, ncomp, compstride, 382 lsize, offsets, curl_orients, memtype=lib.CEED_MEM_HOST, cmode=lib.CEED_COPY_VALUES): 383 """Ceed Blocked Curl Oriented ElemRestriction: blocked curl-oriented 384 restriction from local vectors to elements. 385 386 Args: 387 nelem: number of elements described by the restriction 388 elemsize: size (number of nodes) per element 389 blksize: number of elements in a block 390 ncomp: number of field components per interpolation node 391 (1 for scalar fields) 392 compstride: Stride between components for the same L-vector "node". 393 Data for node i, component k can be found in the 394 L-vector at index [offsets[i] + k*compstride]. 395 lsize: The size of the L-vector. This vector may be larger than 396 the elements and fields given by this restriction. 397 *offsets: Numpy array of shape [nelem, elemsize]. Row i holds the 398 ordered list of the offsets (into the input Ceed Vector) 399 for the unknowns corresponding to element i, where 400 0 <= i < nelem. All offsets must be in the range 401 [0, lsize - 1]. The backend will permute and pad this 402 array to the desired ordering for the blocksize, which is 403 typically given by the backend. The default reordering is 404 to interlace elements. 405 *curl_orients: Numpy array of shape [nelem, 3 * elemsize]. Row i holds 406 a row-major tridiagonal matrix (curl_orients[i, 0] = 407 curl_orients[i, 3 * elemsize - 1] = 0, where 0 <= i < nelem) 408 which is applied to the element unknowns upon restriction. 409 **memtype: memory type of the offsets array, default CEED_MEM_HOST 410 **cmode: copy mode for the offsets array, default CEED_COPY_VALUES 411 412 Returns: 413 elemrestriction: Ceed Blocked Curl Oriented ElemRestiction""" 414 415 return BlockedCurlOrientedElemRestriction( 416 self, nelem, elemsize, blksize, ncomp, compstride, lsize, offsets, curl_orients, memtype=memtype, cmode=cmode) 417 418 def BlockedStridedElemRestriction(self, nelem, elemsize, blksize, ncomp, 419 lsize, strides): 420 """Ceed Blocked Strided ElemRestriction: blocked and strided restriction 421 from local vectors to elements. 422 423 Args: 424 nelem: number of elements described in the restriction 425 elemsize: size (number of nodes) per element 426 blksize: number of elements in a block 427 ncomp: number of field components per interpolation node 428 (1 for scalar fields) 429 lsize: The size of the L-vector. This vector may be larger than 430 the elements and fields given by this restriction. 431 *strides: Array for strides between [nodes, components, elements]. 432 The data for node i, component j, element k in the 433 L-vector is given by 434 i*strides[0] + j*strides[1] + k*strides[2] 435 436 Returns: 437 elemrestriction: Ceed Blocked Strided ElemRestiction""" 438 439 return BlockedStridedElemRestriction(self, nelem, elemsize, blksize, 440 ncomp, lsize, strides) 441 442 # CeedBasis 443 def BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 444 qref1d, qweight1d): 445 """Ceed Tensor H1 Basis: finite element tensor-product basis objects for 446 H^1 discretizations. 447 448 Args: 449 dim: topological dimension 450 ncomp: number of field components (1 for scalar fields) 451 P1d: number of nodes in one dimension 452 Q1d: number of quadrature points in one dimension 453 *interp1d: Numpy array holding the row-major (Q1d * P1d) matrix 454 expressing the values of nodal basis functions at 455 quadrature points 456 *grad1d: Numpy array holding the row-major (Q1d * P1d) matrix 457 expressing the derivatives of nodal basis functions at 458 quadrature points 459 *qref1d: Numpy array of length Q1d holding the locations of 460 quadrature points on the 1D reference element [-1, 1] 461 *qweight1d: Numpy array of length Q1d holding the quadrature 462 weights on the reference element 463 464 Returns: 465 basis: Ceed Basis""" 466 467 return BasisTensorH1(self, dim, ncomp, P1d, Q1d, interp1d, grad1d, 468 qref1d, qweight1d) 469 470 def BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode): 471 """Ceed Tensor H1 Lagrange Basis: finite element tensor-product Lagrange 472 basis objects for H^1 discretizations. 473 474 Args: 475 dim: topological dimension 476 ncomp: number of field components (1 for scalar fields) 477 P: number of Gauss-Lobatto nodes in one dimension. The 478 polynomial degree of the resulting Q_k element is k=P-1. 479 Q: number of quadrature points in one dimension 480 qmode: distribution of the Q quadrature points (affects order of 481 accuracy for the quadrature) 482 483 Returns: 484 basis: Ceed Basis""" 485 486 return BasisTensorH1Lagrange(self, dim, ncomp, P, Q, qmode) 487 488 def BasisH1(self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight): 489 """Ceed H1 Basis: finite element non tensor-product basis for H^1 490 discretizations. 491 492 Args: 493 topo: topology of the element, e.g. hypercube, simplex, etc 494 ncomp: number of field components (1 for scalar fields) 495 nnodes: total number of nodes 496 nqpts: total number of quadrature points 497 *interp: Numpy array holding the row-major (nqpts * nnodes) matrix 498 expressing the values of nodal basis functions at 499 quadrature points 500 *grad: Numpy array holding the row-major (dim * nqpts * nnodes) 501 matrix expressing the derivatives of nodal basis functions 502 at quadrature points 503 *qref: Numpy array of length (nqpts * dim) holding the locations of 504 quadrature points on the reference element [-1, 1] 505 *qweight: Numpy array of length nnodes holding the quadrature 506 weights on the reference element 507 508 Returns: 509 basis: Ceed Basis""" 510 511 return BasisH1(self, topo, ncomp, nnodes, nqpts, 512 interp, grad, qref, qweight) 513 514 def BasisHdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight): 515 """Ceed Hdiv Basis: finite element non tensor-product basis for H(div) 516 discretizations. 517 518 Args: 519 topo: topology of the element, e.g. hypercube, simplex, etc 520 ncomp: number of field components (1 for scalar fields) 521 nnodes: total number of nodes 522 nqpts: total number of quadrature points 523 *interp: Numpy array holding the row-major (dim * nqpts * nnodes) 524 matrix expressing the values of basis functions at 525 quadrature points 526 *div: Numpy array holding the row-major (nqpts * nnodes) matrix 527 expressing the divergence of basis functions at 528 quadrature points 529 *qref: Numpy array of length (nqpts * dim) holding the locations of 530 quadrature points on the reference element [-1, 1] 531 *qweight: Numpy array of length nnodes holding the quadrature 532 weights on the reference element 533 534 Returns: 535 basis: Ceed Basis""" 536 537 return BasisHdiv(self, topo, ncomp, nnodes, nqpts, 538 interp, div, qref, qweight) 539 540 def BasisHcurl(self, topo, ncomp, nnodes, nqpts, 541 interp, curl, qref, qweight): 542 """Ceed Hcurl Basis: finite element non tensor-product basis for H(curl) 543 discretizations. 544 545 Args: 546 topo: topology of the element, e.g. hypercube, simplex, etc 547 ncomp: number of field components (1 for scalar fields) 548 nnodes: total number of nodes 549 nqpts: total number of quadrature points 550 *interp: Numpy array holding the row-major (dim * nqpts * nnodes) 551 matrix expressing the values of basis functions at 552 quadrature points 553 *curl: Numpy array holding the row-major (curlcomp * nqpts * nnodes), 554 curlcomp = 1 if dim < 3 else dim, matrix expressing the curl 555 of basis functions at quadrature points 556 *qref: Numpy array of length (nqpts * dim) holding the locations of 557 quadrature points on the reference element [-1, 1] 558 *qweight: Numpy array of length nnodes holding the quadrature 559 weights on the reference element 560 561 Returns: 562 basis: Ceed Basis""" 563 564 return BasisHcurl(self, topo, ncomp, nnodes, nqpts, 565 interp, curl, qref, qweight) 566 567 # CeedQFunction 568 def QFunction(self, vlength, f, source): 569 """Ceed QFunction: point-wise operation at quadrature points for 570 evaluating volumetric terms. 571 572 Args: 573 vlength: vector length. Caller must ensure that number of quadrature 574 points is a multiple of vlength 575 f: ctypes function pointer to evaluate action at quadrature points 576 source: absolute path to source of QFunction, 577 "\\abs_path\\file.h:function_name 578 579 Returns: 580 qfunction: Ceed QFunction""" 581 582 return QFunction(self, vlength, f, source) 583 584 def QFunctionByName(self, name): 585 """Ceed QFunction By Name: point-wise operation at quadrature points 586 from a given gallery, for evaluating volumetric terms. 587 588 Args: 589 name: name of QFunction to use from gallery 590 591 Returns: 592 qfunction: Ceed QFunction By Name""" 593 594 return QFunctionByName(self, name) 595 596 def IdentityQFunction(self, size, inmode, outmode): 597 """Ceed Idenity QFunction: identity qfunction operation. 598 599 Args: 600 size: size of the qfunction fields 601 **inmode: CeedEvalMode for input to Ceed QFunction 602 **outmode: CeedEvalMode for output to Ceed QFunction 603 604 Returns: 605 qfunction: Ceed Identity QFunction""" 606 607 return IdentityQFunction(self, size, inmode, outmode) 608 609 def QFunctionContext(self): 610 """Ceed QFunction Context: stores Ceed QFunction user context data. 611 612 Returns: 613 userContext: Ceed QFunction Context""" 614 615 return QFunctionContext(self) 616 617 # CeedOperator 618 def Operator(self, qf, dqf=None, qdfT=None): 619 """Ceed Operator: composed FE-type operations on vectors. 620 621 Args: 622 qf: QFunction defining the action of the operator at quadrature 623 points 624 **dqf: QFunction defining the action of the Jacobian, default None 625 **dqfT: QFunction defining the action of the transpose of the 626 Jacobian, default None 627 628 Returns: 629 operator: Ceed Operator""" 630 631 return Operator(self, qf, dqf, qdfT) 632 633 def CompositeOperator(self): 634 """Composite Ceed Operator: composition of multiple CeedOperators. 635 636 Returns: 637 operator: Ceed Composite Operator""" 638 639 return CompositeOperator(self) 640 641 # Destructor 642 def __del__(self): 643 # libCEED call 644 lib.CeedDestroy(self._pointer) 645 646# ------------------------------------------------------------------------------ 647