LORENE
et_rot_diff.C
1 /*
2  * Methods for class Et_rot_diff.
3  *
4  * (see file et_rot_diff.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2001 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 as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 
31 
32 /*
33  * $Id: et_rot_diff.C,v 1.5 2016/12/05 16:17:53 j_novak Exp $
34  * $Log: et_rot_diff.C,v $
35  * Revision 1.5 2016/12/05 16:17:53 j_novak
36  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
37  *
38  * Revision 1.4 2014/10/13 08:52:57 j_novak
39  * Lorene classes and functions now belong to the namespace Lorene.
40  *
41  * Revision 1.3 2004/03/25 10:29:05 j_novak
42  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
43  *
44  * Revision 1.2 2001/12/04 21:27:53 e_gourgoulhon
45  *
46  * All writing/reading to a binary file are now performed according to
47  * the big endian convention, whatever the system is big endian or
48  * small endian, thanks to the functions fwrite_be and fread_be
49  *
50  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
51  * LORENE
52  *
53  * Revision 1.3 2001/10/25 09:20:54 eric
54  * Ajout de la fonction virtuelle display_poly.
55  * Affichage de Max nbar, ener et press.
56  *
57  * Revision 1.2 2001/10/24 16:23:01 eric
58  * *** empty log message ***
59  *
60  * Revision 1.1 2001/10/19 08:18:10 eric
61  * Initial revision
62  *
63  *
64  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_diff.C,v 1.5 2016/12/05 16:17:53 j_novak Exp $
65  *
66  */
67 
68 // Headers C
69 #include "math.h"
70 
71 // Headers Lorene
72 #include "et_rot_diff.h"
73 #include "eos.h"
74 #include "nbr_spx.h"
75 #include "utilitaires.h"
76 #include "unites.h"
77 
78  //--------------//
79  // Constructors //
80  //--------------//
81 
82 // Standard constructor
83 // --------------------
84 namespace Lorene {
85 Et_rot_diff::Et_rot_diff(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i,
86  double (*frot_i)(double, const Tbl&),
87  double (*primfrot_i)(double, const Tbl&),
88  const Tbl& par_frot_i)
89  : Etoile_rot(mp_i, nzet_i, relat, eos_i),
90  frot(frot_i),
91  primfrot(primfrot_i),
92  par_frot(par_frot_i),
93  omega_field(mp_i),
94  prim_field(mp_i)
95  {
96 
97  // To make sure that omega is not used
98  omega = __infinity ;
99 
100  // Initialization to a static state :
101  omega_field = 0 ;
102  prim_field = 0 ;
103  omega_min = 0 ;
104  omega_max = 0 ;
105 
106 }
107 
108 // Copy constructor
109 // ----------------
111  : Etoile_rot(et),
112  frot(et.frot),
113  primfrot(et.primfrot),
114  par_frot(et.par_frot),
115  omega_field(et.omega_field),
116  omega_min(et.omega_min),
117  omega_max(et.omega_max),
118  prim_field(et.prim_field)
119  {}
120 
121 
122 // Constructor from a file
123 // -----------------------
124 Et_rot_diff::Et_rot_diff(Map& mp_i, const Eos& eos_i, FILE* fich,
125  double (*frot_i)(double, const Tbl&),
126  double (*primfrot_i)(double, const Tbl&) )
127  : Etoile_rot(mp_i, eos_i, fich),
128  frot(frot_i),
129  primfrot(primfrot_i),
130  par_frot(fich),
131  omega_field(mp_i),
132  prim_field(mp_i)
133  {
134 
135  Tenseur omega_field_file(mp, fich) ;
136  omega_field = omega_field_file ;
137  fait_prim_field() ;
138 
139  // omega_min and omega_max are read in the file:
140  fread_be(&omega_min, sizeof(double), 1, fich) ;
141  fread_be(&omega_max, sizeof(double), 1, fich) ;
142 
143 }
144 
145 
146  //------------//
147  // Destructor //
148  //------------//
149 
151 
152 
153  //--------------//
154  // Assignment //
155  //--------------//
156 
157 // Assignment to another Et_rot_diff
158 // ---------------------------------
160 
161  // Assignment of quantities common to all the derived classes of Etoile_rot
162  Etoile_rot::operator=(et) ;
163 
164  // Assignment of proper quantities of class Etoile_rot
165  frot = et.frot ;
166  primfrot = et.primfrot ;
167  par_frot = et.par_frot ;
168  omega_field = et.omega_field ;
169  prim_field = et.prim_field ;
170  omega_min = et.omega_min ;
171  omega_max = et.omega_max ;
172 
173 }
174 
175  //--------------//
176  // Outputs //
177  //--------------//
178 
179 // Save in a file
180 // --------------
181 
182 void Et_rot_diff::sauve(FILE* fich) const {
183 
184  Etoile_rot::sauve(fich) ;
185 
186  par_frot.sauve(fich) ;
187 
188  omega_field.sauve(fich) ;
189 
190  fwrite_be(&omega_min, sizeof(double), 1, fich) ;
191  fwrite_be(&omega_max, sizeof(double), 1, fich) ;
192 
193 }
194 
195 
196 // Printing
197 // --------
198 
199 ostream& Et_rot_diff::operator>>(ostream& ost) const {
200 
201  using namespace Unites ;
202 
203  Etoile_rot::operator>>(ost) ;
204 
205  ost << endl ;
206  ost << "Differentially rotating star" << endl ;
207  ost << "-----------------------------" << endl ;
208 
209  ost << endl << "Parameters of F(Omega) : " << endl ;
210  ost << par_frot << endl ;
211 
212  ost << "Min, Max of Omega/(2pi) : " << omega_min / (2*M_PI) * f_unit
213  << " Hz, " << omega_max / (2*M_PI) * f_unit << " Hz" << endl ;
214  int lsurf = nzet - 1;
215  int nt = mp.get_mg()->get_nt(lsurf) ;
216  int nr = mp.get_mg()->get_nr(lsurf) ;
217  ost << "Central, equatorial value of Omega: "
218  << omega_field()(0, 0, 0, 0) * f_unit << " rad/s, "
219  << omega_field()(nzet-1, 0, nt-1, nr-1) * f_unit << " rad/s" << endl ;
220 
221  ost << "Central, equatorial value of Omega/(2 Pi): "
222  << omega_field()(0, 0, 0, 0) * f_unit / (2*M_PI) << " Hz, "
223  << omega_field()(nzet-1, 0, nt-1, nr-1) * f_unit / (2*M_PI)
224  << " Hz" << endl ;
225 
226  double nbar_max = max( max( nbar() ) ) ;
227  double ener_max = max( max( ener() ) ) ;
228  double press_max = max( max( press() ) ) ;
229  ost << "Max prop. bar. dens. : " << nbar_max
230  << " x 0.1 fm^-3 = " << nbar_max / nbar()(0, 0, 0, 0) << " central"
231  << endl ;
232  ost << "Max prop. ener. dens. (e_max) : " << ener_max
233  << " rho_nuc c^2 = " << ener_max / ener()(0, 0, 0, 0) << " central"
234  << endl ;
235  ost << "Max pressure : " << press_max
236  << " rho_nuc c^2 = " << press_max / press()(0, 0, 0, 0) << " central"
237  << endl ;
238 
239  // Length scale set by the maximum energy density:
240  double len_rho = 1. / sqrt( ggrav * ener_max ) ;
241  ost << endl << "Value of A = par_frot(1) in units of e_max^{1/2} : " <<
242  par_frot(1) / len_rho << endl ;
243 
244  ost << "Value of A / r_eq : " <<
245  par_frot(1) / ray_eq() << endl ;
246 
247  ost << "r_p/r_eq : " << aplat() << endl ;
248  ost << "KEH l^2 = (c/G^2)^{2/3} J^2 e_max^{1/3} M_B^{-10/3} : " <<
249  angu_mom() * angu_mom() / pow(len_rho, 0.6666666666666666)
250  / pow(mass_b(), 3.3333333333333333)
251  / pow(ggrav, 1.3333333333333333) << endl ;
252 
253  ost << "M e_max^{1/2} : " << ggrav * mass_g() / len_rho << endl ;
254 
255  ost << "r_eq e_max^{1/2} : " << ray_eq() / len_rho << endl ;
256 
257  ost << "T/W : " << tsw() << endl ;
258 
259  ost << "Omega_c / e_max^{1/2} : " << get_omega_c() * len_rho << endl ;
260 
261  display_poly(ost) ;
262 
263  return ost ;
264 
265 
266 }
267 
268 // display_poly
269 // ------------
270 
271 void Et_rot_diff::display_poly(ostream& ost) const {
272 
273  using namespace Unites ;
274 
275  Etoile_rot::display_poly( ost ) ;
276 
277  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
278 
279  if (p_eos_poly != 0x0) {
280 
281  double kappa = p_eos_poly->get_kap() ;
282  double gamma = p_eos_poly->get_gam() ; ;
283 
284  // kappa^{n/2}
285  double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
286 
287  // Polytropic unit of length in terms of r_unit :
288  double r_poly = kap_ns2 / sqrt(ggrav) ;
289 
290  // Polytropic unit of time in terms of t_unit :
291  double t_poly = r_poly ;
292 
293  // Polytropic unit of density in terms of rho_unit :
294  double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
295 
296  ost.precision(10) ;
297  ost << " n_max : " << max( max( nbar() ) ) / rho_poly << endl ;
298  ost << " e_max : " << max( max( ener() ) ) / rho_poly << endl ;
299  ost << " P_min : " << 2.*M_PI / omega_max / t_poly << endl ;
300  ost << " P_max : " << 2.*M_PI / omega_min / t_poly << endl ;
301 
302  int lsurf = nzet - 1;
303  int nt = mp.get_mg()->get_nt(lsurf) ;
304  int nr = mp.get_mg()->get_nr(lsurf) ;
305  ost << " P_eq : " << 2.*M_PI /
306  omega_field()(nzet-1, 0, nt-1, nr-1) / t_poly << endl ;
307 
308  }
309 
310 }
311 
312 
313 
314  //-----------------------//
315  // Miscellaneous //
316  //-----------------------//
317 
318 double Et_rot_diff::funct_omega(double omeg) const {
319 
320  return frot(omeg, par_frot) ;
321 
322 }
323 
324 double Et_rot_diff::prim_funct_omega(double omeg) const {
325 
326  return primfrot(omeg, par_frot) ;
327 
328 }
329 
330 double Et_rot_diff::get_omega_c() const {
331 
332  return omega_field()(0, 0, 0, 0) ;
333 
334 }
335 
336 
337 
338 }
Et_rot_diff(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, double(*frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
Definition: et_rot_diff.C:85
void operator=(const Et_rot_diff &)
Assignment to another Et_rot_diff.
Definition: et_rot_diff.C:159
Tenseur prim_field
Field .
Definition: et_rot_diff.h:114
double prim_funct_omega(double omeg) const
Evaluates the primitive of , where F is the function defining the rotation profile.
Definition: et_rot_diff.C:324
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: etoile_rot.C:708
Lorene prototypes.
Definition: app_hor.h:67
double(* frot)(double, const Tbl &)
Function defining the rotation profile.
Definition: et_rot_diff.h:88
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:777
double ray_eq() const
Coordinate radius at , [r_unit].
Base class for coordinate mappings.
Definition: map.h:682
Class for isolated rotating stars *** DEPRECATED : use class Star_rot instead ***.
Definition: etoile.h:1499
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Definition: et_rot_diff.h:96
void operator=(const Etoile_rot &)
Assignment to another Etoile_rot.
Definition: etoile_rot.C:415
Tenseur press
Fluid pressure.
Definition: etoile.h:464
double funct_omega(double omeg) const
Evaluates , where F is the function defining the rotation profile.
Definition: et_rot_diff.C:318
double omega_max
Maximum value of .
Definition: et_rot_diff.h:111
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Definition: eos_poly.C:271
Class for differentially rotating stars.
Definition: et_rot_diff.h:73
void fait_prim_field()
Computes the member prim_field from omega_field .
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition: et_rot_diff.C:271
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: etoile_rot.C:477
virtual double angu_mom() const
Angular momentum.
void sauve(FILE *) const
Save in a file.
Definition: tenseur.C:1341
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:462
virtual void sauve(FILE *) const
Save in a file.
Definition: etoile_rot.C:452
Tbl par_frot
Parameters of the function .
Definition: et_rot_diff.h:105
virtual double mass_b() const
Baryon mass.
Polytropic equation of state (relativistic case).
Definition: eos.h:812
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
const Eos & eos
Equation of state of the stellar matter.
Definition: etoile.h:454
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
double omega
Rotation angular velocity ([f_unit] )
Definition: etoile.h:1504
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:435
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:469
virtual void sauve(FILE *) const
Save in a file.
Definition: et_rot_diff.C:182
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:463
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: et_rot_diff.C:199
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:474
double omega_min
Minimum value of .
Definition: et_rot_diff.h:110
virtual double mass_g() const
Gravitational mass.
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] )
Definition: et_rot_diff.C:330
virtual ~Et_rot_diff()
Destructor.
Definition: et_rot_diff.C:150
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
Tenseur omega_field
Field .
Definition: et_rot_diff.h:108
virtual double aplat() const
Flatening r_pole/r_eq.
virtual double tsw() const
Ratio T/W.