dolfin.cpp.mesh

Mesh library module

Functions

cells(*args, **kwargs) Overloaded function.
edges(*args, **kwargs) Overloaded function.
entities(*args, **kwargs) Overloaded function.
faces(*args, **kwargs) Overloaded function.
facets(*args, **kwargs) Overloaded function.
make_dolfin_subdomain(arg0)
vertices(*args, **kwargs) Overloaded function.

Classes

BoundaryMesh DOLFIN BoundaryMesh object
Cell DOLFIN Cell object
CellIterator DOLFIN CellIterator object
CellType
DomainBoundary
Edge DOLFIN Edge object
EdgeIterator DOLFIN EdgeIterator object
Face DOLFIN Face object
FaceIterator DOLFIN FaceIterator object
Facet DOLFIN Facet object
FacetIterator DOLFIN FacetIterator object
Mesh DOLFIN Mesh object
MeshColoring
MeshConnectivity DOLFIN MeshConnectivity object
MeshData Mesh data object
MeshDomains Mesh domains object
MeshEditor DOLFIN MeshEditor object
MeshEntity DOLFIN MeshEntity object
MeshEntityIterator DOLFIN MeshEntityIterator object
MeshFunctionBool DOLFIN MeshFunction object
MeshFunctionDouble DOLFIN MeshFunction object
MeshFunctionInt DOLFIN MeshFunction object
MeshFunctionSizet DOLFIN MeshFunction object
MeshGeometry DOLFIN MeshGeometry object
MeshQuality DOLFIN MeshQuality class
MeshTopology DOLFIN MeshTopology object
MeshTransformation
MeshValueCollection_bool DOLFIN MeshValueCollection object
MeshValueCollection_double DOLFIN MeshValueCollection object
MeshValueCollection_int DOLFIN MeshValueCollection object
MeshValueCollection_sizet DOLFIN MeshValueCollection object
MultiMesh
PeriodicBoundaryComputation
SubDomain DOLFIN SubDomain object
SubMesh DOLFIN SubMesh
SubsetIterator
Vertex DOLFIN Vertex object
VertexIterator DOLFIN VertexIterator object
class dolfin.cpp.mesh.BoundaryMesh

Bases: dolfin.cpp.mesh.Mesh

DOLFIN BoundaryMesh object

entity_map(*args, **kwargs)

Overloaded function.

  1. entity_map(self: dolfin.cpp.mesh.BoundaryMesh, arg0: int) -> dolfin::MeshFunction<unsigned long>
  2. entity_map(self: dolfin.cpp.mesh.BoundaryMesh, arg0: int) -> dolfin::MeshFunction<unsigned long>
class dolfin.cpp.mesh.Cell

Bases: dolfin.cpp.mesh.MeshEntity

DOLFIN Cell object

cell_normal(self: dolfin.cpp.mesh.Cell) → dolfin::Point
circumradius(self: dolfin.cpp.mesh.Cell) → float
collides(*args, **kwargs)

Overloaded function.

  1. collides(self: dolfin.cpp.mesh.Cell, arg0: dolfin::Point) -> bool
  2. collides(self: dolfin.cpp.mesh.Cell, arg0: dolfin.cpp.mesh.MeshEntity) -> bool
contains(self: dolfin.cpp.mesh.Cell, arg0: dolfin::Point) → bool
distance(self: dolfin.cpp.mesh.Cell, arg0: dolfin::Point) → float
facet_area(self: dolfin.cpp.mesh.Cell, arg0: int) → float
get_vertex_coordinates(self: dolfin.cpp.mesh.Cell) → List[float]

Get cell vertex coordinates

h(self: dolfin.cpp.mesh.Cell) → float
inradius(self: dolfin.cpp.mesh.Cell) → float
normal(self: dolfin.cpp.mesh.Cell, arg0: int) → dolfin::Point
orientation(*args, **kwargs)

Overloaded function.

  1. orientation(self: dolfin.cpp.mesh.Cell) -> int
  2. orientation(self: dolfin.cpp.mesh.Cell, arg0: dolfin::Point) -> int
radius_ratio(self: dolfin.cpp.mesh.Cell) → float
volume(self: dolfin.cpp.mesh.Cell) → float
class dolfin.cpp.mesh.CellIterator

Bases: pybind11_builtins.pybind11_object

DOLFIN CellIterator object

class dolfin.cpp.mesh.DomainBoundary

Bases: dolfin.cpp.mesh.SubDomain

class dolfin.cpp.mesh.Edge

Bases: dolfin.cpp.mesh.MeshEntity

DOLFIN Edge object

dot(self: dolfin.cpp.mesh.Edge, arg0: dolfin.cpp.mesh.Edge) → float
length(self: dolfin.cpp.mesh.Edge) → float
class dolfin.cpp.mesh.EdgeIterator

Bases: pybind11_builtins.pybind11_object

DOLFIN EdgeIterator object

class dolfin.cpp.mesh.Face

Bases: dolfin.cpp.mesh.MeshEntity

DOLFIN Face object

area(self: dolfin.cpp.mesh.Face) → float
normal(*args, **kwargs)

Overloaded function.

  1. normal(self: dolfin.cpp.mesh.Face) -> dolfin::Point
  2. normal(self: dolfin.cpp.mesh.Face, arg0: int) -> float
class dolfin.cpp.mesh.FaceIterator

Bases: pybind11_builtins.pybind11_object

DOLFIN FaceIterator object

class dolfin.cpp.mesh.Facet

Bases: dolfin.cpp.mesh.MeshEntity

DOLFIN Facet object

exterior(self: dolfin.cpp.mesh.Facet) → bool
normal(*args, **kwargs)

Overloaded function.

  1. normal(self: dolfin.cpp.mesh.Facet) -> dolfin::Point
  2. normal(self: dolfin.cpp.mesh.Facet, arg0: int) -> float
class dolfin.cpp.mesh.FacetIterator

Bases: pybind11_builtins.pybind11_object

DOLFIN FacetIterator object

class dolfin.cpp.mesh.Mesh

Bases: dolfin.cpp.common.Variable

DOLFIN Mesh object

bounding_box_tree(self: dolfin.cpp.mesh.Mesh) → dolfin::BoundingBoxTree
cell_name(self: dolfin.cpp.mesh.Mesh) → str
cell_orientations(self: dolfin.cpp.mesh.Mesh) → List[int]
cells(self: dolfin.cpp.mesh.Mesh) → array
color(*args, **kwargs)

Overloaded function.

  1. color(self: dolfin.cpp.mesh.Mesh, arg0: str) -> List[int]
  2. color(self: dolfin.cpp.mesh.Mesh, arg0: List[int]) -> List[int]
coordinates(self: dolfin.cpp.mesh.Mesh) → numpy.ndarray[float64[m, n], flags.writeable, flags.c_contiguous]
data(self: dolfin.cpp.mesh.Mesh) → dolfin::MeshData

Data associated with a mesh

domains(self: dolfin.cpp.mesh.Mesh) → dolfin::MeshDomains
geometric_dimension()

Returns geometric dimension for ufl interface

geometry(self: dolfin.cpp.mesh.Mesh) → dolfin.cpp.mesh.MeshGeometry

Mesh geometry

hash(self: dolfin.cpp.mesh.Mesh) → int
hmax(self: dolfin.cpp.mesh.Mesh) → float
hmin(self: dolfin.cpp.mesh.Mesh) → float
id(self: dolfin.cpp.mesh.Mesh) → int
init(*args, **kwargs)

Overloaded function.

  1. init(self: dolfin.cpp.mesh.Mesh) -> None
  2. init(self: dolfin.cpp.mesh.Mesh, arg0: int) -> int
  3. init(self: dolfin.cpp.mesh.Mesh, arg0: int, arg1: int) -> None
init_cell_orientations(*args, **kwargs)

Overloaded function.

  1. init_cell_orientations(self: dolfin.cpp.mesh.Mesh, arg0: dolfin.cpp.function.Expression) -> None
  2. init_cell_orientations(self: dolfin.cpp.mesh.Mesh, arg0: object) -> None
init_global(self: dolfin.cpp.mesh.Mesh, arg0: int) → None
mpi_comm(self: dolfin.cpp.mesh.Mesh) → MPICommWrapper
num_cells(self: dolfin.cpp.mesh.Mesh) → int

Number of cells

num_edges(self: dolfin.cpp.mesh.Mesh) → int

Number of edges

num_entities(self: dolfin.cpp.mesh.Mesh, arg0: int) → int

Number of mesh entities

num_entities_global(self: dolfin.cpp.mesh.Mesh, arg0: int) → int
num_faces(self: dolfin.cpp.mesh.Mesh) → int

Number of faces

num_facets(self: dolfin.cpp.mesh.Mesh) → int

Number of facets

num_vertices(self: dolfin.cpp.mesh.Mesh) → int

Number of vertices

order(self: dolfin.cpp.mesh.Mesh) → None
ordered(self: dolfin.cpp.mesh.Mesh) → bool
rmax(self: dolfin.cpp.mesh.Mesh) → float
rmin(self: dolfin.cpp.mesh.Mesh) → float
rotate(*args, **kwargs)

Overloaded function.

  1. rotate(self: dolfin.cpp.mesh.Mesh, arg0: float, arg1: int, arg2: dolfin::Point) -> None
  2. rotate(self: dolfin.cpp.mesh.Mesh, angle: float, axis: int=2) -> None
scale(self: dolfin.cpp.mesh.Mesh, arg0: float) → None
smooth(self: dolfin.cpp.mesh.Mesh, num_iterations: int=1) → None
smooth_boundary(self: dolfin.cpp.mesh.Mesh, arg0: int, arg1: bool) → None
snap_boundary(self: dolfin.cpp.mesh.Mesh, subdomain: dolfin::SubDomain, harmonic_smoothing: bool=True) → None
topology(self: dolfin.cpp.mesh.Mesh) → dolfin.cpp.mesh.MeshTopology

Mesh topology

translate(self: dolfin.cpp.mesh.Mesh, arg0: dolfin::Point) → None
type(self: dolfin.cpp.mesh.Mesh) → dolfin.cpp.mesh.CellType
ufl_coordinate_element()

Return the finite element of the coordinate vector field of this domain.

ufl_domain()

Returns the ufl domain corresponding to the mesh.

ufl_id(self: dolfin.cpp.mesh.Mesh) → int
class dolfin.cpp.mesh.MeshConnectivity

Bases: pybind11_builtins.pybind11_object

DOLFIN MeshConnectivity object

size(*args, **kwargs)

Overloaded function.

  1. size(self: dolfin.cpp.mesh.MeshConnectivity) -> int
  2. size(self: dolfin.cpp.mesh.MeshConnectivity, arg0: int) -> int
class dolfin.cpp.mesh.MeshData

Bases: dolfin.cpp.common.Variable

Mesh data object

array(self: dolfin.cpp.mesh.MeshData, arg0: str, arg1: int) → numpy.ndarray[uint64]
create_array(self: dolfin.cpp.mesh.MeshData, arg0: str, arg1: int) → List[int]
class dolfin.cpp.mesh.MeshDomains

Bases: pybind11_builtins.pybind11_object

Mesh domains object

get_marker(self: dolfin.cpp.mesh.MeshDomains, arg0: int, arg1: int) → int
init(self: dolfin.cpp.mesh.MeshDomains, arg0: int) → None
markers(self: dolfin.cpp.mesh.MeshDomains, arg0: int) → Dict[int, int]
set_marker(self: dolfin.cpp.mesh.MeshDomains, arg0: Tuple[int, int], arg1: int) → bool
class dolfin.cpp.mesh.MeshEditor

Bases: pybind11_builtins.pybind11_object

DOLFIN MeshEditor object

add_cell(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: List[int]) → None
add_vertex(*args, **kwargs)

Overloaded function.

  1. add_vertex(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: dolfin::Point) -> None
  2. add_vertex(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: List[float]) -> None
add_vertex_global(*args, **kwargs)

Overloaded function.

  1. add_vertex_global(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: int, arg2: dolfin::Point) -> None
  2. add_vertex_global(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: int, arg2: List[float]) -> None
close(self: dolfin.cpp.mesh.MeshEditor, order: bool=True) → None
init_cells(self: dolfin.cpp.mesh.MeshEditor, arg0: int) → None
init_cells_global(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: int) → None
init_vertices(self: dolfin.cpp.mesh.MeshEditor, arg0: int) → None
init_vertices_global(self: dolfin.cpp.mesh.MeshEditor, arg0: int, arg1: int) → None
open(self: dolfin.cpp.mesh.MeshEditor, mesh: dolfin.cpp.mesh.Mesh, type: str, tdim: int, gdim: int, degree: int=1) → None
class dolfin.cpp.mesh.MeshEntity

Bases: pybind11_builtins.pybind11_object

DOLFIN MeshEntity object

dim(self: dolfin.cpp.mesh.MeshEntity) → int

Topological dimension

entities(self: dolfin.cpp.mesh.MeshEntity, arg0: int) → numpy.ndarray[uint32[m, 1]]
global_index(self: dolfin.cpp.mesh.MeshEntity) → int

Global index

index(self: dolfin.cpp.mesh.MeshEntity) → int

Index

is_ghost(self: dolfin.cpp.mesh.MeshEntity) → bool
is_shared(self: dolfin.cpp.mesh.MeshEntity) → bool
mesh(self: dolfin.cpp.mesh.MeshEntity) → dolfin.cpp.mesh.Mesh

Associated mesh

midpoint(self: dolfin.cpp.mesh.MeshEntity) → dolfin::Point

Midpoint of Entity

num_entities(self: dolfin.cpp.mesh.MeshEntity, arg0: int) → int

Number of incident entities of given dimension

num_global_entities(self: dolfin.cpp.mesh.MeshEntity, arg0: int) → int

Global number of incident entities of given dimension

sharing_processes(self: dolfin.cpp.mesh.MeshEntity) → Set[int]
class dolfin.cpp.mesh.MeshEntityIterator

Bases: pybind11_builtins.pybind11_object

DOLFIN MeshEntityIterator object

class dolfin.cpp.mesh.MeshFunctionBool

Bases: dolfin.cpp.common.Variable

DOLFIN MeshFunction object

array(self: dolfin.cpp.mesh.MeshFunctionBool) → numpy.ndarray[bool[m, 1], flags.writeable]
dim(self: dolfin.cpp.mesh.MeshFunctionBool) → int
id(self: dolfin.cpp.mesh.MeshFunctionBool) → int
mesh(self: dolfin.cpp.mesh.MeshFunctionBool) → dolfin.cpp.mesh.Mesh
set_all(self: dolfin.cpp.mesh.MeshFunctionBool, arg0: bool) → None
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshFunctionBool, arg0: int, arg1: bool) -> None
  2. set_value(self: dolfin.cpp.mesh.MeshFunctionBool, arg0: int, arg1: bool, arg2: dolfin.cpp.mesh.Mesh) -> None
set_values(self: dolfin.cpp.mesh.MeshFunctionBool, arg0: List[bool]) → None
size(self: dolfin.cpp.mesh.MeshFunctionBool) → int
ufl_id(self: dolfin.cpp.mesh.MeshFunctionBool) → int
where_equal(self: dolfin.cpp.mesh.MeshFunctionBool, arg0: bool) → List[int]
class dolfin.cpp.mesh.MeshFunctionDouble

Bases: dolfin.cpp.common.Variable

DOLFIN MeshFunction object

array(self: dolfin.cpp.mesh.MeshFunctionDouble) → numpy.ndarray[float64[m, 1], flags.writeable]
dim(self: dolfin.cpp.mesh.MeshFunctionDouble) → int
id(self: dolfin.cpp.mesh.MeshFunctionDouble) → int
mesh(self: dolfin.cpp.mesh.MeshFunctionDouble) → dolfin.cpp.mesh.Mesh
set_all(self: dolfin.cpp.mesh.MeshFunctionDouble, arg0: float) → None
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshFunctionDouble, arg0: int, arg1: float) -> None
  2. set_value(self: dolfin.cpp.mesh.MeshFunctionDouble, arg0: int, arg1: float, arg2: dolfin.cpp.mesh.Mesh) -> None
set_values(self: dolfin.cpp.mesh.MeshFunctionDouble, arg0: List[float]) → None
size(self: dolfin.cpp.mesh.MeshFunctionDouble) → int
ufl_id(self: dolfin.cpp.mesh.MeshFunctionDouble) → int
where_equal(self: dolfin.cpp.mesh.MeshFunctionDouble, arg0: float) → List[int]
class dolfin.cpp.mesh.MeshFunctionInt

Bases: dolfin.cpp.common.Variable

DOLFIN MeshFunction object

array(self: dolfin.cpp.mesh.MeshFunctionInt) → numpy.ndarray[int32[m, 1], flags.writeable]
dim(self: dolfin.cpp.mesh.MeshFunctionInt) → int
id(self: dolfin.cpp.mesh.MeshFunctionInt) → int
mesh(self: dolfin.cpp.mesh.MeshFunctionInt) → dolfin.cpp.mesh.Mesh
set_all(self: dolfin.cpp.mesh.MeshFunctionInt, arg0: int) → None
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshFunctionInt, arg0: int, arg1: int) -> None
  2. set_value(self: dolfin.cpp.mesh.MeshFunctionInt, arg0: int, arg1: int, arg2: dolfin.cpp.mesh.Mesh) -> None
set_values(self: dolfin.cpp.mesh.MeshFunctionInt, arg0: List[int]) → None
size(self: dolfin.cpp.mesh.MeshFunctionInt) → int
ufl_id(self: dolfin.cpp.mesh.MeshFunctionInt) → int
where_equal(self: dolfin.cpp.mesh.MeshFunctionInt, arg0: int) → List[int]
class dolfin.cpp.mesh.MeshFunctionSizet

Bases: dolfin.cpp.common.Variable

DOLFIN MeshFunction object

array(self: dolfin.cpp.mesh.MeshFunctionSizet) → numpy.ndarray[uint64[m, 1], flags.writeable]
dim(self: dolfin.cpp.mesh.MeshFunctionSizet) → int
id(self: dolfin.cpp.mesh.MeshFunctionSizet) → int
mesh(self: dolfin.cpp.mesh.MeshFunctionSizet) → dolfin.cpp.mesh.Mesh
set_all(self: dolfin.cpp.mesh.MeshFunctionSizet, arg0: int) → None
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshFunctionSizet, arg0: int, arg1: int) -> None
  2. set_value(self: dolfin.cpp.mesh.MeshFunctionSizet, arg0: int, arg1: int, arg2: dolfin.cpp.mesh.Mesh) -> None
set_values(self: dolfin.cpp.mesh.MeshFunctionSizet, arg0: List[int]) → None
size(self: dolfin.cpp.mesh.MeshFunctionSizet) → int
ufl_id(self: dolfin.cpp.mesh.MeshFunctionSizet) → int
where_equal(self: dolfin.cpp.mesh.MeshFunctionSizet, arg0: int) → List[int]
class dolfin.cpp.mesh.MeshGeometry

Bases: pybind11_builtins.pybind11_object

DOLFIN MeshGeometry object

degree(self: dolfin.cpp.mesh.MeshGeometry) → int

Degree

dim(self: dolfin.cpp.mesh.MeshGeometry) → int

Geometrical dimension

get_entity_index(self: dolfin.cpp.mesh.MeshGeometry, arg0: int, arg1: int, arg2: int) → int
num_entity_coordinates(self: dolfin.cpp.mesh.MeshGeometry, arg0: int) → int
class dolfin.cpp.mesh.MeshQuality

Bases: pybind11_builtins.pybind11_object

DOLFIN MeshQuality class

aspect_ratio_gamma(arg0: dolfin.cpp.mesh.Mesh) → dolfin.cpp.mesh.MeshFunctionDouble
dihedral_angles_matplotlib_histogram(arg0: dolfin.cpp.mesh.Mesh, arg1: int) → str
dihedral_angles_min_max(arg0: dolfin.cpp.mesh.Mesh) → Tuple[float, float]
radius_ratio_histogram_data(arg0: dolfin.cpp.mesh.Mesh, arg1: int) → Tuple[List[float], List[float]]
radius_ratio_matplotlib_histogram(mesh: dolfin.cpp.mesh.Mesh, num_bins: int=50) → str
radius_ratio_min_max(arg0: dolfin.cpp.mesh.Mesh) → Tuple[float, float]
radius_ratios(arg0: dolfin.cpp.mesh.Mesh) → dolfin.cpp.mesh.MeshFunctionDouble
class dolfin.cpp.mesh.MeshTopology

Bases: dolfin.cpp.common.Variable

DOLFIN MeshTopology object

cell_owner(self: dolfin.cpp.mesh.MeshTopology) → List[int]
dim(self: dolfin.cpp.mesh.MeshTopology) → int

Topological dimension

ghost_offset(self: dolfin.cpp.mesh.MeshTopology, arg0: int) → int
global_indices(self: dolfin.cpp.mesh.MeshTopology, arg0: int) → numpy.ndarray[int64]
hash(self: dolfin.cpp.mesh.MeshTopology) → int
have_global_indices(self: dolfin.cpp.mesh.MeshTopology, arg0: int) → bool
have_shared_entities(self: dolfin.cpp.mesh.MeshTopology, arg0: int) → bool
init(*args, **kwargs)

Overloaded function.

  1. init(self: dolfin.cpp.mesh.MeshTopology, arg0: int) -> None
  2. init(self: dolfin.cpp.mesh.MeshTopology, arg0: int, arg1: int, arg2: int) -> None
init_global_indices(self: dolfin.cpp.mesh.MeshTopology, arg0: int, arg1: int) → None
set_global_index(self: dolfin.cpp.mesh.MeshTopology, arg0: int, arg1: int, arg2: int) → None
shared_entities(self: dolfin.cpp.mesh.MeshTopology, arg0: int) → Dict[int, Set[int]]
size(self: dolfin.cpp.mesh.MeshTopology, arg0: int) → int
str(self: dolfin.cpp.mesh.MeshTopology, arg0: bool) → str
class dolfin.cpp.mesh.MeshValueCollection_bool

Bases: dolfin.cpp.common.Variable

DOLFIN MeshValueCollection object

assign(*args, **kwargs)

Overloaded function.

  1. assign(self: dolfin.cpp.mesh.MeshValueCollection_bool, arg0: dolfin.cpp.mesh.MeshFunctionBool) -> None
  2. assign(self: dolfin.cpp.mesh.MeshValueCollection_bool, arg0: dolfin.cpp.mesh.MeshValueCollection_bool) -> None
dim(self: dolfin.cpp.mesh.MeshValueCollection_bool) → int
get_value(self: dolfin.cpp.mesh.MeshValueCollection_bool, arg0: int, arg1: int) → bool
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshValueCollection_bool, arg0: int, arg1: bool) -> bool
  2. set_value(self: dolfin.cpp.mesh.MeshValueCollection_bool, arg0: int, arg1: int, arg2: bool) -> bool
size(self: dolfin.cpp.mesh.MeshValueCollection_bool) → int
values(self: dolfin.cpp.mesh.MeshValueCollection_bool) → Dict[Tuple[int, int], bool]
class dolfin.cpp.mesh.MeshValueCollection_double

Bases: dolfin.cpp.common.Variable

DOLFIN MeshValueCollection object

assign(*args, **kwargs)

Overloaded function.

  1. assign(self: dolfin.cpp.mesh.MeshValueCollection_double, arg0: dolfin.cpp.mesh.MeshFunctionDouble) -> None
  2. assign(self: dolfin.cpp.mesh.MeshValueCollection_double, arg0: dolfin.cpp.mesh.MeshValueCollection_double) -> None
dim(self: dolfin.cpp.mesh.MeshValueCollection_double) → int
get_value(self: dolfin.cpp.mesh.MeshValueCollection_double, arg0: int, arg1: int) → float
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshValueCollection_double, arg0: int, arg1: float) -> bool
  2. set_value(self: dolfin.cpp.mesh.MeshValueCollection_double, arg0: int, arg1: int, arg2: float) -> bool
size(self: dolfin.cpp.mesh.MeshValueCollection_double) → int
values(self: dolfin.cpp.mesh.MeshValueCollection_double) → Dict[Tuple[int, int], float]
class dolfin.cpp.mesh.MeshValueCollection_int

Bases: dolfin.cpp.common.Variable

DOLFIN MeshValueCollection object

assign(*args, **kwargs)

Overloaded function.

  1. assign(self: dolfin.cpp.mesh.MeshValueCollection_int, arg0: dolfin.cpp.mesh.MeshFunctionInt) -> None
  2. assign(self: dolfin.cpp.mesh.MeshValueCollection_int, arg0: dolfin.cpp.mesh.MeshValueCollection_int) -> None
dim(self: dolfin.cpp.mesh.MeshValueCollection_int) → int
get_value(self: dolfin.cpp.mesh.MeshValueCollection_int, arg0: int, arg1: int) → int
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshValueCollection_int, arg0: int, arg1: int) -> bool
  2. set_value(self: dolfin.cpp.mesh.MeshValueCollection_int, arg0: int, arg1: int, arg2: int) -> bool
size(self: dolfin.cpp.mesh.MeshValueCollection_int) → int
values(self: dolfin.cpp.mesh.MeshValueCollection_int) → Dict[Tuple[int, int], int]
class dolfin.cpp.mesh.MeshValueCollection_sizet

Bases: dolfin.cpp.common.Variable

DOLFIN MeshValueCollection object

assign(*args, **kwargs)

Overloaded function.

  1. assign(self: dolfin.cpp.mesh.MeshValueCollection_sizet, arg0: dolfin.cpp.mesh.MeshFunctionSizet) -> None
  2. assign(self: dolfin.cpp.mesh.MeshValueCollection_sizet, arg0: dolfin.cpp.mesh.MeshValueCollection_sizet) -> None
dim(self: dolfin.cpp.mesh.MeshValueCollection_sizet) → int
get_value(self: dolfin.cpp.mesh.MeshValueCollection_sizet, arg0: int, arg1: int) → int
set_value(*args, **kwargs)

Overloaded function.

  1. set_value(self: dolfin.cpp.mesh.MeshValueCollection_sizet, arg0: int, arg1: int) -> bool
  2. set_value(self: dolfin.cpp.mesh.MeshValueCollection_sizet, arg0: int, arg1: int, arg2: int) -> bool
size(self: dolfin.cpp.mesh.MeshValueCollection_sizet) → int
values(self: dolfin.cpp.mesh.MeshValueCollection_sizet) → Dict[Tuple[int, int], int]
class dolfin.cpp.mesh.MultiMesh

Bases: dolfin.cpp.common.Variable

add(self: dolfin.cpp.mesh.MultiMesh, arg0: dolfin.cpp.mesh.Mesh) → None
auto_cover(self: dolfin.cpp.mesh.MultiMesh, arg0: int, arg1: dolfin::Point) → None

Marks all uncut and cut cells connected to the given point as covered.

build(self: dolfin.cpp.mesh.MultiMesh, quadrature_order: int=2) → None
compute_area(self: dolfin.cpp.mesh.MultiMesh) → float
compute_volume(self: dolfin.cpp.mesh.MultiMesh) → float
covered_cells(self: dolfin.cpp.mesh.MultiMesh, arg0: int) → List[int]
cut_cells(self: dolfin.cpp.mesh.MultiMesh, arg0: int) → List[int]
mark_covered(self: dolfin.cpp.mesh.MultiMesh, arg0: int, arg1: List[int]) → None

Function that marks a set of cells, given by indicies in a list, as covered.

num_parts(self: dolfin.cpp.mesh.MultiMesh) → int
part(self: dolfin.cpp.mesh.MultiMesh, arg0: int) → dolfin.cpp.mesh.Mesh
quadrature_rules_cut_cells(self: dolfin.cpp.mesh.MultiMesh, arg0: int, arg1: int) → Tuple[List[float], List[float]]
quadrature_rules_interface(*args, **kwargs)

Overloaded function.

  1. quadrature_rules_interface(self: dolfin.cpp.mesh.MultiMesh, arg0: int) -> Dict[int, List[Tuple[List[float], List[float]]]]
  2. quadrature_rules_interface(self: dolfin.cpp.mesh.MultiMesh, arg0: int, arg1: int) -> List[Tuple[List[float], List[float]]]
ufl_cell()

Returns the ufl cell of the mesh.

ufl_coordinate_element()

Return the finite element of the coordinate vector field of this domain.

ufl_domain()

Returns the ufl domain corresponding to the mesh.

ufl_id()

Returns an id that UFL can use to decide if two objects are the same.

uncut_cells(self: dolfin.cpp.mesh.MultiMesh, arg0: int) → List[int]
class dolfin.cpp.mesh.SubDomain

Bases: pybind11_builtins.pybind11_object

DOLFIN SubDomain object

get_property(self: dolfin.cpp.mesh.SubDomain, arg0: str) → float
inside(self: dolfin.cpp.mesh.SubDomain, arg0: numpy.ndarray[float64[m, 1]], arg1: bool) → bool
map(self: dolfin.cpp.mesh.SubDomain, arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, 1], flags.writeable]) → None
mark(*args, **kwargs)

Overloaded function.

  1. mark(self: dolfin.cpp.mesh.SubDomain, meshfunction: dolfin.cpp.mesh.MeshFunctionSizet, marker: int, check_midpoint: bool=True) -> None
  2. mark(self: dolfin.cpp.mesh.SubDomain, meshfunction: dolfin.cpp.mesh.MeshFunctionDouble, marker: float, check_midpoint: bool=True) -> None
  3. mark(self: dolfin.cpp.mesh.SubDomain, meshfunction: dolfin.cpp.mesh.MeshFunctionBool, marker: bool, check_midpoint: bool=True) -> None
mark_cells(self: dolfin.cpp.mesh.SubDomain, mesh: dolfin.cpp.mesh.Mesh, sub_domain: int, check_midpoint: bool=True) → None
mark_facets(self: dolfin.cpp.mesh.SubDomain, mesh: dolfin.cpp.mesh.Mesh, sub_domain: int, check_midpoint: bool=True) → None
set_property(self: dolfin.cpp.mesh.SubDomain, arg0: str, arg1: float) → None
class dolfin.cpp.mesh.SubMesh

Bases: dolfin.cpp.mesh.Mesh

DOLFIN SubMesh

class dolfin.cpp.mesh.Vertex

Bases: dolfin.cpp.mesh.MeshEntity

DOLFIN Vertex object

point(self: dolfin.cpp.mesh.Vertex) → dolfin::Point
x(self: dolfin.cpp.mesh.Vertex, arg0: int) → float
class dolfin.cpp.mesh.VertexIterator

Bases: pybind11_builtins.pybind11_object

DOLFIN VertexIterator object

dolfin.cpp.mesh.cells(*args, **kwargs)

Overloaded function.

  1. cells(mesh: dolfin.cpp.mesh.Mesh, type: str=’regular’) -> dolfin.cpp.mesh.CellIterator
  2. cells(arg0: dolfin.cpp.mesh.MeshEntity) -> dolfin.cpp.mesh.CellIterator
dolfin.cpp.mesh.edges(*args, **kwargs)

Overloaded function.

  1. edges(mesh: dolfin.cpp.mesh.Mesh, type: str=’regular’) -> dolfin.cpp.mesh.EdgeIterator
  2. edges(arg0: dolfin.cpp.mesh.MeshEntity) -> dolfin.cpp.mesh.EdgeIterator
dolfin.cpp.mesh.entities(*args, **kwargs)

Overloaded function.

  1. entities(arg0: dolfin.cpp.mesh.Mesh, arg1: int) -> dolfin.cpp.mesh.MeshEntityIterator
  2. entities(arg0: dolfin.cpp.mesh.MeshEntity, arg1: int) -> dolfin.cpp.mesh.MeshEntityIterator
dolfin.cpp.mesh.faces(*args, **kwargs)

Overloaded function.

  1. faces(mesh: dolfin.cpp.mesh.Mesh, type: str=’regular’) -> dolfin.cpp.mesh.FaceIterator
  2. faces(arg0: dolfin.cpp.mesh.MeshEntity) -> dolfin.cpp.mesh.FaceIterator
dolfin.cpp.mesh.facets(*args, **kwargs)

Overloaded function.

  1. facets(mesh: dolfin.cpp.mesh.Mesh, type: str=’regular’) -> dolfin.cpp.mesh.FacetIterator
  2. facets(arg0: dolfin.cpp.mesh.MeshEntity) -> dolfin.cpp.mesh.FacetIterator
dolfin.cpp.mesh.make_dolfin_subdomain(arg0: int) → dolfin::SubDomain
dolfin.cpp.mesh.vertices(*args, **kwargs)

Overloaded function.

  1. vertices(mesh: dolfin.cpp.mesh.Mesh, type: str=’regular’) -> dolfin.cpp.mesh.VertexIterator
  2. vertices(arg0: dolfin.cpp.mesh.MeshEntity) -> dolfin.cpp.mesh.VertexIterator