dolfin.cpp.fem

FEM module

Functions

add_diagonal(arg0, arg1, arg2, arg3)

apply_lifting(*args, **kwargs)

Overloaded function.

assemble_matrix(*args, **kwargs)

Overloaded function.

assemble_scalar(arg0)

Assemble functional over mesh

assemble_vector(*args, **kwargs)

Overloaded function.

bcs_cols(arg0, arg1)

bcs_rows(arg0, arg1)

block_function_spaces(arg0)

build_dofmap(arg0, arg1)

Build and dofmap on a mesh.

create_dofmap(arg0, arg1)

Create DOLFIN DofMap object from a ufc dofmap.

create_element_dof_layout(arg0, arg1, arg2)

Create ElementDofLayout object from a ufc dofmap.

create_form(arg0, arg1)

Create DOLFIN form from a ufc form.

create_matrix(arg0)

Create a PETSc Mat for bilinear form.

create_matrix_block(arg0)

Create monolithic sparse matrix for stacked bilinear forms.

create_matrix_nest(arg0)

Create nested sparse matrix for bilinear forms.

create_vector_block(arg0)

Create a monolithic vector for multiple (stacked) linear forms.

create_vector_nest(arg0)

Create nested vector for multiple (stacked) linear forms.

make_coordinate_mapping(arg0)

Create a CoordinateElement object from a pointer to a ufc_coordinate_map.

make_ufc_dofmap(arg0)

Create a ufc_dofmap object from a pointer.

make_ufc_finite_element(arg0)

Create a ufc_finite_element object from a pointer.

make_ufc_form(arg0)

Create a ufc_form object from a pointer.

set_bc(*args, **kwargs)

Overloaded function.

Classes

CoordinateElement

Coordinate mapping object

DirichletBC(*args, **kwargs)

Object for representing Dirichlet (essential) boundary conditions

DiscreteOperators

DofMap

DofMap object

ElementDofLayout

Object describing the layout of dofs on a cell

FiniteElement(self, arg0)

Finite element object

Form(self, arg0)

Variational form object

FormIntegrals

Holder for integral kernels and domains

PETScDMCollection(self, arg0)

ufc_coordinate_mapping

UFC coordinate_mapping object

ufc_dofmap

UFC dofmap object

ufc_finite_element

UFC finite element object

ufc_form

UFC form object

class dolfin.cpp.fem.CoordinateElement

Bases: pybind11_builtins.pybind11_object

Coordinate mapping object

push_forward(self: dolfin.cpp.fem.CoordinateElement, arg0: numpy.ndarray[float64[m, n], flags.writeable, flags.c_contiguous], arg1: numpy.ndarray[float64[m, n], flags.c_contiguous], arg2: numpy.ndarray[float64[m, n], flags.c_contiguous]) → None
class dolfin.cpp.fem.DirichletBC(*args, **kwargs)

Bases: pybind11_builtins.pybind11_object

Object for representing Dirichlet (essential) boundary conditions

Overloaded function.

  1. __init__(self: dolfin.cpp.fem.DirichletBC, V: dolfin.cpp.function.FunctionSpace, g: dolfin.cpp.function.Function, mark: Callable[[numpy.ndarray[float64[3, n], flags.c_contiguous]], numpy.ndarray[bool[m, 1]]], method: dolfin.cpp.fem.DirichletBC.Method) -> None

  2. __init__(self: dolfin.cpp.fem.DirichletBC, V: dolfin.cpp.function.FunctionSpace, g: dolfin.cpp.function.Function, facets: List[int], method: dolfin.cpp.fem.DirichletBC.Method) -> None

class Method(self: dolfin.cpp.fem.DirichletBC.Method, arg0: int) → None

Bases: pybind11_builtins.pybind11_object

Members:

topological

geometric

pointwise

property name

handle) -> str

Type

(self

class dolfin.cpp.fem.DofMap

Bases: pybind11_builtins.pybind11_object

DofMap object

cell_dofs(self: dolfin.cpp.fem.DofMap, arg0: int) → numpy.ndarray[int32[m, 1]]
dof_array(self: dolfin.cpp.fem.DofMap) → numpy.ndarray[int32[m, 1]]
dofs(self: dolfin.cpp.fem.DofMap, arg0: dolfin.cpp.mesh.Mesh, arg1: int) → numpy.ndarray[int32[m, 1]]
set(self: dolfin.cpp.fem.DofMap, arg0: numpy.ndarray[float64[m, 1], flags.writeable], arg1: float) → None
class dolfin.cpp.fem.ElementDofLayout

Bases: pybind11_builtins.pybind11_object

Object describing the layout of dofs on a cell

entity_closure_dofs(self: dolfin.cpp.fem.ElementDofLayout, arg0: int, arg1: int) → numpy.ndarray[int32[m, 1]]
entity_dofs(self: dolfin.cpp.fem.ElementDofLayout, arg0: int, arg1: int) → numpy.ndarray[int32[m, 1]]
num_entity_closure_dofs(self: dolfin.cpp.fem.ElementDofLayout, arg0: int) → int
num_entity_dofs(self: dolfin.cpp.fem.ElementDofLayout, arg0: int) → int
class dolfin.cpp.fem.FiniteElement(self: dolfin.cpp.fem.FiniteElement, arg0: dolfin.cpp.fem.ufc_finite_element) → None

Bases: pybind11_builtins.pybind11_object

Finite element object

dof_reference_coordinates(self: dolfin.cpp.fem.FiniteElement) → numpy.ndarray[float64[m, n]]
num_sub_elements(self: dolfin.cpp.fem.FiniteElement) → int
signature(self: dolfin.cpp.fem.FiniteElement) → str
space_dimension(self: dolfin.cpp.fem.FiniteElement) → int
value_dimension(self: dolfin.cpp.fem.FiniteElement, arg0: int) → int
class dolfin.cpp.fem.Form(self: dolfin.cpp.fem.Form, arg0: List[dolfin.cpp.function.FunctionSpace]) → None

Bases: pybind11_builtins.pybind11_object

Variational form object

coordinate_mapping(self: dolfin.cpp.fem.Form) → dolfin.cpp.fem.CoordinateElement
function_space(self: dolfin.cpp.fem.Form, arg0: int) → dolfin.cpp.function.FunctionSpace
mesh(self: dolfin.cpp.fem.Form) → dolfin.cpp.mesh.Mesh
num_coefficients(self: dolfin.cpp.fem.Form) → int

Return number of coefficients in form

original_coefficient_position(self: dolfin.cpp.fem.Form, arg0: int) → int
set_cell_domains(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.MeshFunctionSizet) → None
set_coefficient(self: dolfin.cpp.fem.Form, arg0: int, arg1: dolfin.cpp.function.Function) → None
set_constants(self: dolfin.cpp.fem.Form, arg0: List[dolfin.cpp.function.Constant]) → None
set_exterior_facet_domains(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.MeshFunctionSizet) → None
set_interior_facet_domains(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.MeshFunctionSizet) → None
set_mesh(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.Mesh) → None
set_tabulate_tensor(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.fem.FormIntegrals.Type, arg1: int, arg2: int) → None
set_vertex_domains(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.MeshFunctionSizet) → None
class dolfin.cpp.fem.FormIntegrals

Bases: pybind11_builtins.pybind11_object

Holder for integral kernels and domains

class Type(self: dolfin.cpp.fem.FormIntegrals.Type, arg0: int) → None

Bases: pybind11_builtins.pybind11_object

Members:

cell

exterior_facet

interior_facet

property name

handle) -> str

Type

(self

dolfin.cpp.fem.add_diagonal(arg0: mat, arg1: dolfin.cpp.function.FunctionSpace, arg2: List[dolfin.cpp.fem.DirichletBC], arg3: float) → None
dolfin.cpp.fem.apply_lifting(*args, **kwargs)

Overloaded function.

  1. apply_lifting(arg0: vec, arg1: List[dolfinx::fem::Form], arg2: List[List[dolfin.cpp.fem.DirichletBC]], arg3: List[vec], arg4: float) -> None

Modify vector for lifted boundary conditions

  1. apply_lifting(arg0: numpy.ndarray[float64[m, 1], flags.writeable], arg1: List[dolfinx::fem::Form], arg2: List[List[dolfin.cpp.fem.DirichletBC]], arg3: List[numpy.ndarray[float64[m, 1]]], arg4: float) -> None

Modify vector for lifted boundary conditions

dolfin.cpp.fem.assemble_matrix(*args, **kwargs)

Overloaded function.

  1. assemble_matrix(arg0: mat, arg1: dolfinx::fem::Form, arg2: List[dolfin.cpp.fem.DirichletBC]) -> None

  2. assemble_matrix(arg0: mat, arg1: dolfinx::fem::Form, arg2: List[bool], arg3: List[bool]) -> None

dolfin.cpp.fem.assemble_scalar(arg0: dolfinx::fem::Form) → float

Assemble functional over mesh

dolfin.cpp.fem.assemble_vector(*args, **kwargs)

Overloaded function.

  1. assemble_vector(b: vec, L: dolfinx::fem::Form) -> None

Assemble linear form into an existing vector

  1. assemble_vector(b: numpy.ndarray[float64[m, 1], flags.writeable], L: dolfinx::fem::Form) -> None

Assemble linear form into an existing Eigen vector

dolfin.cpp.fem.bcs_cols(arg0: List[List[dolfinx::fem::Form]], arg1: List[dolfin.cpp.fem.DirichletBC]) → List[List[List[dolfin.cpp.fem.DirichletBC]]]
dolfin.cpp.fem.bcs_rows(arg0: List[dolfinx::fem::Form], arg1: List[dolfin.cpp.fem.DirichletBC]) → List[List[dolfin.cpp.fem.DirichletBC]]
dolfin.cpp.fem.block_function_spaces(arg0: List[List[dolfinx::fem::Form]]) → List[List[dolfin.cpp.function.FunctionSpace][2]]
dolfin.cpp.fem.build_dofmap(arg0: dolfin.cpp.mesh.Mesh, arg1: dolfinx::fem::ElementDofLayout) → dolfinx::fem::DofMap

Build and dofmap on a mesh.

dolfin.cpp.fem.create_dofmap(arg0: dolfin.cpp.fem.ufc_dofmap, arg1: dolfin.cpp.mesh.Mesh) → dolfinx::fem::DofMap

Create DOLFIN DofMap object from a ufc dofmap.

dolfin.cpp.fem.create_element_dof_layout(arg0: dolfin.cpp.fem.ufc_dofmap, arg1: dolfin.cpp.mesh.CellType, arg2: List[int]) → dolfinx::fem::ElementDofLayout

Create ElementDofLayout object from a ufc dofmap.

dolfin.cpp.fem.create_form(arg0: dolfin.cpp.fem.ufc_form, arg1: List[dolfin.cpp.function.FunctionSpace]) → dolfinx::fem::Form

Create DOLFIN form from a ufc form.

dolfin.cpp.fem.create_matrix(arg0: dolfinx::fem::Form) → mat

Create a PETSc Mat for bilinear form.

dolfin.cpp.fem.create_matrix_block(arg0: List[List[dolfinx::fem::Form]]) → mat

Create monolithic sparse matrix for stacked bilinear forms.

dolfin.cpp.fem.create_matrix_nest(arg0: List[List[dolfinx::fem::Form]]) → mat

Create nested sparse matrix for bilinear forms.

dolfin.cpp.fem.create_vector_block(arg0: List[dolfin.cpp.common.IndexMap]) → vec

Create a monolithic vector for multiple (stacked) linear forms.

dolfin.cpp.fem.create_vector_nest(arg0: List[dolfin.cpp.common.IndexMap]) → vec

Create nested vector for multiple (stacked) linear forms.

dolfin.cpp.fem.make_coordinate_mapping(arg0: int) → dolfinx::fem::CoordinateElement

Create a CoordinateElement object from a pointer to a ufc_coordinate_map.

dolfin.cpp.fem.make_ufc_dofmap(arg0: int) → dolfin.cpp.fem.ufc_dofmap

Create a ufc_dofmap object from a pointer.

dolfin.cpp.fem.make_ufc_finite_element(arg0: int) → dolfin.cpp.fem.ufc_finite_element

Create a ufc_finite_element object from a pointer.

dolfin.cpp.fem.make_ufc_form(arg0: int) → dolfin.cpp.fem.ufc_form

Create a ufc_form object from a pointer.

dolfin.cpp.fem.set_bc(*args, **kwargs)

Overloaded function.

  1. set_bc(arg0: vec, arg1: List[dolfin.cpp.fem.DirichletBC], arg2: vec, arg3: float) -> None

Insert boundary condition values into vector

  1. set_bc(b: numpy.ndarray[float64[m, 1], flags.writeable], bcs: List[dolfin.cpp.fem.DirichletBC], x0: numpy.ndarray[float64] = None, scale: float = 1.0) -> None

class dolfin.cpp.fem.ufc_coordinate_mapping

Bases: pybind11_builtins.pybind11_object

UFC coordinate_mapping object

class dolfin.cpp.fem.ufc_dofmap

Bases: pybind11_builtins.pybind11_object

UFC dofmap object

class dolfin.cpp.fem.ufc_finite_element

Bases: pybind11_builtins.pybind11_object

UFC finite element object

class dolfin.cpp.fem.ufc_form

Bases: pybind11_builtins.pybind11_object

UFC form object