HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Classes | Typedefs | Functions
Mesh

Classes to construct and describe a 2D mesh. More...

Classes

class  Mesh2D::Mesh
 
class  Mesh2D::MeshBuilder
 The MeshBuilder class provides build tools to create a full mesh with all connectivities. More...
 
class  Mesh2D::MeshReaderTyp2
 The MeshReaderTyp2 class provides functions to read a typ2 mesh file. More...
 
class  Mesh2D::Polytope< object_dim >
 A Vertex is a Polytope with object_dim = 0. More...
 

Typedefs

typedef std::function< std::array< double, 2 >const std::array< double, 2 > &)> Mesh2D::MeshBuilder::TransformationType
 
using Mesh2D::Vertex = Polytope< 0 >
 An Edge is a Polytope with object_dim = 1. More...
 
using Mesh2D::Edge = Polytope< 1 >
 A Face is a Polytope with object_dim = DIMENSION - 1. More...
 
using Mesh2D::Face = Polytope< DIMENSION - 1 >
 A Cell is a Polytope with object_dim = DIMENSION. More...
 
using Mesh2D::Cell = Polytope< DIMENSION >
 

Functions

 Mesh2D::Mesh::Mesh ()
 
 Mesh2D::Mesh::~Mesh ()
 
void Mesh2D::Mesh::set_name (std::string name)
 set the name of the mesh More...
 
std::string Mesh2D::Mesh::get_name ()
 getter for the mesh name More...
 
double Mesh2D::Mesh::h_max () const
 < max diameter of cells More...
 
std::size_t Mesh2D::Mesh::dim () const
 dimension of the mesh More...
 
std::size_t Mesh2D::Mesh::n_vertices () const
 number of vertices in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_edges () const
 number of edges in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_faces () const
 number of faces in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_cells () const
 number of cells in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_b_vertices () const
 number of boundary vertices in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_b_edges () const
 number of boundary edges in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_b_faces () const
 number of boundary faces in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_b_cells () const
 number of boundary cells in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_i_vertices () const
 number of internal vertices in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_i_edges () const
 number of internal edges in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_i_faces () const
 number of internal faces in the mesh. More...
 
std::size_t Mesh2D::Mesh::n_i_cells () const
 number of internal cells in the mesh. More...
 
std::vector< Vertex * > Mesh2D::Mesh::get_vertices () const
 lists the vertices in the mesh. More...
 
std::vector< Edge * > Mesh2D::Mesh::get_edges () const
 lists the edges in the mesh. More...
 
std::vector< Face * > Mesh2D::Mesh::get_faces () const
 lists the faces in the mesh. More...
 
std::vector< Cell * > Mesh2D::Mesh::get_cells () const
 lists the cells in the mesh. More...
 
std::vector< Vertex * > Mesh2D::Mesh::get_b_vertices () const
 lists the boundary vertices in the mesh. More...
 
std::vector< Edge * > Mesh2D::Mesh::get_b_edges () const
 lists the boundary edges in the mesh. More...
 
std::vector< Face * > Mesh2D::Mesh::get_b_faces () const
 lists the boundary faces in the mesh. More...
 
std::vector< Cell * > Mesh2D::Mesh::get_b_cells () const
 lists the boundary cells in the mesh. More...
 
std::vector< Vertex * > Mesh2D::Mesh::get_i_vertices () const
 lists the internal vertices in the mesh. More...
 
std::vector< Edge * > Mesh2D::Mesh::get_i_edges () const
 lists the internal edges in the mesh. More...
 
std::vector< Face * > Mesh2D::Mesh::get_i_faces () const
 lists the internal faces in the mesh. More...
 
std::vector< Cell * > Mesh2D::Mesh::get_i_cells () const
 lists the internal cells in the mesh. More...
 
void Mesh2D::Mesh::add_vertex (Vertex *vertex)
 
void Mesh2D::Mesh::add_edge (Edge *edge)
 
void Mesh2D::Mesh::add_face (Face *face)
 
void Mesh2D::Mesh::add_cell (Cell *cell)
 
void Mesh2D::Mesh::add_b_vertex (Vertex *vertex)
 
void Mesh2D::Mesh::add_b_edge (Edge *edge)
 
void Mesh2D::Mesh::add_b_face (Face *face)
 
void Mesh2D::Mesh::add_b_cell (Cell *cell)
 
void Mesh2D::Mesh::add_i_vertex (Vertex *vertex)
 
void Mesh2D::Mesh::add_i_edge (Edge *edge)
 
void Mesh2D::Mesh::add_i_face (Face *face)
 
void Mesh2D::Mesh::add_i_cell (Cell *cell)
 
VertexMesh2D::Mesh::vertex (std::size_t index) const
 get a constant pointer to a vertex using its global index More...
 
EdgeMesh2D::Mesh::edge (std::size_t index) const
 get a constant pointer to a edge using its global index More...
 
FaceMesh2D::Mesh::face (std::size_t index) const
 get a constant pointer to a face using its global index More...
 
CellMesh2D::Mesh::cell (std::size_t index) const
 get a constant pointer to a cell using its global index More...
 
VertexMesh2D::Mesh::b_vertex (std::size_t index) const
 get a constant pointer to a boundary vertex using an index More...
 
EdgeMesh2D::Mesh::b_edge (std::size_t index) const
 get a constant pointer to boundary a edge using an index More...
 
FaceMesh2D::Mesh::b_face (std::size_t index) const
 get a constant pointer to boundary a face using an index More...
 
CellMesh2D::Mesh::b_cell (std::size_t index) const
 get a constant pointer to boundary a cell using an index More...
 
VertexMesh2D::Mesh::i_vertex (std::size_t index) const
 get a constant pointer to an internal vertex using an index More...
 
EdgeMesh2D::Mesh::i_edge (std::size_t index) const
 get a constant pointer to an internal edge using an index More...
 
FaceMesh2D::Mesh::i_face (std::size_t index) const
 get a constant pointer to an internal face using an index More...
 
CellMesh2D::Mesh::i_cell (std::size_t index) const
 get a constant pointer to an internal cell using an index More...
 
std::vector< double > Mesh2D::Mesh::regularity ()
 
void Mesh2D::Mesh::renum (const char B, const std::vector< size_t > new_to_old)
 
size_t Mesh2D::Mesh::find_cell (const VectorRd x) const
 < returns the index of the cell containing the point x More...
 
void Mesh2D::Mesh::plot_simplices (std::ofstream *out)
 
void Mesh2D::Mesh::plot_mesh (std::ofstream *out)
 
 Mesh2D::MeshBuilder::MeshBuilder ()
 
 Mesh2D::MeshBuilder::MeshBuilder (const std::string mesh_file)
 
std::unique_ptr< MeshMesh2D::MeshBuilder::build_the_mesh (const TransformationType &transformation=[](const std::array< double, 2 > &x) { return x;})
 
 Mesh2D::MeshReaderTyp2::MeshReaderTyp2 (std::string file_name)
 
 Mesh2D::MeshReaderTyp2::~MeshReaderTyp2 ()
 
void Mesh2D::MeshReaderTyp2::read_mesh (std::vector< std::array< double, 2 >> &vertices, std::vector< std::vector< std::size_t >> &cells)
 
 Mesh2D::Polytope< object_dim >::Polytope (size_t index, Simplices< object_dim > simplices)
 < Constructor for a Polytope defined by its simplices More...
 
 Mesh2D::Polytope< object_dim >::Polytope (size_t index, Simplex< object_dim > simplex)
 Constructor for a Polytope defined by a single coordinate (Vertex) More...
 
 Mesh2D::Polytope< object_dim >::Polytope (size_t index, VectorRd vertex)
 Null constructor. More...
 
 Mesh2D::Polytope< object_dim >::Polytope ()
 Destructor. More...
 
 Mesh2D::Polytope< object_dim >::~Polytope ()
 
size_t Mesh2D::Polytope< object_dim >::global_index () const
 Return the global index of the Polytope. More...
 
double Mesh2D::Polytope< object_dim >::diam () const
 Return the diameter of the Polytope. More...
 
VectorRd Mesh2D::Polytope< object_dim >::center_mass () const
 Return the center mass of the Polytope. More...
 
double Mesh2D::Polytope< object_dim >::measure () const
 Return the Lebesgue measure of the Polytope. More...
 
Simplices< object_dim > Mesh2D::Polytope< object_dim >::get_simplices () const
 Return the simplices making up the Polytope. More...
 
void Mesh2D::Polytope< object_dim >::set_global_index (const size_t idx)
 Set the global index. More...
 
bool Mesh2D::Polytope< object_dim >::is_boundary () const
 Return true if Polytope is a boundary object, false otherwise. More...
 
void Mesh2D::Polytope< object_dim >::set_boundary (bool val)
 Set the boundary value of the Polytope. More...
 
std::vector< Polytope< 0 > * > Mesh2D::Polytope< object_dim >::get_vertices () const
 Return the vertices of the Polytope. More...
 
std::vector< Polytope< 1 > * > Mesh2D::Polytope< object_dim >::get_edges () const
 Return the edges of the Polytope. More...
 
std::vector< Polytope< DIMENSION - 1 > * > Mesh2D::Polytope< object_dim >::get_faces () const
 Return the faces of the Polytope. More...
 
std::vector< Polytope< DIMENSION > * > Mesh2D::Polytope< object_dim >::get_cells () const
 Return the cells of the Polytope. More...
 
size_t Mesh2D::Polytope< object_dim >::n_vertices () const
 Return the number of vertices of the Polytope. More...
 
size_t Mesh2D::Polytope< object_dim >::n_edges () const
 Return the number of edges of the Polytope. More...
 
size_t Mesh2D::Polytope< object_dim >::n_faces () const
 Return the number of faces of the Polytope. More...
 
size_t Mesh2D::Polytope< object_dim >::n_cells () const
 Return the number of cells of the Polytope. More...
 
Polytope< 0 > * Mesh2D::Polytope< object_dim >::vertex (const size_t i) const
 Return the i-th vertex of the Polytope. More...
 
Polytope< 1 > * Mesh2D::Polytope< object_dim >::edge (const size_t i) const
 Return the i-th edge of the Polytope. More...
 
Polytope< DIMENSION - 1 > * Mesh2D::Polytope< object_dim >::face (const size_t i) const
 Return the i-th face of the Polytope. More...
 
Polytope< DIMENSION > * Mesh2D::Polytope< object_dim >::cell (const size_t i) const
 Return the i-th cell of the Polytope. More...
 
void Mesh2D::Polytope< object_dim >::add_vertex (Polytope< 0 > *vertex)
 Add a vertex to the Polytope. More...
 
void Mesh2D::Polytope< object_dim >::add_edge (Polytope< 1 > *edge)
 Add an edge to the Polytope. More...
 
void Mesh2D::Polytope< object_dim >::add_face (Polytope< DIMENSION - 1 > *face)
 Add a face to the Polytope. More...
 
void Mesh2D::Polytope< object_dim >::add_cell (Polytope< DIMENSION > *cell)
 Add a cell to the Polytope. More...
 
int Mesh2D::Polytope< object_dim >::index_vertex (const Polytope< 0 > *vertex) const
 Returns the local index of a vertex. More...
 
int Mesh2D::Polytope< object_dim >::index_edge (const Polytope< 1 > *edge) const
 Returns the local index of an edge. More...
 
int Mesh2D::Polytope< object_dim >::index_face (const Polytope< DIMENSION - 1 > *face) const
 Returns the local index of a face. More...
 
int Mesh2D::Polytope< object_dim >::index_cell (const Polytope< DIMENSION > *cell) const
 Returns the local index of a cell. More...
 
VectorRd Mesh2D::Polytope< object_dim >::coords () const
 Return the coordinates of a Vertex. More...
 
VectorRd Mesh2D::Polytope< object_dim >::face_normal (const size_t face_index) const
 Return the outer normal of a Cell towards the Face located at face_index. More...
 
VectorRd Mesh2D::Polytope< object_dim >::edge_normal (const size_t edge_index) const
 Return the edge normal of a 2D object. More...
 
int Mesh2D::Polytope< object_dim >::face_orientation (const size_t face_index) const
 Return the orientation of a Face. More...
 
int Mesh2D::Polytope< object_dim >::edge_orientation (const size_t edge_index) const
 Return the orientation of a Edge. More...
 
int Mesh2D::Polytope< object_dim >::vertex_orientation (const size_t vertex_index) const
 Return the orientation of a Vertex. More...
 
VectorRd Mesh2D::Polytope< object_dim >::normal () const
 Return the normal of a Face. More...
 
VectorRd Mesh2D::Polytope< object_dim >::tangent () const
 Return the tangent of a Edge. More...
 
void Mesh2D::Polytope< object_dim >::construct_face_normals ()
 Set the directions of the face normals of a cell. More...
 
void Mesh2D::Polytope< object_dim >::plot_simplices (std::ofstream *out) const
 Plot the simplices to out. More...
 
void Mesh2D::Polytope< object_dim >::plot (std::ofstream *out) const
 Plot the polytope to out. More...
 
template<size_t object_dim>
VectorRd Mesh2D::simplex_center_mass (Simplex< object_dim > simplex)
 Method to find the Lebesgue measure of an arbitrary simplex in arbitrary space. More...
 
template<size_t object_dim>
double Mesh2D::simplex_measure (Simplex< object_dim > simplex)
 

Detailed Description

Classes to construct and describe a 2D mesh.

Class providing a template for all objects (Vertex, Edge, Face, Cell) that appear in a Mesh.

Typedef Documentation

◆ Cell

using Mesh2D::Cell = typedef Polytope<DIMENSION>

◆ Edge

using Mesh2D::Edge = typedef Polytope<1>

A Face is a Polytope with object_dim = DIMENSION - 1.

◆ Face

using Mesh2D::Face = typedef Polytope<DIMENSION - 1>

A Cell is a Polytope with object_dim = DIMENSION.

◆ TransformationType

◆ Vertex

using Mesh2D::Vertex = typedef Polytope<0>

An Edge is a Polytope with object_dim = 1.

Function Documentation

◆ add_b_cell()

void Mesh2D::Mesh::add_b_cell ( Cell cell)
inline
Parameters
celladds a boundary cell to the mesh

◆ add_b_edge()

void Mesh2D::Mesh::add_b_edge ( Edge edge)
inline
Parameters
edgeadds a boundary edge to the mesh

◆ add_b_face()

void Mesh2D::Mesh::add_b_face ( Face face)
inline
Parameters
faceadds a boundary face to the mesh

◆ add_b_vertex()

void Mesh2D::Mesh::add_b_vertex ( Vertex vertex)
inline
Parameters
vertexadds a boundary vertex to the mesh

◆ add_cell() [1/2]

void Mesh2D::Mesh::add_cell ( Cell cell)
inline
Parameters
celladds a cell to the mesh

◆ add_cell() [2/2]

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::add_cell ( Polytope< DIMENSION > *  cell)

Add a cell to the Polytope.

◆ add_edge() [1/2]

void Mesh2D::Mesh::add_edge ( Edge edge)
inline
Parameters
edgeadds a edge to the mesh

◆ add_edge() [2/2]

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::add_edge ( Polytope< 1 > *  edge)

Add an edge to the Polytope.

◆ add_face() [1/2]

void Mesh2D::Mesh::add_face ( Face face)
inline
Parameters
faceadds a face to the mesh

◆ add_face() [2/2]

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::add_face ( Polytope< DIMENSION - 1 > *  face)

Add a face to the Polytope.

◆ add_i_cell()

void Mesh2D::Mesh::add_i_cell ( Cell cell)
inline
Parameters
celladds an internal cell to the mesh

◆ add_i_edge()

void Mesh2D::Mesh::add_i_edge ( Edge edge)
inline
Parameters
edgeadds an internal edge to the mesh

◆ add_i_face()

void Mesh2D::Mesh::add_i_face ( Face face)
inline
Parameters
faceadds an internal face to the mesh

◆ add_i_vertex()

void Mesh2D::Mesh::add_i_vertex ( Vertex vertex)
inline
Parameters
vertexadds an internal vertex to the mesh

◆ add_vertex() [1/2]

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::add_vertex ( Polytope< 0 > *  vertex)

Add a vertex to the Polytope.

◆ add_vertex() [2/2]

void Mesh2D::Mesh::add_vertex ( Vertex vertex)
inline
Parameters
vertexadds a vertex to the mesh

◆ b_cell()

Cell* Mesh2D::Mesh::b_cell ( std::size_t  index) const
inline

get a constant pointer to boundary a cell using an index

◆ b_edge()

Edge* Mesh2D::Mesh::b_edge ( std::size_t  index) const
inline

get a constant pointer to boundary a edge using an index

◆ b_face()

Face* Mesh2D::Mesh::b_face ( std::size_t  index) const
inline

get a constant pointer to boundary a face using an index

◆ b_vertex()

Vertex* Mesh2D::Mesh::b_vertex ( std::size_t  index) const
inline

get a constant pointer to a boundary vertex using an index

◆ build_the_mesh()

std::unique_ptr<Mesh> Mesh2D::MeshBuilder::build_the_mesh ( const TransformationType transformation = [](const std::array<double, 2> & x) { return x; })
inline

Build mesh

◆ cell() [1/2]

template<size_t object_dim>
Polytope< DIMENSION > * Mesh2D::Polytope< object_dim >::cell ( const size_t  i) const

Return the i-th cell of the Polytope.

◆ cell() [2/2]

Cell* Mesh2D::Mesh::cell ( std::size_t  index) const
inline

get a constant pointer to a cell using its global index

◆ center_mass()

template<size_t object_dim>
VectorRd Mesh2D::Polytope< object_dim >::center_mass ( ) const
inline

Return the center mass of the Polytope.

◆ construct_face_normals()

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::construct_face_normals

Set the directions of the face normals of a cell.

◆ coords()

template<size_t object_dim>
VectorRd Mesh2D::Polytope< object_dim >::coords

Return the coordinates of a Vertex.

◆ diam()

template<size_t object_dim>
double Mesh2D::Polytope< object_dim >::diam ( ) const
inline

Return the diameter of the Polytope.

◆ dim()

std::size_t Mesh2D::Mesh::dim ( ) const
inline

dimension of the mesh

◆ edge() [1/2]

template<size_t object_dim>
Polytope< 1 > * Mesh2D::Polytope< object_dim >::edge ( const size_t  i) const

Return the i-th edge of the Polytope.

◆ edge() [2/2]

Edge* Mesh2D::Mesh::edge ( std::size_t  index) const
inline

get a constant pointer to a edge using its global index

◆ edge_normal()

template<size_t object_dim>
VectorRd Mesh2D::Polytope< object_dim >::edge_normal ( const size_t  edge_index) const

Return the edge normal of a 2D object.

◆ edge_orientation()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::edge_orientation ( const size_t  edge_index) const

Return the orientation of a Edge.

◆ face() [1/2]

template<size_t object_dim>
Polytope< DIMENSION - 1 > * Mesh2D::Polytope< object_dim >::face ( const size_t  i) const

Return the i-th face of the Polytope.

◆ face() [2/2]

Face* Mesh2D::Mesh::face ( std::size_t  index) const
inline

get a constant pointer to a face using its global index

◆ face_normal()

template<size_t object_dim>
VectorRd Mesh2D::Polytope< object_dim >::face_normal ( const size_t  face_index) const

Return the outer normal of a Cell towards the Face located at face_index.

◆ face_orientation()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::face_orientation ( const size_t  face_index) const

Return the orientation of a Face.

◆ find_cell()

size_t Mesh2D::Mesh::find_cell ( const VectorRd  x) const
inline

< returns the index of the cell containing the point x

◆ get_b_cells()

std::vector<Cell *> Mesh2D::Mesh::get_b_cells ( ) const
inline

lists the boundary cells in the mesh.

◆ get_b_edges()

std::vector<Edge *> Mesh2D::Mesh::get_b_edges ( ) const
inline

lists the boundary edges in the mesh.

◆ get_b_faces()

std::vector<Face *> Mesh2D::Mesh::get_b_faces ( ) const
inline

lists the boundary faces in the mesh.

◆ get_b_vertices()

std::vector<Vertex *> Mesh2D::Mesh::get_b_vertices ( ) const
inline

lists the boundary vertices in the mesh.

◆ get_cells() [1/2]

std::vector<Cell *> Mesh2D::Mesh::get_cells ( ) const
inline

lists the cells in the mesh.

◆ get_cells() [2/2]

template<size_t object_dim>
std::vector<Polytope<DIMENSION> *> Mesh2D::Polytope< object_dim >::get_cells ( ) const
inline

Return the cells of the Polytope.

◆ get_edges() [1/2]

std::vector<Edge *> Mesh2D::Mesh::get_edges ( ) const
inline

lists the edges in the mesh.

◆ get_edges() [2/2]

template<size_t object_dim>
std::vector<Polytope<1> *> Mesh2D::Polytope< object_dim >::get_edges ( ) const
inline

Return the edges of the Polytope.

◆ get_faces() [1/2]

std::vector<Face *> Mesh2D::Mesh::get_faces ( ) const
inline

lists the faces in the mesh.

◆ get_faces() [2/2]

template<size_t object_dim>
std::vector<Polytope<DIMENSION - 1> *> Mesh2D::Polytope< object_dim >::get_faces ( ) const
inline

Return the faces of the Polytope.

◆ get_i_cells()

std::vector<Cell *> Mesh2D::Mesh::get_i_cells ( ) const
inline

lists the internal cells in the mesh.

◆ get_i_edges()

std::vector<Edge *> Mesh2D::Mesh::get_i_edges ( ) const
inline

lists the internal edges in the mesh.

◆ get_i_faces()

std::vector<Face *> Mesh2D::Mesh::get_i_faces ( ) const
inline

lists the internal faces in the mesh.

◆ get_i_vertices()

std::vector<Vertex *> Mesh2D::Mesh::get_i_vertices ( ) const
inline

lists the internal vertices in the mesh.

◆ get_name()

std::string Mesh2D::Mesh::get_name ( )
inline

getter for the mesh name

◆ get_simplices()

template<size_t object_dim>
Simplices<object_dim> Mesh2D::Polytope< object_dim >::get_simplices ( ) const
inline

Return the simplices making up the Polytope.

◆ get_vertices() [1/2]

std::vector<Vertex *> Mesh2D::Mesh::get_vertices ( ) const
inline

lists the vertices in the mesh.

◆ get_vertices() [2/2]

template<size_t object_dim>
std::vector<Polytope<0> *> Mesh2D::Polytope< object_dim >::get_vertices ( ) const
inline

Return the vertices of the Polytope.

◆ global_index()

template<size_t object_dim>
size_t Mesh2D::Polytope< object_dim >::global_index ( ) const
inline

Return the global index of the Polytope.

◆ h_max()

double Mesh2D::Mesh::h_max ( ) const
inline

< max diameter of cells

◆ i_cell()

Cell* Mesh2D::Mesh::i_cell ( std::size_t  index) const
inline

get a constant pointer to an internal cell using an index

◆ i_edge()

Edge* Mesh2D::Mesh::i_edge ( std::size_t  index) const
inline

get a constant pointer to an internal edge using an index

◆ i_face()

Face* Mesh2D::Mesh::i_face ( std::size_t  index) const
inline

get a constant pointer to an internal face using an index

◆ i_vertex()

Vertex* Mesh2D::Mesh::i_vertex ( std::size_t  index) const
inline

get a constant pointer to an internal vertex using an index

◆ index_cell()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::index_cell ( const Polytope< DIMENSION > *  cell) const

Returns the local index of a cell.

◆ index_edge()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::index_edge ( const Polytope< 1 > *  edge) const

Returns the local index of an edge.

◆ index_face()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::index_face ( const Polytope< DIMENSION - 1 > *  face) const

Returns the local index of a face.

◆ index_vertex()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::index_vertex ( const Polytope< 0 > *  vertex) const

Returns the local index of a vertex.

◆ is_boundary()

template<size_t object_dim>
bool Mesh2D::Polytope< object_dim >::is_boundary ( ) const
inline

Return true if Polytope is a boundary object, false otherwise.

◆ measure()

template<size_t object_dim>
double Mesh2D::Polytope< object_dim >::measure ( ) const
inline

Return the Lebesgue measure of the Polytope.

◆ Mesh()

Mesh2D::Mesh::Mesh ( )
inline

◆ MeshBuilder() [1/2]

Mesh2D::MeshBuilder::MeshBuilder ( )
inline

Constructor for MeshBuilder.

◆ MeshBuilder() [2/2]

Mesh2D::MeshBuilder::MeshBuilder ( const std::string  mesh_file)
inline

Overloaded constructor for MeshBuilder so read_mesh can be called from build_the_mesh().

◆ MeshReaderTyp2()

Mesh2D::MeshReaderTyp2::MeshReaderTyp2 ( std::string  file_name)
inline

Constructor for mesh reader

Parameters
file_namename of the file name, needs to include the full path

◆ n_b_cells()

std::size_t Mesh2D::Mesh::n_b_cells ( ) const
inline

number of boundary cells in the mesh.

◆ n_b_edges()

std::size_t Mesh2D::Mesh::n_b_edges ( ) const
inline

number of boundary edges in the mesh.

◆ n_b_faces()

std::size_t Mesh2D::Mesh::n_b_faces ( ) const
inline

number of boundary faces in the mesh.

◆ n_b_vertices()

std::size_t Mesh2D::Mesh::n_b_vertices ( ) const
inline

number of boundary vertices in the mesh.

◆ n_cells() [1/2]

std::size_t Mesh2D::Mesh::n_cells ( ) const
inline

number of cells in the mesh.

◆ n_cells() [2/2]

template<size_t object_dim>
size_t Mesh2D::Polytope< object_dim >::n_cells ( ) const
inline

Return the number of cells of the Polytope.

◆ n_edges() [1/2]

std::size_t Mesh2D::Mesh::n_edges ( ) const
inline

number of edges in the mesh.

◆ n_edges() [2/2]

template<size_t object_dim>
size_t Mesh2D::Polytope< object_dim >::n_edges ( ) const
inline

Return the number of edges of the Polytope.

◆ n_faces() [1/2]

std::size_t Mesh2D::Mesh::n_faces ( ) const
inline

number of faces in the mesh.

◆ n_faces() [2/2]

template<size_t object_dim>
size_t Mesh2D::Polytope< object_dim >::n_faces ( ) const
inline

Return the number of faces of the Polytope.

◆ n_i_cells()

std::size_t Mesh2D::Mesh::n_i_cells ( ) const
inline

number of internal cells in the mesh.

◆ n_i_edges()

std::size_t Mesh2D::Mesh::n_i_edges ( ) const
inline

number of internal edges in the mesh.

◆ n_i_faces()

std::size_t Mesh2D::Mesh::n_i_faces ( ) const
inline

number of internal faces in the mesh.

◆ n_i_vertices()

std::size_t Mesh2D::Mesh::n_i_vertices ( ) const
inline

number of internal vertices in the mesh.

◆ n_vertices() [1/2]

std::size_t Mesh2D::Mesh::n_vertices ( ) const
inline

number of vertices in the mesh.

◆ n_vertices() [2/2]

template<size_t object_dim>
size_t Mesh2D::Polytope< object_dim >::n_vertices ( ) const
inline

Return the number of vertices of the Polytope.

◆ normal()

template<size_t object_dim>
VectorRd Mesh2D::Polytope< object_dim >::normal

Return the normal of a Face.

◆ plot()

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::plot ( std::ofstream *  out) const

Plot the polytope to out.

◆ plot_mesh()

void Mesh2D::Mesh::plot_mesh ( std::ofstream *  out)
inline

◆ plot_simplices() [1/2]

void Mesh2D::Mesh::plot_simplices ( std::ofstream *  out)
inline

◆ plot_simplices() [2/2]

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::plot_simplices ( std::ofstream *  out) const

Plot the simplices to out.

◆ Polytope() [1/4]

template<size_t object_dim>
Mesh2D::Polytope< object_dim >::Polytope

Destructor.

◆ Polytope() [2/4]

template<size_t object_dim>
Mesh2D::Polytope< object_dim >::Polytope ( size_t  index,
Simplex< object_dim >  simplex 
)

Constructor for a Polytope defined by a single coordinate (Vertex)

◆ Polytope() [3/4]

template<size_t object_dim>
Mesh2D::Polytope< object_dim >::Polytope ( size_t  index,
Simplices< object_dim >  simplices 
)

< Constructor for a Polytope defined by its simplices

Constructor for a Simplex

◆ Polytope() [4/4]

template<size_t object_dim>
Mesh2D::Polytope< object_dim >::Polytope ( size_t  index,
VectorRd  vertex 
)

Null constructor.

◆ read_mesh()

void Mesh2D::MeshReaderTyp2::read_mesh ( std::vector< std::array< double, 2 >> &  vertices,
std::vector< std::vector< std::size_t >> &  cells 
)
inline

Reads the file into the specified containers

Parameters
verticesreference to a vector to hold the vertices coordinates
cellsreference to a vector to hold the cell vertices

◆ regularity()

std::vector<double> Mesh2D::Mesh::regularity ( )
inline

Regularity factor = 1st component: maximum of

  • diameter of cell / (measure of cell)^{1/dim}
  • diameter of cell / diameter of face [for each face of the cell]

2nd component: evaluation of max of ratio "diam of cell / radius ball inscribed in cell"

◆ renum()

void Mesh2D::Mesh::renum ( const char  B,
const std::vector< size_t >  new_to_old 
)
inline

◆ set_boundary()

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::set_boundary ( bool  val)
inline

Set the boundary value of the Polytope.

◆ set_global_index()

template<size_t object_dim>
void Mesh2D::Polytope< object_dim >::set_global_index ( const size_t  idx)
inline

Set the global index.

◆ set_name()

void Mesh2D::Mesh::set_name ( std::string  name)
inline

set the name of the mesh

◆ simplex_center_mass()

template<size_t object_dim>
VectorRd Mesh2D::simplex_center_mass ( Simplex< object_dim >  simplex)

Method to find the Lebesgue measure of an arbitrary simplex in arbitrary space.

◆ simplex_measure()

template<size_t object_dim>
double Mesh2D::simplex_measure ( Simplex< object_dim >  simplex)

◆ tangent()

template<size_t object_dim>
VectorRd Mesh2D::Polytope< object_dim >::tangent

Return the tangent of a Edge.

◆ vertex() [1/2]

template<size_t object_dim>
Polytope< 0 > * Mesh2D::Polytope< object_dim >::vertex ( const size_t  i) const

Return the i-th vertex of the Polytope.

◆ vertex() [2/2]

Vertex* Mesh2D::Mesh::vertex ( std::size_t  index) const
inline

get a constant pointer to a vertex using its global index

◆ vertex_orientation()

template<size_t object_dim>
int Mesh2D::Polytope< object_dim >::vertex_orientation ( const size_t  vertex_index) const

Return the orientation of a Vertex.

◆ ~Mesh()

Mesh2D::Mesh::~Mesh ( )
inline

◆ ~MeshReaderTyp2()

Mesh2D::MeshReaderTyp2::~MeshReaderTyp2 ( )
inline

◆ ~Polytope()

template<size_t object_dim>
Mesh2D::Polytope< object_dim >::~Polytope