Equations.hpp
1 // Distributed under the MIT License.
2 // See LICENSE.txt for details.
3 
4 #pragma once
5 
6 #include <cstddef>
7 
9 #include "Utilities/Gsl.hpp"
10 
11 /// \cond
12 class DataVector;
13 /// \endcond
14 
15 namespace Xcts {
16 namespace detail {
17 // Tensor-contraction helper functions that should be replaced by tensor
18 // expressions once those work
19 void fully_contract_flat_cartesian(
21  const tnsr::II<DataVector, 3>& tensor1,
22  const tnsr::II<DataVector, 3>& tensor2) noexcept;
23 void fully_contract(gsl::not_null<Scalar<DataVector>*> result,
26  const tnsr::II<DataVector, 3>& tensor1,
27  const tnsr::II<DataVector, 3>& tensor2,
28  const tnsr::ii<DataVector, 3>& metric) noexcept;
29 } // namespace detail
30 
31 /*!
32  * \brief Add the nonlinear source to the Hamiltonian constraint on a flat
33  * conformal background in Cartesian coordinates and with
34  * \f$\bar{u}_{ij}=0=\beta^i\f$.
35  *
36  * Adds \f$\frac{1}{12}\psi^5 K^2 - 2\pi\psi^{5-n}\bar{\rho}\f$ where \f$n\f$ is
37  * the `ConformalMatterScale` and \f$\bar{\rho}=\psi^n\rho\f$ is the
38  * conformally-scaled energy density. Additional sources can be added with
39  * `add_distortion_hamiltonian_sources` and
40  * `add_curved_hamiltonian_or_lapse_sources`.
41  *
42  * \see Xcts
43  */
44 template <int ConformalMatterScale>
46  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
47  const Scalar<DataVector>& conformal_energy_density,
48  const Scalar<DataVector>& extrinsic_curvature_trace,
49  const Scalar<DataVector>& conformal_factor) noexcept;
50 
51 /// The linearization of `add_hamiltonian_sources`
52 template <int ConformalMatterScale>
54  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
55  const Scalar<DataVector>& conformal_energy_density,
56  const Scalar<DataVector>& extrinsic_curvature_trace,
57  const Scalar<DataVector>& conformal_factor,
58  const Scalar<DataVector>& conformal_factor_correction) noexcept;
59 
60 /*!
61  * \brief Add the "distortion" source term to the Hamiltonian constraint.
62  *
63  * Adds \f$-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\f$.
64  *
65  * \see Xcts
66  * \see Xcts::Tags::LongitudinalShiftMinusDtConformalMetricOverLapseSquare
67  */
69  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
70  const Scalar<DataVector>&
71  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
72  const Scalar<DataVector>& conformal_factor) noexcept;
73 
74 /*!
75  * \brief The linearization of `add_distortion_hamiltonian_sources`
76  *
77  * Note that this linearization is only w.r.t. \f$\psi\f$.
78  */
80  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
81  const Scalar<DataVector>&
82  longitudinal_shift_minus_dt_conformal_metric_over_lapse_square,
83  const Scalar<DataVector>& conformal_factor,
84  const Scalar<DataVector>& conformal_factor_correction) noexcept;
85 
86 /*!
87  * \brief Add the contributions from a curved background geometry to the
88  * Hamiltonian constraint or lapse equation
89  *
90  * Adds \f$\frac{1}{8}\psi\bar{R}\f$. This term appears both in the Hamiltonian
91  * constraint and the lapse equation (where in the latter \f$\psi\f$ is replaced
92  * by \f$\alpha\psi\f$).
93  *
94  * This term is linear.
95  *
96  * \see Xcts
97  */
99  gsl::not_null<Scalar<DataVector>*> hamiltonian_or_lapse_equation,
100  const Scalar<DataVector>& conformal_ricci_scalar,
101  const Scalar<DataVector>& field) noexcept;
102 
103 /*!
104  * \brief Add the nonlinear source to the lapse equation on a flat conformal
105  * background in Cartesian coordinates and with \f$\bar{u}_{ij}=0=\beta^i\f$.
106  *
107  * Adds \f$(\alpha\psi)\left(\frac{5}{12}\psi^4 K^2 + 2\pi\psi^{4-n}
108  * \left(\bar{\rho} + 2\bar{S}\right)\right) + \psi^5
109  * \left(\beta^i\partial_i K - \partial_t K\right)\f$ where \f$n\f$ is the
110  * `ConformalMatterScale` and matter quantities are conformally-scaled.
111  * Additional sources can be added with
112  * `add_distortion_hamiltonian_and_lapse_sources` and
113  * `add_curved_hamiltonian_or_lapse_sources`.
114  *
115  * \see Xcts
116  */
117 template <int ConformalMatterScale>
118 void add_lapse_sources(
119  gsl::not_null<Scalar<DataVector>*> lapse_equation,
120  const Scalar<DataVector>& conformal_energy_density,
121  const Scalar<DataVector>& conformal_stress_trace,
122  const Scalar<DataVector>& extrinsic_curvature_trace,
123  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
124  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
125  const Scalar<DataVector>& conformal_factor,
126  const Scalar<DataVector>& lapse_times_conformal_factor) noexcept;
127 
128 /*!
129  * \brief The linearization of `add_lapse_sources`
130  *
131  * Note that this linearization is only w.r.t. \f$\psi\f$ and \f$\alpha\psi\f$.
132  * The linearization w.r.t. \f$\beta^i\f$ is added in
133  * `add_curved_linearized_momentum_sources` or
134  * `add_flat_cartesian_linearized_momentum_sources`.
135  */
136 template <int ConformalMatterScale>
138  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
139  const Scalar<DataVector>& conformal_energy_density,
140  const Scalar<DataVector>& conformal_stress_trace,
141  const Scalar<DataVector>& extrinsic_curvature_trace,
142  const Scalar<DataVector>& dt_extrinsic_curvature_trace,
143  const Scalar<DataVector>& shift_dot_deriv_extrinsic_curvature_trace,
144  const Scalar<DataVector>& conformal_factor,
145  const Scalar<DataVector>& lapse_times_conformal_factor,
146  const Scalar<DataVector>& conformal_factor_correction,
147  const Scalar<DataVector>& lapse_times_conformal_factor_correction) noexcept;
148 
149 /*!
150  * \brief Add the "distortion" source term to the Hamiltonian constraint and the
151  * lapse equation.
152  *
153  * Adds \f$-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\f$ to the Hamiltonian
154  * constraint and \f$\frac{7}{8}\alpha\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\f$ to
155  * the lapse equation.
156  *
157  * \see Xcts
158  * \see Xcts::Tags::LongitudinalShiftMinusDtConformalMetricSquare
159  */
161  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
162  gsl::not_null<Scalar<DataVector>*> lapse_equation,
163  const Scalar<DataVector>&
164  longitudinal_shift_minus_dt_conformal_metric_square,
165  const Scalar<DataVector>& conformal_factor,
166  const Scalar<DataVector>& lapse_times_conformal_factor) noexcept;
167 
168 /*!
169  * \brief The linearization of `add_distortion_hamiltonian_and_lapse_sources`
170  *
171  * Note that this linearization is only w.r.t. \f$\psi\f$ and \f$\alpha\psi\f$.
172  */
174  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
175  gsl::not_null<Scalar<DataVector>*> lapse_equation,
176  const Scalar<DataVector>&
177  longitudinal_shift_minus_dt_conformal_metric_square,
178  const Scalar<DataVector>& conformal_factor,
179  const Scalar<DataVector>& lapse_times_conformal_factor,
180  const Scalar<DataVector>& conformal_factor_correction,
181  const Scalar<DataVector>& lapse_times_conformal_factor_correction) noexcept;
182 
183 /*!
184  * \brief Add the nonlinear source to the momentum constraint and add the
185  * "distortion" source term to the Hamiltonian constraint and lapse equation.
186  *
187  * Adds \f$\left((\bar{L}\beta)^{ij} - \bar{u}^{ij}\right)
188  * \left(\frac{\partial_j (\alpha\psi)}{\alpha\psi}
189  * - 7 \frac{\partial_j \psi}{\psi}\right)
190  * + \partial_j\bar{u}^{ij}
191  * + \frac{4}{3}\frac{\alpha\psi}{\psi}\bar{\gamma}^{ij}\partial_j K
192  * + 16\pi\left(\alpha\psi\right)\psi^{3-n} \bar{S}^i\f$ to the momentum
193  * constraint, where \f$n\f$ is the `ConformalMatterScale` and
194  * \f$\bar{S}^i=\psi^n S^i\f$ is the conformally-scaled momentum density.
195  *
196  * Note that the \f$\partial_j\bar{u}^{ij}\f$ term is not the full covariant
197  * divergence, but only the partial-derivatives part of it. The curved
198  * contribution to this term can be added together with the curved contribution
199  * to the flux divergence of the dynamic shift variable with the
200  * `Elasticity::add_curved_sources` function.
201  *
202  * Also adds \f$-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\f$ to the
203  * Hamiltonian constraint and
204  * \f$\frac{7}{8}\alpha\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\f$ to the lapse
205  * equation.
206  *
207  * \see Xcts
208  */
209 /// @{
210 template <int ConformalMatterScale>
212  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
213  gsl::not_null<Scalar<DataVector>*> lapse_equation,
214  gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
215  const tnsr::I<DataVector, 3>& conformal_momentum_density,
216  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
217  const tnsr::ii<DataVector, 3>& conformal_metric,
218  const tnsr::II<DataVector, 3>& inv_conformal_metric,
219  const tnsr::I<DataVector, 3>& minus_div_dt_conformal_metric,
220  const Scalar<DataVector>& conformal_factor,
221  const Scalar<DataVector>& lapse_times_conformal_factor,
222  const tnsr::I<DataVector, 3>& conformal_factor_flux,
223  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
224  const tnsr::II<DataVector, 3>&
225  longitudinal_shift_minus_dt_conformal_metric) noexcept;
226 
227 template <int ConformalMatterScale>
229  gsl::not_null<Scalar<DataVector>*> hamiltonian_constraint,
230  gsl::not_null<Scalar<DataVector>*> lapse_equation,
231  gsl::not_null<tnsr::I<DataVector, 3>*> momentum_constraint,
232  const tnsr::I<DataVector, 3>& conformal_momentum_density,
233  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
234  const tnsr::I<DataVector, 3>& minus_div_dt_conformal_metric,
235  const Scalar<DataVector>& conformal_factor,
236  const Scalar<DataVector>& lapse_times_conformal_factor,
237  const tnsr::I<DataVector, 3>& conformal_factor_flux,
238  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
239  const tnsr::II<DataVector, 3>&
240  longitudinal_shift_minus_dt_conformal_metric) noexcept;
241 /// @}
242 
243 /// The linearization of `add_curved_momentum_sources`
244 template <int ConformalMatterScale>
246  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
247  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
248  gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
249  const tnsr::I<DataVector, 3>& conformal_momentum_density,
250  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
251  const tnsr::ii<DataVector, 3>& conformal_metric,
252  const tnsr::II<DataVector, 3>& inv_conformal_metric,
253  const Scalar<DataVector>& conformal_factor,
254  const Scalar<DataVector>& lapse_times_conformal_factor,
255  const tnsr::I<DataVector, 3>& conformal_factor_flux,
256  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
257  const tnsr::II<DataVector, 3>& longitudinal_shift_minus_dt_conformal_metric,
258  const Scalar<DataVector>& conformal_factor_correction,
259  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
260  const tnsr::I<DataVector, 3>& shift_correction,
261  const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
262  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux_correction,
263  const tnsr::II<DataVector, 3>& longitudinal_shift_correction) noexcept;
264 
265 /// The linearization of `add_flat_cartesian_momentum_sources`
266 template <int ConformalMatterScale>
268  gsl::not_null<Scalar<DataVector>*> linearized_hamiltonian_constraint,
269  gsl::not_null<Scalar<DataVector>*> linearized_lapse_equation,
270  gsl::not_null<tnsr::I<DataVector, 3>*> linearized_momentum_constraint,
271  const tnsr::I<DataVector, 3>& conformal_momentum_density,
272  const tnsr::i<DataVector, 3>& extrinsic_curvature_trace_gradient,
273  const Scalar<DataVector>& conformal_factor,
274  const Scalar<DataVector>& lapse_times_conformal_factor,
275  const tnsr::I<DataVector, 3>& conformal_factor_flux,
276  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux,
277  const tnsr::II<DataVector, 3>& longitudinal_shift_minus_dt_conformal_metric,
278  const Scalar<DataVector>& conformal_factor_correction,
279  const Scalar<DataVector>& lapse_times_conformal_factor_correction,
280  const tnsr::I<DataVector, 3>& shift_correction,
281  const tnsr::I<DataVector, 3>& conformal_factor_flux_correction,
282  const tnsr::I<DataVector, 3>& lapse_times_conformal_factor_flux_correction,
283  const tnsr::II<DataVector, 3>& longitudinal_shift_correction) noexcept;
284 } // namespace Xcts
Xcts
Items related to solving the Extended Conformal Thin Sandwich (XCTS) decomposition of the Einstein co...
Definition: Flatness.hpp:22
Xcts::add_linearized_distortion_hamiltonian_sources
void add_linearized_distortion_hamiltonian_sources(gsl::not_null< Scalar< DataVector > * > linearized_hamiltonian_constraint, const Scalar< DataVector > &longitudinal_shift_minus_dt_conformal_metric_over_lapse_square, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &conformal_factor_correction) noexcept
The linearization of add_distortion_hamiltonian_sources
Xcts::add_distortion_hamiltonian_and_lapse_sources
void add_distortion_hamiltonian_and_lapse_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_constraint, gsl::not_null< Scalar< DataVector > * > lapse_equation, const Scalar< DataVector > &longitudinal_shift_minus_dt_conformal_metric_square, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor) noexcept
Add the "distortion" source term to the Hamiltonian constraint and the lapse equation.
Xcts::add_distortion_hamiltonian_sources
void add_distortion_hamiltonian_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_constraint, const Scalar< DataVector > &longitudinal_shift_minus_dt_conformal_metric_over_lapse_square, const Scalar< DataVector > &conformal_factor) noexcept
Add the "distortion" source term to the Hamiltonian constraint.
Xcts::add_curved_linearized_momentum_sources
void add_curved_linearized_momentum_sources(gsl::not_null< Scalar< DataVector > * > linearized_hamiltonian_constraint, gsl::not_null< Scalar< DataVector > * > linearized_lapse_equation, gsl::not_null< tnsr::I< DataVector, 3 > * > linearized_momentum_constraint, const tnsr::I< DataVector, 3 > &conformal_momentum_density, const tnsr::i< DataVector, 3 > &extrinsic_curvature_trace_gradient, const tnsr::ii< DataVector, 3 > &conformal_metric, const tnsr::II< DataVector, 3 > &inv_conformal_metric, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor, const tnsr::I< DataVector, 3 > &conformal_factor_flux, const tnsr::I< DataVector, 3 > &lapse_times_conformal_factor_flux, const tnsr::II< DataVector, 3 > &longitudinal_shift_minus_dt_conformal_metric, const Scalar< DataVector > &conformal_factor_correction, const Scalar< DataVector > &lapse_times_conformal_factor_correction, const tnsr::I< DataVector, 3 > &shift_correction, const tnsr::I< DataVector, 3 > &conformal_factor_flux_correction, const tnsr::I< DataVector, 3 > &lapse_times_conformal_factor_flux_correction, const tnsr::II< DataVector, 3 > &longitudinal_shift_correction) noexcept
The linearization of add_curved_momentum_sources
Xcts::add_flat_cartesian_linearized_momentum_sources
void add_flat_cartesian_linearized_momentum_sources(gsl::not_null< Scalar< DataVector > * > linearized_hamiltonian_constraint, gsl::not_null< Scalar< DataVector > * > linearized_lapse_equation, gsl::not_null< tnsr::I< DataVector, 3 > * > linearized_momentum_constraint, const tnsr::I< DataVector, 3 > &conformal_momentum_density, const tnsr::i< DataVector, 3 > &extrinsic_curvature_trace_gradient, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor, const tnsr::I< DataVector, 3 > &conformal_factor_flux, const tnsr::I< DataVector, 3 > &lapse_times_conformal_factor_flux, const tnsr::II< DataVector, 3 > &longitudinal_shift_minus_dt_conformal_metric, const Scalar< DataVector > &conformal_factor_correction, const Scalar< DataVector > &lapse_times_conformal_factor_correction, const tnsr::I< DataVector, 3 > &shift_correction, const tnsr::I< DataVector, 3 > &conformal_factor_flux_correction, const tnsr::I< DataVector, 3 > &lapse_times_conformal_factor_flux_correction, const tnsr::II< DataVector, 3 > &longitudinal_shift_correction) noexcept
The linearization of add_flat_cartesian_momentum_sources
Xcts::add_flat_cartesian_momentum_sources
void add_flat_cartesian_momentum_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_constraint, gsl::not_null< Scalar< DataVector > * > lapse_equation, gsl::not_null< tnsr::I< DataVector, 3 > * > momentum_constraint, const tnsr::I< DataVector, 3 > &conformal_momentum_density, const tnsr::i< DataVector, 3 > &extrinsic_curvature_trace_gradient, const tnsr::I< DataVector, 3 > &minus_div_dt_conformal_metric, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor, const tnsr::I< DataVector, 3 > &conformal_factor_flux, const tnsr::I< DataVector, 3 > &lapse_times_conformal_factor_flux, const tnsr::II< DataVector, 3 > &longitudinal_shift_minus_dt_conformal_metric) noexcept
Add the nonlinear source to the momentum constraint and add the "distortion" source term to the Hamil...
cstddef
Xcts::add_curved_hamiltonian_or_lapse_sources
void add_curved_hamiltonian_or_lapse_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_or_lapse_equation, const Scalar< DataVector > &conformal_ricci_scalar, const Scalar< DataVector > &field) noexcept
Add the contributions from a curved background geometry to the Hamiltonian constraint or lapse equati...
DataVector
Stores a collection of function values.
Definition: DataVector.hpp:46
Xcts::add_linearized_distortion_hamiltonian_and_lapse_sources
void add_linearized_distortion_hamiltonian_and_lapse_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_constraint, gsl::not_null< Scalar< DataVector > * > lapse_equation, const Scalar< DataVector > &longitudinal_shift_minus_dt_conformal_metric_square, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor, const Scalar< DataVector > &conformal_factor_correction, const Scalar< DataVector > &lapse_times_conformal_factor_correction) noexcept
The linearization of add_distortion_hamiltonian_and_lapse_sources
Xcts::add_curved_momentum_sources
void add_curved_momentum_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_constraint, gsl::not_null< Scalar< DataVector > * > lapse_equation, gsl::not_null< tnsr::I< DataVector, 3 > * > momentum_constraint, const tnsr::I< DataVector, 3 > &conformal_momentum_density, const tnsr::i< DataVector, 3 > &extrinsic_curvature_trace_gradient, const tnsr::ii< DataVector, 3 > &conformal_metric, const tnsr::II< DataVector, 3 > &inv_conformal_metric, const tnsr::I< DataVector, 3 > &minus_div_dt_conformal_metric, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor, const tnsr::I< DataVector, 3 > &conformal_factor_flux, const tnsr::I< DataVector, 3 > &lapse_times_conformal_factor_flux, const tnsr::II< DataVector, 3 > &longitudinal_shift_minus_dt_conformal_metric) noexcept
Add the nonlinear source to the momentum constraint and add the "distortion" source term to the Hamil...
Xcts::add_lapse_sources
void add_lapse_sources(gsl::not_null< Scalar< DataVector > * > lapse_equation, const Scalar< DataVector > &conformal_energy_density, const Scalar< DataVector > &conformal_stress_trace, const Scalar< DataVector > &extrinsic_curvature_trace, const Scalar< DataVector > &dt_extrinsic_curvature_trace, const Scalar< DataVector > &shift_dot_deriv_extrinsic_curvature_trace, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor) noexcept
Add the nonlinear source to the lapse equation on a flat conformal background in Cartesian coordinate...
Xcts::add_linearized_lapse_sources
void add_linearized_lapse_sources(gsl::not_null< Scalar< DataVector > * > linearized_lapse_equation, const Scalar< DataVector > &conformal_energy_density, const Scalar< DataVector > &conformal_stress_trace, const Scalar< DataVector > &extrinsic_curvature_trace, const Scalar< DataVector > &dt_extrinsic_curvature_trace, const Scalar< DataVector > &shift_dot_deriv_extrinsic_curvature_trace, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &lapse_times_conformal_factor, const Scalar< DataVector > &conformal_factor_correction, const Scalar< DataVector > &lapse_times_conformal_factor_correction) noexcept
The linearization of add_lapse_sources
Scalar
Tensor< T, Symmetry<>, index_list<> > Scalar
Definition: TypeAliases.hpp:21
Gsl.hpp
TypeAliases.hpp
Xcts::add_hamiltonian_sources
void add_hamiltonian_sources(gsl::not_null< Scalar< DataVector > * > hamiltonian_constraint, const Scalar< DataVector > &conformal_energy_density, const Scalar< DataVector > &extrinsic_curvature_trace, const Scalar< DataVector > &conformal_factor) noexcept
Add the nonlinear source to the Hamiltonian constraint on a flat conformal background in Cartesian co...
Xcts::add_linearized_hamiltonian_sources
void add_linearized_hamiltonian_sources(gsl::not_null< Scalar< DataVector > * > linearized_hamiltonian_constraint, const Scalar< DataVector > &conformal_energy_density, const Scalar< DataVector > &extrinsic_curvature_trace, const Scalar< DataVector > &conformal_factor, const Scalar< DataVector > &conformal_factor_correction) noexcept
The linearization of add_hamiltonian_sources
gsl::not_null
Require a pointer to not be a nullptr
Definition: ReadSpecPiecewisePolynomial.hpp:13