Lorene::App_hor Class Reference

Class describing an apparent horizon. More...

#include <spheroid.h>

Inheritance diagram for Lorene::App_hor:
Lorene::Spheroid

List of all members.

Public Member Functions

 App_hor (const Mg3d &grid_angu, double radius)
 Standard constructor.
 App_hor (const Scalar &h_in, const Metric &gamij, const Sym_tensor &Kij)
 Constructor of an apparent horizon embedded in a 3-slice (Time_slice ) of 3+1 formalism.
 App_hor (const App_hor &)
 Copy constructor.
 App_hor (FILE *)
 Constructor from a file (see sauve(FILE*) ).
virtual ~App_hor ()
 Destructor.
void operator= (const App_hor &)
 Assignment to another App_hor.
bool check_expansion (double thres=1.e-7) const
const Sym_tensorlie_derive_shear (const Scalar &bb, const Scalar &lapse)
 Lie derivative of shear tensor with respect to the evolution vector field.
const Sym_tensorlie_derive_theta_plus (const Scalar &bb, const Scalar &lapse)
 Lie derivative of the null outgoing expansion rate with respect to the evolution vector field.
const Sym_tensorlie_derive_theta_minus (const Scalar &bb, const Scalar &lapse)
 Lie derivative of the null ingoing expansion rate with respect to the evolution vector field.
const Sym_tensorlie_derive_q_ab (const Scalar &bb, const Scalar &lapse)
 Lie derivative of 2-metric with respect to the evolution vector field.
const Scalarl_non_affinity (const Scalar &bb, const Scalar &lapse)
 non-affinity (or surface gravity) with respect to the outgoing null vector field
void set_ephi (const Scalar &)
 Assigns the conformal Killing vector field to phi.
const Scalarget_hsurf () const
 Returns the field h_surf.
const Metricget_qab () const
 Returns the metric $ q_{ab} $.
const Scalarget_ricci () const
 Returns the 2-ricci scalar $ R = q^{ab}R_{ab} $.
const Sym_tensorget_hh () const
 Returns the symmetric tensor $ H_{ab} $.
const Sym_tensorget_qq () const
 returns the 3-d degenerate 2-metric $ Q_{ab}$
const Tensorget_proj () const
 returns the 3-d projector on 2-surface $ Pi $
const Tensorget_jac2d () const
 returns the 2-d jacobian of coordinate transformation $ J $
const Scalarget_trk () const
 Returns the trace K on the 2-surface.
const Vectorget_ll () const
 Returns the vector $ L_a $.
const Vectorget_ss () const
 Returns the vector $ S_a $.
const Vectorget_ephi () const
 Returns the conformal Killing symmetry vector on the 2-surface.
const Sym_tensorget_jj () const
 Returns the symmetric tensor $ J_{ab} $.
const Scalarget_fff () const
 Returns the normalization scalar F.
const Scalarget_ggg () const
 Returns the normalization scalar G.
bool get_issphere () const
 Returns the flag saying whether or not the horizon is geometrically round.
Scalarset_hsurf ()
 Sets the field h_surf.
Metricset_qab ()
 Sets the modified metric (non degenerated) $ q_{ab} $.
Scalarset_ricci ()
 Sets the 2-Ricci scalar $ R = q^{ab}R_{ab} $.
Sym_tensorset_qq ()
 Sets the degenerated metric $ Q_{ab} $.
Tensorset_proj ()
 Sets the projector $ Pi $.
Sym_tensorset_hh ()
 Sets the symmetric tensor $ H_{ab} $.
Scalarset_trk ()
 Sets the trace K on the 2-surface.
Vectorset_ll ()
 Sets the vector $ L_a $.
Vectorset_ss ()
 Sets the vector $ s_a $.
Sym_tensorset_jj ()
 Sets the symmetric tensor $ J_{ab} $.
Scalarset_fff ()
 Sets the normalization factor F.
Scalarset_ggg ()
 Sets the normalization factor G.
bool set_issphere ()
 Sets the boolean linked to geometrical shape of the horizon.
void update_from_tslice (const Metric &gamij, const Sym_tensor &Kij)
 Updates from the 3-slice data.
const Scalarsqrt_q () const
 Computes the normal vector field to the 2-surface.
double area () const
 Computes the area of the 2-surface.
double angu_mom () const
 Computes the angular momentum with respect to a divergence-free vector field tangent to the 2-surface.
double mass () const
 Computes the mass as defined from the calculus of angular momentum, done with respect to a divergence free tangent vector field $ phi $.
double multipole_mass (const int order) const
 Computes the mass multipole of a given order for the spheroid, assumed to be spherical.
double multipole_angu (const int order) const
 Computes the angular multipole of a given order for the spheroid, assumed to be spherical.
double epsilon_A_minus_one () const
 Computation of the refined Penrose parameter for axisymmetric spacetimes, and its difference wrt one.
double epsilon_P_minus_one () const
 Computation of the classical Penrose parameter, and its difference wrt one.
const Scalartheta_plus () const
 Computes the outgoing null expansion $ \theta_+ $.
const Scalartheta_minus () const
 Computes the ingoing null expansion $ \theta_- $.
const Sym_tensorshear () const
 Computes the shear of the 2-surface $ \sigma_{ab} $.
Tensor derive_cov2dflat (const Tensor &uu) const
 Computes the round covariant derivative on the spheroid.
const Tensordelta () const
 Computes the delta coefficients for covariant derivative.
Tensor derive_cov2d (const Tensor &uu) const
 Computes the total covariant derivative on the spheroid.
virtual void sauve (FILE *) const
 Save in a file.

Protected Member Functions

virtual 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

Scalar h_surf
 The location of the 2-surface as r = h_surf$(\theta, \varphi)$.
Tensor jac2d
 The jacobian of the adaptation of coordinates (contravariant/covariant representation).
Tensor proj
 The 3-d projector on the 2-surface (contravariant-covariant form).
Sym_tensor qq
 The 3-d covariant degenerated 2-metric on the surface.
Vector ss
 The adapted normal vector field to spheroid in the 3-slice.
Vector ephi
 The conformal Killing vector field on the 2-surface (set to by default to the axial vector associated with coordinate $ \varphi $).
Metric qab
Scalar ricci
 Induced metric on the 2-surface $ q_{ab} $.
Sym_tensor hh
 The ricci scalar on the 2-surface.
Scalar trk
 Trace K of the extrinsic curvature of the 3-slice.
Vector ll
 Normal-tangent component of the extrinsic curvature of the 3-slice.
Sym_tensor jj
 Tangent components of the extrinsic curvature of the 3-slice.
Scalar fff
 Normalization function for the outgoing null vector l.
Scalar ggg
 Normalization function for the ingoing null vector k.
Scalar zeta
bool issphere
 Flag to know whether the horizon is geometrically round or distorted.
Scalarp_sqrt_q
 Surface element $\sqrt{\det q_{ab}} $.
double * p_area
 The area of the 2-surface.
double * p_angu_mom
 The angular momentum.
double * p_mass
 Mass defined from angular momentum.
double * p_multipole_mass
 Mass multipole for the spheroid.
double * p_multipole_angu
 Angular momentum multipole for the spheroid.
double * p_epsilon_A_minus_one
double * p_epsilon_P_minus_one
 Refined Penrose parameter, difference wrt one.
Scalarp_theta_plus
 Classical Penrose parameter, difference wrt one.
Scalarp_theta_minus
 Null ingoing expansion.
Sym_tensorp_shear
 The shear tensor.
Tensorp_delta

Friends

ostream & operator<< (ostream &, const Spheroid &)
 Display.

Detailed Description

Class describing an apparent horizon.

It is derived from the class Spheroid, with the property that the outgoing null expansion $ \theta_+ = 0$.

Definition at line 408 of file spheroid.h.


Constructor & Destructor Documentation

Lorene::App_hor::App_hor ( const Mg3d grid_angu,
double  radius 
)

Standard constructor.

Lorene::App_hor::App_hor ( const Scalar h_in,
const Metric gamij,
const Sym_tensor Kij 
)

Constructor of an apparent horizon embedded in a 3-slice (Time_slice ) of 3+1 formalism.

This is done from the Time_slice data.

Parameters:
h_in : the location of the surface r = h_in
gamij : the 3-metric on the 3-slice
Kij : the extrinsic curvature of the 3-slice (covariant representation)
Lorene::App_hor::App_hor ( const App_hor  ) 

Copy constructor.

Lorene::App_hor::App_hor ( FILE *   ) 

Constructor from a file (see sauve(FILE*) ).

virtual Lorene::App_hor::~App_hor (  )  [virtual]

Destructor.


Member Function Documentation

double Lorene::Spheroid::angu_mom (  )  const [inherited]

Computes the angular momentum with respect to a divergence-free vector field tangent to the 2-surface.

This is defined as

\[ {\cal J} = \int h^2 L_i \phi^i \sqrt{\det q_{ab}} \sin \theta {\rm d}\theta {\rm d}\varphi \]

Parameters:
phi : the divergence-free vector field $ \phi $

Definition at line 724 of file spheroid.C.

References Lorene::contract(), Lorene::Map::get_bvect_spher(), Lorene::Spheroid::get_ephi(), Lorene::Tensor::get_mp(), Lorene::Spheroid::h_surf, Lorene::Map_af::integrale_surface(), Lorene::Spheroid::jac2d, Lorene::Spheroid::ll, Lorene::Spheroid::p_angu_mom, and Lorene::Spheroid::sqrt_q().

double Lorene::Spheroid::area (  )  const [inherited]

Computes the area of the 2-surface.

This is defined as

\[ {\cal A} = \int h^2 \sqrt{\det q_{ab}} \sin \theta {\rm d}\theta {\rm d}\varphi \]

Definition at line 711 of file spheroid.C.

References Lorene::Tensor::get_mp(), Lorene::Spheroid::h_surf, Lorene::Map_af::integrale_surface(), Lorene::Spheroid::p_area, and Lorene::Spheroid::sqrt_q().

void Lorene::Spheroid::del_deriv (  )  const [protected, virtual, inherited]
const Tensor & Lorene::Spheroid::delta (  )  const [inherited]
Tensor Lorene::Spheroid::derive_cov2d ( const Tensor uu  )  const [inherited]

Computes the total covariant derivative on the spheroid.

Definition at line 1243 of file spheroid.C.

References Lorene::contract(), Lorene::Spheroid::delta(), Lorene::Spheroid::derive_cov2dflat(), Lorene::Tensor::get_index_type(), and Lorene::Tensor::get_valence().

Tensor Lorene::Spheroid::derive_cov2dflat ( const Tensor uu  )  const [inherited]
double Lorene::Spheroid::epsilon_A_minus_one (  )  const [inherited]

Computation of the refined Penrose parameter for axisymmetric spacetimes, and its difference wrt one.

Definition at line 867 of file spheroid.C.

References Lorene::Spheroid::angu_mom(), Lorene::Spheroid::area(), Lorene::Spheroid::mass(), Lorene::pow(), and Lorene::sqrt().

double Lorene::Spheroid::epsilon_P_minus_one (  )  const [inherited]

Computation of the classical Penrose parameter, and its difference wrt one.

To use in replacement of epsilon_A_minus_one when the computed spacetime is not axisymmetric.

Definition at line 879 of file spheroid.C.

References Lorene::Spheroid::angu_mom(), Lorene::Spheroid::area(), Lorene::Spheroid::mass(), Lorene::Spheroid::p_epsilon_P_minus_one, and Lorene::pow().

const Vector& Lorene::Spheroid::get_ephi (  )  const [inline, inherited]

Returns the conformal Killing symmetry vector on the 2-surface.

Definition at line 253 of file spheroid.h.

References Lorene::Spheroid::ephi.

const Scalar& Lorene::Spheroid::get_fff (  )  const [inline, inherited]

Returns the normalization scalar F.

Definition at line 259 of file spheroid.h.

References Lorene::Spheroid::fff.

const Scalar& Lorene::Spheroid::get_ggg (  )  const [inline, inherited]

Returns the normalization scalar G.

Definition at line 262 of file spheroid.h.

References Lorene::Spheroid::ggg.

const Sym_tensor& Lorene::Spheroid::get_hh (  )  const [inline, inherited]

Returns the symmetric tensor $ H_{ab} $.

Definition at line 232 of file spheroid.h.

References Lorene::Spheroid::hh.

const Scalar& Lorene::Spheroid::get_hsurf (  )  const [inline, inherited]

Returns the field h_surf.

Definition at line 223 of file spheroid.h.

References Lorene::Spheroid::h_surf.

bool Lorene::Spheroid::get_issphere (  )  const [inline, inherited]

Returns the flag saying whether or not the horizon is geometrically round.

Definition at line 265 of file spheroid.h.

References Lorene::Spheroid::issphere.

const Tensor& Lorene::Spheroid::get_jac2d (  )  const [inline, inherited]

returns the 2-d jacobian of coordinate transformation $ J $

Definition at line 241 of file spheroid.h.

References Lorene::Spheroid::jac2d.

const Sym_tensor& Lorene::Spheroid::get_jj (  )  const [inline, inherited]

Returns the symmetric tensor $ J_{ab} $.

Definition at line 256 of file spheroid.h.

References Lorene::Spheroid::jj.

const Vector& Lorene::Spheroid::get_ll (  )  const [inline, inherited]

Returns the vector $ L_a $.

Definition at line 247 of file spheroid.h.

References Lorene::Spheroid::ll.

const Tensor& Lorene::Spheroid::get_proj (  )  const [inline, inherited]

returns the 3-d projector on 2-surface $ Pi $

Definition at line 238 of file spheroid.h.

References Lorene::Spheroid::proj.

const Metric& Lorene::Spheroid::get_qab (  )  const [inline, inherited]

Returns the metric $ q_{ab} $.

Definition at line 226 of file spheroid.h.

const Sym_tensor& Lorene::Spheroid::get_qq (  )  const [inline, inherited]

returns the 3-d degenerate 2-metric $ Q_{ab}$

Definition at line 235 of file spheroid.h.

References Lorene::Spheroid::qq.

const Scalar& Lorene::Spheroid::get_ricci (  )  const [inline, inherited]

Returns the 2-ricci scalar $ R = q^{ab}R_{ab} $.

Definition at line 229 of file spheroid.h.

References Lorene::Spheroid::ricci.

const Vector& Lorene::Spheroid::get_ss (  )  const [inline, inherited]

Returns the vector $ S_a $.

Definition at line 250 of file spheroid.h.

References Lorene::Spheroid::ss.

const Scalar& Lorene::Spheroid::get_trk (  )  const [inline, inherited]

Returns the trace K on the 2-surface.

Definition at line 244 of file spheroid.h.

References Lorene::Spheroid::trk.

const Scalar& Lorene::App_hor::l_non_affinity ( const Scalar bb,
const Scalar lapse 
)

non-affinity (or surface gravity) with respect to the outgoing null vector field

const Sym_tensor& Lorene::App_hor::lie_derive_q_ab ( const Scalar bb,
const Scalar lapse 
)

Lie derivative of 2-metric with respect to the evolution vector field.

const Sym_tensor& Lorene::App_hor::lie_derive_shear ( const Scalar bb,
const Scalar lapse 
)

Lie derivative of shear tensor with respect to the evolution vector field.

const Sym_tensor& Lorene::App_hor::lie_derive_theta_minus ( const Scalar bb,
const Scalar lapse 
)

Lie derivative of the null ingoing expansion rate with respect to the evolution vector field.

const Sym_tensor& Lorene::App_hor::lie_derive_theta_plus ( const Scalar bb,
const Scalar lapse 
)

Lie derivative of the null outgoing expansion rate with respect to the evolution vector field.

double Lorene::Spheroid::mass (  )  const [inherited]

Computes the mass as defined from the calculus of angular momentum, done with respect to a divergence free tangent vector field $ phi $.

Spheroid has to be a real sphere (flag issphere true), of constant radius

\[R_{s} \]

. defined as

\[ M = \frac{1}{2 R_{s}} \sqrt{R_{s}^{4} + 4{\cal J}^{2}} \]

Definition at line 739 of file spheroid.C.

References Lorene::Spheroid::angu_mom(), Lorene::Spheroid::area(), Lorene::Spheroid::p_mass, and Lorene::sqrt().

double Lorene::Spheroid::multipole_angu ( const int  order  )  const [inherited]

Computes the angular multipole of a given order for the spheroid, assumed to be spherical.

$ phi $ is a divergence free tangent vector field. WARNING:order has to be strictly higher than zero (no topological defects here...), and an odd number for technical reasons.

Definition at line 804 of file spheroid.C.

References Lorene::Spheroid::area(), Lorene::contract(), Lorene::Map::get_bvect_spher(), Lorene::Spheroid::get_ephi(), Lorene::Tensor::get_mp(), Lorene::Spheroid::h_surf, Lorene::Map_af::integrale_surface(), Lorene::Spheroid::jac2d, Lorene::Spheroid::ll, Lorene::Spheroid::p_multipole_angu, Lorene::Scalar::set_spectral_va(), Lorene::sqrt(), Lorene::Spheroid::sqrt_q(), Lorene::Scalar::std_spectral_base(), and Lorene::Valeur::ylm().

double Lorene::Spheroid::multipole_mass ( const int  order  )  const [inherited]

Computes the mass multipole of a given order for the spheroid, assumed to be spherical.

WARNING: For technical reasons, only even orders are supported by the code.

Definition at line 751 of file spheroid.C.

References Lorene::Spheroid::area(), Lorene::Tensor::get_mp(), Lorene::Spheroid::get_ricci(), Lorene::Spheroid::h_surf, Lorene::Map_af::integrale_surface(), Lorene::Spheroid::mass(), Lorene::Spheroid::p_multipole_mass, Lorene::Scalar::set_spectral_va(), Lorene::sqrt(), Lorene::Spheroid::sqrt_q(), Lorene::Scalar::std_spectral_base(), and Lorene::Valeur::ylm().

void Lorene::App_hor::operator= ( const App_hor  ) 

Assignment to another App_hor.

Reimplemented from Lorene::Spheroid.

void Lorene::Spheroid::sauve ( FILE *   )  const [virtual, inherited]

Save in a file.

Definition at line 1189 of file spheroid.C.

void Lorene::Spheroid::set_der_0x0 (  )  const [protected, inherited]
void Lorene::Spheroid::set_ephi ( const Scalar  )  [inherited]

Assigns the conformal Killing vector field to phi.

Scalar& Lorene::Spheroid::set_fff (  )  [inline, inherited]

Sets the normalization factor F.

Definition at line 298 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::fff.

Scalar& Lorene::Spheroid::set_ggg (  )  [inline, inherited]

Sets the normalization factor G.

Definition at line 301 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::ggg.

Sym_tensor& Lorene::Spheroid::set_hh (  )  [inline, inherited]

Sets the symmetric tensor $ H_{ab} $.

Definition at line 283 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::hh.

Scalar& Lorene::Spheroid::set_hsurf (  )  [inline, inherited]

Sets the field h_surf.

Definition at line 268 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::h_surf.

bool Lorene::Spheroid::set_issphere (  )  [inline, inherited]

Sets the boolean linked to geometrical shape of the horizon.

Definition at line 304 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::issphere.

Sym_tensor& Lorene::Spheroid::set_jj (  )  [inline, inherited]

Sets the symmetric tensor $ J_{ab} $.

Definition at line 295 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::jj.

Vector& Lorene::Spheroid::set_ll (  )  [inline, inherited]

Sets the vector $ L_a $.

Definition at line 289 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::ll.

Tensor& Lorene::Spheroid::set_proj (  )  [inline, inherited]

Sets the projector $ Pi $.

Definition at line 280 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::proj.

Metric& Lorene::Spheroid::set_qab (  )  [inline, inherited]

Sets the modified metric (non degenerated) $ q_{ab} $.

Definition at line 271 of file spheroid.h.

References Lorene::Spheroid::del_deriv().

Sym_tensor& Lorene::Spheroid::set_qq (  )  [inline, inherited]

Sets the degenerated metric $ Q_{ab} $.

Definition at line 277 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::qq.

Scalar& Lorene::Spheroid::set_ricci (  )  [inline, inherited]

Sets the 2-Ricci scalar $ R = q^{ab}R_{ab} $.

Definition at line 274 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::ricci.

Vector& Lorene::Spheroid::set_ss (  )  [inline, inherited]

Sets the vector $ s_a $.

Definition at line 292 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::ss.

Scalar& Lorene::Spheroid::set_trk (  )  [inline, inherited]

Sets the trace K on the 2-surface.

Definition at line 286 of file spheroid.h.

References Lorene::Spheroid::del_deriv(), and Lorene::Spheroid::trk.

const Sym_tensor & Lorene::Spheroid::shear (  )  const [inherited]
const Scalar & Lorene::Spheroid::sqrt_q (  )  const [inherited]

Computes the normal vector field to the 2-surface.

Computes the square root of the determinant of $ q_{ab} $.

Definition at line 696 of file spheroid.C.

References Lorene::Spheroid::get_qq(), Lorene::Spheroid::p_sqrt_q, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().

const Scalar & Lorene::Spheroid::theta_minus (  )  const [inherited]
const Scalar & Lorene::Spheroid::theta_plus (  )  const [inherited]
void Lorene::Spheroid::update_from_tslice ( const Metric gamij,
const Sym_tensor Kij 
) [inherited]

Updates from the 3-slice data.


Friends And Related Function Documentation

ostream& operator<< ( ostream &  ,
const Spheroid  
) [friend, inherited]

Display.


Member Data Documentation

Vector Lorene::Spheroid::ephi [protected, inherited]

The conformal Killing vector field on the 2-surface (set to by default to the axial vector associated with coordinate $ \varphi $).

Definition at line 113 of file spheroid.h.

Scalar Lorene::Spheroid::fff [protected, inherited]

Normalization function for the outgoing null vector l.

Definition at line 138 of file spheroid.h.

Scalar Lorene::Spheroid::ggg [protected, inherited]

Normalization function for the ingoing null vector k.

Definition at line 143 of file spheroid.h.

Scalar Lorene::Spheroid::h_surf [protected, inherited]

The location of the 2-surface as r = h_surf$(\theta, \varphi)$.

Definition at line 91 of file spheroid.h.

Sym_tensor Lorene::Spheroid::hh [protected, inherited]

The ricci scalar on the 2-surface.

Extrinsic curvature of the 2-surface in the 3-slice. $ H_{ab} $ (covariant representation)

Definition at line 122 of file spheroid.h.

bool Lorene::Spheroid::issphere [protected, inherited]

Flag to know whether the horizon is geometrically round or distorted.

Definition at line 151 of file spheroid.h.

Tensor Lorene::Spheroid::jac2d [protected, inherited]

The jacobian of the adaptation of coordinates (contravariant/covariant representation).

Definition at line 96 of file spheroid.h.

Sym_tensor Lorene::Spheroid::jj [protected, inherited]

Tangent components of the extrinsic curvature of the 3-slice.

$ J_{ab} $ (covariant representation)

Definition at line 134 of file spheroid.h.

Vector Lorene::Spheroid::ll [protected, inherited]

Normal-tangent component of the extrinsic curvature of the 3-slice.

$ L_a $ (covariant representation)

Definition at line 129 of file spheroid.h.

double* Lorene::Spheroid::p_angu_mom [mutable, protected, inherited]

The angular momentum.

Definition at line 159 of file spheroid.h.

double* Lorene::Spheroid::p_area [mutable, protected, inherited]

The area of the 2-surface.

Definition at line 158 of file spheroid.h.

double* Lorene::Spheroid::p_epsilon_P_minus_one [mutable, protected, inherited]

Refined Penrose parameter, difference wrt one.

Definition at line 164 of file spheroid.h.

double* Lorene::Spheroid::p_mass [mutable, protected, inherited]

Mass defined from angular momentum.

Definition at line 160 of file spheroid.h.

double* Lorene::Spheroid::p_multipole_angu [mutable, protected, inherited]

Angular momentum multipole for the spheroid.

Definition at line 162 of file spheroid.h.

double* Lorene::Spheroid::p_multipole_mass [mutable, protected, inherited]

Mass multipole for the spheroid.

Definition at line 161 of file spheroid.h.

Sym_tensor* Lorene::Spheroid::p_shear [mutable, protected, inherited]

The shear tensor.

Definition at line 167 of file spheroid.h.

Scalar* Lorene::Spheroid::p_sqrt_q [mutable, protected, inherited]

Surface element $\sqrt{\det q_{ab}} $.

Definition at line 157 of file spheroid.h.

Scalar* Lorene::Spheroid::p_theta_minus [mutable, protected, inherited]

Null ingoing expansion.

Definition at line 166 of file spheroid.h.

Scalar* Lorene::Spheroid::p_theta_plus [mutable, protected, inherited]

Classical Penrose parameter, difference wrt one.

Null outgoing expansion

Definition at line 165 of file spheroid.h.

Tensor Lorene::Spheroid::proj [protected, inherited]

The 3-d projector on the 2-surface (contravariant-covariant form).

Definition at line 100 of file spheroid.h.

Sym_tensor Lorene::Spheroid::qq [protected, inherited]

The 3-d covariant degenerated 2-metric on the surface.

Definition at line 104 of file spheroid.h.

Scalar Lorene::Spheroid::ricci [protected, inherited]

Induced metric on the 2-surface $ q_{ab} $.

Definition at line 117 of file spheroid.h.

Vector Lorene::Spheroid::ss [protected, inherited]

The adapted normal vector field to spheroid in the 3-slice.

Definition at line 108 of file spheroid.h.

Scalar Lorene::Spheroid::trk [protected, inherited]

Trace K of the extrinsic curvature of the 3-slice.

Definition at line 124 of file spheroid.h.


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

Generated on 7 Dec 2019 for LORENE by  doxygen 1.6.1