petsc4py.PETSc.FE#
- class petsc4py.PETSc.FE#
Bases:
Object
A PETSc object that manages a finite element space.
Enumerations
Methods Summary
create
([comm])Create an empty
FE
object.createDefault
(dim, nc, isSimplex[, qorder, ...])Create a
FE
for basic FEM computation.createLagrange
(dim, nc, isSimplex, k[, ...])Create a
FE
for the basic Lagrange space of degree k.destroy
()Destroy the
FE
object.Return the
Space
used for the approximation of theFE
solution.Return the dimension of the finite element space on a cell.
Return the
DualSpace
used to define the inner product for theFE
.Return the
Quad
used to calculate inner products on faces.Return the number of components in the element.
Return the number of DOFs.
Return the
Quad
used to calculate inner products.Return the spatial dimension of the element.
Return the tile sizes for evaluation.
setBasisSpace
(sp)Set the
Space
used for the approximation of the solution.setDualSpace
(dspace)Set the
DualSpace
used to define the inner product.setFaceQuadrature
(quad)Set the
Quad
used to calculate inner products on faces.Set parameters in a
FE
from the options database.setNumComponents
(comp)Set the number of field components in the element.
setQuadrature
(quad)Set the
Quad
used to calculate inner products.setTileSizes
(blockSize, numBlocks, ...)Set the tile sizes for evaluation.
setType
(fe_type)Build a particular
FE
.setUp
()Construct data structures for the
FE
after theType
has been set.view
([viewer])View a
FE
object.viewFromOptions
(name[, obj])View from a
FE
based on values in the options database.Methods Documentation
- create(comm=None)#
Create an empty
FE
object.Collective.
The type can then be set with
setType
.- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- createDefault(dim, nc, isSimplex, qorder=DETERMINE, prefix=None, comm=None)#
Create a
FE
for basic FEM computation.Collective.
- Parameters:
dim (int) – The spatial dimension.
nc (int) – The number of components.
isSimplex (bool) – Flag for simplex reference cell, otherwise it’s a tensor product.
qorder (int) – The quadrature order or
DETERMINE
to useSpace
polynomial degree.comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
See also
- createLagrange(dim, nc, isSimplex, k, qorder=DETERMINE, comm=None)#
Create a
FE
for the basic Lagrange space of degree k.Collective.
- Parameters:
dim (int) – The spatial dimension.
nc (int) – The number of components.
isSimplex (bool) – Flag for simplex reference cell, otherwise it’s a tensor product.
k (int) – The degree of the space.
qorder (int) – The quadrature order or
DETERMINE
to useSpace
polynomial degree.comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
See also
- destroy()#
Destroy the
FE
object.Collective.
See also
Source code at petsc4py/PETSc/FE.pyx:38
- Return type:
- getBasisSpace()#
Return the
Space
used for the approximation of theFE
solution.Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:386
- Return type:
- getDimension()#
Return the dimension of the finite element space on a cell.
Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:180
- Return type:
- getDualSpace()#
Return the
DualSpace
used to define the inner product for theFE
.Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:441
- Return type:
- getFaceQuadrature()#
Return the
Quad
used to calculate inner products on faces.Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:316
- Return type:
- getNumComponents()#
Return the number of components in the element.
Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:208
- Return type:
- getNumDof()#
Return the number of DOFs.
Not collective.
Return the number of DOFs (dual basis vectors) associated with mesh points on the reference cell of a given dimension.
See also
Source code at petsc4py/PETSc/FE.pyx:240
- Return type:
- getQuadrature()#
Return the
Quad
used to calculate inner products.Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:166
- Return type:
- getSpatialDimension()#
Return the spatial dimension of the element.
Not collective.
See also
Source code at petsc4py/PETSc/FE.pyx:194
- Return type:
- getTileSizes()#
Return the tile sizes for evaluation.
Not collective.
- Returns:
- Return type:
See also
- setBasisSpace(sp)#
Set the
Space
used for the approximation of the solution.Not collective.
See also
- setFaceQuadrature(quad)#
Set the
Quad
used to calculate inner products on faces.Not collective.
See also
- setFromOptions()#
Set parameters in a
FE
from the options database.Collective.
Source code at petsc4py/PETSc/FE.pyx:417
- Return type:
- setNumComponents(comp)#
Set the number of field components in the element.
Not collective.
See also
- setTileSizes(blockSize, numBlocks, batchSize, numBatches)#
Set the tile sizes for evaluation.
Not collective.
- Parameters:
- Return type:
See also
- setUp()#
Construct data structures for the
FE
after theType
has been set.Collective.
See also
Source code at petsc4py/PETSc/FE.pyx:429
- Return type: