dolfin/generation

Documentation for C++ code found in dolfin/generation/*.h

Classes

BoxMesh

C++ documentation for BoxMesh from dolfin/generation/BoxMesh.h:

class dolfin::BoxMesh : public dolfin::Mesh

Tetrahedral mesh of the 3D rectangular prism spanned by two points p0 and p1. Given the number of cells (nx, ny, nz) in each direction, the total number of tetrahedra will be 6*nx*ny*nz and the total number of vertices will be (nx + 1)*(ny + 1)*(nz + 1).

dolfin::BoxMesh::BoxMesh(MPI_Comm comm, const Point &p0, const Point &p1, std::size_t nx, std::size_t ny, std::size_t nz)

Create a uniform finite element Mesh over the rectangular prism spanned by the two _Point_s p0 and p1. The order of the two points is not important in terms of minimum and maximum coordinates.

// Mesh with 8 cells in each direction on the
// set [-1,2] x [-1,2] x [-1,2].
Point p0(-1, -1, -1);
Point p1(2, 2, 2);
BoxMesh mesh(MPI_COMM_WORLD, p0, p1, 8, 8, 8);
Parameters:
  • comm – (MPI_Comm) MPI communicator
  • p0 – (Point ) First point.
  • p1 – (Point ) Second point.
  • nx – (double) Number of cells in x-direction.
  • ny – (double) Number of cells in y-direction.
  • nz – (double) Number of cells in z-direction.
dolfin::BoxMesh::BoxMesh(const Point &p0, const Point &p1, std::size_t nx, std::size_t ny, std::size_t nz)

Create a uniform finite element Mesh over the rectangular prism spanned by the two _Point_s p0 and p1. The order of the two points is not important in terms of minimum and maximum coordinates.

// Mesh with 8 cells in each direction on the
// set [-1,2] x [-1,2] x [-1,2].
Point p0(-1, -1, -1);
Point p1(2, 2, 2);
BoxMesh mesh(p0, p1, 8, 8, 8);
Parameters:
  • p0 – (Point ) First point.
  • p1 – (Point ) Second point.
  • nx – (double) Number of cells in x-direction.
  • ny – (double) Number of cells in y-direction.
  • nz – (double) Number of cells in z-direction.
void dolfin::BoxMesh::build(const Point &p0, const Point &p1, std::size_t nx, std::size_t ny, std::size_t nz)
Parameters:
  • p0
  • p1
  • nx
  • ny
  • nz

IntervalMesh

C++ documentation for IntervalMesh from dolfin/generation/IntervalMesh.h:

class dolfin::IntervalMesh : public dolfin::Mesh

Interval mesh of the 1D line [a,b]. Given the number of cells (nx) in the axial direction, the total number of intervals will be nx and the total number of vertices will be (nx + 1).

dolfin::IntervalMesh::IntervalMesh(MPI_Comm comm, std::size_t nx, double a, double b)

Constructor

// Create a mesh of 25 cells in the interval [-1,1]
IntervalMesh mesh(MPI_COMM_WORLD, 25, -1.0, 1.0);
Parameters:
  • comm – (MPI_Comm) MPI communicator
  • nx – (std::size_t) The number of cells.
  • a – (double) The minimum point (inclusive).
  • b – (double) The maximum point (inclusive).
dolfin::IntervalMesh::IntervalMesh(std::size_t nx, double a, double b)

Constructor

// Create a mesh of 25 cells in the interval [-1,1]
IntervalMesh mesh(25, -1.0, 1.0);
Parameters:
  • nx – (std::size_t) The number of cells.
  • a – (double) The minimum point (inclusive).
  • b – (double) The maximum point (inclusive).
void dolfin::IntervalMesh::build(std::size_t nx, double a, double b)
Parameters:
  • nx
  • a
  • b

RectangleMesh

C++ documentation for RectangleMesh from dolfin/generation/RectangleMesh.h:

class dolfin::RectangleMesh : public dolfin::Mesh

Triangular mesh of the 2D rectangle spanned by two points p0 and p1. Given the number of cells (nx, ny) in each direction, the total number of triangles will be 2*nx*ny and the total number of vertices will be (nx + 1)*(ny + 1).

dolfin::RectangleMesh::RectangleMesh(MPI_Comm comm, const Point &p0, const Point &p1, std::size_t nx, std::size_t ny, std::string diagonal = "right")
// Mesh with 8 cells in each direction on the
// set [-1,2] x [-1,2]
Point p0(-1, -1);
Point p1(2, 2);
RectangleMesh mesh(MPI_COMM_WORLD, p0, p1, 8, 8);
Parameters:
  • comm – (MPI_Comm) MPI communicator
  • p0 – (Point ) First point.
  • p1 – (Point ) Second point.
  • nx – (double) Number of cells in \(x\) -direction.
  • ny – (double) Number of cells in \(y\) -direction.
  • diagonal – (string) Direction of diagonals: “left”, “right”, “left/right”, “crossed”
dolfin::RectangleMesh::RectangleMesh(const Point &p0, const Point &p1, std::size_t nx, std::size_t ny, std::string diagonal = "right")
// Mesh with 8 cells in each direction on the
// set [-1,2] x [-1,2]
Point p0(-1, -1);
Point p1(2, 2);
RectangleMesh mesh(p0, p1, 8, 8);
Parameters:
  • p0 – (Point ) First point.
  • p1 – (Point ) Second point.
  • nx – (double) Number of cells in \(x\) -direction.
  • ny – (double) Number of cells in \(y\) -direction.
  • diagonal – (string) Direction of diagonals: “left”, “right”, “left/right”, “crossed”
void dolfin::RectangleMesh::build(const Point &p0, const Point &p1, std::size_t nx, std::size_t ny, std::string diagonal = "right")
Parameters:
  • p0
  • p1
  • nx
  • ny
  • diagonal

SphericalShellMesh

C++ documentation for SphericalShellMesh from dolfin/generation/SphericalShellMesh.h:

class dolfin::SphericalShellMesh : public dolfin::Mesh

Spherical shell approximation, icosahedral mesh, with degree=1 or degree=2.

dolfin::SphericalShellMesh::SphericalShellMesh(MPI_Comm comm, std::size_t degree)

Create a spherical shell manifold for testing.

Parameters:
  • comm
  • degree

UnitCubeMesh

C++ documentation for UnitCubeMesh from dolfin/generation/UnitCubeMesh.h:

class dolfin::UnitCubeMesh : public dolfin::BoxMesh

Tetrahedral mesh of the 3D unit cube [0,1] x [0,1] x [0,1]. Given the number of cells (nx, ny, nz) in each direction, the total number of tetrahedra will be 6*nx*ny*nz and the total number of vertices will be (nx + 1)*(ny + 1)*(nz + 1).

dolfin::UnitCubeMesh::UnitCubeMesh(MPI_Comm comm, std::size_t nx, std::size_t ny, std::size_t nz)

Create a uniform finite element Mesh over the unit cube [0,1] x [0,1] x [0,1].

UnitCubeMesh mesh(MPI_COMM_WORLD, 32, 32, 32);
Parameters:
  • comm – (MPI_Comm) MPI communicator
  • nx – (std::size_t) Number of cells in \(x\) direction.
  • ny – (std::size_t) Number of cells in \(y\) direction.
  • nz – (std::size_t) Number of cells in \(z\) direction.
dolfin::UnitCubeMesh::UnitCubeMesh(std::size_t nx, std::size_t ny, std::size_t nz)

Create a uniform finite element Mesh over the unit cube [0,1] x [0,1] x [0,1].

UnitCubeMesh mesh(32, 32, 32);
Parameters:
  • nx – (std::size_t) Number of cells in \(x\) direction.
  • ny – (std::size_t) Number of cells in \(y\) direction.
  • nz – (std::size_t) Number of cells in \(z\) direction.

UnitDiscMesh

C++ documentation for UnitDiscMesh from dolfin/generation/UnitDiscMesh.h:

class dolfin::UnitDiscMesh : public dolfin::Mesh

A unit disc mesh in 2D or 3D geometry.

dolfin::UnitDiscMesh::UnitDiscMesh(MPI_Comm comm, std::size_t n, std::size_t degree, std::size_t gdim)

Create a unit disc mesh in 2D or 3D geometry with n steps, and given degree polynomial mesh

Parameters:
  • comm
  • n
  • degree
  • gdim

UnitHexMesh

C++ documentation for UnitHexMesh from dolfin/generation/UnitHexMesh.h:

class dolfin::UnitHexMesh : public dolfin::Mesh

NB: this code is experimental, just for testing, and will generally not work with anything else

dolfin::UnitHexMesh::UnitHexMesh(MPI_Comm comm, std::size_t nx, std::size_t ny, std::size_t nz)

NB: this code is experimental, just for testing, and will generally not work with anything else

Parameters:
  • comm
  • nx
  • ny
  • nz

UnitIntervalMesh

C++ documentation for UnitIntervalMesh from dolfin/generation/UnitIntervalMesh.h:

class dolfin::UnitIntervalMesh : public dolfin::IntervalMesh

A mesh of the unit interval (0, 1) with a given number of cells (nx) in the axial direction. The total number of intervals will be nx and the total number of vertices will be (nx + 1).

dolfin::UnitIntervalMesh::UnitIntervalMesh(MPI_Comm comm, std::size_t nx)

Constructor

// Create a mesh of 25 cells in the interval [0,1]
UnitIntervalMesh mesh(MPI_COMM_WORLD, 25);
Parameters:
  • comm – (MPI_Comm) MPI communicator
  • nx – (std::size_t) The number of cells.
dolfin::UnitIntervalMesh::UnitIntervalMesh(std::size_t nx)

Constructor

// Create a mesh of 25 cells in the interval [0,1]
UnitIntervalMesh mesh(25);
Parameters:nx – (std::size_t) The number of cells.

UnitQuadMesh

C++ documentation for UnitQuadMesh from dolfin/generation/UnitQuadMesh.h:

class dolfin::UnitQuadMesh : public dolfin::Mesh

NB: this code is experimental, just for testing, and will generally not work with anything else

dolfin::UnitQuadMesh::UnitQuadMesh(MPI_Comm comm, std::size_t nx, std::size_t ny)

NB: this code is experimental, just for testing, and will generally not work with anything else

Parameters:
  • comm
  • nx
  • ny

UnitSquareMesh

C++ documentation for UnitSquareMesh from dolfin/generation/UnitSquareMesh.h:

class dolfin::UnitSquareMesh : public dolfin::RectangleMesh

Triangular mesh of the 2D unit square [0,1] x [0,1]. Given the number of cells (nx, ny) in each direction, the total number of triangles will be 2*nx*ny and the total number of vertices will be (nx + 1)*(ny + 1). std::string diagonal (“left”, “right”, “right/left”, “left/right”, or “crossed”) indicates the direction of the diagonals.

dolfin::UnitSquareMesh::UnitSquareMesh(MPI_Comm comm, std::size_t nx, std::size_t ny, std::string diagonal = "right")

Create a uniform finite element Mesh over the unit square [0,1] x [0,1].

UnitSquareMesh mesh1(MPI_COMM_WORLD, 32, 32);
UnitSquareMesh mesh2(MPI_COMM_WORLD, 32, 32, "crossed");
Parameters:
  • comm – (MPI_Comm) MPI communicator
  • nx – (std::size_t) Number of cells in horizontal direction.
  • ny – (std::size_t) Number of cells in vertical direction.
  • diagonal – (std::string) Optional argument: A std::string indicating the direction of the diagonals.
dolfin::UnitSquareMesh::UnitSquareMesh(std::size_t nx, std::size_t ny, std::string diagonal = "right")

Create a uniform finite element Mesh over the unit square [0,1] x [0,1].

UnitSquareMesh mesh1(32, 32);
UnitSquareMesh mesh2(32, 32, "crossed");
Parameters:
  • nx – (std::size_t) Number of cells in horizontal direction.
  • ny – (std::size_t) Number of cells in vertical direction.
  • diagonal – (std::string) Optional argument: A std::string indicating the direction of the diagonals.

UnitTetrahedronMesh

C++ documentation for UnitTetrahedronMesh from dolfin/generation/UnitTetrahedronMesh.h:

class dolfin::UnitTetrahedronMesh : public dolfin::Mesh

A mesh consisting of a single tetrahedron with vertices at (0, 0, 0) (1, 0, 0) (0, 1, 0) (0, 0, 1) This class is useful for testing.

dolfin::UnitTetrahedronMesh::UnitTetrahedronMesh()

Create mesh of unit tetrahedron.

UnitTriangleMesh

C++ documentation for UnitTriangleMesh from dolfin/generation/UnitTriangleMesh.h:

class dolfin::UnitTriangleMesh : public dolfin::Mesh

A mesh consisting of a single triangle with vertices at (0, 0) (1, 0) (0, 1) This class is useful for testing.

dolfin::UnitTriangleMesh::UnitTriangleMesh()

Create mesh of unit triangle.