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 <limits>
8 : #include <utility>
9 :
10 : #include "DataStructures/DataBox/PrefixHelpers.hpp"
11 : #include "DataStructures/DataBox/Prefixes.hpp"
12 : #include "DataStructures/DataVector.hpp"
13 : #include "DataStructures/TaggedContainers.hpp"
14 : #include "DataStructures/TaggedTuple.hpp"
15 : #include "DataStructures/Tensor/Tensor.hpp"
16 : #include "Evolution/DiscontinuousGalerkin/TimeDerivativeDecisions.hpp"
17 : #include "Evolution/PassVariables.hpp"
18 : #include "Evolution/Systems/GeneralizedHarmonic/GaugeSourceFunctions/Harmonic.hpp"
19 : #include "Evolution/Systems/GeneralizedHarmonic/System.hpp"
20 : #include "Evolution/Systems/GeneralizedHarmonic/TimeDerivative.hpp"
21 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/AllSolutions.hpp"
22 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/StressEnergy.hpp"
23 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/System.hpp"
24 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Tags.hpp"
25 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/System.hpp"
26 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/TimeDerivativeTerms.hpp"
27 : #include "Evolution/VariableFixing/FixToAtmosphere.hpp"
28 : #include "Evolution/VariableFixing/Tags.hpp"
29 : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp" // For reference tag
30 : #include "PointwiseFunctions/Hydro/Tags.hpp"
31 : #include "Time/Tags/Time.hpp"
32 : #include "Utilities/Gsl.hpp"
33 : #include "Utilities/Literals.hpp"
34 : #include "Utilities/TMPL.hpp"
35 :
36 : namespace grmhd::GhValenciaDivClean {
37 : namespace detail {
38 : // Some temporary tags appear both in the GRMHD temporary list and the
39 : // GeneralizedHarmonic temporary list, so we wrap the GRMHD temporaries in this
40 : // prefix tag to avoid collisions in data structures used to store the
41 : // temporaries.
42 : template <typename Tag>
43 : struct ValenciaTempTag : db::SimpleTag, db::PrefixTag {
44 : using tag = Tag;
45 : using type = typename Tag::type;
46 : };
47 :
48 : namespace dt_type_aliases {
49 : using gh_dt_tags =
50 : db::wrap_tags_in<::Tags::dt,
51 : typename gh::System<3_st>::variables_tag::tags_list>;
52 : using valencia_dt_tags = db::wrap_tags_in<
53 : ::Tags::dt,
54 : typename grmhd::ValenciaDivClean::System::variables_tag::tags_list>;
55 :
56 : using dt_tags = tmpl::append<gh_dt_tags, valencia_dt_tags>;
57 : } // namespace dt_type_aliases
58 :
59 : struct TimeDerivativeTermsImpl;
60 : struct TimeDerivativeTermsImpl {
61 : // Tags related to GeneralizedHarmonic system
62 : using gh_arg_tags = typename gh::TimeDerivative<
63 : ghmhd::GhValenciaDivClean::InitialData::analytic_solutions_and_data_list,
64 : 3_st>::argument_tags;
65 :
66 : using gh_temp_tags = typename gh::TimeDerivative<
67 : ghmhd::GhValenciaDivClean::InitialData::analytic_solutions_and_data_list,
68 : 3_st>::temporary_tags;
69 :
70 : using valencia_flux_tags = db::wrap_tags_in<
71 : ::Tags::Flux,
72 : typename grmhd::ValenciaDivClean::System::variables_tag::tags_list,
73 : tmpl::size_t<3>, Frame::Inertial>;
74 :
75 : using valencia_temp_tags =
76 : typename grmhd::ValenciaDivClean::TimeDerivativeTerms::temporary_tags;
77 :
78 : using valencia_time_derivative_arg_tags =
79 : typename grmhd::ValenciaDivClean::TimeDerivativeTerms::argument_tags;
80 :
81 : using trace_reversed_stress_result_tags =
82 : tmpl::list<Tags::TraceReversedStressEnergy, Tags::FourVelocityOneForm,
83 : Tags::ComovingMagneticFieldOneForm>;
84 : using trace_reversed_stress_argument_tags = tmpl::list<
85 : hydro::Tags::RestMassDensity<DataVector>,
86 : hydro::Tags::SpatialVelocityOneForm<DataVector, 3_st, Frame::Inertial>,
87 : hydro::Tags::MagneticFieldOneForm<DataVector, 3_st, Frame::Inertial>,
88 : hydro::Tags::MagneticFieldSquared<DataVector>,
89 : hydro::Tags::MagneticFieldDotSpatialVelocity<DataVector>,
90 : hydro::Tags::LorentzFactor<DataVector>,
91 : grmhd::ValenciaDivClean::TimeDerivativeTerms::OneOverLorentzFactorSquared,
92 : hydro::Tags::Pressure<DataVector>,
93 : hydro::Tags::SpecificInternalEnergy<DataVector>,
94 : gr::Tags::SpacetimeMetric<DataVector, 3>,
95 : gr::Tags::Shift<DataVector, 3_st>, gr::Tags::Lapse<DataVector>>;
96 :
97 : template <typename TemporaryTagsList, typename... ExtraTags>
98 : static evolution::dg::TimeDerivativeDecisions<3> apply(
99 : const gsl::not_null<Variables<dt_type_aliases::dt_tags>*> dt_vars_ptr,
100 : const gsl::not_null<Variables<db::wrap_tags_in<
101 : ::Tags::Flux, typename ValenciaDivClean::System::flux_variables,
102 : tmpl::size_t<3>, Frame::Inertial>>*>
103 : fluxes_ptr,
104 : const gsl::not_null<Variables<TemporaryTagsList>*> temps_ptr,
105 :
106 : const tnsr::iaa<DataVector, 3>& d_spacetime_metric,
107 : const tnsr::iaa<DataVector, 3>& d_pi,
108 : const tnsr::ijaa<DataVector, 3>& d_phi,
109 :
110 : const tuples::TaggedTuple<ExtraTags...>& arguments) {
111 : generalized_harmonic_time_derivative(
112 : dt_vars_ptr, temps_ptr, d_spacetime_metric, d_pi, d_phi, arguments,
113 : dt_type_aliases::gh_dt_tags{}, gh_temp_tags{}, gh_arg_tags{});
114 :
115 : // Track whether or not we extracted the 3+1 quantities
116 : bool three_plus_one_extracted = false;
117 :
118 : // If we are in the atmosphere, then we can skip the evolution
119 : // of the GRMHD system completely. This inevitably depends on parameters
120 : // of the FixToAtmosphere code, so we grab those directly.
121 : if (const auto& fix_to_atmosphere = get<Tags::detail::TemporaryReference<
122 : ::Tags::VariableFixer<::VariableFixing::FixToAtmosphere<3>>>>(
123 : arguments);
124 : max(get(get<Tags::detail::TemporaryReference<
125 : hydro::Tags::RestMassDensity<DataVector>>>(arguments))) <=
126 : std::min({fix_to_atmosphere.density_of_atmosphere(),
127 : (fix_to_atmosphere.velocity_limiting().has_value()
128 : ? fix_to_atmosphere.velocity_limiting()
129 : ->atmosphere_density_cutoff
130 : : std::numeric_limits<double>::infinity()),
131 : (fix_to_atmosphere.kappa_limiting().has_value()
132 : ? fix_to_atmosphere.kappa_limiting()->density_lower_bound
133 : : std::numeric_limits<double>::infinity())}) *
134 : (1.0 + 10.0 * std::numeric_limits<double>::epsilon())) {
135 : // Point into the right memory, then set it to zero.
136 :
137 : ASSERT(
138 : max(get(get<tmpl::front<dt_type_aliases::valencia_dt_tags>>(
139 : *dt_vars_ptr))) == 0.0,
140 : "GH+GRMHD assumes the time derivatives are set to zero in general."
141 : " If this is no longer the case, please set them to zero in "
142 : "atmosphere by changing the code where this ASSERT was triggered.");
143 : // Code that we could use to set the sources to zero if needed.
144 : // Variables<tmpl::list<ValenciaDtTags...>> dt_div_clean(
145 : // get<tmpl::front<tmpl::list<ValenciaDtTags...>>>(*dt_vars_ptr)[0]
146 : // .data(),
147 : // 0.0);
148 : fluxes_ptr->initialize(fluxes_ptr->number_of_grid_points(), 0.0);
149 : return evolution::dg::TimeDerivativeDecisions<3>{false};
150 : } else {
151 : // MHD calls
152 :
153 : // extracting 3+1 quantities to assign spacetime spatial derivatives
154 : extract_three_plus_one_quantities(temps_ptr, arguments);
155 : three_plus_one_extracted = true;
156 :
157 : // MHD computation
158 : aggregate_time_derivative_terms(
159 : dt_vars_ptr, fluxes_ptr, temps_ptr, arguments,
160 : dt_type_aliases::valencia_dt_tags{}, valencia_flux_tags{},
161 : valencia_temp_tags{}, valencia_time_derivative_arg_tags{},
162 : trace_reversed_stress_result_tags{},
163 : trace_reversed_stress_argument_tags{});
164 : }
165 :
166 : if (!three_plus_one_extracted) {
167 : // If in atmosphere, extract 3+1 quantities for neutrino
168 : // evolution
169 : extract_three_plus_one_quantities(temps_ptr, arguments);
170 : three_plus_one_extracted = true;
171 : // Neutrino evolution will be called below
172 : }
173 :
174 : return evolution::dg::TimeDerivativeDecisions<3>{true};
175 : }
176 :
177 : template <typename OutputTags, typename TemporaryTagsList,
178 : typename... ExtraTags, typename... GhDtTags, typename... GhTempTags,
179 : typename... GhArgTags>
180 : static void generalized_harmonic_time_derivative(
181 : const gsl::not_null<Variables<OutputTags>*> dt_vars_ptr,
182 : const gsl::not_null<Variables<TemporaryTagsList>*> temps_ptr,
183 : const tnsr::iaa<DataVector, 3>& d_spacetime_metric,
184 : const tnsr::iaa<DataVector, 3>& d_pi,
185 : const tnsr::ijaa<DataVector, 3>& d_phi,
186 : const tuples::TaggedTuple<ExtraTags...>& arguments,
187 : tmpl::list<GhDtTags...> /*meta*/, tmpl::list<GhTempTags...> /*meta*/,
188 : tmpl::list<GhArgTags...> /*meta*/) {
189 : gh::TimeDerivative<
190 : ghmhd::GhValenciaDivClean::InitialData::
191 : analytic_solutions_and_data_list,
192 : 3_st>::apply(get<GhDtTags>(dt_vars_ptr)...,
193 : get<GhTempTags>(temps_ptr)..., d_spacetime_metric, d_pi,
194 : d_phi,
195 : get<Tags::detail::TemporaryReference<GhArgTags>>(
196 : arguments)...);
197 : }
198 :
199 : template <typename TemporaryTagsList, typename... ExtraTags>
200 : static void extract_three_plus_one_quantities(
201 : const gsl::not_null<Variables<TemporaryTagsList>*> temps_ptr,
202 : const tuples::TaggedTuple<ExtraTags...>& arguments) {
203 : // Extract the 3+1 quantities from the spacetime metric
204 :
205 : // Check whether or not sqrt(det(g)) exists based on gauge condition
206 : if (get<Tags::detail::TemporaryReference<gh::gauges::Tags::GaugeCondition>>(
207 : arguments)
208 : .is_harmonic()) {
209 : get(get<gr::Tags::SqrtDetSpatialMetric<DataVector>>(*temps_ptr)) =
210 : sqrt(get(get<gr::Tags::DetSpatialMetric<DataVector>>(*temps_ptr)));
211 : }
212 :
213 : // extract spatial derivative of lapse
214 : for (size_t i = 0; i < 3; ++i) {
215 : get<::Tags::deriv<gr::Tags::Lapse<DataVector>, tmpl::size_t<3>,
216 : Frame::Inertial>>(*temps_ptr)
217 : .get(i) = -get(get<gr::Tags::Lapse<DataVector>>(*temps_ptr)) *
218 : get<gh::Tags::HalfPhiTwoNormals<3>>(*temps_ptr).get(i);
219 : }
220 :
221 : const auto& phi =
222 : get<Tags::detail::TemporaryReference<gh::Tags::Phi<DataVector, 3>>>(
223 : arguments);
224 : const auto& inv_spatial_metric =
225 : get<gr::Tags::InverseSpatialMetric<DataVector, 3>>(*temps_ptr);
226 : const auto& shift = get<gr::Tags::Shift<DataVector, 3>>(*temps_ptr);
227 :
228 : // extract spatial derivative of shift
229 : for (size_t i = 0; i < 3; ++i) {
230 : for (size_t j = 0; j < 3; ++j) {
231 : get<::Tags::deriv<gr::Tags::Shift<DataVector, 3>, tmpl::size_t<3>,
232 : Frame::Inertial>>(*temps_ptr)
233 : .get(i, j) = inv_spatial_metric.get(j, 0) * phi.get(i, 0, 1);
234 : for (size_t k = 1; k < 3; ++k) {
235 : get<::Tags::deriv<gr::Tags::Shift<DataVector, 3>, tmpl::size_t<3>,
236 : Frame::Inertial>>(*temps_ptr)
237 : .get(i, j) += inv_spatial_metric.get(j, k) * phi.get(i, 0, k + 1);
238 : }
239 : for (size_t k = 0; k < 3; ++k) {
240 : for (size_t l = 0; l < 3; ++l) {
241 : get<::Tags::deriv<gr::Tags::Shift<DataVector, 3>, tmpl::size_t<3>,
242 : Frame::Inertial>>(*temps_ptr)
243 : .get(i, j) -= shift.get(k) * inv_spatial_metric.get(j, l) *
244 : phi.get(i, l + 1, k + 1);
245 : }
246 : }
247 : }
248 : }
249 :
250 : const auto& pi =
251 : get<Tags::detail::TemporaryReference<gh::Tags::Pi<DataVector, 3>>>(
252 : arguments);
253 : for (size_t i = 0; i < 3; ++i) {
254 : for (size_t j = i; j < 3; ++j) {
255 : get<gr::Tags::ExtrinsicCurvature<DataVector, 3>>(*temps_ptr).get(i, j) =
256 : 0.5 * (pi.get(i + 1, j + 1) +
257 : get<gh::Tags::PhiOneNormal<3>>(*temps_ptr).get(i, j + 1) +
258 : get<gh::Tags::PhiOneNormal<3>>(*temps_ptr).get(j, i + 1));
259 : }
260 : }
261 :
262 : } // end extract3plus1
263 :
264 : template <typename OutputTags, typename TemporaryTagsList,
265 : typename... ExtraTags, typename... ValenciaDtTags,
266 : typename... ValenciaFluxTags, typename... ValenciaTempTags,
267 : typename... ValenciaTimeDerivativeArgTags,
268 : typename... TraceReversedStressResultTags,
269 : typename... TraceReversedStressArgumentTags>
270 : static void aggregate_time_derivative_terms(
271 : const gsl::not_null<Variables<OutputTags>*> dt_vars_ptr,
272 : const gsl::not_null<Variables<db::wrap_tags_in<
273 : ::Tags::Flux, typename ValenciaDivClean::System::flux_variables,
274 : tmpl::size_t<3>, Frame::Inertial>>*>
275 : fluxes_ptr,
276 : const gsl::not_null<Variables<TemporaryTagsList>*> temps_ptr,
277 : const tuples::TaggedTuple<ExtraTags...>& arguments,
278 : tmpl::list<ValenciaDtTags...> /*meta*/,
279 : tmpl::list<ValenciaFluxTags...> /*meta*/,
280 : tmpl::list<ValenciaTempTags...> /*meta*/,
281 : tmpl::list<ValenciaTimeDerivativeArgTags...> /*meta*/,
282 : tmpl::list<TraceReversedStressResultTags...> /*meta*/,
283 : tmpl::list<TraceReversedStressArgumentTags...> /*meta*/) {
284 : using extra_tags_list = tmpl::list<ExtraTags...>;
285 :
286 : grmhd::ValenciaDivClean::TimeDerivativeTerms::apply(
287 : get<ValenciaDtTags>(dt_vars_ptr)...,
288 : get<ValenciaFluxTags>(fluxes_ptr)...,
289 : get<ValenciaTempTags>(temps_ptr)...,
290 :
291 : get<tmpl::conditional_t<
292 : tmpl::list_contains_v<extra_tags_list,
293 : Tags::detail::TemporaryReference<
294 : ValenciaTimeDerivativeArgTags>>,
295 : Tags::detail::TemporaryReference<ValenciaTimeDerivativeArgTags>,
296 : ValenciaTimeDerivativeArgTags>>(arguments, *temps_ptr)...);
297 :
298 : trace_reversed_stress_energy(
299 : get<TraceReversedStressResultTags>(temps_ptr)...,
300 : get<tmpl::conditional_t<
301 : tmpl::list_contains_v<extra_tags_list,
302 : Tags::detail::TemporaryReference<
303 : TraceReversedStressArgumentTags>>,
304 : Tags::detail::TemporaryReference<TraceReversedStressArgumentTags>,
305 : TraceReversedStressArgumentTags>>(*temps_ptr, arguments)...);
306 :
307 : add_stress_energy_term_to_dt_pi(
308 : get<::Tags::dt<gh::Tags::Pi<DataVector, 3>>>(dt_vars_ptr),
309 : get<grmhd::GhValenciaDivClean::Tags::TraceReversedStressEnergy>(
310 : *temps_ptr),
311 : get<gr::Tags::Lapse<DataVector>>(*temps_ptr));
312 : }
313 : }; // namespace detail
314 : } // namespace detail
315 :
316 : /*!
317 : * \brief Compute the RHS terms and flux values for both the Generalized
318 : * Harmonic formulation of Einstein's equations and the Valencia formulation of
319 : * the GRMHD equations with divergence cleaning.
320 : *
321 : * \details The bulk of the computations in this class dispatch to
322 : * `gh::TimeDerivative` and
323 : * `grmhd::ValenciaDivClean::TimeDerivativeTerms` as a 'product system' -- each
324 : * independently operating on its own subset of the supplied variable
325 : * collections.
326 : * The additional step is taken to compute the trace-reversed stress energy
327 : * tensor associated with the GRMHD part of the system and add its contribution
328 : * to the \f$\partial_t \Pi_{a b}\f$ variable in the Generalized Harmonic
329 : * system, which is the only explicit coupling required to back-react the effect
330 : * of matter on the spacetime solution.
331 : *
332 : * \note The MHD calculation reuses any spacetime quantities in its
333 : * argument_tags that are computed by the GH time derivative. However, other
334 : * quantities that aren't computed by the GH time derivative like the extrinsic
335 : * curvature are currently still retrieved from the DataBox. Those calculations
336 : * can be explicitly inlined here to reduce memory pressure and the number of
337 : * compute tags.
338 : */
339 1 : struct TimeDerivativeTerms : evolution::PassVariables {
340 0 : using gh_dt_tags =
341 : db::wrap_tags_in<::Tags::dt,
342 : typename gh::System<3_st>::variables_tag::tags_list>;
343 0 : using valencia_dt_tags = db::wrap_tags_in<
344 : ::Tags::dt,
345 : typename grmhd::ValenciaDivClean::System::variables_tag::tags_list>;
346 :
347 0 : using dt_tags = tmpl::append<gh_dt_tags, valencia_dt_tags>;
348 :
349 0 : using d_spatial_metric = ::Tags::deriv<gr::Tags::SpatialMetric<DataVector, 3>,
350 : tmpl::size_t<3>, Frame::Inertial>;
351 :
352 0 : using gh_temp_tags = typename gh::TimeDerivative<
353 : ghmhd::GhValenciaDivClean::InitialData::analytic_solutions_and_data_list,
354 : 3_st>::temporary_tags;
355 0 : using gh_gradient_tags = typename gh::System<3_st>::gradients_tags;
356 0 : using gh_arg_tags = typename gh::TimeDerivative<
357 : ghmhd::GhValenciaDivClean::InitialData::analytic_solutions_and_data_list,
358 : 3_st>::argument_tags;
359 :
360 0 : using valencia_temp_tags =
361 : typename grmhd::ValenciaDivClean::TimeDerivativeTerms::temporary_tags;
362 : // Additional temp tags are the derivatives of the metric since GH doesn't
363 : // explicitly calculate those.
364 0 : using valencia_extra_temp_tags =
365 : tmpl::list<::Tags::deriv<gr::Tags::Lapse<DataVector>, tmpl::size_t<3>,
366 : Frame::Inertial>,
367 : ::Tags::deriv<gr::Tags::Shift<DataVector, 3>, tmpl::size_t<3>,
368 : Frame::Inertial>,
369 : gr::Tags::ExtrinsicCurvature<DataVector, 3>>;
370 0 : using valencia_arg_tags = tmpl::list_difference<
371 : typename grmhd::ValenciaDivClean::TimeDerivativeTerms::argument_tags,
372 : tmpl::append<gh_temp_tags, valencia_extra_temp_tags>>;
373 :
374 0 : using trace_reversed_stress_result_tags =
375 : tmpl::list<Tags::TraceReversedStressEnergy, Tags::FourVelocityOneForm,
376 : Tags::ComovingMagneticFieldOneForm>;
377 0 : using extra_temp_tags = tmpl::list<gr::Tags::SpatialMetric<DataVector, 3>>;
378 :
379 0 : using temporary_tags = tmpl::remove<
380 : tmpl::remove_duplicates<tmpl::append<
381 : gh_temp_tags, valencia_temp_tags, valencia_extra_temp_tags,
382 : trace_reversed_stress_result_tags, extra_temp_tags>>,
383 : gr::Tags::SpatialMetric<DataVector, 3>>;
384 0 : using argument_tags = tmpl::remove<
385 : tmpl::remove<tmpl::append<gh_arg_tags,
386 :
387 : valencia_arg_tags,
388 :
389 : tmpl::list<::Tags::VariableFixer<
390 : ::VariableFixing::FixToAtmosphere<3>>>>,
391 : gr::Tags::SpatialMetric<DataVector, 3>>,
392 : d_spatial_metric>;
393 :
394 : template <typename... Args>
395 0 : static evolution::dg::TimeDerivativeDecisions<3> apply(
396 : const gsl::not_null<Variables<dt_tags>*> dt_vars_ptr,
397 : const gsl::not_null<Variables<db::wrap_tags_in<
398 : ::Tags::Flux, typename ValenciaDivClean::System::flux_variables,
399 : tmpl::size_t<3>, Frame::Inertial>>*>
400 : fluxes_ptr,
401 : const gsl::not_null<Variables<temporary_tags>*> temps_ptr,
402 : const tnsr::iaa<DataVector, 3>& d_spacetime_metric,
403 : const tnsr::iaa<DataVector, 3>& d_pi,
404 : const tnsr::ijaa<DataVector, 3>& d_phi, const Args&... args) {
405 : using args_list = tmpl::push_back<
406 : db::wrap_tags_in<Tags::detail::TemporaryReference, argument_tags>,
407 : gr::Tags::SpatialMetric<DataVector, 3>, d_spatial_metric>;
408 : tuples::tagged_tuple_from_typelist<args_list> arguments{
409 : args..., typename gr::Tags::SpatialMetric<DataVector, 3>::type{},
410 : typename d_spatial_metric::type{}};
411 : const size_t number_of_points = get<Tags::detail::TemporaryReference<
412 : gr::Tags::SpacetimeMetric<DataVector, 3>>>(arguments)[0]
413 : .size();
414 : for (size_t i = 0; i < 3; ++i) {
415 : for (size_t j = i; j < 3; ++j) {
416 : make_const_view(
417 : make_not_null(
418 : &std::as_const(
419 : get<gr::Tags::SpatialMetric<DataVector, 3>>(arguments))
420 : .get(i, j)),
421 : get<Tags::detail::TemporaryReference<
422 : gr::Tags::SpacetimeMetric<DataVector, 3>>>(arguments)
423 : .get(i + 1, j + 1),
424 : 0, number_of_points);
425 : }
426 : }
427 : for (size_t i = 0; i < 3; ++i) {
428 : for (size_t j = 0; j < 3; ++j) {
429 : for (size_t k = j; k < 3; ++k) {
430 : make_const_view(
431 : make_not_null(&std::as_const(get<d_spatial_metric>(arguments))
432 : .get(i, j, k)),
433 : get<Tags::detail::TemporaryReference<
434 : gh::Tags::Phi<DataVector, 3>>>(arguments)
435 : .get(i, j + 1, k + 1),
436 : 0, number_of_points);
437 : }
438 : }
439 : }
440 :
441 : return detail::TimeDerivativeTermsImpl::apply(dt_vars_ptr, fluxes_ptr,
442 : temps_ptr, d_spacetime_metric,
443 : d_pi, d_phi, arguments);
444 : }
445 : };
446 : } // namespace grmhd::GhValenciaDivClean
|