SpECTRE  v2024.04.12
Xcts Namespace Reference

Items related to solving the Extended Conformal Thin Sandwich (XCTS) decomposition of the Einstein constraint equations. More...

Namespaces

namespace  AnalyticData
 Analytic data for the XCTS equations, i.e. field configurations that do not solve the equations but are used as background or initial guess.
 
namespace  Solutions
 Analytic solutions of the XCTS equations.
 
namespace  Tags
 Tags related to the XCTS equations.
 

Classes

struct  FirstOrderSystem
 The Extended Conformal Thin Sandwich (XCTS) decomposition of the Einstein constraint equations, formulated as a set of coupled first-order partial differential equations. More...
 
struct  Fluxes
 The fluxes \(F^i\) for the first-order formulation of the XCTS equations. More...
 
struct  LinearizedSources
 The linearization of the sources \(S\) for the first-order formulation of the XCTS equations. More...
 
struct  Sources
 The sources \(S\) for the first-order formulation of the XCTS equations. More...
 
struct  SpacetimeQuantitiesComputer
 CachedTempBuffer computer class for 3+1 quantities from XCTS variables. See Xcts::SpacetimeQuantities. More...
 

Typedefs

using SpacetimeQuantities = CachedTempBuffer< Tags::ConformalFactor< DataVector >, Tags::LapseTimesConformalFactor< DataVector >, ::Tags::deriv< Tags::ConformalFactor< DataVector >, tmpl::size_t< 3 >, Frame::Inertial >, detail::ConformalLaplacianOfConformalFactor< DataVector >, ::Tags::deriv< Tags::LapseTimesConformalFactor< DataVector >, tmpl::size_t< 3 >, Frame::Inertial >, ::Tags::deriv< Tags::ShiftExcess< DataVector, 3, Frame::Inertial >, tmpl::size_t< 3 >, Frame::Inertial >, Xcts::Tags::ShiftStrain< DataVector, 3, Frame::Inertial >, Xcts::Tags::LongitudinalShiftExcess< DataVector, 3, Frame::Inertial >, ::Tags::div< Tags::LongitudinalShiftExcess< DataVector, 3, Frame::Inertial > >, detail::LongitudinalShiftMinusDtConformalMetric< DataVector >, gr::Tags::SpatialMetric< DataVector, 3 >, gr::Tags::InverseSpatialMetric< DataVector, 3 >, gr::Tags::Lapse< DataVector >, gr::Tags::Shift< DataVector, 3 >, gr::Tags::ExtrinsicCurvature< DataVector, 3 >, gr::Tags::SpatialChristoffelSecondKind< DataVector, 3 >, gr::Tags::HamiltonianConstraint< DataVector >, gr::Tags::MomentumConstraint< DataVector, 3 > >
 General-relativistic 3+1 quantities computed from XCTS variables.
 

Enumerations

enum class  Equations { Hamiltonian , HamiltonianAndLapse , HamiltonianLapseAndShift }
 Indicates a subset of the XCTS equations. More...
 
enum class  Geometry { FlatCartesian , Curved }
 Types of conformal geometries for the XCTS equations. More...
 

Functions

template<int ConformalMatterScale>
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_minus_one)
 Add the nonlinear source to the Hamiltonian constraint on a flat conformal background in Cartesian coordinates and with \(\bar{u}_{ij}=0=\beta^i\). More...
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &conformal_factor_correction)
 The linearization of add_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_minus_one)
 Add the "distortion" source term to the Hamiltonian constraint. More...
 
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_minus_one, const Scalar< DataVector > &conformal_factor_correction)
 The linearization of add_distortion_hamiltonian_sources More...
 
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, double add_constant=0.)
 Add the contributions from a curved background geometry to the Hamiltonian constraint or lapse equation. More...
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one)
 Add the nonlinear source to the lapse equation on a flat conformal background in Cartesian coordinates and with \(\bar{u}_{ij}=0=\beta^i\). More...
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one, const Scalar< DataVector > &conformal_factor_correction, const Scalar< DataVector > &lapse_times_conformal_factor_correction)
 The linearization of add_lapse_sources More...
 
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one)
 Add the "distortion" source term to the Hamiltonian constraint and the lapse equation. More...
 
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one, const Scalar< DataVector > &conformal_factor_correction, const Scalar< DataVector > &lapse_times_conformal_factor_correction)
 The linearization of add_distortion_hamiltonian_and_lapse_sources More...
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one, 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)
 The linearization of add_curved_momentum_sources
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one, 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)
 The linearization of add_flat_cartesian_momentum_sources
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one, 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)
 Add the nonlinear source to the momentum constraint and add the "distortion" source term to the Hamiltonian constraint and lapse equation. More...
 
template<int ConformalMatterScale>
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_minus_one, const Scalar< DataVector > &lapse_times_conformal_factor_minus_one, 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)
 Add the nonlinear source to the momentum constraint and add the "distortion" source term to the Hamiltonian constraint and lapse equation. More...
 
template<typename DataType >
void extrinsic_curvature (const gsl::not_null< tnsr::ii< DataType, 3 > * > result, const Scalar< DataType > &conformal_factor, const Scalar< DataType > &lapse, const tnsr::ii< DataType, 3 > &conformal_metric, const tnsr::II< DataType, 3 > &longitudinal_shift_minus_dt_conformal_metric, const Scalar< DataType > &trace_extrinsic_curvature)
 Extrinsic curvature computed from the conformal decomposition used in the XCTS system. More...
 
template<typename DataType >
tnsr::ii< DataType, 3 > extrinsic_curvature (const Scalar< DataType > &conformal_factor, const Scalar< DataType > &lapse, const tnsr::ii< DataType, 3 > &conformal_metric, const tnsr::II< DataType, 3 > &longitudinal_shift_minus_dt_conformal_metric, const Scalar< DataType > &trace_extrinsic_curvature)
 Return-by-value overload.
 
template<typename DataType >
void longitudinal_operator (gsl::not_null< tnsr::II< DataType, 3 > * > result, const tnsr::ii< DataType, 3 > &strain, const tnsr::II< DataType, 3 > &inv_metric)
 The longitudinal operator, or vector gradient, \((L\beta)^{ij}\). More...
 
template<typename DataType >
void longitudinal_operator (gsl::not_null< tnsr::II< DataType, 3 > * > result, const tnsr::I< DataType, 3 > &shift, const tnsr::iJ< DataType, 3 > &deriv_shift, const tnsr::II< DataType, 3 > &inv_metric, const tnsr::Ijj< DataType, 3 > &christoffel_second_kind)
 The longitudinal operator, or vector gradient, \((L\beta)^{ij}\). More...
 
template<typename DataType >
void longitudinal_operator_flat_cartesian (gsl::not_null< tnsr::II< DataType, 3 > * > result, const tnsr::ii< DataType, 3 > &strain)
 The conformal longitudinal operator \((L\beta)^{ij}\) on a flat conformal metric in Cartesian coordinates \(\gamma_{ij}=\delta_{ij}\). More...
 

Detailed Description

Items related to solving the Extended Conformal Thin Sandwich (XCTS) decomposition of the Einstein constraint equations.

The XCTS equations

\begin{align} \bar{D}^2 \psi - \frac{1}{8}\psi\bar{R} - \frac{1}{12}\psi^5 K^2 + \frac{1}{8}\psi^{-7}\bar{A}_{ij}\bar{A}^{ij} &= -2\pi\psi^5\rho \\ \bar{D}_i(\bar{L}\beta)^{ij} - (\bar{L}\beta)^{ij}\bar{D}_i \ln(\bar{\alpha}) &= \bar{\alpha}\bar{D}_i\left(\bar{\alpha}^{-1}\bar{u}^{ij} \right) + \frac{4}{3}\bar{\alpha}\psi^6\bar{D}^j K + 16\pi\bar{\alpha} \psi^{10}S^j \\ \bar{D}^2\left(\alpha\psi\right) &= \alpha\psi\left(\frac{7}{8}\psi^{-8}\bar{A}_{ij}\bar{A}^{ij} + \frac{5}{12}\psi^4 K^2 + \frac{1}{8}\bar{R} + 2\pi\psi^4\left(\rho + 2S\right)\right) - \psi^5\partial_t K + \psi^5\beta^i\bar{D}_i K \\ \text{with} \quad \bar{A} &= \frac{1}{2\bar{\alpha}} \left((\bar{L}\beta)^{ij} - \bar{u}^{ij}\right) \\ \quad \text{and} \quad \bar{\alpha} &= \alpha \psi^{-6} \end{align}

are a set of nonlinear elliptic equations that the spacetime metric in general relativity must satisfy at all times. For an introduction see e.g. [13], in particular Box 3.3 which is largely mirrored here. We solve the XCTS equations for the conformal factor \(\psi\), the product of lapse times conformal factor \(\alpha\psi\) and the shift vector \(\beta^j\). The remaining quantities in the equations, i.e. the conformal metric \(\bar{\gamma}_{ij}\), the trace of the extrinsic curvature \(K\), their respective time derivatives \(\bar{u}_{ij}\) and \(\partial_t K\), the energy density \(\rho\), the stress-energy trace \(S\) and the momentum density \(S^i\), are freely specifyable fields that define the physical scenario at hand. Of particular importance is the conformal metric, which defines the background geometry, the covariant derivative \(\bar{D}\), the Ricci scalar \(\bar{R}\) and the longitudinal operator

\begin{equation} \left(\bar{L}\beta\right)^{ij} = \bar{D}^i\beta^j + \bar{D}^j\beta^i - \frac{2}{3}\bar{\gamma}^{ij}\bar{D}_k\beta^k \text{.} \end{equation}

Note that the XCTS equations are essentially two Poisson equations and one Elasticity equation with nonlinear sources on a curved geometry. In this analogy, the longitudinal operator plays the role of the elastic constitutive relation that connects the symmetric "shift strain" \(\bar{D}_{(i}\beta_{j)}\) with the "stress" \((\bar{L}\beta)^{ij}\) of which we take the divergence in the momentum constraint. This particular constitutive relation is equivalent to an isotropic and homogeneous material with bulk modulus \(K=0\) (not to be confused with the extrinsic curvature trace \(K\) in this context) and shear modulus \(\mu=1\) (see Elasticity::ConstitutiveRelations::IsotropicHomogeneous).

Once the XCTS equations are solved we can construct the spatial metric and extrinsic curvature as

\begin{align} \gamma_{ij} &= \psi^4\bar{\gamma}_{ij} \\ K_{ij} &= \psi^{-2}\bar{A}_{ij} + \frac{1}{3}\gamma_{ij} K \end{align}

from which we can compose the full spacetime metric.

Enumeration Type Documentation

◆ Equations

enum class Xcts::Equations
strong

Indicates a subset of the XCTS equations.

Enumerator
Hamiltonian 

Only the Hamiltonian constraint, solved for \(\psi\).

HamiltonianAndLapse 

Both the Hamiltonian constraint and the lapse equation, solved for \(\psi\) and \(\alpha\psi\).

HamiltonianLapseAndShift 

The full XCTS equations, solved for \(\psi\), \(\alpha\psi\) and \(\beta_\mathrm{excess}\).

◆ Geometry

enum class Xcts::Geometry
strong

Types of conformal geometries for the XCTS equations.

Enumerator
FlatCartesian 

Euclidean (flat) manifold with Cartesian coordinates, i.e. the conformal metric has components \(\bar{\gamma}_{ij} = \delta_{ij}\) in these coordinates and thus all Christoffel symbols vanish: \(\bar{\Gamma}^i_{jk}=0\).

Curved 

The conformal geometry is either curved or employs curved coordinates, so non-vanishing Christoffel symbols must be taken into account.

Function Documentation

◆ add_curved_hamiltonian_or_lapse_sources()

void Xcts::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,
double  add_constant = 0. 
)

Add the contributions from a curved background geometry to the Hamiltonian constraint or lapse equation.

Adds \(\frac{1}{8}\psi\bar{R}\). This term appears both in the Hamiltonian constraint and the lapse equation (where in the latter \(\psi\) is replaced by \(\alpha\psi\)).

This term is linear.

See also
Xcts

◆ add_curved_momentum_sources()

template<int ConformalMatterScale>
void Xcts::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_minus_one,
const Scalar< DataVector > &  lapse_times_conformal_factor_minus_one,
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 
)

Add the nonlinear source to the momentum constraint and add the "distortion" source term to the Hamiltonian constraint and lapse equation.

Adds \(\left((\bar{L}\beta)^{ij} - \bar{u}^{ij}\right) \left(\frac{\partial_j (\alpha\psi)}{\alpha\psi} - 7 \frac{\partial_j \psi}{\psi}\right) + \partial_j\bar{u}^{ij} + \frac{4}{3}\frac{\alpha\psi}{\psi}\bar{\gamma}^{ij}\partial_j K + 16\pi\left(\alpha\psi\right)\psi^{3-n} \bar{S}^i\) to the momentum constraint, where \(n\) is the ConformalMatterScale and \(\bar{S}^i=\psi^n S^i\) is the conformally-scaled momentum density.

Note that the \(\partial_j\bar{u}^{ij}\) term is not the full covariant divergence, but only the partial-derivatives part of it. The curved contribution to this term can be added together with the curved contribution to the flux divergence of the dynamic shift variable with the Elasticity::add_curved_sources function.

Also adds \(-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\) to the Hamiltonian constraint and \(\frac{7}{8}\alpha\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\) to the lapse equation.

See also
Xcts

◆ add_distortion_hamiltonian_and_lapse_sources()

void Xcts::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_minus_one,
const Scalar< DataVector > &  lapse_times_conformal_factor_minus_one 
)

Add the "distortion" source term to the Hamiltonian constraint and the lapse equation.

Adds \(-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\) to the Hamiltonian constraint and \(\frac{7}{8}\alpha\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\) to the lapse equation.

See also
Xcts
Xcts::Tags::LongitudinalShiftMinusDtConformalMetricSquare

◆ add_distortion_hamiltonian_sources()

void Xcts::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_minus_one 
)

Add the "distortion" source term to the Hamiltonian constraint.

Adds \(-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\).

See also
Xcts
Xcts::Tags::LongitudinalShiftMinusDtConformalMetricOverLapseSquare

◆ add_flat_cartesian_momentum_sources()

template<int ConformalMatterScale>
void Xcts::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_minus_one,
const Scalar< DataVector > &  lapse_times_conformal_factor_minus_one,
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 
)

Add the nonlinear source to the momentum constraint and add the "distortion" source term to the Hamiltonian constraint and lapse equation.

Adds \(\left((\bar{L}\beta)^{ij} - \bar{u}^{ij}\right) \left(\frac{\partial_j (\alpha\psi)}{\alpha\psi} - 7 \frac{\partial_j \psi}{\psi}\right) + \partial_j\bar{u}^{ij} + \frac{4}{3}\frac{\alpha\psi}{\psi}\bar{\gamma}^{ij}\partial_j K + 16\pi\left(\alpha\psi\right)\psi^{3-n} \bar{S}^i\) to the momentum constraint, where \(n\) is the ConformalMatterScale and \(\bar{S}^i=\psi^n S^i\) is the conformally-scaled momentum density.

Note that the \(\partial_j\bar{u}^{ij}\) term is not the full covariant divergence, but only the partial-derivatives part of it. The curved contribution to this term can be added together with the curved contribution to the flux divergence of the dynamic shift variable with the Elasticity::add_curved_sources function.

Also adds \(-\frac{1}{8}\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\) to the Hamiltonian constraint and \(\frac{7}{8}\alpha\psi^{-7}\bar{A}^{ij}\bar{A}_{ij}\) to the lapse equation.

See also
Xcts

◆ add_hamiltonian_sources()

template<int ConformalMatterScale>
void Xcts::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_minus_one 
)

Add the nonlinear source to the Hamiltonian constraint on a flat conformal background in Cartesian coordinates and with \(\bar{u}_{ij}=0=\beta^i\).

Adds \(\frac{1}{12}\psi^5 K^2 - 2\pi\psi^{5-n}\bar{\rho}\) where \(n\) is the ConformalMatterScale and \(\bar{\rho}=\psi^n\rho\) is the conformally-scaled energy density. Additional sources can be added with add_distortion_hamiltonian_sources and add_curved_hamiltonian_or_lapse_sources.

See also
Xcts

◆ add_lapse_sources()

template<int ConformalMatterScale>
void Xcts::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_minus_one,
const Scalar< DataVector > &  lapse_times_conformal_factor_minus_one 
)

Add the nonlinear source to the lapse equation on a flat conformal background in Cartesian coordinates and with \(\bar{u}_{ij}=0=\beta^i\).

Adds \((\alpha\psi)\left(\frac{5}{12}\psi^4 K^2 + 2\pi\psi^{4-n} \left(\bar{\rho} + 2\bar{S}\right)\right) + \psi^5 \left(\beta^i\partial_i K - \partial_t K\right)\) where \(n\) is the ConformalMatterScale and matter quantities are conformally-scaled. Additional sources can be added with add_distortion_hamiltonian_and_lapse_sources and add_curved_hamiltonian_or_lapse_sources.

See also
Xcts

◆ add_linearized_distortion_hamiltonian_and_lapse_sources()

void Xcts::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_minus_one,
const Scalar< DataVector > &  lapse_times_conformal_factor_minus_one,
const Scalar< DataVector > &  conformal_factor_correction,
const Scalar< DataVector > &  lapse_times_conformal_factor_correction 
)

The linearization of add_distortion_hamiltonian_and_lapse_sources

Note that this linearization is only w.r.t. \(\psi\) and \(\alpha\psi\).

◆ add_linearized_distortion_hamiltonian_sources()

void Xcts::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_minus_one,
const Scalar< DataVector > &  conformal_factor_correction 
)

The linearization of add_distortion_hamiltonian_sources

Note that this linearization is only w.r.t. \(\psi\).

◆ add_linearized_lapse_sources()

template<int ConformalMatterScale>
void Xcts::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_minus_one,
const Scalar< DataVector > &  lapse_times_conformal_factor_minus_one,
const Scalar< DataVector > &  conformal_factor_correction,
const Scalar< DataVector > &  lapse_times_conformal_factor_correction 
)

The linearization of add_lapse_sources

Note that this linearization is only w.r.t. \(\psi\) and \(\alpha\psi\). The linearization w.r.t. \(\beta^i\) is added in add_curved_linearized_momentum_sources or add_flat_cartesian_linearized_momentum_sources.

◆ extrinsic_curvature()

template<typename DataType >
void Xcts::extrinsic_curvature ( const gsl::not_null< tnsr::ii< DataType, 3 > * >  result,
const Scalar< DataType > &  conformal_factor,
const Scalar< DataType > &  lapse,
const tnsr::ii< DataType, 3 > &  conformal_metric,
const tnsr::II< DataType, 3 > &  longitudinal_shift_minus_dt_conformal_metric,
const Scalar< DataType > &  trace_extrinsic_curvature 
)

Extrinsic curvature computed from the conformal decomposition used in the XCTS system.

The extrinsic curvature decomposition is (see Eq. 3.113 in [13]):

\begin{equation} K_{ij} = A_{ij} + \frac{1}{3}\gamma_{ij}K = \frac{\psi^4}{2\lapse}\left((\bar{L}\beta)_{ij} - \bar{u}_{ij}\right) + \frac{\psi^4}{3} \bar{\gamma}_{ij} K \end{equation}

Parameters
resultoutput buffer for the extrinsic curvature
conformal_factorthe conformal factor \(\psi\)
lapsethe lapse \(\alpha\)
conformal_metricthe conformal metric \(\bar{\gamma}_{ij}\)
longitudinal_shift_minus_dt_conformal_metricthe term \((\bar{L}\beta)^{ij} - \bar{u}^{ij}\). Note that \((\bar{L}\beta)^{ij}\) is the conformal longitudinal shift, and \((\bar{L}\beta)^{ij}=\psi^4(L\beta)^{ij}\) (Eq. 3.98 in [13]). See also Xcts::longitudinal_operator.
trace_extrinsic_curvaturethe trace of the extrinsic curvature, \(K=\gamma^{ij}K_{ij}\). Note that it is a conformal invariant, \(K=\bar{K}\) (by choice).

◆ longitudinal_operator() [1/2]

template<typename DataType >
void Xcts::longitudinal_operator ( gsl::not_null< tnsr::II< DataType, 3 > * >  result,
const tnsr::I< DataType, 3 > &  shift,
const tnsr::iJ< DataType, 3 > &  deriv_shift,
const tnsr::II< DataType, 3 > &  inv_metric,
const tnsr::Ijj< DataType, 3 > &  christoffel_second_kind 
)

The longitudinal operator, or vector gradient, \((L\beta)^{ij}\).

Computes the longitudinal operator

\begin{equation} (L\beta)^{ij} = \nabla^i \beta^j + \nabla^j \beta^i - \frac{2}{3}\gamma^{ij}\nabla_k\beta^k \end{equation}

of a vector field \(\beta^i\), where \(\nabla\) denotes the covariant derivative w.r.t. the metric \(\gamma\) (see e.g. Eq. (3.50) in [13]). Note that in the XCTS equations the longitudinal operator is typically applied to conformal quantities and w.r.t. the conformal metric \(\bar{\gamma}\).

In terms of the symmetric "strain" quantity \(B_{ij}=\nabla_{(i}\gamma_{j)k}\beta^k\) the longitudinal operator is:

\begin{equation} (L\beta)^{ij} = 2\left(\gamma^{ik}\gamma^{jl} - \frac{1}{3} \gamma^{ij}\gamma^{kl}\right) B_{kl} \end{equation}

Note that the strain can be computed with Elasticity::strain.

◆ longitudinal_operator() [2/2]

template<typename DataType >
void Xcts::longitudinal_operator ( gsl::not_null< tnsr::II< DataType, 3 > * >  result,
const tnsr::ii< DataType, 3 > &  strain,
const tnsr::II< DataType, 3 > &  inv_metric 
)

The longitudinal operator, or vector gradient, \((L\beta)^{ij}\).

Computes the longitudinal operator

\begin{equation} (L\beta)^{ij} = \nabla^i \beta^j + \nabla^j \beta^i - \frac{2}{3}\gamma^{ij}\nabla_k\beta^k \end{equation}

of a vector field \(\beta^i\), where \(\nabla\) denotes the covariant derivative w.r.t. the metric \(\gamma\) (see e.g. Eq. (3.50) in [13]). Note that in the XCTS equations the longitudinal operator is typically applied to conformal quantities and w.r.t. the conformal metric \(\bar{\gamma}\).

In terms of the symmetric "strain" quantity \(B_{ij}=\nabla_{(i}\gamma_{j)k}\beta^k\) the longitudinal operator is:

\begin{equation} (L\beta)^{ij} = 2\left(\gamma^{ik}\gamma^{jl} - \frac{1}{3} \gamma^{ij}\gamma^{kl}\right) B_{kl} \end{equation}

Note that the strain can be computed with Elasticity::strain.

◆ longitudinal_operator_flat_cartesian()

template<typename DataType >
void Xcts::longitudinal_operator_flat_cartesian ( gsl::not_null< tnsr::II< DataType, 3 > * >  result,
const tnsr::ii< DataType, 3 > &  strain 
)

The conformal longitudinal operator \((L\beta)^{ij}\) on a flat conformal metric in Cartesian coordinates \(\gamma_{ij}=\delta_{ij}\).

See also
Xcts::longitudinal_operator