Tensorium
Loading...
Searching...
No Matches
tensorium_RG::bssn Namespace Reference

Namespaces

namespace  constraint_monitor_detail
namespace  detail

Classes

struct  ConstraintMonitorStats
 Snapshot of max/L2 Hamiltonian error plus algebraic drifts. More...
struct  GaugeParameters
 Tunable coefficients for the Gamma-driver system. More...
struct  ADMVariables
 ADM \((3+1)\) variables evaluated at a spatial point. More...
struct  BSSNVariables
 Conformal variables \(\{\chi,\tilde{\gamma}_{ij},\tilde{\gamma}^{ij}\}\) derived from ADM data. More...
struct  SampleLocation
 Index-space and physical coordinates where the worst invariant violation occurred. More...
struct  InvariantStats
 Aggregated statistics computed over the interior domain. More...
struct  InvariantTolerances
 Thresholds against which assert_invariants compares the measured drifts. More...
struct  ProjectionConfig
 Fine-grained control over which invariants get enforced during projection. More...
struct  BSSNRHSWorkspace
 Stores RHS buffers for every evolved variable. More...
class  BSSNRKStepper
 Classical RK4 stepper that manages halo exchanges, geometry refresh, and diagnostics. More...

Functions

template<typename T>
ConstraintMonitorStats compute_constraint_monitor (const BSSNGridSoA< T > &G, const Field3D< T > &H, size_t padding=4)
 Aggregate Linf/L2 norms of the Hamiltonian constraint and track det/trace violations.
void print_constraint_monitor (const ConstraintMonitorStats &stats, const char *label="constraints")
 Convenience printf helper for debugging constraint convergence.
template<typename T>
static void compute_bssn_constraints (BSSNGridSoA< T > &G, const Field3D< T > *Ricci6, Field3D< T > &H_out, Field3D< T > M_out[3], Field3D< T > C_out[3], double r_min, double r_max, double xc, double yc, double zc)
 Fill the Hamiltonian \(H\), momentum \(M_i\), and Gamma \(C_i\) constraint fields.
static void print_constraint_norms (BSSNGridSoA< double > &G, const Field3D< double > &H, const Field3D< double > M[3], const Field3D< double > C[3], double r_min, double r_max, double xc, double yc, double zc)
 Convenience logger that prints Linf/L2 norms of \(H\), \(|\vec{M}|\), and \(|\vec{C}|\) to stdout.
size_t clamped_lower (size_t lower, size_t guard, size_t upper)
size_t clamped_upper (size_t upper, size_t guard, size_t lower)
template<typename T>
void compute_rhs_A_tilde (const BSSNGridSoA< T > &G, Field3D< T > rhs[6], size_t padding=4)
 Assemble \(\partial_t\tilde{A}_{ij}\) for each symmetric component.
template<typename T>
void compute_rhs_chi (const BSSNGridSoA< T > &G, Field3D< T > &rhs_chi, size_t padding=4)
 Fill the RHS buffer for \(\chi\) inside the padded interior region.
template<typename T>
void compute_rhs_Gamma (const BSSNGridSoA< T > &G, Field3D< T > rhs[3], size_t padding=4)
 Assemble \(\partial_t\tilde{\Gamma}^i\) inside the interior region.
template<typename T>
void compute_rhs_gamma_tilde (const BSSNGridSoA< T > &G, Field3D< T > rhs[6], size_t padding=4)
 Fill \(\partial_t\tilde{\gamma}_{ij}\) for all six symmetric components.
template<typename T>
void compute_rhs_alpha (const BSSNGridSoA< T > &G, Field3D< T > &rhs_alpha, size_t padding=4)
 Build \(\partial_t\alpha\) including advection and KO6 dissipation.
template<typename T>
void compute_rhs_beta (const BSSNGridSoA< T > &G, Field3D< T > rhs[3], const GaugeParameters< T > &params={}, size_t padding=4)
 Assemble \(\partial_t\beta^i = \beta^k\partial_k\beta^i + (3/4)B^i\) plus dissipation.
template<typename T>
void compute_rhs_B (const BSSNGridSoA< T > &G, const Field3D< T > rhs_Gamma[3], Field3D< T > rhs_B[3], const GaugeParameters< T > &params={}, size_t padding=4)
 Assemble \(\partial_t B^i\) given the previously computed \(\partial_t\tilde{\Gamma}^i\).
template<typename T>
void compute_rhs_K (const BSSNGridSoA< T > &G, Field3D< T > &rhs_K, size_t padding=4)
 Compute \(\partial_t K\) throughout the interior region.
template<typename T>
void compute_tildeGamma_contracted (BSSNGridSoA< T > &G)
 Refresh the evolved \(\tilde{\Gamma}^i\) from the metric inverse divergence.
template<typename T>
void metric_inverse_divergence (const BSSNGridSoA< T > &G, size_t i, size_t j, size_t k, T out[3])
 Evaluate \(\partial_j\tilde{\gamma}^{ij}\) using the 4th-order derivative operators.
template<typename T>
void compute_contracted_gamma_from_metric (BSSNGridSoA< T > &G, size_t i, size_t j, size_t k, T out[3])
 Convenience wrapper that produces \(\tilde{\Gamma}^i\) from the metric and stores the negative divergence as required by the BSSN definition.
double & chi_tolerance_override_storage ()
void set_chi_tolerance_override (double tol)
double chi_tolerance_override ()
template<typename T>
InvariantStats compute_invariant_stats (const BSSNGridSoA< T > &G, size_t padding=4, double r_min=0.0, double r_max=std::numeric_limits< double >::max(), double xc=0.0, double yc=0.0, double zc=0.0)
 Measure det, trace, metric-identity, Gamma-constraint, and minimum- \(\chi\) over the grid.
template<typename T>
InvariantTolerances compute_invariant_tolerances (const BSSNGridSoA< T > &G, double scheme_order=4.0, double safety_factor=4.0)
 Build tolerances that reflect the expected truncation error of the chosen stencil order.
template<typename T>
void assert_invariants (const BSSNGridSoA< T > &G, const char *stage, size_t padding, const InvariantTolerances &tol, bool throw_on_violation=true)
 Compare measured invariants against tolerances and throw/log when drifts exceed bounds.
template<typename T>
void assert_invariants (const BSSNGridSoA< T > &G, const char *stage, size_t padding=4, bool throw_on_violation=true)
template<typename T>
void project_bssn_state (BSSNGridSoA< T > &G, const ProjectionConfig &cfg={})
 Apply determinant renormalization, \(\tilde{A}\) trace-free projection, and \(\tilde{\Gamma}^i\) resynchronization.
template<typename T>
void project_bssn_after_update (BSSNGridSoA< T > &G, size_t padding=4)
 Convenience wrapper that enforces all invariants after every RK stage.
template<typename T>
ConstraintMonitorStats project_and_monitor (BSSNGridSoA< T > &G, Field3D< T > &H, Field3D< T > M[3], Field3D< T > C[3], double r_min, double r_max, double xc, double yc, double zc, size_t padding=4)
 Project the BSSN state and immediately recompute Hamiltonian, momentum, and Gamma constraints.
template<typename T>
void compute_ricci_bssn (BSSNGridSoA< T > &G, Field3D< T > *Ricci6, bool throw_on_violation=true)
 Populate the Ricci tensor cache \(R_{ij}\) throughout the interior grid.
template<typename Boundary, typename T>
void apply_halos_grid (BSSNGridSoA< T > &G)
 Apply the boundary functor to every scalar, vector, and tensor field in the grid.
template<typename T>
ADMVariables< T > split_3p1 (const tensorium::Vector< T > &X, const Metric< T > &metric)
 Evaluate the ADM 3+1 split of metric at spatial point X.
template<typename T>
BSSNVariables< T > adm_to_bssn (const tensorium::Tensor< T, 2 > &gamma_ij)
 Convert ADM spatial metric \(\gamma_{ij}\) into conformal variables \(\chi,\tilde{\gamma}_{ij}\).
template<typename T>
compute_dt_cfl (const BSSNGridSoA< T > &grid, T cfl_factor, size_t padding=4)
 Compute a timestep satisfying \(dt = \,\text{CFL}\times \min(\Delta x)/(|\vec{\beta}|+\alpha_{max})\).

Function Documentation

◆ adm_to_bssn()

template<typename T>
BSSNVariables< T > tensorium_RG::bssn::adm_to_bssn ( const tensorium::Tensor< T, 2 > & gamma_ij)
inline

Convert ADM spatial metric \(\gamma_{ij}\) into conformal variables \(\chi,\tilde{\gamma}_{ij}\).

◆ apply_halos_grid()

template<typename Boundary, typename T>
void tensorium_RG::bssn::apply_halos_grid ( BSSNGridSoA< T > & G)
inline

Apply the boundary functor to every scalar, vector, and tensor field in the grid.

Template Parameters
BoundaryFunctor type that exposes static void apply(Field3D<T>&, const GridDims&).

Mapping to code:

  • The first loop handles scalar gauges \(\alpha,\chi,K\).
  • The second loop applies the functor to each component of \(\beta^i, B^i, \tilde{\Gamma}^i\).
  • Symmetric tensors (metric, inverse, \f\tilde{A}_{ij}\f) reuse the packed 6-component arrays.
  • Cached geometry (Gamma_tilde) follows so the Ricci builder can be invoked immediately after halo application.
Note
Guard updates must precede any derivative or KO6 evaluation because the latter access up to three off-domain samples. Refer to BSSN_Grid for halo width requirements.

Referenced by tensorium_RG::init::binary_bowen_york_puncture_init(), tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), tensorium_RG::init::minkowski(), tensorium_RG::init::schwarzschild_isotropic(), and tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::step().

Here is the caller graph for this function:

◆ assert_invariants() [1/2]

template<typename T>
void tensorium_RG::bssn::assert_invariants ( const BSSNGridSoA< T > & G,
const char * stage,
size_t padding,
const InvariantTolerances & tol,
bool throw_on_violation = true )
inline

Compare measured invariants against tolerances and throw/log when drifts exceed bounds.

Parameters
stageText label identifying the caller (used in exception messages).
paddingSame meaning as in compute_invariant_stats.
tolTunable tolerances; use compute_invariant_tolerances for defaults.
throw_on_violationToggle between exceptions and stderr warnings.

References tensorium_RG::bssn::InvariantTolerances::chi_tol, compute_invariant_stats(), tensorium_RG::bssn::InvariantTolerances::det_tol, tensorium_RG::bssn::InvariantTolerances::gamma_tol, tensorium_RG::bssn::InvariantTolerances::metric_tol, and tensorium_RG::bssn::InvariantTolerances::trace_tol.

Referenced by assert_invariants(), tensorium_RG::init::binary_bowen_york_puncture_init(), tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), compute_bssn_constraints(), compute_ricci_bssn(), tensorium_RG::init::minkowski(), tensorium_RG::init::schwarzschild_isotropic(), and tensorium_RG::init::solve_lichnerowicz_u_SOR().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ assert_invariants() [2/2]

template<typename T>
void tensorium_RG::bssn::assert_invariants ( const BSSNGridSoA< T > & G,
const char * stage,
size_t padding = 4,
bool throw_on_violation = true )
inline

References assert_invariants(), and compute_invariant_tolerances().

Here is the call graph for this function:

◆ chi_tolerance_override()

double tensorium_RG::bssn::chi_tolerance_override ( )
inline

References chi_tolerance_override_storage().

Referenced by compute_invariant_tolerances().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ chi_tolerance_override_storage()

double & tensorium_RG::bssn::chi_tolerance_override_storage ( )
inline

Referenced by chi_tolerance_override(), and set_chi_tolerance_override().

Here is the caller graph for this function:

◆ clamped_lower()

size_t tensorium_RG::bssn::clamped_lower ( size_t lower,
size_t guard,
size_t upper )
inline

Referenced by compute_rhs_A_tilde(), compute_rhs_alpha(), compute_rhs_chi(), compute_rhs_Gamma(), compute_rhs_gamma_tilde(), and compute_rhs_K().

Here is the caller graph for this function:

◆ clamped_upper()

size_t tensorium_RG::bssn::clamped_upper ( size_t upper,
size_t guard,
size_t lower )
inline

Referenced by compute_rhs_A_tilde(), compute_rhs_alpha(), compute_rhs_chi(), compute_rhs_Gamma(), compute_rhs_gamma_tilde(), and compute_rhs_K().

Here is the caller graph for this function:

◆ compute_bssn_constraints()

template<typename T>
void tensorium_RG::bssn::compute_bssn_constraints ( BSSNGridSoA< T > & G,
const Field3D< T > * Ricci6,
Field3D< T > & H_out,
Field3D< T > M_out[3],
Field3D< T > C_out[3],
double r_min,
double r_max,
double xc,
double yc,
double zc )
inlinestatic

Fill the Hamiltonian \(H\), momentum \(M_i\), and Gamma \(C_i\) constraint fields.

Parameters
Ricci6Packed \(R_{ij}\) produced by compute_ricci_bssn.
H_out,M_out,C_outDestinations that must be allocated like the grid.
r_min,r_maxOptional radial filter to skip interior/exterior regions (useful around punctures).
  • Hamiltonian block maps gI_* multiplications to \(R=\gamma^{ij}R_{ij}\).
  • Momentum block forms \(M_i = \tilde{\gamma}^{jk}\tilde{D}_j\tilde{A}_{ki} + 6\tilde{A}_{ij}\partial^j\phi - \tfrac{2}{3}\partial_i K\) using the helper grad_6phi_from_chi.
  • Gamma constraint reuses metric_inverse_divergence to compare the evolved \(\tilde{\Gamma}^i\) against \(-\partial_j\tilde{\gamma}^{ij}\).

References assert_invariants(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dz(), metric_inverse_divergence(), tensorium_RG::Field3D< T >::ptr(), tensorium_RG::sym6_get(), tensorium_RG::sym6_index(), tensorium_RG::XX, tensorium_RG::XY, tensorium_RG::XZ, tensorium_RG::YY, tensorium_RG::YZ, and tensorium_RG::ZZ.

Referenced by tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), tensorium_RG::init::minkowski(), project_and_monitor(), tensorium_RG::init::schwarzschild_isotropic(), and tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::step().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_constraint_monitor()

template<typename T>
ConstraintMonitorStats tensorium_RG::bssn::compute_constraint_monitor ( const BSSNGridSoA< T > & G,
const Field3D< T > & H,
size_t padding = 4 )
inline

Aggregate Linf/L2 norms of the Hamiltonian constraint and track det/trace violations.

Parameters
paddingGuard cells excluded from the reduction to avoid halo artifacts.

References tensorium_RG::bssn::constraint_monitor_detail::det3(), tensorium_RG::Field3D< T >::idx(), tensorium_RG::bssn::ConstraintMonitorStats::l2_H, tensorium_RG::bssn::ConstraintMonitorStats::max_det_drift, tensorium_RG::bssn::ConstraintMonitorStats::max_H, tensorium_RG::bssn::ConstraintMonitorStats::max_trace_A, tensorium_RG::Field3D< T >::ptr(), tensorium_RG::bssn::ConstraintMonitorStats::samples, and tensorium_RG::sym6_get().

Referenced by project_and_monitor(), and tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::step().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_contracted_gamma_from_metric()

template<typename T>
void tensorium_RG::bssn::compute_contracted_gamma_from_metric ( BSSNGridSoA< T > & G,
size_t i,
size_t j,
size_t k,
T out[3] )
inline

Convenience wrapper that produces \(\tilde{\Gamma}^i\) from the metric and stores the negative divergence as required by the BSSN definition.

Note
Avoid calling right before constraint evaluation, otherwise \(C^i\) will trivially vanish.

References metric_inverse_divergence().

Here is the call graph for this function:

◆ compute_dt_cfl()

template<typename T>
T tensorium_RG::bssn::compute_dt_cfl ( const BSSNGridSoA< T > & grid,
T cfl_factor,
size_t padding = 4 )
inline

Compute a timestep satisfying \(dt = \,\text{CFL}\times \min(\Delta x)/(|\vec{\beta}|+\alpha_{max})\).

◆ compute_invariant_stats()

template<typename T>
InvariantStats tensorium_RG::bssn::compute_invariant_stats ( const BSSNGridSoA< T > & G,
size_t padding = 4,
double r_min = 0.0,
double r_max = std::numeric_limits<double>::max(),
double xc = 0.0,
double yc = 0.0,
double zc = 0.0 )
inline

Measure det, trace, metric-identity, Gamma-constraint, and minimum- \(\chi\) over the grid.

Parameters
paddingNumber of cells excluded near physical boundaries to avoid halo contamination.
r_min,r_maxOptional spherical mask to focus on physically relevant regions.
  • Determinant deviation: \(|\det\tilde{\gamma}-1|\).
  • Trace of \(\tilde{A}\): \(|\tilde{\gamma}^{ij}\tilde{A}_{ij}|\).
  • Metric identity: maximum norm of \(\tilde{\gamma}_{ik}\tilde{\gamma}^{kj}-\delta_i^{\ j}\).
  • Gamma constraint: \(|\tilde{\Gamma}^i + \partial_j\tilde{\gamma}^{ij}|\).
  • Minimum \(\chi\): used to detect slice stretching approaching the excision floor.

References tensorium_RG::bssn::InvariantStats::chi_location, tensorium_RG::bssn::detail::det3(), tensorium_RG::bssn::InvariantStats::det_location, tensorium_RG::bssn::InvariantStats::gamma_location, tensorium_RG::bssn::InvariantStats::l2_det_deviation, tensorium_RG::bssn::InvariantStats::l2_gamma_constraint, tensorium_RG::bssn::InvariantStats::l2_metric_identity, tensorium_RG::bssn::InvariantStats::l2_trace_A, tensorium_RG::bssn::InvariantStats::max_det_deviation, tensorium_RG::bssn::InvariantStats::max_gamma_constraint, tensorium_RG::bssn::InvariantStats::max_metric_identity, tensorium_RG::bssn::InvariantStats::max_trace_A, tensorium_RG::bssn::InvariantStats::mean_det_deviation, tensorium_RG::bssn::InvariantStats::mean_gamma_constraint, tensorium_RG::bssn::InvariantStats::mean_metric_identity, tensorium_RG::bssn::InvariantStats::mean_trace_A, metric_inverse_divergence(), tensorium_RG::bssn::InvariantStats::metric_location, tensorium_RG::bssn::InvariantStats::min_chi, tensorium_RG::bssn::InvariantStats::samples, tensorium_RG::sym6_get(), tensorium_RG::bssn::InvariantStats::trace_location, X(), and tensorium_RG::XX.

Referenced by assert_invariants().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_invariant_tolerances()

template<typename T>
InvariantTolerances tensorium_RG::bssn::compute_invariant_tolerances ( const BSSNGridSoA< T > & G,
double scheme_order = 4.0,
double safety_factor = 4.0 )
inline

Build tolerances that reflect the expected truncation error of the chosen stencil order.

The leading truncation error of a \(p\)th-order scheme scales like \((h/L)^p\). This routine uses \(h=\max(\Delta x,\Delta y,\Delta z)\) and \(L\) as the largest domain extent to produce resolution-aware tolerances, blended with \f\sqrt{\epsilon_{\text{machine}}}\f to respect floating-point limits.

References tensorium_RG::bssn::InvariantTolerances::chi_tol, chi_tolerance_override(), tensorium_RG::bssn::InvariantTolerances::det_tol, tensorium_RG::bssn::InvariantTolerances::gamma_tol, tensorium_RG::bssn::InvariantTolerances::metric_tol, and tensorium_RG::bssn::InvariantTolerances::trace_tol.

Referenced by assert_invariants().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_rhs_A_tilde()

template<typename T>
void tensorium_RG::bssn::compute_rhs_A_tilde ( const BSSNGridSoA< T > & G,
Field3D< T > rhs[6],
size_t padding = 4 )
inline

Assemble \(\partial_t\tilde{A}_{ij}\) for each symmetric component.

  • lie implements the \(\mathcal{L}_\beta\tilde{A}_{ij}\) term.
  • term_geom corresponds to \(\chi[-D_iD_j\alpha + \alpha R_{ij}]^{TF}\).
  • term_quad implements \(\alpha(K\tilde{A}_{ij} - 2\tilde{A}_{ik}\tilde{A}^k_{\ j})\).

References clamped_lower(), clamped_upper(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dx_upwind(), tensorium_RG::fd::Dxx(), tensorium_RG::fd::Dxy4(), tensorium_RG::fd::Dxz4(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dy_upwind(), tensorium_RG::fd::Dyy(), tensorium_RG::fd::Dyz4(), tensorium_RG::fd::Dz(), tensorium_RG::fd::Dz_upwind(), tensorium_RG::fd::Dzz(), tensorium_RG::fd::KO6(), tensorium_RG::Field3D< T >::ptr(), tensorium_RG::sym6_get(), and tensorium_RG::sym6_index().

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_rhs_alpha()

template<typename T>
void tensorium_RG::bssn::compute_rhs_alpha ( const BSSNGridSoA< T > & G,
Field3D< T > & rhs_alpha,
size_t padding = 4 )
inline

Build \(\partial_t\alpha\) including advection and KO6 dissipation.

References clamped_lower(), clamped_upper(), tensorium_RG::fd::Dx_upwind(), tensorium_RG::fd::Dy_upwind(), tensorium_RG::fd::Dz_upwind(), tensorium_RG::fd::KO6(), and tensorium_RG::Field3D< T >::ptr().

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_rhs_B()

template<typename T>
void tensorium_RG::bssn::compute_rhs_B ( const BSSNGridSoA< T > & G,
const Field3D< T > rhs_Gamma[3],
Field3D< T > rhs_B[3],
const GaugeParameters< T > & params = {},
size_t padding = 4 )
inline

Assemble \(\partial_t B^i\) given the previously computed \(\partial_t\tilde{\Gamma}^i\).

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the caller graph for this function:

◆ compute_rhs_beta()

template<typename T>
void tensorium_RG::bssn::compute_rhs_beta ( const BSSNGridSoA< T > & G,
Field3D< T > rhs[3],
const GaugeParameters< T > & params = {},
size_t padding = 4 )
inline

Assemble \(\partial_t\beta^i = \beta^k\partial_k\beta^i + (3/4)B^i\) plus dissipation.

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the caller graph for this function:

◆ compute_rhs_chi()

template<typename T>
void tensorium_RG::bssn::compute_rhs_chi ( const BSSNGridSoA< T > & G,
Field3D< T > & rhs_chi,
size_t padding = 4 )
inline

Fill the RHS buffer for \(\chi\) inside the padded interior region.

Parameters
paddingNumber of guard cells skipped on each side before looping.

Mapping to code.

  • d_chi corresponds to the advective term \(\beta^i\partial_i\chi\).
  • div_beta uses centered derivatives to compute \(\partial_i\beta^i\).
  • The scalar multiplier T(2/3)*chi*(alpha*K - div_beta) encodes the source term.
  • KO6 call adds dissipation with \(\sigma=0.1\).

References clamped_lower(), clamped_upper(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dx_upwind(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dy_upwind(), tensorium_RG::fd::Dz(), tensorium_RG::fd::Dz_upwind(), tensorium_RG::fd::KO6(), and tensorium_RG::Field3D< T >::ptr().

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_rhs_Gamma()

template<typename T>
void tensorium_RG::bssn::compute_rhs_Gamma ( const BSSNGridSoA< T > & G,
Field3D< T > rhs[3],
size_t padding = 4 )
inline

Assemble \(\partial_t\tilde{\Gamma}^i\) inside the interior region.

Parameters
rhsArray of 3 fields that will store the RHS.
  • adv = \(\beta^k\partial_k\tilde{\Gamma}^i\) uses upwind stencils.
  • gamma_beta and stretch implement the Lie derivative terms \(-\tilde{\Gamma}^k\partial_k\beta^i\) and \(\tfrac{2}{3}\tilde{\Gamma}^i\partial_k\beta^k\).
  • hess_term and grad_div_term build the shift Laplacian pieces.
  • A_grad_alpha and source encode \(-2\tilde{A}^{ij}\partial_j\alpha\) and the Ricci coupling.

References clamped_lower(), clamped_upper(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dx_upwind(), tensorium_RG::fd::Dxx(), tensorium_RG::fd::Dxy4(), tensorium_RG::fd::Dxz4(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dy_upwind(), tensorium_RG::fd::Dyy(), tensorium_RG::fd::Dyz4(), tensorium_RG::fd::Dz(), tensorium_RG::fd::Dz_upwind(), tensorium_RG::fd::Dzz(), tensorium_RG::fd::KO6(), tensorium_RG::Field3D< T >::ptr(), tensorium_RG::sym6_get(), and tensorium_RG::sym6_index().

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_rhs_gamma_tilde()

template<typename T>
void tensorium_RG::bssn::compute_rhs_gamma_tilde ( const BSSNGridSoA< T > & G,
Field3D< T > rhs[6],
size_t padding = 4 )
inline

Fill \(\partial_t\tilde{\gamma}_{ij}\) for all six symmetric components.

References clamped_lower(), clamped_upper(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dx_upwind(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dy_upwind(), tensorium_RG::fd::Dz(), tensorium_RG::fd::Dz_upwind(), tensorium_RG::fd::KO6(), tensorium_RG::Field3D< T >::ptr(), and tensorium_RG::sym6_index().

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_rhs_K()

template<typename T>
void tensorium_RG::bssn::compute_rhs_K ( const BSSNGridSoA< T > & G,
Field3D< T > & rhs_K,
size_t padding = 4 )
inline

Compute \(\partial_t K\) throughout the interior region.

References clamped_lower(), clamped_upper(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dx_upwind(), tensorium_RG::fd::Dxx(), tensorium_RG::fd::Dxy4(), tensorium_RG::fd::Dxz4(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dy_upwind(), tensorium_RG::fd::Dyy(), tensorium_RG::fd::Dyz4(), tensorium_RG::fd::Dz(), tensorium_RG::fd::Dz_upwind(), tensorium_RG::fd::Dzz(), tensorium_RG::fd::KO6(), tensorium_RG::Field3D< T >::ptr(), tensorium_RG::sym6_get(), and tensorium_RG::sym6_index().

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::evaluate_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_ricci_bssn()

template<typename T>
void tensorium_RG::bssn::compute_ricci_bssn ( BSSNGridSoA< T > & G,
Field3D< T > * Ricci6,
bool throw_on_violation = true )

Populate the Ricci tensor cache \(R_{ij}\) throughout the interior grid.

Parameters
GGrid supplying \(\chi\), \(\tilde{\gamma}_{ij}\), \(\tilde{\gamma}^{ij}\), and cached \(\tilde{\Gamma}^i\).
[out]Ricci6Packed symmetric destination fields.
throw_on_violationIf true, assert_invariants raises when the projector detects determinant or trace drifts beyond tolerance.

Mapping to code.

  • dg lambda → \(\partial_m\tilde{\gamma}_{ab}\) samples.
  • Gijk array → \(\tilde{\Gamma}^i_{\ jk}\) assembled from inverse metric contractions.
  • cov_dchi block → \(\tilde{D}_i\tilde{D}_j\chi\) computed via subtracting connection terms from second derivatives.
  • R_chi / R_tilde blocks implement the formulas quoted above and store the sum.
  • dGt corresponds to \(\partial_i\tilde{\Gamma}^k\) entering the \(\tilde{R}_{ij}\) expression.
Warning
Ricci evaluation assumes halos contain valid data for ±2 offsets. Call apply_halos_grid before invoking this routine.

References assert_invariants(), tensorium_RG::fd::Dx(), tensorium_RG::fd::Dxx(), tensorium_RG::fd::Dxy4(), tensorium_RG::fd::Dxz4(), tensorium_RG::fd::Dy(), tensorium_RG::fd::Dyy(), tensorium_RG::fd::Dyz4(), tensorium_RG::fd::Dz(), tensorium_RG::fd::Dzz(), tensorium_RG::Field3D< T >::ptr(), tensorium_RG::sym6_get(), tensorium_RG::sym6_index(), tensorium_RG::XX, tensorium_RG::XY, tensorium_RG::XZ, tensorium_RG::YY, tensorium_RG::YZ, and tensorium_RG::ZZ.

Referenced by tensorium_RG::init::binary_bowen_york_puncture_init(), tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::ensure_geometry(), tensorium_RG::init::minkowski(), and tensorium_RG::init::schwarzschild_isotropic().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_tildeGamma_contracted()

template<typename T>
void tensorium_RG::bssn::compute_tildeGamma_contracted ( BSSNGridSoA< T > & G)
inline

Refresh the evolved \(\tilde{\Gamma}^i\) from the metric inverse divergence.

Mapping to code:

  1. Sample metric_inverse_divergence (see BSSNGamma.hpp).
  2. Store \(\tilde{\Gamma}^i = -\partial_j\tilde{\gamma}^{ij}\) in the SoA vector fields.
    See also
    metric_inverse_divergence

References metric_inverse_divergence(), and tensorium_RG::XX.

Referenced by tensorium_RG::init::binary_bowen_york_puncture_init(), tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), tensorium_RG::init::minkowski(), tensorium_RG::init::schwarzschild_isotropic(), and tensorium_RG::init::solve_lichnerowicz_u_SOR().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ metric_inverse_divergence()

template<typename T>
void tensorium_RG::bssn::metric_inverse_divergence ( const BSSNGridSoA< T > & G,
size_t i,
size_t j,
size_t k,
T out[3] )
inline

Evaluate \(\partial_j\tilde{\gamma}^{ij}\) using the 4th-order derivative operators.

Called during Gamma-constraint monitoring as well as projection steps to re-synchronize the evolved \(\tilde{\Gamma}^i\) with the metric after renormalization.

References tensorium_RG::bssn::detail::diff_gamma_tilde_inv_component().

Referenced by compute_bssn_constraints(), compute_contracted_gamma_from_metric(), compute_invariant_stats(), and compute_tildeGamma_contracted().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ print_constraint_monitor()

void tensorium_RG::bssn::print_constraint_monitor ( const ConstraintMonitorStats & stats,
const char * label = "constraints" )
inline

◆ print_constraint_norms()

void tensorium_RG::bssn::print_constraint_norms ( BSSNGridSoA< double > & G,
const Field3D< double > & H,
const Field3D< double > M[3],
const Field3D< double > C[3],
double r_min,
double r_max,
double xc,
double yc,
double zc )
inlinestatic

Convenience logger that prints Linf/L2 norms of \(H\), \(|\vec{M}|\), and \(|\vec{C}|\) to stdout.

References tensorium_RG::Field3D< T >::ptr().

Referenced by tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), tensorium_RG::init::minkowski(), and tensorium_RG::init::schwarzschild_isotropic().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ project_and_monitor()

template<typename T>
ConstraintMonitorStats tensorium_RG::bssn::project_and_monitor ( BSSNGridSoA< T > & G,
Field3D< T > & H,
Field3D< T > M[3],
Field3D< T > C[3],
double r_min,
double r_max,
double xc,
double yc,
double zc,
size_t padding = 4 )
inline

Project the BSSN state and immediately recompute Hamiltonian, momentum, and Gamma constraints.

Parameters
r_min,r_maxRadial window for constraint evaluation (useful for puncture data).
paddingInterior padding forwarded to project_bssn_after_update and monitor kernels.
Returns
Aggregated statistics of \(H\) gathered by compute_constraint_monitor.

References compute_bssn_constraints(), compute_constraint_monitor(), and project_bssn_after_update().

Here is the call graph for this function:

◆ project_bssn_after_update()

template<typename T>
void tensorium_RG::bssn::project_bssn_after_update ( BSSNGridSoA< T > & G,
size_t padding = 4 )
inline

Convenience wrapper that enforces all invariants after every RK stage.

References tensorium_RG::bssn::ProjectionConfig::padding, tensorium_RG::bssn::ProjectionConfig::project_A_tilde, project_bssn_state(), tensorium_RG::bssn::ProjectionConfig::recompute_inverse, tensorium_RG::bssn::ProjectionConfig::renormalize_metric, and tensorium_RG::bssn::ProjectionConfig::resync_contracted_gamma.

Referenced by tensorium_RG::bssn::BSSNRKStepper< T, Boundary >::ensure_geometry(), and project_and_monitor().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ project_bssn_state()

template<typename T>
void tensorium_RG::bssn::project_bssn_state ( BSSNGridSoA< T > & G,
const ProjectionConfig & cfg = {} )
inline

Apply determinant renormalization, \(\tilde{A}\) trace-free projection, and \(\tilde{\Gamma}^i\) resynchronization.

Parameters
cfgControls which projections are executed; defaults enforce all invariants inside the interior padding region.
  • Metric renormalization rescales \(\tilde{\gamma}_{ij}\) by \(\det(\tilde{\gamma})^{-1/3}\) so the new determinant equals 1.
  • Trace-free projection subtracts \(\tfrac{1}{3}\tilde{\gamma}_{ij}\tilde{\gamma}^{mn}\tilde{A}_{mn}\) from \(\tilde{A}_{ij}\).
  • Optionally recompute \(\tilde{\gamma}^{ij}\) from the renormalized metric and resynchronize \(\tilde{\Gamma}^i\) via metric_inverse_divergence.

Referenced by tensorium_RG::init::binary_bowen_york_puncture_init(), tensorium_RG::init::binary_schwarzschild_isotropic_2centers(), tensorium_RG::init::minkowski(), project_bssn_after_update(), tensorium_RG::init::schwarzschild_isotropic(), and tensorium_RG::init::solve_lichnerowicz_u_SOR().

Here is the caller graph for this function:

◆ set_chi_tolerance_override()

void tensorium_RG::bssn::set_chi_tolerance_override ( double tol)
inline

References chi_tolerance_override_storage().

Here is the call graph for this function:

◆ split_3p1()

template<typename T>
ADMVariables< T > tensorium_RG::bssn::split_3p1 ( const tensorium::Vector< T > & X,
const Metric< T > & metric )
inline

Evaluate the ADM 3+1 split of metric at spatial point X.

References tensorium_RG::bssn::ADMVariables< T >::alpha, tensorium_RG::bssn::ADMVariables< T >::beta, tensorium_RG::Metric< T >::BSSN(), tensorium_RG::bssn::ADMVariables< T >::gamma_ij, and X().

Here is the call graph for this function: