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  // Overwriting the call on construction in Star_rot_CFC, which uses omega instead of omega_field
101  hydro_euler() ;
102 
103  // Recomputing \hat{A}^{ij} and its norm, necessary for correct GRV2 and GRV3
104  //------------------------------------------
105  Vector sou_Xi = 2*Unites::qpig*psi4*psi4*psi2*j_euler ;
106  double lambda = 1./3. ;
107  Vector Xi_new = sou_Xi.poisson(lambda) ;
108 
109  hatA = Xi_new.ope_killing_conf(flat) ;
110  hatA_quad = contract(hatA, 0, 1, hatA.up_down(flat), 0, 1) ;
111 
112 }
113 
114  //------------//
115  // Destructor //
116  //------------//
117 
119 
120  //---------------//
121  // Assignment //
122  //---------------//
123 
124 // Assignment to another Star_rot_Dirac_diff
125 // ------------------------------------
126 
128 
129  // Assignment of quantities common to all the derived classes of
130  // Star_rot_CFC
132 
133  // Assignment of proper quantities of class Star_rot
134  omega_frot = star.omega_frot ;
135  primfrot = star.primfrot ;
136  par_frot = star.par_frot ;
137  omega_field = star.omega_field ;
138  prim_field = star.prim_field ;
139  omega_min = star.omega_min ;
140  omega_max = star.omega_max ;
141 
142 }
143 
144 
145  //--------------//
146  // Outputs //
147  //--------------//
148 
149 // Save in a file
150 // --------------
151 
152 void Star_rot_diff_CFC::sauve(FILE* fich) const {
153 
154  Star_rot_CFC::sauve(fich) ;
155 
156  par_frot.sauve(fich) ;
157 
158  omega_field.sauve(fich) ;
159 
160  fwrite_be(&omega_min, sizeof(double), 1, fich) ;
161  fwrite_be(&omega_max, sizeof(double), 1, fich) ;
162 
163 }
164 
165 
166 // Printing
167 // --------
168 
169 ostream& Star_rot_diff_CFC::operator>>(ostream& ost) const {
170 
171  using namespace Unites ;
172 
174 
175  // TODO: check this section and figure out what to output
176 
177  ost << endl ;
178  ost << "Differentially rotating star" << endl ;
179  ost << "-----------------------------" << endl ;
180 
181  ost << endl << "Parameters of Omega(F) : " << endl ;
182  ost << par_frot << endl ;
183 
184  ost << "Min, Max of Omega/(2pi) : " << omega_min / (2*M_PI) * f_unit
185  << " Hz, " << omega_max / (2*M_PI) * f_unit << " Hz" << endl ;
186  int lsurf = nzet - 1;
187  int nt = mp.get_mg()->get_nt(lsurf) ;
188  int nr = mp.get_mg()->get_nr(lsurf) ;
189  ost << "Central, equatorial value of Omega: "
190  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit << " rad/s, "
191  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit << " rad/s" << endl ;
192 
193  ost << "Central, equatorial value of Omega/(2 Pi): "
194  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit / (2*M_PI) << " Hz, "
195  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit / (2*M_PI)
196  << " Hz" << endl ;
197 
198  double nbar_max = max( max( nbar ) ) ;
199  double ener_max = max( max( ener ) ) ;
200  double press_max = max( max( press ) ) ;
201  ost << "Max prop. bar. dens. : " << nbar_max
202  << " x 0.1 fm^-3 = " << nbar_max / nbar.val_grid_point(0, 0, 0, 0) << " central"
203  << endl ;
204  ost << "Max prop. ener. dens. (e_max) : " << ener_max
205  << " rho_nuc c^2 = " << ener_max / ener.val_grid_point(0, 0, 0, 0) << " central"
206  << endl ;
207  ost << "Max pressure : " << press_max
208  << " rho_nuc c^2 = " << press_max / press.val_grid_point(0, 0, 0, 0) << " central"
209  << endl ;
210 
211  // Length scale set by the maximum energy density:
212  double len_rho = 1. / sqrt( ggrav * ener_max ) ;
213  ost << endl << "Value of A = par_frot(1) in units of e_max^{1/2} : " <<
214  par_frot(1) / len_rho << endl ;
215 
216  ost << "Value of A / r_eq : " <<
217  par_frot(1) / ray_eq() << endl ;
218 
219  ost << "r_p/r_eq : " << aplat() << endl ;
220  ost << "KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
221  angu_mom() * angu_mom() / pow(len_rho, 0.6666666666666666)
222  / pow(mass_b(), 3.3333333333333333)
223  / pow(ggrav, 1.3333333333333333) << endl ;
224 
225  ost << "M e_max^{1/2} : " << ggrav * mass_g() / len_rho << endl ;
226 
227  ost << "r_eq e_max^{1/2} : " << ray_eq() / len_rho << endl ;
228 
229  ost << "T/W : " << tsw() << endl ;
230 
231  ost << "Omega_c / e_max^{1/2} : " << get_omega_c() * len_rho << endl ;
232 
233  display_poly(ost) ;
234 
235  return ost ;
236 
237 }
238 
239 // display_poly
240 // ------------
241 
242 void Star_rot_diff_CFC::display_poly(ostream& ost) const {
243 
244  using namespace Unites ;
245 
246  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
247 
248  if (p_eos_poly != 0x0) {
249 
250  double kappa = p_eos_poly->get_kap() ;
251 
252  // kappa^{n/2}
253  double kap_ns2 = pow( kappa, 0.5 /(p_eos_poly->get_gam()-1) ) ;
254 
255  // Polytropic unit of length in terms of r_unit :
256  double r_poly = kap_ns2 / sqrt(ggrav) ;
257 
258  // Polytropic unit of time in terms of t_unit :
259  double t_poly = r_poly ;
260 
261  // Polytropic unit of mass in terms of m_unit :
262  double m_poly = r_poly / ggrav ;
263 
264  // Polytropic unit of angular momentum in terms of j_unit :
265  double j_poly = r_poly * r_poly / ggrav ;
266 
267  // Polytropic unit of density in terms of rho_unit :
268  double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
269 
270  ost.precision(10) ;
271  ost << endl << "Quantities in polytropic units : " << endl ;
272  ost << "==============================" << endl ;
273  ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
274  ost << " n_c : " << nbar.val_grid_point(0, 0, 0, 0) / rho_poly << endl ;
275  ost << " e_c : " << ener.val_grid_point(0, 0, 0, 0) / rho_poly << endl ;
276  ost << " Omega_c : " << get_omega_c() * t_poly << endl ;
277  ost << " P_c : " << 2.*M_PI / get_omega_c() / t_poly << endl ;
278  ost << " M_bar : " << mass_b() / m_poly << endl ;
279  ost << " M : " << mass_g() / m_poly << endl ;
280  ost << " J : " << angu_mom() / j_poly << endl ;
281  ost << " r_eq : " << ray_eq() / r_poly << endl ;
282  ost << " R_circ_eq : " << r_circ() / r_poly << endl ;
283 
284  ost << " n_max : " << max( max( nbar ) ) / rho_poly << endl ;
285  ost << " e_max : " << max( max( ener ) ) / rho_poly << endl ;
286  ost << " P_min : " << 2.*M_PI / omega_max / t_poly << endl ;
287  ost << " P_max : " << 2.*M_PI / omega_min / t_poly << endl ;
288 
289  int lsurf = nzet - 1;
290  int nt = mp.get_mg()->get_nt(lsurf) ;
291  int nr = mp.get_mg()->get_nr(lsurf) ;
292  ost << " P_eq : " << 2.*M_PI /
293  omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) / t_poly << endl ;
294 
295  }
296  }
297 
298 
299 
300  //-----------------------//
301  // Miscellaneous //
302  //-----------------------//
303 
304 double Star_rot_diff_CFC::omega_funct(double F) const {
305 
306  return omega_frot(F, par_frot) ;
307 
308 }
309 
310 double Star_rot_diff_CFC::prim_funct_omega(double F) const {
311 
312  return primfrot(F, par_frot) ;
313 
314 }
315 
317 
318  return omega_field.val_grid_point(0, 0, 0, 0) ;
319 
320 }
321 
322 void Star_rot_diff_CFC::partial_display(ostream& ost) const {
323 
324  using namespace Unites ;
325 
326  double omega_c = get_omega_c() ;
327  double freq = omega_c / (2.*M_PI) ;
328  ost << "Central Omega : " << omega_c * f_unit
329  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
330  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
331  << endl ;
332  ost << endl << "Central enthalpy : " << ent.val_grid_point(0,0,0,0) << " c^2" << endl ;
333  ost << "Central proper baryon density : " << nbar.val_grid_point(0,0,0,0)
334  << " x 0.1 fm^-3" << endl ;
335  ost << "Central proper energy density : " << ener.val_grid_point(0,0,0,0)
336  << " rho_nuc c^2" << endl ;
337  ost << "Central pressure : " << press.val_grid_point(0,0,0,0)
338  << " rho_nuc c^2" << endl ;
339 
340  ost << "Central value of log(N) : "
341  << logn.val_grid_point(0, 0, 0, 0) << endl ;
342  ost << "Central lapse N : " << nn.val_grid_point(0,0,0,0) << endl ;
343  ost << "Central value of psi : " << psi.val_grid_point(0,0,0,0) << endl ;
344 
345  double nphi_c = beta(3).val_grid_point(0, 0, 0, 0) ;
346  if ( (omega_c==0) && (nphi_c==0) ) {
347  ost << "Central azimuthal shift : " << nphi_c << endl ;
348  }
349  else{
350  ost << "Central azimuthal shift /Omega : " << nphi_c / omega_c
351  << endl ;
352  }
353 
354 
355  int lsurf = nzet - 1;
356  int nt = mp.get_mg()->get_nt(lsurf) ;
357  int nr = mp.get_mg()->get_nr(lsurf) ;
358  ost << "Equatorial value of the velocity U: "
359  << u_euler(3).val_grid_point(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
360 
361  ost << endl
362  << "Coordinate equatorial radius r_eq = "
363  << ray_eq()/km << " km" << endl ;
364  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
365 
366 }
367 
368 }
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.
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
void operator=(const Star_rot_CFC &)
Assignment to another Star_rot_CFC.
Definition: star_rot_cfc.C:212
Scalar psi4
Conformal factor .
Definition: star_rot_cfc.h:81
Class for differentially rotating stars in Conformal Flatness Condition and maximal slicing...
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
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
virtual void sauve(FILE *) const
Save in a file.
Definition: star_rot_cfc.C:243
Vector beta
Shift vector.
Definition: star.h:228
Tensor field of valence 1.
Definition: vector.h:188
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
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
Definition: vector.C:468
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
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
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
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
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 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 .
Vector j_euler
Momentum density 3-vector with respect to the Eulerian observer.
Definition: star_rot_cfc.h:75
void fait_prim_field()
Computes the member prim_field from omega_field .
virtual double mass_g() const
Gravitational mass.
const Metric_flat & flat
flat metric (spherical components)
Definition: star_rot_cfc.h:85