Line data Source code
1 0 : // Distributed under the MIT License.
2 : // See LICENSE.txt for details.
3 :
4 : #pragma once
5 :
6 : #include <array>
7 : #include <cstddef>
8 : #include <utility>
9 :
10 : #include "DataStructures/DataBox/Prefixes.hpp"
11 : #include "DataStructures/Tensor/TypeAliases.hpp"
12 : #include "Domain/Structure/DirectionalIdMap.hpp"
13 : #include "Evolution/Systems/GeneralizedHarmonic/Tags.hpp"
14 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Tags.hpp"
15 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/FiniteDifference/ReconstructWork.hpp"
16 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/TagsDeclarations.hpp"
17 : #include "Evolution/VariableFixing/FixToAtmosphere.hpp"
18 : #include "PointwiseFunctions/GeneralRelativity/TagsDeclarations.hpp"
19 : #include "PointwiseFunctions/Hydro/TagsDeclarations.hpp"
20 : #include "Utilities/TMPL.hpp"
21 :
22 : /// \cond
23 : class DataVector;
24 : template <typename TagsList>
25 : class Variables;
26 : namespace gsl {
27 : template <typename>
28 : class not_null;
29 : } // namespace gsl
30 : template <size_t Dim>
31 : class Direction;
32 : template <size_t Dim>
33 : class ElementId;
34 : template <size_t Dim>
35 : class Element;
36 : template <size_t Dim>
37 : class Mesh;
38 : namespace EquationsOfState {
39 : template <bool IsRelativistic, size_t ThermodynamicDim>
40 : class EquationOfState;
41 : } // namespace EquationsOfState
42 : namespace evolution::dg::subcell {
43 : class GhostData;
44 : } // namespace evolution::dg::subcell
45 : namespace evolution::dg::Actions::detail {
46 : template <size_t Dim>
47 : struct NormalVector;
48 : } // namespace evolution::dg::Actions::detail
49 : namespace gh::Tags {
50 : struct ConstraintGamma1;
51 : struct ConstraintGamma2;
52 : } // namespace gh::Tags
53 : /// \endcond
54 :
55 : namespace grmhd::GhValenciaDivClean::fd {
56 0 : using tags_list_for_reconstruct =
57 : tmpl::push_front<grmhd::ValenciaDivClean::fd::tags_list_for_reconstruct,
58 : gr::Tags::SpacetimeMetric<DataVector, 3>,
59 : gh::Tags::Pi<DataVector, 3>, gh::Tags::Phi<DataVector, 3>>;
60 :
61 : namespace detail {
62 : using tags_list_for_reconstruct_split_lapse =
63 : tmpl::split<tags_list_for_reconstruct, gr::Tags::Lapse<DataVector>>;
64 : } // namespace detail
65 :
66 0 : using tags_list_for_reconstruct_fd_neighbor = tmpl::append<
67 : tmpl::front<detail::tags_list_for_reconstruct_split_lapse>,
68 : tmpl::push_front<tmpl::back<detail::tags_list_for_reconstruct_split_lapse>,
69 : gh::Tags::ConstraintGamma1, gh::Tags::ConstraintGamma2,
70 : gr::Tags::Lapse<DataVector>>>;
71 :
72 : /*!
73 : * \brief Reconstructs \f$\rho, p, Wv^i, B^i\f$, \f$\Phi\f$, and the spacetime
74 : * metric, then computes the Lorentz factor, upper spatial velocity, specific
75 : * internal energy, and the conserved variables. All results are written into
76 : * `vars_on_lower_face` and `vars_on_upper_face`.
77 : */
78 : template <typename SpacetimeTagsToReconstruct,
79 : typename PrimTagsForReconstruction, typename PrimsTags,
80 : typename SpacetimeAndConsTags, size_t ThermodynamicDim,
81 : typename HydroReconstructor, typename SpacetimeReconstructor,
82 : typename ComputeGrmhdSpacetimeVarsFromReconstructedSpacetimeTags,
83 : typename PrimsTagsSentByNeighbor>
84 1 : void reconstruct_prims_work(
85 : gsl::not_null<std::array<Variables<tags_list_for_reconstruct>, 3>*>
86 : vars_on_lower_face,
87 : gsl::not_null<std::array<Variables<tags_list_for_reconstruct>, 3>*>
88 : vars_on_upper_face,
89 : const HydroReconstructor& hydro_reconstructor,
90 : const SpacetimeReconstructor& spacetime_reconstructor,
91 : const ComputeGrmhdSpacetimeVarsFromReconstructedSpacetimeTags&
92 : spacetime_vars_for_grmhd,
93 : const Variables<PrimsTags>& volume_prims,
94 : const Variables<SpacetimeAndConsTags>& volume_spacetime_and_cons_vars,
95 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>& eos,
96 : const Element<3>& element,
97 : const DirectionalIdMap<3, Variables<PrimsTagsSentByNeighbor>>&
98 : neighbor_data,
99 : const Mesh<3>& subcell_mesh, size_t ghost_zone_size,
100 : bool compute_conservatives, bool reconstruct_density_times_temperature,
101 : const VariableFixing::FixToAtmosphere<3>* fix_to_atmosphere);
102 :
103 : /*!
104 : * \brief Reconstructs \f$\rho, p, Wv^i, B^i\f$, \f$\Phi\f$, the spacetime
105 : * metric, \f$\Phi_{iab}\f$, and \f$\Pi_{ab}\f$, then computes the Lorentz
106 : * factor, upper spatial velocity, specific internal energy, and the conserved
107 : * variables. All results are written into `vars_on_face`.
108 : *
109 : * This is used on DG elements to reconstruct their subcell neighbors' solution
110 : * on the shared faces.
111 : */
112 : template <
113 : typename SpacetimeTagsToReconstruct, typename PrimTagsForReconstruction,
114 : typename PrimsTagsSentByNeighbor, typename PrimsTags,
115 : size_t ThermodynamicDim, typename LowerHydroReconstructor,
116 : typename LowerSpacetimeReconstructor, typename UpperHydroReconstructor,
117 : typename UpperSpacetimeReconstructor,
118 : typename ComputeGrmhdSpacetimeVarsFromReconstructedSpacetimeTags>
119 1 : void reconstruct_fd_neighbor_work(
120 : gsl::not_null<Variables<tags_list_for_reconstruct_fd_neighbor>*>
121 : vars_on_face,
122 : const LowerHydroReconstructor& reconstruct_lower_neighbor_hydro,
123 : const LowerSpacetimeReconstructor& reconstruct_lower_neighbor_spacetime,
124 : const UpperHydroReconstructor& reconstruct_upper_neighbor_hydro,
125 : const UpperSpacetimeReconstructor& reconstruct_upper_neighbor_spacetime,
126 : const ComputeGrmhdSpacetimeVarsFromReconstructedSpacetimeTags&
127 : spacetime_vars_for_grmhd,
128 : const Variables<PrimsTags>& subcell_volume_prims,
129 : const Variables<
130 : grmhd::GhValenciaDivClean::Tags::spacetime_reconstruction_tags>&
131 : subcell_volume_spacetime_vars,
132 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>& eos,
133 : const Element<3>& element,
134 : const DirectionalIdMap<3, evolution::dg::subcell::GhostData>& ghost_data,
135 : const Mesh<3>& subcell_mesh, const Direction<3>& direction_to_reconstruct,
136 : size_t ghost_zone_size, bool compute_conservatives,
137 : bool reconstruct_density_times_temperature,
138 : const VariableFixing::FixToAtmosphere<3>* fix_to_atmosphere);
139 : } // namespace grmhd::GhValenciaDivClean::fd
|