The Reconstruction class allows to build and apply a polynomial reconstructions. More...
Public Member Functions | |
Reconstruction (uint8_t degree, uint8_t dimensions) | |
void | assemble () |
void | clear (bool release=true) |
uint16_t | getCoefficientCount () const |
uint8_t | getDegree () const |
uint8_t | getDimensions () const |
void | initialize (uint8_t degree, uint8_t dimensions, bool release=true) |
void | swap (Reconstruction &other) noexcept |
Public Member Functions inherited from bitpit::ReconstructionKernel | |
ReconstructionKernel () | |
ReconstructionKernel (const ReconstructionKernel &other) | |
ReconstructionKernel (ReconstructionKernel &&other)=default | |
ReconstructionKernel (uint8_t degree, uint8_t dimensions, int nEquations) | |
void | assemblePolynomial (const std::array< double, 3 > &origin, const double *values, ReconstructionPolynomial *polynomial) const |
void | assemblePolynomial (const std::array< double, 3 > &origin, int nFields, const double **values, ReconstructionPolynomial *polynomial) const |
void | assemblePolynomial (uint8_t degree, const std::array< double, 3 > &origin, const double *values, ReconstructionPolynomial *polynomial) const |
void | assemblePolynomial (uint8_t degree, const std::array< double, 3 > &origin, int nFields, const double **values, ReconstructionPolynomial *polynomial) const |
void | computeDerivativeLimitedWeights (const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, const double *limiters, double *derivativeWeights) const |
void | computeDerivativeLimitedWeights (int nEquations, uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, const double *limiters, double *derivativeWeights) const |
void | computeDerivativeLimitedWeights (uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, const double *limiters, double *derivativeWeights) const |
void | computeDerivativeWeights (const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double *derivativeWeights) const |
void | computeDerivativeWeights (int nEquations, uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double *derivativeWeights) const |
void | computeDerivativeWeights (uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double *derivativeWeights) const |
void | computeGradientLimitedWeights (const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, std::array< double, 3 > *gradientWeights) const |
void | computeGradientLimitedWeights (int nEquations, uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, std::array< double, 3 > *gradientWeights) const |
void | computeGradientLimitedWeights (uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, std::array< double, 3 > *gradientWeights) const |
void | computeGradientWeights (const std::array< double, 3 > &origin, const std::array< double, 3 > &point, std::array< double, 3 > *gradientWeights) const |
void | computeGradientWeights (int nEquations, uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, std::array< double, 3 > *gradientWeights) const |
void | computeGradientWeights (uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, std::array< double, 3 > *gradientWeights) const |
void | computeValueLimitedWeights (const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, double *valueWeights) const |
void | computeValueLimitedWeights (int nEquations, uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, double *valueWeights) const |
void | computeValueLimitedWeights (uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, double *valueWeights) const |
void | computeValueWeights (const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double *valueWeights) const |
void | computeValueWeights (int nEquations, uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double *valueWeights) const |
void | computeValueWeights (uint8_t degree, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double *valueWeights) const |
void | display (std::ostream &out, double tolerance=1.e-10) const |
uint16_t | getCoefficientCount () const |
uint8_t | getDegree () const |
uint8_t | getDimensions () const |
int | getEquationCount () const |
double * | getPolynomialWeights () |
const double * | getPolynomialWeights () const |
void | initialize (uint8_t degree, uint8_t dimensions, int nEquations, bool release=true) |
ReconstructionKernel & | operator= (const ReconstructionKernel &other) |
ReconstructionKernel & | operator= (ReconstructionKernel &&other)=default |
void | swap (ReconstructionKernel &other) noexcept |
void | updatePolynomial (const double *values, ReconstructionPolynomial *polynomial) const |
void | updatePolynomial (int nFields, const double **values, ReconstructionPolynomial *polynomial) const |
void | updatePolynomial (uint8_t degree, const double *values, ReconstructionPolynomial *polynomial) const |
void | updatePolynomial (uint8_t degree, int nFields, const double **values, ReconstructionPolynomial *polynomial) const |
Public Member Functions inherited from bitpit::ReconstructionAssembler | |
ReconstructionAssembler () | |
ReconstructionAssembler (uint8_t degree, uint8_t dimensions) | |
void | addCellAverageEquation (ReconstructionType type, const Cell &cell, const std::array< double, 3 > &origin, const std::array< double, 3 > *vertexCoords, double scaleFactor=1.) |
void | addPointDerivativeEquation (ReconstructionType type, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double scaleFactor=1.) |
void | addPointValueEquation (ReconstructionType type, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double scaleFactor=1.) |
void | assembleKernel (ReconstructionKernel *kernel) const |
void | clear (bool release=true) |
int | countConstraints () const |
int | countEquations () const |
int | countLeastSquares () const |
uint16_t | getCoefficientCount () const |
uint8_t | getDegree () const |
uint8_t | getDimensions () const |
void | initialize (uint8_t degree, uint8_t dimensions, bool release=true) |
void | swap (ReconstructionAssembler &other) noexcept |
void | updateKernel (ReconstructionKernel *kernel) const |
Additional Inherited Members | |
Public Types inherited from bitpit::ReconstructionAssembler | |
enum | ReconstructionType { TYPE_CONSTRAINT , TYPE_LEAST_SQUARE } |
Protected Member Functions inherited from bitpit::ReconstructionKernel | |
void | applyLimiter (uint8_t degree, const double *limiters, double *coeffs) const |
The Reconstruction class allows to build and apply a polynomial reconstructions.
Definition at line 307 of file reconstruction.hpp.
bitpit::Reconstruction::Reconstruction | ( | uint8_t | degree, |
uint8_t | dimensions ) |
Default constructor.
degree | is the degree of the polynomial |
dimensions | is the number of space dimensions |
Definition at line 3620 of file reconstruction.cpp.
void bitpit::Reconstruction::assemble | ( | ) |
Computes the weights to be used in order to calculate the coefficients of the polynomial. The coefficients are such that the conditions decoded in equations are enforced.
Definition at line 3673 of file reconstruction.cpp.
void bitpit::Reconstruction::clear | ( | bool | release = true | ) |
Clear the reconstruction.
release | if true, the memory hold by the reconstruction will be released, otherwise the reconstruction will be cleared but its memory will not be released |
Definition at line 3659 of file reconstruction.cpp.
uint16_t bitpit::ReconstructionAssembler::getCoefficientCount | ( | ) | const |
Get the number of coefficients of the reconstruction polynomial.
Definition at line 264 of file reconstruction.cpp.
uint8_t bitpit::ReconstructionAssembler::getDegree | ( | ) | const |
Get the degree of the polynomial.
Definition at line 261 of file reconstruction.cpp.
uint8_t bitpit::ReconstructionAssembler::getDimensions | ( | ) | const |
Get the number of space dimensions.
Definition at line 262 of file reconstruction.cpp.
void bitpit::Reconstruction::initialize | ( | uint8_t | degree, |
uint8_t | dimensions, | ||
bool | release = true ) |
Initialize the reconstruction.
degree | is the degree of the polynomial |
dimensions | is the number of space dimensions |
release | if true, possible unneeded memory hold by the reconstruction will be released, otherwise the reconstruction will be initialized but possible unneeded memory will not be released |
Definition at line 3647 of file reconstruction.cpp.
|
noexcept |
Exchanges the content of the reconstruction by the content the specified other reconstruction.
other | is another reconstruction whose content is swapped with that of this reconstruction |
Definition at line 3632 of file reconstruction.cpp.