SpECTRE Documentation Coverage Report
Current view: top level - Evolution/Systems/GrMhd/ValenciaDivClean - System.hpp Hit Total Coverage
Commit: 35a1e98cd3e4fdea528eb8100f99c2f707894fda Lines: 1 19 5.3 %
Date: 2024-04-19 00:10:48
Legend: Lines: hit not hit

          Line data    Source code
       1           0 : // Distributed under the MIT License.
       2             : // See LICENSE.txt for details.
       3             : 
       4             : #pragma once
       5             : 
       6             : #include <cstddef>
       7             : 
       8             : #include "DataStructures/VariablesTag.hpp"
       9             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/BoundaryConditions/BoundaryCondition.hpp"
      10             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/BoundaryCorrections/BoundaryCorrection.hpp"
      11             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Characteristics.hpp"
      12             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/ConservativeFromPrimitive.hpp"
      13             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/NewmanHamlin.hpp"
      14             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/PrimitiveFromConservative.hpp"
      15             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Tags.hpp"
      16             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/TimeDerivativeTerms.hpp"
      17             : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
      18             : #include "PointwiseFunctions/Hydro/Tags.hpp"
      19             : #include "Utilities/TMPL.hpp"
      20             : 
      21             : /// \ingroup EvolutionSystemsGroup
      22             : /// \brief Items related to general relativistic magnetohydrodynamics (GRMHD)
      23             : namespace grmhd {
      24             : /// Items related to the Valencia formulation of ideal GRMHD with divergence
      25             : /// cleaning coupled with electron fraction.
      26             : ///
      27             : /// References:
      28             : /// - Numerical 3+1 General Relativistic Magnetohydrodynamics: A Local
      29             : /// Characteristic Approach \cite Anton2006
      30             : /// - GRHydro: a new open-source general-relativistic magnetohydrodynamics code
      31             : /// for the Einstein toolkit \cite Moesta2014
      32             : /// - Black hole-neutron star mergers with a hot nuclear equation of state:
      33             : /// outflow and neutrino-cooled disk for a low-mass, high-spin case
      34             : /// \cite Deaton2013
      35             : ///
      36             : namespace ValenciaDivClean {
      37             : 
      38             : /*!
      39             :  * \brief Ideal general relativistic magnetohydrodynamics (GRMHD) system with
      40             :  * divergence cleaning coupled with electron fraction.
      41             :  *
      42             :  * We adopt the standard 3+1 form of metric
      43             :  * \f{align*}
      44             :  *  ds^2 = -\alpha^2 dt^2 + \gamma_{ij}(dx^i + \beta^i dt)(dx^j + \beta^j dt)
      45             :  * \f}
      46             :  * where \f$\alpha\f$ is the lapse, \f$\beta^i\f$ is the shift vector, and
      47             :  * \f$\gamma_{ij}\f$ is the spatial metric.
      48             :  *
      49             :  * Primitive variables of the system are
      50             :  *  - rest mass density \f$\rho\f$
      51             :  *  - electron fraction \f$Y_e\f$
      52             :  *  - the spatial velocity \f{align*} v^i =
      53             :  * \frac{1}{\alpha}\left(\frac{u^i}{u^0} + \beta^i\right) \f} measured by the
      54             :  * Eulerian observer
      55             :  *  - specific internal energy density \f$\epsilon\f$
      56             :  *  - magnetic field \f$B^i = -^*F^{ia}n_a = -\alpha (^*F^{0i})\f$
      57             :  *  - the divergence cleaning field \f$\Phi\f$
      58             :  *
      59             :  * with corresponding derived physical quantities which frequently appear in
      60             :  * equations:
      61             :  * - The transport velocity \f$v^i_{tr} = \alpha v^i - \beta^i\f$
      62             :  * - The Lorentz factor \f{align*} W = -u^an_a = \alpha u^0 =
      63             :  * \frac{1}{\sqrt{1-\gamma_{ij}v^iv^j}} = \sqrt{1+\gamma^{ij}u_iu_j} =
      64             :  * \sqrt{1+\gamma^{ij}W^2v_iv_j} \f}
      65             :  * - The specific enthalpy \f$h = 1 + \epsilon + p/\rho\f$ where \f$p\f$ is the
      66             :  * pressure specified by a particular equation of state (EoS)
      67             :  * - The comoving magnetic field \f$b^b = -^*F^{ba} u_a\f$ in the ideal MHD
      68             :  * limit
      69             :  * \f{align*}
      70             :  *  b^0 & = W B^i v_i / \alpha \\
      71             :  *  b^i & = (B^i + \alpha b^0 u^i) / W \\
      72             :  *  b^2 & = b^ab_a = B^2/W^2 + (B^iv_i)^2
      73             :  * \f}
      74             :  * - Augmented enthalpy density \f$(\rho h)^* = \rho h + b^2\f$ and augmented
      75             :  * pressure \f$p^* = p + b^2/2\f$ which include contributions from magnetic
      76             :  * field
      77             :  *
      78             :  * \note We are using the electromagnetic variables with the scaling convention
      79             :  * that the factor \f$4\pi\f$ does not appear in Maxwell equations and the
      80             :  * stress-energy tensor of the EM fields (geometrized Heaviside-Lorentz units).
      81             :  * To recover the physical value of magnetic field in the usual CGS Gaussian
      82             :  * unit, the conversion factor is
      83             :  * \f{align*}
      84             :  *  \sqrt{4\pi}\frac{c^4}{G^{3/2}M_\odot} \approx 8.35 \times 10^{19}
      85             :  *  \,\text{Gauss}
      86             :  * \f}
      87             :  * For example, magnetic field $10^{-5}$ with the code unit corresponds to the
      88             :  * $8.35\times 10^{14}\,\text{G}$ in the CGS Gaussian unit. See also
      89             :  * documentation of hydro::units::cgs::gauss_unit for details.
      90             :  *
      91             :  * The GRMHD equations can be written in a flux-balanced form
      92             :  * \f[
      93             :  *  \partial_t U+ \partial_i F^i(U) = S(U).
      94             :  * \f]
      95             :  *
      96             :  * Evolved (conserved) variables \f$U\f$ are
      97             :  * \f{align*}{
      98             :  * U = \sqrt{\gamma}\left[\,\begin{matrix}
      99             :  *      D \\
     100             :  *      D Y_e \\
     101             :  *      S_j \\
     102             :  *      \tau \\
     103             :  *      B^j \\
     104             :  *      \Phi \\
     105             :  * \end{matrix}\,\right] \equiv \left[\,\,\begin{matrix}
     106             :  *      \tilde{D} \\
     107             :  *      \tilde{Y_e} \\
     108             :  *      \tilde{S}_j \\
     109             :  *      \tilde{\tau} \\
     110             :  *      \tilde{B}^j \\
     111             :  *      \tilde{\Phi} \\
     112             :  * \end{matrix}\,\,\right] = \sqrt{\gamma} \left[\,\,\begin{matrix}
     113             :  *      \rho W \\
     114             :  *      \rho W Y_e \\
     115             :  *      (\rho h)^* W^2 v_j - \alpha b^0 b_j \\
     116             :  *      (\rho h)^* W^2 - p^* - (\alpha b^0)^2 - \rho W \\
     117             :  *      B^j \\
     118             :  *      \Phi \\
     119             :  * \end{matrix}\,\,\right]
     120             :  * \f}
     121             :  *
     122             :  * where \f${\tilde D}\f$, \f${\tilde Y}_e\f$,\f${\tilde S}_i\f$, \f${\tilde
     123             :  * \tau}\f$, \f${\tilde B}^i\f$, and \f${\tilde \Phi}\f$ are a generalized
     124             :  * mass-energy density, electron fraction, momentum density, specific internal
     125             :  * energy density, magnetic field, and divergence cleaning field. Also,
     126             :  * \f$\gamma\f$ is the determinant of the spatial metric \f$\gamma_{ij}\f$.
     127             :  *
     128             :  * Corresponding fluxes \f$F^i(U)\f$ are
     129             :  *
     130             :  * \f{align*}
     131             :  * F^i({\tilde D}) &= {\tilde D} v^i_{tr} \\
     132             :  * F^i({\tilde Y}_e) &= {\tilde Y}_e v^i_{tr} \\
     133             :  * F^i({\tilde S}_j) &= {\tilde S}_j v^i_{tr} + \alpha \sqrt{\gamma} p^*
     134             :  *      \delta^i_j - \alpha b_j \tilde{B}^i / W \\
     135             :  * F^i({\tilde \tau}) &= {\tilde \tau} v^i_{tr} + \alpha \sqrt{\gamma} p^* v^i
     136             :  *      - \alpha^2 b^0 \tilde{B}^i / W \\
     137             :  * F^i({\tilde B}^j) &= {\tilde B}^j v^i_{tr} - \alpha v^j {\tilde B}^i
     138             :  *      + \alpha \gamma^{ij} {\tilde \Phi} \\
     139             :  * F^i({\tilde \Phi}) &= \alpha {\tilde B^i} - {\tilde \Phi} \beta^i
     140             :  * \f}
     141             :  *
     142             :  * and source terms \f$S(U)\f$ are
     143             :  *
     144             :  * \f{align*}
     145             :  * S({\tilde D}) &= 0 \\
     146             :  * S({\tilde Y}_e) &= 0 \\
     147             :  * S({\tilde S}_j) &= \frac{\alpha}{2} {\tilde S}^{mn} \partial_j \gamma_{mn}
     148             :  *      + {\tilde S}_k \partial_j \beta^k - ({\tilde D}
     149             :  *      + {\tilde \tau}) \partial_j \alpha \\
     150             :  * S({\tilde \tau}) &= \alpha {\tilde S}^{mn} K_{mn}
     151             :  *      - {\tilde S}^m \partial_m \alpha \\
     152             :  * S({\tilde B}^j) &= -{\tilde B}^m \partial_m \beta^j
     153             :  *      + \Phi \partial_k (\alpha \sqrt{\gamma}\gamma^{jk}) \\
     154             :  * S({\tilde \Phi}) &= {\tilde B}^k \partial_k \alpha - \alpha K
     155             :  *      {\tilde \Phi} - \alpha \kappa {\tilde \Phi}
     156             :  * \f}
     157             :  *
     158             :  * with
     159             :  * \f{align*}
     160             :  * {\tilde S}^{ij} = & \sqrt{\gamma} \left[ (\rho h)^* W^2 v^i v^j + p^*
     161             :  *      \gamma^{ij} - \gamma^{ik}\gamma^{jl}b_kb_l \right]
     162             :  * \f}
     163             :  *
     164             :  * where \f$K\f$ is the trace of the extrinsic curvature \f$K_{ij}\f$, and
     165             :  * \f$\kappa\f$ is a damping parameter that damps violations of the
     166             :  * divergence-free (no-monopole) condition \f$\Phi = \partial_i {\tilde B}^i =
     167             :  * 0\f$.
     168             :  *
     169             :  * \note On the electron fraction side, the source term is currently set to
     170             :  * \f$S(\tilde{Y}_e) = 0\f$. Implementing the source term using neutrino scheme
     171             :  * is in progress (Last update : Oct 2022).
     172             :  *
     173             :  */
     174           1 : struct System {
     175           0 :   static constexpr bool is_in_flux_conservative_form = true;
     176           0 :   static constexpr bool has_primitive_and_conservative_vars = true;
     177           0 :   static constexpr size_t volume_dim = 3;
     178             : 
     179           0 :   using boundary_conditions_base = BoundaryConditions::BoundaryCondition;
     180           0 :   using boundary_correction_base = BoundaryCorrections::BoundaryCorrection;
     181             : 
     182           0 :   using variables_tag = ::Tags::Variables<
     183             :       tmpl::list<Tags::TildeD, Tags::TildeYe, Tags::TildeTau, Tags::TildeS<>,
     184             :                  Tags::TildeB<>, Tags::TildePhi>>;
     185           0 :   using flux_variables =
     186             :       tmpl::list<Tags::TildeD, Tags::TildeYe, Tags::TildeTau, Tags::TildeS<>,
     187             :                  Tags::TildeB<>, Tags::TildePhi>;
     188           0 :   using non_conservative_variables = tmpl::list<>;
     189           0 :   using gradient_variables = tmpl::list<>;
     190           0 :   using primitive_variables_tag =
     191             :       ::Tags::Variables<hydro::grmhd_tags<DataVector>>;
     192           0 :   using spacetime_variables_tag =
     193             :       ::Tags::Variables<gr::tags_for_hydro<volume_dim, DataVector>>;
     194           0 :   using flux_spacetime_variables_tag = ::Tags::Variables<
     195             :       tmpl::list<gr::Tags::Lapse<DataVector>, gr::Tags::Shift<DataVector, 3>,
     196             :                  gr::Tags::SpatialMetric<DataVector, 3>,
     197             :                  gr::Tags::SqrtDetSpatialMetric<DataVector>,
     198             :                  gr::Tags::InverseSpatialMetric<DataVector, 3>>>;
     199             : 
     200           0 :   using compute_volume_time_derivative_terms = TimeDerivativeTerms;
     201             : 
     202           0 :   using conservative_from_primitive = ConservativeFromPrimitive;
     203             :   template <typename OrderedListOfPrimitiveRecoverySchemes>
     204           0 :   using primitive_from_conservative =
     205             :       PrimitiveFromConservative<OrderedListOfPrimitiveRecoverySchemes>;
     206             : 
     207           0 :   using compute_largest_characteristic_speed =
     208             :       Tags::ComputeLargestCharacteristicSpeed;
     209             : 
     210           0 :   using inverse_spatial_metric_tag =
     211             :       gr::Tags::InverseSpatialMetric<DataVector, volume_dim>;
     212             : };
     213             : }  // namespace ValenciaDivClean
     214             : }  // namespace grmhd

Generated by: LCOV version 1.14