LORENE
Lorene::Single_hor Class Reference

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_afget_mp () const
 Returns the mapping (readonly). More...
 
Map_afset_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 Scalarget_n_auto () const
 Lapse function $ N_{auto} $. More...
 
const Scalarget_n_comp () const
 Lapse function $ N_{comp} $. More...
 
const Scalarget_nn () const
 Lapse function $ N$. More...
 
const Scalarget_psi_auto () const
 Conformal factor $ \Psi_{auto} $. More...
 
const Scalarget_psi_comp () const
 Conformal factor $ \Psi_{comp} $. More...
 
const Scalarget_psi () const
 Conformal factor $ \Psi $. More...
 
const Scalarget_psi4 () const
 Conformal factor $ \Psi^4 $. More...
 
const Vectorget_dn () const
 Covariant derivative of the lapse function $ \overline\nabla_i N $. More...
 
const Vectorget_dpsi () const
 Covariant derivative with respect to the flat metric of the conformal factor $ \overline\nabla_i \Psi $. More...
 
const Vectorget_beta_auto () const
 Shift function $ \beta^i_{auto} $. More...
 
const Vectorget_beta_comp () const
 Shift function $ \beta^i_{comp} $. More...
 
const Vectorget_beta () const
 Shift function $ \beta^i $. More...
 
const Sym_tensorget_aa_auto () const
 Conformal representation $ A^{ij}_{auto} $ of the traceless part of the extrinsic curvature: More...
 
const Sym_tensorget_aa_comp () const
 Conformal representation $ A^{ij}_{comp} $ of the traceless part of the extrinsic curvature: More...
 
const Sym_tensorget_aa () const
 Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: More...
 
const Metricget_tgam () const
 Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $. More...
 
const Metricget_gam () const
 metric $ \gamma_{ij} $ More...
 
const Sym_tensorget_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 $\Psi$ due to the companion hole comp . More...
 
void beta_comp_import (const Single_hor &comp)
 Imports the part of $\beta^i$ 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 $\Psi$ 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 $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $. 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 ( $ \bf n + \bf s $) 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 $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$. More...
 
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. 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_afmp
 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 $ N_{auto} $. More...
 
Scalar n_comp
 Lapse function $ N_{comp} $. More...
 
Scalar nn
 Lapse function $ N $. More...
 
Scalar psi_auto
 Conformal factor $ \Psi_{auto} $. More...
 
Scalar psi_comp
 Conformal factor $ \Psi_{comp} $. More...
 
Scalar psi
 Conformal factor $ \Psi $. More...
 
Scalarp_psi4
 Conformal factor $ \Psi^4 $. More...
 
Vector dn
 Covariant derivative of the lapse with respect to the flat metric $ \overline\nabla_i N $. More...
 
Vector dpsi
 Covariant derivative of the conformal factor $ \overline\nabla_i \Psi $. More...
 
Vector beta_auto
 Shift function $ \beta^i_{auto} $. More...
 
Vector beta_comp
 Shift function $ \beta^i_{comp} $. More...
 
Vector beta
 Shift function $ \beta^i $. More...
 
Metricp_gam
 Spatial metric $ \gamma $. More...
 
Sym_tensor aa_auto
 Components $ A^{ij}_{auto} $ of the conformal representation of the traceless part of the extrinsic curvature: More...
 
Sym_tensor aa_comp
 Components $ A^{ij}_{comp} $ of the conformal representation of the traceless part of the extrinsic curvature: More...
 
Sym_tensor aa
 Components $ A^{ij} $ of the conformal representation of the traceless part of the extrinsic curvature: More...
 
Sym_tensorp_k_dd
 Components $ K^{ij} $ 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 $ A^{ij}_{auto} $ from the total $A^{ij}$. More...
 

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

Constructor & Destructor Documentation

◆ Single_hor() [1/3]

Lorene::Single_hor::Single_hor ( Map_af mpi)

Standard constructor.

Parameters
mpiaffine mapping

Definition at line 73 of file single_hor.C.

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

◆ Single_hor() [2/3]

Lorene::Single_hor::Single_hor ( const Single_hor singlehor_in)

Copy constructor.

Definition at line 99 of file single_hor.C.

References set_der_0x0().

◆ Single_hor() [3/3]

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

Constructor from a binary file.

Parameters
mpiaffine mapping
fichfile containing the saved isol_hor
partial_readindicates 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.

◆ ~Single_hor()

Lorene::Single_hor::~Single_hor ( )
virtual

Destructor.

Definition at line 179 of file single_hor.C.

References del_deriv().

Member Function Documentation

◆ ang_mom_adm()

double Lorene::Single_hor::ang_mom_adm ( ) const

◆ ang_mom_hor()

◆ area_hor()

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.

◆ b_tilde()

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.

◆ beta_comp_import()

void Lorene::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 468 of file single_hor.C.

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

◆ boundary_psi_app_hor()

◆ darea_hor()

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().

◆ del_deriv()

void Lorene::Single_hor::del_deriv ( ) const
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().

◆ expansion()

Scalar Lorene::Single_hor::expansion ( ) const

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

Definition at line 193 of file single_param.C.

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

◆ get_aa()

const Sym_tensor & Lorene::Single_hor::get_aa ( ) const

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

Definition at line 337 of file single_hor.C.

References aa.

◆ get_aa_auto()

const Sym_tensor & Lorene::Single_hor::get_aa_auto ( ) const

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

Definition at line 327 of file single_hor.C.

References aa_auto.

◆ get_aa_comp()

const Sym_tensor & Lorene::Single_hor::get_aa_comp ( ) const

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

Definition at line 332 of file single_hor.C.

References aa_comp.

◆ get_beta()

const Vector & Lorene::Single_hor::get_beta ( ) const

Shift function $ \beta^i $.

Definition at line 322 of file single_hor.C.

References beta.

◆ get_beta_auto()

const Vector & Lorene::Single_hor::get_beta_auto ( ) const

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

Definition at line 312 of file single_hor.C.

References beta_auto.

◆ get_beta_comp()

const Vector & Lorene::Single_hor::get_beta_comp ( ) const

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

Definition at line 317 of file single_hor.C.

References beta_comp.

◆ get_decouple()

const Scalar Lorene::Single_hor::get_decouple ( ) const
inline

Returns the function used to construct tkij_auto from tkij_tot .

Definition at line 1146 of file isol_hor.h.

References decouple.

◆ get_dn()

const Vector & Lorene::Single_hor::get_dn ( ) const

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

Definition at line 302 of file single_hor.C.

References dn.

◆ get_dpsi()

const Vector & Lorene::Single_hor::get_dpsi ( ) const

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

Definition at line 307 of file single_hor.C.

References dpsi.

◆ get_gam()

const Metric & Lorene::Single_hor::get_gam ( ) const

metric $ \gamma_{ij} $

Definition at line 342 of file single_hor.C.

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

◆ get_k_dd()

const Sym_tensor & Lorene::Single_hor::get_k_dd ( ) const

◆ get_mp()

const Map_af& Lorene::Single_hor::get_mp ( ) const
inline

Returns the mapping (readonly).

Definition at line 1038 of file isol_hor.h.

References mp.

◆ get_n_auto()

const Scalar & Lorene::Single_hor::get_n_auto ( ) const

Lapse function $ N_{auto} $.

Definition at line 264 of file single_hor.C.

References n_auto.

◆ get_n_comp()

const Scalar & Lorene::Single_hor::get_n_comp ( ) const

Lapse function $ N_{comp} $.

Definition at line 269 of file single_hor.C.

References n_comp.

◆ get_nn()

const Scalar & Lorene::Single_hor::get_nn ( ) const

Lapse function $ N$.

Definition at line 273 of file single_hor.C.

References nn.

◆ get_omega()

double Lorene::Single_hor::get_omega ( ) const
inline

Returns the angular velocity.

Definition at line 1056 of file isol_hor.h.

References omega.

◆ get_psi()

const Scalar & Lorene::Single_hor::get_psi ( ) const

Conformal factor $ \Psi $.

Definition at line 287 of file single_hor.C.

References psi.

◆ get_psi4()

const Scalar & Lorene::Single_hor::get_psi4 ( ) const

Conformal factor $ \Psi^4 $.

Definition at line 291 of file single_hor.C.

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

◆ get_psi_auto()

const Scalar & Lorene::Single_hor::get_psi_auto ( ) const

Conformal factor $ \Psi_{auto} $.

Definition at line 278 of file single_hor.C.

References psi_auto.

◆ get_psi_comp()

const Scalar & Lorene::Single_hor::get_psi_comp ( ) const

Conformal factor $ \Psi_{comp} $.

Definition at line 283 of file single_hor.C.

References psi_comp.

◆ get_radius()

double Lorene::Single_hor::get_radius ( ) const
inline

Returns the radius of the horizon.

Definition at line 1046 of file isol_hor.h.

References radius.

◆ get_tgam()

const Metric& Lorene::Single_hor::get_tgam ( ) const
inline

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

Definition at line 1133 of file isol_hor.h.

References tgam.

◆ init_bhole()

void Lorene::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 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().

◆ init_bhole_seul()

void Lorene::Single_hor::init_bhole_seul ( )

◆ init_met_trK()

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.

◆ kappa_hor()

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().

◆ mass_hor()

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().

◆ n_comp_import()

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

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().

◆ omega_hor()

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().

◆ operator=()

void Lorene::Single_hor::operator= ( const Single_hor singlehor_in)

Assignment to another Single_hor.

Definition at line 189 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.

◆ psi_comp_import()

void Lorene::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 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().

◆ radius_hor()

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().

◆ regularisation()

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

◆ regularise_one()

double Lorene::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 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.

◆ sauve()

void Lorene::Single_hor::sauve ( FILE *  fich) const
virtual

Total or partial saves in a binary file.

Parameters
fichbinary file
partial_saveindicates 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.

◆ set_aa()

void Lorene::Single_hor::set_aa ( const Scalar aa_in)

Sets aa.

Definition at line 385 of file single_hor.C.

References aa.

◆ set_aa_auto()

void Lorene::Single_hor::set_aa_auto ( const Scalar aa_auto_in)

Sets aa_auto.

Definition at line 377 of file single_hor.C.

References aa_auto.

◆ set_aa_comp()

void Lorene::Single_hor::set_aa_comp ( const Scalar aa_comp_in)

Sets aa_comp.

Definition at line 381 of file single_hor.C.

References aa_comp.

◆ set_beta_auto()

void Lorene::Single_hor::set_beta_auto ( const Scalar shift_in)

Sets the shift.

Definition at line 373 of file single_hor.C.

References beta_auto.

◆ set_der_0x0()

void Lorene::Single_hor::set_der_0x0 ( ) const
protected

Sets to 0x0 all the pointers on derived quantities.

Definition at line 234 of file single_hor.C.

References p_gam, p_k_dd, and p_psi4.

◆ set_mp()

Map_af& Lorene::Single_hor::set_mp ( )
inline

Read/write of the mapping.

Definition at line 1041 of file isol_hor.h.

References mp.

◆ set_n_auto()

void Lorene::Single_hor::set_n_auto ( const Scalar nn_in)

Sets the lapse.

Definition at line 369 of file single_hor.C.

References n_auto.

◆ set_omega()

void Lorene::Single_hor::set_omega ( double  ome)
inline

Sets the angular velocity to ome .

Definition at line 1060 of file isol_hor.h.

References omega.

◆ set_psi_auto()

void Lorene::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 365 of file single_hor.C.

References psi_auto.

◆ set_radius()

void Lorene::Single_hor::set_radius ( double  rad)
inline

Sets the radius of the horizon to rad .

Definition at line 1051 of file isol_hor.h.

References radius.

◆ viriel_seul()

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

Member Data Documentation

◆ aa

Sym_tensor Lorene::Single_hor::aa
protected

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

Definition at line 971 of file isol_hor.h.

◆ aa_auto

Sym_tensor Lorene::Single_hor::aa_auto
protected

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

Definition at line 959 of file isol_hor.h.

◆ aa_comp

Sym_tensor Lorene::Single_hor::aa_comp
protected

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

Definition at line 965 of file isol_hor.h.

◆ beta

Vector Lorene::Single_hor::beta
protected

Shift function $ \beta^i $.

Definition at line 950 of file isol_hor.h.

◆ beta_auto

Vector Lorene::Single_hor::beta_auto
protected

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

Definition at line 944 of file isol_hor.h.

◆ beta_comp

Vector Lorene::Single_hor::beta_comp
protected

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

Definition at line 947 of file isol_hor.h.

◆ decouple

Scalar Lorene::Single_hor::decouple
protected

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

◆ dn

Vector Lorene::Single_hor::dn
protected

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

Definition at line 937 of file isol_hor.h.

◆ dpsi

Vector Lorene::Single_hor::dpsi
protected

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

Definition at line 941 of file isol_hor.h.

◆ ff

Metric_flat Lorene::Single_hor::ff
protected

3 metric flat

Definition at line 980 of file isol_hor.h.

◆ gamt_point

Sym_tensor Lorene::Single_hor::gamt_point
protected

Time derivative of the 3-metric tilde.

Definition at line 986 of file isol_hor.h.

◆ hh

Sym_tensor Lorene::Single_hor::hh
protected

Deviation metric.

Definition at line 983 of file isol_hor.h.

◆ mp

Map_af& Lorene::Single_hor::mp
protected

Affine mapping.

Definition at line 900 of file isol_hor.h.

◆ N [1/2]

Dirichlet boundary condition for c N f partial_r N a Lorene::Single_hor::N
Initial value:
= 0 \f$
const Valeur boundary_nn_Dir(double aa) const

Definition at line 1282 of file isol_hor.h.

◆ N [2/2]

Neumann boundary condition on nn f partial_r N a Lorene::Single_hor::N
Initial value:
= 0 \f$
const Valeur boundary_nn_Neu(double aa) const

Definition at line 1286 of file isol_hor.h.

◆ n_auto

Scalar Lorene::Single_hor::n_auto
protected

Lapse function $ N_{auto} $.

Definition at line 915 of file isol_hor.h.

◆ n_comp

Scalar Lorene::Single_hor::n_comp
protected

Lapse function $ N_{comp} $.

Definition at line 918 of file isol_hor.h.

◆ nn

Scalar Lorene::Single_hor::nn
protected

Lapse function $ N $.

Definition at line 921 of file isol_hor.h.

◆ nz

int Lorene::Single_hor::nz
protected

Number of zones.

Definition at line 903 of file isol_hor.h.

◆ omega

double Lorene::Single_hor::omega
protected

Angular velocity in LORENE's units.

Definition at line 909 of file isol_hor.h.

◆ p_gam

Metric* Lorene::Single_hor::p_gam
mutableprotected

Spatial metric $ \gamma $.

Definition at line 953 of file isol_hor.h.

◆ p_k_dd

Sym_tensor* Lorene::Single_hor::p_k_dd
mutableprotected

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

Definition at line 974 of file isol_hor.h.

◆ p_psi4

Scalar* Lorene::Single_hor::p_psi4
mutableprotected

Conformal factor $ \Psi^4 $.

Definition at line 933 of file isol_hor.h.

◆ psi

Scalar Lorene::Single_hor::psi
protected

Conformal factor $ \Psi $.

Definition at line 930 of file isol_hor.h.

◆ psi_auto

Scalar Lorene::Single_hor::psi_auto
protected

Conformal factor $ \Psi_{auto} $.

Definition at line 924 of file isol_hor.h.

◆ psi_comp

Scalar Lorene::Single_hor::psi_comp
protected

Conformal factor $ \Psi_{comp} $.

Definition at line 927 of file isol_hor.h.

◆ radius

double Lorene::Single_hor::radius
protected

Radius of the horizon in LORENE's units.

Definition at line 906 of file isol_hor.h.

◆ regul

double Lorene::Single_hor::regul
protected

Intensity of the correction on the shift vector.

Definition at line 912 of file isol_hor.h.

◆ tgam

Metric Lorene::Single_hor::tgam
protected

3 metric tilde

Definition at line 977 of file isol_hor.h.

◆ trK

Scalar Lorene::Single_hor::trK
protected

Trace of the extrinsic curvature.

Definition at line 989 of file isol_hor.h.

◆ trK_point

Scalar Lorene::Single_hor::trK_point
protected

Time derivative of the trace of the extrinsic curvature.

Definition at line 992 of file isol_hor.h.


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