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_gradient*/,
229 : const tnsr::I<DataVector, 3>& conformal_factor_flux);
230 : };
231 :
232 : template <int ConformalMatterScale>
233 : struct Sources<Equations::Hamiltonian, Geometry::Curved, ConformalMatterScale> {
234 : using argument_tags = tmpl::push_back<
235 : typename Sources<Equations::Hamiltonian, Geometry::FlatCartesian,
236 : ConformalMatterScale>::argument_tags,
237 : Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
238 : Tags::ConformalRicciScalar<DataVector>>;
239 : using const_global_cache_tags = tmpl::list<>;
240 : static void apply(
241 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
242 : const Scalar<DataVector>& conformal_energy_density,
243 : const Scalar<DataVector>& extrinsic_curvature_trace,
244 : const Scalar<DataVector>&
245 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
246 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
247 : const Scalar<DataVector>& conformal_ricci_scalar,
248 : const Scalar<DataVector>& conformal_factor_minus_one,
249 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient*/,
250 : const tnsr::I<DataVector, 3>& conformal_factor_flux);
251 : };
252 :
253 : template <int ConformalMatterScale>
254 : struct Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
255 : ConformalMatterScale> {
256 : using argument_tags = tmpl::list<
257 : gr::Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
258 : ConformalMatterScale>,
259 : gr::Tags::Conformal<gr::Tags::StressTrace<DataVector>,
260 : ConformalMatterScale>,
261 : gr::Tags::TraceExtrinsicCurvature<DataVector>,
262 : ::Tags::dt<gr::Tags::TraceExtrinsicCurvature<DataVector>>,
263 : Tags::LongitudinalShiftMinusDtConformalMetricSquare<DataVector>,
264 : Tags::ShiftDotDerivExtrinsicCurvatureTrace<DataVector>>;
265 : using const_global_cache_tags = tmpl::list<>;
266 : static void apply(
267 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
268 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
269 : const Scalar<DataVector>& conformal_energy_density,
270 : const Scalar<DataVector>& conformal_stress_trace,
271 : const Scalar<DataVector>& extrinsic_curvature_trace,
272 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
273 : const Scalar<DataVector>&
274 : longitudinal_shift_minus_dt_conformal_metric_square,
275 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
276 : const Scalar<DataVector>& conformal_factor_minus_one,
277 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
278 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient*/,
279 : const tnsr::i<DataVector, 3>&
280 : /*lapse_times_conformal_factor_gradient*/,
281 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
282 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux);
283 : };
284 :
285 : template <int ConformalMatterScale>
286 : struct Sources<Equations::HamiltonianAndLapse, Geometry::Curved,
287 : ConformalMatterScale> {
288 : using argument_tags = tmpl::push_back<
289 : typename Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
290 : ConformalMatterScale>::argument_tags,
291 : Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
292 : Tags::ConformalRicciScalar<DataVector>>;
293 : using const_global_cache_tags = tmpl::list<>;
294 : static void apply(
295 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
296 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
297 : const Scalar<DataVector>& conformal_energy_density,
298 : const Scalar<DataVector>& conformal_stress_trace,
299 : const Scalar<DataVector>& extrinsic_curvature_trace,
300 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
301 : const Scalar<DataVector>&
302 : longitudinal_shift_minus_dt_conformal_metric_square,
303 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
304 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
305 : const Scalar<DataVector>& conformal_ricci_scalar,
306 : const Scalar<DataVector>& conformal_factor_minus_one,
307 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
308 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient*/,
309 : const tnsr::i<DataVector, 3>&
310 : /*lapse_times_conformal_factor_gradient*/,
311 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
312 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux);
313 : };
314 :
315 : template <int ConformalMatterScale>
316 : struct Sources<Equations::HamiltonianLapseAndShift, Geometry::FlatCartesian,
317 : ConformalMatterScale> {
318 : using argument_tags = tmpl::list<
319 : gr::Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
320 : ConformalMatterScale>,
321 : gr::Tags::Conformal<gr::Tags::StressTrace<DataVector>,
322 : ConformalMatterScale>,
323 : gr::Tags::Conformal<gr::Tags::MomentumDensity<DataVector, 3>,
324 : ConformalMatterScale>,
325 : gr::Tags::TraceExtrinsicCurvature<DataVector>,
326 : ::Tags::dt<gr::Tags::TraceExtrinsicCurvature<DataVector>>,
327 : ::Tags::deriv<gr::Tags::TraceExtrinsicCurvature<DataVector>,
328 : tmpl::size_t<3>, Frame::Inertial>,
329 : Tags::ShiftBackground<DataVector, 3, Frame::Inertial>,
330 : Tags::LongitudinalShiftBackgroundMinusDtConformalMetric<DataVector, 3,
331 : Frame::Inertial>,
332 : ::Tags::div<Tags::LongitudinalShiftBackgroundMinusDtConformalMetric<
333 : DataVector, 3, Frame::Inertial>>>;
334 : using const_global_cache_tags = tmpl::list<>;
335 : static void apply(
336 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
337 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
338 : gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
339 : const Scalar<DataVector>& conformal_energy_density,
340 : const Scalar<DataVector>& conformal_stress_trace,
341 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
342 : const Scalar<DataVector>& extrinsic_curvature_trace,
343 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
344 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
345 : const tnsr::I<DataVector, 3>& shift_background,
346 : const tnsr::II<DataVector, 3>&
347 : longitudinal_shift_background_minus_dt_conformal_metric,
348 : const tnsr::I<DataVector, 3>&
349 : div_longitudinal_shift_background_minus_dt_conformal_metric,
350 : const Scalar<DataVector>& conformal_factor_minus_one,
351 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
352 : const tnsr::I<DataVector, 3>& shift_excess,
353 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient*/,
354 : const tnsr::i<DataVector, 3>&
355 : /*lapse_times_conformal_factor_gradient*/,
356 : const tnsr::iJ<DataVector, 3>& /*deriv_shift_excess*/,
357 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
358 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
359 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess);
360 : };
361 :
362 : template <int ConformalMatterScale>
363 : struct Sources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
364 : ConformalMatterScale> {
365 : using argument_tags = tmpl::push_back<
366 : typename Sources<Equations::HamiltonianLapseAndShift,
367 : Geometry::FlatCartesian,
368 : ConformalMatterScale>::argument_tags,
369 : Tags::ConformalMetric<DataVector, 3, Frame::Inertial>,
370 : Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>,
371 : Tags::ConformalChristoffelFirstKind<DataVector, 3, Frame::Inertial>,
372 : Tags::ConformalChristoffelSecondKind<DataVector, 3, Frame::Inertial>,
373 : Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
374 : Tags::ConformalRicciScalar<DataVector>>;
375 : using const_global_cache_tags = tmpl::list<>;
376 : static void apply(
377 : gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
378 : gsl::not_null<Scalar<DataVector>*> lapse_equation,
379 : gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
380 : const Scalar<DataVector>& conformal_energy_density,
381 : const Scalar<DataVector>& conformal_stress_trace,
382 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
383 : const Scalar<DataVector>& extrinsic_curvature_trace,
384 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
385 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
386 : const tnsr::I<DataVector, 3>& shift_background,
387 : const tnsr::II<DataVector, 3>&
388 : longitudinal_shift_background_minus_dt_conformal_metric,
389 : const tnsr::I<DataVector, 3>&
390 : div_longitudinal_shift_background_minus_dt_conformal_metric,
391 : const tnsr::ii<DataVector, 3>& conformal_metric,
392 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
393 : const tnsr::ijj<DataVector, 3>& /*conformal_christoffel_first_kind*/,
394 : const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
395 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
396 : const Scalar<DataVector>& conformal_ricci_scalar,
397 : const Scalar<DataVector>& conformal_factor_minus_one,
398 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
399 : const tnsr::I<DataVector, 3>& shift_excess,
400 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient*/,
401 : const tnsr::i<DataVector, 3>&
402 : /*lapse_times_conformal_factor_gradient*/,
403 : const tnsr::iJ<DataVector, 3>& /*deriv_shift_excess*/,
404 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
405 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
406 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess);
407 : };
408 : /// \endcond
409 :
410 : /// The linearization of the sources \f$S\f$ for the first-order formulation of
411 : /// the XCTS equations.
412 : ///
413 : /// \see Xcts::FirstOrderSystem
414 : template <Equations EnabledEquations, Geometry ConformalGeometry,
415 : int ConformalMatterScale>
416 1 : struct LinearizedSources;
417 :
418 : /// \cond
419 : template <int ConformalMatterScale>
420 : struct LinearizedSources<Equations::Hamiltonian, Geometry::FlatCartesian,
421 : ConformalMatterScale> {
422 : using argument_tags = tmpl::push_back<
423 : typename Sources<Equations::Hamiltonian, Geometry::FlatCartesian,
424 : ConformalMatterScale>::argument_tags,
425 : Tags::ConformalFactorMinusOne<DataVector>>;
426 : using const_global_cache_tags = tmpl::list<>;
427 : static void apply(
428 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
429 : const Scalar<DataVector>& conformal_energy_density,
430 : const Scalar<DataVector>& extrinsic_curvature_trace,
431 : const Scalar<DataVector>&
432 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
433 : const Scalar<DataVector>& conformal_factor_minus_one,
434 : const Scalar<DataVector>& conformal_factor_correction,
435 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient_correction*/,
436 : const tnsr::I<DataVector, 3>&
437 : /*conformal_factor_flux_correction*/);
438 : };
439 :
440 : template <int ConformalMatterScale>
441 : struct LinearizedSources<Equations::Hamiltonian, Geometry::Curved,
442 : ConformalMatterScale> {
443 : using argument_tags =
444 : tmpl::push_back<typename Sources<Equations::Hamiltonian, Geometry::Curved,
445 : ConformalMatterScale>::argument_tags,
446 : Tags::ConformalFactorMinusOne<DataVector>>;
447 : using const_global_cache_tags = tmpl::list<>;
448 : static void apply(
449 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
450 : const Scalar<DataVector>& conformal_energy_density,
451 : const Scalar<DataVector>& extrinsic_curvature_trace,
452 : const Scalar<DataVector>&
453 : longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
454 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
455 : const Scalar<DataVector>& conformal_ricci_scalar,
456 : const Scalar<DataVector>& conformal_factor_minus_one,
457 : const Scalar<DataVector>& conformal_factor_correction,
458 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient_correction*/,
459 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction);
460 : };
461 :
462 : template <int ConformalMatterScale>
463 : struct LinearizedSources<Equations::HamiltonianAndLapse,
464 : Geometry::FlatCartesian, ConformalMatterScale> {
465 : using argument_tags = tmpl::push_back<
466 : typename Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
467 : ConformalMatterScale>::argument_tags,
468 : Tags::ConformalFactorMinusOne<DataVector>,
469 : Tags::LapseTimesConformalFactorMinusOne<DataVector>>;
470 : using const_global_cache_tags = tmpl::list<>;
471 : static void apply(
472 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
473 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
474 : const Scalar<DataVector>& conformal_energy_density,
475 : const Scalar<DataVector>& conformal_stress_trace,
476 : const Scalar<DataVector>& extrinsic_curvature_trace,
477 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
478 : const Scalar<DataVector>&
479 : longitudinal_shift_minus_dt_conformal_metric_square,
480 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
481 : const Scalar<DataVector>& conformal_factor_minus_one,
482 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
483 : const Scalar<DataVector>& conformal_factor_correction,
484 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
485 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient_correction*/,
486 : const tnsr::i<DataVector, 3>&
487 : /*lapse_times_conformal_factor_gradient_correction*/,
488 : const tnsr::I<DataVector, 3>& /*conformal_factor_flux_correction*/,
489 : const tnsr::I<DataVector, 3>&
490 : lapse_times_conformal_factor_flux_correction);
491 : };
492 :
493 : template <int ConformalMatterScale>
494 : struct LinearizedSources<Equations::HamiltonianAndLapse, Geometry::Curved,
495 : ConformalMatterScale> {
496 : using argument_tags = tmpl::push_back<
497 : typename Sources<Equations::HamiltonianAndLapse, Geometry::Curved,
498 : ConformalMatterScale>::argument_tags,
499 : Tags::ConformalFactorMinusOne<DataVector>,
500 : Tags::LapseTimesConformalFactorMinusOne<DataVector>>;
501 : using const_global_cache_tags = tmpl::list<>;
502 : static void apply(
503 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
504 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
505 : const Scalar<DataVector>& conformal_energy_density,
506 : const Scalar<DataVector>& conformal_stress_trace,
507 : const Scalar<DataVector>& extrinsic_curvature_trace,
508 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
509 : const Scalar<DataVector>&
510 : longitudinal_shift_minus_dt_conformal_metric_square,
511 : const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
512 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
513 : const Scalar<DataVector>& conformal_ricci_scalar,
514 : const Scalar<DataVector>& conformal_factor_minus_one,
515 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
516 : const Scalar<DataVector>& conformal_factor_correction,
517 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
518 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient_correction*/,
519 : const tnsr::i<DataVector, 3>&
520 : /*lapse_times_conformal_factor_gradient_correction*/,
521 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
522 : const tnsr::I<DataVector, 3>&
523 : lapse_times_conformal_factor_flux_correction);
524 : };
525 :
526 : template <int ConformalMatterScale>
527 : struct LinearizedSources<Equations::HamiltonianLapseAndShift,
528 : Geometry::FlatCartesian, ConformalMatterScale> {
529 : using argument_tags = tmpl::push_back<
530 : typename Sources<Equations::HamiltonianLapseAndShift,
531 : Geometry::FlatCartesian,
532 : ConformalMatterScale>::argument_tags,
533 : Tags::ConformalFactorMinusOne<DataVector>,
534 : Tags::LapseTimesConformalFactorMinusOne<DataVector>,
535 : Tags::ShiftExcess<DataVector, 3, Frame::Inertial>,
536 : ::Tags::Flux<Tags::ConformalFactorMinusOne<DataVector>, tmpl::size_t<3>,
537 : Frame::Inertial>,
538 : ::Tags::Flux<Tags::LapseTimesConformalFactorMinusOne<DataVector>,
539 : tmpl::size_t<3>, Frame::Inertial>,
540 : Tags::LongitudinalShiftExcess<DataVector, 3, Frame::Inertial>>;
541 : using const_global_cache_tags = tmpl::list<>;
542 : static void apply(
543 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
544 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
545 : gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
546 : const Scalar<DataVector>& conformal_energy_density,
547 : const Scalar<DataVector>& conformal_stress_trace,
548 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
549 : const Scalar<DataVector>& extrinsic_curvature_trace,
550 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
551 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
552 : const tnsr::I<DataVector, 3>& shift_background,
553 : const tnsr::II<DataVector, 3>&
554 : longitudinal_shift_background_minus_dt_conformal_metric,
555 : const tnsr::I<DataVector, 3>&
556 : div_longitudinal_shift_background_minus_dt_conformal_metric,
557 : const Scalar<DataVector>& conformal_factor_minus_one,
558 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
559 : const tnsr::I<DataVector, 3>& shift_excess,
560 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
561 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
562 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
563 : const Scalar<DataVector>& conformal_factor_correction,
564 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
565 : const tnsr::I<DataVector, 3>& shift_excess_correction,
566 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient_correction*/,
567 : const tnsr::i<DataVector, 3>&
568 : /*lapse_times_conformal_factor_gradient_correction*/,
569 : const tnsr::iJ<DataVector, 3>& /*deriv_shift_excess_correction*/,
570 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
571 : const tnsr::I<DataVector, 3>&
572 : lapse_times_conformal_factor_flux_correction,
573 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess_correction);
574 : };
575 :
576 : template <int ConformalMatterScale>
577 : struct LinearizedSources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
578 : ConformalMatterScale> {
579 : using argument_tags = tmpl::push_back<
580 : typename Sources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
581 : ConformalMatterScale>::argument_tags,
582 : Tags::ConformalFactorMinusOne<DataVector>,
583 : Tags::LapseTimesConformalFactorMinusOne<DataVector>,
584 : Tags::ShiftExcess<DataVector, 3, Frame::Inertial>,
585 : ::Tags::Flux<Tags::ConformalFactorMinusOne<DataVector>, tmpl::size_t<3>,
586 : Frame::Inertial>,
587 : ::Tags::Flux<Tags::LapseTimesConformalFactorMinusOne<DataVector>,
588 : tmpl::size_t<3>, Frame::Inertial>,
589 : Tags::LongitudinalShiftExcess<DataVector, 3, Frame::Inertial>>;
590 : using const_global_cache_tags = tmpl::list<>;
591 : static void apply(
592 : gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
593 : gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
594 : gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
595 : const Scalar<DataVector>& conformal_energy_density,
596 : const Scalar<DataVector>& conformal_stress_trace,
597 : const tnsr::I<DataVector, 3>& conformal_momentum_density,
598 : const Scalar<DataVector>& extrinsic_curvature_trace,
599 : const Scalar<DataVector>& dt_extrinsic_curvature_trace,
600 : const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
601 : const tnsr::I<DataVector, 3>& shift_background,
602 : const tnsr::II<DataVector, 3>&
603 : longitudinal_shift_background_minus_dt_conformal_metric,
604 : const tnsr::I<DataVector, 3>&
605 : /*div_longitudinal_shift_background_minus_dt_conformal_metric*/,
606 : const tnsr::ii<DataVector, 3>& conformal_metric,
607 : const tnsr::II<DataVector, 3>& inv_conformal_metric,
608 : const tnsr::ijj<DataVector, 3>& /*conformal_christoffel_first_kind*/,
609 : const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
610 : const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
611 : const Scalar<DataVector>& conformal_ricci_scalar,
612 : const Scalar<DataVector>& conformal_factor_minus_one,
613 : const Scalar<DataVector>& lapse_times_conformal_factor_minus_one,
614 : const tnsr::I<DataVector, 3>& shift_excess,
615 : const tnsr::I<DataVector, 3>& conformal_factor_flux,
616 : const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
617 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
618 : const Scalar<DataVector>& conformal_factor_correction,
619 : const Scalar<DataVector>& lapse_times_conformal_factor_correction,
620 : const tnsr::I<DataVector, 3>& shift_excess_correction,
621 : const tnsr::i<DataVector, 3>& /*conformal_factor_gradient_correction*/,
622 : const tnsr::i<DataVector, 3>&
623 : /*lapse_times_conformal_factor_gradient_correction*/,
624 : const tnsr::iJ<DataVector, 3>& /*deriv_shift_excess_correction*/,
625 : const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
626 : const tnsr::I<DataVector, 3>&
627 : lapse_times_conformal_factor_flux_correction,
628 : const tnsr::II<DataVector, 3>& longitudinal_shift_excess_correction);
629 : };
630 : /// \endcond
631 :
632 : } // namespace Xcts
|