SpECTRE  v2024.04.12
grmhd::ValenciaDivClean Namespace Reference

Items related to the Valencia formulation of ideal GRMHD with divergence cleaning coupled with electron fraction. More...

Namespaces

namespace  BoundaryConditions
 Boundary conditions for the GRMHD Valencia Divergence Cleaning system.
 
namespace  BoundaryCorrections
 Boundary corrections/numerical fluxes.
 
namespace  fd
 Finite difference functionality for the ValenciaDivClean form of the GRMHD equations.
 
namespace  PrimitiveRecoverySchemes
 Schemes for recovering primitive variables from conservative variables.
 
namespace  subcell
 Code required by the DG-subcell/FD hybrid solver.
 
namespace  Tags
 Tags for the Valencia formulation of the ideal GRMHD equations with divergence cleaning.
 

Classes

struct  ComputeFluxes
 The fluxes of the conservative variables. More...
 
struct  ComputeSources
 Compute the source terms for the flux-conservative Valencia formulation of GRMHD with divergence cleaning, coupled with electron fraction. More...
 
struct  ConservativeFromPrimitive
 Compute the conservative variables from primitive variables. More...
 
class  FixConservatives
 Fix conservative variables using method developed by Foucart. More...
 
class  Flattener
 Reduces oscillations inside an element in an attempt to guarantee a physical solution of the conserved variables for which the primitive variables can be recovered. More...
 
class  NumericInitialData
 Numeric initial data loaded from volume data files. More...
 
struct  PrimitiveFromConservative
 Compute the primitive variables from the conservative variables. More...
 
class  PrimitiveFromConservativeOptions
 Options to be passed to the Con2Prim algorithm. Currently, we simply set a threshold for tildeD below which the inversion is not performed and the density is set to atmosphere values. More...
 
struct  SetVariablesNeededFixingToFalse
 Mutator used with Initialization::Actions::AddSimpleTags to initialize the VariablesNeededFixing to false More...
 
struct  System
 Ideal general relativistic magnetohydrodynamics (GRMHD) system with divergence cleaning coupled with electron fraction. More...
 
struct  TimeDerivativeTerms
 Compute the time derivative of the conserved variables for the Valencia formulation of the GRMHD equations with divergence cleaning. More...
 

Functions

bool operator!= (const FixConservatives &lhs, const FixConservatives &rhs)
 
template<typename RecoverySchemesList >
bool operator!= (const Flattener< RecoverySchemesList > &lhs, const Flattener< RecoverySchemesList > &rhs)
 
bool operator!= (const PrimitiveFromConservativeOptions &lhs, const PrimitiveFromConservativeOptions &rhs)
 
template<size_t ThermodynamicDim>
std::array< DataVector, 9 > characteristic_speeds (const Scalar< DataVector > &rest_mass_density, const Scalar< DataVector > &electron_fraction, const Scalar< DataVector > &specific_internal_energy, const Scalar< DataVector > &specific_enthalpy, const tnsr::I< DataVector, 3, Frame::Inertial > &spatial_velocity, const Scalar< DataVector > &lorentz_factor, const tnsr::I< DataVector, 3, Frame::Inertial > &magnetic_field, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3 > &shift, const tnsr::ii< DataVector, 3, Frame::Inertial > &spatial_metric, const tnsr::i< DataVector, 3 > &unit_normal, const EquationsOfState::EquationOfState< true, ThermodynamicDim > &equation_of_state)
 Compute the characteristic speeds for the Valencia formulation of GRMHD with divergence cleaning. More...
 
template<size_t ThermodynamicDim>
void characteristic_speeds (gsl::not_null< std::array< DataVector, 9 > * > char_speeds, const Scalar< DataVector > &rest_mass_density, const Scalar< DataVector > &electron_fraction, const Scalar< DataVector > &specific_internal_energy, const Scalar< DataVector > &specific_enthalpy, const tnsr::I< DataVector, 3, Frame::Inertial > &spatial_velocity, const Scalar< DataVector > &lorentz_factor, const tnsr::I< DataVector, 3, Frame::Inertial > &magnetic_field, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3 > &shift, const tnsr::ii< DataVector, 3, Frame::Inertial > &spatial_metric, const tnsr::i< DataVector, 3 > &unit_normal, const EquationsOfState::EquationOfState< true, ThermodynamicDim > &equation_of_state)
 Compute the characteristic speeds for the Valencia formulation of GRMHD with divergence cleaning. More...
 

Detailed Description

Items related to the Valencia formulation of ideal GRMHD with divergence cleaning coupled with electron fraction.

References:

  • Numerical 3+1 General Relativistic Magnetohydrodynamics: A Local Characteristic Approach [4]
  • GRHydro: a new open-source general-relativistic magnetohydrodynamics code for the Einstein toolkit [136]
  • Black hole-neutron star mergers with a hot nuclear equation of state: outflow and neutrino-cooled disk for a low-mass, high-spin case [45]

Function Documentation

◆ characteristic_speeds() [1/2]

template<size_t ThermodynamicDim>
std::array< DataVector, 9 > grmhd::ValenciaDivClean::characteristic_speeds ( const Scalar< DataVector > &  rest_mass_density,
const Scalar< DataVector > &  electron_fraction,
const Scalar< DataVector > &  specific_internal_energy,
const Scalar< DataVector > &  specific_enthalpy,
const tnsr::I< DataVector, 3, Frame::Inertial > &  spatial_velocity,
const Scalar< DataVector > &  lorentz_factor,
const tnsr::I< DataVector, 3, Frame::Inertial > &  magnetic_field,
const Scalar< DataVector > &  lapse,
const tnsr::I< DataVector, 3 > &  shift,
const tnsr::ii< DataVector, 3, Frame::Inertial > &  spatial_metric,
const tnsr::i< DataVector, 3 > &  unit_normal,
const EquationsOfState::EquationOfState< true, ThermodynamicDim > &  equation_of_state 
)

Compute the characteristic speeds for the Valencia formulation of GRMHD with divergence cleaning.

Obtaining the exact form of the characteristic speeds involves the solution of a nontrivial quartic equation for the fast and slow modes. Here we make use of a common approximation in the literature (e.g. [73]) where the resulting characteristic speeds are analogous to those of the Valencia formulation of the 3-D relativistic Euler system (see RelativisticEuler::Valencia::characteristic_speeds),

\begin{align*} \lambda_2 &= \alpha \Lambda^- - \beta_n,\\ \lambda_{3, 4, 5, 6, 7} &= \alpha v_n - \beta_n,\\ \lambda_{8} &= \alpha \Lambda^+ - \beta_n, \end{align*}

with the substitution

\begin{align*} c_s^2 \longrightarrow c_s^2 + v_A^2(1 - c_s^2) \end{align*}

in the definition of \(\Lambda^\pm\). Here \(v_A\) is the Alfvén speed. In addition, two more speeds corresponding to the divergence cleaning mode and the longitudinal magnetic field are added,

\begin{align*} \lambda_1 = -\alpha - \beta_n,\\ \lambda_9 = \alpha - \beta_n. \end{align*}

Note
The ordering assumed here is such that, in the Newtonian limit, the exact expressions for \(\lambda_{2, 8}\), \(\lambda_{3, 7}\), and \(\lambda_{4, 6}\) should reduce to the corresponding fast modes, Alfvén modes, and slow modes, respectively. See [46] for a detailed description of the hyperbolic characterization of Newtonian MHD. In terms of the primitive variables:

\begin{align*} v^2 &= \gamma_{mn} v^m v^n \\ c_s^2 &= \frac{1}{h} \left[ \left( \frac{\partial p}{\partial \rho} \right)_\epsilon + \frac{p}{\rho^2} \left(\frac{\partial p}{\partial \epsilon} \right)_\rho \right] \\ v_A^2 &= \frac{b^2}{b^2 + \rho h} \\ b^2 &= \frac{1}{W^2} \gamma_{mn} B^m B^n + \left( \gamma_{mn} B^m v^n \right)^2 \end{align*}

where \(\gamma_{mn}\) is the spatial metric, \(\rho\) is the rest mass density, \(W = 1/\sqrt{1-v_i v^i}\) is the Lorentz factor, \(h = 1 + \epsilon + \frac{p}{\rho}\) is the specific enthalpy, \(v^i\) is the spatial velocity, \(\epsilon\) is the specific internal energy, \(p\) is the pressure, and \(B^i\) is the spatial magnetic field measured by an Eulerian observer.

◆ characteristic_speeds() [2/2]

template<size_t ThermodynamicDim>
void grmhd::ValenciaDivClean::characteristic_speeds ( gsl::not_null< std::array< DataVector, 9 > * >  char_speeds,
const Scalar< DataVector > &  rest_mass_density,
const Scalar< DataVector > &  electron_fraction,
const Scalar< DataVector > &  specific_internal_energy,
const Scalar< DataVector > &  specific_enthalpy,
const tnsr::I< DataVector, 3, Frame::Inertial > &  spatial_velocity,
const Scalar< DataVector > &  lorentz_factor,
const tnsr::I< DataVector, 3, Frame::Inertial > &  magnetic_field,
const Scalar< DataVector > &  lapse,
const tnsr::I< DataVector, 3 > &  shift,
const tnsr::ii< DataVector, 3, Frame::Inertial > &  spatial_metric,
const tnsr::i< DataVector, 3 > &  unit_normal,
const EquationsOfState::EquationOfState< true, ThermodynamicDim > &  equation_of_state 
)

Compute the characteristic speeds for the Valencia formulation of GRMHD with divergence cleaning.

Obtaining the exact form of the characteristic speeds involves the solution of a nontrivial quartic equation for the fast and slow modes. Here we make use of a common approximation in the literature (e.g. [73]) where the resulting characteristic speeds are analogous to those of the Valencia formulation of the 3-D relativistic Euler system (see RelativisticEuler::Valencia::characteristic_speeds),

\begin{align*} \lambda_2 &= \alpha \Lambda^- - \beta_n,\\ \lambda_{3, 4, 5, 6, 7} &= \alpha v_n - \beta_n,\\ \lambda_{8} &= \alpha \Lambda^+ - \beta_n, \end{align*}

with the substitution

\begin{align*} c_s^2 \longrightarrow c_s^2 + v_A^2(1 - c_s^2) \end{align*}

in the definition of \(\Lambda^\pm\). Here \(v_A\) is the Alfvén speed. In addition, two more speeds corresponding to the divergence cleaning mode and the longitudinal magnetic field are added,

\begin{align*} \lambda_1 = -\alpha - \beta_n,\\ \lambda_9 = \alpha - \beta_n. \end{align*}

Note
The ordering assumed here is such that, in the Newtonian limit, the exact expressions for \(\lambda_{2, 8}\), \(\lambda_{3, 7}\), and \(\lambda_{4, 6}\) should reduce to the corresponding fast modes, Alfvén modes, and slow modes, respectively. See [46] for a detailed description of the hyperbolic characterization of Newtonian MHD. In terms of the primitive variables:

\begin{align*} v^2 &= \gamma_{mn} v^m v^n \\ c_s^2 &= \frac{1}{h} \left[ \left( \frac{\partial p}{\partial \rho} \right)_\epsilon + \frac{p}{\rho^2} \left(\frac{\partial p}{\partial \epsilon} \right)_\rho \right] \\ v_A^2 &= \frac{b^2}{b^2 + \rho h} \\ b^2 &= \frac{1}{W^2} \gamma_{mn} B^m B^n + \left( \gamma_{mn} B^m v^n \right)^2 \end{align*}

where \(\gamma_{mn}\) is the spatial metric, \(\rho\) is the rest mass density, \(W = 1/\sqrt{1-v_i v^i}\) is the Lorentz factor, \(h = 1 + \epsilon + \frac{p}{\rho}\) is the specific enthalpy, \(v^i\) is the spatial velocity, \(\epsilon\) is the specific internal energy, \(p\) is the pressure, and \(B^i\) is the spatial magnetic field measured by an Eulerian observer.