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

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

#include <ddr_spaces.hpp>

Inheritance diagram for Manicore::DDR_Spaces< dimension >:
Inheritance graph
[legend]

Public Types

template<size_t k>
using FunctionType = std::function< Eigen::Vector< double, Dimension::ExtDim(k, dimension)>(size_t map_id, const Eigen::Vector< double, dimension > &)>
 Signature of the function to use in interpolate. More...
 

Public Member Functions

 DDR_Spaces (Mesh< dimension > const &mesh, int r, bool use_threads=true, std::array< int, dimension > const *dqr_p=nullptr, std::ostream &output=std::cout)
 
template<size_t k>
Eigen::VectorXd interpolate (FunctionType< k > const &func, std::array< int, dimension > const *dqr_p=nullptr) const
 Interpolate the given function on the discrete spaces. More...
 
const Eigen::MatrixXd & full_diff (size_t k, size_t d, size_t i_cell) const
 Compute the full differential operator. More...
 
Eigen::MatrixXd compose_diff (size_t k, size_t d, size_t i_cell) const
 Compute the projected differential operator. More...
 
const Eigen::MatrixXd & potential (size_t k, size_t d, size_t i_cell) const
 Compute the potential operator. More...
 
DOFSpace< dimension > const & dofspace (size_t k) const
 Return DOFSpace associated to k forms. More...
 
const Mesh< dimension > * mesh () const
 Return the Mesh associated with the class. More...
 
int degree () const
 Return the polynomial degree associated with the class. More...
 
Eigen::MatrixXd computeL2Product (size_t k, size_t i_cell) const
 Compute the local \(L^2\) product (including the cell boundary) More...
 
double hmax () const
 Return the mesh size. More...
 

Detailed Description

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

Implement the discrete operators of DDR-PEC.

Interpolator, discrete differential and potential reconstruction.

This class compute store a matrix for each operator on each cell. It also create and store the corresponding PEC object (holding the masses and traces of the cells)

Template Parameters
dimensionDimension of the manifold
Warning
This does not take ownership of the mesh but keep a pointer of it. Ensure that the mesh survives this class.

Member Typedef Documentation

◆ FunctionType

template<size_t dimension>
template<size_t k>
using Manicore::DDR_Spaces< dimension >::FunctionType = std::function<Eigen::Vector<double,Dimension::ExtDim(k,dimension)>(size_t map_id, const Eigen::Vector<double,dimension> &)>

Signature of the function to use in interpolate.

Template Parameters
kForm degree

Constructor & Destructor Documentation

◆ DDR_Spaces()

template<size_t dimension>
DDR_Spaces::DDR_Spaces ( Mesh< dimension > const &  mesh,
int  r,
bool  use_threads = true,
std::array< int, dimension > const *  dqr_p = nullptr,
std::ostream &  output = std::cout 
)
Warning
Ensure that the mesh survives this class
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

◆ compose_diff()

template<size_t dimension>
Eigen::MatrixXd DDR_Spaces::compose_diff ( size_t  k,
size_t  d,
size_t  i_cell 
) const

Compute the projected differential operator.

Returns
\(\star \underline{d}_r^k\) on a cell including its boundary in \(\cup_{d' = k}^{d} \mathcal{P}_r\Lambda^{d-k-1}(\mathbb{R}^{d'})\)
Parameters
kForm degree
dCell dimension
i_cellCell index

◆ computeL2Product()

template<size_t dimension>
Eigen::MatrixXd DDR_Spaces::computeL2Product ( size_t  k,
size_t  i_cell 
) const

Compute the local \(L^2\) product (including the cell boundary)

The contribution on the cell is \(\int_f \langle P^k , P^k \rangle \text{vol}_f \). The contribution from a cell \(f' \in \partial f \) of the boundary is \(\int_f \langle \text{tr}_{f'} P^k_f - P^k_{f'} , \text{tr}_{f'} P^k_f - P^k_{f'} \rangle \text{vol}_{f'} \) multiplied by a scaling factor.

Parameters
kForm degree
i_cellCell index

◆ degree()

template<size_t dimension>
int Manicore::DDR_Spaces< dimension >::degree ( ) const
inline

Return the polynomial degree associated with the class.

◆ dofspace()

template<size_t dimension>
DOFSpace<dimension> const& Manicore::DDR_Spaces< dimension >::dofspace ( size_t  k) const
inline

Return DOFSpace associated to k forms.

Parameters
kForm degree

◆ full_diff()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::DDR_Spaces< dimension >::full_diff ( size_t  k,
size_t  d,
size_t  i_cell 
) const
inline

Compute the full differential operator.

Returns
\(\star d_r^k \) in \(\mathcal{P}_r\Lambda^{d-k-1}(\mathbb{R}^d)\)
Parameters
kForm degree
dCell dimension
i_cellCell index

◆ hmax()

template<size_t dimension>
double DDR_Spaces::hmax

Return the mesh size.

Returns
\(\max_{f\in \delta_{\text{dimension}}} h_f\)

◆ interpolate()

template<size_t dimension>
template<size_t k>
Eigen::VectorXd DDR_Spaces::interpolate ( FunctionType< k > const &  func,
std::array< int, dimension > const *  dqr_p = nullptr 
) const

Interpolate the given function on the discrete spaces.

Template Parameters
kForm degree
Parameters
funcFunction to interpolate
dqr_pQuadrature degree to use

◆ mesh()

template<size_t dimension>
const Mesh<dimension>* Manicore::DDR_Spaces< dimension >::mesh ( ) const
inline

Return the Mesh associated with the class.

◆ potential()

template<size_t dimension>
const Eigen::MatrixXd& Manicore::DDR_Spaces< dimension >::potential ( size_t  k,
size_t  d,
size_t  i_cell 
) const
inline

Compute the potential operator.

Returns
\(\star P_r^k\) in \(\mathcal{P}_r\Lambda^{d-k}(\mathbb{R}^d)\)
Parameters
kForm degree
dCell dimension
i_cellCell index

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