LORENE
spheroid.h
1 /*
2  * Definition of Lorene class Spheroid
3  *
4  */
5 
6 /*
7  * Copyright (c) 2006 Nicolas Vasset, Jerome Novak & Jose-Luis Jaramillo
8  *
9  * This file is part of LORENE.
10  *
11  * LORENE is free software; you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License version 2
13  * as published by the Free Software Foundation.
14  *
15  * LORENE is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the GNU General Public License
21  * along with LORENE; if not, write to the Free Software
22  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23  *
24  */
25 
26 #ifndef __SPHEROID_H_
27 #define __SPHEROID_H_
28 
29 /*
30  * $Id: spheroid.h,v 1.13 2014/10/13 08:52:36 j_novak Exp $
31  * $Log: spheroid.h,v $
32  * Revision 1.13 2014/10/13 08:52:36 j_novak
33  * Lorene classes and functions now belong to the namespace Lorene.
34  *
35  * Revision 1.12 2008/12/10 13:53:56 jl_jaramillo
36  * versions developed at Meudon in Novemver 2008
37  *
38  * Revision 1.11 2008/11/12 15:18:32 n_vasset
39  * New definition for the computation of Ricci scalar (instead of Ricci
40  * tensor previously)
41  *
42  * Revision 1.10 2008/07/09 08:48:23 n_vasset
43  * protected member zeta added.
44  *
45  * Revision 1.9 2008/06/04 12:31:56 n_vasset
46  * new functions multipole_mass and multipole_angu.
47  *
48  * Revision 1.8 2008/06/03 14:32:23 n_vasset
49  * dzpuis corrections (provisory). New function mass implemented.
50  *
51  * Revision 1.7 2006/09/07 08:39:43 j_novak
52  * Minor changes.
53  *
54  * Revision 1.6 2006/07/03 10:12:51 n_vasset
55  * adding of flag isphere
56  *
57  * Revision 1.5 2006/06/13 08:07:55 n_vasset
58  * Addition of 2 members
59  *
60  * Revision 1.4 2006/06/08 09:50:49 n_vasset
61  * corrected version for coordinates transformation
62  *
63  * Revision 1.3 2006/06/07 14:32:24 n_vasset
64  * modified constructor for 3d projector (Sym_tensor to Tensor)
65  *
66  * Revision 1.2 2006/05/26 13:25:42 j_novak
67  * Removed a 'const'
68  *
69  * Revision 1.1 2006/05/26 13:20:42 j_novak
70  * New class spheroid.
71  *
72  *
73  *
74  * $Header: /cvsroot/Lorene/C++/Include/spheroid.h,v 1.13 2014/10/13 08:52:36 j_novak Exp $
75  *
76  */
77 #include "metric.h"
78 
79 namespace Lorene {
84 class Spheroid {
85 
86 
87  // Data :
88  // -----
89  protected:
92 
97 
101 
105 
109 
114 
115  Metric qab ;
116 
118 
123 
125 
130 
135 
139 
144 
145  /*** */
146 
147  Scalar zeta;
148 
151  bool issphere ;
152 
153 
154  // Derived data :
155  // ------------
156  protected:
157  mutable Scalar* p_sqrt_q ;
158  mutable double* p_area ;
159  mutable double* p_angu_mom ;
160  mutable double* p_mass ;
161  mutable double* p_multipole_mass ;
162  mutable double* p_multipole_angu ;
163  mutable double* p_epsilon_A_minus_one ;
164  mutable double* p_epsilon_P_minus_one ;
165  mutable Scalar* p_theta_plus ;
166  mutable Scalar* p_theta_minus ;
167  mutable Sym_tensor* p_shear ;
168  mutable Tensor* p_delta ;
169 
170 
171 
172  // Constructors - Destructor
173  // -------------------------
174  public:
179  Spheroid(const Map_af& map, double radius) ;
189  Spheroid(const Scalar& h_in, const Metric& gamij, const Sym_tensor& Kij) ;
190  Spheroid(const Spheroid& ) ;
191 
193  Spheroid(FILE* ) ;
194 
195  virtual ~Spheroid() ;
196 
197 
198  // Memory management
199  // -----------------
200  protected:
202  virtual void del_deriv() const ;
203 
205  void set_der_0x0() const ;
206 
207 
208  // Mutators / assignment
209  // ---------------------
210  public:
212  void operator=(const Spheroid&) ;
213 
214 
216  void set_ephi(const Scalar&) ;
217 
218 
219  // Accessors
220  // ---------
221  public:
223  const Scalar& get_hsurf() const {return h_surf ; } ;
224 
226  const Metric& get_qab() const {return qab ; } ;
227 
229  const Scalar& get_ricci() const {return ricci ; } ;
230 
232  const Sym_tensor& get_hh() const {return hh ; } ;
233 
235  const Sym_tensor& get_qq() const {return qq ; } ;
236 
238  const Tensor& get_proj() const {return proj ; } ;
239 
241  const Tensor& get_jac2d() const {return jac2d ; } ;
242 
244  const Scalar& get_trk() const {return trk; } ;
245 
247  const Vector& get_ll() const {return ll;} ;
248 
250  const Vector& get_ss() const {return ss;} ;
251 
253  const Vector& get_ephi() const {return ephi;};
254 
256  const Sym_tensor& get_jj() const {return jj;} ;
257 
259  const Scalar& get_fff() const {return fff;} ;
260 
262  const Scalar& get_ggg() const {return ggg;} ;
263 
265  bool get_issphere() const{return issphere;};
266 
268  Scalar& set_hsurf() {del_deriv() ; return h_surf ; } ;
269 
271  Metric& set_qab() {del_deriv() ; return qab ; } ;
272 
274  Scalar& set_ricci() {del_deriv() ; return ricci ; } ;
275 
277  Sym_tensor& set_qq() {del_deriv() ; return qq ; } ;
278 
280  Tensor& set_proj() {del_deriv() ; return proj ; } ;
281 
283  Sym_tensor& set_hh() {del_deriv() ; return hh ; } ;
284 
286  Scalar& set_trk() {del_deriv() ; return trk; } ;
287 
289  Vector& set_ll() {del_deriv() ; return ll;} ;
290 
292  Vector& set_ss() {del_deriv() ; return ss;} ;
293 
295  Sym_tensor& set_jj() {del_deriv() ; return jj;} ;
296 
298  Scalar& set_fff() {del_deriv() ; return fff;} ;
299 
301  Scalar& set_ggg() {del_deriv() ; return ggg;} ;
302 
304  bool set_issphere() {del_deriv() ; return issphere; };
305 
307  void update_from_tslice(const Metric& gamij, const Sym_tensor& Kij) ;
308 
309  // Computational functions
310  // -----------------------
311  public:
312 
314 // Vector compute_normal(const Metric& gamij) ;
315 
317  const Scalar& sqrt_q() const ;
318 
323  double area() const ;
324 
332  double angu_mom() const ;
333 
340  double mass() const;
341 
347  double multipole_mass(const int order) const;
348 
349 
355  double multipole_angu(const int order) const;
356 
357 
358 
360  double epsilon_A_minus_one() const;
361 
365  double epsilon_P_minus_one() const;
366 
368  const Scalar& theta_plus() const ;
369 
371  const Scalar& theta_minus() const ;
372 
374  const Sym_tensor& shear () const ;
375 
377 
378  Tensor derive_cov2dflat(const Tensor& uu) const ;
379 
381 
382  const Tensor& delta() const;
383 
385 
386  Tensor derive_cov2d(const Tensor& uu) const ;
387 
388  // Outputs
389  // -------
390  public:
391  virtual void sauve(FILE *) const ;
392 
394  friend ostream& operator<<(ostream& , const Spheroid& ) ;
395 
396 };
397 
398 ostream& operator<<(ostream& , const Spheroid& ) ;
399 
400  //------------------
401  // Class App_hor
402  //------------------
403 
408 class App_hor : public Spheroid {
409 
410  // Constructors - Destructor
411  // -------------------------
412  public:
413  App_hor(const Mg3d& grid_angu, double radius) ;
414 
422  App_hor(const Scalar& h_in, const Metric& gamij, const Sym_tensor& Kij) ;
423  App_hor(const App_hor& ) ;
424 
426  App_hor(FILE* ) ;
427 
428  virtual ~App_hor() ;
429 
430  // Mutators / assignment
431  // ---------------------
432  public:
434  void operator=(const App_hor&) ;
435 
436  bool check_expansion(double thres = 1.e-7) const ;
437 
439  const Sym_tensor& lie_derive_shear(const Scalar& bb, const Scalar& lapse) ;
440 
444  const Sym_tensor& lie_derive_theta_plus(const Scalar& bb, const Scalar& lapse);
445 
449  const Sym_tensor& lie_derive_theta_minus(const Scalar& bb, const Scalar& lapse);
450 
452  const Sym_tensor& lie_derive_q_ab(const Scalar& bb, const Scalar& lapse) ;
453 
457  const Scalar& l_non_affinity(const Scalar& bb, const Scalar& lapse) ;
458 
459 } ;
460 
461 
462 }
463 #endif
Sym_tensor * p_shear
The shear tensor.
Definition: spheroid.h:167
Scalar h_surf
The location of the 2-surface as r = h_surf .
Definition: spheroid.h:91
const Scalar & theta_plus() const
Computes the outgoing null expansion .
Definition: spheroid.C:892
Tensor & set_proj()
Sets the projector .
Definition: spheroid.h:280
Metric for tensor calculation.
Definition: metric.h:90
double * p_epsilon_P_minus_one
Refined Penrose parameter, difference wrt one.
Definition: spheroid.h:164
double multipole_mass(const int order) const
Computes the mass multipole of a given order for the spheroid, assumed to be spherical.
Definition: spheroid.C:751
double * p_area
The area of the 2-surface.
Definition: spheroid.h:158
double angu_mom() const
Computes the angular momentum with respect to a divergence-free vector field tangent to the 2-surface...
Definition: spheroid.C:724
Tensor jac2d
The jacobian of the adaptation of coordinates (contravariant/covariant representation) ...
Definition: spheroid.h:96
const Scalar & sqrt_q() const
Computes the normal vector field to the 2-surface.
Definition: spheroid.C:696
Scalar & set_hsurf()
Sets the field h_surf.
Definition: spheroid.h:268
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: spheroid.C:653
double * p_angu_mom
The angular momentum.
Definition: spheroid.h:159
const Scalar & get_fff() const
Returns the normalization scalar F.
Definition: spheroid.h:259
double multipole_angu(const int order) const
Computes the angular multipole of a given order for the spheroid, assumed to be spherical.
Definition: spheroid.C:804
Lorene prototypes.
Definition: app_hor.h:67
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
bool get_issphere() const
Returns the flag saying whether or not the horizon is geometrically round.
Definition: spheroid.h:265
Vector ephi
The conformal Killing vector field on the 2-surface (set to by default to the axial vector associated...
Definition: spheroid.h:113
Sym_tensor & set_jj()
Sets the symmetric tensor .
Definition: spheroid.h:295
Sym_tensor & set_hh()
Sets the symmetric tensor .
Definition: spheroid.h:283
const Sym_tensor & shear() const
Computes the shear of the 2-surface .
Definition: spheroid.C:928
const Scalar & l_non_affinity(const Scalar &bb, const Scalar &lapse)
non-affinity (or surface gravity) with respect to the outgoing null vector field
Tensor field of valence 1.
Definition: vector.h:188
Scalar * p_theta_minus
Null ingoing expansion.
Definition: spheroid.h:166
Sym_tensor & set_qq()
Sets the degenerated metric .
Definition: spheroid.h:277
Scalar & set_trk()
Sets the trace K on the 2-surface.
Definition: spheroid.h:286
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Definition: spheroid.h:84
Spheroid(const Map_af &map, double radius)
The delta tensorial fields linked to Christoffel symbols.
Definition: spheroid.C:97
const Vector & get_ephi() const
Returns the conformal Killing symmetry vector on the 2-surface.
Definition: spheroid.h:253
Tensor derive_cov2dflat(const Tensor &uu) const
Computes the round covariant derivative on the spheroid.
Definition: spheroid.C:957
Tensor derive_cov2d(const Tensor &uu) const
Computes the total covariant derivative on the spheroid.
Definition: spheroid.C:1243
const Tensor & get_jac2d() const
returns the 2-d jacobian of coordinate transformation
Definition: spheroid.h:241
bool issphere
Flag to know whether the horizon is geometrically round or distorted.
Definition: spheroid.h:151
const Sym_tensor & get_hh() const
Returns the symmetric tensor .
Definition: spheroid.h:232
Scalar trk
Trace K of the extrinsic curvature of the 3-slice.
Definition: spheroid.h:124
Vector ll
Normal-tangent component of the extrinsic curvature of the 3-slice.
Definition: spheroid.h:129
double * p_multipole_angu
Angular momentum multipole for the spheroid.
Definition: spheroid.h:162
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: spheroid.C:669
void operator=(const Spheroid &)
Assignment to another Spheroid.
Definition: spheroid.C:617
virtual void sauve(FILE *) const
Save in a file.
Definition: spheroid.C:1189
double epsilon_P_minus_one() const
Computation of the classical Penrose parameter, and its difference wrt one.
Definition: spheroid.C:879
const Scalar & get_ricci() const
Returns the 2-ricci scalar .
Definition: spheroid.h:229
const Scalar & get_ggg() const
Returns the normalization scalar G.
Definition: spheroid.h:262
double area() const
Computes the area of the 2-surface.
Definition: spheroid.C:711
Tensor handling.
Definition: tensor.h:294
Sym_tensor jj
Tangent components of the extrinsic curvature of the 3-slice.
Definition: spheroid.h:134
App_hor(const Mg3d &grid_angu, double radius)
Standard constructor.
const Sym_tensor & 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...
virtual ~App_hor()
Destructor.
Vector & set_ss()
Sets the vector .
Definition: spheroid.h:292
Scalar ggg
Normalization function for the ingoing null vector k.
Definition: spheroid.h:143
Multi-domain grid.
Definition: grilles.h:279
Affine radial mapping.
Definition: map.h:2048
const Sym_tensor & lie_derive_shear(const Scalar &bb, const Scalar &lapse)
Lie derivative of shear tensor with respect to the evolution vector field.
const Scalar & theta_minus() const
Computes the ingoing null expansion .
Definition: spheroid.C:912
Sym_tensor hh
The ricci scalar on the 2-surface.
Definition: spheroid.h:122
Vector & set_ll()
Sets the vector .
Definition: spheroid.h:289
const Sym_tensor & get_qq() const
returns the 3-d degenerate 2-metric
Definition: spheroid.h:235
const Scalar & get_trk() const
Returns the trace K on the 2-surface.
Definition: spheroid.h:244
const Vector & get_ss() const
Returns the vector .
Definition: spheroid.h:250
friend ostream & operator<<(ostream &, const Spheroid &)
Display.
void update_from_tslice(const Metric &gamij, const Sym_tensor &Kij)
Updates from the 3-slice data.
const Vector & get_ll() const
Returns the vector .
Definition: spheroid.h:247
Scalar fff
Normalization function for the outgoing null vector l.
Definition: spheroid.h:138
const Metric & get_qab() const
Returns the metric .
Definition: spheroid.h:226
const Sym_tensor & 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 & 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...
Scalar & set_ricci()
Sets the 2-Ricci scalar .
Definition: spheroid.h:274
Metric & set_qab()
Sets the modified metric (non degenerated) .
Definition: spheroid.h:271
double * p_mass
Mass defined from angular momentum.
Definition: spheroid.h:160
const Tensor & get_proj() const
returns the 3-d projector on 2-surface
Definition: spheroid.h:238
Vector ss
The adapted normal vector field to spheroid in the 3-slice.
Definition: spheroid.h:108
double * p_multipole_mass
Mass multipole for the spheroid.
Definition: spheroid.h:161
Scalar & set_ggg()
Sets the normalization factor G.
Definition: spheroid.h:301
Tensor proj
The 3-d projector on the 2-surface (contravariant-covariant form).
Definition: spheroid.h:100
bool set_issphere()
Sets the boolean linked to geometrical shape of the horizon.
Definition: spheroid.h:304
Scalar & set_fff()
Sets the normalization factor F.
Definition: spheroid.h:298
void operator=(const App_hor &)
Assignment to another App_hor.
const Sym_tensor & get_jj() const
Returns the symmetric tensor .
Definition: spheroid.h:256
Scalar * p_theta_plus
Classical Penrose parameter, difference wrt one.
Definition: spheroid.h:165
Scalar ricci
Induced metric on the 2-surface .
Definition: spheroid.h:117
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
double epsilon_A_minus_one() const
Computation of the refined Penrose parameter for axisymmetric spacetimes, and its difference wrt one...
Definition: spheroid.C:867
Class describing an apparent horizon.
Definition: spheroid.h:408
Scalar * p_sqrt_q
Surface element .
Definition: spheroid.h:157
Sym_tensor qq
The 3-d covariant degenerated 2-metric on the surface.
Definition: spheroid.h:104
double mass() const
Computes the mass as defined from the calculus of angular momentum, done with respect to a divergence...
Definition: spheroid.C:739
const Scalar & get_hsurf() const
Returns the field h_surf.
Definition: spheroid.h:223
void set_ephi(const Scalar &)
Assigns the conformal Killing vector field to phi.
virtual ~Spheroid()
Destructor.
Definition: spheroid.C:645
const Tensor & delta() const
Computes the delta coefficients for covariant derivative.
Definition: spheroid.C:1206