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