Manicore
Library to implement schemes on n-dimensionnal manifolds.
|
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... | |
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
dimension | Dimension of the manifold |
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 |
||
) |
mesh | Mesh to use |
r | Polynomial degree |
use_threads | Enable pthreads parallelism |
dqr_p | Degree 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. |
output | Output stream for status messages. |
|
inline |
Return the image of the differential operator.
l | Form degree |
d | Dimension |
|
inline |
Return the image of the differential operator.
l | Form degree |
d | Dimension |
|
inline |
Return the image of the Koszul operator.
l | Form degree |
d | 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)
k | Form degree |
d | Cell dimension |
i_cell | Cell index |
|
inline |
Return the image of the Koszul operator.
l | Form degree |
d | 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.
k | Form degree |
d | Cell dimension |
i_cell | Cell index |
j_bd | Relative index of the boundary element (e.g. between 0 and 2 for a triangle) |
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.
k | Form degree |
d | Cell dimension |
i_cell | Cell index |
j_bd | Relative index of the boundary element (e.g. between 0 and 2 for a triangle) |
|
inline |
Return the image of the trimmed polynomial basis.
l | Form degree |
d | Dimension |
|
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
i_cell | Cell index |