Line data Source code
1 0 : // Distributed under the MIT License.
2 : // See LICENSE.txt for details.
3 :
4 : #pragma once
5 :
6 : #include "DataStructures/Tensor/TypeAliases.hpp"
7 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/TagsDeclarations.hpp" // IWYU pragma: keep
8 : #include "PointwiseFunctions/GeneralRelativity/TagsDeclarations.hpp" // IWYU pragma: keep
9 : #include "PointwiseFunctions/Hydro/TagsDeclarations.hpp" // IWYU pragma: keep
10 : #include "Utilities/TMPL.hpp"
11 :
12 : // IWYU pragma: no_include "Evolution/Systems/GrMhd/ValenciaDivClean/Tags.hpp"
13 : // IWYU pragma: no_include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
14 : // IWYU pragma: no_include "PointwiseFunctions/Hydro/Tags.hpp"
15 :
16 : /// \cond
17 : namespace gsl {
18 : template <typename T>
19 : class not_null;
20 : } // namespace gsl
21 :
22 : class DataVector;
23 : /// \endcond
24 :
25 : // IWYU pragma: no_forward_declare Tensor
26 :
27 : namespace grmhd {
28 : namespace ValenciaDivClean {
29 :
30 : /*!
31 : * \brief Compute the conservative variables from primitive variables
32 : *
33 : * \f{align*}
34 : * {\tilde D} = & \sqrt{\gamma} \rho W \\
35 : * {\tilde Y}_e = & \sqrt{\gamma} \rho Y_e W \\
36 : * {\tilde S}_i = & \sqrt{\gamma} \left( \rho h W^2 v_i + B^m B_m v_i - B^m v_m
37 : * B_i \right) \\
38 : * {\tilde \tau} = & \sqrt{\gamma} \left[ \rho h W^2 - p - \rho W - \frac{1}{2}
39 : * (B^m v_m)^2 + \frac{1}{2} B^m B_m \left( 1 + v^m v_m \right) \right] \\
40 : * {\tilde B}^i = & \sqrt{\gamma} B^i \\
41 : * {\tilde \Phi} = & \sqrt{\gamma} \Phi
42 : * \f}
43 : *
44 : * where the conserved variables \f${\tilde D}\f$, \f$\tilde{Y}_e\f$, \f${\tilde
45 : * S}_i\f$, \f${\tilde \tau}\f$, \f${\tilde B}^i\f$, and \f${\tilde \Phi}\f$ are
46 : * a generalized mass-energy density, electron fraction, momentum density,
47 : * specific internal energy density, magnetic field, and divergence cleaning
48 : * field. Furthermore \f$\gamma\f$ is the determinant of the spatial metric,
49 : * \f$\rho\f$ is the rest mass density, \f$Y_e\f$ is the electron fraction,
50 : * \f$W = 1/\sqrt{1-v_i v^i}\f$ is the Lorentz factor, \f$h = 1 + \epsilon +
51 : * \frac{p}{\rho}\f$ is the specific enthalpy, \f$v^i\f$ is the spatial
52 : * velocity, \f$\epsilon\f$ is the specific internal energy, \f$p\f$ is the
53 : * pressure, \f$B^i\f$ is the spatial magnetic field measured by an Eulerian
54 : * observer, and \f$\Phi\f$ is a divergence cleaning field.
55 : *
56 : * The quantity \f${\tilde \tau}\f$ is rewritten as in `RelativisticEuler`
57 : * to avoid cancellation error in the non-relativistic limit:
58 : * \f[
59 : * \left( \rho h W^2 - p - \rho W \right) \longrightarrow
60 : * W^2 \left[ \rho \left( \epsilon + v^2
61 : * \frac{W}{W + 1} \right) + p v^2 \right] .\f]
62 : */
63 1 : struct ConservativeFromPrimitive {
64 0 : using return_tags = tmpl::list<grmhd::ValenciaDivClean::Tags::TildeD,
65 : grmhd::ValenciaDivClean::Tags::TildeYe,
66 : grmhd::ValenciaDivClean::Tags::TildeTau,
67 : grmhd::ValenciaDivClean::Tags::TildeS<>,
68 : grmhd::ValenciaDivClean::Tags::TildeB<>,
69 : grmhd::ValenciaDivClean::Tags::TildePhi>;
70 :
71 0 : using argument_tags =
72 : tmpl::list<hydro::Tags::RestMassDensity<DataVector>,
73 : hydro::Tags::ElectronFraction<DataVector>,
74 : hydro::Tags::SpecificInternalEnergy<DataVector>,
75 : hydro::Tags::Pressure<DataVector>,
76 : hydro::Tags::SpatialVelocity<DataVector, 3>,
77 : hydro::Tags::LorentzFactor<DataVector>,
78 : hydro::Tags::MagneticField<DataVector, 3>,
79 : gr::Tags::SqrtDetSpatialMetric<DataVector>,
80 : gr::Tags::SpatialMetric<DataVector, 3>,
81 : hydro::Tags::DivergenceCleaningField<DataVector>>;
82 :
83 0 : static void apply(
84 : gsl::not_null<Scalar<DataVector>*> tilde_d,
85 : gsl::not_null<Scalar<DataVector>*> tilde_ye,
86 : gsl::not_null<Scalar<DataVector>*> tilde_tau,
87 : gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*> tilde_s,
88 : gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*> tilde_b,
89 : gsl::not_null<Scalar<DataVector>*> tilde_phi,
90 : const Scalar<DataVector>& rest_mass_density,
91 : const Scalar<DataVector>& electron_fraction,
92 : const Scalar<DataVector>& specific_internal_energy,
93 : const Scalar<DataVector>& pressure,
94 : const tnsr::I<DataVector, 3, Frame::Inertial>& spatial_velocity,
95 : const Scalar<DataVector>& lorentz_factor,
96 : const tnsr::I<DataVector, 3, Frame::Inertial>& magnetic_field,
97 : const Scalar<DataVector>& sqrt_det_spatial_metric,
98 : const tnsr::ii<DataVector, 3, Frame::Inertial>& spatial_metric,
99 : const Scalar<DataVector>& divergence_cleaning_field);
100 : };
101 : } // namespace ValenciaDivClean
102 : } // namespace grmhd
|