LORENE
star_rot_diff.C
1 /*
2  * Methods of class Star_rot_diff
3  *
4  * (see file star_rot_diff.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 
29 // Lorene headers
30 #include "star_rot_diff.h"
31 #include "eos.h"
32 #include "nbr_spx.h"
33 #include "utilitaires.h"
34 #include "unites.h"
35 
36 
37  //--------------//
38  // Constructors //
39  //--------------//
40 
41 // Standard constructor
42 //-------------------------
43 namespace Lorene {
44 Star_rot_diff::Star_rot_diff(Map& mpi, int nzet_i, bool relat, const Eos& eos_i,
45  double (*omega_frot_i)(double, const Tbl&),
46  double (*primfrot_i)(double, const Tbl&),
47  const Tbl& par_frot_i)
48  :Star_rot(mpi, nzet_i, relat, eos_i),
49  omega_frot(omega_frot_i),
50  primfrot(primfrot_i),
51  par_frot(par_frot_i),
52  omega_field(mpi),
53  prim_field(mpi)
54  {
55 
56  // To make sure that omega is not used
57  omega = __infinity ;
58 
59  // Initialization to a static state :
60  omega_field = 0 ;
61  prim_field = 0 ;
62  omega_min = 0 ;
63  omega_max = 0 ;
64 
65 }
66 
67 // Copy constructor
68 //-----------------
70  : Star_rot(star),
71  omega_frot(star.omega_frot),
72  primfrot(star.primfrot),
73  par_frot(star.par_frot),
74  omega_field(star.omega_field),
75  omega_min(star.omega_min),
76  omega_max(star.omega_max),
77  prim_field(star.prim_field)
78  {}
79 
80 
81 //Constructor from a file
82 //------------------------
83 Star_rot_diff::Star_rot_diff(Map& mpi, const Eos& eos_i, FILE* fich,
84  double (*omega_frot_i)(double, const Tbl&),
85  double (*primfrot_i)(double, const Tbl&) )
86  : Star_rot(mpi, eos_i, fich),
87  omega_frot(omega_frot_i),
88  primfrot(primfrot_i),
89  par_frot(fich),
90  omega_field(mpi),
91  prim_field(mpi)
92  {
93 
94  Scalar omega_field_file(mp, *(mp.get_mg()), fich) ;
95  omega_field = omega_field_file ;
96  fait_prim_field() ;
97 
98  // omega_min and omega_max are read in the file:
99  fread_be(&omega_min, sizeof(double), 1, fich) ;
100  fread_be(&omega_max, sizeof(double), 1, fich) ;
101 
102 }
103 
104 
105  //------------//
106  // Destructor //
107  //------------//
108 
110 
111 
112  //---------------//
113  // Assignment //
114  //---------------//
115 
116 // Assignment to another Star_rot_Dirac_diff
117 // ------------------------------------
118 
120 
121  // Assignment of quantities common to all the derived classes of
122  // Star_rot
123  Star_rot::operator=(star) ;
124 
125  // Assignment of proper quantities of class Star_rot
126  omega_frot = star.omega_frot ;
127  primfrot = star.primfrot ;
128  par_frot = star.par_frot ;
129  omega_field = star.omega_field ;
130  prim_field = star.prim_field ;
131  omega_min = star.omega_min ;
132  omega_max = star.omega_max ;
133 
134 }
135 
136 
137  //--------------//
138  // Outputs //
139  //--------------//
140 
141 // Save in a file
142 // --------------
143 
144 void Star_rot_diff::sauve(FILE* fich) const {
145 
146  Star::sauve(fich) ;
147 
148  par_frot.sauve(fich) ;
149 
150  omega_field.sauve(fich) ;
151 
152  fwrite_be(&omega_min, sizeof(double), 1, fich) ;
153  fwrite_be(&omega_max, sizeof(double), 1, fich) ;
154 
155 }
156 
157 
158 // Printing
159 // --------
160 
161 ostream& Star_rot_diff::operator>>(ostream& ost) const {
162 
163  using namespace Unites ;
164 
165  Star_rot::operator>>(ost) ;
166 
167  ost << endl ;
168  ost << "Differentially rotating star" << endl ;
169  ost << "-----------------------------" << endl ;
170 
171  ost << endl << "Parameters of Omega(F) : " << endl ;
172  ost << par_frot << endl ;
173 
174  ost << "Min, Max of Omega/(2pi) : " << omega_min / (2*M_PI) * f_unit
175  << " Hz, " << omega_max / (2*M_PI) * f_unit << " Hz" << endl ;
176  int lsurf = nzet - 1;
177  int nt = mp.get_mg()->get_nt(lsurf) ;
178  int nr = mp.get_mg()->get_nr(lsurf) ;
179  ost << "Central, equatorial value of Omega: "
180  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit << " rad/s, "
181  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit << " rad/s" << endl ;
182 
183  ost << "Central, equatorial value of Omega/(2 Pi): "
184  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit / (2*M_PI) << " Hz, "
185  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit / (2*M_PI)
186  << " Hz" << endl ;
187 
188  double nbar_max = max( max( nbar ) ) ;
189  double ener_max = max( max( ener ) ) ;
190  double press_max = max( max( press ) ) ;
191  ost << "Max prop. bar. dens. : " << nbar_max
192  << " x 0.1 fm^-3 = " << nbar_max / nbar.val_grid_point(0, 0, 0, 0) << " central"
193  << endl ;
194  ost << "Max prop. ener. dens. (e_max) : " << ener_max
195  << " rho_nuc c^2 = " << ener_max / ener.val_grid_point(0, 0, 0, 0) << " central"
196  << endl ;
197  ost << "Max pressure : " << press_max
198  << " rho_nuc c^2 = " << press_max / press.val_grid_point(0, 0, 0, 0) << " central"
199  << endl ;
200 
201  // Length scale set by the maximum energy density:
202  double len_rho = 1. / sqrt( ggrav * ener_max ) ;
203  ost << endl << "Value of A = par_frot(1) in units of e_max^{1/2} : " <<
204  par_frot(1) / len_rho << endl ;
205 
206  ost << "Value of A / r_eq : " <<
207  par_frot(1) / ray_eq() << endl ;
208 
209  ost << "r_p/r_eq : " << aplat() << endl ;
210  ost << "KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
211  angu_mom() * angu_mom() / pow(len_rho, 0.6666666666666666)
212  / pow(mass_b(), 3.3333333333333333)
213  / pow(ggrav, 1.3333333333333333) << endl ;
214 
215  ost << "M e_max^{1/2} : " << ggrav * mass_g() / len_rho << endl ;
216 
217  ost << "r_eq e_max^{1/2} : " << ray_eq() / len_rho << endl ;
218 
219  ost << "T/W : " << tsw() << endl ;
220 
221  ost << "Omega_c / e_max^{1/2} : " << get_omega_c() * len_rho << endl ;
222 
223  display_poly(ost) ;
224 
225  return ost ;
226 
227 }
228 
229 // display_poly
230 // ------------
231 
232 void Star_rot_diff::display_poly(ostream& ost) const {
233 
234  using namespace Unites ;
235 
236  Star_rot::display_poly( ost ) ;
237 
238  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
239 
240  if (p_eos_poly != 0x0) {
241 
242  double kappa = p_eos_poly->get_kap() ;
243  double gam = p_eos_poly->get_gam() ; ;
244 
245  // kappa^{n/2}
246  double kap_ns2 = pow( kappa, 0.5 /(gam-1) ) ;
247 
248  // Polytropic unit of length in terms of r_unit :
249  double r_poly = kap_ns2 / sqrt(ggrav) ;
250 
251  // Polytropic unit of time in terms of t_unit :
252  double t_poly = r_poly ;
253 
254  // Polytropic unit of density in terms of rho_unit :
255  double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
256 
257  ost.precision(10) ;
258  ost << " n_max : " << max( max( nbar ) ) / rho_poly << endl ;
259  ost << " e_max : " << max( max( ener ) ) / rho_poly << endl ;
260  ost << " P_min : " << 2.*M_PI / omega_max / t_poly << endl ;
261  ost << " P_max : " << 2.*M_PI / omega_min / t_poly << endl ;
262 
263  int lsurf = nzet - 1;
264  int nt = mp.get_mg()->get_nt(lsurf) ;
265  int nr = mp.get_mg()->get_nr(lsurf) ;
266  ost << " P_eq : " << 2.*M_PI /
267  omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) / t_poly << endl ;
268 
269  }
270 }
271 
272 
273  //-----------------------//
274  // Miscellaneous //
275  //-----------------------//
276 
277 double Star_rot_diff::omega_funct(double F) const {
278 
279  return omega_frot(F, par_frot) ;
280 
281 }
282 
283 double Star_rot_diff::prim_funct_omega(double F) const {
284 
285  return primfrot(F, par_frot) ;
286 
287 }
288 
290 
291  return omega_field.val_grid_point(0, 0, 0, 0) ;
292 
293 }
294 
295 }
virtual double mass_g() const
Gravitational mass.
virtual ~Star_rot_diff()
Destructor.
void operator=(const Star_rot_diff &)
Assignment to another Star_rot_diff.
Class for differentially rotating stars in quasi-isotropic gauge and maximal slicing.
Definition: star_rot_diff.h:39
Scalar prim_field
Field .
Definition: star_rot_diff.h:79
Map & mp
Mapping associated with the star.
Definition: star.h:180
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
Tbl par_frot
Parameters of the function .
Definition: star_rot_diff.h:70
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
Base class for coordinate mappings.
Definition: map.h:697
void fait_prim_field()
Computes the member prim_field from omega_field .
double omega_funct(double F) const
Evaluates , where F is linked to the components of the fluid 4-velocity by .
double omega
Rotation angular velocity ([f_unit] )
Definition: star_rot.h:113
Scalar nbar
Baryon density in the fluid frame.
Definition: star.h:192
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:98
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
void operator=(const Star_rot &)
Assignment to another Star_rot.
Definition: star_rot.C:376
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star_rot.C:445
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
virtual void sauve(FILE *) const
Save in a file.
double prim_funct_omega(double F) const
Evaluates the primitive of .
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
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Scalar press
Fluid pressure.
Definition: star.h:194
virtual double mass_b() const
Baryon mass.
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
virtual double angu_mom() const
Angular momentum.
Polytropic equation of state (relativistic case).
Definition: eos.h:812
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: star_rot.C:675
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
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
virtual double tsw() const
Ratio T/W.
double omega_max
Maximum value of .
Definition: star_rot_diff.h:76
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
Star_rot_diff(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, double(*omega_frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
Definition: star_rot_diff.C:44
virtual void display_poly(ostream &) const
Display in polytropic units.
double(* omega_frot)(double, const Tbl &)
Function defining the rotation profile.
Definition: star_rot_diff.h:53
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
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
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_min
Minimum value of .
Definition: star_rot_diff.h:75
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Definition: star_rot_diff.h:61
virtual double aplat() const
Flatening r_pole/r_eq.
Scalar omega_field
Field .
Definition: star_rot_diff.h:73