LORENE
star_rot_diff_cfc.C
1 /*
2  * Methods of class Star_rot_diff_CFC
3  *
4  * (see file star_rot_diff_cfc.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2025 Santiago Jaraba
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  // Lorene headers
29 #include "star_rot_diff_cfc.h"
30 #include "eos.h"
31 #include "nbr_spx.h"
32 #include "utilitaires.h"
33 #include "unites.h"
34 
35 
36  //--------------//
37  // Constructors //
38  //--------------//
39 
40 // Standard constructor
41 //-------------------------
42 namespace Lorene {
43 Star_rot_diff_CFC::Star_rot_diff_CFC(Map& mpi, int nzet_i, const Eos& eos_i, int relat_i, int filter,
44  double (*omega_frot_i)(double, const Tbl&),
45  double (*primfrot_i)(double, const Tbl&),
46  const Tbl& par_frot_i)
47  :Star_rot_CFC(mpi, nzet_i, eos_i, relat_i, filter),
48  omega_frot(omega_frot_i),
49  primfrot(primfrot_i),
50  par_frot(par_frot_i),
51  omega_field(mpi),
52  prim_field(mpi)
53  {
54 
55  // To make sure that omega is not used
56  omega = __infinity ;
57 
58  // Initialization to a static state :
59  omega_field = 0 ;
60  prim_field = 0 ;
61  omega_min = 0 ;
62  omega_max = 0 ;
63 
64 }
65 
66 // Copy constructor
67 //-----------------
69  : Star_rot_CFC(star),
70  omega_frot(star.omega_frot),
71  primfrot(star.primfrot),
72  par_frot(star.par_frot),
73  omega_field(star.omega_field),
74  omega_min(star.omega_min),
75  omega_max(star.omega_max),
76  prim_field(star.prim_field)
77  {}
78 
79 //Constructor from a file
80 //------------------------
81 Star_rot_diff_CFC::Star_rot_diff_CFC(Map& mpi, const Eos& eos_i, FILE* fich,
82  double (*omega_frot_i)(double, const Tbl&),
83  double (*primfrot_i)(double, const Tbl&) )
84  : Star_rot_CFC(mpi, eos_i, fich),
85  omega_frot(omega_frot_i),
86  primfrot(primfrot_i),
87  par_frot(fich),
88  omega_field(mpi),
89  prim_field(mpi)
90  {
91 
92  Scalar omega_field_file(mp, *(mp.get_mg()), fich) ;
93  omega_field = omega_field_file ;
94  fait_prim_field() ;
95 
96  // omega_min and omega_max are read in the file:
97  fread_be(&omega_min, sizeof(double), 1, fich) ;
98  fread_be(&omega_max, sizeof(double), 1, fich) ;
99 
100 }
101 
102  //------------//
103  // Destructor //
104  //------------//
105 
107 
108  //---------------//
109  // Assignment //
110  //---------------//
111 
112 // Assignment to another Star_rot_Dirac_diff
113 // ------------------------------------
114 
116 
117  // Assignment of quantities common to all the derived classes of
118  // Star_rot_CFC
120 
121  // Assignment of proper quantities of class Star_rot
122  omega_frot = star.omega_frot ;
123  primfrot = star.primfrot ;
124  par_frot = star.par_frot ;
125  omega_field = star.omega_field ;
126  prim_field = star.prim_field ;
127  omega_min = star.omega_min ;
128  omega_max = star.omega_max ;
129 
130 }
131 
132 
133  //--------------//
134  // Outputs //
135  //--------------//
136 
137 // Save in a file
138 // --------------
139 
140 void Star_rot_diff_CFC::sauve(FILE* fich) const {
141 
142  Star::sauve(fich) ;
143 
144  par_frot.sauve(fich) ;
145 
146  omega_field.sauve(fich) ;
147 
148  fwrite_be(&omega_min, sizeof(double), 1, fich) ;
149  fwrite_be(&omega_max, sizeof(double), 1, fich) ;
150 
151 }
152 
153 
154 // Printing
155 // --------
156 
157 ostream& Star_rot_diff_CFC::operator>>(ostream& ost) const {
158 
159  using namespace Unites ;
160 
162 
163  // TODO: check this section and figure out what to output
164 
165  ost << endl ;
166  ost << "Differentially rotating star" << endl ;
167  ost << "-----------------------------" << endl ;
168 
169  ost << endl << "Parameters of Omega(F) : " << endl ;
170  ost << par_frot << endl ;
171 
172  ost << "Min, Max of Omega/(2pi) : " << omega_min / (2*M_PI) * f_unit
173  << " Hz, " << omega_max / (2*M_PI) * f_unit << " Hz" << endl ;
174  int lsurf = nzet - 1;
175  int nt = mp.get_mg()->get_nt(lsurf) ;
176  int nr = mp.get_mg()->get_nr(lsurf) ;
177  ost << "Central, equatorial value of Omega: "
178  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit << " rad/s, "
179  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit << " rad/s" << endl ;
180 
181  ost << "Central, equatorial value of Omega/(2 Pi): "
182  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit / (2*M_PI) << " Hz, "
183  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit / (2*M_PI)
184  << " Hz" << endl ;
185 
186  double nbar_max = max( max( nbar ) ) ;
187  double ener_max = max( max( ener ) ) ;
188  double press_max = max( max( press ) ) ;
189  ost << "Max prop. bar. dens. : " << nbar_max
190  << " x 0.1 fm^-3 = " << nbar_max / nbar.val_grid_point(0, 0, 0, 0) << " central"
191  << endl ;
192  ost << "Max prop. ener. dens. (e_max) : " << ener_max
193  << " rho_nuc c^2 = " << ener_max / ener.val_grid_point(0, 0, 0, 0) << " central"
194  << endl ;
195  ost << "Max pressure : " << press_max
196  << " rho_nuc c^2 = " << press_max / press.val_grid_point(0, 0, 0, 0) << " central"
197  << endl ;
198 
199  // Length scale set by the maximum energy density:
200  double len_rho = 1. / sqrt( ggrav * ener_max ) ;
201  ost << endl << "Value of A = par_frot(1) in units of e_max^{1/2} : " <<
202  par_frot(1) / len_rho << endl ;
203 
204  ost << "Value of A / r_eq : " <<
205  par_frot(1) / ray_eq() << endl ;
206 
207  ost << "r_p/r_eq : " << aplat() << endl ;
208  ost << "KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
209  angu_mom() * angu_mom() / pow(len_rho, 0.6666666666666666)
210  / pow(mass_b(), 3.3333333333333333)
211  / pow(ggrav, 1.3333333333333333) << endl ;
212 
213  ost << "M e_max^{1/2} : " << ggrav * mass_g() / len_rho << endl ;
214 
215  ost << "r_eq e_max^{1/2} : " << ray_eq() / len_rho << endl ;
216 
217  ost << "T/W : " << tsw() << endl ;
218 
219  ost << "Omega_c / e_max^{1/2} : " << get_omega_c() * len_rho << endl ;
220 
221  display_poly(ost) ;
222 
223  return ost ;
224 
225 }
226 
227 // display_poly
228 // ------------
229 
230 void Star_rot_diff_CFC::display_poly(ostream& ost) const {
231 
232  using namespace Unites ;
233 
234  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
235 
236  if (p_eos_poly != 0x0) {
237 
238  double kappa = p_eos_poly->get_kap() ;
239 
240  // kappa^{n/2}
241  double kap_ns2 = pow( kappa, 0.5 /(p_eos_poly->get_gam()-1) ) ;
242 
243  // Polytropic unit of length in terms of r_unit :
244  double r_poly = kap_ns2 / sqrt(ggrav) ;
245 
246  // Polytropic unit of time in terms of t_unit :
247  double t_poly = r_poly ;
248 
249  // Polytropic unit of mass in terms of m_unit :
250  double m_poly = r_poly / ggrav ;
251 
252  // Polytropic unit of angular momentum in terms of j_unit :
253  double j_poly = r_poly * r_poly / ggrav ;
254 
255  // Polytropic unit of density in terms of rho_unit :
256  double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
257 
258  ost.precision(10) ;
259  ost << endl << "Quantities in polytropic units : " << endl ;
260  ost << "==============================" << endl ;
261  ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
262  ost << " n_c : " << nbar.val_grid_point(0, 0, 0, 0) / rho_poly << endl ;
263  ost << " e_c : " << ener.val_grid_point(0, 0, 0, 0) / rho_poly << endl ;
264  ost << " Omega_c : " << get_omega_c() * t_poly << endl ;
265  ost << " P_c : " << 2.*M_PI / get_omega_c() / t_poly << endl ;
266  ost << " M_bar : " << mass_b() / m_poly << endl ;
267  ost << " M : " << mass_g() / m_poly << endl ;
268  ost << " J : " << angu_mom() / j_poly << endl ;
269  ost << " r_eq : " << ray_eq() / r_poly << endl ;
270  ost << " R_circ_eq : " << r_circ() / r_poly << endl ;
271 
272  ost << " n_max : " << max( max( nbar ) ) / rho_poly << endl ;
273  ost << " e_max : " << max( max( ener ) ) / rho_poly << endl ;
274  ost << " P_min : " << 2.*M_PI / omega_max / t_poly << endl ;
275  ost << " P_max : " << 2.*M_PI / omega_min / t_poly << endl ;
276 
277  int lsurf = nzet - 1;
278  int nt = mp.get_mg()->get_nt(lsurf) ;
279  int nr = mp.get_mg()->get_nr(lsurf) ;
280  ost << " P_eq : " << 2.*M_PI /
281  omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) / t_poly << endl ;
282 
283  }
284  }
285 
286 
287 
288  //-----------------------//
289  // Miscellaneous //
290  //-----------------------//
291 
292 double Star_rot_diff_CFC::omega_funct(double F) const {
293 
294  return omega_frot(F, par_frot) ;
295 
296 }
297 
298 double Star_rot_diff_CFC::prim_funct_omega(double F) const {
299 
300  return primfrot(F, par_frot) ;
301 
302 }
303 
305 
306  return omega_field.val_grid_point(0, 0, 0, 0) ;
307 
308 }
309 
310 void Star_rot_diff_CFC::partial_display(ostream& ost) const {
311 
312  using namespace Unites ;
313 
314  double omega_c = get_omega_c() ;
315  double freq = omega_c / (2.*M_PI) ;
316  ost << "Central Omega : " << omega_c * f_unit
317  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
318  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
319  << endl ;
320  ost << endl << "Central enthalpy : " << ent.val_grid_point(0,0,0,0) << " c^2" << endl ;
321  ost << "Central proper baryon density : " << nbar.val_grid_point(0,0,0,0)
322  << " x 0.1 fm^-3" << endl ;
323  ost << "Central proper energy density : " << ener.val_grid_point(0,0,0,0)
324  << " rho_nuc c^2" << endl ;
325  ost << "Central pressure : " << press.val_grid_point(0,0,0,0)
326  << " rho_nuc c^2" << endl ;
327 
328  ost << "Central value of log(N) : "
329  << logn.val_grid_point(0, 0, 0, 0) << endl ;
330  ost << "Central lapse N : " << nn.val_grid_point(0,0,0,0) << endl ;
331  ost << "Central value of psi : " << psi.val_grid_point(0,0,0,0) << endl ;
332 
333  double nphi_c = beta(3).val_grid_point(0, 0, 0, 0) ;
334  if ( (omega_c==0) && (nphi_c==0) ) {
335  ost << "Central azimuthal shift : " << nphi_c << endl ;
336  }
337  else{
338  ost << "Central azimuthal shift /Omega : " << nphi_c / omega_c
339  << endl ;
340  }
341 
342 
343  int lsurf = nzet - 1;
344  int nt = mp.get_mg()->get_nt(lsurf) ;
345  int nr = mp.get_mg()->get_nr(lsurf) ;
346  ost << "Equatorial value of the velocity U: "
347  << u_euler(3).val_grid_point(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
348 
349  ost << endl
350  << "Coordinate equatorial radius r_eq = "
351  << ray_eq()/km << " km" << endl ;
352  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
353 
354 }
355 
356 }
void operator=(const Star_rot_diff_CFC &)
Assignment to another Star_rot_diff.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star_rot_cfc.C:260
Map & mp
Mapping associated with the star.
Definition: star.h:180
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
double omega_min
Minimum value of .
Scalar psi
Conformal factor .
Definition: star_rot_cfc.h:65
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
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
virtual double angu_mom() const
Angular momentum.
Base class for coordinate mappings.
Definition: map.h:697
virtual double r_circ() const
Circumferential equatorial radius.
void operator=(const Star_rot_CFC &)
Assignment to another Star_rot_CFC.
Definition: star_rot_cfc.C:212
Class for differentially rotating stars in Conformal Flatness Condition and maximal slicing...
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Scalar nbar
Baryon density in the fluid frame.
Definition: star.h:192
Scalar ent
Log-enthalpy.
Definition: star.h:190
Vector beta
Shift vector.
Definition: star.h:228
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Definition: eos_poly.C:271
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 ~Star_rot_diff_CFC()
Destructor.
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
Class for relativistic rotating stars in Conformal Flatness Condition and maximal slicing...
Definition: star_rot_cfc.h:39
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
double prim_funct_omega(double F) const
Evaluates the primitive of .
Scalar ener
Total energy density in the fluid frame.
Definition: star.h:193
virtual double aplat() const
Flattening r_pole/r_eq.
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot_cfc.h:60
Scalar press
Fluid pressure.
Definition: star.h:194
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
virtual void sauve(FILE *) const
Save in a file.
Polytropic equation of state (relativistic case).
Definition: eos.h:812
Tbl par_frot
Parameters of the function .
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
virtual double mass_b() const
Baryonic mass.
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: star.h:207
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
Scalar logn
Logarithm of the lapse N .
Definition: star.h:222
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
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
double(* omega_frot)(double, const Tbl &)
Function defining the rotation profile.
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Scalar nn
Lapse function N .
Definition: star.h:225
double omega_max
Maximum value of .
Star_rot_diff_CFC(Map &mp_i, int nzet_i, const Eos &eos_i, int relat_i, int filter, double(*omega_frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
void sauve(FILE *) const
Save in a file.
Definition: tbl.C:329
virtual void sauve(FILE *) const
Save in a file.
Definition: star.C:400
virtual double tsw() const
Ratio T/W.
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
double omega_funct(double F) const
Evaluates , where F is linked to the components of the fluid 4-velocity by .
void fait_prim_field()
Computes the member prim_field from omega_field .
virtual double mass_g() const
Gravitational mass.