Select Git revision
MathModule.cpp
FiniteVolumesDiffusion.hpp 9.83 KiB
#ifndef FINITE_VOLUMES_DIFFUSION_HPP
#define FINITE_VOLUMES_DIFFUSION_HPP
// --- INCLUSION fichiers headers ---
#include <Kokkos_Core.hpp>
#include <rang.hpp>
#include <BlockPerfectGas.hpp>
#include <TinyVector.hpp>
#include <TinyMatrix.hpp>
#include <Mesh.hpp>
#include <MeshData.hpp>
#include <FiniteVolumesEulerUnknowns.hpp>
// ---------------------------------
// Creation classe FiniteVolumesDiffusion
template<typename MeshData>
class FiniteVolumesDiffusion
{
typedef typename MeshData::MeshType MeshType;
typedef FiniteVolumesEulerUnknowns<MeshData> UnknownsType;
MeshData& m_mesh_data;
const MeshType& m_mesh;
const typename MeshType::Connectivity& m_connectivity;
constexpr static size_t dimension = MeshType::dimension;
typedef TinyVector<dimension> Rd;
typedef TinyMatrix<dimension> Rdd;
private:
// Sert a calculer les reductions (en gros calculer le min sur des
// vecteurs en parallele)
struct ReduceMin
{
private:
const Kokkos::View<const double*> x_;
public:
typedef Kokkos::View<const double*>::non_const_value_type value_type;
ReduceMin(const Kokkos::View<const double*>& x) : x_ (x) {}
typedef Kokkos::View<const 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();
}
};
// Calcule un Fl
Kokkos::View<Rdd*>
computeFl(const Kokkos::View<const Rd**>& Cjr,
const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const double*>& kj) {
const Kokkos::View<const unsigned int**>& face_cells = m_connectivity.faceCells();
const Kokkos::View<const unsigned short*> face_nb_cells
= m_connectivity.faceNbCells();
const Kokkos::View<const unsigned short**> face_cell_local_face
= m_mesh.connectivity().faceCellLocalFace();
const Kokkos::View<const double*>& Vl = m_mesh_data.Vl();
Kokkos::parallel_for(m_mesh.numberOfFaces(), KOKKOS_LAMBDA(const int& l) {
Rdd sum = zero;
double sum2 = 0.;
for (int j=0; j<face_nb_cells(l); ++j) {
int cell_here = face_cells(l,j);
int local_face_number_in_cell = face_cell_local_face(l,j);
sum -= tensorProduct(uj(cell_here, local_face_number_in_cell), Cjr(cell_here, local_face_number_in_cell));
sum2 += kj(cell_here);
}
// k = x
m_Fl(l) = ((sum2*0.5)/Vl(l))*sum;
// k = 2
//m_Fl(l)= (2./Vl(l))*sum;
});
return m_Fl ;
}
// Calcule un Gl
Kokkos::View<Rd*>
computeGl(const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const Rdd*>& Fl) {
const Kokkos::View<const unsigned int**>& face_cells = m_connectivity.faceCells();
const Kokkos::View<const unsigned short*> face_nb_cells
= m_connectivity.faceNbCells();
const Kokkos::View<const unsigned short**> face_cell_local_face
= m_mesh.connectivity().faceCellLocalFace();
Kokkos::parallel_for(m_mesh.numberOfFaces(), KOKKOS_LAMBDA(const int& l) {
Rd sum = zero;
for (int j=0; j<face_nb_cells(l); ++j) {
int cell_here = face_cells(l,j);
int local_face_number_in_cell = face_cell_local_face(l,j);
sum += uj(cell_here, local_face_number_in_cell);
}
m_Gl(l) = 0.5*Fl(l)*sum;
});
return m_Gl ;
}
// Calcule la liste des inverse d'une liste de matrices
// (pour l'instant, juste 1x1)
void inverse(const Kokkos::View<const Rdd*>& A,
Kokkos::View<Rdd*>& inv_A) const {
Kokkos::parallel_for(A.size(), KOKKOS_LAMBDA(const int& r) {
inv_A(r) = Rdd{1./(A(r)(0,0))};
});
}
// Calcule la liste des inverses d'une liste de reels
void 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);
});
}
// Enchaine les operations pour calculer les flux (Fjr et Gjr) pour
// pouvoir derouler le schema
KOKKOS_INLINE_FUNCTION
void computeExplicitFluxes(const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const Rd**>& Cjr,
const Kokkos::View<const double*>& kj) {
Kokkos::View<Rdd*> Fl = m_Fl ;
Fl = computeFl (Cjr, uj, kj);
Kokkos::View<Rd*> Gl = m_Gl ;
Gl = computeGl (uj, Fl );
}
Kokkos::View<Rdd*> m_Fl;
Kokkos::View<Rd*> m_Gl;
public:
FiniteVolumesDiffusion(MeshData& mesh_data,
UnknownsType& unknowns)
: m_mesh_data(mesh_data),
m_mesh(mesh_data.mesh()),
m_connectivity(m_mesh.connectivity()),
m_Fl("Fl", m_mesh.numberOfFaces()),
m_Gl("Gl", m_mesh.numberOfFaces())
{
;
}
// Calcule une evaluation du pas de temps verifiant le CFL parabolique
// Utilise la reduction definie dans la structure ReduceMin.
KOKKOS_INLINE_FUNCTION
double diffusion_dt(const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const double*>& kj) const {
Kokkos::View<double*> dt_j("dt_j", m_mesh.numberOfCells());
const Kokkos::View<const Rd**> Cjr = m_mesh_data.Cjr();
const Kokkos::View<const double*>& Vl = m_mesh_data.Vl();
const Kokkos::View<const double*>& Vj = m_mesh_data.Vj();
const Kokkos::View<const unsigned int**>& cell_faces = m_connectivity.cellFaces();
const Kokkos::View<const unsigned short*> cell_nb_faces
= m_connectivity.cellNbFaces();
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_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
double minVl = std::numeric_limits<double>::max();
for (int ll=0; ll<cell_nb_faces(j); ++ll) {
minVl = std::min(minVl, Vl(cell_faces(j, ll)));
}
//k=2 => (kj(j+1) + 2*kj(j) + kj(j-1)) = 8
// dt_j[j]= 0.5*rhoj(j)*Vj(j)*(2./8.)*minVl;
// k=x
double sum = 0.;
for (int m = 0; m < cell_nb_nodes(j); ++m) {
sum += kj(cell_nodes(j,m));
}
dt_j[j]= 0.5*rhoj(j)*Vj(j)*(1./(2.*kj(j) + sum))*minVl;
});
// for (int j=0; j<m_mesh.numberOfCells(); ++j) {
// std::cout << "dt_j[" << j << "]=" << dt_j[j] << '\n';
// }
double dt = std::numeric_limits<double>::max();
Kokkos::parallel_reduce(m_mesh.numberOfCells(), ReduceMin(dt_j), dt);
// std::cout << dt << std::endl;
// std::exit(0);
return dt;
}
// Avance la valeur des inconnues pendant un pas de temps dt
void computeNextStep(const double& t, const double& dt,
UnknownsType& unknowns)
{
Kokkos::View<double*> rhoj = unknowns.rhoj();
Kokkos::View<Rd*> uj = unknowns.uj();
Kokkos::View<Rd*> Sj = unknowns.Sj();
Kokkos::View<double*> Ej = unknowns.Ej();
Kokkos::View<double*> ej = unknowns.ej();
Kokkos::View<double*> gammaj = unknowns.gammaj();
const Kokkos::View<const double*> kj = unknowns.kj();
const Kokkos::View<const Rd*> xj = m_mesh_data.xj();
const Kokkos::View<const double*> Vj = m_mesh_data.Vj();
const Kokkos::View<const Rd**> Cjr = m_mesh_data.Cjr();
// Calcule les flux
computeExplicitFluxes(uj, Cjr, kj);
const Kokkos::View<const Rdd*> Fl = m_Fl ;
const Kokkos::View<const Rd *> Gl = m_Gl ;
const Kokkos::View<const unsigned short*>& cell_nb_faces
= m_connectivity.cellNbFaces();
const Kokkos::View<const unsigned int*[2]>& cell_faces
= m_connectivity.cellFaces();
// Mise a jour de la vitesse et de l'energie totale specifique
const Kokkos::View<const double*> inv_mj = unknowns.invMj();
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
Rd momentum_fluxes = zero;
double energy_fluxes = 0;
int l = 0.;
for (int R=0; R<cell_nb_faces(j); ++R) {
l = cell_faces(j,R);
momentum_fluxes += Fl(l)*Cjr(j,R);
energy_fluxes += (Gl(l), Cjr(j,R));
}
uj[j] += std::exp(-t)*(dt*inv_mj[j])*Vj(j)*Sj(j) + (dt*inv_mj[j]) * momentum_fluxes;
//uj[j] += (dt*inv_mj[j]) * momentum_fluxes;
Ej[j] += (dt*inv_mj[j]) * energy_fluxes;
});
// Calcul de e par la formule e = E-0.5 u^2
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j) {
ej[j] = Ej[j] - 0.5 * (uj[j],uj[j]);
});
// met a jour les quantites (geometriques) associees au maillage
m_mesh_data.updateAllData();
// Calcul de rho avec la formule Mj = Vj rhoj
const Kokkos::View<const double*> mj = unknowns.mj();
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
rhoj[j] = mj[j]/Vj[j];
});
}
// Calcul erreur entre solution analytique et solution numerique en norme L2
// (quand la solution exacte est connue)
double error_L2(UnknownsType& unknowns) {
Kokkos::View<Rd*> uj = unknowns.uj();
const Kokkos::View<const Rd*> xj = m_mesh_data.xj();
double pi = 4.*std::atan(1.);
double erreur = 0.;
double exacte = 0.;
for (size_t j=0; j<m_mesh.numberOfCells(); ++j) {
exacte = std::sin(pi*xj[j][0])*std::exp(-0.2); // solution exacte
erreur += (exacte - uj[j][0])*(exacte - uj[j][0]);
}
erreur = std::sqrt(erreur);
return erreur;
}
double error_Linf(UnknownsType& unknowns) {
Kokkos::View<Rd*> uj = unknowns.uj();
const Kokkos::View<const Rd*> xj = m_mesh_data.xj();
double pi = 4.*std::atan(1.);
double exacte = std::sin(pi*xj[0][0])*std::exp(-0.2);
double erreur = std::abs(exacte - uj[0][0]);
for (size_t j=1; j<m_mesh.numberOfCells(); ++j) {
exacte = std::sin(pi*xj[j][0])*std::exp(-0.2);
if (std::abs(exacte - uj[j][0]) > erreur) {
erreur = std::abs(exacte - uj[j][0]);
}
}
return erreur;
}
};
#endif // FINITE_VOLUMES_DIFFUSION_HPP