LORENE
excised_slice.h
1 /*
2  * Definition of Lorene class Excised slice
3  *
4  */
5 
6 /*
7  * Copyright (c) 2010 Nicolas Vasset
8  */
9 
10 
11 #ifndef __EXCISEDSLICE_H_
12 #define __EXCISEDSLICE_H_
13 
14 
15 // Headers Lorene
16 #include "spheroid.h"
17 
18 //---------------------------//
19 // Class Excised_slice //
20 //---------------------------//
21 
22 namespace Lorene {
42 
43  // Data :
44  // -----
45  protected:
46  const Map& mp ;
47 
48  // Metric data
49  // -----------------
50 
52  int type_hor;
53 
55  int field_set;
56 
59 
60  // Conformal factor
61  Scalar conf_fact;
62 
65 
70 
75 
80 
81 
82  // Derived data :
83  // ------------
84  protected:
85 
86 
88 
89  mutable Spheroid* p_hor ;
91  mutable double* p_adm_mass ;
92 
96  mutable double* p_komar_angmom ;
97 
101  mutable double* p_virial_residue ;
102 
103 
104 
105 
106 
107  // Constructors - Destructor
108  // -------------------------
109  public:
110 
125  Excised_slice(const Map& mp_i, int hor_type, int scheme_type) ;
126 
127 
128  Excised_slice(const Excised_slice& ) ;
129 
145  Excised_slice(const Map& mp_i, int hor_type, int scheme_type, FILE* fich) ;
146 
147  virtual ~Excised_slice() ;
148 
149 
150  // Memory management
151  // -----------------
152  protected:
154  virtual void del_deriv() const ;
155 
157  virtual void set_der_0x0() const ;
158 
159 
160  // Mutators / assignment
161  // ---------------------
162  public:
164  void operator=(const Excised_slice&) ;
165 
183  void compute_stat_metric(double precis, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i=true, double relax=0.2, int mer_max = 2000, int mer_max2=200, bool isvoid = true) ;
184 
190  void secmembre_kerr(Sym_tensor& source_hh);
191 
192 
193  // Accessors
194  // ---------
195  public:
197  const Map& get_mp() const {return mp; } ;
198 
200  int get_type_hor() const {return type_hor; };
201 
203 
204  int get_field_set() const {return field_set;};
205 
207  const Scalar& get_lapse() const {return lapse;} ;
208 
210  const Scalar& get_conf_fact() const {return conf_fact;};
211 
213  const Vector& get_shift() const {return shift;} ;
214 
216  const Sym_tensor& get_hij() const {return hij;} ;
217 
219  const Sym_tensor& get_hatA() const {return hatA;} ;
220 
222  const Vector& get_Xx() const {
223  if (field_set==2) return Xx;
224  else if (field_set==1) cout <<"Error: the scheme used is the original FCF; this variable is irrelevant" << endl;
225  else cout <<"error in the scheme definition; please check the class consistency" << endl;}
226 
228  int set_type_hor() {del_deriv() ; return type_hor ; } ;
229 
231  Scalar& set_lapse() {del_deriv() ; return lapse ; } ;
232 
234  Scalar& set_conf_fact() {del_deriv() ; return conf_fact ; } ;
235 
237  Vector& set_shift() {del_deriv() ; return shift ; } ;
238 
240  Sym_tensor& set_hij() {del_deriv() ; return hij ; } ;
241 
243  Sym_tensor& set_hatA() {del_deriv() ; return hatA ; } ;
244 
246  Vector& set_Xx() {del_deriv() ; return Xx ; } ; // Include an assert here!!!!
247 
248 
249  // Outputs
250  // -------
251  public:
252  virtual void sauve(FILE* ) const ;
253 
254  void Einstein_errors();
255 
256 
257  // Global quantities
258  // -----------------
259  public:
260 
261 
266  Spheroid hor() ;
267 
269  double adm_mass() ;
270 
274  double komar_angmom() ;
275 
279  double virial_residue() ;
280 
281 };
282 
283 }
284 #endif
Vector shift
Shift vector.
Definition: excised_slice.h:64
Sym_tensor hatA
Rescaled tracefree extrinsic curvature tensor: rescaled same way as Cordero et al.
Definition: excised_slice.h:74
void Einstein_errors()
Prints out errors in Einstein equations for the data obtained.
const Scalar & get_conf_fact() const
Returns the conformal factor.
Vector & set_shift()
Sets the shift vector.
Lorene prototypes.
Definition: app_hor.h:67
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
Base class for coordinate mappings.
Definition: map.h:682
const Map & get_mp() const
Returns the mapping.
virtual void sauve(FILE *) const
Save in a file.
Class to compute single black hole spacetime excised slices.
Definition: excised_slice.h:41
Tensor field of valence 1.
Definition: vector.h:188
int field_set
Chose field set type.
Definition: excised_slice.h:55
void secmembre_kerr(Sym_tensor &source_hh)
If hor_type=1, computes the rhs of hyperbolic equation for conformal metric assuming stationarity; WA...
double * p_adm_mass
Computation of the ADM mass of the BH spacetime.
Definition: excised_slice.h:91
void operator=(const Excised_slice &)
Assignment to another Excised_slice.
virtual void del_deriv() const
Deletes all the derived quantities.
const Vector & get_shift() const
Returns the shift vector .
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Definition: spheroid.h:84
double * p_virial_residue
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
double komar_angmom()
Computation of the Komar angular momentum w.r.t.
Vector & set_Xx()
Sets the longitudinal part of Einstein Equations if the scheme is apropriate.
void compute_stat_metric(double precis, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i=true, double relax=0.2, int mer_max=2000, int mer_max2=200, bool isvoid=true)
If hor_type=1, computes a quasi-stationary 3-slice from the chosen parameters.
Vector Xx
Longitudinal part of the extrinsic curvature.
Definition: excised_slice.h:79
Sym_tensor hij
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al...
Definition: excised_slice.h:69
Spheroid * p_hor
Computation of the spheroid associated with the black hole horizon.
Definition: excised_slice.h:89
const Vector & get_Xx() const
Return the longitudinal part of Einstein Equations if the scheme is appropriate.
const Map & mp
Mapping associated with the slice.
Definition: excised_slice.h:46
Spheroid hor()
Spheroid at the excised boundary associated with the black hole MOTS on the slice.
int type_hor
Chosen horizon type.
Definition: excised_slice.h:52
Scalar lapse
Lapse function.
Definition: excised_slice.h:58
Sym_tensor & set_hij()
Sets the deviation tensor.
int get_type_hor() const
Returns the type of horizon that performs the excision.
Sym_tensor & set_hatA()
Sets the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2) ...
double virial_residue()
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
Excised_slice(const Map &mp_i, int hor_type, int scheme_type)
Standard constructor.
Definition: excised_slice.C:58
double adm_mass()
Computation of the ADM mass of the BH spacetime.
Scalar & set_conf_fact()
Sets the conformal factor.
int get_field_set() const
Returns the field set chosen for the data.
const Scalar & get_lapse() const
Returns the lapse function N.
virtual ~Excised_slice()
Destructor.
Scalar & set_lapse()
Sets the lapse.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
int set_type_hor()
Sets the horizon type.
const Sym_tensor & get_hatA() const
Returns the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2) ...
double * p_komar_angmom
Computation of the Komar angular momentum w.r.t.
Definition: excised_slice.h:96
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
const Sym_tensor & get_hij() const
Returns the deviation tensor .