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