dolfin.cpp.fem

FEM module

Functions

assemble(*args, **kwargs) Overloaded function.
assemble_local((arg0: dolfin.cpp.fem.Form, …)
assemble_system(*args, **kwargs) Overloaded function.
create_mesh(*args, **kwargs) Overloaded function.
dof_to_vertex_map(*args, **kwargs) Overloaded function.
get_coordinates(*args, **kwargs) Overloaded function.
make_ufc_dofmap(…)
make_ufc_finite_element(…)
make_ufc_form(…)
set_coordinates(*args, **kwargs) Overloaded function.
vertex_to_dof_map(*args, **kwargs) Overloaded function.

Classes

Assembler DOLFIN Assembler object
AssemblerBase
DirichletBC DOLFIN DirichletBC object
DiscreteOperators
DofMap DOLFIN DofMap object
FiniteElement DOLFIN FiniteElement object
Form DOLFIN Form object
GenericDofMap DOLFIN DofMap object
LinearVariationalProblem
LinearVariationalSolver
LocalSolver
NonlinearVariationalProblem
NonlinearVariationalSolver
PETScDMCollection
PointSource
SparsityPatternBuilder
SystemAssembler DOLFIN SystemAssembler object
ufc_dofmap UFC dofmap object
ufc_finite_element UFC finite element object
ufc_form UFC form object
class dolfin.cpp.fem.Assembler

Bases: dolfin.cpp.fem.AssemblerBase

DOLFIN Assembler object

assemble(self: dolfin.cpp.fem.Assembler, arg0: dolfin::GenericTensor, arg1: dolfin::Form) → None
class dolfin.cpp.fem.DirichletBC

Bases: dolfin.cpp.common.Variable

DOLFIN DirichletBC object

apply(*args, **kwargs)

Overloaded function.

  1. apply(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin::GenericVector) -> None
  2. apply(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin::GenericMatrix) -> None
  3. apply(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin::GenericMatrix, arg1: dolfin::GenericVector) -> None
  4. apply(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin::GenericVector, arg1: dolfin::GenericVector) -> None
  5. apply(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin::GenericMatrix, arg1: dolfin::GenericVector, arg2: dolfin::GenericVector) -> None
function_space(self: dolfin.cpp.fem.DirichletBC) → dolfin.cpp.function.FunctionSpace
get_boundary_values(self: dolfin.cpp.fem.DirichletBC) → Dict[int, float]
homogenize(self: dolfin.cpp.fem.DirichletBC) → None
method(self: dolfin.cpp.fem.DirichletBC) → str
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin.cpp.function.GenericFunction) -> None
  2. set_value(self: dolfin.cpp.fem.DirichletBC, arg0: object) -> None
user_subdomain(self: dolfin.cpp.fem.DirichletBC) → dolfin.cpp.mesh.SubDomain
value(self: dolfin.cpp.fem.DirichletBC) → dolfin.cpp.function.GenericFunction
zero(self: dolfin.cpp.fem.DirichletBC, arg0: dolfin::GenericMatrix) → None
zero_columns(self: dolfin.cpp.fem.DirichletBC, A: dolfin::GenericMatrix, b: dolfin::GenericVector, diagonal_value: float = 0.0) → None
class dolfin.cpp.fem.DofMap

Bases: dolfin.cpp.fem.GenericDofMap

DOLFIN DofMap object

cell_dofs(self: dolfin.cpp.fem.DofMap, arg0: int) → numpy.ndarray[int32[m, 1]]
ownership_range(self: dolfin.cpp.fem.DofMap) → Tuple[int, int]
class dolfin.cpp.fem.FiniteElement

Bases: pybind11_builtins.pybind11_object

DOLFIN FiniteElement object

evaluate_basis(self: dolfin.cpp.fem.FiniteElement, arg0: int, arg1: numpy.ndarray[float64], arg2: numpy.ndarray[float64], arg3: int) → numpy.ndarray[float64]
evaluate_basis_derivatives(self: dolfin.cpp.fem.FiniteElement, arg0: int, arg1: int, arg2: numpy.ndarray[float64], arg3: numpy.ndarray[float64], arg4: int) → numpy.ndarray[float64]
evaluate_dofs(self: dolfin.cpp.fem.FiniteElement, arg0: object, arg1: numpy.ndarray[float64], arg2: int, arg3: dolfin.cpp.mesh.Cell) → numpy.ndarray[float64]

Evaluate degrees of freedom on element for a given function

geometric_dimension(self: dolfin.cpp.fem.FiniteElement) → int
num_sub_elements(self: dolfin.cpp.fem.FiniteElement) → int
signature(self: dolfin.cpp.fem.FiniteElement) → str
space_dimension(self: dolfin.cpp.fem.FiniteElement) → int
tabulate_dof_coordinates(self: dolfin.cpp.fem.FiniteElement, arg0: dolfin.cpp.mesh.Cell) → numpy.ndarray[float64[m, n]]

Tabulate coordinates of dofs on cell

value_dimension(self: dolfin.cpp.fem.FiniteElement, arg0: int) → int
class dolfin.cpp.fem.Form

Bases: pybind11_builtins.pybind11_object

DOLFIN Form object

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
rank(self: dolfin.cpp.fem.Form) → int
set_cell_domains(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.MeshFunctionSizet) → None
set_coefficient(*args, **kwargs)

Overloaded function.

  1. set_coefficient(self: dolfin.cpp.fem.Form, arg0: int, arg1: dolfin.cpp.function.GenericFunction) -> None

Doc

  1. set_coefficient(self: dolfin.cpp.fem.Form, arg0: str, arg1: dolfin.cpp.function.GenericFunction) -> None

Doc

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_vertex_domains(self: dolfin.cpp.fem.Form, arg0: dolfin.cpp.mesh.MeshFunctionSizet) → None
class dolfin.cpp.fem.GenericDofMap

Bases: dolfin.cpp.common.Variable

DOLFIN DofMap object

block_size(self: dolfin.cpp.fem.GenericDofMap) → int
cell_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: int) → numpy.ndarray[int32[m, 1]]
clear_sub_map_data(self: dolfin.cpp.fem.GenericDofMap) → None
dofs(*args, **kwargs)

Overloaded function.

  1. dofs(self: dolfin.cpp.fem.GenericDofMap) -> List[int]
  2. dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: dolfin.cpp.mesh.Mesh, arg1: int) -> List[int]
entity_closure_dofs(*args, **kwargs)

Overloaded function.

  1. entity_closure_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: dolfin.cpp.mesh.Mesh, arg1: int) -> List[int]
  2. entity_closure_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: dolfin.cpp.mesh.Mesh, arg1: int, arg2: List[int]) -> List[int]
entity_dofs(*args, **kwargs)

Overloaded function.

  1. entity_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: dolfin.cpp.mesh.Mesh, arg1: int) -> List[int]
  2. entity_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: dolfin.cpp.mesh.Mesh, arg1: int, arg2: List[int]) -> List[int]
global_dimension(self: dolfin.cpp.fem.GenericDofMap) → int

The dimension of the global finite element function space

index_map(self: dolfin.cpp.fem.GenericDofMap) → dolfin::IndexMap
local_to_global_index(self: dolfin.cpp.fem.GenericDofMap, arg0: int) → int
local_to_global_unowned(self: dolfin.cpp.fem.GenericDofMap) → numpy.ndarray[uint64[m, 1]]

Return view into unowned part of local-to-global map

neighbours(self: dolfin.cpp.fem.GenericDofMap) → Set[int]
num_entity_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: int) → int
off_process_owner(self: dolfin.cpp.fem.GenericDofMap) → List[int]
set(self: dolfin.cpp.fem.GenericDofMap, arg0: dolfin::GenericVector, arg1: float) → None
shared_nodes(self: dolfin.cpp.fem.GenericDofMap) → Dict[int, List[int]]
tabulate_entity_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: int, arg1: int) → numpy.ndarray[uint64]
tabulate_local_to_global_dofs(*args, **kwargs)

Overloaded function.

  1. tabulate_local_to_global_dofs(self: dolfin.cpp.fem.GenericDofMap, arg0: List[int]) -> None
  2. tabulate_local_to_global_dofs(self: dolfin.cpp.fem.GenericDofMap) -> numpy.ndarray[uint64]
class dolfin.cpp.fem.SystemAssembler

Bases: dolfin.cpp.fem.AssemblerBase

DOLFIN SystemAssembler object

assemble(*args, **kwargs)

Overloaded function.

  1. assemble(self: dolfin.cpp.fem.SystemAssembler, arg0: dolfin::GenericMatrix, arg1: dolfin::GenericVector) -> None
  2. assemble(self: dolfin.cpp.fem.SystemAssembler, arg0: dolfin::GenericMatrix) -> None
  3. assemble(self: dolfin.cpp.fem.SystemAssembler, arg0: dolfin::GenericVector) -> None
  4. assemble(self: dolfin.cpp.fem.SystemAssembler, arg0: dolfin::GenericMatrix, arg1: dolfin::GenericVector, arg2: dolfin::GenericVector) -> None
  5. assemble(self: dolfin.cpp.fem.SystemAssembler, arg0: dolfin::GenericVector, arg1: dolfin::GenericVector) -> None
dolfin.cpp.fem.assemble(*args, **kwargs)

Overloaded function.

  1. assemble(arg0: dolfin::GenericTensor, arg1: dolfin.cpp.fem.Form) -> None
  2. assemble(arg0: dolfin.cpp.fem.Form) -> float
dolfin.cpp.fem.assemble_local(arg0: dolfin.cpp.fem.Form, arg1: dolfin.cpp.mesh.Cell) → numpy.ndarray[float64[m, n]]
dolfin.cpp.fem.assemble_system(*args, **kwargs)

Overloaded function.

  1. assemble_system(arg0: dolfin::GenericMatrix, arg1: dolfin::GenericVector, arg2: dolfin.cpp.fem.Form, arg3: dolfin.cpp.fem.Form, arg4: List[dolfin.cpp.fem.DirichletBC]) -> None
  2. assemble_system(arg0: dolfin::GenericMatrix, arg1: dolfin::GenericVector, arg2: dolfin.cpp.fem.Form, arg3: dolfin.cpp.fem.Form, arg4: List[dolfin.cpp.fem.DirichletBC], arg5: dolfin::GenericVector) -> None
dolfin.cpp.fem.create_mesh(*args, **kwargs)

Overloaded function.

  1. create_mesh(arg0: dolfin.cpp.function.Function) -> dolfin.cpp.mesh.Mesh
  2. create_mesh(arg0: object) -> dolfin.cpp.mesh.Mesh
dolfin.cpp.fem.dof_to_vertex_map(*args, **kwargs)

Overloaded function.

  1. dof_to_vertex_map(arg0: dolfin.cpp.function.FunctionSpace) -> List[int]
  2. dof_to_vertex_map(arg0: object) -> numpy.ndarray[uint64]
dolfin.cpp.fem.get_coordinates(*args, **kwargs)

Overloaded function.

  1. get_coordinates(arg0: dolfin.cpp.function.Function, arg1: dolfin.cpp.mesh.MeshGeometry) -> None
  2. get_coordinates(arg0: object, arg1: dolfin.cpp.mesh.MeshGeometry) -> None
dolfin.cpp.fem.make_ufc_dofmap(arg0: int) → dolfin.cpp.fem.ufc_dofmap
dolfin.cpp.fem.make_ufc_finite_element(arg0: int) → dolfin.cpp.fem.ufc_finite_element
dolfin.cpp.fem.make_ufc_form(arg0: int) → dolfin.cpp.fem.ufc_form
dolfin.cpp.fem.set_coordinates(*args, **kwargs)

Overloaded function.

  1. set_coordinates(arg0: dolfin.cpp.mesh.MeshGeometry, arg1: dolfin.cpp.function.Function) -> None
  2. set_coordinates(arg0: dolfin.cpp.mesh.MeshGeometry, arg1: object) -> None
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

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

Overloaded function.

  1. vertex_to_dof_map(arg0: dolfin.cpp.function.FunctionSpace) -> numpy.ndarray[int32]
  2. vertex_to_dof_map(arg0: object) -> numpy.ndarray[int32]