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 :
8 : #include "DataStructures/Tensor/Tensor.hpp"
9 : #include "Elliptic/Systems/Xcts/Geometry.hpp"
10 : #include "Elliptic/Systems/Xcts/Tags.hpp"
11 : #include "NumericalAlgorithms/LinearOperators/Divergence.hpp"
12 : #include "NumericalAlgorithms/LinearOperators/PartialDerivatives.hpp"
13 : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
14 : #include "PointwiseFunctions/GeneralRelativity/Tags/Conformal.hpp"
15 : #include "Utilities/Gsl.hpp"
16 : #include "Utilities/TMPL.hpp"
17 :
18 : /// \cond
19 : class DataVector;
20 : /// \endcond
21 :
22 : namespace Xcts {
23 :
24 : /// Indicates a subset of the XCTS equations
25 1 : enum class Equations {
26 : /// Only the Hamiltonian constraint, solved for \f$\psi\f$
27 : Hamiltonian,
28 : /// Both the Hamiltonian constraint and the lapse equation, solved for
29 : /// \f$\psi\f$ and \f$\alpha\psi\f$
30 : HamiltonianAndLapse,
31 : /// The full XCTS equations, solved for \f$\psi\f$, \f$\alpha\psi\f$ and
32 : /// \f$\beta_\mathrm{excess}\f$
33 : HamiltonianLapseAndShift
34 : };
35 :
36 : /// The fluxes \f$F^i\f$ for the first-order formulation of the XCTS equations.
37 : ///
38 : /// \see Xcts::FirstOrderSystem
39 : template <Equations EnabledEquations, Geometry ConformalGeometry>
40 1 : struct Fluxes;
41 :
42 : /// \cond
43 : template <>
44 : struct Fluxes<Equations::Hamiltonian, Geometry::FlatCartesian> {
45 : using argument_tags = tmpl::list<>;
46 : using volume_tags = tmpl::list<>;
47 : using const_global_cache_tags = tmpl::list<>;
48 : static constexpr bool is_trivial = true;
49 : static constexpr bool is_discontinuous = false;
50 : static void apply(
51 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
52 : const Scalar<DataVector>& conformal_factor_minus_one,
53 : const tnsr::i<DataVector, 3>& conformal_factor_gradient);
54 : static void apply(
55 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
56 : const tnsr::i<DataVector, 3>& face_normal,
57 : const tnsr::I<DataVector, 3>& face_normal_vector,
58 : const Scalar<DataVector>& conformal_factor_minus_one);
59 : };
60 :
61 : template <>
62 : struct Fluxes<Equations::Hamiltonian, Geometry::Curved> {
63 : using argument_tags =
64 : tmpl::list<Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>>;
65 : using volume_tags = tmpl::list<>;
66 : using const_global_cache_tags = tmpl::list<>;
67 : static constexpr bool is_trivial = true;
68 : static constexpr bool is_discontinuous = false;
69 : static void apply(
70 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
71 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
72 : const Scalar<DataVector>& conformal_factor_minus_one,
73 : const tnsr::i<DataVector, 3>& conformal_factor_gradient);
74 : static void apply(
75 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
76 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
77 : const tnsr::i<DataVector, 3>& face_normal,
78 : const tnsr::I<DataVector, 3>& face_normal_vector,
79 : const Scalar<DataVector>& conformal_factor_minus_one);
80 : };
81 :
82 : template <>
83 : struct Fluxes<Equations::HamiltonianAndLapse, Geometry::FlatCartesian> {
84 : using argument_tags = tmpl::list<>;
85 : using volume_tags = tmpl::list<>;
86 : using const_global_cache_tags = tmpl::list<>;
87 : static constexpr bool is_trivial = true;
88 : static constexpr bool is_discontinuous = false;
89 : static void apply(
90 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
91 : gsl::not_null<tnsr::I<DataVector, 3>*>
92 : flux_for_lapse_times_conformal_factor,
93 : const Scalar<DataVector>& conformal_factor_minus_one,
94 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
95 : const tnsr::i<DataVector, 3>& conformal_factor_gradient,
96 : const tnsr::i<DataVector, 3>& lapse_times_conformal_factor_gradient);
97 : static void apply(
98 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
99 : gsl::not_null<tnsr::I<DataVector, 3>*>
100 : flux_for_lapse_times_conformal_factor,
101 : const tnsr::i<DataVector, 3>& face_normal,
102 : const tnsr::I<DataVector, 3>& face_normal_vector,
103 : const Scalar<DataVector>& conformal_factor_minus_one,
104 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one);
105 : };
106 :
107 : template <>
108 : struct Fluxes<Equations::HamiltonianAndLapse, Geometry::Curved> {
109 : using argument_tags =
110 : tmpl::list<Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>>;
111 : using volume_tags = tmpl::list<>;
112 : using const_global_cache_tags = tmpl::list<>;
113 : static constexpr bool is_trivial = true;
114 : static constexpr bool is_discontinuous = false;
115 : static void apply(
116 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
117 : gsl::not_null<tnsr::I<DataVector, 3>*>
118 : flux_for_lapse_times_conformal_factor,
119 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
120 : const Scalar<DataVector>& conformal_factor_minus_one,
121 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
122 : const tnsr::i<DataVector, 3>& conformal_factor_gradient,
123 : const tnsr::i<DataVector, 3>& lapse_times_conformal_factor_gradient);
124 : static void apply(
125 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
126 : gsl::not_null<tnsr::I<DataVector, 3>*>
127 : flux_for_lapse_times_conformal_factor,
128 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
129 : const tnsr::i<DataVector, 3>& face_normal,
130 : const tnsr::I<DataVector, 3>& face_normal_vector,
131 : const Scalar<DataVector>& conformal_factor_minus_one,
132 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one);
133 : };
134 :
135 : template <>
136 : struct Fluxes<Equations::HamiltonianLapseAndShift, Geometry::FlatCartesian> {
137 : using argument_tags = tmpl::list<>;
138 : using volume_tags = tmpl::list<>;
139 : using const_global_cache_tags = tmpl::list<>;
140 : static constexpr bool is_trivial = false;
141 : static constexpr bool is_discontinuous = false;
142 : static void apply(
143 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
144 : gsl::not_null<tnsr::I<DataVector, 3>*>
145 : flux_for_lapse_times_conformal_factor,
146 : gsl::not_null<tnsr::II<DataVector, 3>*> longitudinal_shift_excess,
147 : const Scalar<DataVector>& conformal_factor_minus_one,
148 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
149 : const tnsr::I<DataVector, 3>& shift_excess,
150 : const tnsr::i<DataVector, 3>& conformal_factor_gradient,
151 : const tnsr::i<DataVector, 3>& lapse_times_conformal_factor_gradient,
152 : const tnsr::iJ<DataVector, 3>& deriv_shift_excess);
153 : static void apply(
154 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
155 : gsl::not_null<tnsr::I<DataVector, 3>*>
156 : flux_for_lapse_times_conformal_factor,
157 : gsl::not_null<tnsr::II<DataVector, 3>*> longitudinal_shift_excess,
158 : const tnsr::i<DataVector, 3>& face_normal,
159 : const tnsr::I<DataVector, 3>& face_normal_vector,
160 : const Scalar<DataVector>& conformal_factor_minus_one,
161 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
162 : const tnsr::I<DataVector, 3>& shift_excess);
163 : };
164 :
165 : template <>
166 : struct Fluxes<Equations::HamiltonianLapseAndShift, Geometry::Curved> {
167 : using argument_tags = tmpl::list<
168 : Tags::ConformalMetric<DataVector, 3, Frame::Inertial>,
169 : Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>,
170 : Tags::ConformalChristoffelSecondKind<DataVector, 3, Frame::Inertial>>;
171 : using volume_tags = tmpl::list<>;
172 : using const_global_cache_tags = tmpl::list<>;
173 : static constexpr bool is_trivial = false;
174 : static constexpr bool is_discontinuous = false;
175 : static void apply(
176 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
177 : gsl::not_null<tnsr::I<DataVector, 3>*>
178 : flux_for_lapse_times_conformal_factor,
179 : gsl::not_null<tnsr::II<DataVector, 3>*> longitudinal_shift_excess,
180 : const tnsr::ii<DataVector, 3>& /*conformal_metric*/,
181 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
182 : const tnsr::Ijj<DataVector, 3>& christoffel_second_kind,
183 : const Scalar<DataVector>& conformal_factor_minus_one,
184 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
185 : const tnsr::I<DataVector, 3>& shift_excess,
186 : const tnsr::i<DataVector, 3>& conformal_factor_gradient,
187 : const tnsr::i<DataVector, 3>& lapse_times_conformal_factor_gradient,
188 : const tnsr::iJ<DataVector, 3>& deriv_shift_excess);
189 : static void apply(
190 : gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
191 : gsl::not_null<tnsr::I<DataVector, 3>*>
192 : flux_for_lapse_times_conformal_factor,
193 : gsl::not_null<tnsr::II<DataVector, 3>*> longitudinal_shift_excess,
194 : const tnsr::ii<DataVector, 3>& /*conformal_metric*/,
195 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
196 : const tnsr::Ijj<DataVector, 3>& christoffel_second_kind,
197 : const tnsr::i<DataVector, 3>& face_normal,
198 : const tnsr::I<DataVector, 3>& face_normal_vector,
199 : const Scalar<DataVector>& conformal_factor_minus_one,
200 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
201 : const tnsr::I<DataVector, 3>& shift_excess);
202 : };
203 : /// \endcond
204 :
205 : /// The sources \f$S\f$ for the first-order formulation of the XCTS equations.
206 : ///
207 : /// \see Xcts::FirstOrderSystem
208 : template <Equations EnabledEquations, Geometry ConformalGeometry,
209 : int ConformalMatterScale>
210 1 : struct Sources;
211 :
212 : /// \cond
213 : template <int ConformalMatterScale>
214 : struct Sources<Equations::Hamiltonian, Geometry::FlatCartesian,
215 : ConformalMatterScale> {
216 : using argument_tags = tmpl::list<
217 : gr::Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
218 : ConformalMatterScale>,
219 : gr::Tags::TraceExtrinsicCurvature<DataVector>,
220 : Tags::LongitudinalShiftMinusDtConformalMetricOverLapseSquare<DataVector>>;
221 : static void apply(
222 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
223 : const Scalar<DataVector>& conformal_energy_density,
224 : const Scalar<DataVector>& extrinsic_curvature_trace,
225 : const Scalar<DataVector>&
226 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
227 : const Scalar<DataVector>& conformal_factor_minus_one,
228 : const tnsr::I<DataVector, 3>& conformal_factor_flux);
229 : };
230 :
231 : template <int ConformalMatterScale>
232 : struct Sources<Equations::Hamiltonian, Geometry::Curved, ConformalMatterScale> {
233 : using argument_tags = tmpl::push_back<
234 : typename Sources<Equations::Hamiltonian, Geometry::FlatCartesian,
235 : ConformalMatterScale>::argument_tags,
236 : Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
237 : Tags::ConformalRicciScalar<DataVector>>;
238 : using const_global_cache_tags = tmpl::list<>;
239 : static void apply(
240 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
241 : const Scalar<DataVector>& conformal_energy_density,
242 : const Scalar<DataVector>& extrinsic_curvature_trace,
243 : const Scalar<DataVector>&
244 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
245 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
246 : const Scalar<DataVector>& conformal_ricci_scalar,
247 : const Scalar<DataVector>& conformal_factor_minus_one,
248 : const tnsr::I<DataVector, 3>& conformal_factor_flux);
249 : };
250 :
251 : template <int ConformalMatterScale>
252 : struct Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
253 : ConformalMatterScale> {
254 : using argument_tags = tmpl::list<
255 : gr::Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
256 : ConformalMatterScale>,
257 : gr::Tags::Conformal<gr::Tags::StressTrace<DataVector>,
258 : ConformalMatterScale>,
259 : gr::Tags::TraceExtrinsicCurvature<DataVector>,
260 : ::Tags::dt<gr::Tags::TraceExtrinsicCurvature<DataVector>>,
261 : Tags::LongitudinalShiftMinusDtConformalMetricSquare<DataVector>,
262 : Tags::ShiftDotDerivExtrinsicCurvatureTrace<DataVector>>;
263 : using const_global_cache_tags = tmpl::list<>;
264 : static void apply(
265 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
266 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
267 : const Scalar<DataVector>& conformal_energy_density,
268 : const Scalar<DataVector>& conformal_stress_trace,
269 : const Scalar<DataVector>& extrinsic_curvature_trace,
270 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
271 : const Scalar<DataVector>&
272 : longitudinal_shift_minus_dt_conformal_metric_square,
273 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
274 : const Scalar<DataVector>& conformal_factor_minus_one,
275 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
276 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
277 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux);
278 : };
279 :
280 : template <int ConformalMatterScale>
281 : struct Sources<Equations::HamiltonianAndLapse, Geometry::Curved,
282 : ConformalMatterScale> {
283 : using argument_tags = tmpl::push_back<
284 : typename Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
285 : ConformalMatterScale>::argument_tags,
286 : Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
287 : Tags::ConformalRicciScalar<DataVector>>;
288 : using const_global_cache_tags = tmpl::list<>;
289 : static void apply(
290 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
291 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
292 : const Scalar<DataVector>& conformal_energy_density,
293 : const Scalar<DataVector>& conformal_stress_trace,
294 : const Scalar<DataVector>& extrinsic_curvature_trace,
295 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
296 : const Scalar<DataVector>&
297 : longitudinal_shift_minus_dt_conformal_metric_square,
298 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
299 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
300 : const Scalar<DataVector>& conformal_ricci_scalar,
301 : const Scalar<DataVector>& conformal_factor_minus_one,
302 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
303 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
304 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux);
305 : };
306 :
307 : template <int ConformalMatterScale>
308 : struct Sources<Equations::HamiltonianLapseAndShift, Geometry::FlatCartesian,
309 : ConformalMatterScale> {
310 : using argument_tags = tmpl::list<
311 : gr::Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
312 : ConformalMatterScale>,
313 : gr::Tags::Conformal<gr::Tags::StressTrace<DataVector>,
314 : ConformalMatterScale>,
315 : gr::Tags::Conformal<gr::Tags::MomentumDensity<DataVector, 3>,
316 : ConformalMatterScale>,
317 : gr::Tags::TraceExtrinsicCurvature<DataVector>,
318 : ::Tags::dt<gr::Tags::TraceExtrinsicCurvature<DataVector>>,
319 : ::Tags::deriv<gr::Tags::TraceExtrinsicCurvature<DataVector>,
320 : tmpl::size_t<3>, Frame::Inertial>,
321 : Tags::ShiftBackground<DataVector, 3, Frame::Inertial>,
322 : Tags::LongitudinalShiftBackgroundMinusDtConformalMetric<DataVector, 3,
323 : Frame::Inertial>,
324 : ::Tags::div<Tags::LongitudinalShiftBackgroundMinusDtConformalMetric<
325 : DataVector, 3, Frame::Inertial>>>;
326 : using const_global_cache_tags = tmpl::list<>;
327 : static void apply(
328 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
329 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
330 : gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
331 : const Scalar<DataVector>& conformal_energy_density,
332 : const Scalar<DataVector>& conformal_stress_trace,
333 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
334 : const Scalar<DataVector>& extrinsic_curvature_trace,
335 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
336 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
337 : const tnsr::I<DataVector, 3>& shift_background,
338 : const tnsr::II<DataVector, 3>&
339 : longitudinal_shift_background_minus_dt_conformal_metric,
340 : const tnsr::I<DataVector, 3>&
341 : div_longitudinal_shift_background_minus_dt_conformal_metric,
342 : const Scalar<DataVector>& conformal_factor_minus_one,
343 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
344 : const tnsr::I<DataVector, 3>& shift_excess,
345 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
346 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
347 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess);
348 : };
349 :
350 : template <int ConformalMatterScale>
351 : struct Sources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
352 : ConformalMatterScale> {
353 : using argument_tags = tmpl::push_back<
354 : typename Sources<Equations::HamiltonianLapseAndShift,
355 : Geometry::FlatCartesian,
356 : ConformalMatterScale>::argument_tags,
357 : Tags::ConformalMetric<DataVector, 3, Frame::Inertial>,
358 : Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>,
359 : Tags::ConformalChristoffelFirstKind<DataVector, 3, Frame::Inertial>,
360 : Tags::ConformalChristoffelSecondKind<DataVector, 3, Frame::Inertial>,
361 : Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
362 : Tags::ConformalRicciScalar<DataVector>>;
363 : using const_global_cache_tags = tmpl::list<>;
364 : static void apply(
365 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
366 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
367 : gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
368 : const Scalar<DataVector>& conformal_energy_density,
369 : const Scalar<DataVector>& conformal_stress_trace,
370 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
371 : const Scalar<DataVector>& extrinsic_curvature_trace,
372 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
373 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
374 : const tnsr::I<DataVector, 3>& shift_background,
375 : const tnsr::II<DataVector, 3>&
376 : longitudinal_shift_background_minus_dt_conformal_metric,
377 : const tnsr::I<DataVector, 3>&
378 : div_longitudinal_shift_background_minus_dt_conformal_metric,
379 : const tnsr::ii<DataVector, 3>& conformal_metric,
380 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
381 : const tnsr::ijj<DataVector, 3>& /*conformal_christoffel_first_kind*/,
382 : const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
383 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
384 : const Scalar<DataVector>& conformal_ricci_scalar,
385 : const Scalar<DataVector>& conformal_factor_minus_one,
386 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
387 : const tnsr::I<DataVector, 3>& shift_excess,
388 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
389 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
390 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess);
391 : };
392 : /// \endcond
393 :
394 : /// The linearization of the sources \f$S\f$ for the first-order formulation of
395 : /// the XCTS equations.
396 : ///
397 : /// \see Xcts::FirstOrderSystem
398 : template <Equations EnabledEquations, Geometry ConformalGeometry,
399 : int ConformalMatterScale>
400 1 : struct LinearizedSources;
401 :
402 : /// \cond
403 : template <int ConformalMatterScale>
404 : struct LinearizedSources<Equations::Hamiltonian, Geometry::FlatCartesian,
405 : ConformalMatterScale> {
406 : using argument_tags = tmpl::push_back<
407 : typename Sources<Equations::Hamiltonian, Geometry::FlatCartesian,
408 : ConformalMatterScale>::argument_tags,
409 : Tags::ConformalFactorMinusOne<DataVector>>;
410 : using const_global_cache_tags = tmpl::list<>;
411 : static void apply(
412 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
413 : const Scalar<DataVector>& conformal_energy_density,
414 : const Scalar<DataVector>& extrinsic_curvature_trace,
415 : const Scalar<DataVector>&
416 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
417 : const Scalar<DataVector>& conformal_factor_minus_one,
418 : const Scalar<DataVector>& conformal_factor_correction,
419 : const tnsr::I<DataVector, 3>&
420 : /*conformal_factor_flux_correction*/);
421 : };
422 :
423 : template <int ConformalMatterScale>
424 : struct LinearizedSources<Equations::Hamiltonian, Geometry::Curved,
425 : ConformalMatterScale> {
426 : using argument_tags =
427 : tmpl::push_back<typename Sources<Equations::Hamiltonian, Geometry::Curved,
428 : ConformalMatterScale>::argument_tags,
429 : Tags::ConformalFactorMinusOne<DataVector>>;
430 : using const_global_cache_tags = tmpl::list<>;
431 : static void apply(
432 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
433 : const Scalar<DataVector>& conformal_energy_density,
434 : const Scalar<DataVector>& extrinsic_curvature_trace,
435 : const Scalar<DataVector>&
436 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
437 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
438 : const Scalar<DataVector>& conformal_ricci_scalar,
439 : const Scalar<DataVector>& conformal_factor_minus_one,
440 : const Scalar<DataVector>& conformal_factor_correction,
441 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction);
442 : };
443 :
444 : template <int ConformalMatterScale>
445 : struct LinearizedSources<Equations::HamiltonianAndLapse,
446 : Geometry::FlatCartesian, ConformalMatterScale> {
447 : using argument_tags = tmpl::push_back<
448 : typename Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
449 : ConformalMatterScale>::argument_tags,
450 : Tags::ConformalFactorMinusOne<DataVector>,
451 : Tags::LapseTimesConformalFactorMinusOne<DataVector>>;
452 : using const_global_cache_tags = tmpl::list<>;
453 : static void apply(
454 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
455 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
456 : const Scalar<DataVector>& conformal_energy_density,
457 : const Scalar<DataVector>& conformal_stress_trace,
458 : const Scalar<DataVector>& extrinsic_curvature_trace,
459 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
460 : const Scalar<DataVector>&
461 : longitudinal_shift_minus_dt_conformal_metric_square,
462 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
463 : const Scalar<DataVector>& conformal_factor_minus_one,
464 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
465 : const Scalar<DataVector>& conformal_factor_correction,
466 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
467 : const tnsr::I<DataVector, 3>& /*conformal_factor_flux_correction*/,
468 : const tnsr::I<DataVector, 3>&
469 : lapse_times_conformal_factor_flux_correction);
470 : };
471 :
472 : template <int ConformalMatterScale>
473 : struct LinearizedSources<Equations::HamiltonianAndLapse, Geometry::Curved,
474 : ConformalMatterScale> {
475 : using argument_tags = tmpl::push_back<
476 : typename Sources<Equations::HamiltonianAndLapse, Geometry::Curved,
477 : ConformalMatterScale>::argument_tags,
478 : Tags::ConformalFactorMinusOne<DataVector>,
479 : Tags::LapseTimesConformalFactorMinusOne<DataVector>>;
480 : using const_global_cache_tags = tmpl::list<>;
481 : static void apply(
482 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
483 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
484 : const Scalar<DataVector>& conformal_energy_density,
485 : const Scalar<DataVector>& conformal_stress_trace,
486 : const Scalar<DataVector>& extrinsic_curvature_trace,
487 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
488 : const Scalar<DataVector>&
489 : longitudinal_shift_minus_dt_conformal_metric_square,
490 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
491 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
492 : const Scalar<DataVector>& conformal_ricci_scalar,
493 : const Scalar<DataVector>& conformal_factor_minus_one,
494 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
495 : const Scalar<DataVector>& conformal_factor_correction,
496 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
497 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
498 : const tnsr::I<DataVector, 3>&
499 : lapse_times_conformal_factor_flux_correction);
500 : };
501 :
502 : template <int ConformalMatterScale>
503 : struct LinearizedSources<Equations::HamiltonianLapseAndShift,
504 : Geometry::FlatCartesian, ConformalMatterScale> {
505 : using argument_tags = tmpl::push_back<
506 : typename Sources<Equations::HamiltonianLapseAndShift,
507 : Geometry::FlatCartesian,
508 : ConformalMatterScale>::argument_tags,
509 : Tags::ConformalFactorMinusOne<DataVector>,
510 : Tags::LapseTimesConformalFactorMinusOne<DataVector>,
511 : Tags::ShiftExcess<DataVector, 3, Frame::Inertial>,
512 : ::Tags::Flux<Tags::ConformalFactorMinusOne<DataVector>, tmpl::size_t<3>,
513 : Frame::Inertial>,
514 : ::Tags::Flux<Tags::LapseTimesConformalFactorMinusOne<DataVector>,
515 : tmpl::size_t<3>, Frame::Inertial>,
516 : Tags::LongitudinalShiftExcess<DataVector, 3, Frame::Inertial>>;
517 : using const_global_cache_tags = tmpl::list<>;
518 : static void apply(
519 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
520 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
521 : gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
522 : const Scalar<DataVector>& conformal_energy_density,
523 : const Scalar<DataVector>& conformal_stress_trace,
524 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
525 : const Scalar<DataVector>& extrinsic_curvature_trace,
526 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
527 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
528 : const tnsr::I<DataVector, 3>& shift_background,
529 : const tnsr::II<DataVector, 3>&
530 : longitudinal_shift_background_minus_dt_conformal_metric,
531 : const tnsr::I<DataVector, 3>&
532 : div_longitudinal_shift_background_minus_dt_conformal_metric,
533 : const Scalar<DataVector>& conformal_factor_minus_one,
534 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
535 : const tnsr::I<DataVector, 3>& shift_excess,
536 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
537 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
538 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
539 : const Scalar<DataVector>& conformal_factor_correction,
540 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
541 : const tnsr::I<DataVector, 3>& shift_excess_correction,
542 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
543 : const tnsr::I<DataVector, 3>&
544 : lapse_times_conformal_factor_flux_correction,
545 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess_correction);
546 : };
547 :
548 : template <int ConformalMatterScale>
549 : struct LinearizedSources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
550 : ConformalMatterScale> {
551 : using argument_tags = tmpl::push_back<
552 : typename Sources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
553 : ConformalMatterScale>::argument_tags,
554 : Tags::ConformalFactorMinusOne<DataVector>,
555 : Tags::LapseTimesConformalFactorMinusOne<DataVector>,
556 : Tags::ShiftExcess<DataVector, 3, Frame::Inertial>,
557 : ::Tags::Flux<Tags::ConformalFactorMinusOne<DataVector>, tmpl::size_t<3>,
558 : Frame::Inertial>,
559 : ::Tags::Flux<Tags::LapseTimesConformalFactorMinusOne<DataVector>,
560 : tmpl::size_t<3>, Frame::Inertial>,
561 : Tags::LongitudinalShiftExcess<DataVector, 3, Frame::Inertial>>;
562 : using const_global_cache_tags = tmpl::list<>;
563 : static void apply(
564 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
565 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
566 : gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
567 : const Scalar<DataVector>& conformal_energy_density,
568 : const Scalar<DataVector>& conformal_stress_trace,
569 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
570 : const Scalar<DataVector>& extrinsic_curvature_trace,
571 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
572 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
573 : const tnsr::I<DataVector, 3>& shift_background,
574 : const tnsr::II<DataVector, 3>&
575 : longitudinal_shift_background_minus_dt_conformal_metric,
576 : const tnsr::I<DataVector, 3>&
577 : /*div_longitudinal_shift_background_minus_dt_conformal_metric*/,
578 : const tnsr::ii<DataVector, 3>& conformal_metric,
579 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
580 : const tnsr::ijj<DataVector, 3>& /*conformal_christoffel_first_kind*/,
581 : const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
582 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
583 : const Scalar<DataVector>& conformal_ricci_scalar,
584 : const Scalar<DataVector>& conformal_factor_minus_one,
585 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
586 : const tnsr::I<DataVector, 3>& shift_excess,
587 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
588 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
589 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
590 : const Scalar<DataVector>& conformal_factor_correction,
591 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
592 : const tnsr::I<DataVector, 3>& shift_excess_correction,
593 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
594 : const tnsr::I<DataVector, 3>&
595 : lapse_times_conformal_factor_flux_correction,
596 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess_correction);
597 : };
598 : /// \endcond
599 :
600 : } // namespace Xcts
|