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 <limits>
9 : #include <memory>
10 : #include <utility>
11 :
12 : #include "DataStructures/Variables.hpp"
13 : #include "DataStructures/VariablesTag.hpp"
14 : #include "Domain/Structure/DirectionalIdMap.hpp"
15 : #include "Domain/Tags.hpp"
16 : #include "Evolution/DgSubcell/Tags/GhostDataForReconstruction.hpp"
17 : #include "Evolution/DgSubcell/Tags/Mesh.hpp"
18 : #include "Evolution/Systems/GeneralizedHarmonic/Tags.hpp"
19 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Reconstructor.hpp"
20 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/System.hpp"
21 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Tags.hpp"
22 : #include "Evolution/VariableFixing/FixToAtmosphere.hpp"
23 : #include "Evolution/VariableFixing/Tags.hpp"
24 : #include "NumericalAlgorithms/FiniteDifference/FallbackReconstructorType.hpp"
25 : #include "Options/Auto.hpp"
26 : #include "Options/Context.hpp"
27 : #include "Options/String.hpp"
28 : #include "PointwiseFunctions/Hydro/Tags.hpp"
29 : #include "Utilities/TMPL.hpp"
30 :
31 : /// \cond
32 : class DataVector;
33 : template <size_t Dim>
34 : class Direction;
35 : template <size_t Dim>
36 : class Element;
37 : template <size_t Dim>
38 : class ElementId;
39 : namespace EquationsOfState {
40 : template <bool IsRelativistic, size_t ThermodynamicDim>
41 : class EquationOfState;
42 : } // namespace EquationsOfState
43 : template <size_t Dim>
44 : class Mesh;
45 : namespace gsl {
46 : template <typename T>
47 : class not_null;
48 : } // namespace gsl
49 : namespace PUP {
50 : class er;
51 : } // namespace PUP
52 : template <typename TagsList>
53 : class Variables;
54 : namespace evolution::dg::subcell {
55 : class GhostData;
56 : } // namespace evolution::dg::subcell
57 : /// \endcond
58 :
59 : namespace grmhd::GhValenciaDivClean::fd {
60 : /*!
61 : * \brief Fifth order weighted nonlinear compact scheme reconstruction using the
62 : * Z oscillation indicator. See ::fd::reconstruction::wcns5z() for details.
63 : *
64 : */
65 : // template on System instead
66 : template <typename System>
67 1 : class Wcns5zPrim : public Reconstructor<System> {
68 : private:
69 0 : using prims_to_reconstruct_tags =
70 : tmpl::list<hydro::Tags::RestMassDensity<DataVector>,
71 : hydro::Tags::ElectronFraction<DataVector>,
72 : hydro::Tags::Temperature<DataVector>,
73 : hydro::Tags::LorentzFactorTimesSpatialVelocity<DataVector, 3>,
74 : hydro::Tags::MagneticField<DataVector, 3>,
75 : hydro::Tags::DivergenceCleaningField<DataVector>>;
76 :
77 0 : using FallbackReconstructorType =
78 : ::fd::reconstruction::FallbackReconstructorType;
79 :
80 : public:
81 0 : static constexpr size_t dim = 3;
82 :
83 0 : struct NonlinearWeightExponent {
84 0 : using type = size_t;
85 0 : static constexpr Options::String help = {
86 : "The exponent q to which the oscillation indicator term is raised"};
87 : };
88 0 : struct Epsilon {
89 0 : using type = double;
90 0 : static constexpr Options::String help = {
91 : "The parameter added to the oscillation indicators to avoid division "
92 : "by zero"};
93 : };
94 0 : struct FallbackReconstructor {
95 0 : using type = FallbackReconstructorType;
96 0 : static constexpr Options::String help = {
97 : "A reconstruction scheme to fallback to adaptively. Finite difference "
98 : "will switch to this reconstruction scheme if there are more extrema "
99 : "in a FD stencil than a specified number. See also the option "
100 : "'MaxNumberOfExtrema' below. Adaptive fallback is disabled if 'None'."};
101 : };
102 0 : struct MaxNumberOfExtrema {
103 0 : using type = size_t;
104 0 : static constexpr Options::String help = {
105 : "The maximum allowed number of extrema in FD stencil for using Wcns5z "
106 : "reconstruction before switching to a low-order reconstruction. If "
107 : "FallbackReconstructor=None, this option is ignored"};
108 : };
109 :
110 0 : struct AtmosphereTreatment {
111 0 : using type = ::VariableFixing::FixReconstructedStateToAtmosphere;
112 0 : static constexpr Options::String help = {
113 : "What reconstructed states to fix to their atmosphere values."};
114 : };
115 :
116 0 : struct ReconstructRhoTimesTemperature {
117 0 : using type = bool;
118 0 : static constexpr Options::String help = {
119 : "If 'true' then we reconstruct the rho*T, if 'false' we reconstruct "
120 : "T."};
121 : };
122 :
123 0 : using options =
124 : tmpl::list<NonlinearWeightExponent, Epsilon, FallbackReconstructor,
125 : MaxNumberOfExtrema, AtmosphereTreatment,
126 : ReconstructRhoTimesTemperature>;
127 :
128 0 : static constexpr Options::String help{
129 : "WCNS 5Z reconstruction scheme using primitive variables."};
130 :
131 0 : Wcns5zPrim() = default;
132 0 : Wcns5zPrim(Wcns5zPrim&&) = default;
133 0 : Wcns5zPrim& operator=(Wcns5zPrim&&) = default;
134 0 : Wcns5zPrim(const Wcns5zPrim&) = default;
135 0 : Wcns5zPrim& operator=(const Wcns5zPrim&) = default;
136 0 : ~Wcns5zPrim() override = default;
137 :
138 0 : Wcns5zPrim(size_t nonlinear_weight_exponent, double epsilon,
139 : FallbackReconstructorType fallback_reconstructor,
140 : size_t max_number_of_extrema,
141 : ::VariableFixing::FixReconstructedStateToAtmosphere
142 : fix_reconstructed_state_to_atmosphere,
143 : bool reconstruct_rho_times_temperature);
144 :
145 0 : explicit Wcns5zPrim(CkMigrateMessage* msg);
146 :
147 0 : WRAPPED_PUPable_decl_base_template(Reconstructor<System>, Wcns5zPrim<System>);
148 :
149 0 : auto get_clone() const -> std::unique_ptr<Reconstructor<System>> override;
150 :
151 0 : static constexpr bool use_adaptive_order = false;
152 :
153 0 : void pup(PUP::er& p) override;
154 :
155 0 : size_t ghost_zone_size() const override { return 3; }
156 :
157 0 : using reconstruction_argument_tags =
158 : tmpl::list<::Tags::Variables<hydro::grmhd_tags<DataVector>>,
159 : typename System::variables_tag,
160 : hydro::Tags::GrmhdEquationOfState, domain::Tags::Element<dim>,
161 : evolution::dg::subcell::Tags::GhostDataForReconstruction<dim>,
162 : evolution::dg::subcell::Tags::Mesh<dim>,
163 : ::Tags::VariableFixer<VariableFixing::FixToAtmosphere<dim>>>;
164 :
165 : template <size_t ThermodynamicDim, typename TagsList>
166 0 : void reconstruct(
167 : gsl::not_null<std::array<Variables<TagsList>, dim>*> vars_on_lower_face,
168 : gsl::not_null<std::array<Variables<TagsList>, dim>*> vars_on_upper_face,
169 : const Variables<hydro::grmhd_tags<DataVector>>& volume_prims,
170 : const Variables<typename System::variables_tag::type::tags_list>&
171 : volume_spacetime_and_cons_vars,
172 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>& eos,
173 : const Element<dim>& element,
174 : const DirectionalIdMap<dim, evolution::dg::subcell::GhostData>&
175 : ghost_data,
176 : const Mesh<dim>& subcell_mesh,
177 : const VariableFixing::FixToAtmosphere<dim>& fix_to_atmosphere) const;
178 :
179 : template <size_t ThermodynamicDim, typename TagsList>
180 0 : void reconstruct_fd_neighbor(
181 : gsl::not_null<Variables<TagsList>*> vars_on_face,
182 : const Variables<hydro::grmhd_tags<DataVector>>& subcell_volume_prims,
183 : const Variables<
184 : grmhd::GhValenciaDivClean::Tags::spacetime_reconstruction_tags>&
185 : subcell_volume_spacetime_metric,
186 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>& eos,
187 : const Element<dim>& element,
188 : const DirectionalIdMap<dim, evolution::dg::subcell::GhostData>&
189 : ghost_data,
190 : const Mesh<dim>& subcell_mesh,
191 : const VariableFixing::FixToAtmosphere<dim>& fix_to_atmosphere,
192 : const Direction<dim>& direction_to_reconstruct) const;
193 :
194 0 : bool reconstruct_rho_times_temperature() const override;
195 :
196 : private:
197 : template <typename LocalSystem>
198 : // NOLINTNEXTLINE(readability-redundant-declaration)
199 0 : friend bool operator==(const Wcns5zPrim<LocalSystem>& lhs,
200 : const Wcns5zPrim<LocalSystem>& rhs);
201 :
202 : template <typename LocalSystem>
203 0 : friend bool operator!=(const Wcns5zPrim<LocalSystem>& lhs,
204 : const Wcns5zPrim<LocalSystem>& rhs);
205 :
206 0 : size_t nonlinear_weight_exponent_ = 0;
207 0 : double epsilon_ = std::numeric_limits<double>::signaling_NaN();
208 0 : FallbackReconstructorType fallback_reconstructor_ =
209 : FallbackReconstructorType::None;
210 0 : size_t max_number_of_extrema_ = 0;
211 : ::VariableFixing::FixReconstructedStateToAtmosphere
212 0 : fix_reconstructed_state_to_atmosphere_{
213 : ::VariableFixing::FixReconstructedStateToAtmosphere::Never};
214 0 : bool reconstruct_rho_times_temperature_{false};
215 :
216 0 : void (*reconstruct_)(gsl::not_null<std::array<gsl::span<double>, dim>*>,
217 : gsl::not_null<std::array<gsl::span<double>, dim>*>,
218 : const gsl::span<const double>&,
219 : const DirectionMap<dim, gsl::span<const double>>&,
220 : const Index<dim>&, size_t, double, size_t) = nullptr;
221 0 : void (*reconstruct_lower_neighbor_)(gsl::not_null<DataVector*>,
222 : const DataVector&, const DataVector&,
223 : const Index<dim>&, const Index<dim>&,
224 : const Direction<dim>&, const double&,
225 : const size_t&) = nullptr;
226 0 : void (*reconstruct_upper_neighbor_)(gsl::not_null<DataVector*>,
227 : const DataVector&, const DataVector&,
228 : const Index<dim>&, const Index<dim>&,
229 : const Direction<dim>&, const double&,
230 : const size_t&) = nullptr;
231 : };
232 : } // namespace grmhd::GhValenciaDivClean::fd
|