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