Manicore
Library to implement schemes on n-dimensionnal manifolds.
Public Member Functions | List of all members
Manicore::Mesh< dimension > Class Template Reference

Main data structure for the mesh. More...

#include <mesh.hpp>

Public Member Functions

size_t n_cells (size_t d) const
 Return the number of cell of a given dimension. More...
 
std::vector< size_t > const & get_map_ids (size_t d, size_t i_cell) const
 Return the map_ids of the i_cell-th d-cell. More...
 
std::vector< size_t > const & get_boundary (size_t d_boundary, size_t d, size_t i_cell) const
 Return a list filled with the global index of the elements on the boundary of the given dimension. More...
 
std::vector< size_t > const & get_relative_map (size_t d_boundary, size_t d, size_t i_cell) const
 The index for the mapping (of index 0) used on the element among the mappings of the boundary. More...
 
size_t global_to_local (size_t d_boundary, size_t d, size_t i_b, size_t i_cell) const
 Return the local index of the i_b d_boundary-cell with respect to the icell-th d-cell. More...
 
int get_boundary_orientation (size_t d, size_t i_cell, size_t j_bd) const
 Return the relative orientation of the j-th (d-1)-boundary of the i-th d-cell. More...
 
Eigen::Vector3d get_3D_embedding (size_t m_id, Eigen::Vector< double, dimension > const &x) const
 Interface to the extrinsic mapping. More...
 
Eigen::Matrix< double, 3, dimension > get_3D_pushforward (size_t m_id, Eigen::Vector< double, dimension > const &x) const
 Interface to the extrinsic mapping. More...
 
template<size_t d>
dCell_map< dimension, d > const & get_cell_map (size_t i) const
 Return the i-th dCell_map. More...
 
Eigen::Matrix< double, dimension, dimension > metric (size_t m_id, Eigen::Vector< double, dimension > const &x) const
 Evaluate the metric (of the tangent space) More...
 
Eigen::Matrix< double, dimension, dimension > metric_inv (size_t m_id, Eigen::Vector< double, dimension > const &x) const
 Evaluate the metric (of the cotangent space) More...
 
double volume_form (size_t m_id, Eigen::Vector< double, dimension > const &x) const
 Evaluate the scaling factor of the volume form. More...
 
int orientationTopCell (size_t i_cell) const
 Return the relative orientation of a top dimensional cell with respect to the manifold. More...
 
template<size_t k, int d>
Eigen::Matrix< double, Dimension::ExtDim(d-k, d), Dimension::ExtDim(k, d)> getHodge (size_t i_cell, Eigen::Vector< double, d > const &x) const
 Evaluate the hodge star operator. More...
 

Detailed Description

template<size_t dimension>
class Manicore::Mesh< dimension >

Main data structure for the mesh.

The class is intended to be build with Mesh_builder, there is no public constructor.

Store a collection of dCell_graph and dCell_map and keep the shared library describing the mesh.

A dcell is a sub manifold of dimension d. It is defined by its parametrization into some reference charts. Some cells appears in several charts, hence, they have several parametrizations. The map_ids property of a cell give the id of the chart corresponding to each parametrizations (e.g. map_ids[0] is the chart corresponding to the first parametrization ...).

Template Parameters
dimensionDimension of the manifold

Member Function Documentation

◆ get_3D_embedding()

template<size_t dimension>
Eigen::Vector3d Manicore::Mesh< dimension >::get_3D_embedding ( size_t  m_id,
Eigen::Vector< double, dimension > const &  x 
) const
inline

Interface to the extrinsic mapping.

Optional function to help viewing the embedding of the manifold. This is not used in any computation and can be set freely in the mesh shared library

Parameters
m_idId of the embedding to use
xLocation in the chart

◆ get_3D_pushforward()

template<size_t dimension>
Eigen::Matrix<double,3,dimension> Manicore::Mesh< dimension >::get_3D_pushforward ( size_t  m_id,
Eigen::Vector< double, dimension > const &  x 
) const
inline

Interface to the extrinsic mapping.

Optional function to help viewing the embedding of the manifold. This is not used in any computation and can be set freely in the mesh shared library

Using the differential of I should always map into the correct tangent space, whereas J could introduce a deviation if it is not constant along the normal component.

Parameters
m_idId of the embedding to use
xLocation in the chart

◆ get_boundary()

template<size_t dimension>
std::vector<size_t> const& Manicore::Mesh< dimension >::get_boundary ( size_t  d_boundary,
size_t  d,
size_t  i_cell 
) const
inline

Return a list filled with the global index of the elements on the boundary of the given dimension.

Parameters
d_boundaryDimension of the boundary cells
dDimension of the cell
i_cellIndex of the cell

◆ get_boundary_orientation()

template<size_t dimension>
int Manicore::Mesh< dimension >::get_boundary_orientation ( size_t  d,
size_t  i_cell,
size_t  j_bd 
) const

Return the relative orientation of the j-th (d-1)-boundary of the i-th d-cell.

Multiply by this factor when doing an integration on the boundary

Parameters
dDimension of the cell
i_cellIndex of the cell
j_bdRelative index of the boundary cell

◆ get_cell_map()

template<size_t dimension>
template<size_t d>
dCell_map<dimension,d> const& Manicore::Mesh< dimension >::get_cell_map ( size_t  i) const

Return the i-th dCell_map.

Template Parameters
dDimension of the cell

◆ get_map_ids()

template<size_t dimension>
std::vector<size_t> const& Manicore::Mesh< dimension >::get_map_ids ( size_t  d,
size_t  i_cell 
) const
inline

Return the map_ids of the i_cell-th d-cell.

Parameters
dDimension of the cell
i_cellIndex of the cell

◆ get_relative_map()

template<size_t dimension>
std::vector<size_t> const& Manicore::Mesh< dimension >::get_relative_map ( size_t  d_boundary,
size_t  d,
size_t  i_cell 
) const
inline

The index for the mapping (of index 0) used on the element among the mappings of the boundary.

The first parametrization of an element correspond to a given chart U, this function return the index of the parametrization of the element on its boundary corresponding to this chart U.

Parameters
d_boundaryDimension of the boundary cell
dDimension of the cell
i_cellIndex of the cell

◆ getHodge()

template<size_t dimension>
template<size_t k, int d>
Eigen::Matrix< double, Dimension::ExtDim(d-k, d), Dimension::ExtDim(k, d)> Manicore::Mesh< dimension >::getHodge ( size_t  i_cell,
Eigen::Vector< double, d > const &  x 
) const

Evaluate the hodge star operator.

Template Parameters
kForm degree
dDimension of the cell
Returns
Matrix from \( \Lambda^k(\mathbb{R}^d) \) to \( \Lambda^{d-k}(\mathbb{R}^d) \)
Parameters
i_cellIndex of the cell
xLocation on the reference element of the cell

◆ global_to_local()

template<size_t dimension>
size_t Manicore::Mesh< dimension >::global_to_local ( size_t  d_boundary,
size_t  d,
size_t  i_b,
size_t  i_cell 
) const
inline

Return the local index of the i_b d_boundary-cell with respect to the icell-th d-cell.

Search among the boundary elements of the i_cell-th d-cell for a d_boundary-cell with global index i_b. Return Mesh::get_boundary (d_boundary,d,i_cell).size() on failure (and abort when debugging).

Parameters
d_boundaryDimension of the boundary cell
dDimension of the cell
i_bGlobal index of the boundary cell
i_cellIndex of the cell

◆ metric()

template<size_t dimension>
Eigen::Matrix<double,dimension,dimension> Manicore::Mesh< dimension >::metric ( size_t  m_id,
Eigen::Vector< double, dimension > const &  x 
) const
inline

Evaluate the metric (of the tangent space)

Parameters
m_idId of the chart to use
xLocation in the chart

◆ metric_inv()

template<size_t dimension>
Eigen::Matrix<double,dimension,dimension> Manicore::Mesh< dimension >::metric_inv ( size_t  m_id,
Eigen::Vector< double, dimension > const &  x 
) const
inline

Evaluate the metric (of the cotangent space)

Parameters
m_idId of the chart to use
xLocation in the chart

◆ n_cells()

template<size_t dimension>
size_t Manicore::Mesh< dimension >::n_cells ( size_t  d) const
inline

Return the number of cell of a given dimension.

Parameters
dDimension of the cell

◆ orientationTopCell()

template<size_t dimension>
int Manicore::Mesh< dimension >::orientationTopCell ( size_t  i_cell) const
inline

Return the relative orientation of a top dimensional cell with respect to the manifold.

◆ volume_form()

template<size_t dimension>
double Manicore::Mesh< dimension >::volume_form ( size_t  m_id,
Eigen::Vector< double, dimension > const &  x 
) const
inline

Evaluate the scaling factor of the volume form.

Parameters
m_idId of the chart to use
xLocation in the chart

The documentation for this class was generated from the following file: