FluxesAndSources.hpp
1 // Distributed under the MIT License.
2 // See LICENSE.txt for details.
3 
4 #pragma once
5 
6 #include <cstddef>
7 
9 #include "Elliptic/Systems/Xcts/Geometry.hpp"
10 #include "Elliptic/Systems/Xcts/Tags.hpp"
13 #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
14 #include "Utilities/Gsl.hpp"
15 #include "Utilities/TMPL.hpp"
16 
17 /// \cond
18 class DataVector;
19 /// \endcond
20 
21 namespace Xcts {
22 
23 /// Indicates a subset of the XCTS equations
24 enum class Equations {
25  /// Only the Hamiltonian constraint, solved for \f$\psi\f$
27  /// Both the Hamiltonian constraint and the lapse equation, solved for
28  /// \f$\psi\f$ and \f$\alpha\psi\f$
30  /// The full XCTS equations, solved for \f$\psi\f$, \f$\alpha\psi\f$ and
31  /// \f$\beta_\mathrm{excess}\f$
33 };
34 
35 /// The fluxes \f$F^i\f$ for the first-order formulation of the XCTS equations.
36 ///
37 /// \see Xcts::FirstOrderSystem
38 template <Equations EnabledEquations, Geometry ConformalGeometry>
39 struct Fluxes;
40 
41 /// \cond
42 template <>
44  using argument_tags = tmpl::list<>;
45  static void apply(
46  gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
47  const tnsr::i<DataVector, 3>& conformal_factor_gradient) noexcept;
48  static void apply(const gsl::not_null<tnsr::Ij<DataVector, 3>*>
49  flux_for_conformal_factor_gradient,
50  const Scalar<DataVector>& conformal_factor) noexcept;
51 };
52 
53 template <>
55  using argument_tags =
56  tmpl::list<Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>>;
57  static void apply(
58  gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
59  const tnsr::II<DataVector, 3>& inv_conformal_metric,
60  const tnsr::i<DataVector, 3>& conformal_factor_gradient) noexcept;
61  static void apply(gsl::not_null<tnsr::Ij<DataVector, 3>*>
62  flux_for_conformal_factor_gradient,
63  const tnsr::II<DataVector, 3>& inv_conformal_metric,
64  const Scalar<DataVector>& conformal_factor) noexcept;
65 };
66 
67 template <>
68 struct Fluxes<Equations::HamiltonianAndLapse, Geometry::FlatCartesian> {
69  using argument_tags = tmpl::list<>;
70  static void apply(
71  gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
72  gsl::not_null<tnsr::I<DataVector, 3>*>
73  flux_for_lapse_times_conformal_factor,
74  const tnsr::i<DataVector, 3>& conformal_factor_gradient,
75  const tnsr::i<DataVector, 3>&
76  lapse_times_conformal_factor_gradient) noexcept;
77  static void apply(
78  gsl::not_null<tnsr::Ij<DataVector, 3>*>
79  flux_for_conformal_factor_gradient,
80  gsl::not_null<tnsr::Ij<DataVector, 3>*>
81  flux_for_lapse_times_conformal_factor_gradient,
82  const Scalar<DataVector>& conformal_factor,
83  const Scalar<DataVector>& lapse_times_conformal_factor) noexcept;
84 };
85 
86 template <>
88  using argument_tags =
89  tmpl::list<Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>>;
90  static void apply(
91  gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
92  gsl::not_null<tnsr::I<DataVector, 3>*>
93  flux_for_lapse_times_conformal_factor,
94  const tnsr::II<DataVector, 3>& inv_conformal_metric,
95  const tnsr::i<DataVector, 3>& conformal_factor_gradient,
96  const tnsr::i<DataVector, 3>&
97  lapse_times_conformal_factor_gradient) noexcept;
98  static void apply(
99  gsl::not_null<tnsr::Ij<DataVector, 3>*>
100  flux_for_conformal_factor_gradient,
101  gsl::not_null<tnsr::Ij<DataVector, 3>*>
102  flux_for_lapse_times_conformal_factor_gradient,
103  const tnsr::II<DataVector, 3>& inv_conformal_metric,
104  const Scalar<DataVector>& conformal_factor,
105  const Scalar<DataVector>& lapse_times_conformal_factor) noexcept;
106 };
107 
108 template <>
110  using argument_tags = tmpl::list<>;
111  static void apply(
112  gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
113  gsl::not_null<tnsr::I<DataVector, 3>*>
114  flux_for_lapse_times_conformal_factor,
115  gsl::not_null<tnsr::II<DataVector, 3>*> longitudinal_shift_excess,
116  const tnsr::i<DataVector, 3>& conformal_factor_gradient,
117  const tnsr::i<DataVector, 3>& lapse_times_conformal_factor_gradient,
118  const tnsr::ii<DataVector, 3>& shift_strain) noexcept;
119  static void apply(
120  gsl::not_null<tnsr::Ij<DataVector, 3>*>
121  flux_for_conformal_factor_gradient,
122  gsl::not_null<tnsr::Ij<DataVector, 3>*>
123  flux_for_lapse_times_conformal_factor_gradient,
124  gsl::not_null<tnsr::Ijj<DataVector, 3>*> flux_for_shift_strain,
125  const Scalar<DataVector>& conformal_factor,
126  const Scalar<DataVector>& lapse_times_conformal_factor,
127  const tnsr::I<DataVector, 3>& shift) noexcept;
128 };
129 
130 template <>
132  using argument_tags =
133  tmpl::list<Tags::ConformalMetric<DataVector, 3, Frame::Inertial>,
134  Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>>;
135  static void apply(
136  gsl::not_null<tnsr::I<DataVector, 3>*> flux_for_conformal_factor,
137  gsl::not_null<tnsr::I<DataVector, 3>*>
138  flux_for_lapse_times_conformal_factor,
139  gsl::not_null<tnsr::II<DataVector, 3>*> longitudinal_shift_excess,
140  const tnsr::ii<DataVector, 3>& /*conformal_metric*/,
141  const tnsr::II<DataVector, 3>& inv_conformal_metric,
142  const tnsr::i<DataVector, 3>& conformal_factor_gradient,
143  const tnsr::i<DataVector, 3>& lapse_times_conformal_factor_gradient,
144  const tnsr::ii<DataVector, 3>& shift_strain) noexcept;
145  static void apply(
146  gsl::not_null<tnsr::Ij<DataVector, 3>*>
147  flux_for_conformal_factor_gradient,
148  gsl::not_null<tnsr::Ij<DataVector, 3>*>
149  flux_for_lapse_times_conformal_factor_gradient,
150  gsl::not_null<tnsr::Ijj<DataVector, 3>*> flux_for_shift_strain,
151  const tnsr::ii<DataVector, 3>& conformal_metric,
152  const tnsr::II<DataVector, 3>& inv_conformal_metric,
153  const Scalar<DataVector>& conformal_factor,
154  const Scalar<DataVector>& lapse_times_conformal_factor,
155  const tnsr::I<DataVector, 3>& shift_excess) noexcept;
156 };
157 /// \endcond
158 
159 /// The sources \f$S\f$ for the first-order formulation of the XCTS equations.
160 ///
161 /// \see Xcts::FirstOrderSystem
162 template <Equations EnabledEquations, Geometry ConformalGeometry,
163  int ConformalMatterScale>
164 struct Sources;
165 
166 /// \cond
167 template <int ConformalMatterScale>
169  ConformalMatterScale> {
170  using argument_tags = tmpl::list<
172  ConformalMatterScale>,
175  static void apply(
176  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
177  const Scalar<DataVector>& conformal_energy_density,
178  const Scalar<DataVector>& extrinsic_curvature_trace,
179  const Scalar<DataVector>&
180  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
181  const Scalar<DataVector>& conformal_factor,
182  const tnsr::I<DataVector, 3>& conformal_factor_flux) noexcept;
183  static void apply(
184  gsl::not_null<tnsr::i<DataVector, 3>*>
185  equation_for_conformal_factor_gradient,
186  const Scalar<DataVector>& conformal_energy_density,
187  const Scalar<DataVector>& extrinsic_curvature_trace,
188  const Scalar<DataVector>&
189  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
190  const Scalar<DataVector>& conformal_factor) noexcept;
191 };
192 
193 template <int ConformalMatterScale>
194 struct Sources<Equations::Hamiltonian, Geometry::Curved, ConformalMatterScale> {
195  using argument_tags = tmpl::push_back<
197  ConformalMatterScale>::argument_tags,
200  static void apply(
201  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
202  const Scalar<DataVector>& conformal_energy_density,
203  const Scalar<DataVector>& extrinsic_curvature_trace,
204  const Scalar<DataVector>&
205  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
206  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
207  const Scalar<DataVector>& conformal_ricci_scalar,
208  const Scalar<DataVector>& conformal_factor,
209  const tnsr::I<DataVector, 3>& conformal_factor_flux) noexcept;
210  static void apply(
211  gsl::not_null<tnsr::i<DataVector, 3>*>
212  equation_for_conformal_factor_gradient,
213  const Scalar<DataVector>& conformal_energy_density,
214  const Scalar<DataVector>& extrinsic_curvature_trace,
215  const Scalar<DataVector>&
216  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
217  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
218  const Scalar<DataVector>& conformal_ricci_scalar,
219  const Scalar<DataVector>& conformal_factor) noexcept;
220 };
221 
222 template <int ConformalMatterScale>
223 struct Sources<Equations::HamiltonianAndLapse, Geometry::FlatCartesian,
224  ConformalMatterScale> {
225  using argument_tags = tmpl::list<
226  Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
227  ConformalMatterScale>,
228  Tags::Conformal<gr::Tags::StressTrace<DataVector>, ConformalMatterScale>,
231  Tags::LongitudinalShiftMinusDtConformalMetricSquare<DataVector>,
232  Tags::ShiftDotDerivExtrinsicCurvatureTrace<DataVector>>;
233  static void apply(
234  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
235  gsl::not_null<Scalar<DataVector>*> lapse_equation,
236  const Scalar<DataVector>& conformal_energy_density,
237  const Scalar<DataVector>& conformal_stress_trace,
238  const Scalar<DataVector>& extrinsic_curvature_trace,
239  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
240  const Scalar<DataVector>&
241  longitudinal_shift_minus_dt_conformal_metric_square,
242  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
243  const Scalar<DataVector>& conformal_factor,
244  const Scalar<DataVector>& lapse_times_conformal_factor,
245  const tnsr::I<DataVector, 3>& conformal_factor_flux,
246  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux) noexcept;
247  static void apply(
248  gsl::not_null<tnsr::i<DataVector, 3>*>
249  equation_for_conformal_factor_gradient,
250  gsl::not_null<tnsr::i<DataVector, 3>*>
251  equation_for_lapse_times_conformal_factor_gradient,
252  const Scalar<DataVector>& conformal_energy_density,
253  const Scalar<DataVector>& conformal_stress_trace,
254  const Scalar<DataVector>& extrinsic_curvature_trace,
255  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
256  const Scalar<DataVector>&
257  longitudinal_shift_minus_dt_conformal_metric_square,
258  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
259  const Scalar<DataVector>& conformal_factor,
260  const Scalar<DataVector>& lapse_times_conformal_factor) noexcept;
261 };
262 
263 template <int ConformalMatterScale>
265  ConformalMatterScale> {
266  using argument_tags = tmpl::push_back<
268  ConformalMatterScale>::argument_tags,
269  Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
270  Tags::ConformalRicciScalar<DataVector>>;
271  static void apply(
272  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
273  gsl::not_null<Scalar<DataVector>*> lapse_equation,
274  const Scalar<DataVector>& conformal_energy_density,
275  const Scalar<DataVector>& conformal_stress_trace,
276  const Scalar<DataVector>& extrinsic_curvature_trace,
277  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
278  const Scalar<DataVector>&
279  longitudinal_shift_minus_dt_conformal_metric_square,
280  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
281  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
282  const Scalar<DataVector>& conformal_ricci_scalar,
283  const Scalar<DataVector>& conformal_factor,
284  const Scalar<DataVector>& lapse_times_conformal_factor,
285  const tnsr::I<DataVector, 3>& conformal_factor_flux,
286  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux) noexcept;
287  static void apply(
288  gsl::not_null<tnsr::i<DataVector, 3>*>
289  equation_for_conformal_factor_gradient,
290  gsl::not_null<tnsr::i<DataVector, 3>*>
291  equation_for_lapse_times_conformal_factor_gradient,
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,
302  const Scalar<DataVector>& lapse_times_conformal_factor) noexcept;
303 };
304 
305 template <int ConformalMatterScale>
307  ConformalMatterScale> {
308  using argument_tags = tmpl::list<
309  Tags::Conformal<gr::Tags::EnergyDensity<DataVector>,
310  ConformalMatterScale>,
311  Tags::Conformal<gr::Tags::StressTrace<DataVector>, ConformalMatterScale>,
312  Tags::Conformal<gr::Tags::MomentumDensity<3, Frame::Inertial, DataVector>,
313  ConformalMatterScale>,
317  tmpl::size_t<3>, Frame::Inertial>,
318  Tags::ShiftBackground<DataVector, 3, Frame::Inertial>,
319  Tags::LongitudinalShiftBackgroundMinusDtConformalMetric<DataVector, 3,
321  ::Tags::div<Tags::LongitudinalShiftBackgroundMinusDtConformalMetric<
323  static void apply(
324  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
325  gsl::not_null<Scalar<DataVector>*> lapse_equation,
326  gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
327  const Scalar<DataVector>& conformal_energy_density,
328  const Scalar<DataVector>& conformal_stress_trace,
329  const tnsr::I<DataVector, 3>& conformal_momentum_density,
330  const Scalar<DataVector>& extrinsic_curvature_trace,
331  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
332  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
333  const tnsr::I<DataVector, 3>& shift_background,
334  const tnsr::II<DataVector, 3>&
335  longitudinal_shift_background_minus_dt_conformal_metric,
336  const tnsr::I<DataVector, 3>&
337  div_longitudinal_shift_background_minus_dt_conformal_metric,
338  const Scalar<DataVector>& conformal_factor,
339  const Scalar<DataVector>& lapse_times_conformal_factor,
340  const tnsr::I<DataVector, 3>& shift_excess,
341  const tnsr::I<DataVector, 3>& conformal_factor_flux,
342  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
343  const tnsr::II<DataVector, 3>& longitudinal_shift_excess) noexcept;
344  static void apply(
345  gsl::not_null<tnsr::i<DataVector, 3>*>
346  equation_for_conformal_factor_gradient,
347  gsl::not_null<tnsr::i<DataVector, 3>*>
348  equation_for_lapse_times_conformal_factor_gradient,
349  gsl::not_null<tnsr::ii<DataVector, 3>*> equation_for_shift_strain,
350  const Scalar<DataVector>& conformal_energy_density,
351  const Scalar<DataVector>& conformal_stress_trace,
352  const tnsr::I<DataVector, 3>& conformal_momentum_density,
353  const Scalar<DataVector>& extrinsic_curvature_trace,
354  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
355  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
356  const tnsr::I<DataVector, 3>& shift_background,
357  const tnsr::II<DataVector, 3>&
358  longitudinal_shift_background_minus_dt_conformal_metric,
359  const tnsr::I<DataVector, 3>&
360  div_longitudinal_shift_background_minus_dt_conformal_metric,
361  const Scalar<DataVector>& conformal_factor,
362  const Scalar<DataVector>& lapse_times_conformal_factor,
363  const tnsr::I<DataVector, 3>& shift_excess) noexcept;
364 };
365 
366 template <int ConformalMatterScale>
368  ConformalMatterScale> {
369  using argument_tags = tmpl::push_back<
370  typename Sources<Equations::HamiltonianLapseAndShift,
372  ConformalMatterScale>::argument_tags,
373  Tags::ConformalMetric<DataVector, 3, Frame::Inertial>,
374  Tags::InverseConformalMetric<DataVector, 3, Frame::Inertial>,
375  Tags::ConformalChristoffelFirstKind<DataVector, 3, Frame::Inertial>,
376  Tags::ConformalChristoffelSecondKind<DataVector, 3, Frame::Inertial>,
377  Tags::ConformalChristoffelContracted<DataVector, 3, Frame::Inertial>,
378  Tags::ConformalRicciScalar<DataVector>>;
379  static void apply(
380  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
381  gsl::not_null<Scalar<DataVector>*> lapse_equation,
382  gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
383  const Scalar<DataVector>& conformal_energy_density,
384  const Scalar<DataVector>& conformal_stress_trace,
385  const tnsr::I<DataVector, 3>& conformal_momentum_density,
386  const Scalar<DataVector>& extrinsic_curvature_trace,
387  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
388  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
389  const tnsr::I<DataVector, 3>& shift_background,
390  const tnsr::II<DataVector, 3>&
391  longitudinal_shift_background_minus_dt_conformal_metric,
392  const tnsr::I<DataVector, 3>&
393  div_longitudinal_shift_background_minus_dt_conformal_metric,
394  const tnsr::ii<DataVector, 3>& conformal_metric,
395  const tnsr::II<DataVector, 3>& inv_conformal_metric,
396  const tnsr::ijj<DataVector, 3>& /*conformal_christoffel_first_kind*/,
397  const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
398  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
399  const Scalar<DataVector>& conformal_ricci_scalar,
400  const Scalar<DataVector>& conformal_factor,
401  const Scalar<DataVector>& lapse_times_conformal_factor,
402  const tnsr::I<DataVector, 3>& shift_excess,
403  const tnsr::I<DataVector, 3>& conformal_factor_flux,
404  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
405  const tnsr::II<DataVector, 3>& longitudinal_shift_excess) noexcept;
406  static void apply(
407  gsl::not_null<tnsr::i<DataVector, 3>*>
408  equation_for_conformal_factor_gradient,
409  gsl::not_null<tnsr::i<DataVector, 3>*>
410  equation_for_lapse_times_conformal_factor_gradient,
411  gsl::not_null<tnsr::ii<DataVector, 3>*> equation_for_shift_strain,
412  const Scalar<DataVector>& conformal_energy_density,
413  const Scalar<DataVector>& conformal_stress_trace,
414  const tnsr::I<DataVector, 3>& conformal_momentum_density,
415  const Scalar<DataVector>& extrinsic_curvature_trace,
416  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
417  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
418  const tnsr::I<DataVector, 3>& shift_background,
419  const tnsr::II<DataVector, 3>&
420  longitudinal_shift_background_minus_dt_conformal_metric,
421  const tnsr::I<DataVector, 3>&
422  div_longitudinal_shift_background_minus_dt_conformal_metric,
423  const tnsr::ii<DataVector, 3>& conformal_metric,
424  const tnsr::II<DataVector, 3>& inv_conformal_metric,
425  const tnsr::ijj<DataVector, 3>& conformal_christoffel_first_kind,
426  const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
427  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
428  const Scalar<DataVector>& conformal_ricci_scalar,
429  const Scalar<DataVector>& conformal_factor,
430  const Scalar<DataVector>& lapse_times_conformal_factor,
431  const tnsr::I<DataVector, 3>& shift_excess) noexcept;
432 };
433 /// \endcond
434 
435 /// The linearization of the sources \f$S\f$ for the first-order formulation of
436 /// the XCTS equations.
437 ///
438 /// \see Xcts::FirstOrderSystem
439 template <Equations EnabledEquations, Geometry ConformalGeometry,
440  int ConformalMatterScale>
442 
443 /// \cond
444 template <int ConformalMatterScale>
446  ConformalMatterScale> {
447  using argument_tags = tmpl::push_back<
449  ConformalMatterScale>::argument_tags,
451  static void apply(
452  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
453  const Scalar<DataVector>& conformal_energy_density,
454  const Scalar<DataVector>& extrinsic_curvature_trace,
455  const Scalar<DataVector>&
456  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
457  const Scalar<DataVector>& conformal_factor,
458  const Scalar<DataVector>& conformal_factor_correction,
459  const tnsr::I<DataVector, 3>&
460  /*conformal_factor_flux_correction*/) noexcept;
461  static void apply(
462  gsl::not_null<tnsr::i<DataVector, 3>*>
463  source_for_conformal_factor_gradient_correction,
464  const Scalar<DataVector>& conformal_energy_density,
465  const Scalar<DataVector>& extrinsic_curvature_trace,
466  const Scalar<DataVector>&
467  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
468  const Scalar<DataVector>& conformal_factor,
469  const Scalar<DataVector>& conformal_factor_correction) noexcept;
470 };
471 
472 template <int ConformalMatterScale>
474  ConformalMatterScale> {
475  using argument_tags =
477  ConformalMatterScale>::argument_tags,
479  static void apply(
480  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
481  const Scalar<DataVector>& conformal_energy_density,
482  const Scalar<DataVector>& extrinsic_curvature_trace,
483  const Scalar<DataVector>&
484  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
485  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
486  const Scalar<DataVector>& conformal_ricci_scalar,
487  const Scalar<DataVector>& conformal_factor,
488  const Scalar<DataVector>& conformal_factor_correction,
489  const tnsr::I<DataVector, 3>& conformal_factor_flux_correction) noexcept;
490  static void apply(
491  gsl::not_null<tnsr::i<DataVector, 3>*>
492  equation_for_conformal_factor_gradient_correction,
493  const Scalar<DataVector>& conformal_energy_density,
494  const Scalar<DataVector>& extrinsic_curvature_trace,
495  const Scalar<DataVector>&
496  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
497  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
498  const Scalar<DataVector>& conformal_ricci_scalar,
499  const Scalar<DataVector>& conformal_factor,
500  const Scalar<DataVector>& conformal_factor_correction) noexcept;
501 };
502 
503 template <int ConformalMatterScale>
504 struct LinearizedSources<Equations::HamiltonianAndLapse,
505  Geometry::FlatCartesian, ConformalMatterScale> {
506  using argument_tags = tmpl::push_back<
508  ConformalMatterScale>::argument_tags,
509  Tags::ConformalFactor<DataVector>,
510  Tags::LapseTimesConformalFactor<DataVector>>;
511  static void apply(
512  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
513  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
514  const Scalar<DataVector>& conformal_energy_density,
515  const Scalar<DataVector>& conformal_stress_trace,
516  const Scalar<DataVector>& extrinsic_curvature_trace,
517  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
518  const Scalar<DataVector>&
519  longitudinal_shift_minus_dt_conformal_metric_square,
520  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
521  const Scalar<DataVector>& conformal_factor,
522  const Scalar<DataVector>& lapse_times_conformal_factor,
523  const Scalar<DataVector>& conformal_factor_correction,
524  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
525  const tnsr::I<DataVector, 3>& /*conformal_factor_flux_correction*/,
526  const tnsr::I<DataVector, 3>&
527  lapse_times_conformal_factor_flux_correction) noexcept;
528  static void apply(
529  gsl::not_null<tnsr::i<DataVector, 3>*>
530  equation_for_conformal_factor_gradient_correction,
531  gsl::not_null<tnsr::i<DataVector, 3>*>
532  equation_for_lapse_times_conformal_factor_gradient_correction,
533  const Scalar<DataVector>& conformal_energy_density,
534  const Scalar<DataVector>& conformal_stress_trace,
535  const Scalar<DataVector>& extrinsic_curvature_trace,
536  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
537  const Scalar<DataVector>&
538  longitudinal_shift_minus_dt_conformal_metric_square,
539  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
540  const Scalar<DataVector>& conformal_factor,
541  const Scalar<DataVector>& lapse_times_conformal_factor,
542  const Scalar<DataVector>& conformal_factor_correction,
543  const Scalar<DataVector>&
544  lapse_times_conformal_factor_correction) noexcept;
545 };
546 
547 template <int ConformalMatterScale>
548 struct LinearizedSources<Equations::HamiltonianAndLapse, Geometry::Curved,
549  ConformalMatterScale> {
550  using argument_tags = tmpl::push_back<
552  ConformalMatterScale>::argument_tags,
553  Tags::ConformalFactor<DataVector>,
554  Tags::LapseTimesConformalFactor<DataVector>>;
555  static void apply(
556  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
557  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
558  const Scalar<DataVector>& conformal_energy_density,
559  const Scalar<DataVector>& conformal_stress_trace,
560  const Scalar<DataVector>& extrinsic_curvature_trace,
561  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
562  const Scalar<DataVector>&
563  longitudinal_shift_minus_dt_conformal_metric_square,
564  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
565  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
566  const Scalar<DataVector>& conformal_ricci_scalar,
567  const Scalar<DataVector>& conformal_factor,
568  const Scalar<DataVector>& lapse_times_conformal_factor,
569  const Scalar<DataVector>& conformal_factor_correction,
570  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
571  const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
572  const tnsr::I<DataVector, 3>&
573  lapse_times_conformal_factor_flux_correction) noexcept;
574  static void apply(
575  gsl::not_null<tnsr::i<DataVector, 3>*>
576  equation_for_conformal_factor_gradient_correction,
577  gsl::not_null<tnsr::i<DataVector, 3>*>
578  equation_for_lapse_times_conformal_factor_gradient_correction,
579  const Scalar<DataVector>& conformal_energy_density,
580  const Scalar<DataVector>& conformal_stress_trace,
581  const Scalar<DataVector>& extrinsic_curvature_trace,
582  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
583  const Scalar<DataVector>&
584  longitudinal_shift_minus_dt_conformal_metric_square,
585  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
586  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
587  const Scalar<DataVector>& conformal_ricci_scalar,
588  const Scalar<DataVector>& conformal_factor,
589  const Scalar<DataVector>& lapse_times_conformal_factor,
590  const Scalar<DataVector>& conformal_factor_correction,
591  const Scalar<DataVector>&
592  lapse_times_conformal_factor_correction) noexcept;
593 };
594 
595 template <int ConformalMatterScale>
596 struct LinearizedSources<Equations::HamiltonianLapseAndShift,
597  Geometry::FlatCartesian, ConformalMatterScale> {
598  using argument_tags = tmpl::push_back<
599  typename Sources<Equations::HamiltonianLapseAndShift,
601  ConformalMatterScale>::argument_tags,
602  Tags::ConformalFactor<DataVector>,
603  Tags::LapseTimesConformalFactor<DataVector>,
604  Tags::ShiftExcess<DataVector, 3, Frame::Inertial>,
609  Tags::LongitudinalShiftExcess<DataVector, 3, Frame::Inertial>>;
610  static void apply(
611  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
612  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
613  gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
614  const Scalar<DataVector>& conformal_energy_density,
615  const Scalar<DataVector>& conformal_stress_trace,
616  const tnsr::I<DataVector, 3>& conformal_momentum_density,
617  const Scalar<DataVector>& extrinsic_curvature_trace,
618  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
619  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
620  const tnsr::I<DataVector, 3>& shift_background,
621  const tnsr::II<DataVector, 3>&
622  longitudinal_shift_background_minus_dt_conformal_metric,
623  const tnsr::I<DataVector, 3>&
624  div_longitudinal_shift_background_minus_dt_conformal_metric,
625  const Scalar<DataVector>& conformal_factor,
626  const Scalar<DataVector>& lapse_times_conformal_factor,
627  const tnsr::I<DataVector, 3>& shift_excess,
628  const tnsr::I<DataVector, 3>& conformal_factor_flux,
629  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
630  const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
631  const Scalar<DataVector>& conformal_factor_correction,
632  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
633  const tnsr::I<DataVector, 3>& shift_excess_correction,
634  const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
635  const tnsr::I<DataVector, 3>&
636  lapse_times_conformal_factor_flux_correction,
637  const tnsr::II<DataVector, 3>&
638  longitudinal_shift_excess_correction) noexcept;
639  static void apply(
640  gsl::not_null<tnsr::i<DataVector, 3>*>
641  equation_for_conformal_factor_gradient_correction,
642  gsl::not_null<tnsr::i<DataVector, 3>*>
643  equation_for_lapse_times_conformal_factor_gradient_correction,
644  gsl::not_null<tnsr::ii<DataVector, 3>*>
645  equation_for_shift_strain_correction,
646  const Scalar<DataVector>& conformal_energy_density,
647  const Scalar<DataVector>& conformal_stress_trace,
648  const tnsr::I<DataVector, 3>& conformal_momentum_density,
649  const Scalar<DataVector>& extrinsic_curvature_trace,
650  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
651  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
652  const tnsr::I<DataVector, 3>& shift_background,
653  const tnsr::II<DataVector, 3>&
654  longitudinal_shift_background_minus_dt_conformal_metric,
655  const tnsr::I<DataVector, 3>&
656  div_longitudinal_shift_background_minus_dt_conformal_metric,
657  const Scalar<DataVector>& conformal_factor,
658  const Scalar<DataVector>& lapse_times_conformal_factor,
659  const tnsr::I<DataVector, 3>& shift_excess,
660  const tnsr::I<DataVector, 3>& conformal_factor_flux,
661  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
662  const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
663  const Scalar<DataVector>& conformal_factor_correction,
664  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
665  const tnsr::I<DataVector, 3>& shift_excess_correction) noexcept;
666 };
667 
668 template <int ConformalMatterScale>
669 struct LinearizedSources<Equations::HamiltonianLapseAndShift, Geometry::Curved,
670  ConformalMatterScale> {
671  using argument_tags = tmpl::push_back<
673  ConformalMatterScale>::argument_tags,
674  Tags::ConformalFactor<DataVector>,
675  Tags::LapseTimesConformalFactor<DataVector>,
676  Tags::ShiftExcess<DataVector, 3, Frame::Inertial>,
681  Tags::LongitudinalShiftExcess<DataVector, 3, Frame::Inertial>>;
682  static void apply(
683  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
684  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
685  gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
686  const Scalar<DataVector>& conformal_energy_density,
687  const Scalar<DataVector>& conformal_stress_trace,
688  const tnsr::I<DataVector, 3>& conformal_momentum_density,
689  const Scalar<DataVector>& extrinsic_curvature_trace,
690  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
691  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
692  const tnsr::I<DataVector, 3>& shift_background,
693  const tnsr::II<DataVector, 3>&
694  longitudinal_shift_background_minus_dt_conformal_metric,
695  const tnsr::I<DataVector, 3>&
696  /*div_longitudinal_shift_background_minus_dt_conformal_metric*/,
697  const tnsr::ii<DataVector, 3>& conformal_metric,
698  const tnsr::II<DataVector, 3>& inv_conformal_metric,
699  const tnsr::ijj<DataVector, 3>& /*conformal_christoffel_first_kind*/,
700  const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
701  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
702  const Scalar<DataVector>& conformal_ricci_scalar,
703  const Scalar<DataVector>& conformal_factor,
704  const Scalar<DataVector>& lapse_times_conformal_factor,
705  const tnsr::I<DataVector, 3>& shift_excess,
706  const tnsr::I<DataVector, 3>& conformal_factor_flux,
707  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
708  const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
709  const Scalar<DataVector>& conformal_factor_correction,
710  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
711  const tnsr::I<DataVector, 3>& shift_excess_correction,
712  const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
713  const tnsr::I<DataVector, 3>&
714  lapse_times_conformal_factor_flux_correction,
715  const tnsr::II<DataVector, 3>&
716  longitudinal_shift_excess_correction) noexcept;
717  static void apply(
718  gsl::not_null<tnsr::i<DataVector, 3>*>
719  equation_for_conformal_factor_gradient_correction,
720  gsl::not_null<tnsr::i<DataVector, 3>*>
721  equation_for_lapse_times_conformal_factor_gradient_correction,
722  gsl::not_null<tnsr::ii<DataVector, 3>*>
723  equation_for_shift_strain_correction,
724  const Scalar<DataVector>& conformal_energy_density,
725  const Scalar<DataVector>& conformal_stress_trace,
726  const tnsr::I<DataVector, 3>& conformal_momentum_density,
727  const Scalar<DataVector>& extrinsic_curvature_trace,
728  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
729  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
730  const tnsr::I<DataVector, 3>& shift_background,
731  const tnsr::II<DataVector, 3>&
732  longitudinal_shift_background_minus_dt_conformal_metric,
733  const tnsr::I<DataVector, 3>&
734  div_longitudinal_shift_background_minus_dt_conformal_metric,
735  const tnsr::ii<DataVector, 3>& conformal_metric,
736  const tnsr::II<DataVector, 3>& inv_conformal_metric,
737  const tnsr::ijj<DataVector, 3>& conformal_christoffel_first_kind,
738  const tnsr::Ijj<DataVector, 3>& conformal_christoffel_second_kind,
739  const tnsr::i<DataVector, 3>& conformal_christoffel_contracted,
740  const Scalar<DataVector>& conformal_ricci_scalar,
741  const Scalar<DataVector>& conformal_factor,
742  const Scalar<DataVector>& lapse_times_conformal_factor,
743  const tnsr::I<DataVector, 3>& shift_excess,
744  const tnsr::I<DataVector, 3>& conformal_factor_flux,
745  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
746  const tnsr::II<DataVector, 3>& longitudinal_shift_excess,
747  const Scalar<DataVector>& conformal_factor_correction,
748  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
749  const tnsr::I<DataVector, 3>& shift_excess_correction) noexcept;
750 };
751 /// \endcond
752 
753 } // namespace Xcts
gr::Tags::TraceExtrinsicCurvature< DataVector >
std::apply
T apply(T... args)
Xcts::Equations::HamiltonianAndLapse
@ HamiltonianAndLapse
Both the Hamiltonian constraint and the lapse equation, solved for and .
Xcts
Items related to solving the Extended Conformal Thin Sandwich (XCTS) decomposition of the Einstein co...
Definition: Flatness.hpp:22
Frame::Inertial
Definition: IndexType.hpp:44
Xcts::Tags::Conformal
The quantity Tag scaled by the Xcts::Tags::ConformalFactor to the given Power
Definition: Tags.hpp:30
Divergence.hpp
Xcts::Tags::ConformalChristoffelContracted
The Christoffel symbols of the second kind (related to the conformal metric ) contracted in their fir...
Definition: Tags.hpp:222
Tags::Flux
Prefix indicating a flux.
Definition: Prefixes.hpp:40
Xcts::Tags::ConformalFactor
The conformal factor that rescales the spatial metric .
Definition: Tags.hpp:21
domain::push_back
CoordinateMap< SourceFrame, TargetFrame, Maps..., NewMap > push_back(CoordinateMap< SourceFrame, TargetFrame, Maps... > old_map, NewMap new_map) noexcept
Creates a CoordinateMap by appending the new map to the end of the old maps.
Xcts::Geometry::Curved
@ Curved
The conformal geometry is either curved or employs curved coordinates, so non-vanishing Christoffel s...
Xcts::Tags::LongitudinalShiftMinusDtConformalMetricOverLapseSquare
The conformal longitudinal operator applied to the shift vector minus the time derivative of the conf...
Definition: Tags.hpp:181
Xcts::Equations::HamiltonianLapseAndShift
@ HamiltonianLapseAndShift
The full XCTS equations, solved for , and .
Xcts::Sources
The sources for the first-order formulation of the XCTS equations.
Definition: FluxesAndSources.hpp:164
Tags::div
Prefix indicating the divergence.
Definition: Divergence.hpp:44
Xcts::LinearizedSources
The linearization of the sources for the first-order formulation of the XCTS equations.
Definition: FluxesAndSources.hpp:441
cstddef
DataVector
Stores a collection of function values.
Definition: DataVector.hpp:46
Xcts::Equations
Equations
Indicates a subset of the XCTS equations.
Definition: FluxesAndSources.hpp:24
Xcts::Fluxes
The fluxes for the first-order formulation of the XCTS equations.
Definition: FluxesAndSources.hpp:39
Xcts::Geometry::FlatCartesian
@ FlatCartesian
Euclidean (flat) manifold with Cartesian coordinates, i.e. the conformal metric has components in th...
gr::shift
tnsr::I< DataType, SpatialDim, Frame > shift(const tnsr::aa< DataType, SpatialDim, Frame > &spacetime_metric, const tnsr::II< DataType, SpatialDim, Frame > &inverse_spatial_metric) noexcept
Compute shift from spacetime metric and inverse spatial metric.
Tags::dt
Prefix indicating a time derivative.
Definition: Prefixes.hpp:29
Scalar
Tensor< T, Symmetry<>, index_list<> > Scalar
Definition: TypeAliases.hpp:21
Gsl.hpp
Tensor.hpp
Xcts::Geometry
Geometry
Types of conformal geometries for the XCTS equations.
Definition: Geometry.hpp:8
Tags::deriv
Prefix indicating spatial derivatives.
Definition: PartialDerivatives.hpp:52
PartialDerivatives.hpp
Xcts::Equations::Hamiltonian
@ Hamiltonian
Only the Hamiltonian constraint, solved for .
TMPL.hpp
Xcts::Tags::ConformalRicciScalar
The Ricci scalar related to the conformal metric .
Definition: Tags.hpp:240
gsl::not_null
Require a pointer to not be a nullptr
Definition: ReadSpecPiecewisePolynomial.hpp:13