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

Added a test copy of the Acoustic solver

parent 713386db
No related branches found
No related tags found
No related merge requests found
...@@ -372,12 +372,12 @@ AcousticSolver::AcousticSolver(const long int& nj) ...@@ -372,12 +372,12 @@ AcousticSolver::AcousticSolver(const long int& nj)
++iteration; ++iteration;
} }
{ // {
std::ofstream fout("rho"); // std::ofstream fout("rho");
for (int j=0; j<nj; ++j) { // for (int j=0; j<nj; ++j) {
fout << xj[j] << ' ' << rhoj[j] << '\n'; // fout << xj[j] << ' ' << rhoj[j] << '\n';
} // }
} // }
std::cout << "* " << rang::style::underline << "Final time" << rang::style::reset std::cout << "* " << rang::style::underline << "Final time" << rang::style::reset
<< ": " << rang::fgB::green << t << rang::fg::reset << " (" << iteration << " iterations)\n"; << ": " << rang::fgB::green << t << rang::fg::reset << " (" << iteration << " iterations)\n";
......
#include <AcousticSolverTest.hpp>
#include <rang.hpp>
#include <memory>
#include <BlockPerfectGas.hpp>
typedef const double my_double;
struct AcousticSolverTest::ReduceMin
{
private:
const Kokkos::View<my_double*> x_;
public:
typedef Kokkos::View<my_double*>::non_const_value_type value_type;
ReduceMin(const Kokkos::View<my_double*>& x) : x_ (x) {}
typedef Kokkos::View<my_double*>::size_type 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
const Kokkos::View<const double*>
AcousticSolverTest::computeRhoCj(const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const double*>& cj)
{
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
m_rhocj[j] = rhoj[j]*cj[j];
});
return m_rhocj;
}
KOKKOS_INLINE_FUNCTION
const Kokkos::View<const double*[2]>
AcousticSolverTest::computeAjr(const Kokkos::View<const double*>& rhocj,
const Kokkos::View<const double*[2]>& Cjr)
{
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
for (int r=0; r<2; ++r) {
m_Ajr(j,r) = rhocj(j)*Cjr(j,r)*Cjr(j,r);
}
});
return m_Ajr;
}
KOKKOS_INLINE_FUNCTION
const Kokkos::View<const double*>
AcousticSolverTest::computeAr(const Kokkos::View<const double*[2]>& Ajr,
const Kokkos::View<const int*[2]>& node_cells,
const Kokkos::View<const int*[2]>& node_cell_local_node,
const Kokkos::View<const int*>& node_nb_cells)
{
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
double sum = 0;
for (int j=0; j<node_nb_cells(r); ++j) {
const int J = node_cells(r,j);
const int R = node_cell_local_node(r,j);
sum += Ajr(J,R);
}
m_Ar(r) = sum;
});
return m_Ar;
}
KOKKOS_INLINE_FUNCTION
const Kokkos::View<const double*>
AcousticSolverTest::computeBr(const Kokkos::View<const double*[2]>& Ajr,
const Kokkos::View<const double*[2]>& Cjr,
const Kokkos::View<const double*>& uj,
const Kokkos::View<const double*>& pj,
const Kokkos::View<const int*[2]>& node_cells,
const Kokkos::View<const int*[2]>& node_cell_local_node,
const Kokkos::View<const int*>& node_nb_cells)
{
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
double sum = 0;
for (int j=0; j<node_nb_cells(r); ++j) {
const int J = node_cells(r,j);
const int R = node_cell_local_node(r,j);
sum += Ajr(J,R)*uj(J) + Cjr(J,R)*pj(J);
}
m_br(r) = sum;
});
return m_br;
}
KOKKOS_INLINE_FUNCTION
Kokkos::View<double*>
AcousticSolverTest::computeUr(const Kokkos::View<const double*>& Ar,
const Kokkos::View<const double*>& br)
{
inverse(Ar, m_inv_Ar);
const Kokkos::View<const double*> invAr = m_inv_Ar;
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r) {
m_ur[r]=br(r)*invAr(r);
});
m_ur[0]=0;
m_ur[m_nr-1]=0;
return m_ur;
}
KOKKOS_INLINE_FUNCTION
Kokkos::View<double*[2]>
AcousticSolverTest::computeFjr(const Kokkos::View<const double*[2]>& Ajr,
const Kokkos::View<const double*>& ur,
const Kokkos::View<const double*[2]>& Cjr,
const Kokkos::View<const double*>& uj,
const Kokkos::View<const double*>& pj,
const Kokkos::View<const int*[2]>& cell_nodes)
{
Kokkos::parallel_for(m_nj, KOKKOS_LAMBDA(const int& j) {
for (int r=0; r<2; ++r) {
m_Fjr(j,r) = Ajr(j,r)*(uj(j)-ur(cell_nodes(j,r)))+Cjr(j,r)*pj(j);
}
});
return m_Fjr;
}
KOKKOS_INLINE_FUNCTION
double AcousticSolverTest::
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){
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);
return dt;
}
KOKKOS_INLINE_FUNCTION
void
AcousticSolverTest::inverse(const Kokkos::View<const double*>& x,
Kokkos::View<double*>& inv_x) const
{
Kokkos::parallel_for(x.size(), KOKKOS_LAMBDA(const int& r) {
inv_x(r) = 1./x(r);
});
}
KOKKOS_INLINE_FUNCTION
void AcousticSolverTest::computeExplicitFluxes(const Kokkos::View<const double*>& xr,
const Kokkos::View<const double*>& xj,
const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const double*>& uj,
const Kokkos::View<const double*>& pj,
const Kokkos::View<const double*>& cj,
const Kokkos::View<const double*>& Vj,
const Kokkos::View<const double*[2]>& Cjr,
const Kokkos::View<const int*[2]>& cell_nodes,
const Kokkos::View<const int*[2]>& node_cells,
const Kokkos::View<const int*>& node_nb_cells,
const Kokkos::View<const int*[2]>& node_cell_local_node,
Kokkos::View<double*>& ur,
Kokkos::View<double*[2]>& Fjr)
{
const Kokkos::View<const double*> rhocj = computeRhoCj(rhoj, cj);
const Kokkos::View<const double*[2]> Ajr = computeAjr(rhocj, Cjr);
const Kokkos::View<const double*> Ar = computeAr(Ajr, node_cells, node_cell_local_node, node_nb_cells);
const Kokkos::View<const double*> br = computeBr(Ajr, Cjr, uj, pj,
node_cells, node_cell_local_node, node_nb_cells);
ur = computeUr(Ar, br);
Fjr = computeFjr(Ajr, ur, Cjr, uj, pj, cell_nodes);
}
AcousticSolverTest::AcousticSolverTest(const long int& nj)
: m_nj(nj),
m_nr(nj+1),
m_br("br", m_nr),
m_Ajr("Ajr", m_nj),
m_Ar("Ar", m_nr),
m_inv_Ar("inv_Ar", m_nr),
m_Fjr("Fjr", m_nj),
m_ur("ur", m_nr),
m_rhocj("rho_c", m_nj),
m_Vj_over_cj("Vj_over_cj", m_nj)
{
Kokkos::View<double*> xj("xj",m_nj);
Kokkos::View<double*> rhoj("rhoj",m_nj);
Kokkos::View<double*> uj("uj",m_nj);
Kokkos::View<double*> Ej("Ej",m_nj);
Kokkos::View<double*> ej("ej",m_nj);
Kokkos::View<double*> pj("pj",m_nj);
Kokkos::View<double*> Vj("Vj",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*> xr("xr", m_nr);
const double delta_x = 1./m_nj;
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r){
xr[r] = r*delta_x;
});
Kokkos::View<int*[2]> cell_nodes("cell_nodes",m_nj,2);
Kokkos::View<int*[2]> node_cells("node_cells",m_nr,2);
Kokkos::View<int*[2]> node_cell_local_node("node_cells",m_nr,2);
Kokkos::View<int*> node_nb_cells("node_cells",m_nr);
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r){
node_nb_cells(r) = 2;
});
node_nb_cells(0) = 1;
node_nb_cells(m_nr-1) = 1;
node_cells(0,0) = 0;
Kokkos::parallel_for(m_nr-2, KOKKOS_LAMBDA(const int& r){
node_cells(r+1,0) = r;
node_cells(r+1,1) = r+1;
});
node_cells(m_nr-1,0) =m_nj-1;
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
cell_nodes(j,0) = j;
cell_nodes(j,1) = j+1;
});
Kokkos::parallel_for(m_nr, KOKKOS_LAMBDA(const int& r){
for (int J=0; J<node_nb_cells(r); ++J) {
int j = node_cells(r,J);
for (int R=0; R<2; ++R) {
if (cell_nodes(j,R) == r) {
node_cell_local_node(r,J) = R;
}
}
}
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
xj[j] = 0.5*(xr[cell_nodes(j,0)]+xr[cell_nodes(j,1)]);
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
Vj[j] = xr[cell_nodes(j,1)]-xr[cell_nodes(j,0)];
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
if (xj[j]<0.5) {
rhoj[j]=1;
} else {
rhoj[j]=0.125;
}
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
if (xj[j]<0.5) {
pj[j]=1;
} else {
pj[j]=0.1;
}
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
gammaj[j] = 1.4;
});
BlockPerfectGas block_eos(rhoj, ej, pj, gammaj, cj);
block_eos.updateEandCFromRhoP();
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
Ej[j] = ej[j]+0.5*uj[j]*uj[j];
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
mj[j] = rhoj[j] * Vj[j];
});
Kokkos::View<double*> inv_mj("inv_mj",m_nj);
inverse(mj, inv_mj);
const double tmax=0.2;
double t=0;
int itermax=std::numeric_limits<int>::max();
int iteration=0;
Kokkos::View<double*[2]> Cjr("Cjr",m_nj);
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j) {
Cjr(j,0)=-1;
Cjr(j,1)= 1;
});
while((t<tmax) and (iteration<itermax)) {
double dt = 0.4*acoustic_dt(Vj, cj);
if (t+dt<tmax) {
t+=dt;
} else {
dt=tmax-t;
t=tmax;
}
computeExplicitFluxes(xr, xj,
rhoj, uj, pj, cj, Vj, Cjr,
cell_nodes, node_cells, node_nb_cells, node_cell_local_node,
m_ur, m_Fjr);
const Kokkos::View<const double*[2]> Fjr = m_Fjr;
const Kokkos::View<const double*> ur = m_ur;
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
double momentum_fluxes = 0;
double energy_fluxes = 0;
for (int R=0; R<2; ++R) {
const int r=cell_nodes(j,R);
momentum_fluxes += Fjr(j,R);
energy_fluxes += Fjr(j,R)*ur[r];
}
uj[j] -= dt*inv_mj[j]*(momentum_fluxes);
Ej[j] -= dt*inv_mj[j]*(energy_fluxes);
});
Kokkos::parallel_for(nj, 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){
xr[r] += dt*ur[r];
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
xj[j] = 0.5*(xr[j]+xr[j+1]);
Vj[j] = xr[j+1]-xr[j];
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
rhoj[j] = mj[j]/Vj[j];
});
block_eos.updatePandCFromRhoE();
++iteration;
}
// {
// std::ofstream fout("rho");
// for (int j=0; j<nj; ++j) {
// fout << xj[j] << ' ' << rhoj[j] << '\n';
// }
// }
std::cout << "* " << rang::style::underline << "Final time" << rang::style::reset
<< ": " << rang::fgB::green << t << rang::fg::reset << " (" << iteration << " iterations)\n";
}
#ifndef ACOUSTIC_SOLVER_TEST_HPP
#define ACOUSTIC_SOLVER_TEST_HPP
#include <Kokkos_Core.hpp>
class AcousticSolverTest
{
private:
inline const Kokkos::View<const double*>
computeRhoCj(const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const double*>& cj);
inline const Kokkos::View<const double*[2]>
computeAjr(const Kokkos::View<const double*>& rhocj,
const Kokkos::View<const double*[2]>& Cjr);
inline const Kokkos::View<const double*>
computeAr(const Kokkos::View<const double*[2]>& Ajr,
const Kokkos::View<const int*[2]>& node_cells,
const Kokkos::View<const int*[2]>& node_cell_local_node,
const Kokkos::View<const int*>& node_nb_cells);
inline const Kokkos::View<const double*>
computeBr(const Kokkos::View<const double*[2]>& Ajr,
const Kokkos::View<const double*[2]>& Cjr,
const Kokkos::View<const double*>& uj,
const Kokkos::View<const double*>& pj,
const Kokkos::View<const int*[2]>& node_cells,
const Kokkos::View<const int*[2]>& node_cell_local_node,
const Kokkos::View<const int*>& node_nb_cells);
Kokkos::View<double*>
computeUr(const Kokkos::View<const double*>& Ar,
const Kokkos::View<const double*>& br);
Kokkos::View<double*[2]>
computeFjr(const Kokkos::View<const double*[2]>& Ajr,
const Kokkos::View<const double*>& ur,
const Kokkos::View<const double*[2]>& Cjr,
const Kokkos::View<const double*>& uj,
const Kokkos::View<const double*>& pj,
const Kokkos::View<const int*[2]>& cell_nodes);
void inverse(const Kokkos::View<const double*>& x,
Kokkos::View<double*>& inv_x) const;
inline double acoustic_dt(const Kokkos::View<const double*>& Vj,
const Kokkos::View<const double*>& cj) const;
inline void computeExplicitFluxes(const Kokkos::View<const double*>& xr,
const Kokkos::View<const double*>& xj,
const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const double*>& uj,
const Kokkos::View<const double*>& pj,
const Kokkos::View<const double*>& cj,
const Kokkos::View<const double*>& Vj,
const Kokkos::View<const double*[2]>& Cjr,
const Kokkos::View<const int*[2]>& cell_nodes,
const Kokkos::View<const int*[2]>& node_cells,
const Kokkos::View<const int*>& node_nb_cells,
const Kokkos::View<const int*[2]>& node_cell_local_node,
Kokkos::View<double*>& ur,
Kokkos::View<double*[2]>& Fjr);
struct ReduceMin;
const long int m_nj;
const long int m_nr;
Kokkos::View<double*> m_br;
Kokkos::View<double*[2]> m_Ajr;
Kokkos::View<double*> m_Ar;
Kokkos::View<double*> m_inv_Ar;
Kokkos::View<double*[2]> m_Fjr;
Kokkos::View<double*> m_ur;
Kokkos::View<double*> m_rhocj;
Kokkos::View<double*> m_Vj_over_cj;
public:
AcousticSolverTest(const long int& nj);
};
#endif // ACOUSTIC_SOLVER_TEST_HPP
...@@ -9,6 +9,7 @@ include_directories(${PASTIS_SOURCE_DIR}/packages/rang/include) ...@@ -9,6 +9,7 @@ include_directories(${PASTIS_SOURCE_DIR}/packages/rang/include)
add_library( add_library(
PastisExperimental PastisExperimental
AcousticSolver.cpp AcousticSolver.cpp
AcousticSolverTest.cpp
RawKokkosAcousticSolver.cpp RawKokkosAcousticSolver.cpp
MeshLessAcousticSolver.cpp MeshLessAcousticSolver.cpp
) )
......
...@@ -9,6 +9,7 @@ ...@@ -9,6 +9,7 @@
#include <RawKokkosAcousticSolver.hpp> #include <RawKokkosAcousticSolver.hpp>
#include <MeshLessAcousticSolver.hpp> #include <MeshLessAcousticSolver.hpp>
#include <AcousticSolver.hpp> #include <AcousticSolver.hpp>
#include <AcousticSolverTest.hpp>
#include <CLI/CLI.hpp> #include <CLI/CLI.hpp>
#include <cassert> #include <cassert>
...@@ -92,6 +93,13 @@ int main(int argc, char *argv[]) ...@@ -92,6 +93,13 @@ int main(int argc, char *argv[])
method_cost_map["MeshLessAcousticSolver"] = timer.seconds(); method_cost_map["MeshLessAcousticSolver"] = timer.seconds();
} }
{ // class for acoustic solver test
Kokkos::Timer timer;
timer.reset();
AcousticSolverTest acoustic_solver(number);
method_cost_map["AcousticSolverTest"] = timer.seconds();
}
{ // class for acoustic solver { // class for acoustic solver
Kokkos::Timer timer; Kokkos::Timer timer;
timer.reset(); timer.reset();
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment