| Manicore
    Library to implement schemes on n-dimensionnal manifolds. | 
Manage the geometry of a cell. More...
#include <dcell.hpp>

| Public Member Functions | |
| dCell_map (Eigen::Vector< double, dimension > const ¢er_mass, Eigen::Matrix< double, dimension, d > const &flat_map, double diam, std::vector< Simplex< d >> const &triangulation) | |
| Flat constructor.  More... | |
| dCell_map (std::vector< std::unique_ptr< ParametrizedMap< dimension, d >>> &maps, std::vector< std::unique_ptr< ParametrizedDerivedMap< dimension, d >>> &pullback_maps, std::vector< Simplex< d >> const &triangulation) | |
| General constructor.  More... | |
| bool | is_flat () const | 
| Is the cell flat.  More... | |
| std::vector< Simplex< d > > const & | get_reference_elem () const | 
| Return the reference element as a collection of simplexes.  More... | |
| Eigen::Vector< double, dimension > | evaluate_I (size_t rel_map_id, Eigen::Vector< double, d > const &x) const | 
| Evaluate the parametrization from the reference element to a chart.  More... | |
| Eigen::Vector< double, d > | evaluate_J (size_t rel_map_id, Eigen::Vector< double, dimension > const &x) const | 
| Evaluate the inverse mapping from the chart to the reference element.  More... | |
| double | evaluate_poly_on_ref (Eigen::Vector< double, d > const &x, size_t i_pbasis, int r) const | 
| Evaluate a scalar polynomial on the reference element.  More... | |
| double | evaluate_poly_pullback (size_t rel_map_id, Eigen::Vector< double, dimension > const &x, size_t i_pbasis, int r) const | 
| Evaluate the pullback of a scalar polynomial on the chart.  More... | |
| Eigen::Matrix< double, dimension, d > | evaluate_DI (size_t rel_map_id, Eigen::Vector< double, d > const &x) const | 
| Evaluate the differential of the parametrization.  More... | |
| Eigen::Matrix< double, d, dimension > | evaluate_DJ (size_t rel_map_id, Eigen::Vector< double, dimension > const &x) const | 
| Evaluate the differential of the inverse mapping.  More... | |
| template<size_t l> | |
| Eigen::Matrix< double, Dimension::ExtDim(l, d), Dimension::ExtDim(l, dimension)> | evaluate_DI_p (size_t rel_map_id, Eigen::Vector< double, d > const &x) const | 
| Compute the action of the pullback of l-forms by the parametrization.  More... | |
| template<size_t l> | |
| Eigen::Matrix< double, Dimension::ExtDim(l, dimension), Dimension::ExtDim(l, d)> | evaluate_DJ_p (size_t rel_map_id, Eigen::Vector< double, dimension > const &x) const | 
| Compute the action of the pullback of l-forms by the inverse mapping.  More... | |
| int | get_orientation (size_t rel_map_id, Eigen::Vector< double, dimension > const &x, Eigen::Matrix< double, dimension, d-1 > const &pM) const requires(d >1) | 
| Return the relative orientation.  More... | |
| Static Public Attributes | |
| static constexpr size_t | cell_dim = d | 
| Dimension of the cell  More... | |
Manage the geometry of a cell.
| dimension | Dimension of the manifold | 
| d | Dimension of the cell Hold the reference element, compute the mappings and pullback between the reference element and the charts. | 
The class is specialize to apply some optimization when dealing with flat cells.
| Manicore::dCell_map< dimension, d >::dCell_map | ( | Eigen::Vector< double, dimension > const & | center_mass, | 
| Eigen::Matrix< double, dimension, d > const & | flat_map, | ||
| double | diam, | ||
| std::vector< Simplex< d >> const & | triangulation | ||
| ) | 
Flat constructor.
Simplify the structure when dealing with flat cells. This must be a flat surface in a chart parametrized by an affine mapping, and all its boundary element must also be flat. As an additional restriction, a flat cell can only be in a single chart
| center_mass | Centroid | 
| flat_map | Basis of the tangent space | 
| diam | Scaling factor to apply | 
| triangulation | Reference element given as a collection of simplexes | 
| Manicore::dCell_map< dimension, d >::dCell_map | ( | std::vector< std::unique_ptr< ParametrizedMap< dimension, d >>> & | maps, | 
| std::vector< std::unique_ptr< ParametrizedDerivedMap< dimension, d >>> & | pullback_maps, | ||
| std::vector< Simplex< d >> const & | triangulation | ||
| ) | 
General constructor.
| maps | List of the parametrization of this cell | 
| pullback_maps | Differentials of the parametrization of this cell | 
| triangulation | Reference element given as a collection of simplexes | 
| Eigen::Matrix< double, dimension, d > Manicore::dCell_map< dimension, d >::evaluate_DI | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, d > const & | x | ||
| ) | const | 
Evaluate the differential of the parametrization.
| rel_map_id | Relative id of the chart to use | 
| x | Location on the reference element | 
| Eigen::Matrix< double, Dimension::ExtDim(l, d), Dimension::ExtDim(l, dimension)> Manicore::dCell_map< dimension, d >::evaluate_DI_p | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, d > const & | x | ||
| ) | const | 
Compute the action of the pullback of l-forms by the parametrization.
| l | Form degree | 
| rel_map_id | Relative id of the chart to use | 
| x | Location on the reference element | 
| Eigen::Matrix< double, d, dimension > Manicore::dCell_map< dimension, d >::evaluate_DJ | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, dimension > const & | x | ||
| ) | const | 
Evaluate the differential of the inverse mapping.
| rel_map_id | Relative id of the chart to use | 
| x | Location on the chart | 
| Eigen::Matrix< double, Dimension::ExtDim(l, dimension), Dimension::ExtDim(l, d)> Manicore::dCell_map< dimension, d >::evaluate_DJ_p | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, dimension > const & | x | ||
| ) | const | 
Compute the action of the pullback of l-forms by the inverse mapping.
| l | Form degree | 
| rel_map_id | Relative id of the chart to use | 
| x | Location on the chart | 
| Eigen::Vector< double, dimension > Manicore::dCell_map< dimension, d >::evaluate_I | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, d > const & | x | ||
| ) | const | 
Evaluate the parametrization from the reference element to a chart.
| rel_map_id | Relative id of the chart to use | 
| x | Location on the reference element | 
| Eigen::Vector< double, d > Manicore::dCell_map< dimension, d >::evaluate_J | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, dimension > const & | x | ||
| ) | const | 
Evaluate the inverse mapping from the chart to the reference element.
| rel_map_id | Relative id of the chart to use | 
| x | Location on the chart | 
| double Manicore::dCell_map< dimension, d >::evaluate_poly_on_ref | ( | Eigen::Vector< double, d > const & | x, | 
| size_t | i_pbasis, | ||
| int | r | ||
| ) | const | 
Evaluate a scalar polynomial on the reference element.
| x | Location on the reference element | 
| i_pbasis | Index of the polynomial to evaluate | 
| r | Polynomial basis | 
| double Manicore::dCell_map< dimension, d >::evaluate_poly_pullback | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, dimension > const & | x, | ||
| size_t | i_pbasis, | ||
| int | r | ||
| ) | const | 
Evaluate the pullback of a scalar polynomial on the chart.
| rel_map_id | Relative id of the chart to use | 
| x | Location on the chart | 
| i_pbasis | Index of the polynomial to evaluate | 
| r | Polynomial basis | 
| int Manicore::dCell_map< dimension, d >::get_orientation | ( | size_t | rel_map_id, | 
| Eigen::Vector< double, dimension > const & | x, | ||
| Eigen::Matrix< double, dimension, d-1 > const & | pM | ||
| ) | const | 
Return the relative orientation.
Given a point on the boundary x, it attempt to find which simplex of the reference element S contains x, then compute a outward normal vector subtracting the center of S to x. The outward normal is then compared with a basis of the tangent space of the boundary to get the orientation of the boundary.
| rel_map_id | Relative id of the chart to use | 
| x | Any point on the boundary in the chart (avoid vertices) | 
| pM | Matrix of a vector basis of the tangent space of the boundary | 
| 
 | inline | 
Return the reference element as a collection of simplexes.
| 
 | inline | 
Is the cell flat.
Only check that the cell was constructed with the flat constructor, does not check if the parametrization is flat
| 
 | staticconstexpr | 
Dimension of the cell