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

Define NodeValuePerCell

NodeValuePerCell<DataType> is an alias to SubItemValuePerItem<Datatype,node,cell>
parent e24dbb69
No related branches found
No related tags found
1 merge request!6Feature/crs
......@@ -23,9 +23,9 @@ public:
private:
const MeshType& m_mesh;
SubItemValuePerItem<Rd> m_Cjr;
SubItemValuePerItem<double> m_ljr;
SubItemValuePerItem<Rd> m_njr;
NodeValuePerCell<Rd> m_Cjr;
NodeValuePerCell<double> m_ljr;
NodeValuePerCell<Rd> m_njr;
Kokkos::View<Rd*> m_xj;
Kokkos::View<double*> m_Vj;
......@@ -88,12 +88,12 @@ private:
}
});
const SubItemValuePerItem<Rd>& Cjr = m_Cjr;
const NodeValuePerCell<Rd>& Cjr = m_Cjr;
Kokkos::parallel_for(m_Cjr.numberOfValues(), KOKKOS_LAMBDA(const int& j){
m_ljr[j] = l2Norm(Cjr[j]);
});
const SubItemValuePerItem<double>& ljr = m_ljr;
const NodeValuePerCell<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 SubItemValuePerItem<Rd>& Cjr = m_Cjr;
const NodeValuePerCell<Rd>& Cjr = m_Cjr;
Kokkos::parallel_for(m_Cjr.numberOfValues(), KOKKOS_LAMBDA(const int& j){
m_ljr[j] = l2Norm(Cjr[j]);
});
const SubItemValuePerItem<double>& ljr = m_ljr;
const NodeValuePerCell<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 SubItemValuePerItem<Rd>& Cjr() const
const NodeValuePerCell<Rd>& Cjr() const
{
return m_Cjr;
}
const SubItemValuePerItem<double>& ljr() const
const NodeValuePerCell<double>& ljr() const
{
return m_ljr;
}
const SubItemValuePerItem<Rd>& njr() const
const NodeValuePerCell<Rd>& njr() const
{
return m_njr;
}
......
......@@ -82,9 +82,9 @@ private:
KOKKOS_INLINE_FUNCTION
void computeAjr(const Kokkos::View<const double*>& rhocj,
const SubItemValuePerItem<Rd>& Cjr,
const SubItemValuePerItem<double>& ljr,
const SubItemValuePerItem<Rd>& njr)
const NodeValuePerCell<Rd>& Cjr,
const NodeValuePerCell<double>& ljr,
const NodeValuePerCell<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 SubItemValuePerItem<Rdd>& Ajr) {
computeAr(const NodeValuePerCell<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 SubItemValuePerItem<Rdd>& Ajr,
const SubItemValuePerItem<Rd>& Cjr,
computeBr(const NodeValuePerCell<Rdd>& Ajr,
const NodeValuePerCell<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 SubItemValuePerItem<Rdd>& Ajr,
computeFjr(const NodeValuePerCell<Rdd>& Ajr,
const Kokkos::View<const Rd*>& ur,
const SubItemValuePerItem<Rd>& Cjr,
const NodeValuePerCell<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 SubItemValuePerItem<Rd>& Cjr,
const SubItemValuePerItem<double>& ljr,
const SubItemValuePerItem<Rd>& njr) {
const NodeValuePerCell<Rd>& Cjr,
const NodeValuePerCell<double>& ljr,
const NodeValuePerCell<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;
SubItemValuePerItem<Rdd> m_Ajr;
NodeValuePerCell<Rdd> m_Ajr;
Kokkos::View<Rdd*> m_Ar;
Kokkos::View<Rdd*> m_inv_Ar;
SubItemValuePerItem<Rd> m_Fjr;
NodeValuePerCell<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 SubItemValuePerItem<double>& ljr = m_mesh_data.ljr();
const NodeValuePerCell<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 SubItemValuePerItem<Rd>& Cjr = m_mesh_data.Cjr();
const SubItemValuePerItem<double>& ljr = m_mesh_data.ljr();
const SubItemValuePerItem<Rd>& njr = m_mesh_data.njr();
const NodeValuePerCell<Rd>& Cjr = m_mesh_data.Cjr();
const NodeValuePerCell<double>& ljr = m_mesh_data.ljr();
const NodeValuePerCell<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 SubItemValuePerItem<Rd>& Fjr = m_Fjr;
const NodeValuePerCell<Rd>& Fjr = m_Fjr;
const Kokkos::View<const Rd*> ur = m_ur;
const Kokkos::View<const double*> inv_mj = unknowns.invMj();
......
......@@ -7,12 +7,21 @@
typedef Kokkos::StaticCrsGraph<unsigned int, Kokkos::HostSpace> ConnectivityMatrix;
template <typename DataType>
enum class TypeOfItem {
node,
edge,
face,
cell
};
template <typename DataType,
TypeOfItem SubItemType,
TypeOfItem ItemType>
class SubItemValuePerItem
{
#warning should eventually filter const from DataType
private:
ConnectivityMatrix m_node_id_per_cell_matrix;
ConnectivityMatrix m_subitem_id_per_item_matrix;
Kokkos::View<DataType*> m_values;
......@@ -101,13 +110,13 @@ class SubItemValuePerItem
KOKKOS_FORCEINLINE_FUNCTION
DataType& operator()(const size_t& j, const size_t& r)
{
return m_values[m_node_id_per_cell_matrix.row_map[j]+r];
return m_values[m_subitem_id_per_item_matrix.row_map[j]+r];
}
KOKKOS_FORCEINLINE_FUNCTION
const DataType& operator()(const size_t& j, const size_t& r) const
{
return m_values[m_node_id_per_cell_matrix.row_map[j]+r];
return m_values[m_subitem_id_per_item_matrix.row_map[j]+r];
}
KOKKOS_INLINE_FUNCTION
......@@ -129,37 +138,37 @@ class SubItemValuePerItem
}
KOKKOS_INLINE_FUNCTION
size_t numberOfEntities() const
size_t numberOfItems() const
{
return m_node_id_per_cell_matrix.numRows();
return m_subitem_id_per_item_matrix.numRows();
}
KOKKOS_INLINE_FUNCTION
size_t numberOfSubValues(const size_t& i_cell) const
{
return m_node_id_per_cell_matrix.row_map[i_cell+1]-m_node_id_per_cell_matrix.row_map[i_cell];
return m_subitem_id_per_item_matrix.row_map[i_cell+1]-m_subitem_id_per_item_matrix.row_map[i_cell];
}
KOKKOS_INLINE_FUNCTION
SubView cellValues(const size_t& i_cell)
SubView itemValues(const size_t& i_cell)
{
const ConnectivityMatrix::size_type& cell_begin = m_node_id_per_cell_matrix.row_map[i_cell];
const ConnectivityMatrix::size_type& cell_end = m_node_id_per_cell_matrix.row_map[i_cell+1];
const ConnectivityMatrix::size_type& cell_begin = m_subitem_id_per_item_matrix.row_map[i_cell];
const ConnectivityMatrix::size_type& cell_end = m_subitem_id_per_item_matrix.row_map[i_cell+1];
return SubView(m_values, cell_begin, cell_end);
}
KOKKOS_INLINE_FUNCTION
SubViewConst cellValues(const size_t& i_cell) const
SubViewConst itemValues(const size_t& i_cell) const
{
const ConnectivityMatrix::size_type& cell_begin = m_node_id_per_cell_matrix.row_map[i_cell];
const ConnectivityMatrix::size_type& cell_end = m_node_id_per_cell_matrix.row_map[i_cell+1];
const ConnectivityMatrix::size_type& cell_begin = m_subitem_id_per_item_matrix.row_map[i_cell];
const ConnectivityMatrix::size_type& cell_end = m_subitem_id_per_item_matrix.row_map[i_cell+1];
return SubViewConst(m_values, cell_begin, cell_end);
}
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))
SubItemValuePerItem(const ConnectivityMatrix& subitem_id_per_item_matrix)
: m_subitem_id_per_item_matrix(subitem_id_per_item_matrix),
m_values("values", m_subitem_id_per_item_matrix.entries.extent(0))
{
;
}
......@@ -167,4 +176,8 @@ class SubItemValuePerItem
~SubItemValuePerItem() = default;
};
template <typename DataType>
using NodeValuePerCell = SubItemValuePerItem<DataType, TypeOfItem::node, TypeOfItem::cell>;
#endif // SUBITEM_VALUE_PER_ITEM_HPP
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment