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

Implement the discrete spaces of DDR-PEC. More...

#include <pec.hpp>

Public Member Functions

 PEC (Mesh< dimension > const &mesh, int r, bool use_threads=true, std::array< int, dimension > const *dqr_p=nullptr, std::ostream &output=std::cout)
 
double getTopScaling (size_t i_cell) const
 Return the characteristic scale \(h_T\) of a top dimensional cell \(T\). More...
 
Eigen::MatrixXd get_mass (size_t k, size_t d, size_t i_cell) const
 Return the mass matrix for the k-forms on the i-th d-cell in the basis PL(r,k,d) More...
 
Eigen::MatrixXd get_trace (size_t k, size_t d, size_t i_cell, size_t j_bd) const
 Return the trace for the k-forms on the i-th d-cell onto its j-th (d-1)-neighbour. More...
 
Eigen::MatrixXd get_starTrace (size_t k, size_t d, size_t i_cell, size_t j_bd) const
 Return the Hodge dual of the trace for the Hodge dual of \((d-k)\)-forms on the i-th d-cell onto its j-th (d-1)-neighbour. More...
 
const Eigen::MatrixXd & get_diff (size_t l, size_t d) const
 Return the image of the differential operator. More...
 
const Eigen::MatrixXd & get_Koszul (size_t l, size_t d) const
 Return the image of the Koszul operator. More...
 
const Eigen::MatrixXd & get_diff_as_degr (size_t l, size_t d) const
 Return the image of the differential operator. More...
 
const Eigen::MatrixXd & get_trimmed (size_t l, size_t d) const
 Return the image of the trimmed polynomial basis. More...
 
const Eigen::MatrixXd & get_reduced_Koszul_m1 (size_t l, size_t d) const
 Return the image of the Koszul operator. More...
 

Detailed Description

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

Implement the discrete spaces of DDR-PEC.

This class compute and store the matrices of the differential operator, the Koszul operator, the trimmed polynomial basis, and the masses and traces of every cell.

The matrices for the differential, Koszul and trimmed basis are global and shared for all cells. The mass matrices and traces are computed and stored for each cell.

The basis used for the operator are given in exterior_algebra.hpp

Template Parameters
dimensionDimension of the manifold

Constructor & Destructor Documentation

◆ PEC()

template<size_t dimension>
PEC::PEC ( Mesh< dimension > const &  mesh,
int  r,
bool  use_threads = true,
std::array< int, dimension > const *  dqr_p = nullptr,
std::ostream &  output = std::cout 
)


Parameters
meshMesh to use
rPolynomial degree
use_threadsEnable pthreads parallelism
dqr_pDegree of quadrature used to generate the mass matrices. It cannot be exact for generic metric and default to a pretty high degree. Use a lower degree if the metric and cells are almost flat.
outputOutput stream for status messages.

Member Function Documentation

◆ get_diff()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::PEC< dimension >::get_diff ( size_t  l,
size_t  d 
) const
inline

Return the image of the differential operator.

Returns
\( \mathcal{P}_r\Lambda^l(\mathbb{R}^d) \rightarrow \mathcal{P}_{r-1}\Lambda^{l+1}(\mathbb{R}^{d}) \)
Parameters
lForm degree
dDimension

◆ get_diff_as_degr()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::PEC< dimension >::get_diff_as_degr ( size_t  l,
size_t  d 
) const
inline

Return the image of the differential operator.

Returns
\( \mathcal{P}_r\Lambda^l(\mathbb{R}^d) \rightarrow \mathcal{P}_{r}\Lambda^{l+1}(\mathbb{R}^{d}) \)
Parameters
lForm degree
dDimension

◆ get_Koszul()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::PEC< dimension >::get_Koszul ( size_t  l,
size_t  d 
) const
inline

Return the image of the Koszul operator.

Returns
\( \mathcal{P}_r\Lambda^l(\mathbb{R}^d) \rightarrow \mathcal{P}_{r+1}\Lambda^{l-1}(\mathbb{R}^{d}) \)
Parameters
lForm degree
dDimension

◆ get_mass()

template<size_t dimension>
Eigen::MatrixXd PEC::get_mass ( size_t  k,
size_t  d,
size_t  i_cell 
) const

Return the mass matrix for the k-forms on the i-th d-cell in the basis PL(r,k,d)

Returns
Mass matrix in the basis \( \mathcal{P}_r\Lambda^k(\mathbb{R}^d)\)
Parameters
kForm degree
dCell dimension
i_cellCell index

◆ get_reduced_Koszul_m1()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::PEC< dimension >::get_reduced_Koszul_m1 ( size_t  l,
size_t  d 
) const
inline

Return the image of the Koszul operator.

Returns
\( \mathcal{P}_{r-1}\Lambda^l(\mathbb{R}^d) \rightarrow \mathcal{P}_{r}\Lambda^{l-1}(\mathbb{R}^{d}) \)
Parameters
lForm degree
dDimension

◆ get_starTrace()

template<size_t dimension>
Eigen::MatrixXd PEC::get_starTrace ( size_t  k,
size_t  d,
size_t  i_cell,
size_t  j_bd 
) const

Return the Hodge dual of the trace for the Hodge dual of \((d-k)\)-forms on the i-th d-cell onto its j-th (d-1)-neighbour.

Returns
The matrix of the operator \(\star \text{tr} \star^{-1}\) in the basis \( \mathcal{P}_r\Lambda^{d-k}(\mathbb{R}^d) \rightarrow \mathcal{P}_r\Lambda^{d-1-k}(\mathbb{R}^{d-1}) \)
Parameters
kForm degree
dCell dimension
i_cellCell index
j_bdRelative index of the boundary element (e.g. between 0 and 2 for a triangle)

◆ get_trace()

template<size_t dimension>
Eigen::MatrixXd PEC::get_trace ( size_t  k,
size_t  d,
size_t  i_cell,
size_t  j_bd 
) const

Return the trace for the k-forms on the i-th d-cell onto its j-th (d-1)-neighbour.

Returns
Trace matrix in the basis \( \mathcal{P}_r\Lambda^k(\mathbb{R}^d) \rightarrow \mathcal{P}_r\Lambda^k(\mathbb{R}^{d-1}) \)
Parameters
kForm degree
dCell dimension
i_cellCell index
j_bdRelative index of the boundary element (e.g. between 0 and 2 for a triangle)

◆ get_trimmed()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::PEC< dimension >::get_trimmed ( size_t  l,
size_t  d 
) const
inline

Return the image of the trimmed polynomial basis.

Returns
\( \mathcal{P}_r^{-}\Lambda^l(\mathbb{R}^d) \rightarrow \mathcal{P}_{r}\Lambda^{l}(\mathbb{R}^{d}) \)
Parameters
lForm degree
dDimension

◆ getTopScaling()

template<size_t dimension>
double Manicore::PEC< dimension >::getTopScaling ( size_t  i_cell) const
inline

Return the characteristic scale \(h_T\) of a top dimensional cell \(T\).

Computed as the \(n\)th root of the volume of the cell

Parameters
i_cellCell index

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