LORENE
etoile_rot.C
1 /*
2  * Methods for the class Etoile_rot
3  *
4  * (see file etoile.h for documentation)
5  */
6 
7 /*
8  * Copyright (c) 2000-2001 Eric Gourgoulhon
9  *
10  * This file is part of LORENE.
11  *
12  * LORENE is free software; you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation; either version 2 of the License, or
15  * (at your option) any later version.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 
30 
31 /*
32  * $Id: etoile_rot.C,v 1.11 2026/02/26 09:53:04 j_novak Exp $
33  * $Log: etoile_rot.C,v $
34  * Revision 1.11 2026/02/26 09:53:04 j_novak
35  * Added method to compute and store circumferential radius as function of latitude (radius_Morsink).
36  *
37  * Revision 1.10 2023/07/04 08:59:53 j_novak
38  * Added method "r_circ_merid()" to compute the circumferential meridional radius.
39  *
40  * Revision 1.9 2021/04/13 11:24:53 j_novak
41  * Corrected a bug in Etoile_rot::fait_shift() which was missing the np=1 case.
42  *
43  * Revision 1.8 2016/12/05 16:17:55 j_novak
44  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
45  *
46  * Revision 1.7 2015/12/03 14:17:24 j_novak
47  * Check added for the computation of area (thanks S. Koeppel).
48  *
49  * Revision 1.6 2015/06/10 14:37:44 a_sourie
50  * Corrected the formula for the quadrupole.
51  *
52  * Revision 1.5 2014/10/13 08:52:59 j_novak
53  * Lorene classes and functions now belong to the namespace Lorene.
54  *
55  * Revision 1.4 2004/03/25 10:29:07 j_novak
56  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
57  *
58  * Revision 1.3 2001/12/06 15:11:43 jl_zdunik
59  * Introduction of the new function f_eq() in the class Etoile_rot
60  *
61  * Revision 1.2 2001/12/04 21:27:53 e_gourgoulhon
62  *
63  * All writing/reading to a binary file are now performed according to
64  * the big endian convention, whatever the system is big endian or
65  * small endian, thanks to the functions fwrite_be and fread_be
66  *
67  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
68  * LORENE
69  *
70  * Revision 2.17 2001/10/24 15:36:20 eric
71  * Ajout de la fonction display_poly.
72  *
73  * Revision 2.16 2001/10/16 14:49:02 eric
74  * Appel de get_omega_c() pour avoir la valeur centrale de Omega.
75  * Affichage different si rotation differentielle.
76  *
77  * Revision 2.15 2001/09/13 08:32:01 eric
78  * Ajout du facteur de compacite M/R dans l'affichage.
79  *
80  * Revision 2.14 2001/06/20 14:20:56 novak
81  * Appel a Etoile_rot::set_der0x0 dans del_deriv (au lieu de set_der0x0
82  * tout court).
83  *
84  * Revision 2.13 2001/03/26 09:30:58 jlz
85  * New members p_espec_isco and p_lspec_isco.
86  *
87  * Revision 2.12 2000/11/20 21:42:02 eric
88  * Appel de fait_nphi() dans le constructeur par lecture de fichier.
89  *
90  * Revision 2.11 2000/11/18 23:18:30 eric
91  * Modifs affichage.
92  *
93  * Revision 2.10 2000/11/18 21:09:57 eric
94  * Ajout des membres p_r_isco et p_f_isco.
95  *
96  * Revision 2.9 2000/11/07 16:33:08 eric
97  * Modif affichage.
98  *
99  * Revision 2.8 2000/10/12 15:37:01 eric
100  * Ajout de la fonction fait_nphi().
101  *
102  * Revision 2.7 2000/09/18 16:15:12 eric
103  * Ajout du membre tkij.
104  *
105  * Revision 2.6 2000/08/31 15:38:00 eric
106  * Bases spectrales standards pour bbb et b_car dans le constructeur
107  * standard (initialisation a la metrique plate).
108  *
109  * Revision 2.5 2000/08/31 11:25:45 eric
110  * Ajout des membres tnphi et ak_car.
111  *
112  * Revision 2.4 2000/08/25 12:28:29 eric
113  * Modif affichage.
114  *
115  * Revision 2.3 2000/08/18 14:01:59 eric
116  * Ajout de partial_display
117  *
118  * Revision 2.2 2000/08/17 12:40:04 eric
119  * *** empty log message ***
120  *
121  * Revision 2.1 2000/07/21 16:31:26 eric
122  * *** empty log message ***
123  *
124  * Revision 1.1 2000/07/20 15:32:37 eric
125  * Initial revision
126  *
127  *
128  * $Header: /cvsroot/Lorene/C++/Source/Etoile/etoile_rot.C,v 1.11 2026/02/26 09:53:04 j_novak Exp $
129  *
130  */
131 
132 // Headers C
133 #include "math.h"
134 
135 // Headers Lorene
136 #include "etoile.h"
137 #include "eos.h"
138 #include "nbr_spx.h"
139 #include "utilitaires.h"
140 #include "unites.h"
141 
142  //--------------//
143  // Constructors //
144  //--------------//
145 
146 // Standard constructor
147 // --------------------
148 namespace Lorene {
149 Etoile_rot::Etoile_rot(Map& mpi, int nzet_i, bool relat, const Eos& eos_i)
150  : Etoile(mpi, nzet_i, relat, eos_i),
151  bbb(mpi),
152  b_car(mpi),
153  nphi(mpi),
154  tnphi(mpi),
155  uuu(mpi),
156  logn(logn_auto),
157  nuf(mpi),
158  nuq(mpi),
159  dzeta(beta_auto),
160  tggg(mpi),
161  w_shift(mpi, 1, CON, mp.get_bvect_cart()),
162  khi_shift(mpi),
163  tkij(mpi, 2, COV, mp.get_bvect_cart()),
164  ak_car(mpi),
165  ssjm1_nuf(mpi),
166  ssjm1_nuq(mpi),
167  ssjm1_dzeta(mpi),
168  ssjm1_tggg(mpi),
169  ssjm1_khi(mpi),
170  ssjm1_wshift(mpi, 1, CON, mp.get_bvect_cart())
171 {
172 
173  // Initialization to a static state :
174  omega = 0 ;
175  uuu = 0 ;
176 
177  // Initialization to a flat metric :
178  bbb = 1 ;
179  bbb.set_std_base() ;
180  b_car = 1 ;
181  b_car.set_std_base() ;
182  nphi = 0 ;
183  tnphi = 0 ;
184  nuf = 0 ;
185  nuq = 0 ;
186  tggg = 0 ;
187 
188  w_shift.set_etat_qcq() ;
189  for (int i=0; i<3; i++) {
190  w_shift.set(i) = 0 ;
191  }
192 
194  khi_shift.set() = 0 ;
195 
196  tkij.set_etat_zero() ;
197 
198  ak_car = 0 ;
199 
200  ssjm1_nuf = 0 ;
201  ssjm1_nuq = 0 ;
202  ssjm1_dzeta = 0 ;
203  ssjm1_tggg = 0 ;
204  ssjm1_khi = 0 ;
205 
207  for (int i=0; i<3; i++) {
208  ssjm1_wshift.set(i) = 0 ;
209  }
210 
211  // Pointers of derived quantities initialized to zero :
212  set_der_0x0() ;
213 
214 }
215 
216 // Copy constructor
217 // ----------------
218 
220  : Etoile(et),
221  bbb(et.bbb),
222  b_car(et.b_car),
223  nphi(et.nphi),
224  tnphi(et.tnphi),
225  uuu(et.uuu),
226  logn(logn_auto),
227  nuf(et.nuf),
228  nuq(et.nuq),
229  dzeta(beta_auto),
230  tggg(et.tggg),
231  w_shift(et.w_shift),
232  khi_shift(et.khi_shift),
233  tkij(et.tkij),
234  ak_car(et.ak_car),
235  ssjm1_nuf(et.ssjm1_nuf),
236  ssjm1_nuq(et.ssjm1_nuq),
237  ssjm1_dzeta(et.ssjm1_dzeta),
238  ssjm1_tggg(et.ssjm1_tggg),
239  ssjm1_khi(et.ssjm1_khi),
240  ssjm1_wshift(et.ssjm1_wshift)
241 {
242  omega = et.omega ;
243 
244  // Pointers of derived quantities initialized to zero :
245  set_der_0x0() ;
246 }
247 
248 
249 // Constructor from a file
250 // -----------------------
251 Etoile_rot::Etoile_rot(Map& mpi, const Eos& eos_i, FILE* fich)
252  : Etoile(mpi, eos_i, fich),
253  bbb(mpi),
254  b_car(mpi),
255  nphi(mpi),
256  tnphi(mpi),
257  uuu(mpi),
258  logn(logn_auto),
259  nuf(mpi),
260  nuq(mpi),
261  dzeta(beta_auto),
262  tggg(mpi),
263  w_shift(mpi, 1, CON, mp.get_bvect_cart()),
264  khi_shift(mpi),
265  tkij(mpi, 2, COV, mp.get_bvect_cart()),
266  ak_car(mpi),
267  ssjm1_nuf(mpi),
268  ssjm1_nuq(mpi),
269  ssjm1_dzeta(mpi),
270  ssjm1_tggg(mpi),
271  ssjm1_khi(mpi),
272  ssjm1_wshift(mpi, 1, CON, mp.get_bvect_cart())
273 {
274 
275  // Etoile parameters
276  // -----------------
277 
278  // omega is read in the file:
279  fread_be(&omega, sizeof(double), 1, fich) ;
280 
281 
282  // Read of the saved fields:
283  // ------------------------
284 
285  Tenseur nuf_file(mp, fich) ;
286  nuf = nuf_file ;
287 
288  Tenseur nuq_file(mp, fich) ;
289  nuq = nuq_file ;
290 
291  Tenseur tggg_file(mp, fich) ;
292  tggg = tggg_file ;
293 
294  Tenseur w_shift_file(mp, mp.get_bvect_cart(), fich) ;
295  w_shift = w_shift_file ;
296 
297  Tenseur khi_shift_file(mp, fich) ;
298  khi_shift = khi_shift_file ;
299 
300  fait_shift() ; // constructs shift from w_shift and khi_shift
301  fait_nphi() ; // constructs N^phi from (N^x,N^y,N^z)
302 
303  Cmp ssjm1_nuf_file(mp, *(mp.get_mg()), fich) ;
304  ssjm1_nuf = ssjm1_nuf_file ;
305 
306  Cmp ssjm1_nuq_file(mp, *(mp.get_mg()), fich) ;
307  ssjm1_nuq = ssjm1_nuq_file ;
308 
309  Cmp ssjm1_dzeta_file(mp, *(mp.get_mg()), fich) ;
310  ssjm1_dzeta = ssjm1_dzeta_file ;
311 
312  Cmp ssjm1_tggg_file(mp, *(mp.get_mg()), fich) ;
313  ssjm1_tggg = ssjm1_tggg_file ;
314 
315  Cmp ssjm1_khi_file(mp, *(mp.get_mg()), fich) ;
316  ssjm1_khi = ssjm1_khi_file ;
317 
318  Tenseur ssjm1_wshift_file(mp, mp.get_bvect_cart(), fich) ;
319  ssjm1_wshift = ssjm1_wshift_file ;
320 
321  // All other fields are initialized to zero :
322  // ----------------------------------------
323  bbb = 0 ;
324  b_car = 0 ;
325  uuu = 0 ;
326 
327  // Pointers of derived quantities initialized to zero
328  // --------------------------------------------------
329  set_der_0x0() ;
330 
331 }
332 
333  //------------//
334  // Destructor //
335  //------------//
336 
338 
339  del_deriv() ;
340 
341 }
342 
343  //----------------------------------//
344  // Management of derived quantities //
345  //----------------------------------//
346 
347 void Etoile_rot::del_deriv() const {
348 
349  Etoile::del_deriv() ;
350 
351  if (p_angu_mom != 0x0) delete p_angu_mom ;
352  if (p_tsw != 0x0) delete p_tsw ;
353  if (p_grv2 != 0x0) delete p_grv2 ;
354  if (p_grv3 != 0x0) delete p_grv3 ;
355  if (p_r_circ != 0x0) delete p_r_circ ;
356  if (p_r_circ_merid != 0x0) delete p_r_circ_merid ;
357  if (p_area != 0x0) delete p_area ;
358  if (p_aplat != 0x0) delete p_aplat ;
359  if (p_z_eqf != 0x0) delete p_z_eqf ;
360  if (p_z_eqb != 0x0) delete p_z_eqb ;
361  if (p_z_pole != 0x0) delete p_z_pole ;
362  if (p_mom_quad != 0x0) delete p_mom_quad ;
363  if (p_mom_quad_old != 0x0) delete p_mom_quad_old ;
364  if (p_mom_quad_Bo != 0x0) delete p_mom_quad_Bo ;
365  if (p_radius_Morsink != 0x0) delete p_radius_Morsink ;
366  if (p_r_isco != 0x0) delete p_r_isco ;
367  if (p_f_isco != 0x0) delete p_f_isco ;
368  if (p_lspec_isco != 0x0) delete p_lspec_isco ;
369  if (p_espec_isco != 0x0) delete p_espec_isco ;
370  if (p_f_eq != 0x0) delete p_f_eq ;
371 
373 }
374 
375 
376 
377 
379 
381 
382  p_angu_mom = 0x0 ;
383  p_tsw = 0x0 ;
384  p_grv2 = 0x0 ;
385  p_grv3 = 0x0 ;
386  p_r_circ = 0x0 ;
387  p_r_circ_merid = 0x0 ;
388  p_area = 0x0 ;
389  p_aplat = 0x0 ;
390  p_z_eqf = 0x0 ;
391  p_z_eqb = 0x0 ;
392  p_z_pole = 0x0 ;
393  p_mom_quad = 0x0 ;
394  p_mom_quad_old = 0x0 ;
395  p_mom_quad_Bo = 0x0 ;
396  p_radius_Morsink = 0x0 ;
397  p_r_isco = 0x0 ;
398  p_f_isco = 0x0 ;
399  p_lspec_isco = 0x0 ;
400  p_espec_isco = 0x0 ;
401  p_f_eq = 0x0 ;
402 
403 }
404 
406 
408 
409  del_deriv() ;
410 
411 }
412 
413 
414  //--------------//
415  // Assignment //
416  //--------------//
417 
418 // Assignment to another Etoile_rot
419 // --------------------------------
421 
422  // Assignment of quantities common to all the derived classes of Etoile
423  Etoile::operator=(et) ;
424 
425  // Assignement of proper quantities of class Etoile_rot
426  omega = et.omega ;
427 
428  bbb = et.bbb ;
429  b_car = et.b_car ;
430  nphi = et.nphi ;
431  tnphi = et.tnphi ;
432  uuu = et.uuu ;
433  nuf = et.nuf ;
434  nuq = et.nuq ;
435  tggg = et.tggg ;
436  w_shift = et.w_shift ;
437  khi_shift = et.khi_shift ;
438  tkij = et.tkij ;
439  ak_car = et.ak_car ;
440  ssjm1_nuf = et.ssjm1_nuf ;
441  ssjm1_nuq = et.ssjm1_nuq ;
442  ssjm1_dzeta = et.ssjm1_dzeta ;
443  ssjm1_tggg = et.ssjm1_tggg ;
444  ssjm1_khi = et.ssjm1_khi ;
445  ssjm1_wshift = et.ssjm1_wshift ;
446 
447  del_deriv() ; // Deletes all derived quantities
448 
449 }
450 
451  //--------------//
452  // Outputs //
453  //--------------//
454 
455 // Save in a file
456 // --------------
457 void Etoile_rot::sauve(FILE* fich) const {
458 
459  Etoile::sauve(fich) ;
460 
461  fwrite_be(&omega, sizeof(double), 1, fich) ;
462 
463  nuf.sauve(fich) ;
464  nuq.sauve(fich) ;
465  tggg.sauve(fich) ;
466  w_shift.sauve(fich) ;
467  khi_shift.sauve(fich) ;
468 
469  ssjm1_nuf.sauve(fich) ;
470  ssjm1_nuq.sauve(fich) ;
471  ssjm1_dzeta.sauve(fich) ;
472  ssjm1_tggg.sauve(fich) ;
473  ssjm1_khi.sauve(fich) ;
474  ssjm1_wshift.sauve(fich) ;
475 
476 
477 }
478 
479 // Printing
480 // --------
481 
482 ostream& Etoile_rot::operator>>(ostream& ost) const {
483 
484  using namespace Unites ;
485 
486  Etoile::operator>>(ost) ;
487 
488  double omega_c = get_omega_c() ;
489 
490  ost << endl ;
491  if (omega != __infinity) {
492  ost << "Uniformly rotating star" << endl ;
493  ost << "-----------------------" << endl ;
494 
495  double freq = omega / (2.*M_PI) ;
496  ost << "Omega : " << omega * f_unit
497  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
498  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
499  << endl ;
500 
501  }
502  else {
503  ost << "Differentially rotating star" << endl ;
504  ost << "----------------------------" << endl ;
505 
506  double freq = omega_c / (2.*M_PI) ;
507  ost << "Central value of Omega : " << omega_c * f_unit
508  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
509  ost << "Central rotation period : " << 1000. / (freq * f_unit) << " ms"
510  << endl ;
511 
512  }
513 
514 
515  double nphi_c = nphi()(0, 0, 0, 0) ;
516  if ( (omega_c==0) && (nphi_c==0) ) {
517  ost << "Central N^phi : " << nphi_c << endl ;
518  }
519  else{
520  ost << "Central N^phi/Omega : " << nphi_c / omega_c << endl ;
521  }
522 
523  ost << "Error on the virial identity GRV2 : " << endl ;
524  ost << "GRV2 = " << grv2() << endl ;
525  ost << "Error on the virial identity GRV3 : " << endl ;
526  double xgrv3 = grv3(&ost) ;
527  ost << "GRV3 = " << xgrv3 << endl ;
528 
529  double mom_quad_38si = mom_quad() * rho_unit * (pow(r_unit, double(5.))
530  / double(1.e38) ) ;
531  ost << "Quadrupole moment Q : " << mom_quad_38si << " 10^38 kg m^2"
532  << endl ;
533  ost << "Q / (M R_circ^2) : "
534  << mom_quad() / ( mass_g() * pow( r_circ(), 2. ) ) << endl ;
535  ost << "c^4 Q / (G^2 M^3) : "
536  << mom_quad() / ( pow(qpig/(4*M_PI), 2.) * pow(mass_g(), 3.) )
537  << endl ;
538 
539  ost << "Angular momentum J : "
540  << angu_mom()/( qpig / (4* M_PI) * msol*msol) << " G M_sol^2 / c"
541  << endl ;
542  ost << "c J / (G M^2) : "
543  << angu_mom()/( qpig / (4* M_PI) * pow(mass_g(), 2.) ) << endl ;
544 
545  if (omega != __infinity) {
546  double mom_iner = angu_mom() / omega ;
547  double mom_iner_38si = mom_iner * rho_unit * (pow(r_unit, double(5.))
548  / double(1.e38) ) ;
549  ost << "Moment of inertia: " << mom_iner_38si << " 10^38 kg m^2"
550  << endl ;
551  }
552 
553  ost << "Ratio T/W : " << tsw() << endl ;
554  ost << "Circumferential equatorial radius R_circ_eq : "
555  << r_circ()/km << " km" << endl ;
556  if (mp.get_mg()->get_np(0) == 1) {
557  ost << "Circumferential meridional radius R_circ_merid : "
558  << r_circ_merid()/km << " km" << endl ;
559  ost << "Flattening r_circ_merid/r_circ_eq : "
560  << r_circ_merid() / r_circ() << endl ;
561  ost << "Surface area : " << area()/(km*km) << " km^2" << endl ;
562  ost << "Mean radius R_mean : "
563  << mean_radius()/km << " km" << endl ;
564  } else {
565  ost <<
566  "Skipping polar radius / surface statements due to number of points in phi direction np > 1"
567  << endl;
568  }
569  ost << "Coordinate equatorial radius r_eq : " << ray_eq()/km << " km"
570  << endl ;
571  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
572 
573  double compact = qpig/(4.*M_PI) * mass_g() / r_circ() ;
574  ost << "Compaction parameter M_g / R_circ : " << compact << endl ;
575 
576  int lsurf = nzet - 1;
577  int nt = mp.get_mg()->get_nt(lsurf) ;
578  int nr = mp.get_mg()->get_nr(lsurf) ;
579  ost << "Equatorial value of the velocity U: "
580  << uuu()(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
581  ost << "Redshift at the equator (forward) : " << z_eqf() << endl ;
582  ost << "Redshift at the equator (backward): " << z_eqb() << endl ;
583  ost << "Redshift at the pole : " << z_pole() << endl ;
584 
585 
586  ost << "Central value of log(N) : "
587  << logn()(0, 0, 0, 0) << endl ;
588 
589  ost << "Central value of dzeta=log(AN) : "
590  << dzeta()(0, 0, 0, 0) << endl ;
591 
592  if ( (omega_c==0) && (nphi_c==0) ) {
593  ost << "Central N^phi : " << nphi_c << endl ;
594  }
595  else{
596  ost << "Central N^phi/Omega : " << nphi_c / omega_c << endl ;
597  }
598 
599  ost << " ... w_shift (NB: components in the star Cartesian frame) [c] : "
600  << endl
601  << w_shift(0)(0, 0, 0, 0) << " "
602  << w_shift(1)(0, 0, 0, 0) << " "
603  << w_shift(2)(0, 0, 0, 0) << endl ;
604 
605  ost << "Central value of khi_shift [km c] : "
606  << khi_shift()(0, 0, 0, 0) / km << endl ;
607 
608  ost << "Central value of B^2 : " << b_car()(0,0,0,0) << endl ;
609 
610  Tbl diff_a_b = diffrel( a_car(), b_car() ) ;
611  ost <<
612  "Relative discrepancy in each domain between the metric coef. A^2 and B^2 : "
613  << endl ;
614  for (int l=0; l<diff_a_b.get_taille(); l++) {
615  ost << diff_a_b(l) << " " ;
616  }
617  ost << endl ;
618 
619  // Approximate formula for R_isco = 6 R_g (1-(2/3)^1.5 j )
620  // up to the first order in j
621  double jdimless = angu_mom() / ( ggrav * pow(mass_g(), 2.) ) ;
622  double r_grav = ggrav * mass_g() ;
623  double r_isco_appr = 6. * r_grav * ( 1. - pow(2./3.,1.5) * jdimless ) ;
624 
625  // Approximate formula for the ISCO frequency
626  // freq_ms = 6^{-1.5}/2pi/R_g (1+11*6^(-1.5) j )
627  // up to the first order in j
628  double f_isco_appr = ( 1. + 11. /6. /sqrt(6.) * jdimless ) / r_grav /
629  (12. * M_PI ) / sqrt(6.) ;
630 
631  ost << endl << "Innermost stable circular orbit (ISCO) : " << endl ;
632  double xr_isco = r_isco(&ost) ;
633  ost <<" circumferential radius r_isco = "
634  << xr_isco / km << " km" << endl ;
635  ost <<" (approx. 6M + 1st order in j : "
636  << r_isco_appr / km << " km)" << endl ;
637  ost <<" (approx. 6M : "
638  << 6. * r_grav / km << " km)" << endl ;
639  ost <<" orbital frequency f_isco = "
640  << f_isco() * f_unit << " Hz" << endl ;
641  ost <<" (approx. 1st order in j : "
642  << f_isco_appr * f_unit << " Hz)" << endl ;
643 
644 
645  return ost ;
646 
647 }
648 
649 
650 void Etoile_rot::partial_display(ostream& ost) const {
651 
652  using namespace Unites ;
653 
654  double omega_c = get_omega_c() ;
655  double freq = omega_c / (2.*M_PI) ;
656  ost << "Central Omega : " << omega_c * f_unit
657  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
658  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
659  << endl ;
660  ost << endl << "Central enthalpy : " << ent()(0,0,0,0) << " c^2" << endl ;
661  ost << "Central proper baryon density : " << nbar()(0,0,0,0)
662  << " x 0.1 fm^-3" << endl ;
663  ost << "Central proper energy density : " << ener()(0,0,0,0)
664  << " rho_nuc c^2" << endl ;
665  ost << "Central pressure : " << press()(0,0,0,0)
666  << " rho_nuc c^2" << endl ;
667 
668  ost << "Central value of log(N) : "
669  << logn()(0, 0, 0, 0) << endl ;
670  ost << "Central lapse N : " << nnn()(0,0,0,0) << endl ;
671  ost << "Central value of dzeta=log(AN) : "
672  << dzeta()(0, 0, 0, 0) << endl ;
673  ost << "Central value of A^2 : " << a_car()(0,0,0,0) << endl ;
674  ost << "Central value of B^2 : " << b_car()(0,0,0,0) << endl ;
675 
676  double nphi_c = nphi()(0, 0, 0, 0) ;
677  if ( (omega_c==0) && (nphi_c==0) ) {
678  ost << "Central N^phi : " << nphi_c << endl ;
679  }
680  else{
681  ost << "Central N^phi/Omega : " << nphi_c / omega_c
682  << endl ;
683  }
684 
685 
686  int lsurf = nzet - 1;
687  int nt = mp.get_mg()->get_nt(lsurf) ;
688  int nr = mp.get_mg()->get_nr(lsurf) ;
689  ost << "Equatorial value of the velocity U: "
690  << uuu()(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
691 
692  ost << endl
693  << "Coordinate equatorial radius r_eq = "
694  << ray_eq()/km << " km" << endl ;
695  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
696  if (mp.get_mg()->get_np(0) == 1)
697  ost << "Flattening r_circ_pole/r_circ_eq : "
698  << r_circ_merid() / r_circ() << endl ;
699 
700 }
701 
702 
703 double Etoile_rot::get_omega_c() const {
704 
705  return omega ;
706 
707 }
708 
709 
710 // display_poly
711 // ------------
712 
713 void Etoile_rot::display_poly(ostream& ost) const {
714 
715  using namespace Unites ;
716 
717  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
718 
719  if (p_eos_poly != 0x0) {
720 
721  double kappa = p_eos_poly->get_kap() ;
722  double gamma = p_eos_poly->get_gam() ; ;
723 
724  // kappa^{n/2}
725  double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
726 
727  // Polytropic unit of length in terms of r_unit :
728  double r_poly = kap_ns2 / sqrt(ggrav) ;
729 
730  // Polytropic unit of time in terms of t_unit :
731  double t_poly = r_poly ;
732 
733  // Polytropic unit of mass in terms of m_unit :
734  double m_poly = r_poly / ggrav ;
735 
736  // Polytropic unit of angular momentum in terms of j_unit :
737  double j_poly = r_poly * r_poly / ggrav ;
738 
739  // Polytropic unit of density in terms of rho_unit :
740  double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
741 
742  ost.precision(10) ;
743  ost << endl << "Quantities in polytropic units : " << endl ;
744  ost << "==============================" << endl ;
745  ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
746  ost << " n_c : " << nbar()(0, 0, 0, 0) / rho_poly << endl ;
747  ost << " e_c : " << ener()(0, 0, 0, 0) / rho_poly << endl ;
748  ost << " Omega_c : " << get_omega_c() * t_poly << endl ;
749  ost << " P_c : " << 2.*M_PI / get_omega_c() / t_poly << endl ;
750  ost << " M_bar : " << mass_b() / m_poly << endl ;
751  ost << " M : " << mass_g() / m_poly << endl ;
752  ost << " J : " << angu_mom() / j_poly << endl ;
753  ost << " r_eq : " << ray_eq() / r_poly << endl ;
754  ost << " R_circ : " << r_circ() / r_poly << endl ;
755 
756 
757  }
758 
759 
760 }
761 
762 
763 
764 
765 
766 
767  //-------------------------//
768  // Computational routines //
769  //-------------------------//
770 
772 
773  Tenseur d_khi = khi_shift.gradient() ;
774 
775  if (d_khi.get_etat() == ETATQCQ) {
776  d_khi.dec2_dzpuis() ; // divide by r^2 in the external compactified
777  // domain
778  }
779 
780  // x_k dW^k/dx_i
781 
782  Tenseur x_d_w = skxk( w_shift.gradient() ) ;
783  x_d_w.dec_dzpuis() ;
784 
785  double lambda = double(1) / double(3) ;
786 
787  if ( mp.get_mg()->get_np(0) == 1 ) {
788  lambda = 0 ;
789  }
790 
791  // The final computation is done component by component because
792  // d_khi and x_d_w are covariant comp. whereas w_shift is
793  // contravariant
794 
795  shift.set_etat_qcq() ;
796 
797  for (int i=0; i<3; i++) {
798  shift.set(i) = (lambda+2)/2./(lambda+1) * w_shift(i)
799  - (lambda/2./(lambda+1)) * (d_khi(i) + x_d_w(i)) ;
800  }
801 
802  shift.set_triad( *(w_shift.get_triad()) ) ;
803 
804 }
805 
806 
807 
809 
810  if ( shift.get_etat() == ETATZERO ) {
811  tnphi = 0 ;
812  nphi = 0 ;
813  return ;
814  }
815 
816  assert( shift.get_etat() == ETATQCQ ) ;
817 
818  // Computation of tnphi
819  // --------------------
820  tnphi.set_etat_qcq() ;
821 
822  mp.comp_p_from_cartesian( shift(0), shift(1), tnphi.set() ) ;
823 
824  // Computation of nphi
825  // -------------------
826 
827  nphi = tnphi ;
828  (nphi.set()).div_rsint() ;
829 
830 }
831 }
double * p_z_eqb
Backward redshift factor at equator.
Definition: etoile.h:1646
virtual double z_eqf() const
Forward redshift factor at equator.
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 r_circ() const
Circumferential equatorial radius.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: etoile.C:514
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
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:449
double * p_r_isco
Circumferential radius of the ISCO.
Definition: etoile.h:1651
Tenseur tnphi
Component of the shift vector.
Definition: etoile.h:1521
void dec2_dzpuis()
dzpuis -= 2 ;
Definition: tenseur.C:1146
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tenseur.C:690
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: etoile_rot.C:713
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
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
void operator=(const Etoile &)
Assignment to another Etoile.
Definition: etoile.C:433
Lorene prototypes.
Definition: app_hor.h:67
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: etoile.h:1556
Standard units of space, time and mass.
Equation of state base class.
Definition: eos.h:209
virtual void comp_p_from_cartesian(const Scalar &v_x, const Scalar &v_y, Scalar &v_p) const =0
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian c...
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: etoile.h:1631
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
Tenseur nnn
Total lapse function.
Definition: etoile.h:515
double * p_z_eqf
Forward redshift factor at equator.
Definition: etoile.h:1645
Tenseur nphi
Metric coefficient .
Definition: etoile.h:1516
double ray_eq() const
Coordinate radius at , [r_unit].
Base class for coordinate mappings.
Definition: map.h:696
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
void operator=(const Etoile_rot &)
Assignment to another Etoile_rot.
Definition: etoile_rot.C:420
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
Definition: et_rot_isco.C:87
Tenseur press
Fluid pressure.
Definition: etoile.h:467
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.
Tenseur shift
Total shift vector.
Definition: etoile.h:518
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
double * p_aplat
Flatening r_pole/r_eq.
Definition: etoile.h:1644
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Definition: eos_poly.C:271
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_rot.C:347
virtual double r_circ_merid() const
Circumferential meridional radius.
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile.C:381
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
Definition: etoile_rot.C:703
Cmp ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta ...
Definition: etoile.h:1609
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
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
virtual double angu_mom() const
Angular momentum.
void sauve(FILE *) const
Save in a file.
Definition: tenseur.C:1341
Tenseur skxk(const Tenseur &)
Contraction of the last index of (*this) with or , depending on the type of S .
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
virtual double z_pole() const
Redshift factor at North pole.
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:465
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tenseur.h:710
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
virtual void sauve(FILE *) const
Save in a file.
Definition: etoile_rot.C:457
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.
double * p_area
Surface area.
Definition: etoile.h:1643
Polytropic equation of state (relativistic case).
Definition: eos.h:812
Map & mp
Mapping associated with the star.
Definition: etoile.h:435
const Eos & eos
Equation of state of the stellar matter.
Definition: etoile.h:457
virtual double z_eqb() const
Backward redshift factor at equator.
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
int get_etat() const
Returns the logical state.
Definition: tenseur.h:713
double * p_f_isco
Orbital frequency of the ISCO.
Definition: etoile.h:1652
virtual double tsw() const
Ratio T/W.
Tenseur bbb
Metric factor B.
Definition: etoile.h:1510
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition: fwrite_be.C:73
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tenseur uuu
Norm of u_euler.
Definition: etoile.h:1524
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
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:438
virtual double mom_quad() const
Quadrupole moment.
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:72
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
Tenseur a_car
Total conformal factor .
Definition: etoile.h:521
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata&#39;s prescription [Prog.
Definition: etoile_rot.C:771
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:466
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:817
void dec_dzpuis()
dzpuis -= 1 ;
Definition: tenseur.C:1120
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:463
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.
Tenseur & logn
Metric potential = logn_auto.
Definition: etoile.h:1527
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition: etoile.h:1656
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: etoile_rot.C:378
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:476
virtual double mass_g() const
Gravitational mass.
double * p_grv2
Error on the virial identity GRV2.
Definition: etoile.h:1639
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
void sauve(FILE *) const
Save in a file.
Definition: cmp.C:564
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: etoile_rot.C:405
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: etoile.h:1537
void set_etat_zero()
Sets the logical state to ETATZERO (zero state).
Definition: tenseur.C:661
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
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
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:307
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition: etoile.h:1654
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates)
Definition: tenseur.C:1558
double * p_angu_mom
Angular momentum.
Definition: etoile.h:1637
virtual double aplat() const
Flatening r_pole/r_eq.
double * p_f_eq
Orbital frequency at the equator.
Definition: etoile.h:1657