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