Lorene::Compobj Class Reference
[Stationary compact objects (under development)]

Base class for stationary compact objects (***under development***). More...

#include <compobj.h>

Inheritance diagram for Lorene::Compobj:
Lorene::Compobj_QI Lorene::HiggsMonopole Lorene::ScalarBH Lorene::AltBH_QI Lorene::Kerr_QI Lorene::Star_QI Lorene::Boson_star

List of all members.

Public Member Functions

 Compobj (Map &map_i)
 Standard constructor.
 Compobj (const Compobj &)
 Copy constructor.
 Compobj (Map &map_i, FILE *)
 Constructor from a file (see sauve(FILE* )).
virtual ~Compobj ()
 Destructor.
void operator= (const Compobj &)
 Assignment to another Compobj.
Mapset_mp ()
 Read/write of the mapping.
const Mapget_mp () const
 Returns the mapping.
const Scalarget_nn () const
 Returns the lapse function N .
const Vectorget_beta () const
 Returns the shift vector $\beta^i$.
const Metricget_gamma () const
 Returns the 3-metric $\gamma_{ij}$.
const Scalarget_ener_euler () const
 Returns the total energy density E in the Eulerian frame.
const Vectorget_mom_euler () const
 Returns the total 3-momentum density $P^i$ in the Eulerian frame.
const Sym_tensorget_stress_euler () const
 Returns the stress tensor $S_{ij}$ with respect to the Eulerian observer.
const Sym_tensorget_kk () const
 Returns the extrinsic curvature tensor $K_{ij}$.
virtual void sauve (FILE *) const
 Save in a file.
void gyoto_data (const char *file_name) const
 Save in a file for GYOTO.
virtual void extrinsic_curvature ()
 Computation of the extrinsic curvature.
virtual double adm_mass () const
 ADM mass (computed as a surface integral at spatial infinity).

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.
virtual ostream & operator>> (ostream &) const
 Operator >> (virtual function called by the operator <<).

Protected Attributes

Mapmp
 Mapping describing the coordinate system (r,theta,phi).
Scalar nn
 Lapse function N .
Vector beta
 Shift vector $\beta^i$.
Metric gamma
 3-metric $\gamma_{ij}$
Scalar ener_euler
 Total energy density E in the Eulerian frame.
Vector mom_euler
 Total 3-momentum density $P^i$ in the Eulerian frame.
Sym_tensor stress_euler
 Stress tensor $S_{ij}$ with respect to the Eulerian observer.
Sym_tensor kk
 Extrinsic curvature tensor $K_{ij}$.
double * p_adm_mass
 ADM mass.

Friends

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

Detailed Description

Base class for stationary compact objects (***under development***).

()

A Compobj describes a single compact object (star or black hole), in a stationary state.

The spacetime metric is written according to the 3+1 formalism :

\[ ds^2 = - N^2 dt^2 + \gamma_{ij} ( dx^i + \beta^i dt ) (dx^j + \beta^j dt ) \]

where $\gamma_{ij}$ is the 3-metric, described by a Lorene object of class Metric.

The total energy-momentum tensor is orthogonally split with respect to the Eulerian observer as follows:

\[ T_{\alpha\beta} = E n_\alpha n_\beta + P_\alpha n_\beta + n_\alpha P_\beta + S_{\alpha\beta} \]

Definition at line 129 of file compobj.h.


Constructor & Destructor Documentation

Lorene::Compobj::Compobj ( Map map_i  ) 

Standard constructor.

Parameters:
mp_i Mapping on which the object is defined

Definition at line 85 of file compobj.C.

References beta, ener_euler, kk, mom_euler, nn, set_der_0x0(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::std_spectral_base(), and stress_euler.

Lorene::Compobj::Compobj ( const Compobj co  ) 

Copy constructor.

Definition at line 112 of file compobj.C.

References set_der_0x0().

Lorene::Compobj::Compobj ( Map map_i,
FILE *  fich 
)

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

Parameters:
mp_i Mapping on which the object is defined
fich input file (must have been created by the function sauve)

Definition at line 129 of file compobj.C.

References set_der_0x0().

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

Destructor.

Definition at line 147 of file compobj.C.

References del_deriv().


Member Function Documentation

double Lorene::Compobj::adm_mass (  )  const [virtual]
void Lorene::Compobj::del_deriv (  )  const [protected, virtual]

Deletes all the derived quantities.

Reimplemented in Lorene::Boson_star, Lorene::Compobj_QI, Lorene::Star_QI, Lorene::Kerr_QI, Lorene::AltBH_QI, and Lorene::ScalarBH.

Definition at line 158 of file compobj.C.

References p_adm_mass, and set_der_0x0().

void Lorene::Compobj::extrinsic_curvature (  )  [virtual]

Computation of the extrinsic curvature.

Reimplemented in Lorene::Compobj_QI, and Lorene::AltBH_QI.

Definition at line 293 of file compobj.C.

References beta, Lorene::Tensor::derive_cov(), Lorene::Tensor::down(), gamma, kk, nn, Lorene::Tensor::set(), and Lorene::Tensor::set_etat_qcq().

const Vector& Lorene::Compobj::get_beta (  )  const [inline]

Returns the shift vector $\beta^i$.

Definition at line 216 of file compobj.h.

References beta.

const Scalar& Lorene::Compobj::get_ener_euler (  )  const [inline]

Returns the total energy density E in the Eulerian frame.

Definition at line 222 of file compobj.h.

References ener_euler.

const Metric& Lorene::Compobj::get_gamma (  )  const [inline]

Returns the 3-metric $\gamma_{ij}$.

Definition at line 219 of file compobj.h.

References gamma.

const Sym_tensor& Lorene::Compobj::get_kk (  )  const [inline]

Returns the extrinsic curvature tensor $K_{ij}$.

Definition at line 231 of file compobj.h.

References kk.

const Vector& Lorene::Compobj::get_mom_euler (  )  const [inline]

Returns the total 3-momentum density $P^i$ in the Eulerian frame.

Definition at line 225 of file compobj.h.

References mom_euler.

const Map& Lorene::Compobj::get_mp (  )  const [inline]

Returns the mapping.

Definition at line 210 of file compobj.h.

References mp.

const Scalar& Lorene::Compobj::get_nn (  )  const [inline]

Returns the lapse function N .

Definition at line 213 of file compobj.h.

References nn.

const Sym_tensor& Lorene::Compobj::get_stress_euler (  )  const [inline]

Returns the stress tensor $S_{ij}$ with respect to the Eulerian observer.

Definition at line 228 of file compobj.h.

References stress_euler.

void Lorene::Compobj::gyoto_data ( const char *  file_name  )  const
void Lorene::Compobj::operator= ( const Compobj co  ) 

Assignment to another Compobj.

Reimplemented in Lorene::Boson_star, Lorene::Compobj_QI, Lorene::Star_QI, Lorene::Kerr_QI, Lorene::AltBH_QI, and Lorene::ScalarBH.

Definition at line 178 of file compobj.C.

References beta, del_deriv(), ener_euler, gamma, kk, mom_euler, mp, nn, and stress_euler.

ostream & Lorene::Compobj::operator>> ( ostream &  ost  )  const [protected, virtual]

Operator >> (virtual function called by the operator <<).

Reimplemented in Lorene::Boson_star, Lorene::Compobj_QI, Lorene::Star_QI, Lorene::Kerr_QI, Lorene::AltBH_QI, Lorene::ScalarBH, and Lorene::HiggsMonopole.

Definition at line 242 of file compobj.C.

References Lorene::Metric::cov(), ener_euler, gamma, kk, mp, nn, stress_euler, and Lorene::Scalar::val_grid_point().

void Lorene::Compobj::sauve ( FILE *  fich  )  const [virtual]
void Lorene::Compobj::set_der_0x0 (  )  const [protected]

Sets to 0x0 all the pointers on derived quantities.

Reimplemented in Lorene::Boson_star, Lorene::Compobj_QI, Lorene::Star_QI, Lorene::Kerr_QI, Lorene::AltBH_QI, and Lorene::ScalarBH.

Definition at line 166 of file compobj.C.

References p_adm_mass.

Map& Lorene::Compobj::set_mp (  )  [inline]

Read/write of the mapping.

Definition at line 203 of file compobj.h.

References mp.


Friends And Related Function Documentation

ostream& operator<< ( ostream &  ,
const Compobj  
) [friend]

Display.


Member Data Documentation

Shift vector $\beta^i$.

Definition at line 141 of file compobj.h.

Total energy density E in the Eulerian frame.

Definition at line 147 of file compobj.h.

3-metric $\gamma_{ij}$

Definition at line 144 of file compobj.h.

Extrinsic curvature tensor $K_{ij}$.

Definition at line 156 of file compobj.h.

Total 3-momentum density $P^i$ in the Eulerian frame.

Definition at line 150 of file compobj.h.

Map& Lorene::Compobj::mp [protected]

Mapping describing the coordinate system (r,theta,phi).

Definition at line 135 of file compobj.h.

Lapse function N .

Definition at line 138 of file compobj.h.

double* Lorene::Compobj::p_adm_mass [mutable, protected]

ADM mass.

Definition at line 161 of file compobj.h.

Stress tensor $S_{ij}$ with respect to the Eulerian observer.

Definition at line 153 of file compobj.h.


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

Generated on 7 Dec 2019 for LORENE by  doxygen 1.6.1