LORENE
star_rot.C
1 /*
2  * Methods of class Star_rot
3  *
4  * (see file star_rot.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2010 Eric Gourgoulhon
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License version 2
15  * as published by the Free Software Foundation.
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  * $Id: star_rot.C,v 1.16 2026/02/05 14:23:25 j_novak Exp $
32  * $Log: star_rot.C,v $
33  * Revision 1.16 2026/02/05 14:23:25 j_novak
34  * New method surf_grav_euler().
35  *
36  * Revision 1.15 2026/01/20 08:35:28 j_novak
37  * New method radius_Morsink() to compute the physical radius along the definition of Morsink et al. (2007) in the oblate Schwarschild approximation.
38  *
39  * Revision 1.14 2023/07/04 09:01:09 j_novak
40  * Changed the name r_pol_circ -> r_circ_merid to explicitly show the radius is a circumferential one, along a meridian.
41  *
42  * Revision 1.13 2023/07/03 14:32:31 j_novak
43  * Added method "r_pol_circ()", to compute polar circumferential radius, with a circumference along a meridian.
44  *
45  * Revision 1.12 2021/10/08 13:10:14 j_novak
46  * Corrected a bug in Star_rot::fait_shift) in the np=1 case
47  *
48  * Revision 1.11 2017/04/11 10:46:55 m_bejger
49  * Star_rot::surf_grav() - surface gravity values along the theta direction
50  *
51  * Revision 1.10 2016/12/05 16:18:15 j_novak
52  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
53  *
54  * Revision 1.9 2015/05/19 09:30:56 j_novak
55  * New methods for computing the area of the star and its mean radius.
56  *
57  * Revision 1.8 2014/10/13 08:53:39 j_novak
58  * Lorene classes and functions now belong to the namespace Lorene.
59  *
60  * Revision 1.7 2014/10/06 15:13:17 j_novak
61  * Modified #include directives to use c++ syntax.
62  *
63  * Revision 1.6 2010/02/08 11:45:58 j_novak
64  * Better computation of fait_shift()
65  *
66  * Revision 1.5 2010/02/08 10:56:30 j_novak
67  * Added a few things missing for the reading from resulting file.
68  *
69  * Revision 1.4 2010/02/02 12:45:16 e_gourgoulhon
70  * Improved the display (operator>>)
71  *
72  * Revision 1.3 2010/01/25 22:33:35 e_gourgoulhon
73  * Debugging...
74  *
75  * Revision 1.2 2010/01/25 18:15:32 e_gourgoulhon
76  * Added member unsurc2
77  *
78  * Revision 1.1 2010/01/24 16:09:39 e_gourgoulhon
79  * New class Star_rot.
80  *
81  *
82  * $Header: /cvsroot/Lorene/C++/Source/Star/star_rot.C,v 1.16 2026/02/05 14:23:25 j_novak Exp $
83  *
84  */
85 
86 
87 // Lorene headers
88 #include "star_rot.h"
89 #include "eos.h"
90 #include "unites.h"
91 #include "utilitaires.h"
92 #include "nbr_spx.h"
93 
94 
95 
96  //--------------//
97  // Constructors //
98  //--------------//
99 
100 // Standard constructor
101 // --------------------
102 namespace Lorene {
103 Star_rot::Star_rot(Map& mpi, int nzet_i, bool relat, const Eos& eos_i)
104  : Star(mpi, nzet_i, eos_i),
105  relativistic(relat),
106  a_car(mpi),
107  bbb(mpi),
108  b_car(mpi),
109  nphi(mpi),
110  tnphi(mpi),
111  uuu(mpi),
112  nuf(mpi),
113  nuq(mpi),
114  dzeta(mpi),
115  tggg(mpi),
116  w_shift(mpi, CON, mp.get_bvect_cart()),
117  khi_shift(mpi),
118  tkij(mpi, COV, mp.get_bvect_cart()),
119  ak_car(mpi),
120  ssjm1_nuf(mpi),
121  ssjm1_nuq(mpi),
122  ssjm1_dzeta(mpi),
123  ssjm1_tggg(mpi),
124  ssjm1_khi(mpi),
125  ssjm1_wshift(mpi, CON, mp.get_bvect_cart())
126 {
127 
128  // Parameter 1/c^2 is deduced from relativistic:
129  unsurc2 = relativistic ? double(1) : double(0) ;
130 
131  // Initialization to a static state :
132  omega = 0 ;
133  uuu = 0 ;
134 
135  // Initialization to a flat metric :
136  a_car = 1 ;
137  bbb = 1 ;
139  b_car = 1 ;
140  nphi = 0 ;
141  tnphi = 0 ;
142  nuf = 0 ;
143  nuq = 0 ;
144  dzeta = 0 ;
145  tggg = 0 ;
146 
148  khi_shift = 0 ;
149 
150  beta.set_etat_zero() ;
152 
153  tkij.set_etat_zero() ;
154 
155  ak_car = 0 ;
156 
157  ssjm1_nuf = 0 ;
158  ssjm1_nuq = 0 ;
159  ssjm1_dzeta = 0 ;
160  ssjm1_tggg = 0 ;
161  ssjm1_khi = 0 ;
162 
164 
165  // Pointers of derived quantities initialized to zero :
166  set_der_0x0() ;
167 
168 }
169 
170 // Copy constructor
171 // ----------------
172 
174  : Star(et),
175  relativistic(et.relativistic),
176  unsurc2(et.unsurc2),
177  omega(et.omega),
178  a_car(et.a_car),
179  bbb(et.bbb),
180  b_car(et.b_car),
181  nphi(et.nphi),
182  tnphi(et.tnphi),
183  uuu(et.uuu),
184  nuf(et.nuf),
185  nuq(et.nuq),
186  dzeta(et.dzeta),
187  tggg(et.tggg),
188  w_shift(et.w_shift),
189  khi_shift(et.khi_shift),
190  tkij(et.tkij),
191  ak_car(et.ak_car),
192  ssjm1_nuf(et.ssjm1_nuf),
193  ssjm1_nuq(et.ssjm1_nuq),
194  ssjm1_dzeta(et.ssjm1_dzeta),
195  ssjm1_tggg(et.ssjm1_tggg),
196  ssjm1_khi(et.ssjm1_khi),
197  ssjm1_wshift(et.ssjm1_wshift)
198 {
199  // Pointers of derived quantities initialized to zero :
200  set_der_0x0() ;
201 }
202 
203 
204 // Constructor from a file
205 // -----------------------
206 Star_rot::Star_rot(Map& mpi, const Eos& eos_i, FILE* fich)
207  : Star(mpi, eos_i, fich),
208  a_car(mpi),
209  bbb(mpi),
210  b_car(mpi),
211  nphi(mpi),
212  tnphi(mpi),
213  uuu(mpi),
214  nuf(mpi),
215  nuq(mpi),
216  dzeta(mpi),
217  tggg(mpi),
218  w_shift(mpi, CON, mp.get_bvect_cart()),
219  khi_shift(mpi),
220  tkij(mpi, COV, mp.get_bvect_cart()),
221  ak_car(mpi),
222  ssjm1_nuf(mpi),
223  ssjm1_nuq(mpi),
224  ssjm1_dzeta(mpi),
225  ssjm1_tggg(mpi),
226  ssjm1_khi(mpi),
227  ssjm1_wshift(mpi, CON, mp.get_bvect_cart())
228 {
229 
230  // Star parameters
231  // -----------------
232 
233  // relativistic is read in the file:
234  size_t nread = fread(&relativistic, sizeof(bool), 1, fich) ;
235  if (nread == 0)
236  cerr << "Star_rot::Star_rot(FILE*): Problem reading the data file." << endl ;
237 
238  // Parameter 1/c^2 is deduced from relativistic:
239  unsurc2 = relativistic ? double(1) : double(0) ;
240 
241  // omega is read in the file:
242  fread_be(&omega, sizeof(double), 1, fich) ;
243 
244 
245  // Read of the saved fields:
246  // ------------------------
247 
248  Scalar nuf_file(mp, *(mp.get_mg()), fich) ;
249  nuf = nuf_file ;
250 
251  Scalar nuq_file(mp, *(mp.get_mg()), fich) ;
252  nuq = nuq_file ;
253 
254  Scalar dzeta_file(mp, *(mp.get_mg()), fich) ;
255  dzeta = dzeta_file ;
256 
257  Scalar tggg_file(mp, *(mp.get_mg()), fich) ;
258  tggg = tggg_file ;
259 
260  Vector w_shift_file(mp, mp.get_bvect_cart(), fich) ;
261  w_shift = w_shift_file ;
262 
263  Scalar khi_shift_file(mp, *(mp.get_mg()), fich) ;
264  khi_shift = khi_shift_file ;
265 
266  fait_shift() ; // constructs shift from w_shift and khi_shift
267  fait_nphi() ; // constructs N^phi from (N^x,N^y,N^z)
268 
269  Scalar ssjm1_nuf_file(mp, *(mp.get_mg()), fich) ;
270  ssjm1_nuf = ssjm1_nuf_file ;
271 
272  Scalar ssjm1_nuq_file(mp, *(mp.get_mg()), fich) ;
273  ssjm1_nuq = ssjm1_nuq_file ;
274 
275  Scalar ssjm1_dzeta_file(mp, *(mp.get_mg()), fich) ;
276  ssjm1_dzeta = ssjm1_dzeta_file ;
277 
278  Scalar ssjm1_tggg_file(mp, *(mp.get_mg()), fich) ;
279  ssjm1_tggg = ssjm1_tggg_file ;
280 
281  Scalar ssjm1_khi_file(mp, *(mp.get_mg()), fich) ;
282  ssjm1_khi = ssjm1_khi_file ;
283 
284  Vector ssjm1_wshift_file(mp, mp.get_bvect_cart(), fich) ;
285  ssjm1_wshift = ssjm1_wshift_file ;
286 
287  // All other fields are initialized to zero :
288  // ----------------------------------------
289  a_car = 0 ;
290  bbb = 0 ;
291  b_car = 0 ;
292  uuu = 0 ;
293  tkij.set_etat_zero() ;
294  ak_car = 0 ;
295 
296  // Pointers of derived quantities initialized to zero
297  // --------------------------------------------------
298  set_der_0x0() ;
299 
300 }
301 
302  //------------//
303  // Destructor //
304  //------------//
305 
307 
308  del_deriv() ;
309 
310 }
311 
312  //----------------------------------//
313  // Management of derived quantities //
314  //----------------------------------//
315 
316 void Star_rot::del_deriv() const {
317 
318  Star::del_deriv() ;
319 
320  if (p_angu_mom != 0x0) delete p_angu_mom ;
321  if (p_tsw != 0x0) delete p_tsw ;
322  if (p_grv2 != 0x0) delete p_grv2 ;
323  if (p_grv3 != 0x0) delete p_grv3 ;
324  if (p_r_circ != 0x0) delete p_r_circ ;
325  if (p_r_circ_merid != 0x0) delete p_r_circ_merid ;
326  if (p_area != 0x0) delete p_area ;
327  if (p_aplat != 0x0) delete p_aplat ;
328  if (p_z_eqf != 0x0) delete p_z_eqf ;
329  if (p_z_eqb != 0x0) delete p_z_eqb ;
330  if (p_z_pole != 0x0) delete p_z_pole ;
331  if (p_mom_quad != 0x0) delete p_mom_quad ;
332  if (p_surf_grav != 0x0) delete p_surf_grav ;
333  if (p_surf_grav_euler != 0x0) delete p_surf_grav_euler ;
334  if (p_radius_Morsink != 0x0) delete p_radius_Morsink ;
335  if (p_r_isco != 0x0) delete p_r_isco ;
336  if (p_f_isco != 0x0) delete p_f_isco ;
337  if (p_lspec_isco != 0x0) delete p_lspec_isco ;
338  if (p_espec_isco != 0x0) delete p_espec_isco ;
339  if (p_f_eq != 0x0) delete p_f_eq ;
340 
342 }
343 
344 
345 void Star_rot::set_der_0x0() const {
346 
348 
349  p_angu_mom = 0x0 ;
350  p_tsw = 0x0 ;
351  p_grv2 = 0x0 ;
352  p_grv3 = 0x0 ;
353  p_r_circ = 0x0 ;
354  p_r_circ_merid = 0x0 ;
355  p_area = 0x0 ;
356  p_aplat = 0x0 ;
357  p_z_eqf = 0x0 ;
358  p_z_eqb = 0x0 ;
359  p_z_pole = 0x0 ;
360  p_mom_quad = 0x0 ;
361  p_surf_grav = 0x0 ;
362  p_surf_grav_euler = 0x0 ;
363  p_radius_Morsink = 0x0 ;
364  p_r_isco = 0x0 ;
365  p_f_isco = 0x0 ;
366  p_lspec_isco = 0x0 ;
367  p_espec_isco = 0x0 ;
368  p_f_eq = 0x0 ;
369 
370 }
371 
373 
375 
376  del_deriv() ;
377 
378 }
379 
380  //--------------//
381  // Assignment //
382  //--------------//
383 
384 // Assignment to another Star_rot
385 // --------------------------------
386 void Star_rot::operator=(const Star_rot& et) {
387 
388  // Assignment of quantities common to all the derived classes of Star
389  Star::operator=(et) ;
390 
391  // Assignement of proper quantities of class Star_rot
392  relativistic = et.relativistic ;
393  unsurc2 = et.unsurc2 ;
394  omega = et.omega ;
395 
396  a_car = et.a_car ;
397  bbb = et.bbb ;
398  b_car = et.b_car ;
399  nphi = et.nphi ;
400  tnphi = et.tnphi ;
401  uuu = et.uuu ;
402  nuf = et.nuf ;
403  nuq = et.nuq ;
404  dzeta = et.dzeta ;
405  tggg = et.tggg ;
406  w_shift = et.w_shift ;
407  khi_shift = et.khi_shift ;
408  tkij = et.tkij ;
409  ak_car = et.ak_car ;
410  ssjm1_nuf = et.ssjm1_nuf ;
411  ssjm1_nuq = et.ssjm1_nuq ;
412  ssjm1_dzeta = et.ssjm1_dzeta ;
413  ssjm1_tggg = et.ssjm1_tggg ;
414  ssjm1_khi = et.ssjm1_khi ;
415  ssjm1_wshift = et.ssjm1_wshift ;
416 
417  del_deriv() ; // Deletes all derived quantities
418 
419 }
420 
421 
422  //--------------//
423  // Outputs //
424  //--------------//
425 
426 // Save in a file
427 // --------------
428 void Star_rot::sauve(FILE* fich) const {
429 
430  Star::sauve(fich) ;
431 
432  fwrite(&relativistic, sizeof(bool), 1, fich) ;
433  fwrite_be(&omega, sizeof(double), 1, fich) ;
434 
435  nuf.sauve(fich) ;
436  nuq.sauve(fich) ;
437  dzeta.sauve(fich) ;
438  tggg.sauve(fich) ;
439  w_shift.sauve(fich) ;
440  khi_shift.sauve(fich) ;
441 
442  ssjm1_nuf.sauve(fich) ;
443  ssjm1_nuq.sauve(fich) ;
444  ssjm1_dzeta.sauve(fich) ;
445  ssjm1_tggg.sauve(fich) ;
446  ssjm1_khi.sauve(fich) ;
447  ssjm1_wshift.sauve(fich) ;
448 
449 
450 }
451 
452 // Printing
453 // --------
454 
455 ostream& Star_rot::operator>>(ostream& ost) const {
456 
457  using namespace Unites ;
458 
459  Star::operator>>(ost) ;
460 
461  double omega_c = get_omega_c() ;
462 
463  ost << endl ;
464 
465  if (omega != __infinity) {
466  ost << "Uniformly rotating star" << endl ;
467  ost << "-----------------------" << endl ;
468 
469  double freq = omega / (2.*M_PI) ;
470  ost << "Omega : " << omega * f_unit
471  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
472  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
473  << endl ;
474 
475  }
476  else {
477  ost << "Differentially rotating star" << endl ;
478  ost << "----------------------------" << endl ;
479 
480  double freq = omega_c / (2.*M_PI) ;
481  ost << "Central value of Omega : " << omega_c * f_unit
482  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
483  ost << "Central rotation period : " << 1000. / (freq * f_unit) << " ms"
484  << endl ;
485  }
486  if (relativistic) {
487  ost << "Relativistic star" << endl ;
488  }
489  else {
490  ost << "Newtonian star" << endl ;
491  }
492  double compact = qpig/(4.*M_PI) * mass_g() / r_circ() ;
493  ost << "Compactness G M_g /(c^2 R_circ) : " << compact << endl ;
494 
495  double nphi_c = nphi.val_grid_point(0, 0, 0, 0) ;
496  if ( (omega_c==0) && (nphi_c==0) ) {
497  ost << "Central N^phi : " << nphi_c << endl ;
498  }
499  else{
500  ost << "Central N^phi/Omega : " << nphi_c / omega_c << endl ;
501  }
502 
503  ost << "Error on the virial identity GRV2 : " << grv2() << endl ;
504  double xgrv3 = grv3(&ost) ;
505  ost << "Error on the virial identity GRV3 : " << xgrv3 << endl ;
506 
507  double mom_quad_38si = mom_quad() * rho_unit * (pow(r_unit, double(5.))
508  / double(1.e38) ) ;
509  ost << "Quadrupole moment Q : " << mom_quad_38si << " 10^38 kg m^2"
510  << endl ;
511  ost << "Q / (M R_circ^2) : "
512  << mom_quad() / ( mass_g() * pow( r_circ(), 2. ) ) << endl ;
513  ost << "c^4 Q / (G^2 M^3) : "
514  << mom_quad() / ( pow(qpig/(4*M_PI), 2.) * pow(mass_g(), 3.) )
515  << endl ;
516 
517  ost << "Angular momentum J : "
518  << angu_mom()/( qpig / (4* M_PI) * msol*msol) << " G M_sol^2 / c"
519  << endl ;
520  ost << "c J / (G M^2) : "
521  << angu_mom()/( qpig / (4* M_PI) * pow(mass_g(), 2.) ) << endl ;
522 
523  if (omega != __infinity) {
524  double mom_iner = angu_mom() / omega ;
525  double mom_iner_38si = mom_iner * rho_unit * (pow(r_unit, double(5.))
526  / double(1.e38) ) ;
527  ost << "Moment of inertia: " << mom_iner_38si << " 10^38 kg m^2"
528  << endl ;
529  }
530 
531  ost << "Ratio T/W : " << tsw() << endl ;
532  ost << "Circumferential equatorial radius R_circ : "
533  << r_circ()/km << " km" << endl ;
534  if (mp.get_mg()->get_np(0) == 1) {
535  ost << "Circumferential polar (meridional) radius R_circ_merid : "
536  << r_circ_merid()/km << " km" << endl ;
537  ost << "Flattening R_circ_merid / R_circ_eq : "
538  << r_circ_merid() / r_circ() << endl ;
539  ost << "Surface area : " << area()/(km*km) << " km^2" << endl ;
540  ost << "Mean radius : " << mean_radius()/km << " km" << endl ;
541  } else {
542  ost <<
543  "Skipping polar radius / surface statements due to number of points in phi direction np > 1"
544  << endl;
545  }
546  ost << "Coordinate equatorial radius r_eq : " << ray_eq()/km << " km"
547  << endl ;
548  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
549  int lsurf = nzet - 1;
550  int nt = mp.get_mg()->get_nt(lsurf) ;
551  int nr = mp.get_mg()->get_nr(lsurf) ;
552  ost << "Equatorial value of the velocity U: "
553  << uuu.val_grid_point(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
554  ost << "Redshift at the equator (forward) : " << z_eqf() << endl ;
555  ost << "Redshift at the equator (backward): " << z_eqb() << endl ;
556  ost << "Redshift at the pole : " << z_pole() << endl ;
557 
558 
559  ost << "Central value of log(N) : "
560  << logn.val_grid_point(0, 0, 0, 0) << endl ;
561 
562  ost << "Central value of dzeta=log(AN) : "
563  << dzeta.val_grid_point(0, 0, 0, 0) << endl ;
564 
565  if ( (omega_c==0) && (nphi_c==0) ) {
566  ost << "Central N^phi : " << nphi_c << endl ;
567  }
568  else{
569  ost << "Central N^phi/Omega : " << nphi_c / omega_c << endl ;
570  }
571 
572  ost << " ... w_shift (NB: components in the star Cartesian frame) [c] : "
573  << endl
574  << w_shift(1).val_grid_point(0, 0, 0, 0) << " "
575  << w_shift(2).val_grid_point(0, 0, 0, 0) << " "
576  << w_shift(3).val_grid_point(0, 0, 0, 0) << endl ;
577 
578  ost << "Central value of khi_shift [km c] : "
579  << khi_shift.val_grid_point(0, 0, 0, 0) / km << endl ;
580 
581  ost << "Central value of B^2 : " << b_car.val_grid_point(0,0,0,0) << endl ;
582 
583  Tbl diff_a_b = diffrel( a_car, b_car ) ;
584  ost <<
585  "Relative discrepancy in each domain between the metric coef. A^2 and B^2 : "
586  << endl ;
587  for (int l=0; l<diff_a_b.get_taille(); l++) {
588  ost << diff_a_b(l) << " " ;
589  }
590  ost << endl ;
591 
592  // Approximate formula for R_isco = 6 R_g (1-(2/3)^1.5 j )
593  // up to the first order in j
594  double jdimless = angu_mom() / ( ggrav * pow(mass_g(), 2.) ) ;
595  double r_grav = ggrav * mass_g() ;
596  double r_isco_appr = 6. * r_grav * ( 1. - pow(2./3.,1.5) * jdimless ) ;
597 
598  // Approximate formula for the ISCO frequency
599  // freq_ms = 6^{-1.5}/2pi/R_g (1+11*6^(-1.5) j )
600  // up to the first order in j
601  double f_isco_appr = ( 1. + 11. /6. /sqrt(6.) * jdimless ) / r_grav /
602  (12. * M_PI ) / sqrt(6.) ;
603 
604  ost << endl << "Innermost stable circular orbit (ISCO) : " << endl ;
605  double xr_isco = r_isco(&ost) ;
606  ost <<" circumferential radius r_isco = "
607  << xr_isco / km << " km" << endl ;
608  ost <<" (approx. 6M + 1st order in j : "
609  << r_isco_appr / km << " km)" << endl ;
610  ost <<" (approx. 6M : "
611  << 6. * r_grav / km << " km)" << endl ;
612  ost <<" orbital frequency f_isco = "
613  << f_isco() * f_unit << " Hz" << endl ;
614  ost <<" (approx. 1st order in j : "
615  << f_isco_appr * f_unit << " Hz)" << endl ;
616 
617 
618  return ost ;
619 
620 }
621 
622 
623 void Star_rot::partial_display(ostream& ost) const {
624 
625  using namespace Unites ;
626 
627  double omega_c = get_omega_c() ;
628  double freq = omega_c / (2.*M_PI) ;
629  ost << "Central Omega : " << omega_c * f_unit
630  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
631  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
632  << endl ;
633  ost << endl << "Central enthalpy : " << ent.val_grid_point(0,0,0,0) << " c^2" << endl ;
634  ost << "Central proper baryon density : " << nbar.val_grid_point(0,0,0,0)
635  << " x 0.1 fm^-3" << endl ;
636  ost << "Central proper energy density : " << ener.val_grid_point(0,0,0,0)
637  << " rho_nuc c^2" << endl ;
638  ost << "Central pressure : " << press.val_grid_point(0,0,0,0)
639  << " rho_nuc c^2" << endl ;
640 
641  ost << "Central value of log(N) : "
642  << logn.val_grid_point(0, 0, 0, 0) << endl ;
643  ost << "Central lapse N : " << nn.val_grid_point(0,0,0,0) << endl ;
644  ost << "Central value of dzeta=log(AN) : "
645  << dzeta.val_grid_point(0, 0, 0, 0) << endl ;
646  ost << "Central value of A^2 : " << a_car.val_grid_point(0,0,0,0) << endl ;
647  ost << "Central value of B^2 : " << b_car.val_grid_point(0,0,0,0) << endl ;
648 
649  double nphi_c = nphi.val_grid_point(0, 0, 0, 0) ;
650  if ( (omega_c==0) && (nphi_c==0) ) {
651  ost << "Central N^phi : " << nphi_c << endl ;
652  }
653  else{
654  ost << "Central N^phi/Omega : " << nphi_c / omega_c
655  << endl ;
656  }
657 
658 
659  int lsurf = nzet - 1;
660  int nt = mp.get_mg()->get_nt(lsurf) ;
661  int nr = mp.get_mg()->get_nr(lsurf) ;
662  ost << "Equatorial value of the velocity U: "
663  << uuu.val_grid_point(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
664 
665  ost << endl
666  << "Coordinate equatorial radius r_eq = "
667  << ray_eq()/km << " km" << endl ;
668  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
669  if (mp.get_mg()->get_np(0) == 1)
670  ost << "Flattening r_circ_pole/r_circ_eq : "
671  << r_circ_merid() / r_circ() << endl ;
672 
673 }
674 
675 
676 double Star_rot::get_omega_c() const {
677 
678  return omega ;
679 
680 }
681 
682 // display_poly
683 // ------------
684 
685 void Star_rot::display_poly(ostream& ost) const {
686 
687  using namespace Unites ;
688 
689  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
690 
691  if (p_eos_poly != 0x0) {
692 
693  double kappa = p_eos_poly->get_kap() ;
694 
695  // kappa^{n/2}
696  double kap_ns2 = pow( kappa, 0.5 /(p_eos_poly->get_gam()-1) ) ;
697 
698  // Polytropic unit of length in terms of r_unit :
699  double r_poly = kap_ns2 / sqrt(ggrav) ;
700 
701  // Polytropic unit of time in terms of t_unit :
702  double t_poly = r_poly ;
703 
704  // Polytropic unit of mass in terms of m_unit :
705  double m_poly = r_poly / ggrav ;
706 
707  // Polytropic unit of angular momentum in terms of j_unit :
708  double j_poly = r_poly * r_poly / ggrav ;
709 
710  // Polytropic unit of density in terms of rho_unit :
711  double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
712 
713  ost.precision(10) ;
714  ost << endl << "Quantities in polytropic units : " << endl ;
715  ost << "==============================" << endl ;
716  ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
717  ost << " n_c : " << nbar.val_grid_point(0, 0, 0, 0) / rho_poly << endl ;
718  ost << " e_c : " << ener.val_grid_point(0, 0, 0, 0) / rho_poly << endl ;
719  ost << " Omega_c : " << get_omega_c() * t_poly << endl ;
720  ost << " P_c : " << 2.*M_PI / get_omega_c() / t_poly << endl ;
721  ost << " M_bar : " << mass_b() / m_poly << endl ;
722  ost << " M : " << mass_g() / m_poly << endl ;
723  ost << " J : " << angu_mom() / j_poly << endl ;
724  ost << " r_eq : " << ray_eq() / r_poly << endl ;
725  ost << " R_circ_eq : " << r_circ() / r_poly << endl ;
726  ost << " R_circ_merid : " << r_circ_merid() / r_poly << endl ;
727  ost << " R_mean : " << mean_radius() / r_poly << endl ;
728 
729  }
730 
731 
732 }
733 
734 
735 
736  //-------------------------//
737  // Computational methods //
738  //-------------------------//
739 
741 
742  Vector d_khi = khi_shift.derive_con( mp.flat_met_cart() ) ;
743 
744  d_khi.dec_dzpuis(2) ; // divide by r^2 in the external compactified domain
745 
746  // x_k dW^k/dx_i
747  Scalar xx(mp) ;
748  Scalar yy(mp) ;
749  Scalar zz(mp) ;
750  Scalar sintcosp(mp) ;
751  Scalar sintsinp(mp) ;
752  Scalar cost(mp) ;
753  xx = mp.x ;
754  yy = mp.y ;
755  zz = mp.z ;
756  sintcosp = mp.sint * mp.cosp ;
757  sintsinp = mp.sint * mp.sinp ;
758  cost = mp.cost ;
759 
760  int nz = mp.get_mg()->get_nzone() ;
761  Vector xk(mp, COV, mp.get_bvect_cart()) ;
762  xk.set(1) = xx ;
763  xk.set(1).set_domain(nz-1) = sintcosp.domain(nz-1) ;
764  xk.set(1).set_dzpuis(-1) ;
765  xk.set(2) = yy ;
766  xk.set(2).set_domain(nz-1) = sintsinp.domain(nz-1) ;
767  xk.set(2).set_dzpuis(-1) ;
768  xk.set(3) = zz ;
769  xk.set(3).set_domain(nz-1) = cost.domain(nz-1) ;
770  xk.set(3).set_dzpuis(-1) ;
771  xk.std_spectral_base() ;
772 
774 
775  Vector x_d_w = contract(xk, 0, d_w, 0) ;
776  x_d_w.dec_dzpuis() ;
777 
778  double lambda = double(1) / double(3) ;
779 
780  if ( mp.get_mg()->get_np(0) == 1 ) {
781  lambda = 0 ;
782  }
783 
784  beta = - (lambda+2)/2./(lambda+1) * w_shift
785  + (lambda/2./(lambda+1)) * (d_khi + x_d_w) ;
786 
787 }
788 
790 
791  if ( (beta(1).get_etat() == ETATZERO) && (beta(2).get_etat() == ETATZERO) ) {
792  tnphi = 0 ;
793  nphi = 0 ;
794  return ;
795  }
796 
797  // Computation of tnphi
798  // --------------------
799 
800  mp.comp_p_from_cartesian( -beta(1), -beta(2), tnphi ) ;
801 
802  // Computation of nphi
803  // -------------------
804 
805  nphi = tnphi ;
806  nphi.div_rsint() ;
807 
808 }
809 
810 }
virtual double mass_g() const
Gravitational mass.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: star.C:319
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: star_rot.h:243
Scalar dzeta
Metric potential .
Definition: star_rot.h:152
virtual double mean_radius() const
Mean star radius from the area .
Scalar ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta ...
Definition: star_rot.h:221
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:122
Map & mp
Mapping associated with the star.
Definition: star.h:180
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition: star_rot.h:267
double * p_r_isco
Circumferential radius of the ISCO.
Definition: star_rot.h:262
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: star_rot.h:144
virtual double r_circ_merid() const
Crcumferential meridional radius.
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
double * p_f_eq
Orbital frequency at the equator.
Definition: star_rot.h:268
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
double * p_r_circ
Circumferential radius (equator)
Definition: star_rot.h:254
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: star_rot.h:178
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: star_rot.h:112
const Tbl & domain(int l) const
Read-only of the value in a given domain.
Definition: scalar.h:637
Scalar tnphi
Component of the shift vector.
Definition: star_rot.h:136
double * p_mom_quad
Quadrupole moment.
Definition: star_rot.h:261
Scalar bbb
Metric factor B.
Definition: star_rot.h:125
Lorene prototypes.
Definition: app_hor.h:67
const Eos & eos
Equation of state of the stellar matter.
Definition: star.h:185
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...
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
virtual double z_pole() const
Redshift factor at North pole.
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
double * p_f_isco
Orbital frequency of the ISCO.
Definition: star_rot.h:263
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
virtual double grv2() const
Error on the virial identity GRV2.
Base class for coordinate mappings.
Definition: map.h:697
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
Definition: star_rot.C:676
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star_rot.C:316
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
virtual void sauve(FILE *) const
Save in a binary file.
Definition: tensor.C:916
Scalar nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: star_rot.h:149
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot.h:119
Scalar nbar
Baryon density in the fluid frame.
Definition: star.h:192
Scalar ent
Log-enthalpy.
Definition: star.h:190
Base class for stars.
Definition: star.h:175
Vector beta
Shift vector.
Definition: star.h:228
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Definition: scalar_deriv.C:401
Tensor field of valence 1.
Definition: vector.h:188
Tbl * p_radius_Morsink
Physical radius close to the definition of Morsink et al. (2007)
Definition: star_rot.h:274
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:337
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:627
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
virtual double area() const
Integrated surface area in .
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Definition: eos_poly.C:271
Class for isolated rotating stars.
Definition: star_rot.h:104
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Definition: star_rot.h:226
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:649
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:818
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: star_rot.C:345
void operator=(const Star_rot &)
Assignment to another Star_rot.
Definition: star_rot.C:386
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star_rot.C:455
void operator=(const Star &)
Assignment to another Star.
Definition: star.C:354
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
virtual ~Star_rot()
Destructor.
Definition: star_rot.C:306
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
Definition: star_rot_isco.C:74
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition: star_rot.h:117
Coord sint
Definition: map.h:748
double * p_tsw
Ratio T/W.
Definition: star_rot.h:251
Scalar nphi
Metric coefficient .
Definition: star_rot.h:131
double * p_z_eqf
Forward redshift factor at equator.
Definition: star_rot.h:258
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
Scalar ener
Total energy density in the fluid frame.
Definition: star.h:193
double * p_grv2
Error on the virial identity GRV2.
Definition: star_rot.h:252
double * p_area
Integrated surface area.
Definition: star_rot.h:257
double * p_r_circ_merid
Circumferential radius (meridian)
Definition: star_rot.h:255
Scalar press
Fluid pressure.
Definition: star.h:194
virtual double mass_b() const
Baryon mass.
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition: star_rot.h:265
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: star_rot.C:372
virtual double f_isco() const
Orbital frequency at the innermost stable circular orbit (ISCO).
virtual double mom_quad() const
Quadrupole moment.
Scalar ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition: star_rot.h:216
virtual double angu_mom() const
Angular momentum.
double * p_z_pole
Redshift factor at North pole.
Definition: star_rot.h:260
Polytropic equation of state (relativistic case).
Definition: eos.h:812
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: star_rot.C:685
Star_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition: star_rot.C:103
Coord sinp
Definition: map.h:750
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
double * p_grv3
Error on the virial identity GRV3.
Definition: star_rot.h:253
Tbl * p_surf_grav_euler
Surface gravity for the Eulerian observer (along the theta direction)
Definition: star_rot.h:272
Sym_tensor tkij
Tensor related to the extrinsic curvature tensor by .
Definition: star_rot.h:185
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
virtual double r_circ() const
Circumferential equatorial radius.
double * p_angu_mom
Angular momentum.
Definition: star_rot.h:250
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tensor handling.
Definition: tensor.h:300
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Scalar logn
Logarithm of the lapse N .
Definition: star.h:222
Scalar b_car
Square of the metric factor B.
Definition: star_rot.h:128
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
double * p_aplat
Flatening r_pole/r_eq.
Definition: star_rot.h:256
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
Definition: tensor.C:1024
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star.C:301
virtual double z_eqb() const
Backward redshift factor at equator.
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition: star_rot.C:623
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: star.C:333
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:818
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:664
Scalar nn
Lapse function N .
Definition: star.h:225
Coord y
y coordinate centered on the grid
Definition: map.h:754
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: star_rot.h:234
int get_taille() const
Gives the total size (ie dim.taille)
Definition: tbl.h:417
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:139
Coord cosp
Definition: map.h:751
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
Scalar ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition: star_rot.h:210
Coord x
x coordinate centered on the grid
Definition: map.h:753
void div_rsint()
Division by everywhere; dzpuis is not changed.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star.C:420
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tensor.C:529
virtual void sauve(FILE *) const
Save in a file.
Definition: star_rot.C:428
virtual void sauve(FILE *) const
Save in a file.
Definition: star.C:400
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:507
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
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
Scalar tggg
Metric potential .
Definition: star_rot.h:155
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Definition: star_rot.C:789
Scalar ak_car
Scalar .
Definition: star_rot.h:204
virtual double aplat() const
Flatening r_pole/r_eq.
Tbl * p_surf_grav
Surface gravity (along the theta direction)
Definition: star_rot.h:270
Coord z
z coordinate centered on the grid
Definition: map.h:755
Vector w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: star_rot.h:168
virtual double z_eqf() const
Forward redshift factor at equator.
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata&#39;s prescription [Prog.
Definition: star_rot.C:740
virtual double tsw() const
Ratio T/W.
double * p_z_eqb
Backward redshift factor at equator.
Definition: star_rot.h:259
Coord cost
Definition: map.h:749