SpECTRE Documentation Coverage Report
Current view: top level - Elliptic/Systems/Xcts - Equations.hpp Hit Total Coverage
Commit: 3c072f0ce967e2e56649d3fa12aa2a0e4fe2a42e Lines: 13 14 92.9 %
Date: 2024-04-23 20:50:18
Legend: Lines: hit not hit

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

Generated by: LCOV version 1.14