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

Add ReduceMin and ReduceMax for any array type

parent b8cbc261
Branches
Tags
1 merge request!7Feature/itemvalue
...@@ -34,7 +34,7 @@ class ItemValue ...@@ -34,7 +34,7 @@ class ItemValue
} }
KOKKOS_INLINE_FUNCTION KOKKOS_INLINE_FUNCTION
size_t numberOfValues() const size_t size() const
{ {
Assert(m_is_built); Assert(m_is_built);
return m_values.size(); return m_values.size();
......
#ifndef ACOUSTIC_SOLVER_HPP #ifndef ACOUSTIC_SOLVER_HPP
#define ACOUSTIC_SOLVER_HPP #define ACOUSTIC_SOLVER_HPP
#include <Kokkos_Core.hpp>
#include <rang.hpp> #include <rang.hpp>
#include <ArrayUtils.hpp>
#include <BlockPerfectGas.hpp> #include <BlockPerfectGas.hpp>
#include <PastisAssert.hpp> #include <PastisAssert.hpp>
...@@ -33,43 +34,6 @@ class AcousticSolver ...@@ -33,43 +34,6 @@ class AcousticSolver
typedef TinyMatrix<dimension> Rdd; typedef TinyMatrix<dimension> Rdd;
private: private:
struct ReduceMin
{
private:
const CellValue<const double> x_;
public:
typedef std::remove_const_t<CellValue<const double>::data_type> value_type;
ReduceMin(const CellValue<const double>& x) : x_ (x) {}
#warning was: "typedef Kokkos::View<const double*>::size_type size_type;""
typedef size_t size_type;
KOKKOS_INLINE_FUNCTION
void operator() (const size_type i, value_type& update) const
{
if (x_[i] < update) {
update = x_[i];
}
}
KOKKOS_INLINE_FUNCTION
void join (volatile value_type& dst,
const volatile value_type& src) const
{
if (src < dst) {
dst = src;
}
}
KOKKOS_INLINE_FUNCTION void
init (value_type& dst) const
{ // The identity under max is -Inf.
dst= Kokkos::reduction_identity<value_type>::min();
}
};
KOKKOS_INLINE_FUNCTION KOKKOS_INLINE_FUNCTION
const CellValue<const double> const CellValue<const double>
computeRhoCj(const CellValue<const double>& rhoj, computeRhoCj(const CellValue<const double>& rhoj,
...@@ -308,10 +272,7 @@ class AcousticSolver ...@@ -308,10 +272,7 @@ class AcousticSolver
m_Vj_over_cj[j] = 2*Vj[j]/(S*cj[j]); m_Vj_over_cj[j] = 2*Vj[j]/(S*cj[j]);
}); });
double dt = std::numeric_limits<double>::max(); return ReduceMin(m_Vj_over_cj);
Kokkos::parallel_reduce(m_mesh.numberOfCells(), ReduceMin(m_Vj_over_cj), dt);
return dt;
} }
......
...@@ -29,7 +29,7 @@ public: ...@@ -29,7 +29,7 @@ public:
void updatePandCFromRhoE() void updatePandCFromRhoE()
{ {
const size_t nj = m_ej.numberOfValues(); const size_t nj = m_ej.size();
const CellValue<const double> rho = m_rhoj; const CellValue<const double> rho = m_rhoj;
const CellValue<const double> e = m_ej; const CellValue<const double> e = m_ej;
const CellValue<const double> gamma = m_gammaj; const CellValue<const double> gamma = m_gammaj;
...@@ -43,7 +43,7 @@ public: ...@@ -43,7 +43,7 @@ public:
void updateEandCFromRhoP() void updateEandCFromRhoP()
{ {
const size_t nj = m_ej.numberOfValues(); const size_t nj = m_ej.size();
const CellValue<const double> rho = m_rhoj; const CellValue<const double> rho = m_rhoj;
const CellValue<const double> p = m_pj; const CellValue<const double> p = m_pj;
const CellValue<const double> gamma = m_gammaj; const CellValue<const double> gamma = m_gammaj;
......
#ifndef ARRAY_UTILS_HPP
#define ARRAY_UTILS_HPP
#include <Kokkos_Core.hpp>
template <typename ArrayType>
class ReduceMin
{
private:
const ArrayType& m_array;
using data_type = std::remove_const_t<typename ArrayType::data_type>;
public:
KOKKOS_INLINE_FUNCTION
operator data_type()
{
data_type reduced_value;
Kokkos::parallel_reduce(m_array.size(), *this, reduced_value);
return reduced_value;
}
KOKKOS_INLINE_FUNCTION
void operator()(const size_t& i, data_type& data) const
{
if (m_array[i] < data) {
data = m_array[i];
}
}
KOKKOS_INLINE_FUNCTION
void join(volatile data_type& dst,
const volatile data_type& src) const
{
if (src < dst) {
dst = src;
}
}
KOKKOS_INLINE_FUNCTION void
init(data_type& value) const
{
value = std::numeric_limits<data_type>::max();
}
KOKKOS_INLINE_FUNCTION
ReduceMin(const ArrayType& array)
: m_array(array)
{
;
}
KOKKOS_INLINE_FUNCTION
~ReduceMin() = default;
};
template <typename ArrayType>
class ReduceMax
{
private:
const ArrayType& m_array;
using data_type = std::remove_const_t<typename ArrayType::data_type>;
public:
KOKKOS_INLINE_FUNCTION
operator data_type()
{
data_type reduced_value;
Kokkos::parallel_reduce(m_array.size(), *this, reduced_value);
return reduced_value;
}
KOKKOS_INLINE_FUNCTION
void operator()(const size_t& i, data_type& data) const
{
if (m_array[i] > data) {
data = m_array[i];
}
}
KOKKOS_INLINE_FUNCTION
void join(volatile data_type& dst,
const volatile data_type& src) const
{
if (src > dst) {
dst = src;
}
}
KOKKOS_INLINE_FUNCTION void
init(data_type& value) const
{
value = -std::numeric_limits<data_type>::max();
}
KOKKOS_INLINE_FUNCTION
ReduceMax(const ArrayType& array)
: m_array(array)
{
;
}
KOKKOS_INLINE_FUNCTION
~ReduceMax() = default;
};
#endif //ARRAY_UTILS_HPP
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment