petsc4py.PETSc.FE#

class petsc4py.PETSc.FE#

Bases: Object

A PETSc object that manages a finite element space.

Enumerations

Type

The finite element types.

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.

getBasisSpace()

Return the Space used for the approximation of the FE solution.

getDimension()

Return the dimension of the finite element space on a cell.

getDualSpace()

Return the DualSpace used to define the inner product for the FE.

getFaceQuadrature()

Return the Quad used to calculate inner products on faces.

getNumComponents()

Return the number of components in the element.

getNumDof()

Return the number of DOFs.

getQuadrature()

Return the Quad used to calculate inner products.

getSpatialDimension()

Return the spatial dimension of the element.

getTileSizes()

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.

setFromOptions()

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 the Type has been set.

view([viewer])

View a FE object.

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:

Self

Source code at petsc4py/PETSc/FE.pyx:53

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 use Space polynomial degree.

  • prefix (str) – The options prefix, or None.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

Source code at petsc4py/PETSc/FE.pyx:76

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 use Space polynomial degree.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

Source code at petsc4py/PETSc/FE.pyx:122

destroy()#

Destroy the FE object.

Collective.

See also

PetscFEDestroy

Source code at petsc4py/PETSc/FE.pyx:40

Return type:

Self

getBasisSpace()#

Return the Space used for the approximation of the FE solution.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:385

Return type:

Space

getDimension()#

Return the dimension of the finite element space on a cell.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:180

Return type:

int

getDualSpace()#

Return the DualSpace used to define the inner product for the FE.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:440

Return type:

DualSpace

getFaceQuadrature()#

Return the Quad used to calculate inner products on faces.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:315

Return type:

Quad

getNumComponents()#

Return the number of components in the element.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:208

Return type:

int

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

PetscFEGetNumDof

Source code at petsc4py/PETSc/FE.pyx:240

Return type:

ndarray

getQuadrature()#

Return the Quad used to calculate inner products.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:166

Return type:

Quad

getSpatialDimension()#

Return the spatial dimension of the element.

Not collective.

Source code at petsc4py/PETSc/FE.pyx:194

Return type:

int

getTileSizes()#

Return the tile sizes for evaluation.

Not collective.

Returns:

  • blockSize (int) – The number of elements in a block.

  • numBlocks (int) – The number of blocks in a batch.

  • batchSize (int) – The number of elements in a batch.

  • numBatches (int) – The number of batches in a chunk.

Return type:

tuple(int, int, int, int)

Source code at petsc4py/PETSc/FE.pyx:259

setBasisSpace(sp)#

Set the Space used for the approximation of the solution.

Not collective.

Parameters:

sp (Space) – The Space object.

Return type:

None

Source code at petsc4py/PETSc/FE.pyx:399

setDualSpace(dspace)#

Set the DualSpace used to define the inner product.

Not collective.

Parameters:

dspace (DualSpace) – The DualSpace object.

Return type:

None

Source code at petsc4py/PETSc/FE.pyx:454

setFaceQuadrature(quad)#

Set the Quad used to calculate inner products on faces.

Not collective.

Parameters:

quad (Quad) – The Quad object.

Return type:

Quad

Source code at petsc4py/PETSc/FE.pyx:347

setFromOptions()#

Set parameters in a FE from the options database.

Collective.

Source code at petsc4py/PETSc/FE.pyx:416

Return type:

None

setNumComponents(comp)#

Set the number of field components in the element.

Not collective.

Parameters:

comp (int) – The number of field components.

Return type:

None

Source code at petsc4py/PETSc/FE.pyx:222

setQuadrature(quad)#

Set the Quad used to calculate inner products.

Not collective.

Parameters:

quad (Quad) – The Quad object.

Return type:

Self

Source code at petsc4py/PETSc/FE.pyx:329

setTileSizes(blockSize, numBlocks, batchSize, numBatches)#

Set the tile sizes for evaluation.

Not collective.

Parameters:
  • blockSize (int) – The number of elements in a block.

  • numBlocks (int) – The number of blocks in a batch.

  • batchSize (int) – The number of elements in a batch.

  • numBatches (int) – The number of batches in a chunk.

Return type:

None

Source code at petsc4py/PETSc/FE.pyx:285

setType(fe_type)#

Build a particular FE.

Collective.

Parameters:

fe_type (Type | str) – The kind of FEM space.

Return type:

Self

See also

PetscFESetType

Source code at petsc4py/PETSc/FE.pyx:365

setUp()#

Construct data structures for the FE after the Type has been set.

Collective.

See also

PetscFESetUp

Source code at petsc4py/PETSc/FE.pyx:428

Return type:

None

view(viewer=None)#

View a FE object.

Collective.

Parameters:

viewer (Viewer | None) – A Viewer to display the graph.

Return type:

None

See also

PetscFEView

Source code at petsc4py/PETSc/FE.pyx:21