Select Git revision
MeshNodeBoundary.hpp
FiniteVolumesDiffusion.hpp 7.50 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> // MeshData est le type generique des donnees (geometriques) attachees a un maillage
class FiniteVolumesDiffusion
{
typedef typename MeshData::MeshType MeshType; // type du maillage
typedef FiniteVolumesEulerUnknowns<MeshData> UnknownsType; // type des inconnues
MeshData& m_mesh_data; // reference vers les donnees attachees du maillage
const MeshType& m_mesh; // reference vers le maillage
const typename MeshType::Connectivity& m_connectivity; // references vers la connectivite
constexpr static size_t dimension = MeshType::dimension; // dimension du maillage (connue a la compilation)
typedef TinyVector<dimension> Rd; // type de petits vecteurs (de dimension MeshType::dimension)
typedef TinyMatrix<dimension> Rdd; // type de petites matrices
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();
}
};
Kokkos::View<Rd**> // Fonction qui calcule F_jr
computeFjr(const Kokkos::View<const Rd**>& Cjr,
const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const Rd*>& xr,
const Kokkos::View<const double*>& kj) {
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) {
for (int r=0; r<cell_nb_nodes[j]; ++r) {
m_Fjr(j,r) = ((kj(cell_nodes(j,r)) + kj(cell_nodes(j,r)-1))/(2*(xj(cell_nodes(j,r))-xj(cell_nodes(j,r)-1))))*uj(j,r)*Cjr(j,r);
}
});
return m_Fjr;
}
Kokkos::View<Rd**> // Fonction qui calcule G_jr
computeGjr(const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const double*>& Fjr) {
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) {
for (int r=0; r<cell_nb_nodes[j]; ++r) {
m_Gjr(j,r) = ((uj(cell_nodes(j,r)) + uj(cell_nodes(j,r)-1))/2)*Fjr(j,r);
}
});
return m_Gjr;
}
// Calcul la liste des inverses d'une liste de matrices (pour
// l'instant seulement $R^{1\times 1}$)
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))};
});
}
// Calcul 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*>& xr,
const Kokkos::View<const Rd*>& xj,
const Kokkos::View<const double*>& kj,
const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const Rd*>& uj,
const Kokkos::View<const Rd**>& Cjr) {
Kokkos::View<Rd**> Fjr = m_Fjr;
Fjr = computeFjr(Cjr, uj, xr, kj);
Kokkos::View<Rd**> Gjr = m_Gjr;
Gjr = computeGjr(uj, Fjr);
}
Kokkos::View<Rd**> m_Fjr;
Kokkos::View<double*> m_CFL;
public:
FiniteVolumesDiffusion(MeshData& mesh_data,
UnknownsType& unknowns)
: m_mesh_data(mesh_data),
m_mesh(mesh_data.mesh()),
m_connectivity(m_mesh.connectivity()),
m_Fjr("Fjr", m_mesh.numberOfCells(), m_connectivity.maxNbNodePerCell()),
m_CFL("CFL", m_mesh.numberOfCells())
{
;
}
// Calcule une evaluation du pas de temps verifiant le CFL parabolique
// Utilise la reduction definie dans la structure ReduceMin. Ici, dx_j=V_j
KOKKOS_INLINE_FUNCTION
double diffusion_dt(const Kokkos::View<const double*>& Vj,
const Kokkos::View<const double*>& rhoj,
const Kokkos::View<const Rd*>& xr,
const Kokkos::View<const Rd*>& kj) const {
Kokkos::View<double*> dt_j("dt_j", m_mesh.numberOfCells());
Kokkos::parallel_for(m_mesh.numberOfCells(), KOKKOS_LAMBDA(const int& j){
m_CFL(j) = rhoj(j)*Vj(j)*min(xj(j+1)-xj(j), xj(j)-xj(j-1))*(2./(kj(j+1) + 2*kj(j) + kj(j-1)));
});
double dt = std::numeric_limits<double>::max();
Kokkos::parallel_reduce(m_mesh.numberOfCells(), ReduceMin(dt_j), dt);
return dt;
}
// Avance la valeur des inconnues pendant un pas de temps dt // A MODIFIER
void computeNextStep(const double& t, const double& dt,
UnknownsType& unknowns)
{
Kokkos::View<double*> rhoj = unknowns.rhoj();
Kokkos::View<Rd*> uj = unknowns.uj();
Kokkos::View<double*> Ej = unknowns.Ej();
Kokkos::View<double*> ej = unknowns.ej();
Kokkos::View<double*> gammaj = unknowns.gammaj();
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();
Kokkos::View<Rd*> xr = m_mesh.xr();
// Calcule les flux
computeExplicitFluxes(xr, xj, rhoj, uj, Cjr);
const Kokkos::View<const Rd**> Fjr = m_Fjr;
const Kokkos::View<const unsigned short*> cell_nb_nodes
= m_connectivity.cellNbNodes();
const Kokkos::View<const unsigned int**>& cell_nodes
= m_connectivity.cellNodes();
// 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;
for (int R=0; R<cell_nb_nodes[j]; ++R) {
const int r=cell_nodes(j,R);
momentum_fluxes += Fjr(j,R);
energy_fluxes += Gjr(j,R);
}
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];
});
}
};
#endif // FINITE_VOLUMES_DIFFUSION_HPP