LORENE
isol_hole.h
1 /*
2  * Definition of Lorene class Isol_hole
3  *
4  */
5 
6 /*
7  * Copyright (c) 2009 Nicolas Vasset
8  */
9 
10 
11 #ifndef __ISOLHOLE_H_
12 #define __ISOLHOLE_H_
13 
14 
15 // Headers Lorene
16 #include "tensor.h"
17 #include "metric.h"
18 #include "spheroid.h"
19 #include "excision_surf.h"
20 #include "excision_hor.h"
21 
22  class Eos ;
23 
24  //---------------------------//
25  // Class Isol_hole //
26  //---------------------------//
27 
28 namespace Lorene {
50 class Isol_hole {
51 
52  // Data :
53  // -----
54  protected:
55  const Map& mp ;
56 
57  double Omega ;
61  bool NorKappa ;
69  bool isCF ;
70 
71 
72  // Metric data
73  // -----------------
74 
77 
78  // Conformal factor
79  Scalar conf_fact;
80 
83 
88 
93 
94  // Derived data :
95  // ------------
96  protected:
97 
98 
100 
101  mutable Spheroid* p_hor ;
103  mutable double* p_adm_mass ;
104 
108  mutable double* p_komar_angmom ;
109 
113  mutable double* p_virial_residue ;
114 
115 
116 
117 
118 
119  // Constructors - Destructor
120  // -------------------------
121  public:
122 
135  Isol_hole(const Map& mp_i, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i = false) ;
136 
137 
138  Isol_hole(const Isol_hole& ) ;
139 
154  Isol_hole(const Map& mp_i, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i, FILE* fich) ;
155 
156  virtual ~Isol_hole() ;
157 
158 
159  // Memory management
160  // -----------------
161  protected:
163  virtual void del_deriv() const ;
164 
166  virtual void set_der_0x0() const ;
167 
168 
169  // Mutators / assignment
170  // ---------------------
171  public:
173  void operator=(const Isol_hole&) ;
174 
181  void compute_stat_metric(double precis, double relax, int mer_max, int mer_max2, bool isvoid = true) ;
182 
188  void secmembre_kerr(Sym_tensor& source_hh);
189 
190 
191 
192  // Accessors
193  // ---------
194  public:
196  const Map& get_mp() const {return mp; } ;
197 
199  double get_Omega() const {return Omega ;} ;
200 
202  const Scalar& get_boundN() const{
203  if (NorKappa == false){
204  return boundNoK;
205  }
206  else cout << "The boundary condition is imposed on the surface gravity!" << endl;
207  }
209 
210  const Scalar& get_Kappa() const{
211  if(NorKappa == true){
212  return boundNoK;
213  }
214  else cout << "The boundary condition is imposed on the lapse!" <<endl;
215  }
216 
217 
219  const Scalar& get_lapse() const {return lapse;} ;
220 
222  const Scalar& get_conf_fact() const {return conf_fact;};
223 
225  const Vector& get_shift() const {return shift;} ;
226 
228  const Sym_tensor& get_hij() const {return hij;} ;
229 
231  const Sym_tensor& get_hatA() const {return hatA;} ;
232 
233 
234 
235 
236  // Outputs
237  // -------
238  public:
239  virtual void sauve(FILE* ) const ;
240 
241  void Einstein_errors();
242 
243 
244  // Global quantities
245  // -----------------
246  public:
247 
248 
253  Spheroid hor() ;
254 
256  double adm_mass() ;
257 
261  double komar_angmom() ;
262 
266  double virial_residue() ;
267 
268 
269 
270 };
271 
272 }
273 #endif
void compute_stat_metric(double precis, double relax, int mer_max, int mer_max2, bool isvoid=true)
Computes a quasi-stationary 3-slice from the chosen parameters.
const Map & get_mp() const
Returns the mapping.
Definition: isol_hole.h:196
double * p_adm_mass
Computation of the ADM mass of the BH spacetime.
Definition: isol_hole.h:103
Class to compute quasistationary single black hole spacetimes in vacuum.
Definition: isol_hole.h:50
double komar_angmom()
Computation of the Komar angular momentum w.r.t.
Definition: isol_hole.C:309
double get_Omega() const
Returns the rotation rate.
Definition: isol_hole.h:199
Lorene prototypes.
Definition: app_hor.h:67
const Scalar & get_boundN() const
Returns the boundary value used for the lapse (if it is the one used)
Definition: isol_hole.h:202
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
Base class for coordinate mappings.
Definition: map.h:688
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: isol_hole.C:145
void secmembre_kerr(Sym_tensor &source_hh)
Computes the rhs of hyperbolic equation for conformal metric assuming statioarity; WARNING; up to now...
Sym_tensor hij
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al...
Definition: isol_hole.h:87
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: isol_hole.C:136
Tensor field of valence 1.
Definition: vector.h:188
double * p_komar_angmom
Computation of the Komar angular momentum w.r.t.
Definition: isol_hole.h:108
const Vector & get_shift() const
Returns the shift vector .
Definition: isol_hole.h:225
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Definition: spheroid.h:84
Scalar boundNoK
Indicates if the boundary value for the lapse or the surface gravity is used.
Definition: isol_hole.h:65
void Einstein_errors()
Prints out errors in Einstein equations for the data obtained.
Definition: isol_hole.C:194
const Scalar & get_Kappa() const
Returns the surface gravity value at the boundary (if it is the one used)
Definition: isol_hole.h:210
Vector shift
Shift vector.
Definition: isol_hole.h:82
const Scalar & get_lapse() const
Returns the lapse function N.
Definition: isol_hole.h:219
const Map & mp
Mapping associated with the star.
Definition: isol_hole.h:55
Spheroid * p_hor
Computation of the spheroid associated with the black hole horizon.
Definition: isol_hole.h:101
double * p_virial_residue
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
Definition: isol_hole.h:113
const Sym_tensor & get_hatA() const
Returns the rescaled tracefree extrinsic curvature .
Definition: isol_hole.h:231
Spheroid hor()
Spheroid at the excised boundary associated with the black hole MOTS on the slice.
Definition: isol_hole.C:247
bool NorKappa
Rotation rate of the horizon in the azimuthal direction.
Definition: isol_hole.h:61
Scalar lapse
Lapse function.
Definition: isol_hole.h:76
double adm_mass()
Computation of the ADM mass of the BH spacetime.
Definition: isol_hole.C:289
double virial_residue()
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
Definition: isol_hole.C:334
void operator=(const Isol_hole &)
Assignment to another Isol_hole.
Definition: isol_hole.C:160
virtual void sauve(FILE *) const
Save in a file.
Definition: isol_hole.C:183
const Sym_tensor & get_hij() const
Returns the deviation tensor .
Definition: isol_hole.h:228
virtual ~Isol_hole()
Destructor.
Definition: isol_hole.C:125
bool isCF
Stores the boundary value of the lapse or surface gravity.
Definition: isol_hole.h:69
Isol_hole(const Map &mp_i, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i=false)
Standard constructor.
Definition: isol_hole.C:51
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
Sym_tensor hatA
Rescaled tracefree extrinsic curvature tensor: see Cordero et al.
Definition: isol_hole.h:92
const Scalar & get_conf_fact() const
Returns the conformal factor.
Definition: isol_hole.h:222