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