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 <memory> 9 : #include <utility> 10 : #include <vector> 11 : 12 : #include "DataStructures/DataBox/PrefixHelpers.hpp" 13 : #include "DataStructures/DataBox/Prefixes.hpp" 14 : #include "DataStructures/Tensor/TypeAliases.hpp" 15 : #include "DataStructures/Variables.hpp" 16 : #include "DataStructures/VariablesTag.hpp" 17 : #include "Domain/Structure/DirectionalIdMap.hpp" 18 : #include "Domain/Tags.hpp" 19 : #include "Evolution/DgSubcell/Tags/GhostDataForReconstruction.hpp" 20 : #include "Evolution/DgSubcell/Tags/Mesh.hpp" 21 : #include "Evolution/Systems/GeneralizedHarmonic/Tags.hpp" 22 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Reconstructor.hpp" 23 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/System.hpp" 24 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Tags.hpp" 25 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Tags.hpp" 26 : #include "Evolution/VariableFixing/FixToAtmosphere.hpp" 27 : #include "Evolution/VariableFixing/Tags.hpp" 28 : #include "Options/String.hpp" 29 : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp" 30 : #include "PointwiseFunctions/Hydro/Tags.hpp" 31 : #include "Utilities/Serialization/CharmPupable.hpp" 32 : #include "Utilities/TMPL.hpp" 33 : 34 : /// \cond 35 : class DataVector; 36 : template <size_t Dim> 37 : class Direction; 38 : template <size_t Dim> 39 : class Element; 40 : template <size_t Dim> 41 : class ElementId; 42 : namespace EquationsOfState { 43 : template <bool IsRelativistic, size_t ThermodynamicDim> 44 : class EquationOfState; 45 : } // namespace EquationsOfState 46 : template <size_t Dim> 47 : class Mesh; 48 : namespace gsl { 49 : template <typename T> 50 : class not_null; 51 : } // namespace gsl 52 : namespace PUP { 53 : class er; 54 : } // namespace PUP 55 : namespace evolution::dg::subcell { 56 : class GhostData; 57 : } // namespace evolution::dg::subcell 58 : /// \endcond 59 : 60 : namespace grmhd::GhValenciaDivClean::fd { 61 : /*! 62 : * \brief Monotonised central reconstruction on the GRMHD primitive variables 63 : * (see 64 : * ::fd::reconstruction::monotonised_central() for details) and unlimited 3rd 65 : * order (degree 2 polynomial) reconstruction on the metric variables. 66 : * 67 : * Only the spacetime metric is reconstructed when we and the neighboring 68 : * element in the direction are doing FD. If we are doing DG and a neighboring 69 : * element is doing FD, then the spacetime metric, \f$\Phi_{iab}\f$, and 70 : * \f$\Pi_{ab}\f$ are all reconstructed since the Riemann solver on the DG 71 : * element also needs to solve for the metric variables. 72 : */ 73 : template <typename System> 74 1 : class MonotonisedCentralPrim : public Reconstructor<System> { 75 : public: 76 0 : static constexpr size_t dim = 3; 77 : 78 0 : struct AtmosphereTreatment { 79 0 : using type = ::VariableFixing::FixReconstructedStateToAtmosphere; 80 0 : static constexpr Options::String help = { 81 : "What reconstructed states to fix to their atmosphere values."}; 82 : }; 83 0 : struct ReconstructRhoTimesTemperature { 84 0 : using type = bool; 85 0 : static constexpr Options::String help = { 86 : "If 'true' then we reconstruct the rho*T, if 'false' we reconstruct " 87 : "T."}; 88 : }; 89 : 90 0 : using options = 91 : tmpl::list<AtmosphereTreatment, ReconstructRhoTimesTemperature>; 92 0 : static constexpr Options::String help{ 93 : "Monotonised central reconstruction scheme using primitive variables and " 94 : "the metric variables."}; 95 : 96 0 : MonotonisedCentralPrim(::VariableFixing::FixReconstructedStateToAtmosphere 97 : fix_reconstructed_state_to_atmosphere, 98 : bool reconstruct_rho_times_temperature); 99 : 100 0 : MonotonisedCentralPrim() = default; 101 0 : MonotonisedCentralPrim(MonotonisedCentralPrim&&) = default; 102 0 : MonotonisedCentralPrim& operator=(MonotonisedCentralPrim&&) = default; 103 0 : MonotonisedCentralPrim(const MonotonisedCentralPrim&) = default; 104 0 : MonotonisedCentralPrim& operator=(const MonotonisedCentralPrim&) = default; 105 0 : ~MonotonisedCentralPrim() override = default; 106 : 107 0 : explicit MonotonisedCentralPrim(CkMigrateMessage* msg); 108 : 109 0 : WRAPPED_PUPable_decl_base_template(Reconstructor<System>, 110 : MonotonisedCentralPrim); 111 : 112 0 : auto get_clone() const -> std::unique_ptr<Reconstructor<System>> override; 113 : 114 0 : static constexpr bool use_adaptive_order = false; 115 : 116 0 : void pup(PUP::er& p) override; 117 : 118 0 : size_t ghost_zone_size() const override { return 2; } 119 : 120 0 : using reconstruction_argument_tags = 121 : tmpl::list<::Tags::Variables<hydro::grmhd_tags<DataVector>>, 122 : typename System::variables_tag, 123 : hydro::Tags::GrmhdEquationOfState, domain::Tags::Element<dim>, 124 : evolution::dg::subcell::Tags::GhostDataForReconstruction<dim>, 125 : evolution::dg::subcell::Tags::Mesh<dim>, 126 : ::Tags::VariableFixer<VariableFixing::FixToAtmosphere<dim>>>; 127 : 128 : template <size_t ThermodynamicDim, typename TagsList> 129 0 : void reconstruct( 130 : gsl::not_null<std::array<Variables<TagsList>, dim>*> vars_on_lower_face, 131 : gsl::not_null<std::array<Variables<TagsList>, dim>*> vars_on_upper_face, 132 : const Variables<hydro::grmhd_tags<DataVector>>& volume_prims, 133 : const Variables<typename System::variables_tag::type::tags_list>& 134 : volume_spacetime_and_cons_vars, 135 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>& eos, 136 : const Element<dim>& element, 137 : const DirectionalIdMap<dim, evolution::dg::subcell::GhostData>& 138 : ghost_data, 139 : const Mesh<dim>& subcell_mesh, 140 : const VariableFixing::FixToAtmosphere<dim>& fix_to_atmosphere) const; 141 : 142 : /// Called by an element doing DG when the neighbor is doing subcell. 143 : template <size_t ThermodynamicDim, typename TagsList> 144 1 : void reconstruct_fd_neighbor( 145 : gsl::not_null<Variables<TagsList>*> vars_on_face, 146 : const Variables<hydro::grmhd_tags<DataVector>>& subcell_volume_prims, 147 : const Variables< 148 : grmhd::GhValenciaDivClean::Tags::spacetime_reconstruction_tags>& 149 : subcell_volume_spacetime_metric, 150 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>& eos, 151 : const Element<dim>& element, 152 : const DirectionalIdMap<dim, evolution::dg::subcell::GhostData>& 153 : ghost_data, 154 : const Mesh<dim>& subcell_mesh, 155 : const VariableFixing::FixToAtmosphere<dim>& fix_to_atmosphere, 156 : Direction<dim> direction_to_reconstruct) const; 157 : 158 0 : bool reconstruct_rho_times_temperature() const override; 159 : 160 : private: 161 : template <typename LocalSystem> 162 0 : friend bool operator==(const MonotonisedCentralPrim<LocalSystem>& lhs, 163 : const MonotonisedCentralPrim<LocalSystem>& rhs); 164 : template <typename LocalSystem> 165 0 : friend bool operator!=(const MonotonisedCentralPrim<LocalSystem>& lhs, 166 : const MonotonisedCentralPrim<LocalSystem>& rhs); 167 : 168 : ::VariableFixing::FixReconstructedStateToAtmosphere 169 0 : fix_reconstructed_state_to_atmosphere_{ 170 : ::VariableFixing::FixReconstructedStateToAtmosphere::Never}; 171 0 : bool reconstruct_rho_times_temperature_{false}; 172 : }; 173 : } // namespace grmhd::GhValenciaDivClean::fd