Select Git revision
UnaryOperatorRegisterForN.cpp
-
Stéphane Del Pino authored
One now use OperatorRepository
Stéphane Del Pino authoredOne now use OperatorRepository
MeshLessAcousticSolver.cpp 5.49 KiB
#include <MeshLessAcousticSolver.hpp>
#include <rang.hpp>
#include <memory>
#include <BlockPerfectGas.hpp>
typedef const double my_double;
struct MeshLessAcousticSolver::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
double MeshLessAcousticSolver::
acoustic_dt(const Kokkos::View<const double*>& Vj,
const Kokkos::View<const double*>& cj) const
{
const size_t nj = Vj.size();
double dt = std::numeric_limits<double>::max();
Kokkos::View<double*> Vj_cj("Vj_cj", nj);
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
Vj_cj[j] = Vj[j]/cj[j];
});
Kokkos::parallel_reduce(nj, ReduceMin(Vj_cj), dt);
return dt;
}
KOKKOS_INLINE_FUNCTION
void MeshLessAcousticSolver::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,
Kokkos::View<double*>& ur,
Kokkos::View<double*>& pr) const
{
// calcul de ur
ur[0]=0;
const size_t nr = ur.size();
const size_t nj = uj.size();
Kokkos::parallel_for(nj-1, KOKKOS_LAMBDA(const int& j) {
const int r = j+1;
const int k = r;
const double ujr = uj[j];
const double ukr = uj[k];
const double pjr = pj[j];
const double pkr = pj[k];
ur[r]=(rhoj[j]*cj[j]*ujr + rhoj[k]*cj[k]*ukr + pjr-pkr)/(rhoj[j]*cj[j]+rhoj[k]*cj[k]);
});
ur[nr-1]=0;
// calcul de pr
pr[0] = pj[0] + rhoj[0]*cj[0]*(ur[0] - uj[0]);
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j) {
const int r = j+1;
const double ujr = uj[j];
const double pjr = pj[j];
pr[r]=pjr+rhoj[j]*cj[j]*(ujr-ur[r]);
});
}
MeshLessAcousticSolver::MeshLessAcousticSolver(const long int& nj)
{
Kokkos::View<double*> xj("xj", nj);
Kokkos::View<double*> rhoj("rhoj", nj);
Kokkos::View<double*> uj("uj", nj);
Kokkos::View<double*> Ej("Ej", nj);
Kokkos::View<double*> ej("ej", nj);
Kokkos::View<double*> pj("pj", nj);
Kokkos::View<double*> Vj("Vj", nj);
Kokkos::View<double*> gammaj("gammaj", nj);
Kokkos::View<double*> cj("cj", nj);
Kokkos::View<double*> mj("mj", nj);
Kokkos::View<double*> inv_mj("inv_mj", nj);
const int nr=nj+1;
Kokkos::View<double*> xr("xr", nr);
const double delta_x = 1./nj;
Kokkos::parallel_for(nr, KOKKOS_LAMBDA(const int& r){
xr[r] = r*delta_x;
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
xj[j] = 0.5*(xr[j]+xr[j+1]);
});
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
Vj[j] = xr[j+1]-xr[j];
});
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];
});
const double tmax=0.2;
double t=0;
int itermax=std::numeric_limits<int>::max();
int iteration=0;
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;
}
Kokkos::View<double*> ur("ur", nr);
Kokkos::View<double*> pr("pr", nr);
computeExplicitFluxes(xr, xj,
rhoj, uj, pj, cj, Vj,
ur, pr);
Kokkos::View<double*> new_uj("new_uj", nj);
Kokkos::View<double*> new_Ej("new_Ej", nj);
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
int rm=j;
int rp=j+1;
new_uj[j] = uj[j] + dt/mj[j]*(pr[rm]-pr[rp]);
new_Ej[j] = Ej[j] + dt/mj[j]*(pr[rm]*ur[rm]-pr[rp]*ur[rp]);
});
uj=new_uj;
Ej=new_Ej;
Kokkos::parallel_for(nj, KOKKOS_LAMBDA(const int& j){
ej[j] = Ej[j] - 0.5 * uj[j]*uj[j];
});
Kokkos::View<double*> xr_new("new_xr", nr);
Kokkos::parallel_for(nr, KOKKOS_LAMBDA(const int& r){
xr_new[r] = xr[r] + dt*ur[r];
});
xr = xr_new;
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::cout << "* " << rang::style::underline << "Final time" << rang::style::reset
<< ": " << rang::fgB::green << t << rang::fg::reset << " (" << iteration << " iterations)\n";
}