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.37 2026/02/26 09:53:04 j_novak Exp $
36  * $Log: etoile.h,v $
37  * Revision 1.37 2026/02/26 09:53:04 j_novak
38  * Added method to compute and store circumferential radius as function of latitude (radius_Morsink).
39  *
40  * Revision 1.36 2023/07/04 08:59:53 j_novak
41  * Added method "r_circ_merid()" to compute the circumferential meridional radius.
42  *
43  * Revision 1.35 2015/06/12 12:38:24 j_novak
44  * Implementation of the corrected formula for the quadrupole momentum.
45  *
46  * Revision 1.34 2015/06/11 13:50:19 j_novak
47  * Minor corrections
48  *
49  * Revision 1.33 2015/06/10 14:36:39 a_sourie
50  * Corrected the formula for the computation of the quadrupole momentum.
51  *
52  * Revision 1.32 2014/10/13 08:52:34 j_novak
53  * Lorene classes and functions now belong to the namespace Lorene.
54  *
55  * Revision 1.31 2011/10/06 14:55:36 j_novak
56  * equation_of_state() is now virtual to be able to call to the magnetized
57  * Eos_mag.
58  *
59  * Revision 1.30 2010/02/02 21:05:49 e_gourgoulhon
60  * Etoile_bin:equilibrium : suppressed the argument method_khi added by
61  * mistake during previous commit.
62  *
63  * Revision 1.29 2010/02/02 13:34:12 e_gourgoulhon
64  * Marked DEPRECATED (in the documentation).
65  *
66  * Revision 1.28 2008/11/14 13:51:08 e_gourgoulhon
67  * Added the parameter ent_limit to Etoile::equilibrium_spher and
68  * Etoile_bin::equilibrium.
69  *
70  * Revision 1.27 2005/10/05 15:14:47 j_novak
71  * Added a Param* as parameter of Etoile_rot::equilibrium
72  *
73  * Revision 1.26 2005/08/29 15:10:12 p_grandclement
74  * Addition of things needed :
75  * 1) For BBH with different masses
76  * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
77  * WORKING YET !!!
78  *
79  * Revision 1.25 2005/01/18 22:35:51 k_taniguchi
80  * Delete a pointer for Etoile::ray_eq(int kk).
81  *
82  * Revision 1.24 2005/01/18 20:34:14 k_taniguchi
83  * Addition of Etoile::ray_eq(int kk).
84  *
85  * Revision 1.23 2005/01/17 20:39:32 k_taniguchi
86  * Addition of Etoile::ray_eq_3pis2().
87  *
88  * Revision 1.22 2004/11/30 20:40:06 k_taniguchi
89  * Addition of the method for calculating a spherical star with falloff
90  * condition at the outer boundary.
91  *
92  * Revision 1.21 2004/05/10 10:18:33 f_limousin
93  * Change to avoid warning in the compilation of Lorene
94  *
95  * Revision 1.20 2004/05/07 12:37:12 f_limousin
96  * Add new member ssjm1_psi
97  *
98  * Revision 1.19 2004/04/08 16:42:31 f_limousin
99  * Add a function velocity_potential with argument ssjm1_psi for the
100  * case of strange stars.
101  *
102  * Revision 1.18 2004/03/22 13:12:41 j_novak
103  * Modification of comments to use doxygen instead of doc++
104  *
105  * Revision 1.17 2003/10/24 12:24:41 k_taniguchi
106  * Suppress the methods of update metric for NS-BH
107  *
108  * Revision 1.16 2003/10/21 11:44:43 k_taniguchi
109  * Delete various things for the Bin_ns_bh project.
110  * They are moved to et_bin_nsbh.h.
111  *
112  * Revision 1.15 2003/10/20 14:50:04 k_taniguchi
113  * Addition of various things for the Bin_ns_bh project
114  * which are related with the part of the neutron star.
115  *
116  * Revision 1.14 2003/10/20 13:11:03 k_taniguchi
117  * Back to version 1.12
118  *
119  * Revision 1.13 2003/10/20 12:15:55 k_taniguchi
120  * Addition of various things for the Bin_ns_bh project
121  * which are related with the part of the neutron star.
122  *
123  * Revision 1.12 2003/06/20 14:13:16 f_limousin
124  * Change to virtual the functions equilibrium_spher() and kinematics().
125  *
126  * Revision 1.11 2003/02/13 16:40:24 p_grandclement
127  * Addition of various things for the Bin_ns_bh project, non of them being
128  * completely tested
129  *
130  * Revision 1.10 2003/02/04 16:20:35 f_limousin
131  * Change to virtual the routine extrinsic_curvature
132  *
133  * Revision 1.9 2003/01/31 16:57:12 p_grandclement
134  * addition of the member Cmp decouple used to compute the K_ij auto, once
135  * the K_ij total is known
136  *
137  * Revision 1.8 2002/12/19 14:48:00 e_gourgoulhon
138  *
139  * Class Etoile_bin: added the new functions:
140  * void update_metric(const Bhole& comp)
141  * void update_metric_der_comp(const Bhole& comp)
142  * to treat the case where the companion is a black hole
143  *
144  * Revision 1.7 2002/12/17 21:17:08 e_gourgoulhon
145  * Class Etoile_bin: suppression of the member p_companion
146  * as well as the associated functions set_companion
147  * and get_companion.
148  *
149  * Revision 1.6 2002/09/13 09:17:33 j_novak
150  * Modif. commentaires
151  *
152  * Revision 1.5 2002/06/17 14:05:16 j_novak
153  * friend functions are now also declared outside the class definition
154  *
155  * Revision 1.4 2002/04/05 09:09:36 j_novak
156  * The inversion of the EOS for 2-fluids polytrope has been modified.
157  * Some errors in the determination of the surface were corrected.
158  *
159  * Revision 1.3 2002/01/11 14:09:34 j_novak
160  * Added newtonian version for 2-fluid stars
161  *
162  * Revision 1.2 2001/12/06 15:11:43 jl_zdunik
163  * Introduction of the new function f_eq() in the class Etoile_rot
164  *
165  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
166  * LORENE
167  *
168  * Revision 2.61 2001/10/25 08:32:57 eric
169  * Etoile_rot::display_poly passe de protected a public.
170  *
171  * Revision 2.60 2001/10/24 15:35:55 eric
172  * Classe Etoile_rot: ajout de la fonction display_poly.
173  *
174  * Revision 2.59 2001/10/16 14:48:00 eric
175  * La fonction Etoile_rot::get_omega() devient
176  * virtual Etoile_rot::get_omega_c()
177  * (retourne la valeur centrale de Omega).
178  *
179  * Revision 2.58 2001/10/11 09:24:00 eric
180  * La fonction Etoile_rot::equilibrium est desormais virtuelle.
181  *
182  * Revision 2.57 2001/08/06 15:39:04 keisuke
183  * Addition of a new argument to Etoile_bin::equilibrium and equil_regular.
184  *
185  * Revision 2.56 2001/06/25 12:52:33 eric
186  * Classe Etoile_bin : ajout du membre p_companion et des fonctions
187  * associees set_companion() et get_companion().
188  *
189  * Revision 2.55 2001/06/13 14:11:55 eric
190  * Modif commentaires (mise en conformite Doc++ 3.4.7)
191  *
192  * Revision 2.54 2001/03/26 09:29:59 jlz
193  * Classe Etoile_rot: new members p_espec_isco and p_lspec_isco.
194  *
195  * Revision 2.53 2001/02/08 15:37:09 eric
196  * *** empty log message ***
197  *
198  * Revision 2.52 2001/02/08 15:12:42 eric
199  * Ajout de la fonction Etoile_rot::f_eccentric.
200  *
201  * Revision 2.51 2000/11/23 15:43:24 eric
202  * Ajout de l'argument ent_limit a Etoile_rot::equilibrium.
203  *
204  * Revision 2.50 2000/11/19 18:51:13 eric
205  * Etoile_rot: ajout de la fonction (static) lambda_grv2
206  *
207  * Revision 2.49 2000/11/18 23:17:32 eric
208  * Ajout de l'argument ost a Etoile_rot::r_isco.
209  *
210  * Revision 2.48 2000/11/18 21:08:33 eric
211  * Classe Etoile_rot: ajout des fonctions r_isco() et f_isco()
212  * ainsi que des membres associes p_r_isco et p_f_isco.
213  *
214  * Revision 2.47 2000/11/10 15:15:38 eric
215  * Modif arguments Etoile_rot::equilibrium.
216  *
217  * Revision 2.46 2000/10/20 13:10:23 eric
218  * Etoile_rot::equilibrium: ajout de l'argument nzadapt.
219  *
220  * Revision 2.45 2000/10/17 15:59:14 eric
221  * Modif commentaires Etoile_rot::equilibrium
222  *
223  * Revision 2.44 2000/10/12 15:22:29 eric
224  * Modif commentaires Etoile_rot.
225  *
226  * Revision 2.43 2000/10/12 10:20:12 eric
227  * Ajout de la fonction Etoile_rot::fait_nphi().
228  *
229  * Revision 2.42 2000/09/22 15:50:13 keisuke
230  * d_logn_auto_div devient desormais un membre de la classe Etoile
231  * et non plus de la classe derivee Etoile_bin.
232  *
233  * Revision 2.41 2000/09/18 16:14:37 eric
234  * Classe Etoile_rot: ajout du membre tkij et de la fonction
235  * extrinsic curvature().
236  *
237  * Revision 2.40 2000/09/07 14:31:09 keisuke
238  * Ajout des membres logn_auto_regu et get_logn_auto_regu() a la classe Etoile.
239  * Ajout des membres d_logn_auto_regu et get_d_logn_auto_regu()
240  * a la classe Etoile_bin.
241  *
242  * Revision 2.39 2000/09/07 10:25:50 keisuke
243  * Ajout du membre get_logn_auto_div() a la classe Etoile.
244  * Ajout du membre get_d_logn_auto_div() a la classe Etoile_bin.
245  *
246  * Revision 2.38 2000/08/31 11:25:24 eric
247  * Classe Etoile_rot: ajout des membres tnphi et ak_car.
248  *
249  * Revision 2.37 2000/08/29 11:37:06 eric
250  * Ajout des membres k_div et logn_auto_div a la classe Etoile.
251  * Ajout du membre d_logn_auto_div a la classe Etoile_bin.
252  *
253  * Revision 2.36 2000/08/25 10:25:57 keisuke
254  * Ajout de Etoile_bin::equil_regular.
255  *
256  * Revision 2.35 2000/08/18 14:01:07 eric
257  * Modif commentaires.
258  *
259  * Revision 2.34 2000/08/17 12:38:30 eric
260  * Modifs classe Etoile_rot : ajout des membres nuf, nuq, ssjm1_nuf et ssjm1_nuq
261  * Modif arguments Etoile_rot::equilibrium.
262  * .\
263  *
264  * Revision 2.33 2000/08/07 12:11:13 keisuke
265  * Ajout de Etoile::equil_spher_regular.
266  *
267  * Revision 2.32 2000/07/21 12:02:55 eric
268  * Suppression de Etoile_rot::relaxation.
269  *
270  * Revision 2.31 2000/07/20 15:32:28 eric
271  * *** empty log message ***
272  *
273  * Revision 2.30 2000/07/06 09:39:12 eric
274  * Ajout du membre p_xa_barycenter a Etoile_bin, ainsi que de la
275  * fonction associee xa_barycenter().
276  *
277  * Revision 2.29 2000/05/25 13:47:38 eric
278  * Modif Etoile_bin::equilibrium: ajout de l'argument thres_adapt
279  *
280  * Revision 2.28 2000/03/22 16:41:45 eric
281  * Ajout de la sortie de nouvelles erreurs dans Etoile_bin::equilibrium.
282  *
283  * Revision 2.27 2000/03/22 12:54:44 eric
284  * Nouveau prototype d'Etoile_bin::velocity_potential : l'erreur est
285  * retournee en double.
286  * Nouveau prototype d'Etoile_bin::equilibrium : diff_ent est remplace
287  * par le Tbl diff.
288  *
289  * Revision 2.26 2000/03/15 11:04:15 eric
290  * Ajout des fonctions Etoile_bin::set_w_shift() et Etoile_bin::set_khi_shift()
291  * Amelioration des commentaires.
292  *
293  * Revision 2.25 2000/03/08 12:12:49 eric
294  * Ajout de la fonction Etoile_bin::is_irrotational().
295  *
296  * Revision 2.24 2000/03/07 14:48:01 eric
297  * Ajout de la fonction Etoile_bin::extrinsic_curvature().
298  *
299  * Revision 2.23 2000/02/21 13:57:57 eric
300  * classe Etoile_bin: suppression du membre psi: remplacement par psi0.
301  *
302  * Revision 2.22 2000/02/17 16:51:22 eric
303  * Ajout de l'argument diff_ent dans Etoile_bin::equilibrium.
304  *
305  * Revision 2.21 2000/02/17 15:29:38 eric
306  * Ajout de la fonction Etoile_bin::relaxation.
307  *
308  * Revision 2.20 2000/02/17 13:54:21 eric
309  * Ajout de la fonction Etoile_bin::velocity_potential.
310  *
311  * Revision 2.19 2000/02/16 15:05:14 eric
312  * Ajout des membres w_shift et khi_shift.
313  * (sauves dans les fichiers a la place de shift_auto).
314  * Ajout de la fonction Etoile_bin::fait_shift_auto.
315  *
316  * Revision 2.18 2000/02/16 13:47:02 eric
317  * Classe Etoile_bin: ajout des membres ssjm1_khi et ssjm1_wshift.
318  *
319  * Revision 2.17 2000/02/16 11:54:13 eric
320  * Classe Etoile_bin : ajout des membres ssjm1_logn et ssjm1_beta.
321  *
322  * Revision 2.16 2000/02/15 15:40:07 eric
323  * Ajout de Etoile_bin::equilibrium.
324  *
325  * Revision 2.15 2000/02/12 18:40:15 eric
326  * Modif commentaires.
327  *
328  * Revision 2.14 2000/02/12 14:44:26 eric
329  * Ajout des fonctions Etoile_bin::set_logn_comp et set_pot_centri.
330  *
331  * Revision 2.13 2000/02/10 20:22:25 eric
332  * Modif commentaires.
333  *
334  * Revision 2.12 2000/02/10 16:11:24 eric
335  * Classe Etoile_bin : ajout des accesseurs get_psi, etc...
336  * ajout de la fonction fait_d_psi
337  *
338  * Revision 2.11 2000/02/08 19:28:29 eric
339  * La fonction Etoile_bin::scal_prod est rebaptisee Etoile_bin::sprod
340  *
341  * Revision 2.10 2000/02/04 17:15:15 eric
342  * Classe Etoile_bin: ajout du membre ref_triad.
343  *
344  * Revision 2.9 2000/02/04 16:36:48 eric
345  * Ajout des fonctions update_metric* et kinematics.
346  *
347  * Revision 2.8 2000/02/02 10:12:37 eric
348  * Ajout des fonctions de lecture/ecriture mp, nzet, eos, etc...
349  *
350  * Revision 2.7 2000/02/01 15:59:43 eric
351  * Ajout de la fonction Etoile_bin::scal_prod.
352  *
353  * Revision 2.6 2000/01/31 15:56:45 eric
354  * Introduction de la classe derivee Etoile_bin.
355  *
356  * Revision 2.5 2000/01/28 17:17:45 eric
357  * Ajout des fonctions de calcul des quantites globales.
358  *
359  * Revision 2.4 2000/01/27 16:46:59 eric
360  * Ajout des fonctions get_ent(), etc....
361  *
362  * Revision 2.3 2000/01/24 17:19:48 eric
363  * Modif commentaires.
364  *
365  * Revision 2.2 2000/01/24 17:13:04 eric
366  * Le mapping mp n'est plus constant.
367  * Ajout de la fonction equilibrium_spher.
368  *
369  * Revision 2.1 2000/01/24 13:37:19 eric
370  * *** empty log message ***
371  *
372  * Revision 2.0 2000/01/20 17:04:33 eric
373  * *** empty log message ***
374  *
375  *
376  * $Header: /cvsroot/Lorene/C++/Include/etoile.h,v 1.37 2026/02/26 09:53:04 j_novak Exp $
377  *
378  */
379 
380 // Headers Lorene
381 #include "tenseur.h"
382 namespace Lorene {
383  class Eos ;
384  class Bhole ;
385 
386 
387  //---------------------------//
388  // base class Etoile //
389  //---------------------------//
390 
430  class Etoile {
431 
432  // Data :
433  // -----
434  protected:
435  Map& mp ;
436 
438  int nzet ;
439 
444 
448  double unsurc2 ;
449 
455  int k_div ;
456 
457  const Eos& eos ;
458 
459  // Fluid quantities with respect to the fluid frame
460  // ------------------------------------------------
461 
464 
468 
469  // Fluid quantities with respect to the Eulerian frame
470  // ---------------------------------------------------
472 
475 
478 
481 
482  // Metric potentials
483  // -----------------
484 
491 
498 
504 
508 
513 
516 
519 
522 
523  // Derived data :
524  // ------------
525  protected:
527  mutable double* p_ray_eq ;
528 
530  mutable double* p_ray_eq_pis2 ;
531 
533  mutable double* p_ray_eq_pi ;
534 
536  mutable double* p_ray_eq_3pis2 ;
537 
539  mutable double* p_ray_pole ;
540 
545  mutable Itbl* p_l_surf ;
546 
551  mutable Tbl* p_xi_surf ;
552 
553  mutable double* p_mass_b ;
554  mutable double* p_mass_g ;
555 
556 
557  // Constructors - Destructor
558  // -------------------------
559  public:
560 
570  Etoile(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
571 
572  Etoile(const Etoile& ) ;
573 
581  Etoile(Map& mp_i, const Eos& eos_i, FILE* fich) ;
582 
583  virtual ~Etoile() ;
584 
585  // Memory management
586  // -----------------
587  protected:
589  virtual void del_deriv() const ;
590 
592  virtual void set_der_0x0() const ;
593 
597  virtual void del_hydro_euler() ;
598 
599 
600  // Mutators / assignment
601  // ---------------------
602  public:
604  void operator=(const Etoile&) ;
605 
607  Map& set_mp() {return mp; } ;
608 
610  void set_enthalpy(const Cmp& ) ;
611 
615  virtual void equation_of_state() ;
616 
621  virtual void hydro_euler() ;
622 
635  virtual void equilibrium_spher(double ent_c, double precis = 1.e-14,
636  const Tbl* ent_limit = 0x0 ) ;
637 
648  void equil_spher_regular(double ent_c, double precis = 1.e-14) ;
649 
658  virtual void equil_spher_falloff(double ent_c,
659  double precis = 1.e-14) ;
660 
661  // Accessors
662  // ---------
663  public:
665  const Map& get_mp() const {return mp; } ;
666 
668  int get_nzet() const {return nzet; } ;
669 
673  bool is_relativistic() const {return relativistic; } ;
674 
676  const Eos& get_eos() const {return eos; } ;
677 
679  const Tenseur& get_ent() const {return ent;} ;
680 
682  const Tenseur& get_nbar() const {return nbar;} ;
683 
685  const Tenseur& get_ener() const {return ener;} ;
686 
688  const Tenseur& get_press() const {return press;} ;
689 
691  const Tenseur& get_ener_euler() const {return ener_euler;} ;
692 
694  const Tenseur& get_s_euler() const {return s_euler;} ;
695 
697  const Tenseur& get_gam_euler() const {return gam_euler;} ;
698 
700  const Tenseur& get_u_euler() const {return u_euler;} ;
701 
707  const Tenseur& get_logn_auto() const {return logn_auto;} ;
708 
714  const Tenseur& get_logn_auto_regu() const {return logn_auto_regu;} ;
715 
721  const Tenseur& get_logn_auto_div() const {return logn_auto_div;} ;
722 
725  const Tenseur& get_d_logn_auto_div() const {return d_logn_auto_div;} ;
726 
730  const Tenseur& get_beta_auto() const {return beta_auto;} ;
731 
733  const Tenseur& get_nnn() const {return nnn;} ;
734 
736  const Tenseur& get_shift() const {return shift;} ;
737 
739  const Tenseur& get_a_car() const {return a_car;} ;
740 
741  // Outputs
742  // -------
743  public:
744  virtual void sauve(FILE* ) const ;
745 
747  friend ostream& operator<<(ostream& , const Etoile& ) ;
748 
749  protected:
751  virtual ostream& operator>>(ostream& ) const ;
752 
753  // Global quantities
754  // -----------------
755  public:
757  double ray_eq() const ;
758 
760  double ray_eq_pis2() const ;
761 
763  double ray_eq_pi() const ;
764 
766  double ray_eq_3pis2() const ;
767 
769  double ray_pole() const ;
770 
772  double ray_eq(int kk) const ;
773 
781  virtual const Itbl& l_surf() const ;
782 
790  const Tbl& xi_surf() const ;
791 
793  virtual double mass_b() const ;
794 
796  virtual double mass_g() const ;
797 
798  };
799  ostream& operator<<(ostream& , const Etoile& ) ;
800 
801 
802  //---------------------------//
803  // class Etoile_bin //
804  //---------------------------//
805 
820  class Etoile_bin : public Etoile {
821 
822  // Data :
823  // -----
824  protected:
828  bool irrotational ;
829 
835 
840 
845 
851 
856 
861 
866 
871 
876 
881 
886 
891 
896 
902 
915 
925 
932 
939 
945 
951 
957 
960 
966 
972 
980 
990 
996 
1007 
1008  // Derived data :
1009  // ------------
1010  protected:
1012  mutable double* p_xa_barycenter ;
1013 
1014 
1015  // Constructors - Destructor
1016  // -------------------------
1017  public:
1033  Etoile_bin(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i,
1034  bool irrot, const Base_vect& ref_triad_i) ;
1035 
1036 
1037  Etoile_bin(const Etoile_bin& ) ;
1038 
1051  Etoile_bin(Map& mp_i, const Eos& eos_i, const Base_vect& ref_triad_i,
1052  FILE* fich) ;
1053 
1054  virtual ~Etoile_bin() ;
1055 
1056 
1057  // Memory management
1058  // -----------------
1059  protected:
1061  virtual void del_deriv() const ;
1062 
1064  virtual void set_der_0x0() const ;
1065 
1069  virtual void del_hydro_euler() ;
1070 
1071 
1072  // Mutators / assignment
1073  // ---------------------
1074  public:
1076  void operator=(const Etoile_bin& ) ;
1077 
1081  Tenseur& set_logn_comp() ;
1082 
1084  Tenseur& set_pot_centri() ;
1085 
1087  Tenseur& set_w_shift() ;
1088 
1090  Tenseur& set_khi_shift() ;
1091 
1092  // Accessors
1093  // ---------
1094  public:
1098  bool is_irrotational() const {return irrotational; } ;
1099 
1101  const Tenseur& get_psi0() const {return psi0;} ;
1102 
1106  const Tenseur& get_d_psi() const {return d_psi;} ;
1107 
1112  const Tenseur& get_wit_w() const {return wit_w;} ;
1113 
1117  const Tenseur& get_loggam() const {return loggam;} ;
1118 
1122  const Tenseur& get_logn_comp() const {return logn_comp;} ;
1123 
1127  const Tenseur& get_d_logn_auto() const {return d_logn_auto;} ;
1128 
1132  const Tenseur& get_d_logn_auto_regu() const {return d_logn_auto_regu;} ;
1133 
1137  const Tenseur& get_d_logn_comp() const {return d_logn_comp;} ;
1138 
1142  const Tenseur& get_beta_comp() const {return beta_comp;} ;
1143 
1147  const Tenseur& get_d_beta_auto() const {return d_beta_auto;} ;
1148 
1152  const Tenseur& get_d_beta_comp() const {return d_beta_comp;} ;
1153 
1158  const Tenseur& get_shift_auto() const {return shift_auto;} ;
1159 
1164  const Tenseur& get_shift_comp() const {return shift_comp;} ;
1165 
1178  const Tenseur& get_w_shift() const {return w_shift;} ;
1179 
1192  const Tenseur& get_khi_shift() const {return khi_shift;} ;
1193 
1198  const Tenseur_sym& get_tkij_auto() const {return tkij_auto;} ;
1199 
1204  const Tenseur_sym& get_tkij_comp() const {return tkij_comp;} ;
1205 
1210  const Tenseur& get_akcar_auto() const {return akcar_auto;} ;
1211 
1216  const Tenseur& get_akcar_comp() const {return akcar_comp;} ;
1217 
1222  const Tenseur& get_bsn() const {return bsn;} ;
1223 
1225  const Tenseur& get_pot_centri() const {return pot_centri;} ;
1230  const Cmp get_decouple() const {return decouple ;}
1231 
1232  // Outputs
1233  // -------
1234  public:
1235  virtual void sauve(FILE* ) const ;
1236 
1237  protected:
1239  virtual ostream& operator>>(ostream& ) const ;
1240 
1241  // Global quantities
1242  // -----------------
1243  public:
1245  virtual double mass_b() const ;
1246 
1248  virtual double mass_g() const ;
1249 
1258  virtual double xa_barycenter() const ;
1259 
1260 
1261  // Computational routines
1262  // ----------------------
1263  public:
1271  virtual Tenseur sprod(const Tenseur& t1, const Tenseur& t2) const ;
1272 
1285  virtual void hydro_euler() ;
1286 
1304  void update_metric(const Etoile_bin& comp) ;
1305 
1313  void update_metric(const Etoile_bin& comp, const Etoile_bin& star_prev,
1314  double relax) ;
1315 
1330  void update_metric_der_comp(const Etoile_bin& comp) ;
1331 
1342  virtual void kinematics(double omega, double x_axe) ;
1343 
1347  void fait_d_psi() ;
1348 
1357  void fait_shift_auto() ;
1358 
1362  virtual void extrinsic_curvature() ;
1363 
1404  void equilibrium(double ent_c,
1405  int mermax, int mermax_poisson,
1406  double relax_poisson, int mermax_potvit,
1407  double relax_potvit, double thres_adapt,
1408  const Tbl& fact, Tbl& diff, const Tbl* ent_limit = 0x0) ;
1409 
1449  void equil_regular(double ent_c, int mermax, int mermax_poisson,
1450  double relax_poisson, int mermax_potvit,
1451  double relax_potvit, double thres_adapt,
1452  const Tbl& fact, Tbl& diff) ;
1453 
1467  double velocity_potential(int mermax, double precis, double relax) ;
1468 
1479  void relaxation(const Etoile_bin& star_prev, double relax_ent,
1480  double relax_met, int mer, int fmer_met) ;
1481 
1482  friend class Bin_ns_bh ;
1483  };
1484 
1485 
1486 
1487  //---------------------------//
1488  // class Etoile_rot //
1489  //---------------------------//
1490 
1502 class Etoile_rot : public Etoile {
1503 
1504  // Data :
1505  // -----
1506  protected:
1507  double omega ;
1508 
1511 
1514 
1517 
1522 
1525 
1528 
1533 
1538 
1541 
1544 
1557 
1567 
1574 
1593 
1599 
1605 
1610 
1615 
1623 
1632 
1633  // Derived data :
1634  // ------------
1635  protected:
1636 
1637  mutable double* p_angu_mom ;
1638  mutable double* p_tsw ;
1639  mutable double* p_grv2 ;
1640  mutable double* p_grv3 ;
1641  mutable double* p_r_circ ;
1642  mutable double* p_r_circ_merid ;
1643  mutable double* p_area ;
1644  mutable double* p_aplat ;
1645  mutable double* p_z_eqf ;
1646  mutable double* p_z_eqb ;
1647  mutable double* p_z_pole ;
1648  mutable double* p_mom_quad ;
1649  mutable double* p_mom_quad_old ;
1650  mutable double* p_mom_quad_Bo ;
1651  mutable double* p_r_isco ;
1652  mutable double* p_f_isco ;
1653  mutable double* p_espec_isco ;
1656  mutable double* p_lspec_isco ;
1657  mutable double* p_f_eq ;
1658  mutable Tbl* p_radius_Morsink ;
1660 
1661 
1662 
1663 
1664  // Constructors - Destructor
1665  // -------------------------
1666  public:
1675  Etoile_rot(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i) ;
1676 
1677 
1678  Etoile_rot(const Etoile_rot& ) ;
1679 
1687  Etoile_rot(Map& mp_i, const Eos& eos_i, FILE* fich) ;
1688 
1689  virtual ~Etoile_rot() ;
1690 
1691 
1692  // Memory management
1693  // -----------------
1694  protected:
1696  virtual void del_deriv() const ;
1697 
1699  virtual void set_der_0x0() const ;
1700 
1704  virtual void del_hydro_euler() ;
1705 
1706 
1707  // Mutators / assignment
1708  // ---------------------
1709  public:
1711  void operator=(const Etoile_rot& ) ;
1712 
1713  // Accessors
1714  // ---------
1715  public:
1719  virtual double get_omega_c() const ;
1720 
1722  const Tenseur& get_bbb() const {return bbb;} ;
1723 
1725  const Tenseur& get_b_car() const {return b_car;} ;
1726 
1728  const Tenseur& get_nphi() const {return nphi;} ;
1729 
1733  const Tenseur& get_tnphi() const {return tnphi;} ;
1734 
1736  const Tenseur& get_uuu() const {return uuu;} ;
1737 
1739  const Tenseur& get_logn() const {return logn;} ;
1740 
1744  const Tenseur& get_nuf() const {return nuf;} ;
1745 
1749  const Tenseur& get_nuq() const {return nuq;} ;
1750 
1752  const Tenseur& get_dzeta() const {return dzeta;} ;
1753 
1755  const Tenseur& get_tggg() const {return tggg;} ;
1756 
1769  const Tenseur& get_w_shift() const {return w_shift;} ;
1770 
1783  const Tenseur& get_khi_shift() const {return khi_shift;} ;
1784 
1790  const Tenseur_sym& get_tkij() const {return tkij;} ;
1791 
1809  const Tenseur& get_ak_car() const {return ak_car;} ;
1810 
1811  // Outputs
1812  // -------
1813  public:
1814  virtual void sauve(FILE* ) const ;
1815 
1817  virtual void display_poly(ostream& ) const ;
1818 
1819  protected:
1821  virtual ostream& operator>>(ostream& ) const ;
1822 
1824  virtual void partial_display(ostream& ) const ;
1825 
1826  // Global quantities
1827  // -----------------
1828  public:
1829 
1837  virtual const Itbl& l_surf() const ;
1838 
1839  virtual double mass_b() const ;
1840  virtual double mass_g() const ;
1841  virtual double angu_mom() const ;
1842  virtual double tsw() const ;
1843 
1847  virtual double grv2() const ;
1848 
1860  virtual double grv3(ostream* ost = 0x0) const ;
1861 
1862  virtual double r_circ() const ;
1863  virtual double r_circ_merid() const ;
1864  virtual double area() const ;
1865  virtual double mean_radius() const ;
1866  virtual double aplat() const ;
1867  virtual double z_eqf() const ;
1868  virtual double z_eqb() const ;
1869  virtual double z_pole() const ;
1870 
1885  virtual double mom_quad() const ;
1886 
1894  virtual double mom_quad_old() const ;
1895 
1902  virtual double mom_quad_Bo() const ;
1903 
1907  virtual const Tbl& radius_Morsink() const ;
1908 
1915  virtual double r_isco(ostream* ost = 0x0) const ;
1916 
1918  virtual double f_isco() const ;
1919 
1921  virtual double espec_isco() const ;
1922 
1924  virtual double lspec_isco() const ;
1925 
1926 
1937  virtual double f_eccentric(double ecc, double periast,
1938  ostream* ost = 0x0) const ;
1939 
1941  virtual double f_eq() const ;
1942 
1943 
1944  // Computational routines
1945  // ----------------------
1946  public:
1957  virtual void hydro_euler() ;
1958 
1968  void update_metric() ;
1969 
1978  void fait_shift() ;
1979 
1983  void fait_nphi() ;
1984 
1988  void extrinsic_curvature() ;
1989 
2019  static double lambda_grv2(const Cmp& sou_m, const Cmp& sou_q) ;
2020 
2100  virtual void equilibrium(double ent_c, double omega0, double fact_omega,
2101  int nzadapt, const Tbl& ent_limit,
2102  const Itbl& icontrol, const Tbl& control,
2103  double mbar_wanted, double aexp_mass,
2104  Tbl& diff, Param* = 0x0) ;
2105 
2106 
2107 };
2108 
2109 
2110 
2111 
2112 }
2113 #endif
double * p_z_eqb
Backward redshift factor at equator.
Definition: etoile.h:1646
virtual double z_eqf() const
Forward redshift factor at equator.
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:714
const Tenseur & get_nuq() const
Returns the Part of the Metric potential = logn generated by the quadratic terms.
Definition: etoile.h:1749
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: etoile.C:399
double * p_z_pole
Redshift factor at North pole.
Definition: etoile.h:1647
Base class for stars *** DEPRECATED : use class Star instead ***.
Definition: etoile.h:430
virtual double mom_quad_Bo() const
Part of the quadrupole moment.
int k_div
Index of regularity of the gravitational potential logn_auto .
Definition: etoile.h:455
virtual double r_circ() const
Circumferential equatorial radius.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: etoile.C:514
const Tbl & xi_surf() const
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate...
double * p_mom_quad_old
Part of the quadrupole moment.
Definition: etoile.h:1649
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: etoile.h:1566
Tbl * p_radius_Morsink
Physical radius close to the definition of Morsink et al. (2007)
Definition: etoile.h:1659
Tenseur shift_comp
Part of the shift vector generated principaly by the companion star.
Definition: etoile.h:901
const Tenseur & get_beta_comp() const
Returns the part of the logarithm of AN generated principaly by the companion star.
Definition: etoile.h:1142
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur &#39;s...
Definition: etoile.h:834
double * p_r_circ
Circumferential equatorial radius.
Definition: etoile.h:1641
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: etoile.h:1622
double * p_mom_quad_Bo
Part of the quadrupole moment.
Definition: etoile.h:1650
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
const Cmp get_decouple() const
Returns the function used to construct tkij_auto from tkij_tot .
Definition: etoile.h:1230
virtual double lspec_isco() const
Angular momentum of a particle on the ISCO.
Definition: et_rot_isco.C:287
double * p_r_isco
Circumferential radius of the ISCO.
Definition: etoile.h:1651
virtual void sauve(FILE *) const
Save in a file.
Definition: etoile_bin.C:562
Tenseur tnphi
Component of the shift vector.
Definition: etoile.h:1521
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:1132
const Map & get_mp() const
Returns the mapping.
Definition: etoile.h:665
const Tenseur & get_a_car() const
Returns the total conformal factor .
Definition: etoile.h:739
Neutron star - black hole binary system.
Definition: bin_ns_bh.h:120
void update_metric(const Etoile_bin &comp)
Computes metric coefficients from known potentials, when the companion is another star...
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:1112
const Tenseur_sym & get_tkij_comp() const
Returns the part of the extrinsic curvature tensor generated by shift_comp .
Definition: etoile.h:1204
void fait_shift_auto()
Computes shift_auto from w_shift and khi_shift according to Shibata&#39;s prescription [Prog...
Definition: etoile_bin.C:829
const Tenseur & get_tggg() const
Returns the Metric potential .
Definition: etoile.h:1755
double ray_eq_pis2() const
Coordinate radius at , [r_unit].
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: etoile_rot.C:713
virtual double mass_g() const
Gravitational mass.
Class intended to describe tensors with a symmetry on the two last indices *** DEPRECATED : use class...
Definition: tenseur.h:1256
Tenseur pot_centri
Centrifugal potential.
Definition: etoile.h:959
virtual void sauve(FILE *) const
Save in a file.
Definition: etoile.C:486
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Definition: etoile_rot.C:808
Cmp ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto by mea...
Definition: etoile.h:965
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:497
Lorene prototypes.
Definition: app_hor.h:67
const Tenseur & get_akcar_comp() const
Returns the part of the scalar generated by shift_auto and shift_comp , i.e.
Definition: etoile.h:1216
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: etoile.h:1556
Equation of state base class.
Definition: eos.h:209
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: etoile.h:1631
const Tenseur_sym & get_tkij_auto() const
Returns the part of the extrinsic curvature tensor generated by shift_auto .
Definition: etoile.h:1198
Tenseur nnn
Total lapse function.
Definition: etoile.h:515
virtual void equilibrium_spher(double ent_c, double precis=1.e-14, const Tbl *ent_limit=0x0)
Computes a spherical static configuration.
const Tenseur & get_d_beta_auto() const
Returns the gradient of beta_auto (Cartesian components with respect to ref_triad ) ...
Definition: etoile.h:1147
double * p_z_eqf
Forward redshift factor at equator.
Definition: etoile.h:1645
Tenseur nphi
Metric coefficient .
Definition: etoile.h:1516
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:474
const Tenseur & get_shift() const
Returns the total shift vector .
Definition: etoile.h:736
double ray_eq() const
Coordinate radius at , [r_unit].
Base class for coordinate mappings.
Definition: map.h:697
Class for isolated rotating stars *** DEPRECATED : use class Star_rot instead ***.
Definition: etoile.h:1502
double * p_mom_quad
Quadrupole moment.
Definition: etoile.h:1648
Tenseur b_car
Square of the metric factor B.
Definition: etoile.h:1513
const Tenseur & get_b_car() const
Returns the square of the metric factor B.
Definition: etoile.h:1725
int get_nzet() const
Returns the number of domains occupied by the star.
Definition: etoile.h:668
void operator=(const Etoile_rot &)
Assignment to another Etoile_rot.
Definition: etoile_rot.C:420
Tenseur wit_w
Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Definition: etoile.h:850
const Tenseur & get_logn_auto() const
Returns the logarithm of the part of the lapse N generated principaly by the star.
Definition: etoile.h:707
virtual double mom_quad_old() const
Part of the quadrupole moment.
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: etoile_bin.C:470
Basic integer array class.
Definition: itbl.h:122
const Tenseur & get_d_psi() const
Returns the gradient of the velocity potential (Cartesian components with respect to ref_triad ) ...
Definition: etoile.h:1106
bool is_irrotational() const
Returns true for an irrotational motion, false for a corotating one.
Definition: etoile.h:1098
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition: etoile.h:448
virtual double mass_g() const
Gravitational mass.
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
Definition: et_rot_isco.C:87
Class for stars in binary system.
Definition: etoile.h:820
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:721
virtual double xa_barycenter() const
Absolute coordinate X of the barycenter of the baryon density, defined according to the formula wher...
Tenseur press
Fluid pressure.
Definition: etoile.h:467
const Tenseur & get_uuu() const
Returns the norm of u_euler.
Definition: etoile.h:1736
virtual void extrinsic_curvature()
Computes tkij_auto and akcar_auto from shift_auto , nnn and a_car .
const Tenseur & get_d_logn_auto() const
Returns the gradient of logn_auto (Cartesian components with respect to ref_triad ) ...
Definition: etoile.h:1127
Tenseur d_beta_auto
Gradient of beta_auto (Cartesian components with respect to ref_triad )
Definition: etoile.h:885
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...
Etoile_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition: etoile_rot.C:149
virtual double mean_radius() const
Mean radius.
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition: base_vect.h:105
void set_enthalpy(const Cmp &)
Assignment of the enthalpy field.
Definition: etoile.C:468
Tenseur shift
Total shift vector.
Definition: etoile.h:518
const Tenseur & get_nuf() const
Returns the part of the Metric potential = logn generated by the matter terms.
Definition: etoile.h:1744
double * p_aplat
Flatening r_pole/r_eq.
Definition: etoile.h:1644
const Tenseur & get_nphi() const
Returns the metric coefficient .
Definition: etoile.h:1728
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_rot.C:347
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition: etoile.h:895
const Tenseur & get_d_logn_auto_div() const
Returns the gradient of logn_auto_div.
Definition: etoile.h:725
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
virtual double r_circ_merid() const
Circumferential meridional radius.
double * p_ray_pole
Coordinate radius at .
Definition: etoile.h:539
const Tenseur & get_press() const
Returns the fluid pressure.
Definition: etoile.h:688
const Tenseur & get_logn_comp() const
Returns the part of the lapse logarithm (gravitational potential at the Newtonian limit) generated pr...
Definition: etoile.h:1122
void extrinsic_curvature()
Computes tkij and ak_car from shift , nnn and b_car .
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:503
const Tenseur & get_w_shift() const
Returns the vector used in the decomposition of shift_auto , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:1178
bool irrotational
true for an irrotational star, false for a corotating one
Definition: etoile.h:828
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile.C:381
void update_metric()
Computes metric coefficients from known potentials.
Definition: et_rot_upmetr.C:72
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
Definition: etoile_rot.C:703
Tenseur psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case) ...
Definition: etoile.h:839
Cmp ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta ...
Definition: etoile.h:1609
Tenseur w_shift
Vector used in the decomposition of shift_auto , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:914
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:480
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition: etoile.h:1604
double * p_r_circ_merid
Circumferential meridional radius.
Definition: etoile.h:1642
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: etoile_rot.C:482
virtual ~Etoile_rot()
Destructor.
Definition: etoile_rot.C:337
double * p_grv3
Error on the virial identity GRV3.
Definition: etoile.h:1640
const Tenseur & get_loggam() const
Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer...
Definition: etoile.h:1117
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:890
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_bin.C:450
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: etoile_bin.C:462
virtual double angu_mom() const
Angular momentum.
Tenseur logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principaly by ...
Definition: etoile.h:860
virtual ~Etoile_bin()
Destructor.
Definition: etoile_bin.C:440
Map & set_mp()
Read/write of the mapping.
Definition: etoile.h:607
Tenseur d_logn_auto
Gradient of logn_auto (Cartesian components with respect to ref_triad )
Definition: etoile.h:865
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition: etoile.h:956
virtual double z_pole() const
Redshift factor at North pole.
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:465
void fait_d_psi()
Computes the gradient of the total velocity potential .
Definition: etoile_bin.C:767
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
const Tenseur & get_khi_shift() const
Returns the scalar used in the decomposition of shift_auto following Shibata&#39;s prescription [Prog...
Definition: etoile.h:1192
virtual void sauve(FILE *) const
Save in a file.
Definition: etoile_rot.C:457
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:477
Cmp decouple
Function used to construct the part of generated by the star from the total .
Definition: etoile.h:1006
const Tenseur & get_s_euler() const
Returns the trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:694
Tenseur akcar_comp
Part of the scalar generated by shift_auto and shift_comp , i.e.
Definition: etoile.h:950
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
double * p_tsw
Ratio T/W.
Definition: etoile.h:1638
virtual double grv2() const
Error on the virial identity GRV2.
virtual double mass_b() const
Baryon mass.
const Tenseur & get_logn() const
Returns the metric potential = logn_auto.
Definition: etoile.h:1739
double * p_area
Surface area.
Definition: etoile.h:1643
const Tenseur_sym & get_tkij() const
Returns the tensor related to the extrinsic curvature tensor by .
Definition: etoile.h:1790
Parameter storage.
Definition: param.h:125
Tenseur beta_comp
Part of the logarithm of AN generated principaly by the companion star.
Definition: etoile.h:880
const Tenseur & get_bsn() const
Returns the shift vector, divided by N , of the rotating coordinates, .
Definition: etoile.h:1222
void equil_spher_regular(double ent_c, double precis=1.e-14)
Computes a spherical static configuration.
virtual double mass_b() const
Baryon mass.
Map & mp
Mapping associated with the star.
Definition: etoile.h:435
virtual double espec_isco() const
Energy of a particle on the ISCO.
Definition: et_rot_isco.C:304
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: etoile.C:569
const Eos & eos
Equation of state of the stellar matter.
Definition: etoile.h:457
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.
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:551
const Tenseur & get_khi_shift() const
Returns the scalar used in the decomposition of shift following Shibata&#39;s prescription [Prog...
Definition: etoile.h:1783
const Tenseur & get_nnn() const
Returns the total lapse function N.
Definition: etoile.h:733
virtual double z_eqb() const
Backward redshift factor at equator.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: etoile.C:687
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition: etoile.h:1598
Tenseur ak_car
Scalar .
Definition: etoile.h:1592
static double lambda_grv2(const Cmp &sou_m, const Cmp &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
double * p_f_isco
Orbital frequency of the ISCO.
Definition: etoile.h:1652
double * p_ray_eq
Coordinate radius at , .
Definition: etoile.h:527
virtual double tsw() const
Ratio T/W.
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:924
virtual double f_eq() const
Orbital frequency at the equator.
Definition: et_rot_isco.C:322
Tenseur & set_w_shift()
Read/write of w_shift.
Definition: etoile_bin.C:541
const Tenseur & get_pot_centri() const
Returns the centrifugal potential.
Definition: etoile.h:1225
Tenseur bbb
Metric factor B.
Definition: etoile.h:1510
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: et_rot_hydro.C:86
const Tenseur & get_dzeta() const
Returns the Metric potential = beta_auto.
Definition: etoile.h:1752
const Tenseur & get_gam_euler() const
Returns the Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:697
Etoile(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition: etoile.C:146
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor generated by shift_auto .
Definition: etoile.h:931
const Tenseur & get_d_beta_comp() const
Returns the gradient of beta_comp (Cartesian components with respect to ref_triad ) ...
Definition: etoile.h:1152
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 uuu
Norm of u_euler.
Definition: etoile.h:1524
const Tenseur & get_ener_euler() const
Returns the total energy density with respect to the Eulerian observer.
Definition: etoile.h:691
double * p_mass_g
Gravitational mass.
Definition: etoile.h:554
Tenseur tggg
Metric potential .
Definition: etoile.h:1543
double omega
Rotation angular velocity ([f_unit] )
Definition: etoile.h:1507
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition: etoile_rot.C:650
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: etoile.h:1532
void operator=(const Etoile_bin &)
Assignment to another Etoile_bin.
Definition: etoile_bin.C:485
const Tenseur & get_d_logn_comp() const
Returns the gradient of logn_comp (Cartesian components with respect to ref_triad ) ...
Definition: etoile.h:1137
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:438
virtual double mom_quad() const
Quadrupole moment.
Tenseur a_car
Total conformal factor .
Definition: etoile.h:521
virtual ~Etoile()
Destructor.
Definition: etoile.C:370
double * p_ray_eq_pis2
Coordinate radius at , .
Definition: etoile.h:530
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:443
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata&#39;s prescription [Prog.
Definition: etoile_rot.C:771
double * p_ray_eq_3pis2
Coordinate radius at , .
Definition: etoile.h:536
double ray_pole() const
Coordinate radius at [r_unit].
double * p_xa_barycenter
Absolute coordinate X of the barycenter of the baryon density.
Definition: etoile.h:1012
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:466
const Tenseur & get_w_shift() const
Returns the vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:1769
Cmp ssjm1_beta
Effective source at the previous step for the resolution of the Poisson equation for beta_auto by mea...
Definition: etoile.h:971
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: etoile_bin.C:585
double ray_eq_pi() const
Coordinate radius at , [r_unit].
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition: etoile.h:490
Tenseur d_logn_auto_regu
Gradient of logn_auto_regu (Cartesian components with respect to ref_triad )
Definition: etoile.h:870
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.
const Tenseur & get_akcar_auto() const
Returns the part of the scalar generated by shift_auto , i.e.
Definition: etoile.h:1210
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition: etoile.h:855
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:463
virtual const Itbl & l_surf() const
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l o...
Definition: etoile_global.C:78
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.
Cmp ssjm1_psi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: etoile.h:995
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: etoile.C:413
virtual double area() const
Surface area.
friend ostream & operator<<(ostream &, const Etoile &)
Display.
Definition: etoile.C:509
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
const Tenseur & get_tnphi() const
Returns the component of the shift vector.
Definition: etoile.h:1733
virtual const Tbl & radius_Morsink() const
Physical radius following Morsink et al.
const Tenseur & get_ener() const
Returns the proper total energy density.
Definition: etoile.h:685
virtual void kinematics(double omega, double x_axe)
Computes the quantities bsn and pot_centri .
Definition: et_bin_kinema.C:84
Tenseur & logn
Metric potential = logn_auto.
Definition: etoile.h:1527
double * p_ray_eq_pi
Coordinate radius at , .
Definition: etoile.h:533
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: etoile_rot.C:378
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition: etoile.h:1656
void update_metric_der_comp(const Etoile_bin &comp)
Computes the derivative of metric functions related to the companion star.
Basic array class.
Definition: tbl.h:164
const Tenseur & get_psi0() const
Returns the non-translational part of the velocity potential.
Definition: etoile.h:1101
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...
virtual double mass_g() const
Gravitational mass.
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition: etoile.h:512
double * p_grv2
Error on the virial identity GRV2.
Definition: etoile.h:1639
const Tenseur & get_u_euler() const
Returns the fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:700
Tenseur & set_pot_centri()
Read/write the centrifugal potential.
Definition: etoile_bin.C:534
const Tenseur & get_shift_comp() const
Returns the part of the shift vector generated principaly by the companion star. ...
Definition: etoile.h:1164
const Tenseur & get_shift_auto() const
Returns the part of the shift vector generated principaly by the star.
Definition: etoile.h:1158
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: etoile_rot.C:405
const Tenseur & get_bbb() const
Returns the metric factor B.
Definition: etoile.h:1722
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: etoile.h:1537
Tenseur akcar_auto
Part of the scalar generated by shift_auto , i.e.
Definition: etoile.h:944
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:471
virtual double f_eccentric(double ecc, double periast, ostream *ost=0x0) const
Computation of frequency of eccentric orbits.
virtual double f_isco() const
Orbital frequency at the innermost stable circular orbit (ISCO).
Definition: et_rot_isco.C:270
Tenseur & dzeta
Metric potential = beta_auto.
Definition: etoile.h:1540
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:875
virtual double mass_b() const
Baryon mass.
Tenseur & set_khi_shift()
Read/write of khi_shift.
Definition: etoile_bin.C:548
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 )
Definition: etoile.h:507
double ray_eq_3pis2() const
Coordinate radius at , [r_unit].
const Tenseur & get_nbar() const
Returns the proper baryon density.
Definition: etoile.h:682
bool is_relativistic() const
Returns true for a relativistic star, false for a Newtonian one.
Definition: etoile.h:673
Tenseur_sym tkij_comp
Part of the extrinsic curvature tensor generated by shift_comp .
Definition: etoile.h:938
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Definition: etoile.h:1614
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Definition: etoile.h:1573
double * p_mass_b
Baryon mass.
Definition: etoile.h:553
const Eos & get_eos() const
Returns the equation of state.
Definition: etoile.h:676
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
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
const Tenseur & get_ent() const
Returns the enthalpy field.
Definition: etoile.h:679
const Tenseur & get_beta_auto() const
Returns the logarithm of the part of the product AN generated principaly by the star.
Definition: etoile.h:730
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition: etoile.h:1654
double * p_angu_mom
Angular momentum.
Definition: etoile.h:1637
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: etoile.h:979
Tenseur d_psi
Gradient of (in the irrotational case) (Cartesian components with respect to ref_triad ) ...
Definition: etoile.h:844
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: et_bin_hydro.C:109
virtual double aplat() const
Flatening r_pole/r_eq.
double * p_f_eq
Orbital frequency at the equator.
Definition: etoile.h:1657
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:545
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Definition: etoile.h:989
const Tenseur & get_ak_car() const
Returns the scalar .
Definition: etoile.h:1809