LORENE
|
Binary black holes system. More...
#include <isol_hor.h>
Public Member Functions | |
Single_hor (Map_af &mpi) | |
Standard constructor. More... | |
Single_hor (const Single_hor &) | |
Copy constructor. More... | |
Single_hor (Map_af &mp, FILE *fich) | |
Constructor from a binary file. More... | |
virtual | ~Single_hor () |
Destructor. More... | |
void | operator= (const Single_hor &) |
Assignment to another Single_hor. More... | |
const Map_af & | get_mp () const |
Returns the mapping (readonly). More... | |
Map_af & | set_mp () |
Read/write of the mapping. More... | |
double | get_radius () const |
Returns the radius of the horizon. More... | |
void | set_radius (double rad) |
Sets the radius of the horizon to rad . More... | |
double | get_omega () const |
Returns the angular velocity. More... | |
void | set_omega (double ome) |
Sets the angular velocity to ome . More... | |
const Scalar & | get_n_auto () const |
Lapse function . More... | |
const Scalar & | get_n_comp () const |
Lapse function . More... | |
const Scalar & | get_nn () const |
Lapse function . More... | |
const Scalar & | get_psi_auto () const |
Conformal factor . More... | |
const Scalar & | get_psi_comp () const |
Conformal factor . More... | |
const Scalar & | get_psi () const |
Conformal factor . More... | |
const Scalar & | get_psi4 () const |
Conformal factor . More... | |
const Vector & | get_dn () const |
Covariant derivative of the lapse function . More... | |
const Vector & | get_dpsi () const |
Covariant derivative with respect to the flat metric of the conformal factor . More... | |
const Vector & | get_beta_auto () const |
Shift function . More... | |
const Vector & | get_beta_comp () const |
Shift function . More... | |
const Vector & | get_beta () const |
Shift function . More... | |
const Sym_tensor & | get_aa_auto () const |
Conformal representation of the traceless part of the extrinsic curvature: More... | |
const Sym_tensor & | get_aa_comp () const |
Conformal representation of the traceless part of the extrinsic curvature: More... | |
const Sym_tensor & | get_aa () const |
Conformal representation of the traceless part of the extrinsic curvature: More... | |
const Metric & | get_tgam () const |
Conformal metric . More... | |
const Metric & | get_gam () const |
metric More... | |
const Sym_tensor & | get_k_dd () const |
k_dd More... | |
const Scalar | get_decouple () const |
Returns the function used to construct tkij_auto from tkij_tot . More... | |
void | n_comp_import (const Single_hor &comp) |
Imports the part of N due to the companion hole comp . More... | |
void | psi_comp_import (const Single_hor &comp) |
Imports the part of due to the companion hole comp . More... | |
void | beta_comp_import (const Single_hor &comp) |
Imports the part of due to the companion hole comp . More... | |
double | viriel_seul () const |
Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively and N . More... | |
void | init_bhole () |
Sets the values of the fields to : More... | |
void | init_met_trK () |
Sets the 3-metric tilde to the flat metric and gamt_point, trK and trK_point to zero. More... | |
void | init_bhole_seul () |
Initiates for a single black hole. More... | |
void | set_psi_auto (const Scalar &psi_in) |
Sets the conformal factor relating the physical metric to the conformal one: . More... | |
void | set_n_auto (const Scalar &nn_in) |
Sets the lapse. More... | |
void | set_beta_auto (const Scalar &shift_in) |
Sets the shift. More... | |
void | set_aa_auto (const Scalar &aa_auto_in) |
Sets aa_auto. More... | |
void | set_aa_comp (const Scalar &aa_comp_in) |
Sets aa_comp. More... | |
void | set_aa (const Scalar &aa_in) |
Sets aa. More... | |
const Scalar | b_tilde () const |
Radial component of the shift with respect to the conformal metric. More... | |
const Scalar | darea_hor () const |
Element of area of the horizon. More... | |
double | area_hor () const |
Area of the horizon. More... | |
double | radius_hor () const |
Radius of the horizon. More... | |
double | ang_mom_hor () const |
Angular momentum (modulo) More... | |
double | mass_hor () const |
Mass computed at the horizon. More... | |
double | kappa_hor () const |
Surface gravity. More... | |
double | omega_hor () const |
Orbital velocity. More... | |
double | ang_mom_adm () const |
ADM angular Momentum. More... | |
Scalar | expansion () const |
Expansion of the outgoing null normal ( ) More... | |
const Valeur | boundary_psi_app_hor () const |
Neumann boundary condition for. More... | |
Component x of boundary value of f beta f const Valeur | boundary_beta_x (double om_orb, double om_loc) const |
Component y of boundary value of f beta f const Valeur | boundary_beta_y (double om_orb, double om_loc) const |
Component z of boundary value of f beta f const Valeur | boundary_beta_z () const |
double | regularisation (const Vector &shift_auto, const Vector &shift_comp, double ang_vel) |
Corrects shift_auto in such a way that the total is equal to zero in the horizon, which should ensure the regularity of . More... | |
double | regularise_one () |
Corrects the shift in the innermost shell, so that it remains and that equals zero on the horizon. More... | |
virtual void | sauve (FILE *fich) const |
Total or partial saves in a binary file. More... | |
Public Attributes | |
Dirichlet boundary condition for c N f partial_r N a | N |
Neumann boundary condition on nn f partial_r N a | N |
Protected Member Functions | |
void | del_deriv () const |
Deletes all the derived quantities. More... | |
void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. More... | |
Protected Attributes | |
Map_af & | mp |
Affine mapping. More... | |
int | nz |
Number of zones. More... | |
double | radius |
Radius of the horizon in LORENE's units. More... | |
double | omega |
Angular velocity in LORENE's units. More... | |
double | regul |
Intensity of the correction on the shift vector. More... | |
Scalar | n_auto |
Lapse function . More... | |
Scalar | n_comp |
Lapse function . More... | |
Scalar | nn |
Lapse function . More... | |
Scalar | psi_auto |
Conformal factor . More... | |
Scalar | psi_comp |
Conformal factor . More... | |
Scalar | psi |
Conformal factor . More... | |
Scalar * | p_psi4 |
Conformal factor . More... | |
Vector | dn |
Covariant derivative of the lapse with respect to the flat metric . More... | |
Vector | dpsi |
Covariant derivative of the conformal factor . More... | |
Vector | beta_auto |
Shift function . More... | |
Vector | beta_comp |
Shift function . More... | |
Vector | beta |
Shift function . More... | |
Metric * | p_gam |
Spatial metric . More... | |
Sym_tensor | aa_auto |
Components of the conformal representation of the traceless part of the extrinsic curvature: More... | |
Sym_tensor | aa_comp |
Components of the conformal representation of the traceless part of the extrinsic curvature: More... | |
Sym_tensor | aa |
Components of the conformal representation of the traceless part of the extrinsic curvature: More... | |
Sym_tensor * | p_k_dd |
Components of the extrinsic curvature: More... | |
Metric | tgam |
3 metric tilde More... | |
Metric_flat | ff |
3 metric flat More... | |
Sym_tensor | hh |
Deviation metric. More... | |
Sym_tensor | gamt_point |
Time derivative of the 3-metric tilde. More... | |
Scalar | trK |
Trace of the extrinsic curvature. More... | |
Scalar | trK_point |
Time derivative of the trace of the extrinsic curvature. More... | |
Scalar | decouple |
Function used to construct from the total . More... | |
Friends | |
class | Bin_hor |
Binary black holes system.
()
This class is intended for dealing with binary black holes configurations in the conformaly flat approximation.
Definition at line 894 of file isol_hor.h.
Lorene::Single_hor::Single_hor | ( | Map_af & | mpi | ) |
Standard constructor.
mpi | affine mapping |
Definition at line 73 of file single_hor.C.
References hh, set_der_0x0(), and Lorene::Tensor::set_etat_zero().
Lorene::Single_hor::Single_hor | ( | const Single_hor & | singlehor_in | ) |
Lorene::Single_hor::Single_hor | ( | Map_af & | mp, |
FILE * | fich | ||
) |
Constructor from a binary file.
mpi | affine mapping |
fich | file containing the saved isol_hor |
partial_read | indicates whether the full object must be read in file or whether the final construction is devoted to a constructor of a derived class |
Definition at line 133 of file single_hor.C.
References Lorene::Metric::con(), Lorene::Metric_flat::con(), ff, Lorene::fread_be(), gamt_point, Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), hh, mp, omega, set_der_0x0(), tgam, trK, and trK_point.
|
virtual |
double Lorene::Single_hor::ang_mom_adm | ( | ) | const |
ADM angular Momentum.
Definition at line 177 of file single_param.C.
References Lorene::Metric::cov(), get_gam(), get_k_dd(), Lorene::Map_af::integrale_surface_infini(), mp, Lorene::Scalar::mult_rsint(), and trK.
double Lorene::Single_hor::ang_mom_hor | ( | ) | const |
Angular momentum (modulo)
Definition at line 110 of file single_param.C.
References Lorene::contract(), darea_hor(), ff, get_gam(), get_k_dd(), Lorene::Metric::get_mp(), Lorene::Metric_flat::get_triad(), Lorene::Map_af::integrale_surface(), mp, radius, Lorene::Vector::set(), and Lorene::Scalar::std_spectral_base().
double Lorene::Single_hor::area_hor | ( | ) | const |
Area of the horizon.
Definition at line 91 of file single_param.C.
References darea_hor(), Lorene::Map_af::integrale_surface(), mp, Lorene::Scalar::raccord(), and radius.
const Scalar Lorene::Single_hor::b_tilde | ( | ) | const |
Radial component of the shift with respect to the conformal metric.
Definition at line 71 of file single_param.C.
References beta, Lorene::contract(), Lorene::Tensor::down(), Lorene::Metric::radial_vect(), and tgam.
void Lorene::Single_hor::beta_comp_import | ( | const Single_hor & | comp | ) |
Imports the part of due to the companion hole comp
.
The total is then calculated.
Definition at line 468 of file single_hor.C.
References beta, beta_auto, beta_comp, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), and mp.
const Valeur Lorene::Single_hor::boundary_psi_app_hor | ( | ) | const |
Neumann boundary condition for.
Definition at line 72 of file single_bound.C.
References Lorene::contract(), Lorene::Scalar::derive_cov(), Lorene::Vector::divergence(), ff, Lorene::Map::flat_met_spher(), Lorene::Mg3d::get_angu(), Lorene::Map::get_bvect_spher(), get_k_dd(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, psi, Lorene::Metric::radial_vect(), Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Scalar::std_spectral_base(), tgam, trK, and Lorene::Scalar::val_grid_point().
const Scalar Lorene::Single_hor::darea_hor | ( | ) | const |
Element of area of the horizon.
Definition at line 80 of file single_param.C.
References get_gam(), Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
protected |
Deletes all the derived quantities.
Definition at line 224 of file single_hor.C.
References p_gam, p_k_dd, p_psi4, and set_der_0x0().
Scalar Lorene::Single_hor::expansion | ( | ) | const |
Expansion of the outgoing null normal ( )
Definition at line 193 of file single_param.C.
References Lorene::contract(), get_gam(), get_k_dd(), and trK.
const Sym_tensor & Lorene::Single_hor::get_aa | ( | ) | const |
Conformal representation of the traceless part of the extrinsic curvature:
Definition at line 337 of file single_hor.C.
References aa.
const Sym_tensor & Lorene::Single_hor::get_aa_auto | ( | ) | const |
Conformal representation of the traceless part of the extrinsic curvature:
Definition at line 327 of file single_hor.C.
References aa_auto.
const Sym_tensor & Lorene::Single_hor::get_aa_comp | ( | ) | const |
Conformal representation of the traceless part of the extrinsic curvature:
Definition at line 332 of file single_hor.C.
References aa_comp.
const Vector & Lorene::Single_hor::get_beta | ( | ) | const |
const Vector & Lorene::Single_hor::get_beta_auto | ( | ) | const |
const Vector & Lorene::Single_hor::get_beta_comp | ( | ) | const |
|
inline |
Returns the function used to construct tkij_auto
from tkij_tot
.
Definition at line 1146 of file isol_hor.h.
References decouple.
const Vector & Lorene::Single_hor::get_dn | ( | ) | const |
Covariant derivative of the lapse function .
Definition at line 302 of file single_hor.C.
References dn.
const Vector & Lorene::Single_hor::get_dpsi | ( | ) | const |
Covariant derivative with respect to the flat metric of the conformal factor .
Definition at line 307 of file single_hor.C.
References dpsi.
const Metric & Lorene::Single_hor::get_gam | ( | ) | const |
metric
Definition at line 342 of file single_hor.C.
References Lorene::Metric::cov(), get_psi4(), p_gam, and tgam.
const Sym_tensor & Lorene::Single_hor::get_k_dd | ( | ) | const |
k_dd
Definition at line 351 of file single_hor.C.
References aa, get_gam(), get_psi4(), p_k_dd, Lorene::Tensor::std_spectral_base(), trK, and Lorene::Tensor::up_down().
|
inline |
const Scalar & Lorene::Single_hor::get_n_auto | ( | ) | const |
const Scalar & Lorene::Single_hor::get_n_comp | ( | ) | const |
const Scalar & Lorene::Single_hor::get_nn | ( | ) | const |
|
inline |
const Scalar & Lorene::Single_hor::get_psi | ( | ) | const |
const Scalar & Lorene::Single_hor::get_psi4 | ( | ) | const |
Conformal factor .
Definition at line 291 of file single_hor.C.
References p_psi4, Lorene::pow(), psi, and Lorene::Scalar::std_spectral_base().
const Scalar & Lorene::Single_hor::get_psi_auto | ( | ) | const |
const Scalar & Lorene::Single_hor::get_psi_comp | ( | ) | const |
|
inline |
|
inline |
void Lorene::Single_hor::init_bhole | ( | ) |
Sets the values of the fields to :
n_auto
n_comp
psi_auto
psi_comp
a being the radius of the hole, the other fields being set to zero.
Definition at line 487 of file single_hor.C.
References Lorene::Scalar::annule(), beta, beta_auto, beta_comp, Lorene::Scalar::derive_cov(), dn, dpsi, ff, Lorene::Map::get_bvect_spher(), mp, n_auto, n_comp, nn, psi, psi_auto, psi_comp, Lorene::Map::r, Lorene::Scalar::raccord(), radius, Lorene::Coord::set(), Lorene::Vector::set(), Lorene::Scalar::set_dzpuis(), Lorene::Tensor::set_etat_zero(), and Lorene::Scalar::std_spectral_base().
void Lorene::Single_hor::init_bhole_seul | ( | ) |
Initiates for a single black hole.
WARNING It supposes that the boost is zero and should only be used for an isolated black hole..
Definition at line 551 of file single_hor.C.
References Lorene::Scalar::annule(), beta, beta_auto, beta_comp, Lorene::Scalar::derive_cov(), dn, dpsi, ff, Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), mp, n_auto, n_comp, nn, psi, psi_auto, psi_comp, Lorene::Map::r, Lorene::Scalar::raccord(), radius, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::set_outer_boundary(), and Lorene::Scalar::std_spectral_base().
void Lorene::Single_hor::init_met_trK | ( | ) |
Sets the 3-metric tilde to the flat metric and gamt_point, trK and trK_point to zero.
Definition at line 539 of file single_hor.C.
References Lorene::Map::flat_met_spher(), gamt_point, mp, Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), tgam, trK, and trK_point.
double Lorene::Single_hor::kappa_hor | ( | ) | const |
Surface gravity.
Definition at line 149 of file single_param.C.
References ang_mom_hor(), Lorene::pow(), radius_hor(), and Lorene::sqrt().
double Lorene::Single_hor::mass_hor | ( | ) | const |
Mass computed at the horizon.
Definition at line 138 of file single_param.C.
References ang_mom_hor(), Lorene::pow(), radius_hor(), and Lorene::sqrt().
void Lorene::Single_hor::n_comp_import | ( | const Single_hor & | comp | ) |
Imports the part of N due to the companion hole comp
.
The total N is then calculated.
It also imports the covariant derivative of N and construct the total .
Definition at line 393 of file single_hor.C.
References Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_cov(), ff, Lorene::Map::get_bvect_cart(), Lorene::Scalar::import(), mp, n_auto, n_comp, nn, Lorene::Tensor::set_etat_qcq(), and Lorene::Scalar::std_spectral_base().
double Lorene::Single_hor::omega_hor | ( | ) | const |
Orbital velocity.
Definition at line 163 of file single_param.C.
References ang_mom_hor(), Lorene::pow(), radius_hor(), and Lorene::sqrt().
void Lorene::Single_hor::operator= | ( | const Single_hor & | singlehor_in | ) |
void Lorene::Single_hor::psi_comp_import | ( | const Single_hor & | comp | ) |
Imports the part of due to the companion hole comp
.
The total is then calculated.
It also imports the covariant derivative of and construct the total .
Definition at line 427 of file single_hor.C.
References Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_cov(), ff, Lorene::Map::get_bvect_cart(), Lorene::Scalar::import(), mp, psi, psi_auto, psi_comp, Lorene::Tensor::set_etat_qcq(), and Lorene::Scalar::std_spectral_base().
double Lorene::Single_hor::radius_hor | ( | ) | const |
Radius of the horizon.
Definition at line 100 of file single_param.C.
References area_hor(), and Lorene::pow().
double Lorene::Single_hor::regularisation | ( | const Vector & | shift_auto, |
const Vector & | shift_comp, | ||
double | ang_vel | ||
) |
Corrects shift_auto
in such a way that the total is equal to zero in the horizon, which should ensure the regularity of .
WARNING : this should only be used for a black hole in a binary system Bin_hor
.
comp | [input]: the part of generated by the companion hole. |
Definition at line 62 of file single_regul.C.
References Lorene::Cmp::annule(), Lorene::Tensor::annule_domain(), Lorene::Scalar::annule_hard(), beta_auto, Lorene::Valeur::c, Lorene::Vector::change_triad(), Lorene::Valeur::coef_i(), Lorene::diffrelmax(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Map::get_rot_phi(), Lorene::Scalar::get_spectral_va(), Lorene::Tensor::get_triad(), Lorene::Scalar::import(), Lorene::norme(), nz, Lorene::pow(), Lorene::Vector::set(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Vector::std_spectral_base().
double Lorene::Single_hor::regularise_one | ( | ) |
Corrects the shift in the innermost shell, so that it remains and that equals zero on the horizon.
return the relative difference between the shift before and after the regularisation.
WARNING this should only be used for an isolated black hole.
Definition at line 156 of file single_regul.C.
References beta, Lorene::Vector::change_triad(), Lorene::Valeur::coef_i(), Lorene::Map::get_bvect_cart(), mp, omega, Lorene::Vector::set(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Map::x, and Lorene::Map::y.
|
virtual |
Total or partial saves in a binary file.
fich | binary file |
partial_save | indicates whether the whole object must be saved. |
Definition at line 247 of file single_hor.C.
References beta_auto, Lorene::Metric::con(), Lorene::fwrite_be(), gamt_point, n_auto, omega, psi_auto, Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), Lorene::Tensor_sym::sauve(), tgam, trK, and trK_point.
void Lorene::Single_hor::set_aa | ( | const Scalar & | aa_in | ) |
void Lorene::Single_hor::set_aa_auto | ( | const Scalar & | aa_auto_in | ) |
void Lorene::Single_hor::set_aa_comp | ( | const Scalar & | aa_comp_in | ) |
void Lorene::Single_hor::set_beta_auto | ( | const Scalar & | shift_in | ) |
|
protected |
Sets to 0x0
all the pointers on derived quantities.
Definition at line 234 of file single_hor.C.
|
inline |
void Lorene::Single_hor::set_n_auto | ( | const Scalar & | nn_in | ) |
|
inline |
void Lorene::Single_hor::set_psi_auto | ( | const Scalar & | psi_in | ) |
Sets the conformal factor relating the physical metric to the conformal one: .
is defined by
Sets the value at the current time step (jtime
) and delete all quantities which depend on .
Definition at line 365 of file single_hor.C.
References psi_auto.
|
inline |
Sets the radius of the horizon to rad
.
Definition at line 1051 of file isol_hor.h.
References radius.
double Lorene::Single_hor::viriel_seul | ( | ) | const |
Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively and N .
WARNING this should only be used for an isolated black hole.
Definition at line 60 of file binhor_viriel.C.
References Lorene::Scalar::asymptot(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), mp, n_auto, and psi_auto.
|
protected |
Components of the conformal representation of the traceless part of the extrinsic curvature:
Definition at line 971 of file isol_hor.h.
|
protected |
Components of the conformal representation of the traceless part of the extrinsic curvature:
Definition at line 959 of file isol_hor.h.
|
protected |
Components of the conformal representation of the traceless part of the extrinsic curvature:
Definition at line 965 of file isol_hor.h.
|
protected |
Shift function .
Definition at line 950 of file isol_hor.h.
|
protected |
Shift function .
Definition at line 944 of file isol_hor.h.
|
protected |
Shift function .
Definition at line 947 of file isol_hor.h.
|
protected |
Function used to construct from the total .
Only used for a binary system.
Mainly this Scalar
is 1 around the hole and 0 around the companion and the sum of decouple
for the hole and his companion is 1 everywhere.
Definition at line 1002 of file isol_hor.h.
|
protected |
Covariant derivative of the lapse with respect to the flat metric .
Definition at line 937 of file isol_hor.h.
|
protected |
Covariant derivative of the conformal factor .
Definition at line 941 of file isol_hor.h.
|
protected |
3 metric flat
Definition at line 980 of file isol_hor.h.
|
protected |
Time derivative of the 3-metric tilde.
Definition at line 986 of file isol_hor.h.
|
protected |
Deviation metric.
Definition at line 983 of file isol_hor.h.
|
protected |
Affine mapping.
Definition at line 900 of file isol_hor.h.
Dirichlet boundary condition for c N f partial_r N a Lorene::Single_hor::N |
Definition at line 1282 of file isol_hor.h.
Neumann boundary condition on nn f partial_r N a Lorene::Single_hor::N |
Definition at line 1286 of file isol_hor.h.
|
protected |
Lapse function .
Definition at line 915 of file isol_hor.h.
|
protected |
Lapse function .
Definition at line 918 of file isol_hor.h.
|
protected |
Lapse function .
Definition at line 921 of file isol_hor.h.
|
protected |
Number of zones.
Definition at line 903 of file isol_hor.h.
|
protected |
Angular velocity in LORENE's units.
Definition at line 909 of file isol_hor.h.
|
mutableprotected |
Spatial metric .
Definition at line 953 of file isol_hor.h.
|
mutableprotected |
Components of the extrinsic curvature:
Definition at line 974 of file isol_hor.h.
|
mutableprotected |
Conformal factor .
Definition at line 933 of file isol_hor.h.
|
protected |
Conformal factor .
Definition at line 930 of file isol_hor.h.
|
protected |
Conformal factor .
Definition at line 924 of file isol_hor.h.
|
protected |
Conformal factor .
Definition at line 927 of file isol_hor.h.
|
protected |
Radius of the horizon in LORENE's units.
Definition at line 906 of file isol_hor.h.
|
protected |
Intensity of the correction on the shift vector.
Definition at line 912 of file isol_hor.h.
|
protected |
3 metric tilde
Definition at line 977 of file isol_hor.h.
|
protected |
Trace of the extrinsic curvature.
Definition at line 989 of file isol_hor.h.
|
protected |
Time derivative of the trace of the extrinsic curvature.
Definition at line 992 of file isol_hor.h.