Skip to content
Snippets Groups Projects
Commit fb83f355 authored by Stéphane Del Pino's avatar Stéphane Del Pino
Browse files

Add affine reconstruction for DiscreteFunctionP0 in 1, 2 and 3d

The possible data types are double, TinyVector and TinyMatrix
parent 017a617e
No related branches found
No related tags found
1 merge request!205High-order polynomial reconstruction
......@@ -40,10 +40,10 @@ class Givens
static void
_lift(const MatrixType& A, const RHSVectorType& b, VectorType& x)
{
for (ssize_t i = x.size() - 1; i >= 0; --i) {
for (ssize_t i = x.dimension() - 1; i >= 0; --i) {
const double inv_Aii = 1 / A(i, i);
double sum = 0;
for (size_t j = i + 1; j < x.size(); ++j) {
for (size_t j = i + 1; j < x.dimension(); ++j) {
sum += A(i, j) * x[j];
}
x[i] = (b[i] - sum) * inv_Aii;
......@@ -71,12 +71,13 @@ class Givens
static void
solveInPlace(MatrixType& A, VectorType& x, RHSVectorType& b)
{
Assert(x.size() == A.numberOfColumns(), "The number of columns of A must be the size of x");
Assert(b.size() == A.numberOfRows(), "The number of rows of A must be the size of b");
Assert(x.dimension() == A.numberOfColumns(), "The number of columns of A must be the size of x");
Assert(b.dimension() == A.numberOfRows(), "The number of rows of A must be the size of b");
for (size_t j = 0; j < A.numberOfColumns(); ++j) {
for (size_t i = A.numberOfRows() - 1; i > j; --i) {
double c(0), s(0);
double c;
double s;
_givens(A(i - 1, j), A(i, j), c, s);
for (size_t k = j; k < A.numberOfColumns(); ++k) {
double xi_1(0), xi(0);
......@@ -98,8 +99,8 @@ class Givens
static void
solve(const MatrixType& A, VectorType& x, const RHSVectorType& b)
{
Assert(x.size() == A.numberOfColumns(), "The number of columns of A must be the size of x");
Assert(b.size() == A.numberOfRows(), "The number of rows of A must be the size of b");
Assert(x.dimension() == A.numberOfColumns(), "The number of columns of A must be the size of x");
Assert(b.dimension() == A.numberOfRows(), "The number of rows of A must be the size of b");
MatrixType rotateA = copy(A);
RHSVectorType rotateb = copy(b);
......@@ -109,7 +110,7 @@ class Givens
template <typename MatrixType, typename UnknownMatrixType, typename RHSMatrixType>
static void
solveCollectionInPlace(const MatrixType& A, UnknownMatrixType& X, const RHSMatrixType& B)
solveCollectionInPlace(MatrixType& A, UnknownMatrixType& X, RHSMatrixType& B)
{
Assert(X.numberOfRows() == A.numberOfColumns(), "The number of columns of A must be the number of rows of X");
Assert(B.numberOfRows() == A.numberOfRows(), "The number of rows of A must be the rows of B");
......@@ -117,7 +118,8 @@ class Givens
for (size_t j = 0; j < A.numberOfColumns(); ++j) {
for (size_t i = A.numberOfRows() - 1; i > j; --i) {
double c(0), s(0);
double c;
double s;
_givens(A(i - 1, j), A(i, j), c, s);
for (size_t k = j; k < A.numberOfColumns(); ++k) {
double xi_1(0), xi(0);
......
......@@ -166,6 +166,13 @@ class SmallVector // LCOV_EXCL_LINE
return m_values.size();
}
PUGS_INLINE
size_t
dimension() const noexcept
{
return m_values.size();
}
PUGS_INLINE SmallVector&
fill(const DataType& value) noexcept
{
......
......@@ -10,6 +10,8 @@ add_library(
DiscreteFunctionVectorIntegrator.cpp
DiscreteFunctionVectorInterpoler.cpp
FluxingAdvectionSolver.cpp
test_reconstruction.cpp
)
target_link_libraries(
......
......@@ -20,55 +20,149 @@ class PolynomialReconstruction
{
public:
template <MeshConcept MeshType, typename DataType>
PUGS_INLINE DiscreteFunctionDPk<MeshType::Dimension, std::remove_const_t<DataType>>
[[nodiscard]] PUGS_INLINE DiscreteFunctionDPk<MeshType::Dimension, std::remove_const_t<DataType>>
build(const MeshType& mesh, const DiscreteFunctionP0<DataType> p0_function)
{
using Rd = TinyVector<MeshType::Dimension>;
const size_t degree = 1;
const Stencil stencil = StencilBuilder{}.build(mesh.connectivity());
auto xr = mesh.xr();
auto xj = MeshDataManager::instance().getMeshData(mesh).xj();
auto cell_is_owned = mesh.connectivity().cellIsOwned();
DiscreteFunctionDPk<MeshType::Dimension, std::remove_const_t<DataType>> dpk{mesh.meshVariant(), degree};
if constexpr ((is_polygonal_mesh_v<MeshType>)and(MeshType::Dimension == 1)) {
auto qf = QuadratureManager::instance().getLineFormula(GaussLegendreQuadratureDescriptor(1));
if constexpr (is_polygonal_mesh_v<MeshType>) {
if constexpr (std::is_arithmetic_v<DataType>) {
parallel_for(
mesh.numberOfCells(), PUGS_LAMBDA(const CellId cell_j_id) {
if (cell_is_owned[cell_j_id]) {
auto stencil_cell_list = stencil[cell_j_id];
SmallVector<double> b(stencil_cell_list.size());
const double p_j = p0_function[cell_j_id];
const double pj = p0_function[cell_j_id];
for (size_t i = 0; i < stencil_cell_list.size(); ++i) {
const CellId cell_i_id = stencil_cell_list[i];
b[i] = p0_function[cell_i_id] - p_j;
b[i] = p0_function[cell_i_id] - pj;
}
SmallMatrix<double> A(stencil_cell_list.size(), MeshType::Dimension);
const Rd& Xj = xj[cell_j_id];
for (size_t i = 0; i < stencil_cell_list.size(); ++i) {
const CellId cell_i_id = stencil_cell_list[i];
const Rd Xi_Xj = xj[cell_i_id] - Xj;
for (size_t l = 0; l < MeshType::Dimension; ++l) {
A(i, l) = Xi_Xj[l];
}
}
SmallMatrix<double> A(stencil_cell_list.size(), degree);
SmallVector<double> x(MeshType::Dimension);
Givens::solveInPlace(A, x, b);
using R1 = TinyVector<1>;
auto dpk_j = dpk.coefficients(cell_j_id);
const R1& Xj = xj[cell_j_id];
auto X_Xj = [&](const R1& X) { return X - Xj; };
dpk_j[0] = pj;
for (size_t l = 0; l < MeshType::Dimension; ++l) {
dpk_j[1 + l] = x[l];
}
}
});
} else if constexpr (is_tiny_vector_v<DataType>) {
parallel_for(
mesh.numberOfCells(), PUGS_LAMBDA(const CellId cell_j_id) {
if (cell_is_owned[cell_j_id]) {
auto stencil_cell_list = stencil[cell_j_id];
SmallMatrix<double> B(stencil_cell_list.size(), DataType::Dimension);
const DataType& pj = p0_function[cell_j_id];
for (size_t i = 0; i < stencil_cell_list.size(); ++i) {
const CellId cell_i_id = stencil_cell_list[i];
const DataType& pi_pj = p0_function[cell_i_id] - pj;
for (size_t k = 0; k < DataType::Dimension; ++k) {
B(i, k) = pi_pj[k];
}
}
SmallMatrix<double> A(stencil_cell_list.size(), MeshType::Dimension);
A(i, 0) = X_Xj(xj[cell_i_id])[0];
const Rd& Xj = xj[cell_j_id];
for (size_t i = 0; i < stencil_cell_list.size(); ++i) {
const CellId cell_i_id = stencil_cell_list[i];
const Rd Xi_Xj = xj[cell_i_id] - Xj;
for (size_t l = 0; l < MeshType::Dimension; ++l) {
A(i, l) = Xi_Xj[l];
}
}
SmallVector<double> x(1);
Givens::solveInPlace(A, x, b);
SmallMatrix<double> X(MeshType::Dimension, DataType::Dimension);
Givens::solveCollectionInPlace(A, X, B);
auto dpk_j = dpk.coefficients(cell_j_id);
dpk_j[0] = pj;
for (size_t i = 0; i < MeshType::Dimension; ++i) {
auto& dpk_j_ip1 = dpk_j[i + 1];
for (size_t k = 0; k < DataType::Dimension; ++k) {
dpk_j_ip1[k] = X(i, k);
}
}
}
});
} else if constexpr (is_tiny_matrix_v<DataType>) {
parallel_for(
mesh.numberOfCells(), PUGS_LAMBDA(const CellId cell_j_id) {
if (cell_is_owned[cell_j_id]) {
auto stencil_cell_list = stencil[cell_j_id];
SmallMatrix<double> B(stencil_cell_list.size(), DataType::Dimension);
dpk.coefficients(cell_j_id)[0] = p_j;
dpk.coefficients(cell_j_id)[1] = x[0];
const DataType& pj = p0_function[cell_j_id];
for (size_t i = 0; i < stencil_cell_list.size(); ++i) {
const CellId cell_i_id = stencil_cell_list[i];
const DataType& pi_pj = p0_function[cell_i_id] - pj;
for (size_t k = 0; k < DataType::NumberOfRows; ++k) {
for (size_t l = 0; l < DataType::NumberOfColumns; ++l) {
B(i, k * DataType::NumberOfColumns + l) = pi_pj(k, l);
}
}
}
SmallMatrix<double> A(stencil_cell_list.size(), MeshType::Dimension);
const Rd& Xj = xj[cell_j_id];
for (size_t i = 0; i < stencil_cell_list.size(); ++i) {
const CellId cell_i_id = stencil_cell_list[i];
const Rd Xi_Xj = xj[cell_i_id] - Xj;
for (size_t l = 0; l < MeshType::Dimension; ++l) {
A(i, l) = Xi_Xj[l];
}
}
SmallMatrix<double> X(MeshType::Dimension, DataType::Dimension);
Givens::solveCollectionInPlace(A, X, B);
auto dpk_j = dpk.coefficients(cell_j_id);
dpk_j[0] = pj;
for (size_t i = 0; i < MeshType::Dimension; ++i) {
auto& dpk_j_ip1 = dpk_j[i + 1];
for (size_t k = 0; k < DataType::NumberOfRows; ++k) {
for (size_t l = 0; l < DataType::NumberOfColumns; ++l) {
dpk_j_ip1(k, l) = X(i, k * DataType::NumberOfColumns + l);
}
}
}
}
});
} else {
throw NotImplementedError("dealing with " + demangle<DataType>());
}
}
synchronize(dpk.cellArrays());
......
......@@ -106,6 +106,12 @@ inline constexpr bool is_tiny_vector_v = false;
template <size_t N, typename T>
inline constexpr bool is_tiny_vector_v<TinyVector<N, T>> = true;
template <typename T>
inline constexpr bool is_tiny_vector_v<const T> = is_tiny_vector_v<std::remove_cvref_t<T>>;
template <typename T>
inline constexpr bool is_tiny_vector_v<T&> = is_tiny_vector_v<std::remove_cvref_t<T>>;
// Traits is_tiny_matrix
template <typename T>
......@@ -114,6 +120,12 @@ inline constexpr bool is_tiny_matrix_v = false;
template <size_t M, size_t N, typename T>
inline constexpr bool is_tiny_matrix_v<TinyMatrix<M, N, T>> = true;
template <typename T>
inline constexpr bool is_tiny_matrix_v<const T> = is_tiny_matrix_v<std::remove_cvref_t<T>>;
template <typename T>
inline constexpr bool is_tiny_matrix_v<T&> = is_tiny_matrix_v<std::remove_cvref_t<T>>;
// Trais is ItemValue
template <typename T>
......
This diff is collapsed.
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment