LORENE
etoile.h
1/*
2 * Definition of Lorene classes Etoile
3 * Etoile_bin
4 * Etoile_rot
5 *
6 */
7
8/*
9 * Copyright (c) 2000-2001 Eric Gourgoulhon
10 * Copyright (c) 2000-2001 Keisuke Taniguchi
11 *
12 * This file is part of LORENE.
13 *
14 * LORENE is free software; you can redistribute it and/or modify
15 * it under the terms of the GNU General Public License as published by
16 * the Free Software Foundation; either version 2 of the License, or
17 * (at your option) any later version.
18 *
19 * LORENE is distributed in the hope that it will be useful,
20 * but WITHOUT ANY WARRANTY; without even the implied warranty of
21 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22 * GNU General Public License for more details.
23 *
24 * You should have received a copy of the GNU General Public License
25 * along with LORENE; if not, write to the Free Software
26 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
27 *
28 */
29
30
31#ifndef __ETOILE_H_
32#define __ETOILE_H_
33
34/*
35 * $Id: etoile.h,v 1.36 2023/07/04 08:59:53 j_novak Exp $
36 * $Log: etoile.h,v $
37 * Revision 1.36 2023/07/04 08:59:53 j_novak
38 * Added method "r_circ_merid()" to compute the circumferential meridional radius.
39 *
40 * Revision 1.35 2015/06/12 12:38:24 j_novak
41 * Implementation of the corrected formula for the quadrupole momentum.
42 *
43 * Revision 1.34 2015/06/11 13:50:19 j_novak
44 * Minor corrections
45 *
46 * Revision 1.33 2015/06/10 14:36:39 a_sourie
47 * Corrected the formula for the computation of the quadrupole momentum.
48 *
49 * Revision 1.32 2014/10/13 08:52:34 j_novak
50 * Lorene classes and functions now belong to the namespace Lorene.
51 *
52 * Revision 1.31 2011/10/06 14:55:36 j_novak
53 * equation_of_state() is now virtual to be able to call to the magnetized
54 * Eos_mag.
55 *
56 * Revision 1.30 2010/02/02 21:05:49 e_gourgoulhon
57 * Etoile_bin:equilibrium : suppressed the argument method_khi added by
58 * mistake during previous commit.
59 *
60 * Revision 1.29 2010/02/02 13:34:12 e_gourgoulhon
61 * Marked DEPRECATED (in the documentation).
62 *
63 * Revision 1.28 2008/11/14 13:51:08 e_gourgoulhon
64 * Added the parameter ent_limit to Etoile::equilibrium_spher and
65 * Etoile_bin::equilibrium.
66 *
67 * Revision 1.27 2005/10/05 15:14:47 j_novak
68 * Added a Param* as parameter of Etoile_rot::equilibrium
69 *
70 * Revision 1.26 2005/08/29 15:10:12 p_grandclement
71 * Addition of things needed :
72 * 1) For BBH with different masses
73 * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
74 * WORKING YET !!!
75 *
76 * Revision 1.25 2005/01/18 22:35:51 k_taniguchi
77 * Delete a pointer for Etoile::ray_eq(int kk).
78 *
79 * Revision 1.24 2005/01/18 20:34:14 k_taniguchi
80 * Addition of Etoile::ray_eq(int kk).
81 *
82 * Revision 1.23 2005/01/17 20:39:32 k_taniguchi
83 * Addition of Etoile::ray_eq_3pis2().
84 *
85 * Revision 1.22 2004/11/30 20:40:06 k_taniguchi
86 * Addition of the method for calculating a spherical star with falloff
87 * condition at the outer boundary.
88 *
89 * Revision 1.21 2004/05/10 10:18:33 f_limousin
90 * Change to avoid warning in the compilation of Lorene
91 *
92 * Revision 1.20 2004/05/07 12:37:12 f_limousin
93 * Add new member ssjm1_psi
94 *
95 * Revision 1.19 2004/04/08 16:42:31 f_limousin
96 * Add a function velocity_potential with argument ssjm1_psi for the
97 * case of strange stars.
98 *
99 * Revision 1.18 2004/03/22 13:12:41 j_novak
100 * Modification of comments to use doxygen instead of doc++
101 *
102 * Revision 1.17 2003/10/24 12:24:41 k_taniguchi
103 * Suppress the methods of update metric for NS-BH
104 *
105 * Revision 1.16 2003/10/21 11:44:43 k_taniguchi
106 * Delete various things for the Bin_ns_bh project.
107 * They are moved to et_bin_nsbh.h.
108 *
109 * Revision 1.15 2003/10/20 14:50:04 k_taniguchi
110 * Addition of various things for the Bin_ns_bh project
111 * which are related with the part of the neutron star.
112 *
113 * Revision 1.14 2003/10/20 13:11:03 k_taniguchi
114 * Back to version 1.12
115 *
116 * Revision 1.13 2003/10/20 12:15:55 k_taniguchi
117 * Addition of various things for the Bin_ns_bh project
118 * which are related with the part of the neutron star.
119 *
120 * Revision 1.12 2003/06/20 14:13:16 f_limousin
121 * Change to virtual the functions equilibrium_spher() and kinematics().
122 *
123 * Revision 1.11 2003/02/13 16:40:24 p_grandclement
124 * Addition of various things for the Bin_ns_bh project, non of them being
125 * completely tested
126 *
127 * Revision 1.10 2003/02/04 16:20:35 f_limousin
128 * Change to virtual the routine extrinsic_curvature
129 *
130 * Revision 1.9 2003/01/31 16:57:12 p_grandclement
131 * addition of the member Cmp decouple used to compute the K_ij auto, once
132 * the K_ij total is known
133 *
134 * Revision 1.8 2002/12/19 14:48:00 e_gourgoulhon
135 *
136 * Class Etoile_bin: added the new functions:
137 * void update_metric(const Bhole& comp)
138 * void update_metric_der_comp(const Bhole& comp)
139 * to treat the case where the companion is a black hole
140 *
141 * Revision 1.7 2002/12/17 21:17:08 e_gourgoulhon
142 * Class Etoile_bin: suppression of the member p_companion
143 * as well as the associated functions set_companion
144 * and get_companion.
145 *
146 * Revision 1.6 2002/09/13 09:17:33 j_novak
147 * Modif. commentaires
148 *
149 * Revision 1.5 2002/06/17 14:05:16 j_novak
150 * friend functions are now also declared outside the class definition
151 *
152 * Revision 1.4 2002/04/05 09:09:36 j_novak
153 * The inversion of the EOS for 2-fluids polytrope has been modified.
154 * Some errors in the determination of the surface were corrected.
155 *
156 * Revision 1.3 2002/01/11 14:09:34 j_novak
157 * Added newtonian version for 2-fluid stars
158 *
159 * Revision 1.2 2001/12/06 15:11:43 jl_zdunik
160 * Introduction of the new function f_eq() in the class Etoile_rot
161 *
162 * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
163 * LORENE
164 *
165 * Revision 2.61 2001/10/25 08:32:57 eric
166 * Etoile_rot::display_poly passe de protected a public.
167 *
168 * Revision 2.60 2001/10/24 15:35:55 eric
169 * Classe Etoile_rot: ajout de la fonction display_poly.
170 *
171 * Revision 2.59 2001/10/16 14:48:00 eric
172 * La fonction Etoile_rot::get_omega() devient
173 * virtual Etoile_rot::get_omega_c()
174 * (retourne la valeur centrale de Omega).
175 *
176 * Revision 2.58 2001/10/11 09:24:00 eric
177 * La fonction Etoile_rot::equilibrium est desormais virtuelle.
178 *
179 * Revision 2.57 2001/08/06 15:39:04 keisuke
180 * Addition of a new argument to Etoile_bin::equilibrium and equil_regular.
181 *
182 * Revision 2.56 2001/06/25 12:52:33 eric
183 * Classe Etoile_bin : ajout du membre p_companion et des fonctions
184 * associees set_companion() et get_companion().
185 *
186 * Revision 2.55 2001/06/13 14:11:55 eric
187 * Modif commentaires (mise en conformite Doc++ 3.4.7)
188 *
189 * Revision 2.54 2001/03/26 09:29:59 jlz
190 * Classe Etoile_rot: new members p_espec_isco and p_lspec_isco.
191 *
192 * Revision 2.53 2001/02/08 15:37:09 eric
193 * *** empty log message ***
194 *
195 * Revision 2.52 2001/02/08 15:12:42 eric
196 * Ajout de la fonction Etoile_rot::f_eccentric.
197 *
198 * Revision 2.51 2000/11/23 15:43:24 eric
199 * Ajout de l'argument ent_limit a Etoile_rot::equilibrium.
200 *
201 * Revision 2.50 2000/11/19 18:51:13 eric
202 * Etoile_rot: ajout de la fonction (static) lambda_grv2
203 *
204 * Revision 2.49 2000/11/18 23:17:32 eric
205 * Ajout de l'argument ost a Etoile_rot::r_isco.
206 *
207 * Revision 2.48 2000/11/18 21:08:33 eric
208 * Classe Etoile_rot: ajout des fonctions r_isco() et f_isco()
209 * ainsi que des membres associes p_r_isco et p_f_isco.
210 *
211 * Revision 2.47 2000/11/10 15:15:38 eric
212 * Modif arguments Etoile_rot::equilibrium.
213 *
214 * Revision 2.46 2000/10/20 13:10:23 eric
215 * Etoile_rot::equilibrium: ajout de l'argument nzadapt.
216 *
217 * Revision 2.45 2000/10/17 15:59:14 eric
218 * Modif commentaires Etoile_rot::equilibrium
219 *
220 * Revision 2.44 2000/10/12 15:22:29 eric
221 * Modif commentaires Etoile_rot.
222 *
223 * Revision 2.43 2000/10/12 10:20:12 eric
224 * Ajout de la fonction Etoile_rot::fait_nphi().
225 *
226 * Revision 2.42 2000/09/22 15:50:13 keisuke
227 * d_logn_auto_div devient desormais un membre de la classe Etoile
228 * et non plus de la classe derivee Etoile_bin.
229 *
230 * Revision 2.41 2000/09/18 16:14:37 eric
231 * Classe Etoile_rot: ajout du membre tkij et de la fonction
232 * extrinsic curvature().
233 *
234 * Revision 2.40 2000/09/07 14:31:09 keisuke
235 * Ajout des membres logn_auto_regu et get_logn_auto_regu() a la classe Etoile.
236 * Ajout des membres d_logn_auto_regu et get_d_logn_auto_regu()
237 * a la classe Etoile_bin.
238 *
239 * Revision 2.39 2000/09/07 10:25:50 keisuke
240 * Ajout du membre get_logn_auto_div() a la classe Etoile.
241 * Ajout du membre get_d_logn_auto_div() a la classe Etoile_bin.
242 *
243 * Revision 2.38 2000/08/31 11:25:24 eric
244 * Classe Etoile_rot: ajout des membres tnphi et ak_car.
245 *
246 * Revision 2.37 2000/08/29 11:37:06 eric
247 * Ajout des membres k_div et logn_auto_div a la classe Etoile.
248 * Ajout du membre d_logn_auto_div a la classe Etoile_bin.
249 *
250 * Revision 2.36 2000/08/25 10:25:57 keisuke
251 * Ajout de Etoile_bin::equil_regular.
252 *
253 * Revision 2.35 2000/08/18 14:01:07 eric
254 * Modif commentaires.
255 *
256 * Revision 2.34 2000/08/17 12:38:30 eric
257 * Modifs classe Etoile_rot : ajout des membres nuf, nuq, ssjm1_nuf et ssjm1_nuq
258 * Modif arguments Etoile_rot::equilibrium.
259 * .\
260 *
261 * Revision 2.33 2000/08/07 12:11:13 keisuke
262 * Ajout de Etoile::equil_spher_regular.
263 *
264 * Revision 2.32 2000/07/21 12:02:55 eric
265 * Suppression de Etoile_rot::relaxation.
266 *
267 * Revision 2.31 2000/07/20 15:32:28 eric
268 * *** empty log message ***
269 *
270 * Revision 2.30 2000/07/06 09:39:12 eric
271 * Ajout du membre p_xa_barycenter a Etoile_bin, ainsi que de la
272 * fonction associee xa_barycenter().
273 *
274 * Revision 2.29 2000/05/25 13:47:38 eric
275 * Modif Etoile_bin::equilibrium: ajout de l'argument thres_adapt
276 *
277 * Revision 2.28 2000/03/22 16:41:45 eric
278 * Ajout de la sortie de nouvelles erreurs dans Etoile_bin::equilibrium.
279 *
280 * Revision 2.27 2000/03/22 12:54:44 eric
281 * Nouveau prototype d'Etoile_bin::velocity_potential : l'erreur est
282 * retournee en double.
283 * Nouveau prototype d'Etoile_bin::equilibrium : diff_ent est remplace
284 * par le Tbl diff.
285 *
286 * Revision 2.26 2000/03/15 11:04:15 eric
287 * Ajout des fonctions Etoile_bin::set_w_shift() et Etoile_bin::set_khi_shift()
288 * Amelioration des commentaires.
289 *
290 * Revision 2.25 2000/03/08 12:12:49 eric
291 * Ajout de la fonction Etoile_bin::is_irrotational().
292 *
293 * Revision 2.24 2000/03/07 14:48:01 eric
294 * Ajout de la fonction Etoile_bin::extrinsic_curvature().
295 *
296 * Revision 2.23 2000/02/21 13:57:57 eric
297 * classe Etoile_bin: suppression du membre psi: remplacement par psi0.
298 *
299 * Revision 2.22 2000/02/17 16:51:22 eric
300 * Ajout de l'argument diff_ent dans Etoile_bin::equilibrium.
301 *
302 * Revision 2.21 2000/02/17 15:29:38 eric
303 * Ajout de la fonction Etoile_bin::relaxation.
304 *
305 * Revision 2.20 2000/02/17 13:54:21 eric
306 * Ajout de la fonction Etoile_bin::velocity_potential.
307 *
308 * Revision 2.19 2000/02/16 15:05:14 eric
309 * Ajout des membres w_shift et khi_shift.
310 * (sauves dans les fichiers a la place de shift_auto).
311 * Ajout de la fonction Etoile_bin::fait_shift_auto.
312 *
313 * Revision 2.18 2000/02/16 13:47:02 eric
314 * Classe Etoile_bin: ajout des membres ssjm1_khi et ssjm1_wshift.
315 *
316 * Revision 2.17 2000/02/16 11:54:13 eric
317 * Classe Etoile_bin : ajout des membres ssjm1_logn et ssjm1_beta.
318 *
319 * Revision 2.16 2000/02/15 15:40:07 eric
320 * Ajout de Etoile_bin::equilibrium.
321 *
322 * Revision 2.15 2000/02/12 18:40:15 eric
323 * Modif commentaires.
324 *
325 * Revision 2.14 2000/02/12 14:44:26 eric
326 * Ajout des fonctions Etoile_bin::set_logn_comp et set_pot_centri.
327 *
328 * Revision 2.13 2000/02/10 20:22:25 eric
329 * Modif commentaires.
330 *
331 * Revision 2.12 2000/02/10 16:11:24 eric
332 * Classe Etoile_bin : ajout des accesseurs get_psi, etc...
333 * ajout de la fonction fait_d_psi
334 *
335 * Revision 2.11 2000/02/08 19:28:29 eric
336 * La fonction Etoile_bin::scal_prod est rebaptisee Etoile_bin::sprod
337 *
338 * Revision 2.10 2000/02/04 17:15:15 eric
339 * Classe Etoile_bin: ajout du membre ref_triad.
340 *
341 * Revision 2.9 2000/02/04 16:36:48 eric
342 * Ajout des fonctions update_metric* et kinematics.
343 *
344 * Revision 2.8 2000/02/02 10:12:37 eric
345 * Ajout des fonctions de lecture/ecriture mp, nzet, eos, etc...
346 *
347 * Revision 2.7 2000/02/01 15:59:43 eric
348 * Ajout de la fonction Etoile_bin::scal_prod.
349 *
350 * Revision 2.6 2000/01/31 15:56:45 eric
351 * Introduction de la classe derivee Etoile_bin.
352 *
353 * Revision 2.5 2000/01/28 17:17:45 eric
354 * Ajout des fonctions de calcul des quantites globales.
355 *
356 * Revision 2.4 2000/01/27 16:46:59 eric
357 * Ajout des fonctions get_ent(), etc....
358 *
359 * Revision 2.3 2000/01/24 17:19:48 eric
360 * Modif commentaires.
361 *
362 * Revision 2.2 2000/01/24 17:13:04 eric
363 * Le mapping mp n'est plus constant.
364 * Ajout de la fonction equilibrium_spher.
365 *
366 * Revision 2.1 2000/01/24 13:37:19 eric
367 * *** empty log message ***
368 *
369 * Revision 2.0 2000/01/20 17:04:33 eric
370 * *** empty log message ***
371 *
372 *
373 * $Header: /cvsroot/Lorene/C++/Include/etoile.h,v 1.36 2023/07/04 08:59:53 j_novak Exp $
374 *
375 */
376
377// Headers Lorene
378#include "tenseur.h"
379namespace Lorene {
380 class Eos ;
381 class Bhole ;
382
383
384 //---------------------------//
385 // base class Etoile //
386 //---------------------------//
387
427 class Etoile {
428
429 // Data :
430 // -----
431 protected:
432 Map& mp ;
433
435 int nzet ;
436
441
445 double unsurc2 ;
446
452 int k_div ;
453
454 const Eos& eos ;
455
456 // Fluid quantities with respect to the fluid frame
457 // ------------------------------------------------
458
461
465
466 // Fluid quantities with respect to the Eulerian frame
467 // ---------------------------------------------------
469
472
475
478
479 // Metric potentials
480 // -----------------
481
488
495
501
505
510
513
516
519
520 // Derived data :
521 // ------------
522 protected:
524 mutable double* p_ray_eq ;
525
527 mutable double* p_ray_eq_pis2 ;
528
530 mutable double* p_ray_eq_pi ;
531
533 mutable double* p_ray_eq_3pis2 ;
534
536 mutable double* p_ray_pole ;
537
542 mutable Itbl* p_l_surf ;
543
548 mutable Tbl* p_xi_surf ;
549
550 mutable double* p_mass_b ;
551 mutable double* p_mass_g ;
552
553
554 // Constructors - Destructor
555 // -------------------------
556 public:
557
567 Etoile(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
568
569 Etoile(const Etoile& ) ;
570
578 Etoile(Map& mp_i, const Eos& eos_i, FILE* fich) ;
579
580 virtual ~Etoile() ;
581
582 // Memory management
583 // -----------------
584 protected:
586 virtual void del_deriv() const ;
587
589 virtual void set_der_0x0() const ;
590
594 virtual void del_hydro_euler() ;
595
596
597 // Mutators / assignment
598 // ---------------------
599 public:
601 void operator=(const Etoile&) ;
602
604 Map& set_mp() {return mp; } ;
605
607 void set_enthalpy(const Cmp& ) ;
608
612 virtual void equation_of_state() ;
613
618 virtual void hydro_euler() ;
619
632 virtual void equilibrium_spher(double ent_c, double precis = 1.e-14,
633 const Tbl* ent_limit = 0x0 ) ;
634
645 void equil_spher_regular(double ent_c, double precis = 1.e-14) ;
646
655 virtual void equil_spher_falloff(double ent_c,
656 double precis = 1.e-14) ;
657
658 // Accessors
659 // ---------
660 public:
662 const Map& get_mp() const {return mp; } ;
663
665 int get_nzet() const {return nzet; } ;
666
670 bool is_relativistic() const {return relativistic; } ;
671
673 const Eos& get_eos() const {return eos; } ;
674
676 const Tenseur& get_ent() const {return ent;} ;
677
679 const Tenseur& get_nbar() const {return nbar;} ;
680
682 const Tenseur& get_ener() const {return ener;} ;
683
685 const Tenseur& get_press() const {return press;} ;
686
688 const Tenseur& get_ener_euler() const {return ener_euler;} ;
689
691 const Tenseur& get_s_euler() const {return s_euler;} ;
692
694 const Tenseur& get_gam_euler() const {return gam_euler;} ;
695
697 const Tenseur& get_u_euler() const {return u_euler;} ;
698
704 const Tenseur& get_logn_auto() const {return logn_auto;} ;
705
711 const Tenseur& get_logn_auto_regu() const {return logn_auto_regu;} ;
712
718 const Tenseur& get_logn_auto_div() const {return logn_auto_div;} ;
719
722 const Tenseur& get_d_logn_auto_div() const {return d_logn_auto_div;} ;
723
727 const Tenseur& get_beta_auto() const {return beta_auto;} ;
728
730 const Tenseur& get_nnn() const {return nnn;} ;
731
733 const Tenseur& get_shift() const {return shift;} ;
734
736 const Tenseur& get_a_car() const {return a_car;} ;
737
738 // Outputs
739 // -------
740 public:
741 virtual void sauve(FILE* ) const ;
742
744 friend ostream& operator<<(ostream& , const Etoile& ) ;
745
746 protected:
748 virtual ostream& operator>>(ostream& ) const ;
749
750 // Global quantities
751 // -----------------
752 public:
754 double ray_eq() const ;
755
757 double ray_eq_pis2() const ;
758
760 double ray_eq_pi() const ;
761
763 double ray_eq_3pis2() const ;
764
766 double ray_pole() const ;
767
769 double ray_eq(int kk) const ;
770
778 virtual const Itbl& l_surf() const ;
779
787 const Tbl& xi_surf() const ;
788
790 virtual double mass_b() const ;
791
793 virtual double mass_g() const ;
794
795 };
796 ostream& operator<<(ostream& , const Etoile& ) ;
797
798
799 //---------------------------//
800 // class Etoile_bin //
801 //---------------------------//
802
817 class Etoile_bin : public Etoile {
818
819 // Data :
820 // -----
821 protected:
826
832
837
842
848
853
858
863
868
873
878
883
888
893
899
912
922
929
936
942
948
954
957
963
969
977
987
993
1004
1005 // Derived data :
1006 // ------------
1007 protected:
1009 mutable double* p_xa_barycenter ;
1010
1011
1012 // Constructors - Destructor
1013 // -------------------------
1014 public:
1030 Etoile_bin(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i,
1031 bool irrot, const Base_vect& ref_triad_i) ;
1032
1033
1034 Etoile_bin(const Etoile_bin& ) ;
1035
1048 Etoile_bin(Map& mp_i, const Eos& eos_i, const Base_vect& ref_triad_i,
1049 FILE* fich) ;
1050
1051 virtual ~Etoile_bin() ;
1052
1053
1054 // Memory management
1055 // -----------------
1056 protected:
1058 virtual void del_deriv() const ;
1059
1061 virtual void set_der_0x0() const ;
1062
1066 virtual void del_hydro_euler() ;
1067
1068
1069 // Mutators / assignment
1070 // ---------------------
1071 public:
1073 void operator=(const Etoile_bin& ) ;
1074
1079
1082
1084 Tenseur& set_w_shift() ;
1085
1088
1089 // Accessors
1090 // ---------
1091 public:
1095 bool is_irrotational() const {return irrotational; } ;
1096
1098 const Tenseur& get_psi0() const {return psi0;} ;
1099
1103 const Tenseur& get_d_psi() const {return d_psi;} ;
1104
1109 const Tenseur& get_wit_w() const {return wit_w;} ;
1110
1114 const Tenseur& get_loggam() const {return loggam;} ;
1115
1119 const Tenseur& get_logn_comp() const {return logn_comp;} ;
1120
1124 const Tenseur& get_d_logn_auto() const {return d_logn_auto;} ;
1125
1130
1134 const Tenseur& get_d_logn_comp() const {return d_logn_comp;} ;
1135
1139 const Tenseur& get_beta_comp() const {return beta_comp;} ;
1140
1144 const Tenseur& get_d_beta_auto() const {return d_beta_auto;} ;
1145
1149 const Tenseur& get_d_beta_comp() const {return d_beta_comp;} ;
1150
1155 const Tenseur& get_shift_auto() const {return shift_auto;} ;
1156
1161 const Tenseur& get_shift_comp() const {return shift_comp;} ;
1162
1175 const Tenseur& get_w_shift() const {return w_shift;} ;
1176
1189 const Tenseur& get_khi_shift() const {return khi_shift;} ;
1190
1195 const Tenseur_sym& get_tkij_auto() const {return tkij_auto;} ;
1196
1201 const Tenseur_sym& get_tkij_comp() const {return tkij_comp;} ;
1202
1207 const Tenseur& get_akcar_auto() const {return akcar_auto;} ;
1208
1213 const Tenseur& get_akcar_comp() const {return akcar_comp;} ;
1214
1219 const Tenseur& get_bsn() const {return bsn;} ;
1220
1222 const Tenseur& get_pot_centri() const {return pot_centri;} ;
1227 const Cmp get_decouple() const {return decouple ;}
1228
1229 // Outputs
1230 // -------
1231 public:
1232 virtual void sauve(FILE* ) const ;
1233
1234 protected:
1236 virtual ostream& operator>>(ostream& ) const ;
1237
1238 // Global quantities
1239 // -----------------
1240 public:
1242 virtual double mass_b() const ;
1243
1245 virtual double mass_g() const ;
1246
1255 virtual double xa_barycenter() const ;
1256
1257
1258 // Computational routines
1259 // ----------------------
1260 public:
1268 virtual Tenseur sprod(const Tenseur& t1, const Tenseur& t2) const ;
1269
1282 virtual void hydro_euler() ;
1283
1301 void update_metric(const Etoile_bin& comp) ;
1302
1310 void update_metric(const Etoile_bin& comp, const Etoile_bin& star_prev,
1311 double relax) ;
1312
1327 void update_metric_der_comp(const Etoile_bin& comp) ;
1328
1339 virtual void kinematics(double omega, double x_axe) ;
1340
1344 void fait_d_psi() ;
1345
1354 void fait_shift_auto() ;
1355
1359 virtual void extrinsic_curvature() ;
1360
1401 void equilibrium(double ent_c,
1402 int mermax, int mermax_poisson,
1403 double relax_poisson, int mermax_potvit,
1404 double relax_potvit, double thres_adapt,
1405 const Tbl& fact, Tbl& diff, const Tbl* ent_limit = 0x0) ;
1406
1446 void equil_regular(double ent_c, int mermax, int mermax_poisson,
1447 double relax_poisson, int mermax_potvit,
1448 double relax_potvit, double thres_adapt,
1449 const Tbl& fact, Tbl& diff) ;
1450
1464 double velocity_potential(int mermax, double precis, double relax) ;
1465
1476 void relaxation(const Etoile_bin& star_prev, double relax_ent,
1477 double relax_met, int mer, int fmer_met) ;
1478
1479 friend class Bin_ns_bh ;
1480 };
1481
1482
1483
1484 //---------------------------//
1485 // class Etoile_rot //
1486 //---------------------------//
1487
1499class Etoile_rot : public Etoile {
1500
1501 // Data :
1502 // -----
1503 protected:
1504 double omega ;
1505
1508
1511
1514
1519
1522
1525
1530
1535
1538
1541
1554
1564
1571
1590
1596
1602
1607
1612
1620
1629
1630 // Derived data :
1631 // ------------
1632 protected:
1633
1634 mutable double* p_angu_mom ;
1635 mutable double* p_tsw ;
1636 mutable double* p_grv2 ;
1637 mutable double* p_grv3 ;
1638 mutable double* p_r_circ ;
1639 mutable double* p_r_circ_merid ;
1640 mutable double* p_area ;
1641 mutable double* p_aplat ;
1642 mutable double* p_z_eqf ;
1643 mutable double* p_z_eqb ;
1644 mutable double* p_z_pole ;
1645 mutable double* p_mom_quad ;
1646 mutable double* p_mom_quad_old ;
1647 mutable double* p_mom_quad_Bo ;
1648 mutable double* p_r_isco ;
1649 mutable double* p_f_isco ;
1651 mutable double* p_espec_isco ;
1653 mutable double* p_lspec_isco ;
1654 mutable double* p_f_eq ;
1655
1656
1657
1658 // Constructors - Destructor
1659 // -------------------------
1660 public:
1669 Etoile_rot(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
1670
1671
1672 Etoile_rot(const Etoile_rot& ) ;
1673
1681 Etoile_rot(Map& mp_i, const Eos& eos_i, FILE* fich) ;
1682
1683 virtual ~Etoile_rot() ;
1684
1685
1686 // Memory management
1687 // -----------------
1688 protected:
1690 virtual void del_deriv() const ;
1691
1693 virtual void set_der_0x0() const ;
1694
1698 virtual void del_hydro_euler() ;
1699
1700
1701 // Mutators / assignment
1702 // ---------------------
1703 public:
1705 void operator=(const Etoile_rot& ) ;
1706
1707 // Accessors
1708 // ---------
1709 public:
1713 virtual double get_omega_c() const ;
1714
1716 const Tenseur& get_bbb() const {return bbb;} ;
1717
1719 const Tenseur& get_b_car() const {return b_car;} ;
1720
1722 const Tenseur& get_nphi() const {return nphi;} ;
1723
1727 const Tenseur& get_tnphi() const {return tnphi;} ;
1728
1730 const Tenseur& get_uuu() const {return uuu;} ;
1731
1733 const Tenseur& get_logn() const {return logn;} ;
1734
1738 const Tenseur& get_nuf() const {return nuf;} ;
1739
1743 const Tenseur& get_nuq() const {return nuq;} ;
1744
1746 const Tenseur& get_dzeta() const {return dzeta;} ;
1747
1749 const Tenseur& get_tggg() const {return tggg;} ;
1750
1763 const Tenseur& get_w_shift() const {return w_shift;} ;
1764
1777 const Tenseur& get_khi_shift() const {return khi_shift;} ;
1778
1784 const Tenseur_sym& get_tkij() const {return tkij;} ;
1785
1803 const Tenseur& get_ak_car() const {return ak_car;} ;
1804
1805 // Outputs
1806 // -------
1807 public:
1808 virtual void sauve(FILE* ) const ;
1809
1811 virtual void display_poly(ostream& ) const ;
1812
1813 protected:
1815 virtual ostream& operator>>(ostream& ) const ;
1816
1818 virtual void partial_display(ostream& ) const ;
1819
1820 // Global quantities
1821 // -----------------
1822 public:
1823
1831 virtual const Itbl& l_surf() const ;
1832
1833 virtual double mass_b() const ;
1834 virtual double mass_g() const ;
1835 virtual double angu_mom() const ;
1836 virtual double tsw() const ;
1837
1841 virtual double grv2() const ;
1842
1854 virtual double grv3(ostream* ost = 0x0) const ;
1855
1856 virtual double r_circ() const ;
1857 virtual double r_circ_merid() const ;
1858 virtual double area() const ;
1859 virtual double mean_radius() const ;
1860 virtual double aplat() const ;
1861 virtual double z_eqf() const ;
1862 virtual double z_eqb() const ;
1863 virtual double z_pole() const ;
1864
1879 virtual double mom_quad() const ;
1880
1888 virtual double mom_quad_old() const ;
1889
1896 virtual double mom_quad_Bo() const ;
1897
1904 virtual double r_isco(ostream* ost = 0x0) const ;
1905
1907 virtual double f_isco() const ;
1908
1910 virtual double espec_isco() const ;
1911
1913 virtual double lspec_isco() const ;
1914
1915
1926 virtual double f_eccentric(double ecc, double periast,
1927 ostream* ost = 0x0) const ;
1928
1930 virtual double f_eq() const ;
1931
1932
1933 // Computational routines
1934 // ----------------------
1935 public:
1946 virtual void hydro_euler() ;
1947
1957 void update_metric() ;
1958
1967 void fait_shift() ;
1968
1972 void fait_nphi() ;
1973
1977 void extrinsic_curvature() ;
1978
2008 static double lambda_grv2(const Cmp& sou_m, const Cmp& sou_q) ;
2009
2089 virtual void equilibrium(double ent_c, double omega0, double fact_omega,
2090 int nzadapt, const Tbl& ent_limit,
2091 const Itbl& icontrol, const Tbl& control,
2092 double mbar_wanted, double aexp_mass,
2093 Tbl& diff, Param* = 0x0) ;
2094
2095
2096};
2097
2098
2099
2100
2101}
2102#endif
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition base_vect.h:105
Black hole.
Definition bhole.h:268
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Equation of state base class.
Definition eos.h:206
Class for stars in binary system.
Definition etoile.h:817
Tenseur shift_comp
Part of the shift vector generated principaly by the companion star.
Definition etoile.h:898
void update_metric(const Etoile_bin &comp)
Computes metric coefficients from known potentials, when the companion is another star.
const Tenseur & get_psi0() const
Returns the non-translational part of the velocity potential.
Definition etoile.h:1098
Tenseur psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case).
Definition etoile.h:836
const Tenseur & get_d_logn_auto_regu() const
Returns the gradient of logn_auto_regu (Cartesian components with respect to ref_triad ).
Definition etoile.h:1129
Tenseur d_beta_auto
Gradient of beta_auto (Cartesian components with respect to ref_triad ).
Definition etoile.h:882
Tenseur d_logn_auto
Gradient of logn_auto (Cartesian components with respect to ref_triad ).
Definition etoile.h:862
virtual double mass_b() const
Baryon mass.
double * p_xa_barycenter
Absolute coordinate X of the barycenter of the baryon density.
Definition etoile.h:1009
const Tenseur & get_d_logn_auto() const
Returns the gradient of logn_auto (Cartesian components with respect to ref_triad ).
Definition etoile.h:1124
virtual void sauve(FILE *) const
Save in a file.
Definition etoile_bin.C:562
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition etoile.h:953
Tenseur akcar_comp
Part of the scalar generated by shift_auto and shift_comp , i.e.
Definition etoile.h:947
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Definition etoile_bin.C:470
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:976
const Tenseur & get_wit_w() const
Returns the spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Definition etoile.h:1109
Tenseur d_psi
Gradient of (in the irrotational case) (Cartesian components with respect to ref_triad ).
Definition etoile.h:841
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
void equilibrium(double ent_c, int mermax, int mermax_poisson, double relax_poisson, int mermax_potvit, double relax_potvit, double thres_adapt, const Tbl &fact, Tbl &diff, const Tbl *ent_limit=0x0)
Computes an equilibrium configuration.
void relaxation(const Etoile_bin &star_prev, double relax_ent, double relax_met, int mer, int fmer_met)
Performs a relaxation on ent , logn_auto , beta_auto and shift_auto .
Definition etoile_bin.C:866
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur 's...
Definition etoile.h:831
void fait_shift_auto()
Computes shift_auto from w_shift and khi_shift according to Shibata's prescription [Prog.
Definition etoile_bin.C:829
void update_metric_der_comp(const Etoile_bin &comp)
Computes the derivative of metric functions related to the companion star.
const Tenseur & get_logn_comp() const
Returns the part of the lapse logarithm (gravitational potential at the Newtonian limit) generated pr...
Definition etoile.h:1119
bool irrotational
true for an irrotational star, false for a corotating one
Definition etoile.h:825
const Tenseur & get_shift_auto() const
Returns the part of the shift vector generated principaly by the star.
Definition etoile.h:1155
Tenseur w_shift
Vector used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:911
Tenseur logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principaly by ...
Definition etoile.h:857
const Cmp get_decouple() const
Returns the function used to construct tkij_auto from tkij_tot .
Definition etoile.h:1227
Tenseur & set_w_shift()
Read/write of w_shift.
Definition etoile_bin.C:541
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile_bin.C:462
const Tenseur_sym & get_tkij_auto() const
Returns the part of the extrinsic curvature tensor generated by shift_auto .
Definition etoile.h:1195
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad ).
Definition etoile.h:872
Tenseur beta_comp
Part of the logarithm of AN generated principaly by the companion star.
Definition etoile.h:877
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
const Tenseur & get_w_shift() const
Returns the vector used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:1175
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Definition etoile.h:986
Cmp ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto by mea...
Definition etoile.h:962
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile_bin.C:450
const Tenseur & get_shift_comp() const
Returns the part of the shift vector generated principaly by the companion star.
Definition etoile.h:1161
Cmp ssjm1_psi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:992
const Tenseur & get_loggam() const
Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition etoile.h:1114
Tenseur pot_centri
Centrifugal potential.
Definition etoile.h:956
void equil_regular(double ent_c, int mermax, int mermax_poisson, double relax_poisson, int mermax_potvit, double relax_potvit, double thres_adapt, const Tbl &fact, Tbl &diff)
Computes an equilibrium configuration by regularizing the diverging source.
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition etoile.h:852
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition etoile.h:892
virtual double xa_barycenter() const
Absolute coordinate X of the barycenter of the baryon density, defined according to the formula.
const Tenseur & get_d_beta_auto() const
Returns the gradient of beta_auto (Cartesian components with respect to ref_triad ).
Definition etoile.h:1144
bool is_irrotational() const
Returns true for an irrotational motion, false for a corotating one.
Definition etoile.h:1095
const Tenseur & get_akcar_comp() const
Returns the part of the scalar generated by shift_auto and shift_comp , i.e.
Definition etoile.h:1213
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor generated by shift_auto .
Definition etoile.h:928
Tenseur akcar_auto
Part of the scalar generated by shift_auto , i.e.
Definition etoile.h:941
Cmp ssjm1_beta
Effective source at the previous step for the resolution of the Poisson equation for beta_auto by mea...
Definition etoile.h:968
Cmp decouple
Function used to construct the part of generated by the star from the total .
Definition etoile.h:1003
Tenseur wit_w
Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Definition etoile.h:847
const Tenseur & get_akcar_auto() const
Returns the part of the scalar generated by shift_auto , i.e.
Definition etoile.h:1207
virtual double mass_g() const
Gravitational mass.
Tenseur & set_logn_comp()
Read/write the part of the lapse logarithm (gravitational potential at the Newtonian limit) generated...
Definition etoile_bin.C:527
Tenseur_sym tkij_comp
Part of the extrinsic curvature tensor generated by shift_comp .
Definition etoile.h:935
Tenseur & set_khi_shift()
Read/write of khi_shift.
Definition etoile_bin.C:548
virtual Tenseur sprod(const Tenseur &t1, const Tenseur &t2) const
Performs the scalar product of two tensors by contracting the last index of t1 with the first index o...
Definition etoile_bin.C:751
const Tenseur & get_pot_centri() const
Returns the centrifugal potential.
Definition etoile.h:1222
virtual void extrinsic_curvature()
Computes tkij_auto and akcar_auto from shift_auto , nnn and a_car .
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad ).
Definition etoile.h:887
Etoile_bin(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, bool irrot, const Base_vect &ref_triad_i)
Standard constructor.
Definition etoile_bin.C:210
Tenseur d_logn_auto_regu
Gradient of logn_auto_regu (Cartesian components with respect to ref_triad ).
Definition etoile.h:867
virtual ~Etoile_bin()
Destructor.
Definition etoile_bin.C:440
const Tenseur & get_bsn() const
Returns the shift vector, divided by N , of the rotating coordinates, .
Definition etoile.h:1219
Tenseur & set_pot_centri()
Read/write the centrifugal potential.
Definition etoile_bin.C:534
void fait_d_psi()
Computes the gradient of the total velocity potential .
Definition etoile_bin.C:767
const Tenseur & get_d_logn_comp() const
Returns the gradient of logn_comp (Cartesian components with respect to ref_triad ).
Definition etoile.h:1134
const Tenseur & get_d_beta_comp() const
Returns the gradient of beta_comp (Cartesian components with respect to ref_triad ).
Definition etoile.h:1149
const Tenseur & get_beta_comp() const
Returns the part of the logarithm of AN generated principaly by the companion star.
Definition etoile.h:1139
virtual void kinematics(double omega, double x_axe)
Computes the quantities bsn and pot_centri .
const Tenseur_sym & get_tkij_comp() const
Returns the part of the extrinsic curvature tensor generated by shift_comp .
Definition etoile.h:1201
void operator=(const Etoile_bin &)
Assignment to another Etoile_bin.
Definition etoile_bin.C:485
const Tenseur & get_khi_shift() const
Returns the scalar used in the decomposition of shift_auto following Shibata's prescription [Prog.
Definition etoile.h:1189
const Tenseur & get_d_psi() const
Returns the gradient of the velocity potential (Cartesian components with respect to ref_triad ).
Definition etoile.h:1103
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition etoile_bin.C:585
friend class Bin_ns_bh
Friend class Bin_ns_bh.
Definition etoile.h:1479
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:921
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Definition etoile.h:1628
double * p_r_circ_merid
Circumferential meridional radius.
Definition etoile.h:1639
const Tenseur & get_b_car() const
Returns the square of the metric factor B.
Definition etoile.h:1719
Tenseur uuu
Norm of u_euler.
Definition etoile.h:1521
const Tenseur & get_tnphi() const
Returns the component of the shift vector.
Definition etoile.h:1727
const Tenseur & get_w_shift() const
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1763
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile_rot.C:344
double omega
Rotation angular velocity ([f_unit] ).
Definition etoile.h:1504
double * p_z_pole
Redshift factor at North pole.
Definition etoile.h:1644
virtual double espec_isco() const
Energy of a particle on the ISCO.
virtual double mom_quad() const
Quadrupole moment.
double * p_z_eqf
Forward redshift factor at equator.
Definition etoile.h:1642
void operator=(const Etoile_rot &)
Assignment to another Etoile_rot.
Definition etoile_rot.C:415
virtual double r_circ() const
Circumferential equatorial radius.
virtual double f_eccentric(double ecc, double periast, ostream *ost=0x0) const
Computation of frequency of eccentric orbits.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile_rot.C:374
Tenseur & logn
Metric potential = logn_auto.
Definition etoile.h:1524
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
const Tenseur & get_nphi() const
Returns the metric coefficient .
Definition etoile.h:1722
virtual const Itbl & l_surf() const
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on...
const Tenseur & get_khi_shift() const
Returns the scalar used in the decomposition of shift following Shibata's prescription [Prog.
Definition etoile.h:1777
double * p_mom_quad_old
Part of the quadrupole moment.
Definition etoile.h:1646
const Tenseur & get_tggg() const
Returns the Metric potential .
Definition etoile.h:1749
Cmp ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta .
Definition etoile.h:1606
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition etoile.h:1534
const Tenseur_sym & get_tkij() const
Returns the tensor related to the extrinsic curvature tensor by .
Definition etoile.h:1784
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition etoile_rot.C:477
double * p_z_eqb
Backward redshift factor at equator.
Definition etoile.h:1643
virtual double mass_g() const
Gravitational mass.
virtual double f_isco() const
Orbital frequency at the innermost stable circular orbit (ISCO).
double * p_r_circ
Circumferential equatorial radius.
Definition etoile.h:1638
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1563
virtual double mom_quad_Bo() const
Part of the quadrupole moment.
virtual double aplat() const
Flatening r_pole/r_eq.
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Etoile_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition etoile_rot.C:146
const Tenseur & get_bbb() const
Returns the metric factor B.
Definition etoile.h:1716
virtual double tsw() const
Ratio T/W.
Tenseur tggg
Metric potential .
Definition etoile.h:1540
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Definition etoile.h:1529
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Definition etoile.h:1611
Tenseur nphi
Metric coefficient .
Definition etoile.h:1513
virtual double mass_b() const
Baryon mass.
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
Definition et_rot_isco.C:87
static double lambda_grv2(const Cmp &sou_m, const Cmp &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
const Tenseur & get_uuu() const
Returns the norm of u_euler.
Definition etoile.h:1730
const Tenseur & get_logn() const
Returns the metric potential = logn_auto.
Definition etoile.h:1733
Tenseur bbb
Metric factor B.
Definition etoile.h:1507
double * p_area
Surface area.
Definition etoile.h:1640
void update_metric()
Computes metric coefficients from known potentials.
void extrinsic_curvature()
Computes tkij and ak_car from shift , nnn and b_car .
Tenseur ak_car
Scalar .
Definition etoile.h:1589
virtual ~Etoile_rot()
Destructor.
Definition etoile_rot.C:334
Tenseur & dzeta
Metric potential = beta_auto.
Definition etoile.h:1537
virtual double mean_radius() const
Mean radius.
double * p_grv3
Error on the virial identity GRV3.
Definition etoile.h:1637
double * p_grv2
Error on the virial identity GRV2.
Definition etoile.h:1636
double * p_aplat
Flatening r_pole/r_eq.
Definition etoile.h:1641
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition etoile.h:1595
double * p_mom_quad_Bo
Part of the quadrupole moment.
Definition etoile.h:1647
virtual double lspec_isco() const
Angular momentum of a particle on the ISCO.
const Tenseur & get_nuf() const
Returns the part of the Metric potential = logn generated by the matter terms.
Definition etoile.h:1738
virtual double grv2() const
Error on the virial identity GRV2.
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Definition etoile_rot.C:803
virtual double r_circ_merid() const
Circumferential meridional radius.
double * p_mom_quad
Quadrupole moment.
Definition etoile.h:1645
double * p_f_eq
Orbital frequency at the equator.
Definition etoile.h:1654
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Definition etoile_rot.C:400
virtual double area() const
Surface area.
double * p_angu_mom
Angular momentum.
Definition etoile.h:1634
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:1619
virtual void equilibrium(double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff, Param *=0x0)
Computes an equilibrium configuration.
virtual double angu_mom() const
Angular momentum.
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Definition etoile.h:1570
virtual double z_eqf() const
Forward redshift factor at equator.
const Tenseur & get_dzeta() const
Returns the Metric potential = beta_auto.
Definition etoile.h:1746
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition etoile_rot.C:708
double * p_f_isco
Orbital frequency of the ISCO.
Definition etoile.h:1649
virtual double mom_quad_old() const
Part of the quadrupole moment.
const Tenseur & get_nuq() const
Returns the Part of the Metric potential = logn generated by the quadratic terms.
Definition etoile.h:1743
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition etoile.h:1653
double * p_tsw
Ratio T/W.
Definition etoile.h:1635
Tenseur b_car
Square of the metric factor B.
Definition etoile.h:1510
virtual double z_eqb() const
Backward redshift factor at equator.
Tenseur tnphi
Component of the shift vector.
Definition etoile.h:1518
const Tenseur & get_ak_car() const
Returns the scalar .
Definition etoile.h:1803
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition etoile.h:1651
double * p_r_isco
Circumferential radius of the ISCO.
Definition etoile.h:1648
virtual void sauve(FILE *) const
Save in a file.
Definition etoile_rot.C:452
virtual double f_eq() const
Orbital frequency at the equator.
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] ).
Definition etoile_rot.C:698
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition etoile.h:1601
virtual double z_pole() const
Redshift factor at North pole.
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1553
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition etoile_rot.C:645
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog.
Definition etoile_rot.C:766
Base class for stars *** DEPRECATED : use class Star instead ***.
Definition etoile.h:427
virtual void equilibrium_spher(double ent_c, double precis=1.e-14, const Tbl *ent_limit=0x0)
Computes a spherical static configuration.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile.C:399
Tenseur logn_auto_div
Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by t...
Definition etoile.h:500
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:435
const Tenseur & get_nnn() const
Returns the total lapse function N.
Definition etoile.h:730
double * p_mass_b
Baryon mass.
Definition etoile.h:550
int get_nzet() const
Returns the number of domains occupied by the star.
Definition etoile.h:665
double * p_mass_g
Gravitational mass.
Definition etoile.h:551
Tbl * p_xi_surf
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the su...
Definition etoile.h:548
void operator=(const Etoile &)
Assignment to another Etoile.
Definition etoile.C:433
Tenseur logn_auto_regu
Regular part of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:494
double ray_eq_pi() const
Coordinate radius at , [r_unit].
const Tenseur & get_shift() const
Returns the total shift vector .
Definition etoile.h:733
int k_div
Index of regularity of the gravitational potential logn_auto .
Definition etoile.h:452
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Definition etoile.h:512
const Tenseur & get_logn_auto_regu() const
Returns the regular part of the logarithm of the part of the lapse N generated principaly by the star...
Definition etoile.h:711
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:487
const Map & get_mp() const
Returns the mapping.
Definition etoile.h:662
Tenseur nbar
Baryon density in the fluid frame.
Definition etoile.h:462
const Eos & eos
Equation of state of the stellar matter.
Definition etoile.h:454
bool is_relativistic() const
Returns true for a relativistic star, false for a Newtonian one.
Definition etoile.h:670
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition etoile.C:569
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition etoile.h:477
Itbl * p_l_surf
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surfa...
Definition etoile.h:542
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:474
friend ostream & operator<<(ostream &, const Etoile &)
Display.
Definition etoile.C:509
Map & mp
Mapping associated with the star.
Definition etoile.h:432
double * p_ray_eq_3pis2
Coordinate radius at , .
Definition etoile.h:533
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile.C:381
virtual void equil_spher_falloff(double ent_c, double precis=1.e-14)
Computes a spherical static configuration with the outer boundary condition at a finite radius.
double * p_ray_eq_pi
Coordinate radius at , .
Definition etoile.h:530
const Tenseur & get_ent() const
Returns the enthalpy field.
Definition etoile.h:676
const Eos & get_eos() const
Returns the equation of state.
Definition etoile.h:673
virtual double mass_b() const
Baryon mass.
const Tenseur & get_beta_auto() const
Returns the logarithm of the part of the product AN generated principaly by the star.
Definition etoile.h:727
const Tenseur & get_d_logn_auto_div() const
Returns the gradient of logn_auto_div.
Definition etoile.h:722
Tenseur ener
Total energy density in the fluid frame.
Definition etoile.h:463
Etoile(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition etoile.C:146
void equil_spher_regular(double ent_c, double precis=1.e-14)
Computes a spherical static configuration.
double ray_eq_pis2() const
Coordinate radius at , [r_unit].
const Tenseur & get_a_car() const
Returns the total conformal factor .
Definition etoile.h:736
virtual ~Etoile()
Destructor.
Definition etoile.C:370
double * p_ray_pole
Coordinate radius at .
Definition etoile.h:536
const Tenseur & get_logn_auto_div() const
Returns the divergent part of the logarithm of the part of the lapse N generated principaly by the st...
Definition etoile.h:718
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 ).
Definition etoile.h:504
virtual const Itbl & l_surf() const
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l ...
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition etoile.C:514
double ray_eq_3pis2() const
Coordinate radius at , [r_unit].
Map & set_mp()
Read/write of the mapping.
Definition etoile.h:604
Tenseur press
Fluid pressure.
Definition etoile.h:464
virtual double mass_g() const
Gravitational mass.
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:440
const Tenseur & get_logn_auto() const
Returns the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:704
const Tenseur & get_u_euler() const
Returns the fluid 3-velocity with respect to the Eulerian observer.
Definition etoile.h:697
double * p_ray_eq_pis2
Coordinate radius at , .
Definition etoile.h:527
const Tenseur & get_s_euler() const
Returns the trace of the stress tensor in the Eulerian frame.
Definition etoile.h:691
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition etoile.h:468
virtual void sauve(FILE *) const
Save in a file.
Definition etoile.C:486
Tenseur shift
Total shift vector.
Definition etoile.h:515
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition etoile.C:687
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition etoile.h:471
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Definition etoile.h:460
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Definition etoile.C:413
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition etoile.h:509
const Tenseur & get_ener() const
Returns the proper total energy density.
Definition etoile.h:682
const Tbl & xi_surf() const
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinat...
const Tenseur & get_press() const
Returns the fluid pressure.
Definition etoile.h:685
Tenseur a_car
Total conformal factor .
Definition etoile.h:518
void set_enthalpy(const Cmp &)
Assignment of the enthalpy field.
Definition etoile.C:468
const Tenseur & get_ener_euler() const
Returns the total energy density with respect to the Eulerian observer.
Definition etoile.h:688
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition etoile.h:445
double * p_ray_eq
Coordinate radius at , .
Definition etoile.h:524
const Tenseur & get_nbar() const
Returns the proper baryon density.
Definition etoile.h:679
const Tenseur & get_gam_euler() const
Returns the Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:694
double ray_pole() const
Coordinate radius at [r_unit].
Basic integer array class.
Definition itbl.h:122
Parameter storage.
Definition param.h:125
Basic array class.
Definition tbl.h:161
Class intended to describe tensors with a symmetry on the two last indices *** DEPRECATED : use class...
Definition tenseur.h:1256
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142