SpECTRE  v2025.03.17
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
grmhd::ValenciaDivClean::System Struct Reference

Ideal general relativistic magnetohydrodynamics (GRMHD) system with divergence cleaning coupled with electron fraction. More...

#include <System.hpp>

Public Types

using boundary_conditions_base = BoundaryConditions::BoundaryCondition
 
using boundary_correction_base = BoundaryCorrections::BoundaryCorrection
 
using variables_tag = ::Tags::Variables< tmpl::list< Tags::TildeD, Tags::TildeYe, Tags::TildeTau, Tags::TildeS<>, Tags::TildeB<>, Tags::TildePhi > >
 
using flux_variables = implementation defined
 
using non_conservative_variables = implementation defined
 
using gradient_variables = implementation defined
 
using primitive_variables_tag = ::Tags::Variables< hydro::grmhd_tags< DataVector > >
 
using spacetime_variables_tag = ::Tags::Variables< gr::tags_for_hydro< volume_dim, DataVector > >
 
using flux_spacetime_variables_tag = ::Tags::Variables< tmpl::list< gr::Tags::Lapse< DataVector >, gr::Tags::Shift< DataVector, 3 >, gr::Tags::SpatialMetric< DataVector, 3 >, gr::Tags::SqrtDetSpatialMetric< DataVector >, gr::Tags::InverseSpatialMetric< DataVector, 3 > > >
 
using compute_volume_time_derivative_terms = TimeDerivativeTerms
 
using conservative_from_primitive = ConservativeFromPrimitive
 
template<typename OrderedListOfPrimitiveRecoverySchemes >
using primitive_from_conservative = PrimitiveFromConservative< OrderedListOfPrimitiveRecoverySchemes >
 
using compute_largest_characteristic_speed = Tags::ComputeLargestCharacteristicSpeed
 
using inverse_spatial_metric_tag = gr::Tags::InverseSpatialMetric< DataVector, volume_dim >
 

Static Public Attributes

static constexpr bool is_in_flux_conservative_form = true
 
static constexpr bool has_primitive_and_conservative_vars = true
 
static constexpr size_t volume_dim = 3
 

Detailed Description

Ideal general relativistic magnetohydrodynamics (GRMHD) system with divergence cleaning coupled with electron fraction.

We adopt the standard 3+1 form of metric

ds2=α2dt2+γij(dxi+βidt)(dxj+βjdt)

where α is the lapse, βi is the shift vector, and γij is the spatial metric.

Primitive variables of the system are

  • rest mass density ρ
  • electron fraction Ye
  • the spatial velocity

    vi=1α(uiu0+βi)

    measured by the Eulerian observer
  • specific internal energy density ϵ
  • magnetic field Bi=Fiana=α(F0i)
  • the divergence cleaning field Φ

with corresponding derived physical quantities which frequently appear in equations:

  • The transport velocity vtri=αviβi
  • The Lorentz factor

    W=uana=αu0=11γijvivj=1+γijuiuj=1+γijW2vivj

  • The specific enthalpy h=1+ϵ+p/ρ where p is the pressure specified by a particular equation of state (EoS)
  • The comoving magnetic field bb=Fbaua in the ideal MHD limit

    b0=WBivi/αbi=(Bi+αb0ui)/Wb2=baba=B2/W2+(Bivi)2

  • Augmented enthalpy density (ρh)=ρh+b2 and augmented pressure p=p+b2/2 which include contributions from magnetic field
Note
We are using the electromagnetic variables with the scaling convention that the factor 4π does not appear in Maxwell equations and the stress-energy tensor of the EM fields (geometrized Heaviside-Lorentz units). To recover the physical value of magnetic field in the usual CGS Gaussian unit, the conversion factor is

4πc4G3/2M8.35×1019Gauss

For example, magnetic field 105 with the code unit corresponds to the 8.35×1014G in the CGS Gaussian unit. See also documentation of hydro::units::cgs::gauss_unit for details.

The GRMHD equations can be written in a flux-balanced form

tU+iFi(U)=S(U).

Evolved (conserved) variables U are

U=γ[DDYeSjτBjΦ][D~Ye~S~jτ~B~jΦ~]=γ[ρWρWYe(ρh)W2vjαb0bj(ρh)W2p(αb0)2ρWBjΦ]

where D~, Y~e, S~i, τ~, B~i, and Φ~ are a generalized mass-energy density, electron fraction, momentum density, specific internal energy density, magnetic field, and divergence cleaning field. Also, γ is the determinant of the spatial metric γij.

Corresponding fluxes Fi(U) are

Fi(D~)=D~vtriFi(Y~e)=Y~evtriFi(S~j)=S~jvtri+αγpδjiαbjB~i/WFi(τ~)=τ~vtri+αγpviα2b0B~i/WFi(B~j)=B~jvtrivtrjB~i+αγijΦ~Fi(Φ~)=αB~iΦ~βi

and source terms S(U) are

S(D~)=0S(Y~e)=0S(S~j)=α2S~mnjγmn+S~kjβk(D~+τ~)jαS(τ~)=αS~mnKmnS~mmαS(B~j)=Φk(αγγjk)S(Φ~)=B~kkααKΦ~ακΦ~

with

S~ij=γ[(ρh)W2vivj+pγijγikγjlbkbl]

where K is the trace of the extrinsic curvature Kij, and κ is a damping parameter that damps violations of the divergence-free (no-monopole) condition Φ=iB~i=0.

Note
The fluxes and sources of B~ follow eq. 32 of rather than eq. 87 of
On the electron fraction side, the source term is currently set to S(Y~e)=0. Implementing the source term using neutrino scheme is in progress (Last update : Oct 2022).

The documentation for this struct was generated from the following file: