dolfinx.cpp.fem

FEM module

Classes

Constant_complex128(self, c[, writable, order])

Value constant with respect to integration domain

Constant_complex64(self, c[, writable, order])

Value constant with respect to integration domain

Constant_float32(self, c[, writable, order])

Value constant with respect to integration domain

Constant_float64(self, c[, writable, order])

Value constant with respect to integration domain

CoordinateElement_float32()

Coordinate map element

CoordinateElement_float64()

Coordinate map element

DirichletBC_complex128([order])

Object for representing Dirichlet (essential) boundary conditions

DirichletBC_complex64([order])

Object for representing Dirichlet (essential) boundary conditions

DirichletBC_float32([order])

Object for representing Dirichlet (essential) boundary conditions

DirichletBC_float64([order])

Object for representing Dirichlet (essential) boundary conditions

DofMap(self, element_dof_layout, index_map, ...)

DofMap object

ElementDofLayout(self, block_size, ...)

Object describing the layout of dofs on a cell

Expression_complex128(self, coefficients, ...)

An Expression

Expression_complex64(self, coefficients, ...)

An Expression

Expression_float32(self, coefficients, ...)

An Expression

Expression_float64(self, coefficients, ...)

An Expression

FiniteElement_float32(self, ufcx_element)

Finite element object

FiniteElement_float64(self, ufcx_element)

Finite element object

Form_complex128([order, writable, shape, order])

Variational form object

Form_complex64([order, writable, shape, order])

Variational form object

Form_float32([order, writable, shape, order])

Variational form object

Form_float64([order, writable, shape, order])

Variational form object

FunctionSpace_float32(self, mesh, element, ...)

Finite element function space

FunctionSpace_float64(self, mesh, element, ...)

Finite element function space

Function_complex128()

A finite element function

Function_complex64()

A finite element function

Function_float32()

A finite element function

Function_float64()

A finite element function

IntegralType

class dolfinx.cpp.fem.Constant_complex128(self, c: ndarray[dtype=complex128, writable=False, order='C'])

Bases: object

Value constant with respect to integration domain

Create a constant from a value array

property dtype

(self) -> str

property value

(self) -> numpy.ndarray[dtype=complex128, ]

class dolfinx.cpp.fem.Constant_complex64(self, c: ndarray[dtype=complex64, writable=False, order='C'])

Bases: object

Value constant with respect to integration domain

Create a constant from a value array

property dtype

(self) -> str

property value

(self) -> numpy.ndarray[dtype=complex64, ]

class dolfinx.cpp.fem.Constant_float32(self, c: ndarray[dtype=float32, writable=False, order='C'])

Bases: object

Value constant with respect to integration domain

Create a constant from a value array

property dtype

(self) -> str

property value

(self) -> numpy.ndarray[dtype=float32, ]

class dolfinx.cpp.fem.Constant_float64(self, c: ndarray[dtype=float64, writable=False, order='C'])

Bases: object

Value constant with respect to integration domain

Create a constant from a value array

property dtype

(self) -> str

property value

(self) -> numpy.ndarray[dtype=float64, ]

class dolfinx.cpp.fem.CoordinateElement_float32(self, celltype: dolfinx.cpp.mesh.CellType, degree: int)
class dolfinx.cpp.fem.CoordinateElement_float32(self, element: basix._basixcpp.FiniteElement_float32)
class dolfinx.cpp.fem.CoordinateElement_float32(self, celltype: dolfinx.cpp.mesh.CellType, degree: int, variant: int)

Bases: object

Coordinate map element

create_dof_layout
property degree

(self) -> int

property dim

(self) -> int

property dtype

(self) -> str

pull_back
push_forward
property variant

(self) -> int

class dolfinx.cpp.fem.CoordinateElement_float64(self, celltype: dolfinx.cpp.mesh.CellType, degree: int)
class dolfinx.cpp.fem.CoordinateElement_float64(self, element: basix._basixcpp.FiniteElement_float64)
class dolfinx.cpp.fem.CoordinateElement_float64(self, celltype: dolfinx.cpp.mesh.CellType, degree: int, variant: int)

Bases: object

Coordinate map element

create_dof_layout
property degree

(self) -> int

property dim

(self) -> int

property dtype

(self) -> str

pull_back
push_forward
property variant

(self) -> int

class dolfinx.cpp.fem.DirichletBC_complex128(self, g: ndarray[dtype=complex128, writable=False, order='C'], dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float64)
class dolfinx.cpp.fem.DirichletBC_complex128(self, g: dolfinx.cpp.fem.Constant_complex128, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float64)
class dolfinx.cpp.fem.DirichletBC_complex128(self, g: dolfinx.cpp.fem.Function_complex128, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'])
class dolfinx.cpp.fem.DirichletBC_complex128(self, g: dolfinx.cpp.fem.Function_complex128, dofs: list[ndarray[dtype=int32, writable=False, shape=(*), order='C']], V: dolfinx.cpp.fem.FunctionSpace_float64)

Bases: object

Object for representing Dirichlet (essential) boundary conditions

dof_indices
property dtype

(self) -> str

property function_space

(self) -> dolfinx::fem::FunctionSpace<double>

property value

(self) -> Union[dolfinx::fem::Function<std::complex<double>, double>, dolfinx::fem::Constant<std::complex<double> >]

class dolfinx.cpp.fem.DirichletBC_complex64(self, g: ndarray[dtype=complex64, writable=False, order='C'], dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float32)
class dolfinx.cpp.fem.DirichletBC_complex64(self, g: dolfinx.cpp.fem.Constant_complex64, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float32)
class dolfinx.cpp.fem.DirichletBC_complex64(self, g: dolfinx.cpp.fem.Function_complex64, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'])
class dolfinx.cpp.fem.DirichletBC_complex64(self, g: dolfinx.cpp.fem.Function_complex64, dofs: list[ndarray[dtype=int32, writable=False, shape=(*), order='C']], V: dolfinx.cpp.fem.FunctionSpace_float32)

Bases: object

Object for representing Dirichlet (essential) boundary conditions

dof_indices
property dtype

(self) -> str

property function_space

(self) -> dolfinx::fem::FunctionSpace<float>

property value

(self) -> Union[dolfinx::fem::Function<std::complex<float>, float>, dolfinx::fem::Constant<std::complex<float> >]

class dolfinx.cpp.fem.DirichletBC_float32(self, g: ndarray[dtype=float32, writable=False, order='C'], dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float32)
class dolfinx.cpp.fem.DirichletBC_float32(self, g: dolfinx.cpp.fem.Constant_float32, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float32)
class dolfinx.cpp.fem.DirichletBC_float32(self, g: dolfinx.cpp.fem.Function_float32, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'])
class dolfinx.cpp.fem.DirichletBC_float32(self, g: dolfinx.cpp.fem.Function_float32, dofs: list[ndarray[dtype=int32, writable=False, shape=(*), order='C']], V: dolfinx.cpp.fem.FunctionSpace_float32)

Bases: object

Object for representing Dirichlet (essential) boundary conditions

dof_indices
property dtype

(self) -> str

property function_space

(self) -> dolfinx::fem::FunctionSpace<float>

property value

(self) -> Union[dolfinx::fem::Function<float, float>, dolfinx::fem::Constant<float>]

class dolfinx.cpp.fem.DirichletBC_float64(self, g: ndarray[dtype=float64, writable=False, order='C'], dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float64)
class dolfinx.cpp.fem.DirichletBC_float64(self, g: dolfinx.cpp.fem.Constant_float64, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'], V: dolfinx.cpp.fem.FunctionSpace_float64)
class dolfinx.cpp.fem.DirichletBC_float64(self, g: dolfinx.cpp.fem.Function_float64, dofs: ndarray[dtype=int32, writable=False, shape=(*), order='C'])
class dolfinx.cpp.fem.DirichletBC_float64(self, g: dolfinx.cpp.fem.Function_float64, dofs: list[ndarray[dtype=int32, writable=False, shape=(*), order='C']], V: dolfinx.cpp.fem.FunctionSpace_float64)

Bases: object

Object for representing Dirichlet (essential) boundary conditions

dof_indices
property dtype

(self) -> str

property function_space

(self) -> dolfinx::fem::FunctionSpace<double>

property value

(self) -> Union[dolfinx::fem::Function<double, double>, dolfinx::fem::Constant<double>]

class dolfinx.cpp.fem.DofMap(self, element_dof_layout: dolfinx.cpp.fem.ElementDofLayout, index_map: dolfinx.cpp.common.IndexMap, index_map_bs: int, dofmap: dolfinx.cpp.graph.AdjacencyList_int32, bs: int)

Bases: object

DofMap object

property bs

(self) -> int

cell_dofs
property dof_layout

(self) -> dolfinx.cpp.fem.ElementDofLayout

property index_map

(self) -> dolfinx.cpp.common.IndexMap

property index_map_bs

(self) -> int

map
class dolfinx.cpp.fem.ElementDofLayout(self, block_size: int, endity_dofs: list[list[list[int]]], entity_closure_dofs: list[list[list[int]]], parent_map: list[int], sub_layouts: list[dolfinx.cpp.fem.ElementDofLayout])

Bases: object

Object describing the layout of dofs on a cell

property block_size

(self) -> int

entity_closure_dofs
entity_dofs
property num_dofs

(self) -> int

num_entity_closure_dofs
num_entity_dofs
class dolfinx.cpp.fem.Expression_complex128(self, coefficients: list[dolfinx.cpp.fem.Function_complex128], constants: list[dolfinx.cpp.fem.Constant_complex128], X: ndarray[dtype=float64, writable=False, shape=(*, *), order='C'], fn: int, value_shape: list[int], argument_function_space: dolfinx.cpp.fem.FunctionSpace_float64)

Bases: object

An Expression

X
property dtype

(self) -> str

eval
property value_shape

(self) -> list[int]

property value_size

(self) -> int

class dolfinx.cpp.fem.Expression_complex64(self, coefficients: list[dolfinx.cpp.fem.Function_complex64], constants: list[dolfinx.cpp.fem.Constant_complex64], X: ndarray[dtype=float32, writable=False, shape=(*, *), order='C'], fn: int, value_shape: list[int], argument_function_space: dolfinx.cpp.fem.FunctionSpace_float32)

Bases: object

An Expression

X
property dtype

(self) -> str

eval
property value_shape

(self) -> list[int]

property value_size

(self) -> int

class dolfinx.cpp.fem.Expression_float32(self, coefficients: list[dolfinx.cpp.fem.Function_float32], constants: list[dolfinx.cpp.fem.Constant_float32], X: ndarray[dtype=float32, writable=False, shape=(*, *), order='C'], fn: int, value_shape: list[int], argument_function_space: dolfinx.cpp.fem.FunctionSpace_float32)

Bases: object

An Expression

X
property dtype

(self) -> str

eval
property value_shape

(self) -> list[int]

property value_size

(self) -> int

class dolfinx.cpp.fem.Expression_float64(self, coefficients: list[dolfinx.cpp.fem.Function_float64], constants: list[dolfinx.cpp.fem.Constant_float64], X: ndarray[dtype=float64, writable=False, shape=(*, *), order='C'], fn: int, value_shape: list[int], argument_function_space: dolfinx.cpp.fem.FunctionSpace_float64)

Bases: object

An Expression

X
property dtype

(self) -> str

eval
property value_shape

(self) -> list[int]

property value_size

(self) -> int

class dolfinx.cpp.fem.FiniteElement_float32(self, ufcx_element: int)

Bases: object

Finite element object

T_apply
Tt_apply
Tt_inv_apply
property basix_element

(self) -> basix::FiniteElement<float>

property dtype

(self) -> str

property interpolation_ident

(self) -> bool

interpolation_points
property needs_dof_transformations

(self) -> bool

property num_sub_elements

(self) -> int

signature
property space_dimension

(self) -> int

class dolfinx.cpp.fem.FiniteElement_float64(self, ufcx_element: int)

Bases: object

Finite element object

T_apply
Tt_apply
Tt_inv_apply
property basix_element

(self) -> basix::FiniteElement<double>

property dtype

(self) -> str

property interpolation_ident

(self) -> bool

interpolation_points
property needs_dof_transformations

(self) -> bool

property num_sub_elements

(self) -> int

signature
property space_dimension

(self) -> int

class dolfinx.cpp.fem.Form_complex128(self, spaces: list[dolfinx.cpp.fem.FunctionSpace_float64], integrals: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], coefficients: list[dolfinx.cpp.fem.Function_complex128], constants: list[dolfinx.cpp.fem.Constant_complex128], need_permutation_data: bool, entity_maps: dict[dolfinx.cpp.mesh.Mesh_float64, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float64])
class dolfinx.cpp.fem.Form_complex128(self, form: int, spaces: list[dolfinx.cpp.fem.FunctionSpace_float64], coefficients: list[dolfinx.cpp.fem.Function_complex128], constants: list[dolfinx.cpp.fem.Constant_complex128], subdomains: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], entity_maps: dict[dolfinx.cpp.mesh.Mesh_float64, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float64])

Bases: object

Variational form object

Create a Form from a pointer to a ufcx_form

property coefficients

(self) -> list[dolfinx.cpp.fem.Function_complex128]

domains
property dtype

(self) -> str

property function_spaces

(self) -> list[dolfinx::fem::FunctionSpace<double>]

integral_ids
property integral_types

(self) -> std::set<dolfinx::fem::IntegralType, std::less<dolfinx::fem::IntegralType>, std::allocator<dolfinx::fem::IntegralType> >

property mesh

(self) -> dolfinx.cpp.mesh.Mesh_float64

property needs_facet_permutations

(self) -> bool

property rank

(self) -> int

class dolfinx.cpp.fem.Form_complex64(self, spaces: list[dolfinx.cpp.fem.FunctionSpace_float32], integrals: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], coefficients: list[dolfinx.cpp.fem.Function_complex64], constants: list[dolfinx.cpp.fem.Constant_complex64], need_permutation_data: bool, entity_maps: dict[dolfinx.cpp.mesh.Mesh_float32, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float32])
class dolfinx.cpp.fem.Form_complex64(self, form: int, spaces: list[dolfinx.cpp.fem.FunctionSpace_float32], coefficients: list[dolfinx.cpp.fem.Function_complex64], constants: list[dolfinx.cpp.fem.Constant_complex64], subdomains: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], entity_maps: dict[dolfinx.cpp.mesh.Mesh_float32, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float32])

Bases: object

Variational form object

Create a Form from a pointer to a ufcx_form

property coefficients

(self) -> list[dolfinx.cpp.fem.Function_complex64]

domains
property dtype

(self) -> str

property function_spaces

(self) -> list[dolfinx::fem::FunctionSpace<float>]

integral_ids
property integral_types

(self) -> std::set<dolfinx::fem::IntegralType, std::less<dolfinx::fem::IntegralType>, std::allocator<dolfinx::fem::IntegralType> >

property mesh

(self) -> dolfinx.cpp.mesh.Mesh_float32

property needs_facet_permutations

(self) -> bool

property rank

(self) -> int

class dolfinx.cpp.fem.Form_float32(self, spaces: list[dolfinx.cpp.fem.FunctionSpace_float32], integrals: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], coefficients: list[dolfinx.cpp.fem.Function_float32], constants: list[dolfinx.cpp.fem.Constant_float32], need_permutation_data: bool, entity_maps: dict[dolfinx.cpp.mesh.Mesh_float32, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float32])
class dolfinx.cpp.fem.Form_float32(self, form: int, spaces: list[dolfinx.cpp.fem.FunctionSpace_float32], coefficients: list[dolfinx.cpp.fem.Function_float32], constants: list[dolfinx.cpp.fem.Constant_float32], subdomains: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], entity_maps: dict[dolfinx.cpp.mesh.Mesh_float32, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float32])

Bases: object

Variational form object

Create a Form from a pointer to a ufcx_form

property coefficients

(self) -> list[dolfinx.cpp.fem.Function_float32]

domains
property dtype

(self) -> str

property function_spaces

(self) -> list[dolfinx::fem::FunctionSpace<float>]

integral_ids
property integral_types

(self) -> std::set<dolfinx::fem::IntegralType, std::less<dolfinx::fem::IntegralType>, std::allocator<dolfinx::fem::IntegralType> >

property mesh

(self) -> dolfinx.cpp.mesh.Mesh_float32

property needs_facet_permutations

(self) -> bool

property rank

(self) -> int

class dolfinx.cpp.fem.Form_float64(self, spaces: list[dolfinx.cpp.fem.FunctionSpace_float64], integrals: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], coefficients: list[dolfinx.cpp.fem.Function_float64], constants: list[dolfinx.cpp.fem.Constant_float64], need_permutation_data: bool, entity_maps: dict[dolfinx.cpp.mesh.Mesh_float64, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float64])
class dolfinx.cpp.fem.Form_float64(self, form: int, spaces: list[dolfinx.cpp.fem.FunctionSpace_float64], coefficients: list[dolfinx.cpp.fem.Function_float64], constants: list[dolfinx.cpp.fem.Constant_float64], subdomains: dict[dolfinx.cpp.fem.IntegralType, list[tuple[int, ndarray[dtype=int32, writable=False, shape=(*), order='C']]]], entity_maps: dict[dolfinx.cpp.mesh.Mesh_float64, ndarray[dtype=int32, writable=False, shape=(*), order='C']], mesh: Optional[dolfinx.cpp.mesh.Mesh_float64])

Bases: object

Variational form object

Create a Form from a pointer to a ufcx_form

property coefficients

(self) -> list[dolfinx.cpp.fem.Function_float64]

domains
property dtype

(self) -> str

property function_spaces

(self) -> list[dolfinx::fem::FunctionSpace<double>]

integral_ids
property integral_types

(self) -> std::set<dolfinx::fem::IntegralType, std::less<dolfinx::fem::IntegralType>, std::allocator<dolfinx::fem::IntegralType> >

property mesh

(self) -> dolfinx.cpp.mesh.Mesh_float64

property needs_facet_permutations

(self) -> bool

property rank

(self) -> int

class dolfinx.cpp.fem.FunctionSpace_float32(self, mesh: dolfinx.cpp.mesh.Mesh_float32, element: dolfinx.cpp.fem.FiniteElement_float32, dofmap: dolfinx.cpp.fem.DofMap, value_shape: list[int])

Bases: object

Finite element function space

collapse
component
contains
property dofmap

(self) -> dolfinx.cpp.fem.DofMap

property element

(self) -> dolfinx::fem::FiniteElement<float>

property mesh

(self) -> dolfinx.cpp.mesh.Mesh_float32

sub
tabulate_dof_coordinates
property value_shape

(self) -> numpy.ndarray[dtype=uint64, writable=False, ]

class dolfinx.cpp.fem.FunctionSpace_float64(self, mesh: dolfinx.cpp.mesh.Mesh_float64, element: dolfinx.cpp.fem.FiniteElement_float64, dofmap: dolfinx.cpp.fem.DofMap, value_shape: list[int])

Bases: object

Finite element function space

collapse
component
contains
property dofmap

(self) -> dolfinx.cpp.fem.DofMap

property element

(self) -> dolfinx::fem::FiniteElement<double>

property mesh

(self) -> dolfinx.cpp.mesh.Mesh_float64

sub
tabulate_dof_coordinates
property value_shape

(self) -> numpy.ndarray[dtype=uint64, writable=False, ]

class dolfinx.cpp.fem.Function_complex128(self, arg: dolfinx.cpp.fem.FunctionSpace_float64, /)
class dolfinx.cpp.fem.Function_complex128(self, arg0: dolfinx.cpp.fem.FunctionSpace_float64, arg1: dolfinx.cpp.la.Vector_complex128, /)

Bases: object

A finite element function

Create a function on the given function space

collapse

Collapse sub-function view

eval

Evaluate Function

property function_space

(self) -> dolfinx::fem::FunctionSpace<double>

interpolate

Overloaded function.

  1. interpolate(self, f: ndarray[dtype=complex128, writable=False, shape=(*), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, f: ndarray[dtype=complex128, writable=False, shape=(*, *), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, u: dolfinx.cpp.fem.Function_complex128, cells: ndarray[dtype=int32, writable=False, shape=(*), order='C'], cell_map: ndarray[dtype=int32, writable=False, shape=(*), order='C'], nmm_interpolation_data: tuple[ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=float64, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C']]) -> None

Interpolate a finite element function

  1. interpolate(self, expr: dolfinx.cpp.fem.Expression_complex128, cells: ndarray[dtype=int32, writable=False, order='C'], expr_mesh: dolfinx.cpp.mesh.Mesh_float64, cell_map: ndarray[dtype=int32, writable=False, order='C']) -> None

Interpolate an Expression on a set of cells

interpolate_ptr

Interpolate using a pointer to an expression with a C signature

property name

(self) -> str

sub

Return sub-function (view into parent Function

property x

Return the vector associated with the finite element Function

class dolfinx.cpp.fem.Function_complex64(self, arg: dolfinx.cpp.fem.FunctionSpace_float32, /)
class dolfinx.cpp.fem.Function_complex64(self, arg0: dolfinx.cpp.fem.FunctionSpace_float32, arg1: dolfinx.cpp.la.Vector_complex64, /)

Bases: object

A finite element function

Create a function on the given function space

collapse

Collapse sub-function view

eval

Evaluate Function

property function_space

(self) -> dolfinx::fem::FunctionSpace<float>

interpolate

Overloaded function.

  1. interpolate(self, f: ndarray[dtype=complex64, writable=False, shape=(*), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, f: ndarray[dtype=complex64, writable=False, shape=(*, *), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, u: dolfinx.cpp.fem.Function_complex64, cells: ndarray[dtype=int32, writable=False, shape=(*), order='C'], cell_map: ndarray[dtype=int32, writable=False, shape=(*), order='C'], nmm_interpolation_data: tuple[ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=float32, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C']]) -> None

Interpolate a finite element function

  1. interpolate(self, expr: dolfinx.cpp.fem.Expression_complex64, cells: ndarray[dtype=int32, writable=False, order='C'], expr_mesh: dolfinx.cpp.mesh.Mesh_float32, cell_map: ndarray[dtype=int32, writable=False, order='C']) -> None

Interpolate an Expression on a set of cells

interpolate_ptr

Interpolate using a pointer to an expression with a C signature

property name

(self) -> str

sub

Return sub-function (view into parent Function

property x

Return the vector associated with the finite element Function

class dolfinx.cpp.fem.Function_float32(self, arg: dolfinx.cpp.fem.FunctionSpace_float32, /)
class dolfinx.cpp.fem.Function_float32(self, arg0: dolfinx.cpp.fem.FunctionSpace_float32, arg1: dolfinx.cpp.la.Vector_float32, /)

Bases: object

A finite element function

Create a function on the given function space

collapse

Collapse sub-function view

eval

Evaluate Function

property function_space

(self) -> dolfinx::fem::FunctionSpace<float>

interpolate

Overloaded function.

  1. interpolate(self, f: ndarray[dtype=float32, writable=False, shape=(*), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, f: ndarray[dtype=float32, writable=False, shape=(*, *), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, u: dolfinx.cpp.fem.Function_float32, cells: ndarray[dtype=int32, writable=False, shape=(*), order='C'], cell_map: ndarray[dtype=int32, writable=False, shape=(*), order='C'], nmm_interpolation_data: tuple[ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=float32, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C']]) -> None

Interpolate a finite element function

  1. interpolate(self, expr: dolfinx.cpp.fem.Expression_float32, cells: ndarray[dtype=int32, writable=False, order='C'], expr_mesh: dolfinx.cpp.mesh.Mesh_float32, cell_map: ndarray[dtype=int32, writable=False, order='C']) -> None

Interpolate an Expression on a set of cells

interpolate_ptr

Interpolate using a pointer to an expression with a C signature

property name

(self) -> str

sub

Return sub-function (view into parent Function

property x

Return the vector associated with the finite element Function

class dolfinx.cpp.fem.Function_float64(self, arg: dolfinx.cpp.fem.FunctionSpace_float64, /)
class dolfinx.cpp.fem.Function_float64(self, arg0: dolfinx.cpp.fem.FunctionSpace_float64, arg1: dolfinx.cpp.la.Vector_float64, /)

Bases: object

A finite element function

Create a function on the given function space

collapse

Collapse sub-function view

eval

Evaluate Function

property function_space

(self) -> dolfinx::fem::FunctionSpace<double>

interpolate

Overloaded function.

  1. interpolate(self, f: ndarray[dtype=float64, writable=False, shape=(*), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, f: ndarray[dtype=float64, writable=False, shape=(*, *), order='C'], cells: ndarray[dtype=int32, writable=False, shape=(*), order='C']) -> None

Interpolate an expression function

  1. interpolate(self, u: dolfinx.cpp.fem.Function_float64, cells: ndarray[dtype=int32, writable=False, shape=(*), order='C'], cell_map: ndarray[dtype=int32, writable=False, shape=(*), order='C'], nmm_interpolation_data: tuple[ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C'], ndarray[dtype=float64, writable=False, shape=(*), order='C'], ndarray[dtype=int32, writable=False, shape=(*), order='C']]) -> None

Interpolate a finite element function

  1. interpolate(self, expr: dolfinx.cpp.fem.Expression_float64, cells: ndarray[dtype=int32, writable=False, order='C'], expr_mesh: dolfinx.cpp.mesh.Mesh_float64, cell_map: ndarray[dtype=int32, writable=False, order='C']) -> None

Interpolate an Expression on a set of cells

interpolate_ptr

Interpolate using a pointer to an expression with a C signature

property name

(self) -> str

sub

Return sub-function (view into parent Function

property x

Return the vector associated with the finite element Function

class dolfinx.cpp.fem.IntegralType

Bases: object

cell = dolfinx.cpp.fem.IntegralType.cell
exterior_facet = dolfinx.cpp.fem.IntegralType.exterior_facet
interior_facet = dolfinx.cpp.fem.IntegralType.interior_facet
vertex = dolfinx.cpp.fem.IntegralType.vertex