![]() |
Shamrock 2025.10.0
Astrophysical Code
|
From original version by Thomas Guillet (T.A.G.nosp@m.uill.nosp@m.et@ex.nosp@m.eter.nosp@m..ac.u.nosp@m.k) More...
#include "shambackends/math.hpp"#include "shambackends/typeAliasVec.hpp"#include "shambackends/vec.hpp"#include <algorithm>#include <array>#include <cmath>#include <iostream>
Include dependency graph for riemann.hpp:
This graph shows which files directly or indirectly include this file:Go to the source code of this file.
Classes | |
| struct | shammath::ConsState< Tvec_ > |
| struct | shammath::PrimState< Tvec_ > |
| struct | shammath::Fluxes< Tvec_ > |
Namespaces | |
| namespace | shammath |
| namespace for math utility | |
Functions | |
| template<class Tvec > | |
| const ConsState< Tvec > | shammath::operator+ (const ConsState< Tvec > &lhs, const ConsState< Tvec > &rhs) |
| template<class Tvec > | |
| const ConsState< Tvec > | shammath::operator- (const ConsState< Tvec > &lhs, const ConsState< Tvec > &rhs) |
| template<class Tvec > | |
| const ConsState< Tvec > | shammath::operator* (const typename ConsState< Tvec >::Tscal factor, const ConsState< Tvec > &rhs) |
| template<class Tvec > | |
| const ConsState< Tvec > | shammath::operator* (const ConsState< Tvec > &lhs, const typename ConsState< Tvec >::Tscal factor) |
| template<class Tvec > | |
| constexpr shambase::VecComponent< Tvec > | shammath::rhoekin (shambase::VecComponent< Tvec > rho, Tvec v) |
| template<class Tvec > | |
| constexpr ConsState< Tvec > | shammath::prim_to_cons (const PrimState< Tvec > prim, typename PrimState< Tvec >::Tscal gamma) |
| template<class Tvec > | |
| constexpr PrimState< Tvec > | shammath::cons_to_prim (const ConsState< Tvec > cons, typename ConsState< Tvec >::Tscal gamma) |
| template<class Tvec > | |
| constexpr ConsState< Tvec > | shammath::hydro_flux_x (const ConsState< Tvec > cons, typename ConsState< Tvec >::Tscal gamma) |
| template<class Tvec > | |
| constexpr shambase::VecComponent< Tvec > | shammath::sound_speed (PrimState< Tvec > prim, shambase::VecComponent< Tvec > gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::rusanov_flux_x (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::y_to_x (const Tcons c) |
| template<class Tcons > | |
| constexpr Tcons | shammath::x_to_y (const Tcons c) |
| template<class Tcons > | |
| constexpr Tcons | shammath::z_to_x (const Tcons c) |
| template<class Tcons > | |
| constexpr Tcons | shammath::x_to_z (const Tcons c) |
| template<class Tcons > | |
| constexpr Tcons | shammath::invert_axis (const Tcons c) |
| template<class Tcons > | |
| constexpr Tcons | shammath::rusanov_flux_y (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::rusanov_flux_z (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::rusanov_flux_mx (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::rusanov_flux_my (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::rusanov_flux_mz (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr auto | shammath::hll_flux_x (const Tcons consL, const Tcons consR, const typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::hll_flux_y (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::hll_flux_z (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::hll_flux_mx (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::hll_flux_my (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::hll_flux_mz (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| template<class Tcons > | |
| constexpr Tcons | shammath::hllc_flux_x (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| HLLC solver based on section 10.4 from Toro 3rd Edition , Springer 2009. The wave speeds estimates are based on Bernd Einfeldt (SIAM, 1988), On Godunov-Type Methods for Gas Dynamics. | |
| template<class Tcons > | |
| constexpr Tcons | shammath::hllc_flux_y (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| HLLC flux in the +y direction. | |
| template<class Tcons > | |
| constexpr Tcons | shammath::hllc_flux_z (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| HLLC flux in the +z direction. | |
| template<class Tcons > | |
| constexpr Tcons | shammath::hllc_flux_mx (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| HLLC flux in the -x direction. | |
| template<class Tcons > | |
| constexpr Tcons | shammath::hllc_flux_my (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| HLLC flux in the -y direction. | |
| template<class Tcons > | |
| constexpr Tcons | shammath::hllc_flux_mz (Tcons cL, Tcons cR, typename Tcons::Tscal gamma) |
| HLLC flux in the -z direction. | |
From original version by Thomas Guillet (T.A.G.nosp@m.uill.nosp@m.et@ex.nosp@m.eter.nosp@m..ac.u.nosp@m.k)
Definition in file riemann.hpp.