SpECTRE Documentation Coverage Report
Current view: top level - Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell - TimeDerivative.hpp Hit Total Coverage
Commit: a18e59fda1a195609825c55450f7d61ad20a91a4 Lines: 1 3 33.3 %
Date: 2026-06-11 22:10:41
Legend: Lines: hit not hit

          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 <optional>
       9             : #include <type_traits>
      10             : 
      11             : #include "DataStructures/DataBox/AsAccess.hpp"
      12             : #include "DataStructures/DataBox/DataBox.hpp"
      13             : #include "DataStructures/DataBox/PrefixHelpers.hpp"
      14             : #include "DataStructures/DataBox/Prefixes.hpp"
      15             : #include "DataStructures/DataVector.hpp"
      16             : #include "DataStructures/TaggedContainers.hpp"
      17             : #include "DataStructures/Tensor/Tensor.hpp"
      18             : #include "DataStructures/Variables.hpp"
      19             : #include "DataStructures/VectorImpl.hpp"
      20             : #include "Domain/Structure/Element.hpp"
      21             : #include "Domain/Tags.hpp"
      22             : #include "Domain/TagsTimeDependent.hpp"
      23             : #include "Evolution/BoundaryCorrection.hpp"
      24             : #include "Evolution/BoundaryCorrectionTags.hpp"
      25             : #include "Evolution/DgSubcell/CartesianFluxDivergence.hpp"
      26             : #include "Evolution/DgSubcell/ComputeBoundaryTerms.hpp"
      27             : #include "Evolution/DgSubcell/CorrectPackagedData.hpp"
      28             : #include "Evolution/DgSubcell/Projection.hpp"
      29             : #include "Evolution/DgSubcell/ReconstructionOrder.hpp"
      30             : #include "Evolution/DgSubcell/Tags/Coordinates.hpp"
      31             : #include "Evolution/DgSubcell/Tags/GhostDataForReconstruction.hpp"
      32             : #include "Evolution/DgSubcell/Tags/Jacobians.hpp"
      33             : #include "Evolution/DgSubcell/Tags/Mesh.hpp"
      34             : #include "Evolution/DgSubcell/Tags/OnSubcellFaces.hpp"
      35             : #include "Evolution/DiscontinuousGalerkin/Actions/NormalCovectorAndMagnitude.hpp"
      36             : #include "Evolution/DiscontinuousGalerkin/Actions/PackageDataImpl.hpp"
      37             : #include "Evolution/DiscontinuousGalerkin/MortarTags.hpp"
      38             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/AllSolutions.hpp"
      39             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/BoundaryConditionGhostData.hpp"
      40             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Derivatives.hpp"
      41             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/FilterOptions.hpp"
      42             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Filters.hpp"
      43             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Reconstructor.hpp"
      44             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Tag.hpp"
      45             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/StressEnergy.hpp"
      46             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/System.hpp"
      47             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Tags.hpp"
      48             : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/TimeDerivativeTerms.hpp"
      49             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Fluxes.hpp"
      50             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Sources.hpp"
      51             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Subcell/ComputeFluxes.hpp"
      52             : #include "Evolution/Systems/GrMhd/ValenciaDivClean/TimeDerivativeTerms.hpp"
      53             : #include "NumericalAlgorithms/FiniteDifference/PartialDerivatives.hpp"
      54             : #include "NumericalAlgorithms/Spectral/Mesh.hpp"
      55             : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/DerivSpatialMetric.hpp"
      56             : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/ExtrinsicCurvature.hpp"
      57             : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/SpatialDerivOfLapse.hpp"
      58             : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/SpatialDerivOfShift.hpp"
      59             : #include "PointwiseFunctions/Hydro/Tags.hpp"
      60             : #include "Utilities/CallWithDynamicType.hpp"
      61             : #include "Utilities/ErrorHandling/Assert.hpp"
      62             : #include "Utilities/Gsl.hpp"
      63             : #include "Utilities/TMPL.hpp"
      64             : 
      65             : /// \cond
      66             : namespace Tags {
      67             : struct Time;
      68             : }  // namespace Tags
      69             : /// \endcond
      70             : 
      71             : namespace grmhd::GhValenciaDivClean::subcell {
      72             : namespace detail {
      73             : template <class GhDtTagsList, class GhTemporariesList, class GhGradientTagsList,
      74             :           class GhExtraTagsList, class GrmhdDtTagsList,
      75             :           class GrmhdSourceTagsList, class GrmhdArgumentSourceTagsList,
      76             :           typename System>
      77             : struct ComputeTimeDerivImpl;
      78             : 
      79             : template <class... GhDtTags, class... GhTemporaries, class... GhGradientTags,
      80             :           class... GhExtraTags, class... GrmhdDtTags, class... GrmhdSourceTags,
      81             :           class... GrmhdArgumentSourceTags, typename System>
      82             : struct ComputeTimeDerivImpl<
      83             :     tmpl::list<GhDtTags...>, tmpl::list<GhTemporaries...>,
      84             :     tmpl::list<GhGradientTags...>, tmpl::list<GhExtraTags...>,
      85             :     tmpl::list<GrmhdDtTags...>, tmpl::list<GrmhdSourceTags...>,
      86             :     tmpl::list<GrmhdArgumentSourceTags...>, System> {
      87             :   template <class DbTagsList>
      88             :   static void apply(
      89             :       const gsl::not_null<db::DataBox<DbTagsList>*> box,
      90             :       const tnsr::I<DataVector, 3, Frame::Inertial>& inertial_coords,
      91             :       const Scalar<DataVector>& cell_centered_det_inv_jacobian,
      92             :       const InverseJacobian<DataVector, 3, Frame::ElementLogical,
      93             :                             Frame::Inertial>&
      94             :           cell_centered_logical_to_inertial_inv_jacobian,
      95             :       const std::array<double, 3>& one_over_delta_xi,
      96             :       const std::array<Variables<tmpl::list<GrmhdDtTags...>>, 3>&
      97             :           boundary_corrections,
      98             :       const Variables<
      99             :           db::wrap_tags_in<::Tags::deriv, typename System::gradients_tags,
     100             :                            tmpl::size_t<3>, Frame::Inertial>>& gh_derivs) {
     101             :     const Mesh<3>& subcell_mesh =
     102             :         db::get<evolution::dg::subcell::Tags::Mesh<3>>(*box);
     103             :     const size_t number_of_points = subcell_mesh.number_of_grid_points();
     104             :     // Note: GH+GRMHD tags are always GH,GRMHD
     105             :     using deriv_lapse = ::Tags::deriv<gr::Tags::Lapse<DataVector>,
     106             :                                       tmpl::size_t<3>, Frame::Inertial>;
     107             :     using deriv_shift = ::Tags::deriv<gr::Tags::Shift<DataVector, 3>,
     108             :                                       tmpl::size_t<3>, Frame::Inertial>;
     109             :     using deriv_spatial_metric =
     110             :         ::Tags::deriv<gr::Tags::SpatialMetric<DataVector, 3>, tmpl::size_t<3>,
     111             :                       Frame::Inertial>;
     112             :     using extra_tags_for_grmhd =
     113             :         tmpl::list<deriv_lapse, deriv_shift, deriv_spatial_metric,
     114             :                    gr::Tags::ExtrinsicCurvature<DataVector, 3>>;
     115             :     using temporary_tags = tmpl::remove_duplicates<tmpl::append<
     116             :         typename gh::TimeDerivative<ghmhd::GhValenciaDivClean::InitialData::
     117             :                                         analytic_solutions_and_data_list,
     118             :                                     3_st>::temporary_tags,
     119             :         tmpl::push_front<typename grmhd::ValenciaDivClean::TimeDerivativeTerms::
     120             :                              temporary_tags,
     121             :                          ::gh::Tags::ConstraintGamma0>,
     122             :         extra_tags_for_grmhd,
     123             :         tmpl::list<
     124             :             Tags::TraceReversedStressEnergy, Tags::FourVelocityOneForm,
     125             :             grmhd::ValenciaDivClean::Tags::ComovingMagneticFieldOneForm>>>;
     126             :     Variables<temporary_tags> temp_tags{subcell_mesh.number_of_grid_points()};
     127             :     const auto temp_tags_ptr = make_not_null(&temp_tags);
     128             : 
     129             :     // Compute constraint damping terms.
     130             :     const double time = db::get<::Tags::Time>(*box);
     131             :     const auto& functions_of_time =
     132             :         db::get<::domain::Tags::FunctionsOfTime>(*box);
     133             :     const auto& grid_coords =
     134             :         db::get<evolution::dg::subcell::Tags::Coordinates<3, Frame::Grid>>(
     135             :             *box);
     136             :     db::get<gh::Tags::DampingFunctionGamma0<3, Frame::Grid>> (*box)(
     137             :         get<gh::Tags::ConstraintGamma0>(temp_tags_ptr), grid_coords, time,
     138             :         functions_of_time);
     139             :     db::get<gh::Tags::DampingFunctionGamma1<3, Frame::Grid>> (*box)(
     140             :         get<gh::Tags::ConstraintGamma1>(temp_tags_ptr), grid_coords, time,
     141             :         functions_of_time);
     142             :     db::get<gh::Tags::DampingFunctionGamma2<3, Frame::Grid>> (*box)(
     143             :         get<gh::Tags::ConstraintGamma2>(temp_tags_ptr), grid_coords, time,
     144             :         functions_of_time);
     145             : 
     146             :     using variables_tag = typename System::variables_tag;
     147             :     using dt_variables_tag = db::add_tag_prefix<::Tags::dt, variables_tag>;
     148             :     const gsl::not_null<typename dt_variables_tag::type*> dt_vars_ptr =
     149             :         db::mutate<dt_variables_tag>(
     150             :             [](const auto local_dt_vars_ptr) { return local_dt_vars_ptr; },
     151             :             box);
     152             :     dt_vars_ptr->initialize(subcell_mesh.number_of_grid_points());
     153             : 
     154             :     using primitives_tag = typename System::primitive_variables_tag;
     155             :     using evolved_vars_tag = typename System::variables_tag;
     156             : 
     157             :     const auto& primitive_vars = db::get<primitives_tag>(*box);
     158             :     const auto& evolved_vars = db::get<evolved_vars_tag>(*box);
     159             : 
     160             :     // Velocity of the moving mesh, if applicable. We project the value
     161             :     // stored on the DG grid onto the subcell grid.
     162             :     const Mesh<3>& dg_mesh = db::get<domain::Tags::Mesh<3>>(*box);
     163             :     const std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>&
     164             :         mesh_velocity_dg = db::get<domain::Tags::MeshVelocity<3>>(*box);
     165             :     const std::optional<Scalar<DataVector>>& div_mesh_velocity_dg =
     166             :         db::get<domain::Tags::DivMeshVelocity>(*box);
     167             :     std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>
     168             :         mesh_velocity_subcell = {};
     169             :     if (mesh_velocity_dg.has_value()) {
     170             :       mesh_velocity_subcell = tnsr::I<DataVector, 3, Frame::Inertial>{
     171             :           subcell_mesh.number_of_grid_points()};
     172             :       for (size_t i = 0; i < 3; i++) {
     173             :         mesh_velocity_subcell.value().get(i) =
     174             :             evolution::dg::subcell::fd::project(mesh_velocity_dg.value().get(i),
     175             :                                                 dg_mesh,
     176             :                                                 subcell_mesh.extents());
     177             :       }
     178             :     }
     179             : 
     180             :     gh::TimeDerivative<
     181             :         ghmhd::GhValenciaDivClean::InitialData::
     182             :             analytic_solutions_and_data_list,
     183             :         3_st>::apply(get<::Tags::dt<GhDtTags>>(dt_vars_ptr)...,
     184             :                      get<GhTemporaries>(temp_tags_ptr)...,
     185             :                      get<::Tags::deriv<GhGradientTags, tmpl::size_t<3>,
     186             :                                        Frame::Inertial>>(gh_derivs)...,
     187             :                      get<GhExtraTags>(evolved_vars, temp_tags)...,
     188             : 
     189             :                      db::get<::gh::gauges::Tags::GaugeCondition>(*box),
     190             :                      db::get<evolution::dg::subcell::Tags::Mesh<3>>(*box), time,
     191             :                      inertial_coords,
     192             :                      cell_centered_logical_to_inertial_inv_jacobian,
     193             :                      mesh_velocity_subcell);
     194             :     if (get<gh::gauges::Tags::GaugeCondition>(*box).is_harmonic()) {
     195             :       get(get<gr::Tags::SqrtDetSpatialMetric<DataVector>>(*temp_tags_ptr)) =
     196             :           sqrt(
     197             :               get(get<gr::Tags::DetSpatialMetric<DataVector>>(*temp_tags_ptr)));
     198             :     }
     199             : 
     200             :     // Add source terms from moving mesh
     201             :     if (mesh_velocity_dg.has_value()) {
     202             :       tmpl::for_each<tmpl::list<GhDtTags...>>([&dt_vars_ptr,
     203             :                                                &mesh_velocity_subcell,
     204             :                                                &gh_derivs](
     205             :                                                   auto evolved_var_tag_v) {
     206             :         using evolved_var_tag = tmpl::type_from<decltype(evolved_var_tag_v)>;
     207             :         using dt_tag = ::Tags::dt<evolved_var_tag>;
     208             :         using grad_tag =
     209             :             ::Tags::deriv<evolved_var_tag, tmpl::size_t<3>, Frame::Inertial>;
     210             :         // Flux and gradients use the same indexing conventions,
     211             :         // replacing the direction of the face with the direction
     212             :         // of the derivative.
     213             :         using FluxTensor = typename grad_tag::type;
     214             :         auto& dt_var = get<dt_tag>(*dt_vars_ptr);
     215             :         const auto& grad_var = get<grad_tag>(gh_derivs);
     216             :         for (size_t i = 0; i < dt_var.size(); ++i) {
     217             :           const auto tensor_index = dt_var.get_tensor_index(i);
     218             :           for (size_t j = 0; j < 3; j++) {
     219             :             const auto grad_index =
     220             :                 FluxTensor::get_storage_index(prepend(tensor_index, j));
     221             :             // Add (mesh_velocity)^j grad_j (var[i])
     222             :             dt_var[i] +=
     223             :                 mesh_velocity_subcell.value().get(j) * grad_var[grad_index];
     224             :           }
     225             :         }
     226             :       });
     227             :     }
     228             : 
     229             :     {
     230             :       // Set extra tags needed for GRMHD source terms. We compute these from
     231             :       // quantities already computed inside the GH RHS computation to minimize
     232             :       // FLOPs.
     233             :       const auto& lapse = get<gr::Tags::Lapse<DataVector>>(temp_tags);
     234             :       const auto& half_phi_two_normals =
     235             :           get<gh::Tags::HalfPhiTwoNormals<3>>(temp_tags);
     236             :       const auto& phi = get<gh::Tags::Phi<DataVector, 3>>(evolved_vars);
     237             :       const auto& phi_one_normal = get<gh::Tags::PhiOneNormal<3>>(temp_tags);
     238             :       const auto& spacetime_normal_vector =
     239             :           get<gr::Tags::SpacetimeNormalVector<DataVector, 3>>(temp_tags);
     240             :       const auto& inverse_spacetime_metric =
     241             :           get<gr::Tags::InverseSpacetimeMetric<DataVector, 3>>(temp_tags);
     242             : 
     243             :       auto& spatial_deriv_lapse = get<deriv_lapse>(temp_tags);
     244             :       auto& spatial_deriv_shift = get<deriv_shift>(temp_tags);
     245             :       // Compute d_i beta^i
     246             :       for (size_t i = 0; i < 3; ++i) {
     247             :         // Use spatial_deriv_lapse as temp buffer to reduce number of 2*
     248             :         // operations.
     249             :         const auto& phi_two_normals_i = spatial_deriv_lapse.get(i) =
     250             :             2.0 * half_phi_two_normals.get(i);
     251             :         for (size_t j = 0; j < 3; ++j) {
     252             :           spatial_deriv_shift.get(i, j) =
     253             :               spacetime_normal_vector.get(j + 1) * phi_two_normals_i;
     254             :           for (size_t a = 0; a < 4; ++a) {
     255             :             spatial_deriv_shift.get(i, j) +=
     256             :                 inverse_spacetime_metric.get(j + 1, a) *
     257             :                 phi_one_normal.get(i, a);
     258             :           }
     259             :           spatial_deriv_shift.get(i, j) *= get(lapse);
     260             :         }
     261             :       }
     262             : 
     263             :       // Compute d_i lapse
     264             :       for (size_t i = 0; i < 3; ++i) {
     265             :         spatial_deriv_lapse.get(i) = -get(lapse) * half_phi_two_normals.get(i);
     266             :       }
     267             :       // Extract d_i \gamma_{ij}
     268             :       for (size_t k = 0; k < 3; ++k) {
     269             :         for (size_t i = 0; i < 3; ++i) {
     270             :           for (size_t j = i; j < 3; ++j) {
     271             :             get<deriv_spatial_metric>(temp_tags).get(k, i, j) =
     272             :                 phi.get(k, i + 1, j + 1);
     273             :           }
     274             :         }
     275             :       }
     276             : 
     277             :       // Compute extrinsic curvature
     278             :       const auto& pi = get<gh::Tags::Pi<DataVector, 3>>(evolved_vars);
     279             :       for (size_t i = 0; i < 3; ++i) {
     280             :         for (size_t j = i; j < 3; ++j) {
     281             :           get<gr::Tags::ExtrinsicCurvature<DataVector, 3>>(temp_tags).get(i,
     282             :                                                                           j) =
     283             :               0.5 * (pi.get(i + 1, j + 1) + phi_one_normal.get(i, j + 1) +
     284             :                      phi_one_normal.get(j, i + 1));
     285             :         }
     286             :       }
     287             :     }  // End scope for computing metric terms in GRMHD source terms.
     288             : 
     289             :     grmhd::ValenciaDivClean::ComputeSources::apply(
     290             :         get<::Tags::dt<GrmhdSourceTags>>(dt_vars_ptr)...,
     291             :         get<GrmhdArgumentSourceTags>(temp_tags, primitive_vars, evolved_vars,
     292             :                                      *box)...);
     293             : 
     294             :     // Zero GRMHD tags that don't have sources.
     295             :     tmpl::for_each<tmpl::list<GrmhdDtTags...>>([&dt_vars_ptr](
     296             :                                                    auto evolved_var_tag_v) {
     297             :       using evolved_var_tag = tmpl::type_from<decltype(evolved_var_tag_v)>;
     298             :       using dt_tag = ::Tags::dt<evolved_var_tag>;
     299             :       auto& dt_var = get<dt_tag>(*dt_vars_ptr);
     300             :       for (size_t i = 0; i < dt_var.size(); ++i) {
     301             :         if constexpr (not tmpl::list_contains_v<tmpl::list<GrmhdSourceTags...>,
     302             :                                                 evolved_var_tag>) {
     303             :           // Zero the GRMHD dt(u) for variables that do not have a source term .
     304             :           // This is necessary to avoid `+=` to a `NaN` (debug mode) or random
     305             :           // garbage (release mode) when adding to dt_var below.
     306             :           dt_var[i] = 0.0;
     307             :         }
     308             :       }
     309             :     });
     310             :     // Correction to source terms due to moving mesh
     311             :     if (div_mesh_velocity_dg.has_value()) {
     312             :       const DataVector div_mesh_velocity_subcell =
     313             :           evolution::dg::subcell::fd::project(
     314             :               div_mesh_velocity_dg.value().get(), dg_mesh,
     315             :               subcell_mesh.extents());
     316             :       tmpl::for_each<tmpl::list<GrmhdDtTags...>>(
     317             :           [&dt_vars_ptr, &div_mesh_velocity_subcell,
     318             :            &evolved_vars](auto evolved_var_tag_v) {
     319             :             using evolved_var_tag =
     320             :                 tmpl::type_from<decltype(evolved_var_tag_v)>;
     321             :             using dt_tag = ::Tags::dt<evolved_var_tag>;
     322             :             auto& dt_var = get<dt_tag>(*dt_vars_ptr);
     323             :             const auto& evolved_var = get<evolved_var_tag>(evolved_vars);
     324             :             for (size_t i = 0; i < dt_var.size(); ++i) {
     325             :               dt_var[i] -= div_mesh_velocity_subcell * evolved_var[i];
     326             :             }
     327             :           });
     328             :     }
     329             : 
     330             :     const tnsr::ii<DataVector, 3> spatial_metric{};
     331             :     for (size_t i = 0; i < 3; ++i) {
     332             :       for (size_t j = i; j < 3; ++j) {
     333             :         make_const_view(
     334             :             make_not_null(&spatial_metric.get(i, j)),
     335             :             get<gr::Tags::SpacetimeMetric<DataVector, 3>>(evolved_vars)
     336             :                 .get(i + 1, j + 1),
     337             :             0, number_of_points);
     338             :       }
     339             :     }
     340             : 
     341             :     tenex::evaluate<ti::i>(get<hydro::Tags::SpatialVelocityOneForm<
     342             :                                DataVector, 3, Frame::Inertial>>(temp_tags_ptr),
     343             :                            get<hydro::Tags::SpatialVelocity<DataVector, 3>>(
     344             :                                primitive_vars)(ti::J) *
     345             :                                spatial_metric(ti::i, ti::j));
     346             : 
     347             :     tenex::evaluate<ti::i>(
     348             :         get<hydro::Tags::MagneticFieldOneForm<DataVector, 3, Frame::Inertial>>(
     349             :             temp_tags_ptr),
     350             :         get<hydro::Tags::MagneticField<DataVector, 3>>(primitive_vars)(ti::J) *
     351             :             spatial_metric(ti::i, ti::j));
     352             : 
     353             :     tenex::evaluate(
     354             :         get<hydro::Tags::MagneticFieldSquared<DataVector>>(temp_tags_ptr),
     355             :         get<hydro::Tags::MagneticField<DataVector, 3>>(primitive_vars)(ti::J) *
     356             :             get<hydro::Tags::MagneticFieldOneForm<DataVector, 3>>(temp_tags)(
     357             :                 ti::j));
     358             : 
     359             :     tenex::evaluate(
     360             :         get<hydro::Tags::MagneticFieldDotSpatialVelocity<DataVector>>(
     361             :             temp_tags_ptr),
     362             :         get<hydro::Tags::SpatialVelocity<DataVector, 3>>(primitive_vars)(
     363             :             ti::J) *
     364             :             get<hydro::Tags::MagneticFieldOneForm<DataVector, 3>>(temp_tags)(
     365             :                 ti::j));
     366             : 
     367             :     tenex::evaluate(get<typename ValenciaDivClean::TimeDerivativeTerms::
     368             :                             OneOverLorentzFactorSquared>(temp_tags_ptr),
     369             :                     1.0 / (square(get<hydro::Tags::LorentzFactor<DataVector>>(
     370             :                               primitive_vars)())));
     371             : 
     372             :     trace_reversed_stress_energy(
     373             :         get<Tags::TraceReversedStressEnergy>(temp_tags_ptr),
     374             :         get<Tags::FourVelocityOneForm>(temp_tags_ptr),
     375             :         get<grmhd::ValenciaDivClean::Tags::ComovingMagneticFieldOneForm>(
     376             :             temp_tags_ptr),
     377             : 
     378             :         get<hydro::Tags::RestMassDensity<DataVector>>(evolved_vars, temp_tags,
     379             :                                                       primitive_vars),
     380             :         get<hydro::Tags::SpatialVelocityOneForm<DataVector, 3,
     381             :                                                 Frame::Inertial>>(
     382             :             evolved_vars, temp_tags, primitive_vars),
     383             : 
     384             :         get<hydro::Tags::MagneticFieldOneForm<DataVector, 3, Frame::Inertial>>(
     385             :             evolved_vars, temp_tags, primitive_vars),
     386             : 
     387             :         get<hydro::Tags::MagneticFieldSquared<DataVector>>(
     388             :             evolved_vars, temp_tags, primitive_vars),
     389             : 
     390             :         get<hydro::Tags::MagneticFieldDotSpatialVelocity<DataVector>>(
     391             :             evolved_vars, temp_tags, primitive_vars),
     392             :         get<hydro::Tags::LorentzFactor<DataVector>>(evolved_vars, temp_tags,
     393             :                                                     primitive_vars),
     394             :         get<typename ValenciaDivClean::TimeDerivativeTerms::
     395             :                 OneOverLorentzFactorSquared>(evolved_vars, temp_tags,
     396             :                                              primitive_vars),
     397             :         get<hydro::Tags::Pressure<DataVector>>(evolved_vars, temp_tags,
     398             :                                                primitive_vars),
     399             :         get<hydro::Tags::SpecificInternalEnergy<DataVector>>(
     400             :             evolved_vars, temp_tags, primitive_vars),
     401             :         get<gr::Tags::SpacetimeMetric<DataVector, 3>>(evolved_vars, temp_tags,
     402             :                                                       primitive_vars),
     403             :         get<gr::Tags::Shift<DataVector, 3>>(evolved_vars, temp_tags,
     404             :                                             primitive_vars),
     405             :         get<gr::Tags::Lapse<DataVector>>(evolved_vars, temp_tags,
     406             :                                          primitive_vars));
     407             : 
     408             :     add_stress_energy_term_to_dt_pi(
     409             :         get<::Tags::dt<gh::Tags::Pi<DataVector, 3>>>(dt_vars_ptr),
     410             :         get<Tags::TraceReversedStressEnergy>(temp_tags),
     411             :         get<gr::Tags::Lapse<DataVector>>(temp_tags));
     412             : 
     413             :     for (size_t dim = 0; dim < 3; ++dim) {
     414             :       const auto& boundary_correction_in_axis =
     415             :           gsl::at(boundary_corrections, dim);
     416             :       const double inverse_delta = gsl::at(one_over_delta_xi, dim);
     417             :       EXPAND_PACK_LEFT_TO_RIGHT([&dt_vars_ptr, &boundary_correction_in_axis,
     418             :                                  &cell_centered_det_inv_jacobian, dim,
     419             :                                  inverse_delta, &subcell_mesh]() {
     420             :         auto& dt_var = *get<::Tags::dt<GrmhdDtTags>>(dt_vars_ptr);
     421             :         const auto& var_correction =
     422             :             get<GrmhdDtTags>(boundary_correction_in_axis);
     423             :         for (size_t i = 0; i < dt_var.size(); ++i) {
     424             :           evolution::dg::subcell::add_cartesian_flux_divergence(
     425             :               make_not_null(&dt_var[i]), inverse_delta,
     426             :               get(cell_centered_det_inv_jacobian), var_correction[i],
     427             :               subcell_mesh.extents(), dim);
     428             :         }
     429             :       }());
     430             :     }
     431             :   }
     432             : };
     433             : }  // namespace detail
     434             : 
     435             : /*!
     436             :  * \brief Compute the time derivative on the subcell grid using FD
     437             :  * reconstruction.
     438             :  *
     439             :  * The code makes the following unchecked assumptions:
     440             :  * - Assumes Cartesian coordinates with a diagonal Jacobian matrix
     441             :  * from the logical to the inertial frame
     442             :  */
     443             : template <typename System>
     444           1 : struct TimeDerivative {
     445             :   template <typename DbTagsList>
     446           0 :   static void apply(const gsl::not_null<db::DataBox<DbTagsList>*> box) {
     447             :     using metavariables =
     448             :         typename std::decay_t<decltype(db::get<Parallel::Tags::Metavariables>(
     449             :             *box))>;
     450             :     using evolved_vars_tag = typename System::variables_tag;
     451             :     using evolved_vars_tags = typename evolved_vars_tag::tags_list;
     452             :     using grmhd_evolved_vars_tag =
     453             :         typename grmhd::ValenciaDivClean::System::variables_tag;
     454             :     using grmhd_evolved_vars_tags = typename grmhd_evolved_vars_tag::tags_list;
     455             :     using fluxes_tags =
     456             :         db::wrap_tags_in<::Tags::Flux, typename System::flux_variables,
     457             :                          tmpl::size_t<3>, Frame::Inertial>;
     458             :     using prim_tag = typename System::primitive_variables_tag;
     459             :     using prim_tags = typename prim_tag::tags_list;
     460             :     using recons_prim_tags = tmpl::push_front<tmpl::push_back<
     461             :         prim_tags,
     462             :         hydro::Tags::LorentzFactorTimesSpatialVelocity<DataVector, 3>>>;
     463             :     using gradients_tags = typename System::gradients_tags;
     464             : 
     465             :     const Mesh<3>& dg_mesh = db::get<domain::Tags::Mesh<3>>(*box);
     466             :     const Mesh<3>& subcell_mesh =
     467             :         db::get<evolution::dg::subcell::Tags::Mesh<3>>(*box);
     468             :     ASSERT(
     469             :         subcell_mesh == Mesh<3>(subcell_mesh.extents(0), subcell_mesh.basis(0),
     470             :                                 subcell_mesh.quadrature(0)),
     471             :         "The subcell/FD mesh must be isotropic for the FD time derivative but "
     472             :         "got "
     473             :             << subcell_mesh);
     474             :     const size_t num_pts = subcell_mesh.number_of_grid_points();
     475             :     const size_t reconstructed_num_pts =
     476             :         (subcell_mesh.extents(0) + 1) *
     477             :         subcell_mesh.extents().slice_away(0).product();
     478             : 
     479             :     const tnsr::I<DataVector, 3, Frame::ElementLogical>&
     480             :         cell_centered_logical_coords =
     481             :             db::get<evolution::dg::subcell::Tags::Coordinates<
     482             :                 3, Frame::ElementLogical>>(*box);
     483             :     std::array<double, 3> one_over_delta_xi{};
     484             :     for (size_t i = 0; i < 3; ++i) {
     485             :       // Note: assumes isotropic extents
     486             :       gsl::at(one_over_delta_xi, i) =
     487             :           1.0 / (get<0>(cell_centered_logical_coords)[1] -
     488             :                  get<0>(cell_centered_logical_coords)[0]);
     489             :     }
     490             :     const auto& cell_centered_logical_to_inertial_inv_jacobian = db::get<
     491             :         evolution::dg::subcell::fd::Tags::InverseJacobianLogicalToInertial<3>>(
     492             :         *box);
     493             :     const auto& inertial_coords =
     494             :         db::get<evolution::dg::subcell::Tags::Coordinates<3, Frame::Inertial>>(
     495             :             *box);
     496             : 
     497             :     const Element<3>& element = db::get<domain::Tags::Element<3>>(*box);
     498             :     const bool element_is_interior = element.external_boundaries().empty();
     499             :     constexpr bool subcell_enabled_at_external_boundary =
     500             :         metavariables::SubcellOptions::subcell_enabled_at_external_boundary;
     501             : 
     502             :     ASSERT(element_is_interior or subcell_enabled_at_external_boundary,
     503             :            "Subcell time derivative is called at a boundary element while "
     504             :            "using subcell is disabled at external boundaries."
     505             :            "ElementID "
     506             :                << element.id());
     507             : 
     508             :     const fd::Reconstructor<System>& recons =
     509             :         db::get<fd::Tags::Reconstructor<System>>(*box);
     510             :     // If the element has external boundaries and subcell is enabled for
     511             :     // boundary elements, compute FD ghost data with a given boundary condition.
     512             :     if constexpr (subcell_enabled_at_external_boundary) {
     513             :       if (not element_is_interior) {
     514             :         fd::BoundaryConditionGhostData<System>::apply(box, element, recons);
     515             :       }
     516             :     }
     517             :     std::optional<std::array<gsl::span<std::uint8_t>, 3>>
     518             :         reconstruction_order{};
     519             : 
     520             :     if (const auto& filter_options =
     521             :             db::get<grmhd::GhValenciaDivClean::fd::Tags::FilterOptions>(*box);
     522             :         filter_options.spacetime_dissipation.has_value()) {
     523             :       db::mutate<evolved_vars_tag>(
     524             :           [&filter_options, &recons, &subcell_mesh](const auto evolved_vars_ptr,
     525             :                                                     const auto& ghost_data) {
     526             :             typename evolved_vars_tag::type filtered_vars = *evolved_vars_ptr;
     527             :             // $(recons.ghost_zone_size() - 1) * 2 + 1$ => always use highest
     528             :             // order dissipation filter possible.
     529             :             grmhd::GhValenciaDivClean::fd::spacetime_kreiss_oliger_filter(
     530             :                 make_not_null(&filtered_vars), *evolved_vars_ptr, ghost_data,
     531             :                 subcell_mesh, 2 * recons.ghost_zone_size(),
     532             :                 filter_options.spacetime_dissipation.value());
     533             :             *evolved_vars_ptr = filtered_vars;
     534             :           },
     535             :           box,
     536             :           db::get<evolution::dg::subcell::Tags::GhostDataForReconstruction<3>>(
     537             :               *box));
     538             :     }
     539             : 
     540             :     // Velocity of the moving mesh on the dg grid, if applicable.
     541             :     const std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>&
     542             :         mesh_velocity_dg = db::get<domain::Tags::MeshVelocity<3>>(*box);
     543             :     // Inverse jacobian, to be projected on faces
     544             :     const auto& inv_jacobian_dg =
     545             :         db::get<domain::Tags::InverseJacobian<3, Frame::ElementLogical,
     546             :                                               Frame::Inertial>>(*box);
     547             :     const auto& det_inv_jacobian_dg = db::get<
     548             :         domain::Tags::DetInvJacobian<Frame::ElementLogical, Frame::Inertial>>(
     549             :         *box);
     550             : 
     551             :     // GH+GRMHD is a bit different.
     552             :     // 1. Compute GH time derivative, since this will also give us lapse, shift,
     553             :     //    etc. that we need to reconstruct.
     554             :     // 2. Compute d_t Pi_{ab} source terms from MHD (or do we wait until post
     555             :     //    MHD source terms?)
     556             :     // 3. Reconstruct MHD+spacetime vars to interfaces
     557             :     // 4. Compute MHD time derivatives.
     558             :     //
     559             :     // Compute FD GH derivatives with neighbor data
     560             :     // Use highest possible FD order for number of GZ, 2 * (ghost_zone_size)
     561             :     const auto& evolved_vars = db::get<evolved_vars_tag>(*box);
     562             :     Variables<db::wrap_tags_in<::Tags::deriv, gradients_tags, tmpl::size_t<3>,
     563             :                                Frame::Inertial>>
     564             :         cell_centered_gh_derivs{num_pts};
     565             :     grmhd::GhValenciaDivClean::fd::spacetime_derivatives<System>(
     566             :         make_not_null(&cell_centered_gh_derivs), evolved_vars,
     567             :         db::get<evolution::dg::subcell::Tags::GhostDataForReconstruction<3>>(
     568             :             *box),
     569             :         recons.ghost_zone_size() * 2, subcell_mesh,
     570             :         cell_centered_logical_to_inertial_inv_jacobian);
     571             : 
     572             :     // Now package the data and compute the correction
     573             :     //
     574             :     // Note: Assumes a the GH and GRMHD corrections can be invoked separately.
     575             :     // This is reasonable since the systems are a tensor product system.
     576             :     const auto& base_boundary_correction =
     577             :         db::get<evolution::Tags::BoundaryCorrection>(*box);
     578             :     using derived_boundary_corrections =
     579             :         tmpl::at<typename metavariables::factory_creation::factory_classes,
     580             :                  evolution::BoundaryCorrection>;
     581             :     std::array<Variables<grmhd_evolved_vars_tags>, 3> boundary_corrections{};
     582             :     call_with_dynamic_type<void, derived_boundary_corrections>(
     583             :         &base_boundary_correction, [&](const auto* gh_grmhd_correction) {
     584             :           // Need the GH packaged tags to avoid projecting them.
     585             :           using gh_dg_package_field_tags = typename std::decay_t<
     586             :               decltype(gh_grmhd_correction
     587             :                            ->gh_correction())>::dg_package_field_tags;
     588             :           // Only apply correction to GRMHD variables.
     589             :           const auto& boundary_correction =
     590             :               gh_grmhd_correction->valencia_correction();
     591             :           using DerivedCorrection = std::decay_t<decltype(boundary_correction)>;
     592             :           using dg_package_data_temporary_tags =
     593             :               typename DerivedCorrection::dg_package_data_temporary_tags;
     594             : 
     595             :           using dg_package_data_argument_tags = tmpl::append<
     596             :               evolved_vars_tags, recons_prim_tags, fluxes_tags,
     597             :               tmpl::remove_duplicates<tmpl::push_back<
     598             :                   dg_package_data_temporary_tags,
     599             :                   gr::Tags::SpatialMetric<DataVector, 3>,
     600             :                   gr::Tags::SqrtDetSpatialMetric<DataVector>,
     601             :                   gr::Tags::InverseSpatialMetric<DataVector, 3>,
     602             :                   evolution::dg::Actions::detail::NormalVector<3>>>>;
     603             : 
     604             :           // Computed prims and cons on face via reconstruction
     605             :           auto package_data_argvars_lower_face = make_array<3>(
     606             :               Variables<dg_package_data_argument_tags>(reconstructed_num_pts));
     607             :           auto package_data_argvars_upper_face = make_array<3>(
     608             :               Variables<dg_package_data_argument_tags>(reconstructed_num_pts));
     609             : 
     610             :           // Reconstruct data to the face
     611             :           call_with_dynamic_type<
     612             :               void, typename grmhd::GhValenciaDivClean::fd::Reconstructor<
     613             :                         System>::creatable_classes>(
     614             :               &recons, [&box, &package_data_argvars_lower_face,
     615             :                         &package_data_argvars_upper_face,
     616             :                         &reconstruction_order](const auto& reconstructor) {
     617             :                 using ReconstructorType =
     618             :                     std::decay_t<decltype(*reconstructor)>;
     619             :                 db::apply<
     620             :                     typename ReconstructorType::reconstruction_argument_tags>(
     621             :                     [&package_data_argvars_lower_face,
     622             :                      &package_data_argvars_upper_face, &reconstructor,
     623             :                      &reconstruction_order](const auto&... args) {
     624             :                       if constexpr (ReconstructorType::use_adaptive_order) {
     625             :                         reconstructor->reconstruct(
     626             :                             make_not_null(&package_data_argvars_lower_face),
     627             :                             make_not_null(&package_data_argvars_upper_face),
     628             :                             make_not_null(&reconstruction_order), args...);
     629             :                       } else {
     630             :                         (void)reconstruction_order;
     631             :                         reconstructor->reconstruct(
     632             :                             make_not_null(&package_data_argvars_lower_face),
     633             :                             make_not_null(&package_data_argvars_upper_face),
     634             :                             args...);
     635             :                       }
     636             :                     },
     637             :                     *box);
     638             :               });
     639             : 
     640             :           using dg_package_field_tags =
     641             :               typename DerivedCorrection::dg_package_field_tags;
     642             :           // Allocated outside for loop to reduce allocations
     643             :           Variables<dg_package_field_tags> upper_packaged_data{
     644             :               reconstructed_num_pts};
     645             :           Variables<dg_package_field_tags> lower_packaged_data{
     646             :               reconstructed_num_pts};
     647             : 
     648             :           // Compute fluxes on faces
     649             :           for (size_t i = 0; i < 3; ++i) {
     650             :             auto& vars_upper_face = gsl::at(package_data_argvars_upper_face, i);
     651             :             auto& vars_lower_face = gsl::at(package_data_argvars_lower_face, i);
     652             :             grmhd::ValenciaDivClean::subcell::compute_fluxes(
     653             :                 make_not_null(&vars_upper_face));
     654             :             grmhd::ValenciaDivClean::subcell::compute_fluxes(
     655             :                 make_not_null(&vars_lower_face));
     656             : 
     657             :             // Build extents of mesh shifted by half a grid cell in direction i
     658             :             const unsigned long& num_subcells_1d = subcell_mesh.extents(0);
     659             :             Index<3> face_mesh_extents(std::array<size_t, 3>{
     660             :                 num_subcells_1d, num_subcells_1d, num_subcells_1d});
     661             :             face_mesh_extents[i] = num_subcells_1d + 1;
     662             :             // Add moving mesh corrections to the fluxes, if needed
     663             :             std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>
     664             :                 mesh_velocity_on_face = {};
     665             :             if (mesh_velocity_dg.has_value()) {
     666             :               // Project mesh velocity on face mesh.
     667             :               // Can we get away with only doing the normal component? It
     668             :               // is also used in the packaged data...
     669             :               mesh_velocity_on_face = tnsr::I<DataVector, 3, Frame::Inertial>{
     670             :                   reconstructed_num_pts};
     671             :               for (size_t j = 0; j < 3; j++) {
     672             :                 // j^th component of the velocity on the i^th directed face
     673             :                 mesh_velocity_on_face.value().get(j) =
     674             :                     evolution::dg::subcell::fd::project_to_faces(
     675             :                         mesh_velocity_dg.value().get(j), dg_mesh,
     676             :                         face_mesh_extents, i);
     677             :               }
     678             : 
     679             :               tmpl::for_each<grmhd_evolved_vars_tags>(
     680             :                   [&vars_upper_face, &vars_lower_face,
     681             :                    &mesh_velocity_on_face](auto tag_v) {
     682             :                     using tag = tmpl::type_from<decltype(tag_v)>;
     683             :                     using flux_tag =
     684             :                         ::Tags::Flux<tag, tmpl::size_t<3>, Frame::Inertial>;
     685             :                     using FluxTensor = typename flux_tag::type;
     686             :                     const auto& var_upper = get<tag>(vars_upper_face);
     687             :                     const auto& var_lower = get<tag>(vars_lower_face);
     688             :                     auto& flux_upper = get<flux_tag>(vars_upper_face);
     689             :                     auto& flux_lower = get<flux_tag>(vars_lower_face);
     690             :                     for (size_t storage_index = 0;
     691             :                          storage_index < var_upper.size(); ++storage_index) {
     692             :                       const auto tensor_index =
     693             :                           var_upper.get_tensor_index(storage_index);
     694             :                       for (size_t j = 0; j < 3; j++) {
     695             :                         const auto flux_storage_index =
     696             :                             FluxTensor::get_storage_index(
     697             :                                 prepend(tensor_index, j));
     698             :                         flux_upper[flux_storage_index] -=
     699             :                             mesh_velocity_on_face.value().get(j) *
     700             :                             var_upper[storage_index];
     701             :                         flux_lower[flux_storage_index] -=
     702             :                             mesh_velocity_on_face.value().get(j) *
     703             :                             var_lower[storage_index];
     704             :                       }
     705             :                     }
     706             :                   });
     707             :             }
     708             : 
     709             :             // Normal vectors in curved spacetime normalized by inverse
     710             :             // spatial metric. Since we assume a Cartesian grid, this is
     711             :             // relatively easy. Note that we use the sign convention on
     712             :             // the normal vectors to be compatible with DG.
     713             :             //
     714             :             // Note that these normal vectors are on all faces inside the DG
     715             :             // element since there are a bunch of subcells. We don't use the
     716             :             // NormalCovectorAndMagnitude tag in the DataBox right now to avoid
     717             :             // conflicts with the DG solver. We can explore in the future if
     718             :             // it's possible to reuse that allocation.
     719             :             //
     720             :             // The unnormalized normal vector is
     721             :             // n_j = d \xi^{\hat i}/dx^j
     722             :             // with "i" the current face.
     723             :             tnsr::i<DataVector, 3, Frame::Inertial> lower_outward_conormal{
     724             :                 reconstructed_num_pts, 0.0};
     725             :             for (size_t j = 0; j < 3; j++) {
     726             :               lower_outward_conormal.get(j) =
     727             :                   evolution::dg::subcell::fd::project_to_faces(
     728             :                       inv_jacobian_dg.get(i, j), dg_mesh, face_mesh_extents, i);
     729             :             }
     730             :             const auto det_inv_jacobian_face =
     731             :                 evolution::dg::subcell::fd::project_to_faces(
     732             :                     get(det_inv_jacobian_dg), dg_mesh, face_mesh_extents, i);
     733             : 
     734             :             const Scalar<DataVector> normalization{sqrt(get(
     735             :                 dot_product(lower_outward_conormal, lower_outward_conormal,
     736             :                             get<gr::Tags::InverseSpatialMetric<DataVector, 3>>(
     737             :                                 vars_upper_face))))};
     738             :             for (size_t j = 0; j < 3; j++) {
     739             :               lower_outward_conormal.get(j) =
     740             :                   lower_outward_conormal.get(j) / get(normalization);
     741             :             }
     742             : 
     743             :             tnsr::i<DataVector, 3, Frame::Inertial> upper_outward_conormal{
     744             :                 reconstructed_num_pts, 0.0};
     745             :             for (size_t j = 0; j < 3; j++) {
     746             :               upper_outward_conormal.get(j) = -lower_outward_conormal.get(j);
     747             :             }
     748             :             // Note: we probably should compute the normal vector in addition to
     749             :             // the co-vector. Not a huge issue since we'll get an FPE right now
     750             :             // if it's used by a Riemann solver.
     751             : 
     752             :             // Compute the packaged data
     753             :             using dg_package_data_projected_tags = tmpl::append<
     754             :                 grmhd_evolved_vars_tags, fluxes_tags,
     755             :                 dg_package_data_temporary_tags,
     756             :                 typename DerivedCorrection::dg_package_data_primitive_tags>;
     757             :             evolution::dg::Actions::detail::dg_package_data<System>(
     758             :                 make_not_null(&upper_packaged_data),
     759             :                 dynamic_cast<const DerivedCorrection&>(boundary_correction),
     760             :                 vars_upper_face, upper_outward_conormal, mesh_velocity_on_face,
     761             :                 *box, typename DerivedCorrection::dg_package_data_volume_tags{},
     762             :                 dg_package_data_projected_tags{});
     763             : 
     764             :             evolution::dg::Actions::detail::dg_package_data<System>(
     765             :                 make_not_null(&lower_packaged_data),
     766             :                 dynamic_cast<const DerivedCorrection&>(boundary_correction),
     767             :                 vars_lower_face, lower_outward_conormal, mesh_velocity_on_face,
     768             :                 *box, typename DerivedCorrection::dg_package_data_volume_tags{},
     769             :                 dg_package_data_projected_tags{});
     770             : 
     771             :             // Now need to check if any of our neighbors are doing DG,
     772             :             // because if so then we need to use whatever boundary data
     773             :             // they sent instead of what we computed locally.
     774             :             //
     775             :             // Note: We could check this beforehand to avoid the extra
     776             :             // work of reconstruction and flux computations at the
     777             :             // boundaries.
     778             :             evolution::dg::subcell::correct_package_data<true>(
     779             :                 make_not_null(&lower_packaged_data),
     780             :                 make_not_null(&upper_packaged_data), i, element, subcell_mesh,
     781             :                 db::get<evolution::dg::Tags::MortarData<3>>(*box),
     782             :                 Variables<gh_dg_package_field_tags>::
     783             :                     number_of_independent_components);
     784             : 
     785             :             // Compute the corrections on the faces. We only need to
     786             :             // compute this once because we can just flip the normal
     787             :             // vectors then
     788             :             gsl::at(boundary_corrections, i).initialize(reconstructed_num_pts);
     789             :             evolution::dg::subcell::compute_boundary_terms(
     790             :                 make_not_null(&gsl::at(boundary_corrections, i)),
     791             :                 dynamic_cast<const DerivedCorrection&>(boundary_correction),
     792             :                 upper_packaged_data, lower_packaged_data, db::as_access(*box),
     793             :                 typename DerivedCorrection::dg_boundary_terms_volume_tags{});
     794             :             // We need to multiply by the normal vector normalization
     795             :             gsl::at(boundary_corrections, i) *= get(normalization);
     796             :             // Also multiply by determinant of Jacobian, following Eq.(34)
     797             :             // of 2109.11645
     798             :             gsl::at(boundary_corrections, i) *= 1.0 / det_inv_jacobian_face;
     799             :           }
     800             :         });
     801             : 
     802             :     // Now compute the actual time derivatives.
     803             :     using gh_variables_tags =
     804             :         typename System::gh_system::variables_tag::tags_list;
     805             :     using gh_gradient_tags = typename TimeDerivativeTerms::gh_gradient_tags;
     806             :     using gh_temporary_tags = typename TimeDerivativeTerms::gh_temp_tags;
     807             :     using gh_extra_tags =
     808             :         tmpl::list<gr::Tags::SpacetimeMetric<DataVector, 3>,
     809             :                    gh::Tags::Pi<DataVector, 3>, gh::Tags::Phi<DataVector, 3>,
     810             :                    ::gh::Tags::ConstraintGamma0, ::gh::Tags::ConstraintGamma1,
     811             :                    ::gh::Tags::ConstraintGamma2>;
     812             :     using grmhd_source_tags =
     813             :         tmpl::transform<ValenciaDivClean::ComputeSources::return_tags,
     814             :                         tmpl::bind<db::remove_tag_prefix, tmpl::_1>>;
     815             :     using grmhd_source_argument_tags =
     816             :         ValenciaDivClean::ComputeSources::argument_tags;
     817             :     detail::ComputeTimeDerivImpl<
     818             :         gh_variables_tags, gh_temporary_tags, gh_gradient_tags, gh_extra_tags,
     819             :         grmhd_evolved_vars_tags, grmhd_source_tags, grmhd_source_argument_tags,
     820             :         System>::apply(box, inertial_coords,
     821             :                        db::get<evolution::dg::subcell::fd::Tags::
     822             :                                    DetInverseJacobianLogicalToInertial>(*box),
     823             :                        cell_centered_logical_to_inertial_inv_jacobian,
     824             :                        one_over_delta_xi, boundary_corrections,
     825             :                        cell_centered_gh_derivs);
     826             :     evolution::dg::subcell::store_reconstruction_order_in_databox(
     827             :         box, reconstruction_order);
     828             :   }
     829             : };
     830             : }  // namespace grmhd::GhValenciaDivClean::subcell

Generated by: LCOV version 1.14