Single_hor Class Reference
[Stars and black holes]

Binary black holes system. More...

#include <isol_hor.h>

List of all members.

Public Member Functions

 Single_hor (Map_af &mpi)
 Standard constructor.
 Single_hor (const Single_hor &)
 Copy constructor.
 Single_hor (Map_af &mp, FILE *fich)
 Constructor from a binary file.
virtual ~Single_hor ()
 Destructor.
void operator= (const Single_hor &)
 Assignment to another Single_hor.
const Map_afget_mp () const
 Returns the mapping (readonly).
Map_afset_mp ()
 Read/write of the mapping.
double get_radius () const
 Returns the radius of the horizon.
void set_radius (double rad)
 Sets the radius of the horizon to rad .
double get_omega () const
 Returns the angular velocity.
void set_omega (double ome)
 Sets the angular velocity to ome .
const Scalarget_n_auto () const
 Lapse function $ N_{auto} $.
const Scalarget_n_comp () const
 Lapse function $ N_{comp} $.
const Scalarget_nn () const
 Lapse function $ N$.
const Scalarget_psi_auto () const
 Conformal factor $ \Psi_{auto} $.
const Scalarget_psi_comp () const
 Conformal factor $ \Psi_{comp} $.
const Scalarget_psi () const
 Conformal factor $ \Psi $.
const Scalarget_psi4 () const
 Conformal factor $ \Psi^4 $.
const Vectorget_dn () const
 Covariant derivative of the lapse function $ \overline\nabla_i N $.
const Vectorget_dpsi () const
 Covariant derivative with respect to the flat metric of the conformal factor $ \overline\nabla_i \Psi $.
const Vectorget_beta_auto () const
 Shift function $ \beta^i_{auto} $.
const Vectorget_beta_comp () const
 Shift function $ \beta^i_{comp} $.
const Vectorget_beta () const
 Shift function $ \beta^i $.
const Sym_tensorget_aa_auto () const
 Conformal representation $ A^{ij}_{auto} $ of the traceless part of the extrinsic curvature:.
const Sym_tensorget_aa_comp () const
 Conformal representation $ A^{ij}_{comp} $ of the traceless part of the extrinsic curvature:.
const Sym_tensorget_aa () const
 Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature:.
const Metricget_tgam () const
 Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $.
const Metricget_gam () const
 metric $ \gamma_{ij} $
const Sym_tensorget_k_dd () const
 k_dd
const Scalar get_decouple () const
 Returns the function used to construct tkij_auto from tkij_tot .
void n_comp_import (const Single_hor &comp)
 Imports the part of N due to the companion hole comp .
void psi_comp_import (const Single_hor &comp)
 Imports the part of $\Psi$ due to the companion hole comp .
void beta_comp_import (const Single_hor &comp)
 Imports the part of $\beta^i$ due to the companion hole comp.
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 $\Psi$ and N .
void init_bhole ()
 Sets the values of the fields to :.
void init_met_trK ()
 Sets the 3-metric tilde to the flat metric and gamt_point, trK and trK_point to zero.
void init_bhole_seul ()
 Initiates for a single black hole.
void set_psi_auto (const Scalar &psi_in)
 Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
void set_n_auto (const Scalar &nn_in)
 Sets the lapse.
void set_beta_auto (const Scalar &shift_in)
 Sets the shift.
void set_aa_auto (const Scalar &aa_auto_in)
 Sets aa_auto.
void set_aa_comp (const Scalar &aa_comp_in)
 Sets aa_comp.
void set_aa (const Scalar &aa_in)
 Sets aa.
const Scalar b_tilde () const
 Radial component of the shift with respect to the conformal metric.
const Scalar darea_hor () const
 Element of area of the horizon.
double area_hor () const
 Area of the horizon.
double radius_hor () const
 Radius of the horizon.
double ang_mom_hor () const
 Angular momentum (modulo).
double mass_hor () const
 Mass computed at the horizon.
double kappa_hor () const
 Surface gravity.
double omega_hor () const
 Orbital velocity.
double ang_mom_adm () const
 ADM angular Momentum.
Scalar expansion () const
 Expansion of the outgoing null normal ($ \bf n + \bf s $).
const Valeur boundary_psi_app_hor () const
 Neumann boundary condition for.
const Valeur boundary_nn_Dir (double aa) const
 Dirichlet boundary condition for N $ \partial_r N + a N = 0 $.
const Valeur boundary_nn_Neu (double aa) const
 Neumann boundary condition on nn $ \partial_r N + a N = 0 $.
const Valeur boundary_beta_x (double om_orb, double om_loc) const
 Component x of boundary value of $ \beta $.
const Valeur boundary_beta_y (double om_orb, double om_loc) const
 Component y of boundary value of $ \beta $.
const Valeur boundary_beta_z () const
 Component z of boundary value of $ \beta $.
double regularisation (const Vector &shift_auto, const Vector &shift_comp, double ang_vel)
 Corrects shift_auto in such a way that the total $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$.
double regularise_one ()
 Corrects the shift in the innermost shell, so that it remains $ {\mathcal{C}}^2$ and that $A^{ij}$ equals zero on the horizon.
virtual void sauve (FILE *fich) const
 Total or partial saves in a binary file.

Protected Member Functions

void del_deriv () const
 Deletes all the derived quantities.
void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities.

Protected Attributes

Map_afmp
 Affine mapping.
int nz
 Number of zones.
double radius
 Radius of the horizon in LORENE's units.
double omega
 Angular velocity in LORENE's units.
double regul
 Intensity of the correction on the shift vector.
Scalar n_auto
 Lapse function $ N_{auto} $.
Scalar n_comp
 Lapse function $ N_{comp} $.
Scalar nn
 Lapse function $ N $.
Scalar psi_auto
 Conformal factor $ \Psi_{auto} $.
Scalar psi_comp
 Conformal factor $ \Psi_{comp} $.
Scalar psi
 Conformal factor $ \Psi $.
Scalarp_psi4
 Conformal factor $ \Psi^4 $.
Vector dn
 Covariant derivative of the lapse with respect to the flat metric $ \overline\nabla_i N $.
Vector dpsi
 Covariant derivative of the conformal factor $ \overline\nabla_i \Psi $.
Vector beta_auto
 Shift function $ \beta^i_{auto} $.
Vector beta_comp
 Shift function $ \beta^i_{comp} $.
Vector beta
 Shift function $ \beta^i $.
Metricp_gam
 Spatial metric $ \gamma $.
Sym_tensor aa_auto
 Components $ A^{ij}_{auto} $ of the conformal representation of the traceless part of the extrinsic curvature:.
Sym_tensor aa_comp
 Components $ A^{ij}_{comp} $ of the conformal representation of the traceless part of the extrinsic curvature:.
Sym_tensor aa
 Components $ A^{ij} $ of the conformal representation of the traceless part of the extrinsic curvature:.
Sym_tensorp_k_dd
 Components $ K^{ij} $ of the extrinsic curvature:.
Metric tgam
 3 metric tilde
Metric_flat ff
 3 metric flat
Sym_tensor hh
 Deviation metric.
Sym_tensor gamt_point
 Time derivative of the 3-metric tilde.
Scalar trK
 Trace of the extrinsic curvature.
Scalar trK_point
 Time derivative of the trace of the extrinsic curvature.
Scalar decouple
 Function used to construct $ A^{ij}_{auto} $ from the total $A^{ij}$.

Friends

class Bin_hor

Detailed Description

Binary black holes system.

()

This class is intended for dealing with binary black holes configurations in the conformaly flat approximation.

Definition at line 892 of file isol_hor.h.


Constructor & Destructor Documentation

Single_hor::Single_hor ( Map_af mpi  ) 

Standard constructor.

Parameters:
mpi affine mapping

Definition at line 66 of file single_hor.C.

References hh, set_der_0x0(), and Tensor::set_etat_zero().

Single_hor::Single_hor ( const Single_hor singlehor_in  ) 

Copy constructor.

Definition at line 92 of file single_hor.C.

References set_der_0x0().

Single_hor::Single_hor ( Map_af mp,
FILE *  fich 
)

Constructor from a binary file.

Parameters:
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 126 of file single_hor.C.

References Metric_flat::con(), ff, fread_be(), gamt_point, Map::get_bvect_spher(), Map::get_mg(), hh, mp, omega, set_der_0x0(), tgam, trK, and trK_point.

Single_hor::~Single_hor (  )  [virtual]

Destructor.

Definition at line 172 of file single_hor.C.

References del_deriv().


Member Function Documentation

double Single_hor::ang_mom_adm (  )  const

ADM angular Momentum.

Definition at line 170 of file single_param.C.

References Metric::cov(), get_gam(), get_k_dd(), Map_af::integrale_surface_infini(), mp, Scalar::mult_rsint(), and trK.

double Single_hor::ang_mom_hor (  )  const
double Single_hor::area_hor (  )  const

Area of the horizon.

Definition at line 84 of file single_param.C.

References darea_hor(), Map_af::integrale_surface(), mp, Scalar::raccord(), and radius.

const Scalar Single_hor::b_tilde (  )  const

Radial component of the shift with respect to the conformal metric.

Definition at line 64 of file single_param.C.

References beta, contract(), Tensor::down(), Metric::radial_vect(), and tgam.

void Single_hor::beta_comp_import ( const Single_hor comp  ) 

Imports the part of $\beta^i$ due to the companion hole comp.

The total $\beta^i$ is then calculated.

Definition at line 461 of file single_hor.C.

References beta, beta_auto, beta_comp, Map::get_bvect_cart(), Map::get_bvect_spher(), and mp.

const Valeur Single_hor::boundary_beta_x ( double  om_orb,
double  om_loc 
) const
const Valeur Single_hor::boundary_beta_y ( double  om_orb,
double  om_loc 
) const
const Valeur Single_hor::boundary_beta_z (  )  const
const Valeur Single_hor::boundary_nn_Dir ( double  aa  )  const

Dirichlet boundary condition for N $ \partial_r N + a N = 0 $.

Definition at line 131 of file single_bound.C.

References Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, nn, and Scalar::val_grid_point().

const Valeur Single_hor::boundary_nn_Neu ( double  aa  )  const

Neumann boundary condition on nn $ \partial_r N + a N = 0 $.

Definition at line 102 of file single_bound.C.

References dn, Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, nn, Metric::radial_vect(), tgam, and Scalar::val_grid_point().

const Valeur Single_hor::boundary_psi_app_hor (  )  const
const Scalar Single_hor::darea_hor (  )  const

Element of area of the horizon.

Definition at line 73 of file single_param.C.

References get_gam(), sqrt(), and Scalar::std_spectral_base().

void Single_hor::del_deriv (  )  const [protected]

Deletes all the derived quantities.

Definition at line 217 of file single_hor.C.

References p_gam, p_k_dd, p_psi4, and set_der_0x0().

Scalar Single_hor::expansion (  )  const

Expansion of the outgoing null normal ($ \bf n + \bf s $).

Definition at line 186 of file single_param.C.

References contract(), get_gam(), get_k_dd(), and trK.

const Sym_tensor & Single_hor::get_aa (  )  const

Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature:.

Definition at line 330 of file single_hor.C.

References aa.

const Sym_tensor & Single_hor::get_aa_auto (  )  const

Conformal representation $ A^{ij}_{auto} $ of the traceless part of the extrinsic curvature:.

Definition at line 320 of file single_hor.C.

References aa_auto.

const Sym_tensor & Single_hor::get_aa_comp (  )  const

Conformal representation $ A^{ij}_{comp} $ of the traceless part of the extrinsic curvature:.

Definition at line 325 of file single_hor.C.

References aa_comp.

const Vector & Single_hor::get_beta (  )  const

Shift function $ \beta^i $.

Definition at line 315 of file single_hor.C.

References beta.

const Vector & Single_hor::get_beta_auto (  )  const

Shift function $ \beta^i_{auto} $.

Definition at line 305 of file single_hor.C.

References beta_auto.

const Vector & Single_hor::get_beta_comp (  )  const

Shift function $ \beta^i_{comp} $.

Definition at line 310 of file single_hor.C.

References beta_comp.

const Scalar Single_hor::get_decouple (  )  const [inline]

Returns the function used to construct tkij_auto from tkij_tot .

Definition at line 1144 of file isol_hor.h.

References decouple.

const Vector & Single_hor::get_dn (  )  const

Covariant derivative of the lapse function $ \overline\nabla_i N $.

Definition at line 295 of file single_hor.C.

References dn.

const Vector & Single_hor::get_dpsi (  )  const

Covariant derivative with respect to the flat metric of the conformal factor $ \overline\nabla_i \Psi $.

Definition at line 300 of file single_hor.C.

References dpsi.

const Metric & Single_hor::get_gam (  )  const

metric $ \gamma_{ij} $

Definition at line 335 of file single_hor.C.

References Metric::cov(), get_psi4(), p_gam, and tgam.

const Sym_tensor & Single_hor::get_k_dd (  )  const

k_dd

Definition at line 344 of file single_hor.C.

References aa, get_gam(), get_psi4(), p_k_dd, Tensor::std_spectral_base(), trK, and Tensor::up_down().

const Map_af& Single_hor::get_mp (  )  const [inline]

Returns the mapping (readonly).

Definition at line 1036 of file isol_hor.h.

References mp.

const Scalar & Single_hor::get_n_auto (  )  const

Lapse function $ N_{auto} $.

Definition at line 257 of file single_hor.C.

References n_auto.

const Scalar & Single_hor::get_n_comp (  )  const

Lapse function $ N_{comp} $.

Definition at line 262 of file single_hor.C.

References n_comp.

const Scalar & Single_hor::get_nn (  )  const

Lapse function $ N$.

Definition at line 266 of file single_hor.C.

References nn.

double Single_hor::get_omega (  )  const [inline]

Returns the angular velocity.

Definition at line 1054 of file isol_hor.h.

References omega.

const Scalar & Single_hor::get_psi (  )  const

Conformal factor $ \Psi $.

Definition at line 280 of file single_hor.C.

References psi.

const Scalar & Single_hor::get_psi4 (  )  const

Conformal factor $ \Psi^4 $.

Definition at line 284 of file single_hor.C.

References p_psi4, pow(), psi, and Scalar::std_spectral_base().

const Scalar & Single_hor::get_psi_auto (  )  const

Conformal factor $ \Psi_{auto} $.

Definition at line 271 of file single_hor.C.

References psi_auto.

const Scalar & Single_hor::get_psi_comp (  )  const

Conformal factor $ \Psi_{comp} $.

Definition at line 276 of file single_hor.C.

References psi_comp.

double Single_hor::get_radius (  )  const [inline]

Returns the radius of the horizon.

Definition at line 1044 of file isol_hor.h.

References radius.

const Metric& Single_hor::get_tgam (  )  const [inline]

Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $.

Definition at line 1131 of file isol_hor.h.

References tgam.

void Single_hor::init_bhole (  ) 

Sets the values of the fields to :.

  • n_auto $= \frac{1}{2}-2\frac{a}{r}$
  • n_comp $= \frac{1}{2}$
  • psi_auto $= \frac{1}{2}+\frac{a}{r}$
  • psi_comp $= \frac{1}{2}$

a being the radius of the hole, the other fields being set to zero.

Definition at line 480 of file single_hor.C.

References Scalar::annule(), beta, beta_auto, beta_comp, Scalar::derive_cov(), dn, dpsi, ff, Map::get_bvect_spher(), mp, n_auto, n_comp, nn, psi, psi_auto, psi_comp, Map::r, Scalar::raccord(), radius, Coord::set(), Vector::set(), Scalar::set_dzpuis(), and Scalar::std_spectral_base().

void 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 544 of file single_hor.C.

References Scalar::annule(), beta, beta_auto, beta_comp, Scalar::derive_cov(), dn, dpsi, ff, Map::get_bvect_spher(), Map::get_mg(), mp, n_auto, n_comp, nn, psi, psi_auto, psi_comp, Map::r, Scalar::raccord(), radius, Scalar::set_dzpuis(), Tensor::set_etat_zero(), Scalar::set_etat_zero(), Scalar::set_outer_boundary(), and Scalar::std_spectral_base().

void 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 532 of file single_hor.C.

References Map::flat_met_spher(), gamt_point, mp, Scalar::set_etat_zero(), Tensor::set_etat_zero(), tgam, trK, and trK_point.

double Single_hor::kappa_hor (  )  const

Surface gravity.

Definition at line 142 of file single_param.C.

References ang_mom_hor(), pow(), radius_hor(), and sqrt().

double Single_hor::mass_hor (  )  const

Mass computed at the horizon.

Definition at line 131 of file single_param.C.

References ang_mom_hor(), pow(), radius_hor(), and sqrt().

void 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 $\nabla_i N$.

Definition at line 386 of file single_hor.C.

References Tensor::derive_cov(), Scalar::derive_cov(), dn, ff, Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::import(), mp, n_auto, n_comp, nn, Scalar::set_etat_qcq(), and Scalar::std_spectral_base().

double Single_hor::omega_hor (  )  const

Orbital velocity.

Definition at line 156 of file single_param.C.

References ang_mom_hor(), pow(), radius_hor(), and sqrt().

void Single_hor::operator= ( const Single_hor singlehor_in  ) 

Assignment to another Single_hor.

Definition at line 182 of file single_hor.C.

References aa, aa_auto, aa_comp, beta, beta_auto, beta_comp, decouple, dn, dpsi, ff, gamt_point, hh, mp, n_auto, n_comp, nn, nz, omega, psi, psi_auto, psi_comp, radius, regul, tgam, trK, and trK_point.

void Single_hor::psi_comp_import ( const Single_hor comp  ) 

Imports the part of $\Psi$ due to the companion hole comp .

The total $\Psi$ is then calculated.

It also imports the covariant derivative of $\Psi$ and construct the total $\nabla_i \Psi$.

Definition at line 420 of file single_hor.C.

References Tensor::derive_cov(), Scalar::derive_cov(), dpsi, ff, Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::import(), mp, psi, psi_auto, psi_comp, Scalar::set_etat_qcq(), and Scalar::std_spectral_base().

double Single_hor::radius_hor (  )  const

Radius of the horizon.

Definition at line 93 of file single_param.C.

References area_hor(), and pow().

double Single_hor::regularisation ( const Vector shift_auto,
const Vector shift_comp,
double  ang_vel 
)

Corrects shift_auto in such a way that the total $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$.

WARNING : this should only be used for a black hole in a binary system Bin_hor.

Parameters:
comp [input]: the part of $\beta^i$ generated by the companion hole.

Definition at line 55 of file single_regul.C.

References Cmp::annule(), Tensor::annule_domain(), beta_auto, Valeur::c, Vector::change_triad(), Valeur::coef_i(), diffrelmax(), Map::get_bvect_cart(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Map::get_rot_phi(), Scalar::get_spectral_va(), Tensor::get_triad(), norme(), nz, pow(), Vector::set(), Valeur::set_etat_c_qcq(), Scalar::set_spectral_va(), and Vector::std_spectral_base().

double Single_hor::regularise_one (  ) 

Corrects the shift in the innermost shell, so that it remains $ {\mathcal{C}}^2$ and that $A^{ij}$ 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 149 of file single_regul.C.

References Cmp::annule(), Tensor::annule_domain(), beta, Vector::change_triad(), Valeur::coef_i(), diffrelmax(), Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), mp, norme(), nz, omega, pow(), Map::r, Vector::set(), Valeur::set_etat_c_qcq(), Scalar::set_spectral_va(), Vector::std_spectral_base(), Map_af::val_r(), Map::x, and Map::y.

void Single_hor::sauve ( FILE *  fich  )  const [virtual]

Total or partial saves in a binary file.

Parameters:
fich binary file
partial_save indicates whether the whole object must be saved.

Definition at line 240 of file single_hor.C.

References beta_auto, Metric::con(), fwrite_be(), gamt_point, n_auto, omega, psi_auto, Tensor_sym::sauve(), Tensor::sauve(), Scalar::sauve(), tgam, trK, and trK_point.

void Single_hor::set_aa ( const Scalar aa_in  ) 

Sets aa.

Definition at line 378 of file single_hor.C.

References aa.

void Single_hor::set_aa_auto ( const Scalar aa_auto_in  ) 

Sets aa_auto.

Definition at line 370 of file single_hor.C.

References aa_auto.

void Single_hor::set_aa_comp ( const Scalar aa_comp_in  ) 

Sets aa_comp.

Definition at line 374 of file single_hor.C.

References aa_comp.

void Single_hor::set_beta_auto ( const Scalar shift_in  ) 

Sets the shift.

Definition at line 366 of file single_hor.C.

References beta_auto.

void Single_hor::set_der_0x0 (  )  const [protected]

Sets to 0x0 all the pointers on derived quantities.

Definition at line 227 of file single_hor.C.

References p_gam, p_k_dd, and p_psi4.

Map_af& Single_hor::set_mp (  )  [inline]

Read/write of the mapping.

Definition at line 1039 of file isol_hor.h.

References mp.

void Single_hor::set_n_auto ( const Scalar nn_in  ) 

Sets the lapse.

Definition at line 362 of file single_hor.C.

References n_auto.

void Single_hor::set_omega ( double  ome  )  [inline]

Sets the angular velocity to ome .

Definition at line 1058 of file isol_hor.h.

References omega.

void Single_hor::set_psi_auto ( const Scalar psi_in  ) 

Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.

$ \Psi $ is defined by

\[ \Psi := \left( \frac{\det\gamma_{ij}}{\det f_{ij}} \right) ^{1/12} \]

Sets the value at the current time step (jtime ) and delete all quantities which depend on $ \Psi $.

Definition at line 358 of file single_hor.C.

References psi_auto.

void Single_hor::set_radius ( double  rad  )  [inline]

Sets the radius of the horizon to rad .

Definition at line 1049 of file isol_hor.h.

References radius.

double 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 $\Psi$ and N .

WARNING this should only be used for an isolated black hole.

Definition at line 53 of file binhor_viriel.C.

References Scalar::asymptot(), Map::get_mg(), Mg3d::get_nzone(), mp, n_auto, and psi_auto.


Member Data Documentation

Components $ A^{ij} $ of the conformal representation of the traceless part of the extrinsic curvature:.

Definition at line 969 of file isol_hor.h.

Components $ A^{ij}_{auto} $ of the conformal representation of the traceless part of the extrinsic curvature:.

Definition at line 957 of file isol_hor.h.

Components $ A^{ij}_{comp} $ of the conformal representation of the traceless part of the extrinsic curvature:.

Definition at line 963 of file isol_hor.h.

Vector Single_hor::beta [protected]

Shift function $ \beta^i $.

Definition at line 948 of file isol_hor.h.

Shift function $ \beta^i_{auto} $.

Definition at line 942 of file isol_hor.h.

Shift function $ \beta^i_{comp} $.

Definition at line 945 of file isol_hor.h.

Function used to construct $ A^{ij}_{auto} $ from the total $A^{ij}$.

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 1000 of file isol_hor.h.

Vector Single_hor::dn [protected]

Covariant derivative of the lapse with respect to the flat metric $ \overline\nabla_i N $.

Definition at line 935 of file isol_hor.h.

Vector Single_hor::dpsi [protected]

Covariant derivative of the conformal factor $ \overline\nabla_i \Psi $.

Definition at line 939 of file isol_hor.h.

3 metric flat

Definition at line 978 of file isol_hor.h.

Time derivative of the 3-metric tilde.

Definition at line 984 of file isol_hor.h.

Deviation metric.

Definition at line 981 of file isol_hor.h.

Map_af& Single_hor::mp [protected]

Affine mapping.

Definition at line 898 of file isol_hor.h.

Lapse function $ N_{auto} $.

Definition at line 913 of file isol_hor.h.

Lapse function $ N_{comp} $.

Definition at line 916 of file isol_hor.h.

Scalar Single_hor::nn [protected]

Lapse function $ N $.

Definition at line 919 of file isol_hor.h.

int Single_hor::nz [protected]

Number of zones.

Definition at line 901 of file isol_hor.h.

double Single_hor::omega [protected]

Angular velocity in LORENE's units.

Definition at line 907 of file isol_hor.h.

Metric* Single_hor::p_gam [mutable, protected]

Spatial metric $ \gamma $.

Definition at line 951 of file isol_hor.h.

Sym_tensor* Single_hor::p_k_dd [mutable, protected]

Components $ K^{ij} $ of the extrinsic curvature:.

Definition at line 972 of file isol_hor.h.

Scalar* Single_hor::p_psi4 [mutable, protected]

Conformal factor $ \Psi^4 $.

Definition at line 931 of file isol_hor.h.

Scalar Single_hor::psi [protected]

Conformal factor $ \Psi $.

Definition at line 928 of file isol_hor.h.

Conformal factor $ \Psi_{auto} $.

Definition at line 922 of file isol_hor.h.

Conformal factor $ \Psi_{comp} $.

Definition at line 925 of file isol_hor.h.

double Single_hor::radius [protected]

Radius of the horizon in LORENE's units.

Definition at line 904 of file isol_hor.h.

double Single_hor::regul [protected]

Intensity of the correction on the shift vector.

Definition at line 910 of file isol_hor.h.

Metric Single_hor::tgam [protected]

3 metric tilde

Definition at line 975 of file isol_hor.h.

Scalar Single_hor::trK [protected]

Trace of the extrinsic curvature.

Definition at line 987 of file isol_hor.h.

Time derivative of the trace of the extrinsic curvature.

Definition at line 990 of file isol_hor.h.


The documentation for this class was generated from the following files:

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1