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

starts externalization of variables

parent bfb44870
Branches
Tags
2 merge requests!2Develop,!1Develop
......@@ -10,11 +10,13 @@
#include <TinyMatrix.hpp>
#include <Mesh.hpp>
#include <MeshData.hpp>
#include <FiniteVolumesEulerUnknowns.hpp>
template<typename MeshData>
class AcousticSolverWithMesh
{
typedef typename MeshData::MeshType MeshType;
typedef FiniteVolumesEulerUnknowns<MeshData> UnknownsType;
MeshData& m_mesh_data;
const MeshType& m_mesh;
......@@ -24,8 +26,8 @@ class AcousticSolverWithMesh
typedef TinyVector<dimension> Rd;
typedef TinyMatrix<dimension> Rdd;
private:
private:
struct ReduceMin
{
private:
......@@ -68,7 +70,7 @@ private:
computeRhoCj(const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const double*>& cj)
{
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
m_rhocj[j] = rhoj[j]*cj[j];
});
return m_rhocj;
......@@ -81,7 +83,7 @@ private:
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_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
for (int r=0; r<cell_nb_nodes[j]; ++r) {
m_Ajr(j,r) = tensorProduct(rhocj(j)*Cjr(j,r), Cjr(j,r));
}
......@@ -97,7 +99,7 @@ private:
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();
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
Kokkos::parallel_for(m_mesh.numberOfNodes(), KOKKOS_LAMBDA(const int& r) {
Rdd sum = zero;
for (int j=0; j<node_nb_cells(r); ++j) {
const int J = node_cells(r,j);
......@@ -120,7 +122,7 @@ private:
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();
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
Kokkos::parallel_for(m_mesh.numberOfNodes(), KOKKOS_LAMBDA(const int& r) {
Rd& br = m_br(r);
br = zero;
for (int j=0; j<node_nb_cells(r); ++j) {
......@@ -138,11 +140,11 @@ private:
const Kokkos::View<const Rd*>& br) {
inverse(Ar, m_inv_Ar);
const Kokkos::View<const Rdd*> invAr = m_inv_Ar;
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
Kokkos::parallel_for(m_mesh.numberOfNodes(), KOKKOS_LAMBDA(const int& r) {
m_ur[r]=invAr(r)*br(r);
});
m_ur[0]=zero;
m_ur[m_nr-1]=zero;
m_ur[m_mesh.numberOfNodes()-1]=zero;
return m_ur;
}
......@@ -157,7 +159,7 @@ private:
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_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
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);
}
......@@ -183,12 +185,12 @@ private:
KOKKOS_INLINE_FUNCTION
double acoustic_dt(const Kokkos::View<const double*>& Vj,
const Kokkos::View<const double*>& cj) const {
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
m_Vj_over_cj[j] = Vj[j]/cj[j];
});
double dt = std::numeric_limits<double>::max();
Kokkos::parallel_reduce(m_nj, ReduceMin(m_Vj_over_cj), dt);
Kokkos::parallel_reduce(m_mesh.numberOfCells(), ReduceMin(m_Vj_over_cj), dt);
return dt;
}
......@@ -214,9 +216,6 @@ private:
Fjr = computeFjr(Ajr, ur, Cjr, uj, pj);
}
const size_t m_nj;
const size_t m_nr;
Kokkos::View<Rd*> m_br;
Kokkos::View<Rdd**> m_Ajr;
Kokkos::View<Rdd*> m_Ar;
......@@ -228,31 +227,29 @@ private:
public:
AcousticSolverWithMesh(MeshData& mesh_data)
AcousticSolverWithMesh(MeshData& mesh_data,
UnknownsType& unknowns)
: m_mesh_data(mesh_data),
m_mesh(mesh_data.mesh()),
m_connectivity(m_mesh.connectivity()),
m_nj(m_mesh.numberOfCells()),
m_nr(m_mesh.numberOfNodes()),
m_br("br", m_nr),
m_Ajr("Ajr", m_nj, m_connectivity.maxNbNodePerCell()),
m_Ar("Ar", m_nr),
m_inv_Ar("inv_Ar", m_nr),
m_Fjr("Fjr", m_nj, m_connectivity.maxNbNodePerCell()),
m_ur("ur", m_nr),
m_rhocj("rho_c", m_nj),
m_Vj_over_cj("Vj_over_cj", m_nj)
m_br("br", m_mesh.numberOfNodes()),
m_Ajr("Ajr", m_mesh.numberOfCells(), m_connectivity.maxNbNodePerCell()),
m_Ar("Ar", m_mesh.numberOfNodes()),
m_inv_Ar("inv_Ar", m_mesh.numberOfNodes()),
m_Fjr("Fjr", m_mesh.numberOfCells(), m_connectivity.maxNbNodePerCell()),
m_ur("ur", m_mesh.numberOfNodes()),
m_rhocj("rho_c", m_mesh.numberOfCells()),
m_Vj_over_cj("Vj_over_cj", m_mesh.numberOfCells())
{
Kokkos::View<double*> rhoj("rhoj",m_nj);
Kokkos::View<Rd*> uj("uj",m_nj);
Kokkos::View<double*> rhoj = unknowns.rhoj();
Kokkos::View<Rd*> uj = unknowns.uj();
Kokkos::View<double*> Ej = unknowns.Ej();
Kokkos::View<double*> Ej("Ej",m_nj);
Kokkos::View<double*> ej("ej",m_nj);
Kokkos::View<double*> pj("pj",m_nj);
Kokkos::View<double*> gammaj("gammaj",m_nj);
Kokkos::View<double*> cj("cj",m_nj);
Kokkos::View<double*> mj("mj",m_nj);
Kokkos::View<double*> ej = unknowns.ej();
Kokkos::View<double*> pj = unknowns.pj();
Kokkos::View<double*> gammaj = unknowns.gammaj();
Kokkos::View<double*> cj = unknowns.cj();
const Kokkos::View<const double*> mj = unknowns.mj();
Kokkos::View<Rd*> xr = m_mesh.xr();
......@@ -262,43 +259,7 @@ public:
const Kokkos::View<const unsigned int**>& cell_nodes = m_connectivity.cellNodes();
const Kokkos::View<const Rd**> Cjr = m_mesh_data.Cjr();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
if (xj[j][0]<0.5) {
rhoj[j]=1;
} else {
rhoj[j]=0.125;
}
});
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
if (xj[j][0]<0.5) {
pj[j]=1;
} else {
pj[j]=0.1;
}
});
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
uj[j] = zero;
});
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
gammaj[j] = 1.4;
});
BlockPerfectGas block_eos(rhoj, ej, pj, gammaj, cj);
block_eos.updateEandCFromRhoP();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
Ej[j] = ej[j]+0.5*(uj[j],uj[j]);
});
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
mj[j] = rhoj[j] * Vj[j];
});
Kokkos::View<double*> inv_mj("inv_mj",m_nj);
Kokkos::View<double*> inv_mj("inv_mj",m_mesh.numberOfCells());
inverse(mj, inv_mj);
const double tmax=0.2;
......@@ -307,6 +268,8 @@ public:
int itermax=std::numeric_limits<int>::max();
int iteration=0;
BlockPerfectGas block_eos(rhoj, ej, pj, gammaj, cj);
while((t<tmax) and (iteration<itermax)) {
double dt = 0.4*acoustic_dt(Vj, cj);
if (t+dt<tmax) {
......@@ -323,7 +286,7 @@ public:
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_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
Rd momentum_fluxes = zero;
double energy_fluxes = 0;
for (int R=0; R<cell_nb_nodes[j]; ++R) {
......@@ -335,17 +298,17 @@ public:
Ej[j] -= (dt*inv_mj[j]) * energy_fluxes;
});
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
ej[j] = Ej[j] - 0.5 * (uj[j],uj[j]);
});
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r){
Kokkos::parallel_for(m_mesh.numberOfNodes(), KOKKOS_LAMBDA(const int& r){
xr[r] += dt*ur[r];
});
m_mesh_data.updateAllData();
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j){
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
rhoj[j] = mj[j]/Vj[j];
});
......@@ -354,13 +317,6 @@ public:
++iteration;
}
// {
// std::ofstream fout("rho");
// for (int j=0; j<nj; ++j) {
// fout << xj[j][0] << ' ' << rhoj[j] << '\n';
// }
// }
std::cout << "* " << rang::style::underline << "Final time" << rang::style::reset
<< ": " << rang::fgB::green << t << rang::fg::reset << " (" << iteration << " iterations)\n";
......
#ifndef FINITE_VOLUMES_EULER_UNKNOWNS_HPP
#define FINITE_VOLUMES_EULER_UNKNOWNS_HPP
template <typename TMeshData>
class FiniteVolumesEulerUnknowns
{
public:
typedef TMeshData MeshDataType;
typedef typename MeshDataType::MeshType MeshType;
static constexpr size_t dimension = MeshType::dimension;
typedef TinyVector<dimension> Rd;
private:
const MeshDataType& m_mesh_data;
const MeshType& m_mesh;
Kokkos::View<double*> m_rhoj;
Kokkos::View<Rd*> m_uj;
Kokkos::View<double*> m_Ej;
Kokkos::View<double*> m_ej;
Kokkos::View<double*> m_pj;
Kokkos::View<double*> m_gammaj;
Kokkos::View<double*> m_cj;
Kokkos::View<double*> m_mj;
public:
Kokkos::View<double*> rhoj()
{
return m_rhoj;
}
const Kokkos::View<const double*> rhoj() const
{
return m_rhoj;
}
Kokkos::View<Rd*> uj()
{
return m_uj;
}
const Kokkos::View<const Rd*> uj() const
{
return m_uj;
}
Kokkos::View<double*> Ej()
{
return m_Ej;
}
const Kokkos::View<const double*> Ej() const
{
return m_Ej;
}
Kokkos::View<double*> ej()
{
return m_ej;
}
const Kokkos::View<const double*> ej() const
{
return m_ej;
}
Kokkos::View<double*> pj()
{
return m_pj;
}
const Kokkos::View<const double*> pj() const
{
return m_pj;
}
Kokkos::View<double*> gammaj()
{
return m_gammaj;
}
const Kokkos::View<const double*> gammaj() const
{
return m_gammaj;
}
Kokkos::View<double*> cj()
{
return m_cj;
}
const Kokkos::View<const double*> cj() const
{
return m_cj;
}
Kokkos::View<double*> mj()
{
return m_mj;
}
const Kokkos::View<const double*> mj() const
{
return m_mj;
}
void initializeSod()
{
const Kokkos::View<const Rd*> xj = m_mesh_data.xj();
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
if (xj[j][0]<0.5) {
m_rhoj[j]=1;
} else {
m_rhoj[j]=0.125;
}
});
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
if (xj[j][0]<0.5) {
m_pj[j]=1;
} else {
m_pj[j]=0.1;
}
});
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
m_uj[j] = zero;
});
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
m_gammaj[j] = 1.4;
});
BlockPerfectGas block_eos(m_rhoj, m_ej, m_pj, m_gammaj, m_cj);
block_eos.updateEandCFromRhoP();
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
m_Ej[j] = m_ej[j]+0.5*(m_uj[j],m_uj[j]);
});
const Kokkos::View<const double*> Vj = m_mesh_data.Vj();
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
m_mj[j] = m_rhoj[j] * Vj[j];
});
}
FiniteVolumesEulerUnknowns(const MeshDataType& mesh_data)
: m_mesh_data(mesh_data),
m_mesh(m_mesh_data.mesh()),
m_rhoj("rhoj",m_mesh.numberOfCells()),
m_uj("uj",m_mesh.numberOfCells()),
m_Ej("Ej",m_mesh.numberOfCells()),
m_ej("ej",m_mesh.numberOfCells()),
m_pj("pj",m_mesh.numberOfCells()),
m_gammaj("gammaj",m_mesh.numberOfCells()),
m_cj("cj",m_mesh.numberOfCells()),
m_mj("mj",m_mesh.numberOfCells())
{
;
}
};
#endif // FINITE_VOLUMES_EULER_UNKNOWNS_HPP
......@@ -120,12 +120,27 @@ int main(int argc, char *argv[])
Connectivity1D connectivity(number);
typedef Mesh<Connectivity1D> MeshType;
typedef MeshData<MeshType> MeshDataType;
typedef FiniteVolumesEulerUnknowns<MeshDataType> UnknownsType;
MeshType mesh(connectivity);
MeshDataType mesh_data(mesh);
AcousticSolverWithMesh<MeshDataType> acoustic_solver(mesh_data);
UnknownsType unknowns(mesh_data);
unknowns.initializeSod();
AcousticSolverWithMesh<MeshDataType> acoustic_solver(mesh_data, unknowns);
method_cost_map["AcousticSolverWithMesh"] = timer.seconds();
{
const Kokkos::View<const Rd*> xj = mesh_data.xj();
const Kokkos::View<const double*> rhoj = unknowns.rhoj();
std::ofstream fout("rho");
for (int j=0; j<mesh.numberOfCells(); ++j) {
fout << xj[j][0] << ' ' << rhoj[j] << '\n';
}
}
}
Kokkos::View<TinyVector<2,double>*[2]> test("test", 10);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment