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 : #include <string>
8 : #include <variant>
9 :
10 : #include "DataStructures/DataBox/DataBox.hpp"
11 : #include "DataStructures/DataVector.hpp"
12 : #include "DataStructures/Tensor/EagerMath/DeterminantAndInverse.hpp"
13 : #include "DataStructures/Tensor/Tensor.hpp"
14 : #include "Domain/Structure/ElementId.hpp"
15 : #include "Evolution/DgSubcell/ActiveGrid.hpp"
16 : #include "Evolution/DgSubcell/Tags/ActiveGrid.hpp"
17 : #include "Evolution/DgSubcell/Tags/Coordinates.hpp"
18 : #include "Evolution/DgSubcell/Tags/Jacobians.hpp"
19 : #include "Evolution/DgSubcell/Tags/Mesh.hpp"
20 : #include "Evolution/Systems/GeneralizedHarmonic/Actions/SetInitialData.hpp"
21 : #include "Evolution/Systems/GeneralizedHarmonic/Tags.hpp"
22 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Actions/NumericInitialData.hpp"
23 : #include "IO/Importers/Actions/ReadVolumeData.hpp"
24 : #include "IO/Importers/ElementDataReader.hpp"
25 : #include "IO/Importers/Tags.hpp"
26 : #include "NumericalAlgorithms/Spectral/Mesh.hpp"
27 : #include "Options/String.hpp"
28 : #include "Parallel/AlgorithmExecution.hpp"
29 : #include "Parallel/GlobalCache.hpp"
30 : #include "Parallel/Invoke.hpp"
31 : #include "PointwiseFunctions/GeneralRelativity/SpatialMetric.hpp"
32 : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
33 : #include "PointwiseFunctions/Hydro/EquationsOfState/EquationOfState.hpp"
34 : #include "PointwiseFunctions/Hydro/Tags.hpp"
35 : #include "PointwiseFunctions/InitialDataUtilities/InitialData.hpp"
36 : #include "PointwiseFunctions/InitialDataUtilities/Tags/InitialData.hpp"
37 : #include "Utilities/ErrorHandling/Error.hpp"
38 : #include "Utilities/Gsl.hpp"
39 : #include "Utilities/TMPL.hpp"
40 : #include "Utilities/TaggedTuple.hpp"
41 :
42 : /// \cond
43 : namespace Tags {
44 : struct Time;
45 : } // namespace Tags
46 : /// \endcond
47 :
48 1 : namespace grmhd::GhValenciaDivClean {
49 :
50 : /*!
51 : * \brief Numeric initial data loaded from volume data files
52 : *
53 : * This class can be factory-created in the input file to start an evolution
54 : * from numeric initial data. It selects the set of GH variables to load from
55 : * the volume data file (ADM or GH variables), the set of hydro variables, and
56 : * allows to set constant values for some of the hydro variables. This class
57 : * mostly combines the `gh::NumericInitialData` and
58 : * `grmhd::ValenciaDivClean::NumericInitialData` classes.
59 : */
60 1 : class NumericInitialData : public evolution::initial_data::InitialData {
61 : private:
62 0 : using GhNumericId = gh::NumericInitialData;
63 0 : using HydroNumericId = grmhd::ValenciaDivClean::NumericInitialData;
64 :
65 : public:
66 0 : using all_vars =
67 : tmpl::append<GhNumericId::all_vars, HydroNumericId::all_vars>;
68 :
69 0 : struct GhVariables : GhNumericId::Variables {};
70 0 : struct HydroVariables : HydroNumericId::Variables {};
71 :
72 0 : using options =
73 : tmpl::list<importers::OptionTags::FileGlob,
74 : importers::OptionTags::Subgroup,
75 : importers::OptionTags::ObservationValue,
76 : importers::OptionTags::EnableInterpolation, GhVariables,
77 : HydroVariables, HydroNumericId::DensityCutoff>;
78 :
79 0 : static constexpr Options::String help =
80 : "Numeric initial data loaded from volume data files";
81 :
82 0 : NumericInitialData() = default;
83 0 : NumericInitialData(const NumericInitialData& rhs) = default;
84 0 : NumericInitialData& operator=(const NumericInitialData& rhs) = default;
85 0 : NumericInitialData(NumericInitialData&& /*rhs*/) = default;
86 0 : NumericInitialData& operator=(NumericInitialData&& /*rhs*/) = default;
87 0 : ~NumericInitialData() = default;
88 :
89 : /// \cond
90 : explicit NumericInitialData(CkMigrateMessage* msg);
91 : using PUP::able::register_constructor;
92 : WRAPPED_PUPable_decl_template(NumericInitialData);
93 : /// \endcond
94 :
95 0 : std::unique_ptr<evolution::initial_data::InitialData> get_clone()
96 : const override {
97 : return std::make_unique<NumericInitialData>(*this);
98 : }
99 :
100 0 : NumericInitialData(
101 : std::string file_glob, std::string subfile_name,
102 : std::variant<double, importers::ObservationSelector> observation_value,
103 : bool enable_interpolation,
104 : typename GhNumericId::Variables::type gh_selected_variables,
105 : typename HydroNumericId::Variables::type hydro_selected_variables,
106 : double density_cutoff);
107 :
108 0 : const importers::ImporterOptions& importer_options() const {
109 : return gh_numeric_id_.importer_options();
110 : }
111 :
112 0 : const GhNumericId& gh_numeric_id() const { return gh_numeric_id_; }
113 :
114 0 : const HydroNumericId& hydro_numeric_id() const { return hydro_numeric_id_; }
115 :
116 0 : size_t volume_data_id() const;
117 :
118 : template <typename... AllTags>
119 0 : void select_for_import(
120 : const gsl::not_null<tuples::TaggedTuple<AllTags...>*> fields) const {
121 : gh_numeric_id_.select_for_import(fields);
122 : hydro_numeric_id_.select_for_import(fields);
123 : }
124 :
125 : template <typename... AllTags, size_t ThermodynamicDim>
126 0 : void set_initial_data(
127 : const gsl::not_null<tnsr::aa<DataVector, 3>*> spacetime_metric,
128 : const gsl::not_null<tnsr::aa<DataVector, 3>*> pi,
129 : const gsl::not_null<tnsr::iaa<DataVector, 3>*> phi,
130 : const gsl::not_null<Scalar<DataVector>*> rest_mass_density,
131 : const gsl::not_null<Scalar<DataVector>*> electron_fraction,
132 : const gsl::not_null<Scalar<DataVector>*> specific_internal_energy,
133 : const gsl::not_null<tnsr::I<DataVector, 3>*> spatial_velocity,
134 : const gsl::not_null<tnsr::I<DataVector, 3>*> magnetic_field,
135 : const gsl::not_null<Scalar<DataVector>*> div_cleaning_field,
136 : const gsl::not_null<Scalar<DataVector>*> lorentz_factor,
137 : const gsl::not_null<Scalar<DataVector>*> pressure,
138 : const gsl::not_null<Scalar<DataVector>*> specific_enthalpy,
139 : const gsl::not_null<Scalar<DataVector>*> temperature,
140 : const gsl::not_null<tuples::TaggedTuple<AllTags...>*> numeric_data,
141 : const Mesh<3>& mesh,
142 : const InverseJacobian<DataVector, 3, Frame::ElementLogical,
143 : Frame::Inertial>& inv_jacobian,
144 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>&
145 : equation_of_state) const {
146 : gh_numeric_id_.set_initial_data(spacetime_metric, pi, phi, numeric_data,
147 : mesh, inv_jacobian);
148 : const auto spatial_metric = gr::spatial_metric(*spacetime_metric);
149 : const auto inv_spatial_metric =
150 : determinant_and_inverse(spatial_metric).second;
151 : hydro_numeric_id_.set_initial_data(
152 : rest_mass_density, electron_fraction, specific_internal_energy,
153 : spatial_velocity, magnetic_field, div_cleaning_field, lorentz_factor,
154 : pressure, specific_enthalpy, temperature, numeric_data,
155 : inv_spatial_metric, equation_of_state);
156 : }
157 :
158 0 : void pup(PUP::er& p) override;
159 :
160 0 : friend bool operator==(const NumericInitialData& lhs,
161 : const NumericInitialData& rhs);
162 :
163 : private:
164 0 : GhNumericId gh_numeric_id_{};
165 0 : HydroNumericId hydro_numeric_id_{};
166 : };
167 :
168 0 : namespace Actions {
169 :
170 : /*!
171 : * \brief Dispatch loading numeric initial data from files.
172 : *
173 : * Place this action before
174 : * grmhd::GhValenciaDivClean::Actions::SetNumericInitialData in the action list.
175 : * See importers::Actions::ReadAllVolumeDataAndDistribute for details, which is
176 : * invoked by this action.
177 : */
178 1 : struct SetInitialData {
179 0 : using const_global_cache_tags =
180 : tmpl::list<evolution::initial_data::Tags::InitialData>;
181 :
182 : template <typename DbTagsList, typename... InboxTags, typename Metavariables,
183 : typename ArrayIndex, typename ActionList,
184 : typename ParallelComponent>
185 0 : static Parallel::iterable_action_return_t apply(
186 : db::DataBox<DbTagsList>& box,
187 : const tuples::TaggedTuple<InboxTags...>& /*inboxes*/,
188 : Parallel::GlobalCache<Metavariables>& cache,
189 : const ArrayIndex& /*array_index*/, const ActionList /*meta*/,
190 : const ParallelComponent* const parallel_component) {
191 : // Dispatch to the correct `apply` overload based on type of initial data
192 : using initial_data_classes =
193 : tmpl::at<typename Metavariables::factory_creation::factory_classes,
194 : evolution::initial_data::InitialData>;
195 : return call_with_dynamic_type<Parallel::iterable_action_return_t,
196 : initial_data_classes>(
197 : &db::get<evolution::initial_data::Tags::InitialData>(box),
198 : [&box, &cache, ¶llel_component](const auto* const initial_data) {
199 : return apply(make_not_null(&box), *initial_data, cache,
200 : parallel_component);
201 : });
202 : }
203 :
204 : private:
205 0 : static constexpr size_t Dim = 3;
206 :
207 : // Numeric initial data
208 : template <typename DbTagsList, typename Metavariables,
209 : typename ParallelComponent>
210 0 : static Parallel::iterable_action_return_t apply(
211 : const gsl::not_null<db::DataBox<DbTagsList>*> /*box*/,
212 : const NumericInitialData& initial_data,
213 : Parallel::GlobalCache<Metavariables>& cache,
214 : const ParallelComponent* const /*meta*/) {
215 : // Select the subset of the available variables that we want to read from
216 : // the volume data file
217 : tuples::tagged_tuple_from_typelist<db::wrap_tags_in<
218 : importers::Tags::Selected, NumericInitialData::all_vars>>
219 : selected_fields{};
220 : initial_data.select_for_import(make_not_null(&selected_fields));
221 : // Dispatch loading the variables from the volume data file
222 : // - Not using `ckLocalBranch` here to make sure the simple action
223 : // invocation is asynchronous.
224 : auto& reader_component = Parallel::get_parallel_component<
225 : importers::ElementDataReader<Metavariables>>(cache);
226 : Parallel::simple_action<importers::Actions::ReadAllVolumeDataAndDistribute<
227 : 3, NumericInitialData::all_vars, ParallelComponent>>(
228 : reader_component, initial_data.importer_options(),
229 : initial_data.volume_data_id(), std::move(selected_fields));
230 : return {Parallel::AlgorithmExecution::Continue, std::nullopt};
231 : }
232 :
233 : // "AnalyticData"-type initial data
234 : template <typename DbTagsList, typename InitialData, typename Metavariables,
235 : typename ParallelComponent>
236 0 : static Parallel::iterable_action_return_t apply(
237 : const gsl::not_null<db::DataBox<DbTagsList>*> box,
238 : const InitialData& initial_data,
239 : Parallel::GlobalCache<Metavariables>& /*cache*/,
240 : const ParallelComponent* const /*meta*/) {
241 : // Get ADM + hydro variables from analytic data / solution
242 : const auto& [coords, mesh, inv_jacobian] = [&box]() {
243 : if constexpr (db::tag_is_retrievable_v<
244 : evolution::dg::subcell::Tags::ActiveGrid,
245 : db::DataBox<DbTagsList>>) {
246 : const bool on_subcell =
247 : db::get<evolution::dg::subcell::Tags::ActiveGrid>(*box) ==
248 : evolution::dg::subcell::ActiveGrid::Subcell;
249 : if (on_subcell) {
250 : return std::forward_as_tuple(
251 : db::get<evolution::dg::subcell::Tags::Coordinates<
252 : Dim, Frame::Inertial>>(*box),
253 : db::get<evolution::dg::subcell::Tags::Mesh<Dim>>(*box),
254 : db::get<evolution::dg::subcell::fd::Tags::
255 : InverseJacobianLogicalToInertial<Dim>>(*box));
256 : }
257 : }
258 : return std::forward_as_tuple(
259 : db::get<domain::Tags::Coordinates<Dim, Frame::Inertial>>(*box),
260 : db::get<domain::Tags::Mesh<Dim>>(*box),
261 : db::get<domain::Tags::InverseJacobian<Dim, Frame::ElementLogical,
262 : Frame::Inertial>>(*box));
263 : }();
264 : auto vars = evolution::Initialization::initial_data(
265 : initial_data, coords, db::get<::Tags::Time>(*box),
266 : tmpl::append<tmpl::list<gr::Tags::SpatialMetric<DataVector, 3>,
267 : gr::Tags::Lapse<DataVector>,
268 : gr::Tags::Shift<DataVector, 3>,
269 : gr::Tags::ExtrinsicCurvature<DataVector, 3>>,
270 : hydro::grmhd_tags<DataVector>>{});
271 : const auto& spatial_metric =
272 : get<gr::Tags::SpatialMetric<DataVector, 3>>(vars);
273 : const auto& lapse = get<gr::Tags::Lapse<DataVector>>(vars);
274 : const auto& shift = get<gr::Tags::Shift<DataVector, 3>>(vars);
275 : const auto& extrinsic_curvature =
276 : get<gr::Tags::ExtrinsicCurvature<DataVector, 3>>(vars);
277 :
278 : // Compute GH vars from ADM vars
279 : db::mutate<gr::Tags::SpacetimeMetric<DataVector, 3>,
280 : gh::Tags::Pi<DataVector, 3>, gh::Tags::Phi<DataVector, 3>>(
281 : &gh::initial_gh_variables_from_adm<3>, box, spatial_metric, lapse,
282 : shift, extrinsic_curvature, mesh, inv_jacobian);
283 :
284 : // Move hydro vars directly into the DataBox
285 : tmpl::for_each<hydro::grmhd_tags<DataVector>>(
286 : [&box, &vars](const auto tag_v) {
287 : using tag = tmpl::type_from<std::decay_t<decltype(tag_v)>>;
288 : Initialization::mutate_assign<tmpl::list<tag>>(
289 : box, std::move(get<tag>(vars)));
290 : });
291 :
292 : // No need to import numeric initial data, so we terminate the phase by
293 : // pausing the algorithm on this element
294 : return {Parallel::AlgorithmExecution::Pause, std::nullopt};
295 : }
296 : };
297 :
298 : /*!
299 : * \brief Receive numeric initial data loaded by
300 : * grmhd::GhValenciaDivClean::Actions::SetInitialData.
301 : *
302 : * Place this action in the action list after
303 : * grmhd::GhValenciaDivClean::Actions::SetInitialData to wait until the
304 : * data for this element has arrived, and then compute the GH variables and the
305 : * remaining primitive variables and store them in the DataBox to be used as
306 : * initial data.
307 : *
308 : * This action modifies the GH system tags (spacetime metric, pi, phi) and the
309 : * tags listed in `hydro::grmhd_tags` in the DataBox (i.e., the hydro
310 : * primitives). It does not modify conservative variables, so it relies on a
311 : * primitive-to-conservative update in the action list before the evolution can
312 : * start.
313 : *
314 : * \requires This action requires an equation of state, which is retrieved from
315 : * the DataBox as `hydro::Tags::EquationOfStateBase`.
316 : */
317 1 : struct ReceiveNumericInitialData {
318 0 : static constexpr size_t Dim = 3;
319 0 : using inbox_tags =
320 : tmpl::list<importers::Tags::VolumeData<NumericInitialData::all_vars>>;
321 :
322 : template <typename DbTagsList, typename... InboxTags, typename Metavariables,
323 : typename ActionList, typename ParallelComponent>
324 0 : static Parallel::iterable_action_return_t apply(
325 : db::DataBox<DbTagsList>& box, tuples::TaggedTuple<InboxTags...>& inboxes,
326 : const Parallel::GlobalCache<Metavariables>& /*cache*/,
327 : const ElementId<Dim>& /*element_id*/, const ActionList /*meta*/,
328 : const ParallelComponent* const /*meta*/) {
329 : auto& inbox =
330 : tuples::get<importers::Tags::VolumeData<NumericInitialData::all_vars>>(
331 : inboxes);
332 : const auto& initial_data = dynamic_cast<const NumericInitialData&>(
333 : db::get<evolution::initial_data::Tags::InitialData>(box));
334 : const size_t volume_data_id = initial_data.volume_data_id();
335 : if (inbox.find(volume_data_id) == inbox.end()) {
336 : return {Parallel::AlgorithmExecution::Retry, std::nullopt};
337 : }
338 : auto numeric_data = std::move(inbox.extract(volume_data_id).mapped());
339 :
340 : const auto& [mesh, inv_jacobian] = [&box]() {
341 : if constexpr (db::tag_is_retrievable_v<
342 : evolution::dg::subcell::Tags::ActiveGrid,
343 : db::DataBox<DbTagsList>>) {
344 : const bool on_subcell =
345 : db::get<evolution::dg::subcell::Tags::ActiveGrid>(box) ==
346 : evolution::dg::subcell::ActiveGrid::Subcell;
347 : if (on_subcell) {
348 : return std::forward_as_tuple(
349 : db::get<evolution::dg::subcell::Tags::Mesh<Dim>>(box),
350 : db::get<evolution::dg::subcell::fd::Tags::
351 : InverseJacobianLogicalToInertial<Dim>>(box));
352 : }
353 : }
354 : return std::forward_as_tuple(
355 : db::get<domain::Tags::Mesh<Dim>>(box),
356 : db::get<domain::Tags::InverseJacobian<Dim, Frame::ElementLogical,
357 : Frame::Inertial>>(box));
358 : }();
359 : const auto& equation_of_state =
360 : db::get<hydro::Tags::EquationOfStateBase>(box);
361 :
362 : db::mutate<gr::Tags::SpacetimeMetric<DataVector, 3>,
363 : gh::Tags::Pi<DataVector, 3>, gh::Tags::Phi<DataVector, 3>,
364 : hydro::Tags::RestMassDensity<DataVector>,
365 : hydro::Tags::ElectronFraction<DataVector>,
366 : hydro::Tags::SpecificInternalEnergy<DataVector>,
367 : hydro::Tags::SpatialVelocity<DataVector, 3>,
368 : hydro::Tags::MagneticField<DataVector, 3>,
369 : hydro::Tags::DivergenceCleaningField<DataVector>,
370 : hydro::Tags::LorentzFactor<DataVector>,
371 : hydro::Tags::Pressure<DataVector>,
372 : hydro::Tags::SpecificEnthalpy<DataVector>,
373 : hydro::Tags::Temperature<DataVector>>(
374 : [&initial_data, &numeric_data, &equation_of_state](
375 : const gsl::not_null<tnsr::aa<DataVector, 3>*> spacetime_metric,
376 : const gsl::not_null<tnsr::aa<DataVector, 3>*> pi,
377 : const gsl::not_null<tnsr::iaa<DataVector, 3>*> phi,
378 : const gsl::not_null<Scalar<DataVector>*> rest_mass_density,
379 : const gsl::not_null<Scalar<DataVector>*> electron_fraction,
380 : const gsl::not_null<Scalar<DataVector>*> specific_internal_energy,
381 : const gsl::not_null<tnsr::I<DataVector, 3>*> spatial_velocity,
382 : const gsl::not_null<tnsr::I<DataVector, 3>*> magnetic_field,
383 : const gsl::not_null<Scalar<DataVector>*> div_cleaning_field,
384 : const gsl::not_null<Scalar<DataVector>*> lorentz_factor,
385 : const gsl::not_null<Scalar<DataVector>*> pressure,
386 : const gsl::not_null<Scalar<DataVector>*> specific_enthalpy,
387 : const gsl::not_null<Scalar<DataVector>*> temperature,
388 : const auto& local_mesh, const auto& local_inv_jacobian) {
389 : initial_data.set_initial_data(
390 : spacetime_metric, pi, phi, rest_mass_density, electron_fraction,
391 : specific_internal_energy, spatial_velocity, magnetic_field,
392 : div_cleaning_field, lorentz_factor, pressure, specific_enthalpy,
393 : temperature, make_not_null(&numeric_data), local_mesh,
394 : local_inv_jacobian, equation_of_state);
395 : },
396 : make_not_null(&box), mesh, inv_jacobian);
397 :
398 : return {Parallel::AlgorithmExecution::Continue, std::nullopt};
399 : }
400 : };
401 :
402 : } // namespace Actions
403 :
404 : } // namespace grmhd::GhValenciaDivClean
|