Skip to content

Commit

Permalink
implement full equation for G
Browse files Browse the repository at this point in the history
  • Loading branch information
psharda committed Aug 1, 2023
1 parent f72cce0 commit fe202ce
Show file tree
Hide file tree
Showing 2 changed files with 15 additions and 8 deletions.
15 changes: 11 additions & 4 deletions src/EOS.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -241,8 +241,11 @@ AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE auto EOS<problem_t>::ComputeOtherDeriva
amrex::Real deint_dP = NAN;
// compute derivative of density w/r/t pressure, given density and pressure
amrex::Real dRho_dP = NAN;
// compute derivative of pressure w/r/t density, given density and pressure (needed for the fundamental derivative G)
// compute derivative of pressure w/r/t density at constant entropy, given density and pressure (needed for the fundamental derivative G)
amrex:: Real dP_dRho_s = NAN;
amrex::Real G_gamma = NAN;
// compute second derivative of pressure w/r/t density at constant entropy, given density and pressure (needed for the fundamental derivative G)
amrex::Real d2P_dRho2_s = NAN;

#ifdef PRIMORDIAL_CHEM
eos_t chemstate;
Expand All @@ -264,7 +267,9 @@ AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE auto EOS<problem_t>::ComputeOtherDeriva
deint_dRho = chemstate.dedr;
deint_dP = 1.0 / chemstate.dpde;
dRho_dP = 1.0 / (chemstate.dpdr * C::k_B / boltzmann_constant_);
G_gamma = (chemstate.cs * chemstate.cs) * rho / P;
dP_dRho_s = chemstate.cs * chemstate.cs;
G_gamma = dP_dRho_s * rho / P;
d2P_dRho2_s = chemstate.d2pdr2_s;

#else
if constexpr (gamma_ != 1.0) {
Expand All @@ -276,10 +281,12 @@ AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE auto EOS<problem_t>::ComputeOtherDeriva
deint_dRho = estate.dedr;
deint_dP = 1.0 / estate.dpde;
dRho_dP = 1.0 / (estate.dpdr * C::k_B / boltzmann_constant_);
G_gamma = (estate.cs * estate.cs) * rho / P;
dP_dRho_s = estate.cs * estate.cs;
G_gamma = dP_dRho_s * rho / P;
d2P_dRho2_s = estate.d2pdr2_s;
}
#endif
return std::make_tuple(deint_dRho, deint_dP, dRho_dP, G_gamma);
return std::make_tuple(deint_dRho, deint_dP, dRho_dP, dP_dRho_s, G_gamma, d2P_dRho2_s);
}

template <typename problem_t>
Expand Down
8 changes: 4 additions & 4 deletions src/HLLC.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -44,8 +44,8 @@ AMREX_FORCE_INLINE AMREX_GPU_DEVICE auto HLLC(quokka::HydroState<N_scalars, N_ms
double S_R = NAN;
if (gamma != 1.0) {

auto [dedr_L, dedp_L, drdp_L, G_gamma_L] = quokka::EOS<problem_t>::ComputeOtherDerivatives(sL.rho, sL.P, sL.massScalar);
auto [dedr_R, dedp_R, drdp_R, G_gamma_R] = quokka::EOS<problem_t>::ComputeOtherDerivatives(sR.rho, sR.P, sR.massScalar);
auto [dedr_L, dedp_L, drdp_L, dpdr_s_L, G_gamma_L, d2pdr2_s_L] = quokka::EOS<problem_t>::ComputeOtherDerivatives(sL.rho, sL.P, sL.massScalar);
auto [dedr_R, dedp_R, drdp_R, dpdr_s_R, G_gamma_R, d2pdr2_s_R] = quokka::EOS<problem_t>::ComputeOtherDerivatives(sR.rho, sR.P, sR.massScalar);

// equation A.5a of Kershaw+1998
// need specific internal energy here
Expand All @@ -57,8 +57,8 @@ AMREX_FORCE_INLINE AMREX_GPU_DEVICE auto HLLC(quokka::HydroState<N_scalars, N_ms
// equation 4.12 of Kershaw+1998
cs_tilde = std::sqrt((1.0 / C_tilde_P) * (H_tilde - 0.5 * vsq_tilde - C_tilde_rho));

const double G_L = 0.5 * (G_gamma_L + 1.); // 'fundamental derivative' for ideal gases
const double G_R = 0.5 * (G_gamma_R + 1.);
const double G_L = 0.5 * (G_gamma_L + 1.0 + (std::pow(sL.rho, 2)/(sL.P*G_gamma_L))*d2pdr2_s_L + (sL.rho/(sL.P*G_gamma_L))*dpdr_s_L - (std::pow(sL.rho/sL.P, 2)/G_gamma_L)*std::pow(dpdr_s_L, 2)); // 'fundamental derivative' for ideal gases
const double G_R = 0.5 * (G_gamma_R + 1.0 + (std::pow(sR.rho, 2)/(sR.P*G_gamma_R))*d2pdr2_s_R + (sR.rho/(sR.P*G_gamma_R))*dpdr_s_R - (std::pow(sR.rho/sR.P, 2)/G_gamma_R)*std::pow(dpdr_s_R, 2));

const double s_NL = 0.5 * G_L * std::max(dU, 0.); // second-order wavespeed correction
const double s_NR = 0.5 * G_R * std::max(dU, 0.);
Expand Down

0 comments on commit fe202ce

Please sign in to comment.