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

Changed all static size entries in AcousticSolverWithMesh by dynamic

ones (preparing for multi-d)
parent f0b5a23a
Branches
Tags
2 merge requests!2Develop,!1Develop
...@@ -72,11 +72,14 @@ private: ...@@ -72,11 +72,14 @@ private:
} }
KOKKOS_INLINE_FUNCTION KOKKOS_INLINE_FUNCTION
const Kokkos::View<const Rdd*[2]> const Kokkos::View<const Rdd**>
computeAjr(const Kokkos::View<const double*>& rhocj, computeAjr(const Kokkos::View<const double*>& rhocj,
const Kokkos::View<const Rd*[2]>& Cjr) { const Kokkos::View<const Rd**>& Cjr) {
const Kokkos::View<const unsigned short*> cell_nb_nodes
= m_connectivity.cellNbNodes();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) { Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
for (int r=0; r<2; ++r) { for (int r=0; r<cell_nb_nodes[j]; ++r) {
m_Ajr(j,r) = tensorProduct(rhocj(j)*Cjr(j,r), Cjr(j,r)); m_Ajr(j,r) = tensorProduct(rhocj(j)*Cjr(j,r), Cjr(j,r));
} }
}); });
...@@ -86,9 +89,9 @@ private: ...@@ -86,9 +89,9 @@ private:
KOKKOS_INLINE_FUNCTION KOKKOS_INLINE_FUNCTION
const Kokkos::View<const Rdd*> const Kokkos::View<const Rdd*>
computeAr(const Kokkos::View<const Rdd*[2]>& Ajr) { computeAr(const Kokkos::View<const Rdd**>& Ajr) {
const Kokkos::View<const unsigned int*[2]> node_cells = m_connectivity.nodeCells(); const Kokkos::View<const unsigned int**> node_cells = m_connectivity.nodeCells();
const Kokkos::View<const unsigned short*[2]> node_cell_local_node = m_connectivity.nodeCellLocalNode(); const Kokkos::View<const unsigned short**> node_cell_local_node = m_connectivity.nodeCellLocalNode();
const Kokkos::View<const unsigned short*> node_nb_cells = m_connectivity.nodeNbCells(); const Kokkos::View<const unsigned short*> node_nb_cells = m_connectivity.nodeNbCells();
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) { Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
...@@ -106,12 +109,12 @@ private: ...@@ -106,12 +109,12 @@ private:
KOKKOS_INLINE_FUNCTION KOKKOS_INLINE_FUNCTION
const Kokkos::View<const Rd*> const Kokkos::View<const Rd*>
computeBr(const Kokkos::View<const Rdd*[2]>& Ajr, computeBr(const Kokkos::View<const Rdd**>& Ajr,
const Kokkos::View<const Rd*[2]>& Cjr, const Kokkos::View<const Rd**>& Cjr,
const Kokkos::View<const Rd*>& uj, const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const double*>& pj) { const Kokkos::View<const double*>& pj) {
const Kokkos::View<const unsigned int*[2]>& node_cells = m_connectivity.nodeCells(); const Kokkos::View<const unsigned int**>& node_cells = m_connectivity.nodeCells();
const Kokkos::View<const unsigned short*[2]>& node_cell_local_node = m_connectivity.nodeCellLocalNode(); const Kokkos::View<const unsigned short**>& node_cell_local_node = m_connectivity.nodeCellLocalNode();
const Kokkos::View<const unsigned short*>& node_nb_cells = m_connectivity.nodeNbCells(); const Kokkos::View<const unsigned short*>& node_nb_cells = m_connectivity.nodeNbCells();
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) { Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
...@@ -141,15 +144,18 @@ private: ...@@ -141,15 +144,18 @@ private:
return m_ur; return m_ur;
} }
Kokkos::View<Rd*[2]> Kokkos::View<Rd**>
computeFjr(const Kokkos::View<const Rdd*[2]>& Ajr, computeFjr(const Kokkos::View<const Rdd**>& Ajr,
const Kokkos::View<const Rd*>& ur, const Kokkos::View<const Rd*>& ur,
const Kokkos::View<const Rd*[2]>& Cjr, const Kokkos::View<const Rd**>& Cjr,
const Kokkos::View<const Rd*>& uj, const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const double*>& pj) { const Kokkos::View<const double*>& pj) {
const Kokkos::View<const unsigned int*[2]>& cell_nodes = m_connectivity.cellNodes(); const Kokkos::View<const unsigned int**>& cell_nodes = m_connectivity.cellNodes();
const Kokkos::View<const unsigned short*> cell_nb_nodes
= m_connectivity.cellNbNodes();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) { Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
for (int r=0; r<2; ++r) { for (int r=0; r<cell_nb_nodes[j]; ++r) {
m_Fjr(j,r) = Ajr(j,r)*(uj(j)-ur(cell_nodes(j,r)))+pj(j)*Cjr(j,r); m_Fjr(j,r) = Ajr(j,r)*(uj(j)-ur(cell_nodes(j,r)))+pj(j)*Cjr(j,r);
} }
}); });
...@@ -192,15 +198,15 @@ private: ...@@ -192,15 +198,15 @@ private:
const Kokkos::View<const double*>& pj, const Kokkos::View<const double*>& pj,
const Kokkos::View<const double*>& cj, const Kokkos::View<const double*>& cj,
const Kokkos::View<const double*>& Vj, const Kokkos::View<const double*>& Vj,
const Kokkos::View<const Rd*[2]>& Cjr, const Kokkos::View<const Rd**>& Cjr) {
Kokkos::View<Rd*>& ur,
Kokkos::View<Rd*[2]>& Fjr) {
const Kokkos::View<const double*> rhocj = computeRhoCj(rhoj, cj); const Kokkos::View<const double*> rhocj = computeRhoCj(rhoj, cj);
const Kokkos::View<const Rdd*[2]> Ajr = computeAjr(rhocj, Cjr); const Kokkos::View<const Rdd**> Ajr = computeAjr(rhocj, Cjr);
const Kokkos::View<const Rdd*> Ar = computeAr(Ajr); const Kokkos::View<const Rdd*> Ar = computeAr(Ajr);
const Kokkos::View<const Rd*> br = computeBr(Ajr, Cjr, uj, pj); const Kokkos::View<const Rd*> br = computeBr(Ajr, Cjr, uj, pj);
Kokkos::View<Rd*> ur = m_ur;
Kokkos::View<Rd**> Fjr = m_Fjr;
ur = computeUr(Ar, br); ur = computeUr(Ar, br);
Fjr = computeFjr(Ajr, ur, Cjr, uj, pj); Fjr = computeFjr(Ajr, ur, Cjr, uj, pj);
} }
...@@ -209,10 +215,10 @@ private: ...@@ -209,10 +215,10 @@ private:
const size_t m_nr; const size_t m_nr;
Kokkos::View<Rd*> m_br; Kokkos::View<Rd*> m_br;
Kokkos::View<Rdd*[2]> m_Ajr; Kokkos::View<Rdd**> m_Ajr;
Kokkos::View<Rdd*> m_Ar; Kokkos::View<Rdd*> m_Ar;
Kokkos::View<Rdd*> m_inv_Ar; Kokkos::View<Rdd*> m_inv_Ar;
Kokkos::View<Rd*[2]> m_Fjr; Kokkos::View<Rd**> m_Fjr;
Kokkos::View<Rd*> m_ur; Kokkos::View<Rd*> m_ur;
Kokkos::View<double*> m_rhocj; Kokkos::View<double*> m_rhocj;
Kokkos::View<double*> m_Vj_over_cj; Kokkos::View<double*> m_Vj_over_cj;
...@@ -225,10 +231,10 @@ public: ...@@ -225,10 +231,10 @@ public:
m_nj(mesh.numberOfCells()), m_nj(mesh.numberOfCells()),
m_nr(mesh.numberOfNodes()), m_nr(mesh.numberOfNodes()),
m_br("br", m_nr), m_br("br", m_nr),
m_Ajr("Ajr", m_nj), m_Ajr("Ajr", m_nj, m_connectivity.maxNbNodePerCell()),
m_Ar("Ar", m_nr), m_Ar("Ar", m_nr),
m_inv_Ar("inv_Ar", m_nr), m_inv_Ar("inv_Ar", m_nr),
m_Fjr("Fjr", m_nj), m_Fjr("Fjr", m_nj, m_connectivity.maxNbNodePerCell()),
m_ur("ur", m_nr), m_ur("ur", m_nr),
m_rhocj("rho_c", m_nj), m_rhocj("rho_c", m_nj),
m_Vj_over_cj("Vj_over_cj", m_nj) m_Vj_over_cj("Vj_over_cj", m_nj)
...@@ -248,7 +254,7 @@ public: ...@@ -248,7 +254,7 @@ public:
Kokkos::View<Rd*> xr = mesh.xr(); Kokkos::View<Rd*> xr = mesh.xr();
const Kokkos::View<const unsigned int*[2]>& cell_nodes = m_connectivity.cellNodes(); const Kokkos::View<const unsigned int**>& cell_nodes = m_connectivity.cellNodes();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){ Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
xj[j] = 0.5*(xr[cell_nodes(j,0)]+xr[cell_nodes(j,1)]); xj[j] = 0.5*(xr[cell_nodes(j,0)]+xr[cell_nodes(j,1)]);
}); });
...@@ -318,16 +324,17 @@ public: ...@@ -318,16 +324,17 @@ public:
t=tmax; t=tmax;
} }
computeExplicitFluxes(xr, xj, rhoj, uj, pj, cj, Vj, Cjr, computeExplicitFluxes(xr, xj, rhoj, uj, pj, cj, Vj, Cjr);
m_ur, m_Fjr);
const Kokkos::View<const Rd*[2]> Fjr = m_Fjr; const Kokkos::View<const Rd**> Fjr = m_Fjr;
const Kokkos::View<const Rd*> ur = m_ur; const Kokkos::View<const Rd*> ur = m_ur;
const Kokkos::View<const unsigned short*> cell_nb_nodes
= m_connectivity.cellNbNodes();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) { Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
Rd momentum_fluxes = zero; Rd momentum_fluxes = zero;
double energy_fluxes = 0; double energy_fluxes = 0;
for (int R=0; R<2; ++R) { for (int R=0; R<cell_nb_nodes[j]; ++R) {
const int r=cell_nodes(j,R); const int r=cell_nodes(j,R);
momentum_fluxes += Fjr(j,R); momentum_fluxes += Fjr(j,R);
energy_fluxes += (Fjr(j,R), ur[r]); energy_fluxes += (Fjr(j,R), ur[r]);
......
...@@ -17,6 +17,9 @@ private: ...@@ -17,6 +17,9 @@ private:
Kokkos::View<unsigned int*[2]> m_cell_nodes; Kokkos::View<unsigned int*[2]> m_cell_nodes;
const Kokkos::View<unsigned int*[2]>& m_cell_faces; const Kokkos::View<unsigned int*[2]>& m_cell_faces;
Kokkos::View<unsigned short*> m_cell_nb_nodes;
const Kokkos::View<unsigned short*> m_cell_nb_faces;
Kokkos::View<unsigned short*> m_node_nb_cells; Kokkos::View<unsigned short*> m_node_nb_cells;
const Kokkos::View<unsigned short*>& m_face_nb_cells; const Kokkos::View<unsigned short*>& m_face_nb_cells;
...@@ -26,6 +29,8 @@ private: ...@@ -26,6 +29,8 @@ private:
Kokkos::View<unsigned short*[2]> m_node_cell_local_node; Kokkos::View<unsigned short*[2]> m_node_cell_local_node;
const Kokkos::View<unsigned short*[2]>& m_face_cell_local_face; const Kokkos::View<unsigned short*[2]>& m_face_cell_local_face;
const size_t m_max_nb_node_per_cell;
public: public:
const size_t& numberOfNodes() const const size_t& numberOfNodes() const
{ {
...@@ -42,12 +47,17 @@ public: ...@@ -42,12 +47,17 @@ public:
return m_number_of_cells; return m_number_of_cells;
} }
const Kokkos::View<const unsigned int*[2]> cellNodes() const const size_t& maxNbNodePerCell() const
{
return m_max_nb_node_per_cell;
}
const Kokkos::View<const unsigned int**> cellNodes() const
{ {
return m_cell_nodes; return m_cell_nodes;
} }
const Kokkos::View<const unsigned int*[2]> cellFaces() const const Kokkos::View<const unsigned int**> cellFaces() const
{ {
return m_cell_faces; return m_cell_faces;
} }
...@@ -57,27 +67,37 @@ public: ...@@ -57,27 +67,37 @@ public:
return m_node_nb_cells; return m_node_nb_cells;
} }
const Kokkos::View<const unsigned short*> cellNbNodes() const
{
return m_cell_nb_nodes;
}
const Kokkos::View<const unsigned short*> cellNbFaces() const
{
return m_cell_nb_faces;
}
const Kokkos::View<const unsigned short*> faceNbCells() const const Kokkos::View<const unsigned short*> faceNbCells() const
{ {
return m_face_nb_cells; return m_face_nb_cells;
} }
const Kokkos::View<const unsigned int*[2]> nodeCells() const const Kokkos::View<const unsigned int**> nodeCells() const
{ {
return m_node_cells; return m_node_cells;
} }
const Kokkos::View<const unsigned int*[2]> faceCells() const const Kokkos::View<const unsigned int**> faceCells() const
{ {
return m_face_cells; return m_face_cells;
} }
const Kokkos::View<const unsigned short*[2]> nodeCellLocalNode() const const Kokkos::View<const unsigned short**> nodeCellLocalNode() const
{ {
return m_node_cell_local_node; return m_node_cell_local_node;
} }
const Kokkos::View<const unsigned short*[2]> faceCellLocalFace() const const Kokkos::View<const unsigned short**> faceCellLocalFace() const
{ {
return m_face_cell_local_face; return m_face_cell_local_face;
} }
...@@ -90,18 +110,26 @@ public: ...@@ -90,18 +110,26 @@ public:
m_number_of_nodes (number_of_cells+1), m_number_of_nodes (number_of_cells+1),
m_cell_nodes ("cell_nodes", m_number_of_cells), m_cell_nodes ("cell_nodes", m_number_of_cells),
m_cell_faces (m_cell_nodes), m_cell_faces (m_cell_nodes),
m_cell_nb_nodes ("cell_nb_nodes", m_number_of_cells),
m_cell_nb_faces (m_cell_nb_nodes),
m_node_nb_cells ("node_nb_cells",m_number_of_nodes), m_node_nb_cells ("node_nb_cells",m_number_of_nodes),
m_face_nb_cells (m_node_nb_cells), m_face_nb_cells (m_node_nb_cells),
m_node_cells ("node_cells", m_number_of_nodes), m_node_cells ("node_cells", m_number_of_nodes),
m_face_cells (m_node_cells), m_face_cells (m_node_cells),
m_node_cell_local_node ("node_cell_local_node", m_number_of_nodes), m_node_cell_local_node ("node_cell_local_node", m_number_of_nodes),
m_face_cell_local_face (m_node_cell_local_node) m_face_cell_local_face (m_node_cell_local_node),
m_max_nb_node_per_cell (2)
{ {
assert(number_of_cells>0); assert(number_of_cells>0);
Kokkos::parallel_for(m_number_of_cells, KOKKOS_LAMBDA(const size_t& j) { Kokkos::parallel_for(m_number_of_cells, KOKKOS_LAMBDA(const size_t& j) {
m_cell_nodes(j,0)=j; m_cell_nb_nodes[j] = 2;
m_cell_nodes(j,1)=j+1; });
Kokkos::parallel_for(m_number_of_cells, KOKKOS_LAMBDA(const size_t& j) {
for (int R=0; R<2; ++R) {
m_cell_nodes(j,R)=j+R;
}
}); });
Kokkos::parallel_for(m_number_of_nodes, KOKKOS_LAMBDA(const size_t& r) { Kokkos::parallel_for(m_number_of_nodes, KOKKOS_LAMBDA(const size_t& r) {
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment