diff --git a/src/mesh/MeshData.hpp b/src/mesh/MeshData.hpp index 6b600e4aa4fdf8be1e18a4ff0b46a4f809f29a34..9977c06a18f6dee5ea6443f4bc38f00190613e9f 100644 --- a/src/mesh/MeshData.hpp +++ b/src/mesh/MeshData.hpp @@ -4,7 +4,7 @@ #include <Kokkos_Core.hpp> #include <TinyVector.hpp> -#include <NodeByCellData.hpp> +#include <SubItemValuePerItem.hpp> #include <map> @@ -23,9 +23,9 @@ public: private: const MeshType& m_mesh; - NodeByCellData<Rd> m_Cjr; - NodeByCellData<double> m_ljr; - NodeByCellData<Rd> m_njr; + SubItemValuePerItem<Rd> m_Cjr; + SubItemValuePerItem<double> m_ljr; + SubItemValuePerItem<Rd> m_njr; Kokkos::View<Rd*> m_xj; Kokkos::View<double*> m_Vj; @@ -88,12 +88,12 @@ private: } }); - const NodeByCellData<Rd>& Cjr = m_Cjr; + const SubItemValuePerItem<Rd>& Cjr = m_Cjr; Kokkos::parallel_for(m_Cjr.numberOfValues(), KOKKOS_LAMBDA(const int& j){ m_ljr[j] = l2Norm(Cjr[j]); }); - const NodeByCellData<double>& ljr = m_ljr; + const SubItemValuePerItem<double>& ljr = m_ljr; Kokkos::parallel_for(m_Cjr.numberOfValues(), KOKKOS_LAMBDA(const int& j){ m_njr[j] = (1./ljr[j])*Cjr[j]; }); @@ -150,12 +150,12 @@ private: } }); - const NodeByCellData<Rd>& Cjr = m_Cjr; + const SubItemValuePerItem<Rd>& Cjr = m_Cjr; Kokkos::parallel_for(m_Cjr.numberOfValues(), KOKKOS_LAMBDA(const int& j){ m_ljr[j] = l2Norm(Cjr[j]); }); - const NodeByCellData<double>& ljr = m_ljr; + const SubItemValuePerItem<double>& ljr = m_ljr; Kokkos::parallel_for(m_Cjr.numberOfValues(), KOKKOS_LAMBDA(const int& j){ m_njr[j] = (1./ljr[j])*Cjr[j]; }); @@ -169,17 +169,17 @@ public: return m_mesh; } - const NodeByCellData<Rd>& Cjr() const + const SubItemValuePerItem<Rd>& Cjr() const { return m_Cjr; } - const NodeByCellData<double>& ljr() const + const SubItemValuePerItem<double>& ljr() const { return m_ljr; } - const NodeByCellData<Rd>& njr() const + const SubItemValuePerItem<Rd>& njr() const { return m_njr; } diff --git a/src/scheme/AcousticSolver.hpp b/src/scheme/AcousticSolver.hpp index b0037cb511a068f909a158f4590be7a73b40ef5e..b831ec7d0e5438dead4e9a6b7db69af3c4129074 100644 --- a/src/scheme/AcousticSolver.hpp +++ b/src/scheme/AcousticSolver.hpp @@ -14,7 +14,7 @@ #include <FiniteVolumesEulerUnknowns.hpp> #include <BoundaryCondition.hpp> -#include <NodeByCellData.hpp> +#include <SubItemValuePerItem.hpp> template<typename MeshData> class AcousticSolver @@ -82,9 +82,9 @@ private: KOKKOS_INLINE_FUNCTION void computeAjr(const Kokkos::View<const double*>& rhocj, - const NodeByCellData<Rd>& Cjr, - const NodeByCellData<double>& ljr, - const NodeByCellData<Rd>& njr) + const SubItemValuePerItem<Rd>& Cjr, + const SubItemValuePerItem<double>& ljr, + const SubItemValuePerItem<Rd>& njr) { Kokkos::parallel_for("new nested Ajr", m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) { const size_t& nb_nodes =m_Ajr.numberOfSubValues(j); @@ -97,7 +97,7 @@ private: KOKKOS_INLINE_FUNCTION const Kokkos::View<const Rdd*> - computeAr(const NodeByCellData<Rdd>& Ajr) { + computeAr(const SubItemValuePerItem<Rdd>& Ajr) { Kokkos::parallel_for(m_mesh.numberOfNodes(), KOKKOS_LAMBDA(const int& r) { Rdd sum = zero; const auto& node_to_cell = m_connectivity.m_node_to_cell_matrix.rowConst(r); @@ -115,8 +115,8 @@ private: KOKKOS_INLINE_FUNCTION const Kokkos::View<const Rd*> - computeBr(const NodeByCellData<Rdd>& Ajr, - const NodeByCellData<Rd>& Cjr, + computeBr(const SubItemValuePerItem<Rdd>& Ajr, + const SubItemValuePerItem<Rd>& Cjr, const Kokkos::View<const Rd*>& uj, const Kokkos::View<const double*>& pj) { Kokkos::parallel_for(m_mesh.numberOfNodes(), KOKKOS_LAMBDA(const int& r) { @@ -190,9 +190,9 @@ private: } void - computeFjr(const NodeByCellData<Rdd>& Ajr, + computeFjr(const SubItemValuePerItem<Rdd>& Ajr, const Kokkos::View<const Rd*>& ur, - const NodeByCellData<Rd>& Cjr, + const SubItemValuePerItem<Rd>& Cjr, const Kokkos::View<const Rd*>& uj, const Kokkos::View<const double*>& pj) { @@ -221,9 +221,9 @@ private: const Kokkos::View<const double*>& pj, const Kokkos::View<const double*>& cj, const Kokkos::View<const double*>& Vj, - const NodeByCellData<Rd>& Cjr, - const NodeByCellData<double>& ljr, - const NodeByCellData<Rd>& njr) { + const SubItemValuePerItem<Rd>& Cjr, + const SubItemValuePerItem<double>& ljr, + const SubItemValuePerItem<Rd>& njr) { const Kokkos::View<const double*> rhocj = computeRhoCj(rhoj, cj); computeAjr(rhocj, Cjr, ljr, njr); @@ -238,10 +238,10 @@ private: } Kokkos::View<Rd*> m_br; - NodeByCellData<Rdd> m_Ajr; + SubItemValuePerItem<Rdd> m_Ajr; Kokkos::View<Rdd*> m_Ar; Kokkos::View<Rdd*> m_inv_Ar; - NodeByCellData<Rd> m_Fjr; + SubItemValuePerItem<Rd> m_Fjr; Kokkos::View<Rd*> m_ur; Kokkos::View<double*> m_rhocj; Kokkos::View<double*> m_Vj_over_cj; @@ -270,7 +270,7 @@ public: double acoustic_dt(const Kokkos::View<const double*>& Vj, const Kokkos::View<const double*>& cj) const { - const NodeByCellData<double>& ljr = m_mesh_data.ljr(); + const SubItemValuePerItem<double>& ljr = m_mesh_data.ljr(); Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){ const auto& cell_nodes = m_mesh.connectivity().m_cell_to_node_matrix.rowConst(j); @@ -303,14 +303,14 @@ public: const Kokkos::View<const Rd*> xj = m_mesh_data.xj(); const Kokkos::View<const double*> Vj = m_mesh_data.Vj(); - const NodeByCellData<Rd>& Cjr = m_mesh_data.Cjr(); - const NodeByCellData<double>& ljr = m_mesh_data.ljr(); - const NodeByCellData<Rd>& njr = m_mesh_data.njr(); + const SubItemValuePerItem<Rd>& Cjr = m_mesh_data.Cjr(); + const SubItemValuePerItem<double>& ljr = m_mesh_data.ljr(); + const SubItemValuePerItem<Rd>& njr = m_mesh_data.njr(); Kokkos::View<Rd*> xr = m_mesh.xr(); computeExplicitFluxes(xr, xj, rhoj, uj, pj, cj, Vj, Cjr, ljr, njr); - const NodeByCellData<Rd>& Fjr = m_Fjr; + const SubItemValuePerItem<Rd>& Fjr = m_Fjr; const Kokkos::View<const Rd*> ur = m_ur; const Kokkos::View<const double*> inv_mj = unknowns.invMj(); diff --git a/src/scheme/NodeByCellData.hpp b/src/scheme/SubItemValuePerItem.hpp similarity index 91% rename from src/scheme/NodeByCellData.hpp rename to src/scheme/SubItemValuePerItem.hpp index 313858f881b51354938bc8c61a3767fa8225b3b0..102f046dbe41fc7fc26a6e861207a5c5edd36e14 100644 --- a/src/scheme/NodeByCellData.hpp +++ b/src/scheme/SubItemValuePerItem.hpp @@ -1,5 +1,5 @@ -#ifndef NODE_BY_CELL_DATA_HPP -#define NODE_BY_CELL_DATA_HPP +#ifndef SUBITEM_VALUE_PER_ITEM_HPP +#define SUBITEM_VALUE_PER_ITEM_HPP #include <Kokkos_StaticCrsGraph.hpp> @@ -8,7 +8,7 @@ typedef Kokkos::StaticCrsGraph<unsigned int, Kokkos::HostSpace> ConnectivityMatrix; template <typename DataType> -class NodeByCellData +class SubItemValuePerItem { #warning should eventually filter const from DataType private: @@ -96,7 +96,7 @@ class NodeByCellData } }; - NodeByCellData& operator=(const NodeByCellData&) = default; + SubItemValuePerItem& operator=(const SubItemValuePerItem&) = default; KOKKOS_FORCEINLINE_FUNCTION DataType& operator()(const size_t& j, const size_t& r) @@ -156,15 +156,15 @@ class NodeByCellData return SubViewConst(m_values, cell_begin, cell_end); } - NodeByCellData() = default; - NodeByCellData(const ConnectivityMatrix& node_id_per_cell_matrix) + SubItemValuePerItem() = default; + SubItemValuePerItem(const ConnectivityMatrix& node_id_per_cell_matrix) : m_node_id_per_cell_matrix(node_id_per_cell_matrix), m_values("values", m_node_id_per_cell_matrix.entries.extent(0)) { ; } - ~NodeByCellData() = default; + ~SubItemValuePerItem() = default; }; -#endif // NODE_BY_CELL_DATA_HPP +#endif // SUBITEM_VALUE_PER_ITEM_HPP