LORENE
hot_star_rot_diff_cfc.C
1 /*
2  * Methods of class Hot_star_rot_diff_CFC
3  *
4  * (see file hot_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 "hot_star_rot_diff_cfc.h"
30 #include "hoteos.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 Hot_star_rot_diff_CFC::Hot_star_rot_diff_CFC(Map& mpi, int nzet_i, const Hot_eos& heos_i, int relat_i,
44  int e_type_i, double const_value_i, int filter,
45  double (*omega_frot_i)(double, const Tbl&),
46  double (*primfrot_i)(double, const Tbl&),
47  const Tbl& par_frot_i)
48  :Hot_star_rot_CFC(mpi, nzet_i, heos_i, relat_i, e_type_i, const_value_i, filter),
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  : Hot_star_rot_CFC(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 //Constructor from a file
81 //------------------------
83  double (*omega_frot_i)(double, const Tbl&),
84  double (*primfrot_i)(double, const Tbl&) )
85  : Hot_star_rot_CFC(mpi, heos_i, fich),
86  omega_frot(omega_frot_i),
87  primfrot(primfrot_i),
88  par_frot(fich),
89  omega_field(mpi),
90  prim_field(mpi)
91  {
92 
93  Scalar omega_field_file(mp, *(mp.get_mg()), fich) ;
94  omega_field = omega_field_file ;
95  fait_prim_field() ;
96 
97  // omega_min and omega_max are read in the file:
98  fread_be(&omega_min, sizeof(double), 1, fich) ;
99  fread_be(&omega_max, sizeof(double), 1, fich) ;
100 
101  // Overwriting the call on construction in Hot_star_rot_CFC, which uses omega instead of omega_field
102  hydro_euler() ;
103 
104  // Recomputing \hat{A}^{ij} and its norm, necessary for correct GRV2 and GRV3
105  //------------------------------------------
106  Vector sou_Xi = 2*Unites::qpig*psi4*psi4*psi2*j_euler ;
107  double lambda = 1./3. ;
108  Vector Xi_new = sou_Xi.poisson(lambda) ;
109 
110  hatA = Xi_new.ope_killing_conf(flat) ;
111  hatA_quad = contract(hatA, 0, 1, hatA.up_down(flat), 0, 1) ;
112 
113 }
114 
115  //------------//
116  // Destructor //
117  //------------//
118 
120 
121  //---------------//
122  // Assignment //
123  //---------------//
124 
125 // Assignment to another Hot_star_rot_Dirac_diff
126 // ------------------------------------
127 
129 
130  // Assignment of quantities common to all the derived classes of
131  // Hot_star_rot_CFC
133 
134  // Assignment of proper quantities of class Hot_star_rot_diff_CFC
135  omega_frot = star.omega_frot ;
136  primfrot = star.primfrot ;
137  par_frot = star.par_frot ;
138  omega_field = star.omega_field ;
139  prim_field = star.prim_field ;
140  omega_min = star.omega_min ;
141  omega_max = star.omega_max ;
142 
143 }
144 
145 
146  //--------------//
147  // Outputs //
148  //--------------//
149 
150 // Save in a file
151 // --------------
152 
153 void Hot_star_rot_diff_CFC::sauve(FILE* fich) const {
154 
156 
157  par_frot.sauve(fich) ;
158 
159  omega_field.sauve(fich) ;
160 
161  fwrite_be(&omega_min, sizeof(double), 1, fich) ;
162  fwrite_be(&omega_max, sizeof(double), 1, fich) ;
163 
164 }
165 
166 
167 // Printing
168 // --------
169 
170 ostream& Hot_star_rot_diff_CFC::operator>>(ostream& ost) const {
171 
172  using namespace Unites ;
173 
175 
176  // TODO: check this section and figure out what to output
177 
178  ost << endl ;
179  ost << "Differentially rotating star" << endl ;
180  ost << "-----------------------------" << endl ;
181 
182  ost << endl << "Parameters of Omega(F) : " << endl ;
183  ost << par_frot << endl ;
184 
185  ost << "Min, Max of Omega/(2pi) : " << omega_min / (2*M_PI) * f_unit
186  << " Hz, " << omega_max / (2*M_PI) * f_unit << " Hz" << endl ;
187  int lsurf = nzet - 1;
188  int nt = mp.get_mg()->get_nt(lsurf) ;
189  int nr = mp.get_mg()->get_nr(lsurf) ;
190  ost << "Central, equatorial value of Omega: "
191  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit << " rad/s, "
192  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit << " rad/s" << endl ;
193 
194  ost << "Central, equatorial value of Omega/(2 Pi): "
195  << omega_field.val_grid_point(0, 0, 0, 0) * f_unit / (2*M_PI) << " Hz, "
196  << omega_field.val_grid_point(nzet-1, 0, nt-1, nr-1) * f_unit / (2*M_PI)
197  << " Hz" << endl ;
198 
199  double nbar_max = max( max( nbar ) ) ;
200  double ener_max = max( max( ener ) ) ;
201  double press_max = max( max( press ) ) ;
202  ost << "Max prop. bar. dens. : " << nbar_max
203  << " x 0.1 fm^-3 = " << nbar_max / nbar.val_grid_point(0, 0, 0, 0) << " central"
204  << endl ;
205  ost << "Max prop. ener. dens. (e_max) : " << ener_max
206  << " rho_nuc c^2 = " << ener_max / ener.val_grid_point(0, 0, 0, 0) << " central"
207  << endl ;
208  ost << "Max pressure : " << press_max
209  << " rho_nuc c^2 = " << press_max / press.val_grid_point(0, 0, 0, 0) << " central"
210  << endl ;
211 
212  // Length scale set by the maximum energy density:
213  double len_rho = 1. / sqrt( ggrav * ener_max ) ;
214  ost << endl << "Value of A = par_frot(1) in units of e_max^{1/2} : " <<
215  par_frot(1) / len_rho << endl ;
216 
217  ost << "Value of A / r_eq : " <<
218  par_frot(1) / ray_eq() << endl ;
219 
220  ost << "r_p/r_eq : " << aplat() << endl ;
221  ost << "KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
222  angu_mom() * angu_mom() / pow(len_rho, 0.6666666666666666)
223  / pow(mass_b(), 3.3333333333333333)
224  / pow(ggrav, 1.3333333333333333) << endl ;
225 
226  ost << "M e_max^{1/2} : " << ggrav * mass_g() / len_rho << endl ;
227 
228  ost << "r_eq e_max^{1/2} : " << ray_eq() / len_rho << endl ;
229 
230  ost << "T/W : " << tsw() << endl ;
231 
232  ost << "Omega_c / e_max^{1/2} : " << get_omega_c() * len_rho << endl ;
233 
234  return ost ;
235 
236 }
237 
238 
239  //-----------------------//
240  // Miscellaneous //
241  //-----------------------//
242 
243 double Hot_star_rot_diff_CFC::omega_funct(double F) const {
244 
245  return omega_frot(F, par_frot) ;
246 
247 }
248 
250 
251  return primfrot(F, par_frot) ;
252 
253 }
254 
256 
257  return omega_field.val_grid_point(0, 0, 0, 0) ;
258 
259 }
260 
261 void Hot_star_rot_diff_CFC::partial_display(ostream& ost) const {
262 
263  using namespace Unites ;
264 
265  double omega_c = get_omega_c() ;
266  double freq = omega_c / (2.*M_PI) ;
267  ost << "Central Omega : " << omega_c * f_unit
268  << " rad/s f : " << freq * f_unit << " Hz" << endl ;
269  ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
270  << endl ;
271  ost << endl << "Central enthalpy : " << ent.val_grid_point(0,0,0,0) << " c^2" << endl ;
272  ost << "Central proper baryon density : " << nbar.val_grid_point(0,0,0,0)
273  << " x 0.1 fm^-3" << endl ;
274  ost << "Central proper energy density : " << ener.val_grid_point(0,0,0,0)
275  << " rho_nuc c^2" << endl ;
276  ost << "Central pressure : " << press.val_grid_point(0,0,0,0)
277  << " rho_nuc c^2" << endl ;
278 
279  ost << "Central value of log(N) : "
280  << logn.val_grid_point(0, 0, 0, 0) << endl ;
281  ost << "Central lapse N : " << nn.val_grid_point(0,0,0,0) << endl ;
282  ost << "Central value of psi : " << psi.val_grid_point(0,0,0,0) << endl ;
283 
284  double nphi_c = beta(3).val_grid_point(0, 0, 0, 0) ;
285  if ( (omega_c==0) && (nphi_c==0) ) {
286  ost << "Central azimuthal shift : " << nphi_c << endl ;
287  }
288  else{
289  ost << "Central azimuthal shift /Omega : " << nphi_c / omega_c
290  << endl ;
291  }
292 
293 
294  int lsurf = nzet - 1;
295  int nt = mp.get_mg()->get_nt(lsurf) ;
296  int nr = mp.get_mg()->get_nr(lsurf) ;
297  ost << "Equatorial value of the velocity U: "
298  << u_euler(3).val_grid_point(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
299 
300  ost << endl
301  << "Coordinate equatorial radius r_eq = "
302  << ray_eq()/km << " km" << endl ;
303  ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
304 
305 }
306 
307 }
const Metric_flat & flat
flat metric (spherical components)
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
Scalar nbar
Baryon density in the fluid frame.
Definition: hot_star.h:93
Scalar psi
Conformal factor .
Hot_star_rot_diff_CFC(Map &mp_i, int nzet_i, const Hot_eos &heos_i, int relat_i, int e_type, double const_value_i, int filter, double(*omega_frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
virtual double aplat() const
Flattening r_pole/r_eq.
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:402
Tbl par_frot
Parameters of the function .
Base class for coordinate mappings.
Definition: map.h:696
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Class for hot differentially rotating stars in Conformal Flatness Condition and maximal slicing...
void operator=(const Hot_star_rot_diff_CFC &)
Assignment to another Hot_star_rot_diff_CFC.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Scalar nn
Lapse function N .
Definition: hot_star.h:130
Tensor field of valence 1.
Definition: vector.h:188
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:652
Scalar logn
Logarithm of the lapse N .
Definition: hot_star.h:127
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
Definition: vector.C:468
double ray_eq() const
Coordinate radius at , [r_unit].
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
int nzet
Number of domains of *mp occupied by the star.
Definition: hot_star.h:84
virtual void sauve(FILE *) const
Save in a file.
virtual double mass_b() const
Baryonic mass.
double omega
Rotation angular velocity ([f_unit] )
void operator=(const Hot_star_rot_CFC &)
Assignment to another Hot_star_rot_CFC.
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
double omega_funct(double F) const
Evaluates , where F is linked to the components of the fluid 4-velocity by .
virtual double mass_g() const
Gravitational mass.
Scalar ener
Total energy density in the fluid frame.
Definition: hot_star.h:94
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
Scalar ent
Log-enthalpy.
Definition: hot_star.h:91
Map & mp
Mapping associated with the star.
Definition: hot_star.h:81
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Scalar psi4
Conformal factor .
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
double(* omega_frot)(double, const Tbl &)
Function defining the rotation profile.
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
virtual double angu_mom() const
Angular momentum.
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Class for relativistic rotating stars in Conformal Flatness Condition and maximal slicing...
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
Base class for 2-parameters equations of state (abstract class).
Definition: hoteos.h:85
double omega_max
Maximum value of .
virtual void sauve(FILE *) const
Save in a file.
virtual double tsw() const
Ratio T/W.
void fait_prim_field()
Computes the member prim_field from omega_field .
void sauve(FILE *) const
Save in a file.
Definition: tbl.C:329
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 press
Fluid pressure.
Definition: hot_star.h:95
Vector j_euler
Momentum density 3-vector with respect to the Eulerian observer.
double omega_min
Minimum value of .
virtual ~Hot_star_rot_diff_CFC()
Destructor.
Vector beta
Shift vector.
Definition: hot_star.h:133
double prim_funct_omega(double F) const
Evaluates the primitive of .
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: hot_star.h:112