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 <vector>
8 :
9 : #include "ControlSystem/Actions/InitializeMeasurements.hpp"
10 : #include "ControlSystem/Actions/LimitTimeStep.hpp"
11 : #include "ControlSystem/Component.hpp"
12 : #include "ControlSystem/Measurements/BNSCenterOfMass.hpp"
13 : #include "ControlSystem/Metafunctions.hpp"
14 : #include "ControlSystem/Systems/Expansion.hpp"
15 : #include "ControlSystem/Systems/Rotation.hpp"
16 : #include "ControlSystem/Systems/Translation.hpp"
17 : #include "ControlSystem/Trigger.hpp"
18 : #include "DataStructures/DataBox/PrefixHelpers.hpp"
19 : #include "DataStructures/DataBox/Tag.hpp"
20 : #include "DataStructures/Tensor/EagerMath/RaiseOrLowerIndex.hpp"
21 : #include "DataStructures/Tensor/IndexType.hpp"
22 : #include "Domain/Creators/BinaryCompactObject.hpp"
23 : #include "Domain/Creators/Factory3D.hpp"
24 : #include "Domain/Tags.hpp"
25 : #include "Evolution/Actions/RunEventsAndDenseTriggers.hpp"
26 : #include "Evolution/Actions/RunEventsAndTriggers.hpp"
27 : #include "Evolution/ComputeTags.hpp"
28 : #include "Evolution/Conservative/UpdateConservatives.hpp"
29 : #include "Evolution/Conservative/UpdatePrimitives.hpp"
30 : #include "Evolution/DgSubcell/Actions/Initialize.hpp"
31 : #include "Evolution/DgSubcell/Actions/Labels.hpp"
32 : #include "Evolution/DgSubcell/Actions/ReconstructionCommunication.hpp"
33 : #include "Evolution/DgSubcell/Actions/SelectNumericalMethod.hpp"
34 : #include "Evolution/DgSubcell/Actions/TakeTimeStep.hpp"
35 : #include "Evolution/DgSubcell/Actions/TciAndRollback.hpp"
36 : #include "Evolution/DgSubcell/Actions/TciAndSwitchToDg.hpp"
37 : #include "Evolution/DgSubcell/CartesianFluxDivergence.hpp"
38 : #include "Evolution/DgSubcell/ComputeBoundaryTerms.hpp"
39 : #include "Evolution/DgSubcell/CorrectPackagedData.hpp"
40 : #include "Evolution/DgSubcell/GetActiveTag.hpp"
41 : #include "Evolution/DgSubcell/NeighborReconstructedFaceSolution.hpp"
42 : #include "Evolution/DgSubcell/PerssonTci.hpp"
43 : #include "Evolution/DgSubcell/PrepareNeighborData.hpp"
44 : #include "Evolution/DgSubcell/SetInterpolators.hpp"
45 : #include "Evolution/DgSubcell/Tags/ObserverCoordinates.hpp"
46 : #include "Evolution/DgSubcell/Tags/ObserverMesh.hpp"
47 : #include "Evolution/DgSubcell/Tags/ObserverMeshVelocity.hpp"
48 : #include "Evolution/DgSubcell/Tags/TciStatus.hpp"
49 : #include "Evolution/DgSubcell/TwoMeshRdmpTci.hpp"
50 : #include "Evolution/DiscontinuousGalerkin/Actions/ApplyBoundaryCorrections.hpp"
51 : #include "Evolution/DiscontinuousGalerkin/Actions/ComputeTimeDerivative.hpp"
52 : #include "Evolution/DiscontinuousGalerkin/DgElementArray.hpp"
53 : #include "Evolution/DiscontinuousGalerkin/Initialization/Mortars.hpp"
54 : #include "Evolution/DiscontinuousGalerkin/Initialization/QuadratureTag.hpp"
55 : #include "Evolution/DiscontinuousGalerkin/Limiters/Minmod.hpp"
56 : #include "Evolution/DiscontinuousGalerkin/Limiters/Tags.hpp"
57 : #include "Evolution/DiscontinuousGalerkin/Limiters/Weno.hpp"
58 : #include "Evolution/DiscontinuousGalerkin/UsingSubcell.hpp"
59 : #include "Evolution/Initialization/ConservativeSystem.hpp"
60 : #include "Evolution/Initialization/DgDomain.hpp"
61 : #include "Evolution/Initialization/Evolution.hpp"
62 : #include "Evolution/Initialization/Limiter.hpp"
63 : #include "Evolution/Initialization/NonconservativeSystem.hpp"
64 : #include "Evolution/Initialization/SetVariables.hpp"
65 : #include "Evolution/NumericInitialData.hpp"
66 : #include "Evolution/Systems/Cce/Callbacks/DumpBondiSachsOnWorldtube.hpp"
67 : #include "Evolution/Systems/GeneralizedHarmonic/AllSolutions.hpp"
68 : #include "Evolution/Systems/GeneralizedHarmonic/BoundaryConditions/Factory.hpp"
69 : #include "Evolution/Systems/GeneralizedHarmonic/BoundaryCorrections/Factory.hpp"
70 : #include "Evolution/Systems/GeneralizedHarmonic/Equations.hpp"
71 : #include "Evolution/Systems/GeneralizedHarmonic/GaugeSourceFunctions/Factory.hpp"
72 : #include "Evolution/Systems/GeneralizedHarmonic/GaugeSourceFunctions/Gauges.hpp"
73 : #include "Evolution/Systems/GeneralizedHarmonic/GaugeSourceFunctions/SetPiAndPhiFromConstraints.hpp"
74 : #include "Evolution/Systems/GeneralizedHarmonic/GaugeSourceFunctions/Tags/GaugeCondition.hpp"
75 : #include "Evolution/Systems/GeneralizedHarmonic/Initialize.hpp"
76 : #include "Evolution/Systems/GeneralizedHarmonic/System.hpp"
77 : #include "Evolution/Systems/GeneralizedHarmonic/Tags.hpp"
78 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Actions/SetInitialData.hpp"
79 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/AllSolutions.hpp"
80 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/BoundaryConditions/Factory.hpp"
81 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/BoundaryCorrections/ProductOfCorrections.hpp"
82 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/FiniteDifference/Tag.hpp"
83 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/SetPiAndPhiFromConstraints.hpp"
84 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/FixConservativesAndComputePrims.hpp"
85 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/NeighborPackagedData.hpp"
86 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/PrimitiveGhostData.hpp"
87 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/PrimsAfterRollback.hpp"
88 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/ResizeAndComputePrimitives.hpp"
89 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/TciOnDgGrid.hpp"
90 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/TciOnFdGrid.hpp"
91 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/Subcell/TimeDerivative.hpp"
92 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/System.hpp"
93 : #include "Evolution/Systems/GrMhd/GhValenciaDivClean/TimeDerivativeTerms.hpp"
94 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/AllSolutions.hpp"
95 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/BoundaryConditions/Factory.hpp"
96 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/BoundaryCorrections/Factory.hpp"
97 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/FixConservatives.hpp"
98 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/KastaunEtAl.hpp"
99 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/NewmanHamlin.hpp"
100 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/PalenzuelaEtAl.hpp"
101 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/QuadrupoleFormula.hpp"
102 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/SetVariablesNeededFixingToFalse.hpp"
103 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Subcell/SetInitialRdmpData.hpp"
104 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/System.hpp"
105 : #include "Evolution/Systems/GrMhd/ValenciaDivClean/Tags.hpp"
106 : #include "Evolution/Tags/Filter.hpp"
107 : #include "Evolution/TypeTraits.hpp"
108 : #include "Evolution/VariableFixing/Actions.hpp"
109 : #include "Evolution/VariableFixing/FixToAtmosphere.hpp"
110 : #include "Evolution/VariableFixing/LimitLorentzFactor.hpp"
111 : #include "Evolution/VariableFixing/Tags.hpp"
112 : #include "IO/Importers/Actions/ReadVolumeData.hpp"
113 : #include "IO/Importers/Actions/ReceiveVolumeData.hpp"
114 : #include "IO/Importers/Actions/RegisterWithElementDataReader.hpp"
115 : #include "IO/Importers/ElementDataReader.hpp"
116 : #include "IO/Observer/Actions/ObserverRegistration.hpp"
117 : #include "IO/Observer/Actions/RegisterEvents.hpp"
118 : #include "IO/Observer/Helpers.hpp"
119 : #include "IO/Observer/ObserverComponent.hpp"
120 : #include "IO/Observer/Tags.hpp"
121 : #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp"
122 : #include "NumericalAlgorithms/DiscontinuousGalerkin/Tags.hpp"
123 : #include "NumericalAlgorithms/LinearOperators/ExponentialFilter.hpp"
124 : #include "NumericalAlgorithms/LinearOperators/PartialDerivatives.hpp"
125 : #include "Options/Protocols/FactoryCreation.hpp"
126 : #include "Options/String.hpp"
127 : #include "Parallel/Algorithms/AlgorithmSingleton.hpp"
128 : #include "Parallel/Local.hpp"
129 : #include "Parallel/Phase.hpp"
130 : #include "Parallel/PhaseControl/ExecutePhaseChange.hpp"
131 : #include "Parallel/PhaseControl/Factory.hpp"
132 : #include "Parallel/PhaseControl/VisitAndReturn.hpp"
133 : #include "Parallel/PhaseDependentActionList.hpp"
134 : #include "Parallel/Protocols/RegistrationMetavariables.hpp"
135 : #include "Parallel/Reduction.hpp"
136 : #include "ParallelAlgorithms/Actions/AddComputeTags.hpp"
137 : #include "ParallelAlgorithms/Actions/AddSimpleTags.hpp"
138 : #include "ParallelAlgorithms/Actions/FilterAction.hpp"
139 : #include "ParallelAlgorithms/Actions/FunctionsOfTimeAreReady.hpp"
140 : #include "ParallelAlgorithms/Actions/InitializeItems.hpp"
141 : #include "ParallelAlgorithms/Actions/LimiterActions.hpp"
142 : #include "ParallelAlgorithms/Actions/MutateApply.hpp"
143 : #include "ParallelAlgorithms/Actions/TerminatePhase.hpp"
144 : #include "ParallelAlgorithms/ApparentHorizonFinder/Callbacks/FindApparentHorizon.hpp"
145 : #include "ParallelAlgorithms/ApparentHorizonFinder/InterpolationTarget.hpp"
146 : #include "ParallelAlgorithms/Events/Factory.hpp"
147 : #include "ParallelAlgorithms/Events/ObserveAtExtremum.hpp"
148 : #include "ParallelAlgorithms/EventsAndDenseTriggers/DenseTrigger.hpp"
149 : #include "ParallelAlgorithms/EventsAndDenseTriggers/DenseTriggers/Factory.hpp"
150 : #include "ParallelAlgorithms/EventsAndTriggers/Actions/RunEventsOnFailure.hpp"
151 : #include "ParallelAlgorithms/EventsAndTriggers/Completion.hpp"
152 : #include "ParallelAlgorithms/EventsAndTriggers/Event.hpp"
153 : #include "ParallelAlgorithms/EventsAndTriggers/EventsAndTriggers.hpp"
154 : #include "ParallelAlgorithms/EventsAndTriggers/LogicalTriggers.hpp"
155 : #include "ParallelAlgorithms/EventsAndTriggers/Trigger.hpp"
156 : #include "ParallelAlgorithms/Interpolation/Actions/CleanUpInterpolator.hpp"
157 : #include "ParallelAlgorithms/Interpolation/Actions/ElementInitInterpPoints.hpp"
158 : #include "ParallelAlgorithms/Interpolation/Actions/InitializeInterpolationTarget.hpp"
159 : #include "ParallelAlgorithms/Interpolation/Actions/InterpolationTargetReceiveVars.hpp"
160 : #include "ParallelAlgorithms/Interpolation/Actions/InterpolatorReceivePoints.hpp"
161 : #include "ParallelAlgorithms/Interpolation/Actions/InterpolatorReceiveVolumeData.hpp"
162 : #include "ParallelAlgorithms/Interpolation/Actions/InterpolatorRegisterElement.hpp"
163 : #include "ParallelAlgorithms/Interpolation/Actions/TryToInterpolate.hpp"
164 : #include "ParallelAlgorithms/Interpolation/Callbacks/ObserveTimeSeriesOnSurface.hpp"
165 : #include "ParallelAlgorithms/Interpolation/Events/Interpolate.hpp"
166 : #include "ParallelAlgorithms/Interpolation/Events/InterpolateWithoutInterpComponent.hpp"
167 : #include "ParallelAlgorithms/Interpolation/InterpolationTarget.hpp"
168 : #include "ParallelAlgorithms/Interpolation/Interpolator.hpp"
169 : #include "ParallelAlgorithms/Interpolation/Protocols/InterpolationTargetTag.hpp"
170 : #include "ParallelAlgorithms/Interpolation/Tags.hpp"
171 : #include "ParallelAlgorithms/Interpolation/Targets/KerrHorizon.hpp"
172 : #include "ParallelAlgorithms/Interpolation/Targets/Sphere.hpp"
173 : #include "PointwiseFunctions/AnalyticData/AnalyticData.hpp"
174 : #include "PointwiseFunctions/AnalyticData/GhGrMhd/Factory.hpp"
175 : #include "PointwiseFunctions/AnalyticData/GrMhd/BlastWave.hpp"
176 : #include "PointwiseFunctions/AnalyticData/GrMhd/BondiHoyleAccretion.hpp"
177 : #include "PointwiseFunctions/AnalyticData/GrMhd/InitialMagneticFields/Factory.hpp"
178 : #include "PointwiseFunctions/AnalyticData/GrMhd/MagneticFieldLoop.hpp"
179 : #include "PointwiseFunctions/AnalyticData/GrMhd/MagneticRotor.hpp"
180 : #include "PointwiseFunctions/AnalyticData/GrMhd/MagnetizedFmDisk.hpp"
181 : #include "PointwiseFunctions/AnalyticData/GrMhd/OrszagTangVortex.hpp"
182 : #include "PointwiseFunctions/AnalyticData/Tags.hpp"
183 : #include "PointwiseFunctions/AnalyticSolutions/AnalyticSolution.hpp"
184 : #include "PointwiseFunctions/AnalyticSolutions/GeneralRelativity/Factory.hpp"
185 : #include "PointwiseFunctions/AnalyticSolutions/GeneralRelativity/KerrSchild.hpp"
186 : #include "PointwiseFunctions/AnalyticSolutions/GeneralRelativity/WrappedGr.hpp"
187 : #include "PointwiseFunctions/AnalyticSolutions/GhGrMhd/Factory.hpp"
188 : #include "PointwiseFunctions/AnalyticSolutions/GhRelativisticEuler/Factory.hpp"
189 : #include "PointwiseFunctions/AnalyticSolutions/GrMhd/AlfvenWave.hpp"
190 : #include "PointwiseFunctions/AnalyticSolutions/GrMhd/BondiMichel.hpp"
191 : #include "PointwiseFunctions/AnalyticSolutions/GrMhd/KomissarovShock.hpp"
192 : #include "PointwiseFunctions/AnalyticSolutions/GrMhd/SmoothFlow.hpp"
193 : #include "PointwiseFunctions/AnalyticSolutions/RelativisticEuler/FishboneMoncriefDisk.hpp"
194 : #include "PointwiseFunctions/AnalyticSolutions/RelativisticEuler/SmoothFlow.hpp"
195 : #include "PointwiseFunctions/AnalyticSolutions/RelativisticEuler/TovStar.hpp"
196 : #include "PointwiseFunctions/AnalyticSolutions/Tags.hpp"
197 : #include "PointwiseFunctions/GeneralRelativity/Christoffel.hpp"
198 : #include "PointwiseFunctions/GeneralRelativity/DetAndInverseSpatialMetric.hpp"
199 : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/ConstraintGammas.hpp"
200 : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/ExtrinsicCurvature.hpp"
201 : #include "PointwiseFunctions/GeneralRelativity/GeneralizedHarmonic/SecondTimeDerivOfSpacetimeMetric.hpp"
202 : #include "PointwiseFunctions/GeneralRelativity/Lapse.hpp"
203 : #include "PointwiseFunctions/GeneralRelativity/Psi4Real.hpp"
204 : #include "PointwiseFunctions/GeneralRelativity/Ricci.hpp"
205 : #include "PointwiseFunctions/GeneralRelativity/Shift.hpp"
206 : #include "PointwiseFunctions/GeneralRelativity/SpatialMetric.hpp"
207 : #include "PointwiseFunctions/GeneralRelativity/Surfaces/Tags.hpp"
208 : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
209 : #include "PointwiseFunctions/GeneralRelativity/WeylElectric.hpp"
210 : #include "PointwiseFunctions/Hydro/EquationsOfState/Factory.hpp"
211 : #include "PointwiseFunctions/Hydro/LowerSpatialFourVelocity.hpp"
212 : #include "PointwiseFunctions/Hydro/MassFlux.hpp"
213 : #include "PointwiseFunctions/Hydro/MassWeightedFluidItems.hpp"
214 : #include "PointwiseFunctions/Hydro/Tags.hpp"
215 : #include "PointwiseFunctions/Hydro/TransportVelocity.hpp"
216 : #include "PointwiseFunctions/InitialDataUtilities/Tags/InitialData.hpp"
217 : #include "Time/Actions/AdvanceTime.hpp"
218 : #include "Time/Actions/CleanHistory.hpp"
219 : #include "Time/Actions/RecordTimeStepperData.hpp"
220 : #include "Time/Actions/SelfStartActions.hpp"
221 : #include "Time/Actions/UpdateU.hpp"
222 : #include "Time/ChangeSlabSize/Action.hpp"
223 : #include "Time/StepChoosers/Factory.hpp"
224 : #include "Time/StepChoosers/StepChooser.hpp"
225 : #include "Time/Tags/Time.hpp"
226 : #include "Time/Tags/TimeStepId.hpp"
227 : #include "Time/TimeSequence.hpp"
228 : #include "Time/TimeSteppers/Factory.hpp"
229 : #include "Time/TimeSteppers/LtsTimeStepper.hpp"
230 : #include "Time/TimeSteppers/TimeStepper.hpp"
231 : #include "Time/Triggers/TimeTriggers.hpp"
232 : #include "Utilities/ErrorHandling/Error.hpp"
233 : #include "Utilities/Functional.hpp"
234 : #include "Utilities/NoSuchType.hpp"
235 : #include "Utilities/ProtocolHelpers.hpp"
236 : #include "Utilities/TMPL.hpp"
237 :
238 : /// \cond
239 : namespace Frame {
240 :
241 : struct Inertial;
242 : } // namespace Frame
243 : namespace PUP {
244 : class er;
245 : } // namespace PUP
246 : namespace Parallel {
247 : template <typename Metavariables>
248 : class CProxy_GlobalCache;
249 : } // namespace Parallel
250 : /// \endcond
251 :
252 : template <bool UseDgSubcell>
253 0 : struct GhValenciaDivCleanDefaults {
254 : public:
255 0 : static constexpr size_t volume_dim = 3;
256 0 : using domain_frame = Frame::Inertial;
257 0 : static constexpr bool use_damped_harmonic_rollon = true;
258 0 : using temporal_id = Tags::TimeStepId;
259 0 : using TimeStepperBase = TimeStepper;
260 :
261 0 : static constexpr bool local_time_stepping =
262 : TimeStepperBase::local_time_stepping;
263 :
264 0 : using system = grmhd::GhValenciaDivClean::System;
265 0 : using analytic_variables_tags =
266 : typename system::primitive_variables_tag::tags_list;
267 0 : using analytic_solution_fields =
268 : tmpl::append<typename system::primitive_variables_tag::tags_list,
269 : typename system::gh_system::variables_tag::tags_list>;
270 0 : using ordered_list_of_primitive_recovery_schemes = tmpl::list<
271 : grmhd::ValenciaDivClean::PrimitiveRecoverySchemes::KastaunEtAl,
272 : grmhd::ValenciaDivClean::PrimitiveRecoverySchemes::NewmanHamlin,
273 : grmhd::ValenciaDivClean::PrimitiveRecoverySchemes::PalenzuelaEtAl>;
274 :
275 : // Do not limit the divergence-cleaning field Phi or the GH fields
276 0 : using limiter = Tags::Limiter<
277 : Limiters::Minmod<3, tmpl::list<grmhd::ValenciaDivClean::Tags::TildeD,
278 : grmhd::ValenciaDivClean::Tags::TildeTau,
279 : grmhd::ValenciaDivClean::Tags::TildeS<>,
280 : grmhd::ValenciaDivClean::Tags::TildeB<>>>>;
281 :
282 0 : using initialize_initial_data_dependent_quantities_actions = tmpl::list<
283 : gh::Actions::InitializeGhAnd3Plus1Variables<volume_dim>,
284 : Actions::MutateApply<tmpl::conditional_t<
285 : UseDgSubcell, grmhd::GhValenciaDivClean::SetPiAndPhiFromConstraints,
286 : gh::gauges::SetPiAndPhiFromConstraints<3>>>,
287 : Initialization::Actions::AddComputeTags<
288 : tmpl::list<gr::Tags::SqrtDetSpatialMetricCompute<
289 : DataVector, volume_dim, domain_frame>>>,
290 : VariableFixing::Actions::FixVariables<
291 : VariableFixing::FixToAtmosphere<volume_dim>>,
292 : VariableFixing::Actions::FixVariables<VariableFixing::LimitLorentzFactor>,
293 : Actions::UpdateConservatives,
294 : tmpl::conditional_t<
295 : UseDgSubcell,
296 : tmpl::list<
297 : Initialization::Actions::AddSimpleTags<
298 : grmhd::ValenciaDivClean::SetVariablesNeededFixingToFalse>,
299 :
300 : evolution::dg::subcell::Actions::SetAndCommunicateInitialRdmpData<
301 : volume_dim,
302 : grmhd::ValenciaDivClean::subcell::SetInitialRdmpData>,
303 : evolution::dg::subcell::Actions::ComputeAndSendTciOnInitialGrid<
304 : volume_dim, system,
305 : grmhd::GhValenciaDivClean::subcell::TciOnFdGrid>,
306 : evolution::dg::subcell::Actions::SetInitialGridFromTciData<
307 : volume_dim, system>,
308 : Actions::MutateApply<
309 : grmhd::GhValenciaDivClean::subcell::ResizeAndComputePrims<
310 : ordered_list_of_primitive_recovery_schemes>>,
311 :
312 : VariableFixing::Actions::FixVariables<
313 : VariableFixing::FixToAtmosphere<volume_dim>>,
314 : VariableFixing::Actions::FixVariables<
315 : VariableFixing::LimitLorentzFactor>,
316 : Actions::UpdateConservatives,
317 : Actions::MutateApply<
318 : grmhd::GhValenciaDivClean::SetPiAndPhiFromConstraints>>,
319 : tmpl::list<>>,
320 : Parallel::Actions::TerminatePhase>;
321 :
322 : // NOLINTNEXTLINE(google-runtime-references)
323 0 : void pup(PUP::er& /*p*/) {}
324 : };
325 :
326 : template <typename EvolutionMetavarsDerived, bool UseDgSubcell,
327 : bool UseControlSystems>
328 0 : struct GhValenciaDivCleanTemplateBase;
329 :
330 : template <bool UseDgSubcell, bool UseControlSystems,
331 : template <bool, typename...> class EvolutionMetavarsDerived,
332 : typename... InterpolationTargetTags>
333 0 : struct GhValenciaDivCleanTemplateBase<
334 : EvolutionMetavarsDerived<UseControlSystems, InterpolationTargetTags...>,
335 : UseDgSubcell, UseControlSystems>
336 : : public virtual GhValenciaDivCleanDefaults<UseDgSubcell> {
337 0 : using derived_metavars =
338 : EvolutionMetavarsDerived<UseControlSystems, InterpolationTargetTags...>;
339 0 : using defaults = GhValenciaDivCleanDefaults<UseDgSubcell>;
340 0 : static constexpr size_t volume_dim = defaults::volume_dim;
341 0 : using domain_frame = typename defaults::domain_frame;
342 0 : static constexpr bool use_damped_harmonic_rollon =
343 : defaults::use_damped_harmonic_rollon;
344 0 : using temporal_id = typename defaults::temporal_id;
345 0 : using TimeStepperBase = typename defaults::TimeStepperBase;
346 0 : static constexpr bool local_time_stepping = defaults::local_time_stepping;
347 0 : using system = typename defaults::system;
348 0 : using analytic_variables_tags = typename defaults::analytic_variables_tags;
349 0 : using analytic_solution_fields = typename defaults::analytic_solution_fields;
350 0 : using ordered_list_of_primitive_recovery_schemes =
351 : typename defaults::ordered_list_of_primitive_recovery_schemes;
352 0 : using limiter = typename defaults::limiter;
353 0 : using initialize_initial_data_dependent_quantities_actions =
354 : typename defaults::initialize_initial_data_dependent_quantities_actions;
355 :
356 0 : static constexpr bool use_dg_subcell = UseDgSubcell;
357 0 : static constexpr bool use_control_systems = UseControlSystems;
358 :
359 0 : using initial_data_list =
360 : ghmhd::GhValenciaDivClean::InitialData::initial_data_list;
361 :
362 0 : using initial_data_tag = evolution::initial_data::Tags::InitialData;
363 0 : using equation_of_state_tag = hydro::Tags::GrmhdEquationOfState;
364 0 : using measurement = control_system::measurements::BothNSCenters;
365 0 : using control_systems = tmpl::conditional_t<
366 : use_control_systems,
367 : tmpl::list<control_system::Systems::Rotation<3, measurement>,
368 : control_system::Systems::Expansion<2, measurement>,
369 : control_system::Systems::Translation<2, measurement>>,
370 : tmpl::list<>>;
371 :
372 0 : using interpolator_source_vars =
373 : tmpl::remove_duplicates<tmpl::flatten<tmpl::list<
374 : typename InterpolationTargetTags::vars_to_interpolate_to_target...>>>;
375 :
376 0 : using analytic_compute = evolution::Tags::AnalyticSolutionsCompute<
377 : volume_dim, analytic_solution_fields, use_dg_subcell, initial_data_list>;
378 0 : using error_compute = Tags::ErrorsCompute<analytic_solution_fields>;
379 0 : using error_tags = db::wrap_tags_in<Tags::Error, analytic_solution_fields>;
380 0 : using observe_fields = tmpl::append<
381 : typename system::variables_tag::tags_list,
382 : typename system::primitive_variables_tag::tags_list, error_tags,
383 : tmpl::list<
384 : hydro::Tags::TildeDInHalfPlaneCompute<
385 : DataVector, volume_dim,
386 : ::hydro::HalfPlaneIntegralMask::PositiveXOnly,
387 : Events::Tags::ObserverCoordinates<3, Frame::Grid>>,
388 : hydro::Tags::TildeDInHalfPlaneCompute<
389 : DataVector, volume_dim,
390 : ::hydro::HalfPlaneIntegralMask::NegativeXOnly,
391 : Events::Tags::ObserverCoordinates<3, Frame::Grid>>,
392 : hydro::Tags::MassWeightedCoordsCompute<
393 : DataVector, volume_dim,
394 : ::hydro::HalfPlaneIntegralMask::PositiveXOnly,
395 : Events::Tags::ObserverCoordinates<3, Frame::Grid>,
396 : Events::Tags::ObserverCoordinates<3, Frame::Inertial>,
397 : Frame::Inertial>,
398 : hydro::Tags::MassWeightedCoordsCompute<
399 : DataVector, volume_dim,
400 : ::hydro::HalfPlaneIntegralMask::NegativeXOnly,
401 : Events::Tags::ObserverCoordinates<3, Frame::Grid>,
402 : Events::Tags::ObserverCoordinates<3, Frame::Inertial>,
403 : Frame::Inertial>,
404 :
405 : hydro::Tags::MassWeightedInternalEnergyCompute<DataVector>,
406 : hydro::Tags::MassWeightedKineticEnergyCompute<DataVector>,
407 : hydro::Tags::TildeDUnboundUtCriterionCompute<DataVector, volume_dim,
408 : domain_frame>,
409 : hydro::Tags::LowerSpatialFourVelocityCompute,
410 : hydro::Tags::MassWeightedCoordsCompute<
411 : DataVector, volume_dim, ::hydro::HalfPlaneIntegralMask::None,
412 : Events::Tags::ObserverCoordinates<3, Frame::Grid>,
413 : Events::Tags::ObserverCoordinates<3, Frame::Inertial>,
414 : Frame::Inertial>,
415 : gr::Tags::SpacetimeNormalOneFormCompute<DataVector, volume_dim,
416 : domain_frame>,
417 : gr::Tags::SpacetimeNormalVectorCompute<DataVector, volume_dim,
418 : domain_frame>,
419 : gr::Tags::InverseSpacetimeMetricCompute<DataVector, volume_dim,
420 : domain_frame>,
421 : gr::Tags::Lapse<DataVector>,
422 : gr::Tags::Shift<DataVector, volume_dim, domain_frame>,
423 : gr::Tags::SpatialMetric<DataVector, volume_dim, domain_frame>,
424 : gr::Tags::DetSpatialMetric<DataVector>,
425 : gr::Tags::InverseSpatialMetric<DataVector, volume_dim>,
426 : gh::Tags::ExtrinsicCurvatureCompute<volume_dim, domain_frame>,
427 : gh::Tags::DerivSpatialMetricCompute<volume_dim, ::Frame::Inertial>,
428 : gr::Tags::SpatialChristoffelFirstKindCompute<DataVector, volume_dim,
429 : ::Frame::Inertial>,
430 : gr::Tags::SpatialChristoffelSecondKindCompute<DataVector, volume_dim,
431 : ::Frame::Inertial>,
432 : ::Tags::DerivTensorCompute<
433 : gr::Tags::SpatialChristoffelSecondKind<DataVector, volume_dim>,
434 : ::Events::Tags::ObserverInverseJacobian<
435 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
436 : ::Events::Tags::ObserverMesh<volume_dim>>,
437 : gr::Tags::SpatialRicciCompute<DataVector, volume_dim,
438 : ::Frame::Inertial>,
439 : gr::Tags::SpatialRicciScalarCompute<DataVector, volume_dim,
440 : ::Frame::Inertial>,
441 :
442 : // Constraints
443 : gh::Tags::GaugeConstraintCompute<volume_dim, domain_frame>,
444 : gh::Tags::TwoIndexConstraintCompute<volume_dim, ::Frame::Inertial>,
445 : gh::Tags::ThreeIndexConstraintCompute<volume_dim, ::Frame::Inertial>,
446 : gh::Tags::FourIndexConstraintCompute<3, ::Frame::Inertial>,
447 : gh::Tags::FConstraintCompute<3, ::Frame::Inertial>,
448 :
449 : // L2 norms of constraints
450 : ::Tags::PointwiseL2NormCompute<
451 : gh::Tags::GaugeConstraint<DataVector, volume_dim, domain_frame>>,
452 : ::Tags::PointwiseL2NormCompute<
453 : gh::Tags::TwoIndexConstraint<DataVector, volume_dim>>,
454 : ::Tags::PointwiseL2NormCompute<
455 : gh::Tags::ThreeIndexConstraint<DataVector, volume_dim>>,
456 : ::Tags::PointwiseL2NormCompute<
457 : gh::Tags::FourIndexConstraint<DataVector, 3>>,
458 : gh::Tags::ConstraintEnergyCompute<3, ::Frame::Inertial>,
459 : ::Tags::PointwiseL2NormCompute<gh::Tags::FConstraint<DataVector, 3>>,
460 : ::Tags::PointwiseL2NormCompute<
461 : gh::Tags::GaugeH<DataVector, volume_dim>>,
462 : ::Tags::PointwiseL2NormCompute<gh::Tags::Phi<DataVector, volume_dim>>,
463 : ::Tags::PointwiseL2NormCompute<
464 : ::Tags::deriv<gr::Tags::SpacetimeMetric<DataVector, volume_dim>,
465 : tmpl::size_t<volume_dim>, Frame::Inertial>>,
466 :
467 : // GW Tags
468 : grmhd::ValenciaDivClean::Tags::QuadrupoleMomentCompute<
469 : DataVector, volume_dim,
470 : ::Events::Tags::ObserverCoordinates<volume_dim, Frame::Inertial>>,
471 : grmhd::ValenciaDivClean::Tags::QuadrupoleMomentDerivativeCompute<
472 : DataVector, volume_dim,
473 : ::Events::Tags::ObserverCoordinates<volume_dim, Frame::Inertial>,
474 : hydro::Tags::SpatialVelocity<DataVector, volume_dim,
475 : Frame::Inertial>>,
476 : hydro::Tags::TransportVelocityCompute<DataVector, volume_dim,
477 : Frame::Inertial>,
478 : grmhd::ValenciaDivClean::Tags::QuadrupoleMomentDerivativeCompute<
479 : DataVector, volume_dim,
480 : ::Events::Tags::ObserverCoordinates<volume_dim, Frame::Inertial>,
481 : hydro::Tags::TransportVelocity<DataVector, volume_dim,
482 : Frame::Inertial>>,
483 : gh::Tags::TimeDerivLapseCompute<volume_dim, Frame::Inertial>,
484 : gh::Tags::TimeDerivShiftCompute<volume_dim, Frame::Inertial>,
485 : ::Tags::dt<gh::Tags::Phi<DataVector, volume_dim, Frame::Inertial>>,
486 : ::Tags::dt<gh::Tags::Pi<DataVector, volume_dim, Frame::Inertial>>,
487 : gh::Tags::SecondTimeDerivOfSpacetimeMetricCompute<volume_dim,
488 : Frame::Inertial>,
489 : ::Tags::DerivTensorCompute<
490 : gr::Tags::ExtrinsicCurvature<DataVector, 3>,
491 : ::Events::Tags::ObserverInverseJacobian<
492 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
493 : ::Events::Tags::ObserverMesh<volume_dim>>,
494 : gr::Tags::WeylElectricCompute<DataVector, 3, Frame::Inertial>,
495 : gr::Tags::Psi4RealCompute<Frame::Inertial>
496 :
497 : >,
498 : tmpl::conditional_t<
499 : use_dg_subcell,
500 : tmpl::list<evolution::dg::subcell::Tags::TciStatusCompute<volume_dim>,
501 : evolution::dg::subcell::Tags::ObserverCoordinatesCompute<
502 : volume_dim, Frame::ElementLogical>,
503 : evolution::dg::subcell::Tags::ObserverCoordinatesCompute<
504 : volume_dim, Frame::Grid>,
505 : evolution::dg::subcell::Tags::ObserverCoordinatesCompute<
506 : volume_dim, Frame::Inertial>>,
507 : tmpl::list<::Events::Tags::ObserverCoordinatesCompute<
508 : volume_dim, Frame::ElementLogical>,
509 : ::Events::Tags::ObserverCoordinatesCompute<volume_dim,
510 : Frame::Grid>,
511 : ::Events::Tags::ObserverCoordinatesCompute<
512 : volume_dim, Frame::Inertial>>>>;
513 0 : using integrand_fields = tmpl::append<
514 : typename system::variables_tag::tags_list,
515 : tmpl::list<
516 : // Control system tags
517 : hydro::Tags::TildeDInHalfPlaneCompute<
518 : DataVector, volume_dim,
519 : ::hydro::HalfPlaneIntegralMask::PositiveXOnly,
520 : Events::Tags::ObserverCoordinates<3, Frame::Grid>>,
521 : hydro::Tags::TildeDInHalfPlaneCompute<
522 : DataVector, volume_dim,
523 : ::hydro::HalfPlaneIntegralMask::NegativeXOnly,
524 : Events::Tags::ObserverCoordinates<3, Frame::Grid>>,
525 : hydro::Tags::MassWeightedCoordsCompute<
526 : DataVector, volume_dim,
527 : ::hydro::HalfPlaneIntegralMask::PositiveXOnly,
528 : Events::Tags::ObserverCoordinates<3, Frame::Grid>,
529 : Events::Tags::ObserverCoordinates<3, Frame::Inertial>,
530 : Frame::Inertial>,
531 : hydro::Tags::MassWeightedCoordsCompute<
532 : DataVector, volume_dim,
533 : ::hydro::HalfPlaneIntegralMask::NegativeXOnly,
534 : Events::Tags::ObserverCoordinates<3, Frame::Grid>,
535 : Events::Tags::ObserverCoordinates<3, Frame::Inertial>,
536 : Frame::Inertial>,
537 :
538 : // General tags
539 : hydro::Tags::MassWeightedInternalEnergyCompute<DataVector>,
540 : hydro::Tags::MassWeightedKineticEnergyCompute<DataVector>,
541 : hydro::Tags::TildeDUnboundUtCriterionCompute<DataVector, volume_dim,
542 : domain_frame>,
543 : hydro::Tags::MassWeightedCoordsCompute<
544 : DataVector, volume_dim, ::hydro::HalfPlaneIntegralMask::None,
545 : Events::Tags::ObserverCoordinates<3, Frame::Grid>,
546 : Events::Tags::ObserverCoordinates<3, Frame::Inertial>,
547 : Frame::Inertial>>>;
548 :
549 0 : using non_tensor_compute_tags = tmpl::append<
550 : tmpl::conditional_t<
551 : use_dg_subcell,
552 : tmpl::list<
553 : evolution::dg::subcell::Tags::ObserverMeshCompute<volume_dim>,
554 : evolution::dg::subcell::Tags::ObserverInverseJacobianCompute<
555 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
556 : evolution::dg::subcell::Tags::
557 : ObserverJacobianAndDetInvJacobianCompute<
558 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
559 : evolution::dg::subcell::Tags::ObserverMeshVelocityCompute<
560 : volume_dim>>,
561 : tmpl::list<::Events::Tags::ObserverMeshCompute<volume_dim>,
562 : ::Events::Tags::ObserverInverseJacobianCompute<
563 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
564 : ::Events::Tags::ObserverJacobianCompute<
565 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
566 : ::Events::Tags::ObserverDetInvJacobianCompute<
567 : Frame::ElementLogical, Frame::Inertial>>>,
568 : tmpl::list<analytic_compute, error_compute>,
569 : tmpl::list<::Tags::DerivCompute<
570 : typename system::variables_tag,
571 : ::Events::Tags::ObserverMesh<volume_dim>,
572 : ::Events::Tags::ObserverInverseJacobian<
573 : volume_dim, Frame::ElementLogical, Frame::Inertial>,
574 : typename system::gradient_variables>,
575 : gh::gauges::Tags::GaugeAndDerivativeCompute<volume_dim>>>;
576 :
577 0 : struct factory_creation
578 : : tt::ConformsTo<Options::protocols::FactoryCreation> {
579 : private:
580 0 : using boundary_conditions = grmhd::GhValenciaDivClean::BoundaryConditions::
581 : standard_boundary_conditions;
582 :
583 : public:
584 0 : using factory_classes = tmpl::map<
585 : tmpl::pair<DenseTrigger,
586 : tmpl::append<DenseTriggers::standard_dense_triggers,
587 : control_system::control_system_triggers<
588 : control_systems>>>,
589 : tmpl::pair<DomainCreator<volume_dim>,
590 : // Currently control systems can only be used with BCO
591 : // domains
592 : tmpl::conditional_t<
593 : use_control_systems,
594 : tmpl::list<::domain::creators::BinaryCompactObject>,
595 : domain_creators<volume_dim>>>,
596 : tmpl::pair<Event,
597 : tmpl::flatten<tmpl::list<
598 : Events::Completion,
599 : dg::Events::field_observations<
600 : volume_dim, observe_fields, non_tensor_compute_tags>,
601 : Events::ObserveAtExtremum<observe_fields,
602 : non_tensor_compute_tags>,
603 : Events::time_events<system>,
604 : control_system::metafunctions::control_system_events<
605 : control_systems>,
606 : intrp::Events::InterpolateWithoutInterpComponent<
607 : volume_dim, InterpolationTargetTags,
608 : interpolator_source_vars>...>>>,
609 : tmpl::pair<
610 : grmhd::GhValenciaDivClean::BoundaryConditions::BoundaryCondition,
611 : boundary_conditions>,
612 : tmpl::pair<
613 : grmhd::AnalyticData::InitialMagneticFields::InitialMagneticField,
614 : grmhd::AnalyticData::InitialMagneticFields::
615 : initial_magnetic_fields>,
616 : tmpl::pair<gh::gauges::GaugeCondition, gh::gauges::all_gauges>,
617 : tmpl::pair<evolution::initial_data::InitialData, initial_data_list>,
618 : // Restrict to monotonic time steppers in LTS to avoid control
619 : // systems deadlocking.
620 : tmpl::pair<
621 : LtsTimeStepper,
622 : tmpl::conditional_t<use_control_systems,
623 : TimeSteppers::monotonic_lts_time_steppers,
624 : TimeSteppers::lts_time_steppers>>,
625 : tmpl::pair<PhaseChange, PhaseControl::factory_creatable_classes>,
626 : tmpl::pair<StepChooser<StepChooserUse::LtsStep>,
627 : StepChoosers::standard_step_choosers<system>>,
628 : tmpl::pair<
629 : StepChooser<StepChooserUse::Slab>,
630 : StepChoosers::standard_slab_choosers<system, local_time_stepping>>,
631 : tmpl::pair<TimeSequence<double>,
632 : TimeSequences::all_time_sequences<double>>,
633 : tmpl::pair<TimeSequence<std::uint64_t>,
634 : TimeSequences::all_time_sequences<std::uint64_t>>,
635 : tmpl::pair<TimeStepper, TimeSteppers::time_steppers>,
636 : tmpl::pair<Trigger, tmpl::append<Triggers::logical_triggers,
637 : Triggers::time_triggers>>>;
638 : };
639 :
640 0 : using interpolation_target_tags = tmpl::list<InterpolationTargetTags...>;
641 :
642 0 : using observed_reduction_data_tags = observers::collect_reduction_data_tags<
643 : tmpl::at<typename factory_creation::factory_classes, Event>>;
644 :
645 0 : using const_global_cache_tags = tmpl::flatten<tmpl::list<
646 : tmpl::conditional_t<
647 : use_dg_subcell,
648 : tmpl::list<
649 : grmhd::GhValenciaDivClean::fd::Tags::Reconstructor,
650 : grmhd::GhValenciaDivClean::fd::Tags::FilterOptions,
651 : ::Tags::VariableFixer<grmhd::ValenciaDivClean::FixConservatives>,
652 : grmhd::ValenciaDivClean::subcell::Tags::TciOptions>,
653 : tmpl::list<>>,
654 : grmhd::ValenciaDivClean::Tags::PrimitiveFromConservativeOptions,
655 : gh::gauges::Tags::GaugeCondition, initial_data_tag,
656 : grmhd::ValenciaDivClean::Tags::ConstraintDampingParameter,
657 : equation_of_state_tag,
658 : gh::ConstraintDamping::Tags::DampingFunctionGamma0<volume_dim,
659 : Frame::Grid>,
660 : gh::ConstraintDamping::Tags::DampingFunctionGamma1<volume_dim,
661 : Frame::Grid>,
662 : gh::ConstraintDamping::Tags::DampingFunctionGamma2<volume_dim,
663 : Frame::Grid>>>;
664 :
665 0 : using dg_registration_list =
666 : tmpl::list<intrp::Actions::RegisterElementWithInterpolator,
667 : observers::Actions::RegisterEventsWithObservers>;
668 :
669 0 : static constexpr auto default_phase_order = std::array<Parallel::Phase, 8>{
670 : {Parallel::Phase::Initialization,
671 : Parallel::Phase::RegisterWithElementDataReader,
672 : Parallel::Phase::Register, Parallel::Phase::ImportInitialData,
673 : Parallel::Phase::InitializeInitialDataDependentQuantities,
674 : Parallel::Phase::InitializeTimeStepperHistory, Parallel::Phase::Evolve,
675 : Parallel::Phase::Exit}};
676 :
677 0 : struct SubcellOptions {
678 0 : using evolved_vars_tags = typename system::variables_tag::tags_list;
679 0 : using prim_tags = typename system::primitive_variables_tag::tags_list;
680 :
681 0 : static constexpr bool subcell_enabled = use_dg_subcell;
682 0 : static constexpr bool subcell_enabled_at_external_boundary = true;
683 :
684 : // We send `ghost_zone_size` cell-centered grid points for variable
685 : // reconstruction, of which we need `ghost_zone_size-1` for reconstruction
686 : // to the internal side of the element face, and `ghost_zone_size` for
687 : // reconstruction to the external side of the element face.
688 : template <typename DbTagsList>
689 0 : static constexpr size_t ghost_zone_size(
690 : const db::DataBox<DbTagsList>& box) {
691 : return db::get<grmhd::GhValenciaDivClean::fd::Tags::Reconstructor>(box)
692 : .ghost_zone_size();
693 : }
694 :
695 0 : using DgComputeSubcellNeighborPackagedData =
696 : grmhd::GhValenciaDivClean::subcell::NeighborPackagedData;
697 :
698 0 : using GhostVariables =
699 : grmhd::GhValenciaDivClean::subcell::PrimitiveGhostVariables;
700 : };
701 :
702 0 : using events_and_dense_triggers_subcell_postprocessors = tmpl::list<
703 : ::domain::CheckFunctionsOfTimeAreReadyPostprocessor,
704 : AlwaysReadyPostprocessor<
705 : grmhd::GhValenciaDivClean::subcell::FixConservativesAndComputePrims<
706 : ordered_list_of_primitive_recovery_schemes>>>;
707 :
708 0 : using dg_step_actions = tmpl::flatten<tmpl::list<
709 : evolution::dg::Actions::ComputeTimeDerivative<
710 : volume_dim, system, AllStepChoosers, local_time_stepping>,
711 : tmpl::conditional_t<
712 : local_time_stepping,
713 : tmpl::list<evolution::dg::Actions::ApplyLtsBoundaryCorrections<
714 : system, volume_dim, false>>,
715 : tmpl::list<
716 : evolution::dg::Actions::ApplyBoundaryCorrectionsToTimeDerivative<
717 : system, volume_dim, false>,
718 : Actions::RecordTimeStepperData<system>,
719 : Actions::UpdateU<system>>>,
720 : Actions::CleanHistory<system, local_time_stepping>,
721 : Limiters::Actions::SendData<derived_metavars>,
722 : Limiters::Actions::Limit<derived_metavars>,
723 : VariableFixing::Actions::FixVariables<
724 : grmhd::ValenciaDivClean::FixConservatives>,
725 : Actions::UpdatePrimitives>>;
726 :
727 0 : using dg_subcell_step_actions = tmpl::flatten<tmpl::list<
728 : evolution::dg::subcell::Actions::SelectNumericalMethod,
729 :
730 : Actions::Label<evolution::dg::subcell::Actions::Labels::BeginDg>,
731 : dg::Actions::Filter<::Filters::Exponential<0>,
732 : tmpl::list<gr::Tags::SpacetimeMetric<DataVector, 3>,
733 : gh::Tags::Pi<DataVector, 3>,
734 : gh::Tags::Phi<DataVector, 3>>>,
735 : evolution::dg::Actions::ComputeTimeDerivative<
736 : volume_dim, system, AllStepChoosers, local_time_stepping>,
737 : evolution::dg::Actions::ApplyBoundaryCorrectionsToTimeDerivative<
738 : system, volume_dim, false>,
739 : tmpl::conditional_t<
740 : local_time_stepping, tmpl::list<>,
741 : tmpl::list<Actions::RecordTimeStepperData<system>,
742 : evolution::Actions::RunEventsAndDenseTriggers<
743 : events_and_dense_triggers_subcell_postprocessors>,
744 : control_system::Actions::LimitTimeStep<control_systems>,
745 : Actions::UpdateU<system>>>,
746 : // Note: The primitive variables are computed as part of the TCI.
747 : evolution::dg::subcell::Actions::TciAndRollback<
748 : grmhd::GhValenciaDivClean::subcell::TciOnDgGrid<
749 : tmpl::front<ordered_list_of_primitive_recovery_schemes>>>,
750 : Actions::CleanHistory<system, local_time_stepping>,
751 : VariableFixing::Actions::FixVariables<
752 : VariableFixing::FixToAtmosphere<volume_dim>>,
753 : VariableFixing::Actions::FixVariables<VariableFixing::LimitLorentzFactor>,
754 : Actions::UpdateConservatives,
755 : Actions::Goto<evolution::dg::subcell::Actions::Labels::EndOfSolvers>,
756 :
757 : Actions::Label<evolution::dg::subcell::Actions::Labels::BeginSubcell>,
758 : evolution::dg::subcell::Actions::SendDataForReconstruction<
759 : volume_dim,
760 : grmhd::GhValenciaDivClean::subcell::PrimitiveGhostVariables,
761 : local_time_stepping>,
762 : evolution::dg::subcell::Actions::ReceiveDataForReconstruction<volume_dim>,
763 : Actions::Label<
764 : evolution::dg::subcell::Actions::Labels::BeginSubcellAfterDgRollback>,
765 : Actions::MutateApply<
766 : grmhd::GhValenciaDivClean::subcell::PrimsAfterRollback<
767 : ordered_list_of_primitive_recovery_schemes>>,
768 : evolution::dg::subcell::fd::Actions::TakeTimeStep<
769 : grmhd::GhValenciaDivClean::subcell::TimeDerivative>,
770 : Actions::RecordTimeStepperData<system>,
771 : evolution::Actions::RunEventsAndDenseTriggers<
772 : events_and_dense_triggers_subcell_postprocessors>,
773 : control_system::Actions::LimitTimeStep<control_systems>,
774 : Actions::UpdateU<system>,
775 : Actions::CleanHistory<system, local_time_stepping>,
776 : Actions::MutateApply<
777 : grmhd::GhValenciaDivClean::subcell::FixConservativesAndComputePrims<
778 : ordered_list_of_primitive_recovery_schemes>>,
779 : evolution::dg::subcell::Actions::TciAndSwitchToDg<
780 : grmhd::GhValenciaDivClean::subcell::TciOnFdGrid>,
781 : Actions::MutateApply<
782 : grmhd::GhValenciaDivClean::subcell::ResizeAndComputePrims<
783 : ordered_list_of_primitive_recovery_schemes>>,
784 : VariableFixing::Actions::FixVariables<
785 : VariableFixing::FixToAtmosphere<volume_dim>>,
786 : VariableFixing::Actions::FixVariables<VariableFixing::LimitLorentzFactor>,
787 : Actions::UpdateConservatives,
788 :
789 : Actions::Label<evolution::dg::subcell::Actions::Labels::EndOfSolvers>>>;
790 :
791 0 : using step_actions =
792 : tmpl::conditional_t<use_dg_subcell, dg_subcell_step_actions,
793 : dg_step_actions>;
794 :
795 0 : using initialization_actions = tmpl::list<
796 : Initialization::Actions::InitializeItems<
797 : Initialization::TimeStepping<derived_metavars, TimeStepperBase>,
798 : evolution::dg::Initialization::Domain<3, use_control_systems>,
799 : Initialization::TimeStepperHistory<derived_metavars>>,
800 : Initialization::Actions::ConservativeSystem<system>,
801 : // This conditional is untested and probably doesn't work if
802 : // `use_dg_subcell` is `false`
803 : tmpl::conditional_t<
804 : use_dg_subcell,
805 : tmpl::list<evolution::dg::subcell::Actions::SetSubcellGrid<
806 : volume_dim, system, true>,
807 : Actions::MutateApply<
808 : evolution::dg::subcell::SetInterpolators<volume_dim>>>,
809 : evolution::Initialization::Actions::SetVariables<
810 : ::domain::Tags::Coordinates<volume_dim, Frame::ElementLogical>>>,
811 : Initialization::Actions::AddComputeTags<
812 : StepChoosers::step_chooser_compute_tags<
813 : GhValenciaDivCleanTemplateBase, local_time_stepping>>,
814 : ::evolution::dg::Initialization::Mortars<volume_dim, system>,
815 : Initialization::Actions::Minmod<3>,
816 : evolution::Actions::InitializeRunEventsAndDenseTriggers,
817 : intrp::Actions::ElementInitInterpPoints<
818 : intrp::Tags::InterpPointInfo<derived_metavars>>,
819 : tmpl::conditional_t<
820 : use_control_systems,
821 : control_system::Actions::InitializeMeasurements<control_systems>,
822 : tmpl::list<>>,
823 : Parallel::Actions::TerminatePhase>;
824 :
825 0 : using import_initial_data_action_lists = tmpl::list<
826 : Parallel::PhaseActions<
827 : Parallel::Phase::RegisterWithElementDataReader,
828 : tmpl::list<importers::Actions::RegisterWithElementDataReader,
829 : Parallel::Actions::TerminatePhase>>,
830 : Parallel::PhaseActions<
831 : Parallel::Phase::ImportInitialData,
832 : tmpl::list<
833 : grmhd::GhValenciaDivClean::Actions::SetInitialData,
834 : grmhd::GhValenciaDivClean::Actions::ReceiveNumericInitialData,
835 : Parallel::Actions::TerminatePhase>>>;
836 :
837 0 : using dg_element_array_component = DgElementArray<
838 : derived_metavars,
839 : tmpl::flatten<tmpl::list<
840 : Parallel::PhaseActions<Parallel::Phase::Initialization,
841 : initialization_actions>,
842 : import_initial_data_action_lists,
843 : Parallel::PhaseActions<
844 : Parallel::Phase::InitializeInitialDataDependentQuantities,
845 : initialize_initial_data_dependent_quantities_actions>,
846 : Parallel::PhaseActions<
847 : Parallel::Phase::InitializeTimeStepperHistory,
848 : SelfStart::self_start_procedure<step_actions, system>>,
849 : Parallel::PhaseActions<Parallel::Phase::Register,
850 : tmpl::list<dg_registration_list,
851 : Parallel::Actions::TerminatePhase>>,
852 : Parallel::PhaseActions<
853 : Parallel::Phase::Evolve,
854 : tmpl::list<::domain::Actions::CheckFunctionsOfTimeAreReady,
855 : VariableFixing::Actions::FixVariables<
856 : VariableFixing::FixToAtmosphere<volume_dim>>,
857 : VariableFixing::Actions::FixVariables<
858 : VariableFixing::LimitLorentzFactor>,
859 : Actions::UpdateConservatives,
860 : evolution::Actions::RunEventsAndTriggers,
861 : Actions::ChangeSlabSize, step_actions,
862 : Actions::AdvanceTime,
863 : PhaseControl::Actions::ExecutePhaseChange>>,
864 : Parallel::PhaseActions<
865 : Parallel::Phase::PostFailureCleanup,
866 : tmpl::list<Actions::RunEventsOnFailure<Tags::Time>,
867 : Parallel::Actions::TerminatePhase>>>>>;
868 :
869 0 : struct registration
870 : : tt::ConformsTo<Parallel::protocols::RegistrationMetavariables> {
871 0 : using element_registrars =
872 : tmpl::map<tmpl::pair<dg_element_array_component, dg_registration_list>>;
873 : };
874 :
875 0 : using component_list = tmpl::flatten<tmpl::list<
876 : observers::Observer<derived_metavars>,
877 : observers::ObserverWriter<derived_metavars>,
878 : importers::ElementDataReader<derived_metavars>,
879 : control_system::control_components<derived_metavars, control_systems>,
880 : intrp::Interpolator<derived_metavars>,
881 : intrp::InterpolationTarget<derived_metavars, InterpolationTargetTags>...,
882 : dg_element_array_component>>;
883 : };
884 :
885 0 : struct BondiSachs : tt::ConformsTo<intrp::protocols::InterpolationTargetTag> {
886 0 : static std::string name() { return "BondiSachsInterpolation"; }
887 0 : using temporal_id = ::Tags::Time;
888 0 : using vars_to_interpolate_to_target =
889 : typename gh::System<3>::variables_tag::tags_list;
890 0 : using compute_target_points =
891 : intrp::TargetPoints::Sphere<BondiSachs, ::Frame::Inertial>;
892 0 : using post_interpolation_callbacks =
893 : tmpl::list<intrp::callbacks::DumpBondiSachsOnWorldtube<BondiSachs>>;
894 0 : using compute_items_on_target = tmpl::list<>;
895 : template <typename Metavariables>
896 0 : using interpolating_component =
897 : typename Metavariables::dg_element_array_component;
898 : };
|