LORENE
etoile_rot.C
1/*
2 * Methods for the class Etoile_rot
3 *
4 * (see file etoile.h for documentation)
5 */
6
7/*
8 * Copyright (c) 2000-2001 Eric Gourgoulhon
9 *
10 * This file is part of LORENE.
11 *
12 * LORENE is free software; you can redistribute it and/or modify
13 * it under the terms of the GNU General Public License as published by
14 * the Free Software Foundation; either version 2 of the License, or
15 * (at your option) any later version.
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
30
31/*
32 * $Id: etoile_rot.C,v 1.10 2023/07/04 08:59:53 j_novak Exp $
33 * $Log: etoile_rot.C,v $
34 * Revision 1.10 2023/07/04 08:59:53 j_novak
35 * Added method "r_circ_merid()" to compute the circumferential meridional radius.
36 *
37 * Revision 1.9 2021/04/13 11:24:53 j_novak
38 * Corrected a bug in Etoile_rot::fait_shift() which was missing the np=1 case.
39 *
40 * Revision 1.8 2016/12/05 16:17:55 j_novak
41 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
42 *
43 * Revision 1.7 2015/12/03 14:17:24 j_novak
44 * Check added for the computation of area (thanks S. Koeppel).
45 *
46 * Revision 1.6 2015/06/10 14:37:44 a_sourie
47 * Corrected the formula for the quadrupole.
48 *
49 * Revision 1.5 2014/10/13 08:52:59 j_novak
50 * Lorene classes and functions now belong to the namespace Lorene.
51 *
52 * Revision 1.4 2004/03/25 10:29:07 j_novak
53 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
54 *
55 * Revision 1.3 2001/12/06 15:11:43 jl_zdunik
56 * Introduction of the new function f_eq() in the class Etoile_rot
57 *
58 * Revision 1.2 2001/12/04 21:27:53 e_gourgoulhon
59 *
60 * All writing/reading to a binary file are now performed according to
61 * the big endian convention, whatever the system is big endian or
62 * small endian, thanks to the functions fwrite_be and fread_be
63 *
64 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
65 * LORENE
66 *
67 * Revision 2.17 2001/10/24 15:36:20 eric
68 * Ajout de la fonction display_poly.
69 *
70 * Revision 2.16 2001/10/16 14:49:02 eric
71 * Appel de get_omega_c() pour avoir la valeur centrale de Omega.
72 * Affichage different si rotation differentielle.
73 *
74 * Revision 2.15 2001/09/13 08:32:01 eric
75 * Ajout du facteur de compacite M/R dans l'affichage.
76 *
77 * Revision 2.14 2001/06/20 14:20:56 novak
78 * Appel a Etoile_rot::set_der0x0 dans del_deriv (au lieu de set_der0x0
79 * tout court).
80 *
81 * Revision 2.13 2001/03/26 09:30:58 jlz
82 * New members p_espec_isco and p_lspec_isco.
83 *
84 * Revision 2.12 2000/11/20 21:42:02 eric
85 * Appel de fait_nphi() dans le constructeur par lecture de fichier.
86 *
87 * Revision 2.11 2000/11/18 23:18:30 eric
88 * Modifs affichage.
89 *
90 * Revision 2.10 2000/11/18 21:09:57 eric
91 * Ajout des membres p_r_isco et p_f_isco.
92 *
93 * Revision 2.9 2000/11/07 16:33:08 eric
94 * Modif affichage.
95 *
96 * Revision 2.8 2000/10/12 15:37:01 eric
97 * Ajout de la fonction fait_nphi().
98 *
99 * Revision 2.7 2000/09/18 16:15:12 eric
100 * Ajout du membre tkij.
101 *
102 * Revision 2.6 2000/08/31 15:38:00 eric
103 * Bases spectrales standards pour bbb et b_car dans le constructeur
104 * standard (initialisation a la metrique plate).
105 *
106 * Revision 2.5 2000/08/31 11:25:45 eric
107 * Ajout des membres tnphi et ak_car.
108 *
109 * Revision 2.4 2000/08/25 12:28:29 eric
110 * Modif affichage.
111 *
112 * Revision 2.3 2000/08/18 14:01:59 eric
113 * Ajout de partial_display
114 *
115 * Revision 2.2 2000/08/17 12:40:04 eric
116 * *** empty log message ***
117 *
118 * Revision 2.1 2000/07/21 16:31:26 eric
119 * *** empty log message ***
120 *
121 * Revision 1.1 2000/07/20 15:32:37 eric
122 * Initial revision
123 *
124 *
125 * $Header: /cvsroot/Lorene/C++/Source/Etoile/etoile_rot.C,v 1.10 2023/07/04 08:59:53 j_novak Exp $
126 *
127 */
128
129// Headers C
130#include "math.h"
131
132// Headers Lorene
133#include "etoile.h"
134#include "eos.h"
135#include "nbr_spx.h"
136#include "utilitaires.h"
137#include "unites.h"
138
139 //--------------//
140 // Constructors //
141 //--------------//
142
143// Standard constructor
144// --------------------
145namespace Lorene {
146Etoile_rot::Etoile_rot(Map& mpi, int nzet_i, bool relat, const Eos& eos_i)
147 : Etoile(mpi, nzet_i, relat, eos_i),
148 bbb(mpi),
149 b_car(mpi),
150 nphi(mpi),
151 tnphi(mpi),
152 uuu(mpi),
153 logn(logn_auto),
154 nuf(mpi),
155 nuq(mpi),
157 tggg(mpi),
158 w_shift(mpi, 1, CON, mp.get_bvect_cart()),
159 khi_shift(mpi),
160 tkij(mpi, 2, COV, mp.get_bvect_cart()),
161 ak_car(mpi),
162 ssjm1_nuf(mpi),
163 ssjm1_nuq(mpi),
164 ssjm1_dzeta(mpi),
165 ssjm1_tggg(mpi),
166 ssjm1_khi(mpi),
167 ssjm1_wshift(mpi, 1, CON, mp.get_bvect_cart())
168{
169
170 // Initialization to a static state :
171 omega = 0 ;
172 uuu = 0 ;
173
174 // Initialization to a flat metric :
175 bbb = 1 ;
176 bbb.set_std_base() ;
177 b_car = 1 ;
178 b_car.set_std_base() ;
179 nphi = 0 ;
180 tnphi = 0 ;
181 nuf = 0 ;
182 nuq = 0 ;
183 tggg = 0 ;
184
185 w_shift.set_etat_qcq() ;
186 for (int i=0; i<3; i++) {
187 w_shift.set(i) = 0 ;
188 }
189
190 khi_shift.set_etat_qcq() ;
191 khi_shift.set() = 0 ;
192
193 tkij.set_etat_zero() ;
194
195 ak_car = 0 ;
196
197 ssjm1_nuf = 0 ;
198 ssjm1_nuq = 0 ;
199 ssjm1_dzeta = 0 ;
200 ssjm1_tggg = 0 ;
201 ssjm1_khi = 0 ;
202
203 ssjm1_wshift.set_etat_qcq() ;
204 for (int i=0; i<3; i++) {
205 ssjm1_wshift.set(i) = 0 ;
206 }
207
208 // Pointers of derived quantities initialized to zero :
209 set_der_0x0() ;
210
211}
212
213// Copy constructor
214// ----------------
215
217 : Etoile(et),
218 bbb(et.bbb),
219 b_car(et.b_car),
220 nphi(et.nphi),
221 tnphi(et.tnphi),
222 uuu(et.uuu),
223 logn(logn_auto),
224 nuf(et.nuf),
225 nuq(et.nuq),
227 tggg(et.tggg),
228 w_shift(et.w_shift),
229 khi_shift(et.khi_shift),
230 tkij(et.tkij),
231 ak_car(et.ak_car),
232 ssjm1_nuf(et.ssjm1_nuf),
233 ssjm1_nuq(et.ssjm1_nuq),
236 ssjm1_khi(et.ssjm1_khi),
238{
239 omega = et.omega ;
240
241 // Pointers of derived quantities initialized to zero :
242 set_der_0x0() ;
243}
244
245
246// Constructor from a file
247// -----------------------
248Etoile_rot::Etoile_rot(Map& mpi, const Eos& eos_i, FILE* fich)
249 : Etoile(mpi, eos_i, fich),
250 bbb(mpi),
251 b_car(mpi),
252 nphi(mpi),
253 tnphi(mpi),
254 uuu(mpi),
255 logn(logn_auto),
256 nuf(mpi),
257 nuq(mpi),
259 tggg(mpi),
260 w_shift(mpi, 1, CON, mp.get_bvect_cart()),
261 khi_shift(mpi),
262 tkij(mpi, 2, COV, mp.get_bvect_cart()),
263 ak_car(mpi),
264 ssjm1_nuf(mpi),
265 ssjm1_nuq(mpi),
266 ssjm1_dzeta(mpi),
267 ssjm1_tggg(mpi),
268 ssjm1_khi(mpi),
269 ssjm1_wshift(mpi, 1, CON, mp.get_bvect_cart())
270{
271
272 // Etoile parameters
273 // -----------------
274
275 // omega is read in the file:
276 fread_be(&omega, sizeof(double), 1, fich) ;
277
278
279 // Read of the saved fields:
280 // ------------------------
281
282 Tenseur nuf_file(mp, fich) ;
283 nuf = nuf_file ;
284
285 Tenseur nuq_file(mp, fich) ;
286 nuq = nuq_file ;
287
288 Tenseur tggg_file(mp, fich) ;
289 tggg = tggg_file ;
290
291 Tenseur w_shift_file(mp, mp.get_bvect_cart(), fich) ;
292 w_shift = w_shift_file ;
293
294 Tenseur khi_shift_file(mp, fich) ;
295 khi_shift = khi_shift_file ;
296
297 fait_shift() ; // constructs shift from w_shift and khi_shift
298 fait_nphi() ; // constructs N^phi from (N^x,N^y,N^z)
299
300 Cmp ssjm1_nuf_file(mp, *(mp.get_mg()), fich) ;
301 ssjm1_nuf = ssjm1_nuf_file ;
302
303 Cmp ssjm1_nuq_file(mp, *(mp.get_mg()), fich) ;
304 ssjm1_nuq = ssjm1_nuq_file ;
305
306 Cmp ssjm1_dzeta_file(mp, *(mp.get_mg()), fich) ;
307 ssjm1_dzeta = ssjm1_dzeta_file ;
308
309 Cmp ssjm1_tggg_file(mp, *(mp.get_mg()), fich) ;
310 ssjm1_tggg = ssjm1_tggg_file ;
311
312 Cmp ssjm1_khi_file(mp, *(mp.get_mg()), fich) ;
313 ssjm1_khi = ssjm1_khi_file ;
314
315 Tenseur ssjm1_wshift_file(mp, mp.get_bvect_cart(), fich) ;
316 ssjm1_wshift = ssjm1_wshift_file ;
317
318 // All other fields are initialized to zero :
319 // ----------------------------------------
320 bbb = 0 ;
321 b_car = 0 ;
322 uuu = 0 ;
323
324 // Pointers of derived quantities initialized to zero
325 // --------------------------------------------------
326 set_der_0x0() ;
327
328}
329
330 //------------//
331 // Destructor //
332 //------------//
333
335
336 del_deriv() ;
337
338}
339
340 //----------------------------------//
341 // Management of derived quantities //
342 //----------------------------------//
343
345
347
348 if (p_angu_mom != 0x0) delete p_angu_mom ;
349 if (p_tsw != 0x0) delete p_tsw ;
350 if (p_grv2 != 0x0) delete p_grv2 ;
351 if (p_grv3 != 0x0) delete p_grv3 ;
352 if (p_r_circ != 0x0) delete p_r_circ ;
353 if (p_r_circ_merid != 0x0) delete p_r_circ_merid ;
354 if (p_area != 0x0) delete p_area ;
355 if (p_aplat != 0x0) delete p_aplat ;
356 if (p_z_eqf != 0x0) delete p_z_eqf ;
357 if (p_z_eqb != 0x0) delete p_z_eqb ;
358 if (p_z_pole != 0x0) delete p_z_pole ;
359 if (p_mom_quad != 0x0) delete p_mom_quad ;
360 if (p_mom_quad_old != 0x0) delete p_mom_quad_old ;
361 if (p_mom_quad_Bo != 0x0) delete p_mom_quad_Bo ;
362 if (p_r_isco != 0x0) delete p_r_isco ;
363 if (p_f_isco != 0x0) delete p_f_isco ;
364 if (p_lspec_isco != 0x0) delete p_lspec_isco ;
365 if (p_espec_isco != 0x0) delete p_espec_isco ;
366 if (p_f_eq != 0x0) delete p_f_eq ;
367
369}
370
371
372
373
375
377
378 p_angu_mom = 0x0 ;
379 p_tsw = 0x0 ;
380 p_grv2 = 0x0 ;
381 p_grv3 = 0x0 ;
382 p_r_circ = 0x0 ;
383 p_r_circ_merid = 0x0 ;
384 p_area = 0x0 ;
385 p_aplat = 0x0 ;
386 p_z_eqf = 0x0 ;
387 p_z_eqb = 0x0 ;
388 p_z_pole = 0x0 ;
389 p_mom_quad = 0x0 ;
390 p_mom_quad_old = 0x0 ;
391 p_mom_quad_Bo = 0x0 ;
392 p_r_isco = 0x0 ;
393 p_f_isco = 0x0 ;
394 p_lspec_isco = 0x0 ;
395 p_espec_isco = 0x0 ;
396 p_f_eq = 0x0 ;
397
398}
399
401
403
404 del_deriv() ;
405
406}
407
408
409 //--------------//
410 // Assignment //
411 //--------------//
412
413// Assignment to another Etoile_rot
414// --------------------------------
416
417 // Assignment of quantities common to all the derived classes of Etoile
418 Etoile::operator=(et) ;
419
420 // Assignement of proper quantities of class Etoile_rot
421 omega = et.omega ;
422
423 bbb = et.bbb ;
424 b_car = et.b_car ;
425 nphi = et.nphi ;
426 tnphi = et.tnphi ;
427 uuu = et.uuu ;
428 nuf = et.nuf ;
429 nuq = et.nuq ;
430 tggg = et.tggg ;
431 w_shift = et.w_shift ;
432 khi_shift = et.khi_shift ;
433 tkij = et.tkij ;
434 ak_car = et.ak_car ;
435 ssjm1_nuf = et.ssjm1_nuf ;
436 ssjm1_nuq = et.ssjm1_nuq ;
439 ssjm1_khi = et.ssjm1_khi ;
441
442 del_deriv() ; // Deletes all derived quantities
443
444}
445
446 //--------------//
447 // Outputs //
448 //--------------//
449
450// Save in a file
451// --------------
452void Etoile_rot::sauve(FILE* fich) const {
453
454 Etoile::sauve(fich) ;
455
456 fwrite_be(&omega, sizeof(double), 1, fich) ;
457
458 nuf.sauve(fich) ;
459 nuq.sauve(fich) ;
460 tggg.sauve(fich) ;
461 w_shift.sauve(fich) ;
462 khi_shift.sauve(fich) ;
463
464 ssjm1_nuf.sauve(fich) ;
465 ssjm1_nuq.sauve(fich) ;
466 ssjm1_dzeta.sauve(fich) ;
467 ssjm1_tggg.sauve(fich) ;
468 ssjm1_khi.sauve(fich) ;
469 ssjm1_wshift.sauve(fich) ;
470
471
472}
473
474// Printing
475// --------
476
477ostream& Etoile_rot::operator>>(ostream& ost) const {
478
479 using namespace Unites ;
480
481 Etoile::operator>>(ost) ;
482
483 double omega_c = get_omega_c() ;
484
485 ost << endl ;
486 if (omega != __infinity) {
487 ost << "Uniformly rotating star" << endl ;
488 ost << "-----------------------" << endl ;
489
490 double freq = omega / (2.*M_PI) ;
491 ost << "Omega : " << omega * f_unit
492 << " rad/s f : " << freq * f_unit << " Hz" << endl ;
493 ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
494 << endl ;
495
496 }
497 else {
498 ost << "Differentially rotating star" << endl ;
499 ost << "----------------------------" << endl ;
500
501 double freq = omega_c / (2.*M_PI) ;
502 ost << "Central value of Omega : " << omega_c * f_unit
503 << " rad/s f : " << freq * f_unit << " Hz" << endl ;
504 ost << "Central rotation period : " << 1000. / (freq * f_unit) << " ms"
505 << endl ;
506
507 }
508
509
510 double nphi_c = nphi()(0, 0, 0, 0) ;
511 if ( (omega_c==0) && (nphi_c==0) ) {
512 ost << "Central N^phi : " << nphi_c << endl ;
513 }
514 else{
515 ost << "Central N^phi/Omega : " << nphi_c / omega_c << endl ;
516 }
517
518 ost << "Error on the virial identity GRV2 : " << endl ;
519 ost << "GRV2 = " << grv2() << endl ;
520 ost << "Error on the virial identity GRV3 : " << endl ;
521 double xgrv3 = grv3(&ost) ;
522 ost << "GRV3 = " << xgrv3 << endl ;
523
524 double mom_quad_38si = mom_quad() * rho_unit * (pow(r_unit, double(5.))
525 / double(1.e38) ) ;
526 ost << "Quadrupole moment Q : " << mom_quad_38si << " 10^38 kg m^2"
527 << endl ;
528 ost << "Q / (M R_circ^2) : "
529 << mom_quad() / ( mass_g() * pow( r_circ(), 2. ) ) << endl ;
530 ost << "c^4 Q / (G^2 M^3) : "
531 << mom_quad() / ( pow(qpig/(4*M_PI), 2.) * pow(mass_g(), 3.) )
532 << endl ;
533
534 ost << "Angular momentum J : "
535 << angu_mom()/( qpig / (4* M_PI) * msol*msol) << " G M_sol^2 / c"
536 << endl ;
537 ost << "c J / (G M^2) : "
538 << angu_mom()/( qpig / (4* M_PI) * pow(mass_g(), 2.) ) << endl ;
539
540 if (omega != __infinity) {
541 double mom_iner = angu_mom() / omega ;
542 double mom_iner_38si = mom_iner * rho_unit * (pow(r_unit, double(5.))
543 / double(1.e38) ) ;
544 ost << "Moment of inertia: " << mom_iner_38si << " 10^38 kg m^2"
545 << endl ;
546 }
547
548 ost << "Ratio T/W : " << tsw() << endl ;
549 ost << "Circumferential equatorial radius R_circ_eq : "
550 << r_circ()/km << " km" << endl ;
551 if (mp.get_mg()->get_np(0) == 1) {
552 ost << "Circumferential meridional radius R_circ_merid : "
553 << r_circ_merid()/km << " km" << endl ;
554 ost << "Flattening r_circ_merid/r_circ_eq : "
555 << r_circ_merid() / r_circ() << endl ;
556 ost << "Surface area : " << area()/(km*km) << " km^2" << endl ;
557 ost << "Mean radius R_mean : "
558 << mean_radius()/km << " km" << endl ;
559 } else {
560 ost <<
561 "Skipping polar radius / surface statements due to number of points in phi direction np > 1"
562 << endl;
563 }
564 ost << "Coordinate equatorial radius r_eq : " << ray_eq()/km << " km"
565 << endl ;
566 ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
567
568 double compact = qpig/(4.*M_PI) * mass_g() / r_circ() ;
569 ost << "Compaction parameter M_g / R_circ : " << compact << endl ;
570
571 int lsurf = nzet - 1;
572 int nt = mp.get_mg()->get_nt(lsurf) ;
573 int nr = mp.get_mg()->get_nr(lsurf) ;
574 ost << "Equatorial value of the velocity U: "
575 << uuu()(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
576 ost << "Redshift at the equator (forward) : " << z_eqf() << endl ;
577 ost << "Redshift at the equator (backward): " << z_eqb() << endl ;
578 ost << "Redshift at the pole : " << z_pole() << endl ;
579
580
581 ost << "Central value of log(N) : "
582 << logn()(0, 0, 0, 0) << endl ;
583
584 ost << "Central value of dzeta=log(AN) : "
585 << dzeta()(0, 0, 0, 0) << endl ;
586
587 if ( (omega_c==0) && (nphi_c==0) ) {
588 ost << "Central N^phi : " << nphi_c << endl ;
589 }
590 else{
591 ost << "Central N^phi/Omega : " << nphi_c / omega_c << endl ;
592 }
593
594 ost << " ... w_shift (NB: components in the star Cartesian frame) [c] : "
595 << endl
596 << w_shift(0)(0, 0, 0, 0) << " "
597 << w_shift(1)(0, 0, 0, 0) << " "
598 << w_shift(2)(0, 0, 0, 0) << endl ;
599
600 ost << "Central value of khi_shift [km c] : "
601 << khi_shift()(0, 0, 0, 0) / km << endl ;
602
603 ost << "Central value of B^2 : " << b_car()(0,0,0,0) << endl ;
604
605 Tbl diff_a_b = diffrel( a_car(), b_car() ) ;
606 ost <<
607 "Relative discrepancy in each domain between the metric coef. A^2 and B^2 : "
608 << endl ;
609 for (int l=0; l<diff_a_b.get_taille(); l++) {
610 ost << diff_a_b(l) << " " ;
611 }
612 ost << endl ;
613
614 // Approximate formula for R_isco = 6 R_g (1-(2/3)^1.5 j )
615 // up to the first order in j
616 double jdimless = angu_mom() / ( ggrav * pow(mass_g(), 2.) ) ;
617 double r_grav = ggrav * mass_g() ;
618 double r_isco_appr = 6. * r_grav * ( 1. - pow(2./3.,1.5) * jdimless ) ;
619
620 // Approximate formula for the ISCO frequency
621 // freq_ms = 6^{-1.5}/2pi/R_g (1+11*6^(-1.5) j )
622 // up to the first order in j
623 double f_isco_appr = ( 1. + 11. /6. /sqrt(6.) * jdimless ) / r_grav /
624 (12. * M_PI ) / sqrt(6.) ;
625
626 ost << endl << "Innermost stable circular orbit (ISCO) : " << endl ;
627 double xr_isco = r_isco(&ost) ;
628 ost <<" circumferential radius r_isco = "
629 << xr_isco / km << " km" << endl ;
630 ost <<" (approx. 6M + 1st order in j : "
631 << r_isco_appr / km << " km)" << endl ;
632 ost <<" (approx. 6M : "
633 << 6. * r_grav / km << " km)" << endl ;
634 ost <<" orbital frequency f_isco = "
635 << f_isco() * f_unit << " Hz" << endl ;
636 ost <<" (approx. 1st order in j : "
637 << f_isco_appr * f_unit << " Hz)" << endl ;
638
639
640 return ost ;
641
642}
643
644
645void Etoile_rot::partial_display(ostream& ost) const {
646
647 using namespace Unites ;
648
649 double omega_c = get_omega_c() ;
650 double freq = omega_c / (2.*M_PI) ;
651 ost << "Central Omega : " << omega_c * f_unit
652 << " rad/s f : " << freq * f_unit << " Hz" << endl ;
653 ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
654 << endl ;
655 ost << endl << "Central enthalpy : " << ent()(0,0,0,0) << " c^2" << endl ;
656 ost << "Central proper baryon density : " << nbar()(0,0,0,0)
657 << " x 0.1 fm^-3" << endl ;
658 ost << "Central proper energy density : " << ener()(0,0,0,0)
659 << " rho_nuc c^2" << endl ;
660 ost << "Central pressure : " << press()(0,0,0,0)
661 << " rho_nuc c^2" << endl ;
662
663 ost << "Central value of log(N) : "
664 << logn()(0, 0, 0, 0) << endl ;
665 ost << "Central lapse N : " << nnn()(0,0,0,0) << endl ;
666 ost << "Central value of dzeta=log(AN) : "
667 << dzeta()(0, 0, 0, 0) << endl ;
668 ost << "Central value of A^2 : " << a_car()(0,0,0,0) << endl ;
669 ost << "Central value of B^2 : " << b_car()(0,0,0,0) << endl ;
670
671 double nphi_c = nphi()(0, 0, 0, 0) ;
672 if ( (omega_c==0) && (nphi_c==0) ) {
673 ost << "Central N^phi : " << nphi_c << endl ;
674 }
675 else{
676 ost << "Central N^phi/Omega : " << nphi_c / omega_c
677 << endl ;
678 }
679
680
681 int lsurf = nzet - 1;
682 int nt = mp.get_mg()->get_nt(lsurf) ;
683 int nr = mp.get_mg()->get_nr(lsurf) ;
684 ost << "Equatorial value of the velocity U: "
685 << uuu()(nzet-1, 0, nt-1, nr-1) << " c" << endl ;
686
687 ost << endl
688 << "Coordinate equatorial radius r_eq = "
689 << ray_eq()/km << " km" << endl ;
690 ost << "Flattening r_pole/r_eq : " << aplat() << endl ;
691 if (mp.get_mg()->get_np(0) == 1)
692 ost << "Flattening r_circ_pole/r_circ_eq : "
693 << r_circ_merid() / r_circ() << endl ;
694
695}
696
697
699
700 return omega ;
701
702}
703
704
705// display_poly
706// ------------
707
708void Etoile_rot::display_poly(ostream& ost) const {
709
710 using namespace Unites ;
711
712 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
713
714 if (p_eos_poly != 0x0) {
715
716 double kappa = p_eos_poly->get_kap() ;
717 double gamma = p_eos_poly->get_gam() ; ;
718
719 // kappa^{n/2}
720 double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
721
722 // Polytropic unit of length in terms of r_unit :
723 double r_poly = kap_ns2 / sqrt(ggrav) ;
724
725 // Polytropic unit of time in terms of t_unit :
726 double t_poly = r_poly ;
727
728 // Polytropic unit of mass in terms of m_unit :
729 double m_poly = r_poly / ggrav ;
730
731 // Polytropic unit of angular momentum in terms of j_unit :
732 double j_poly = r_poly * r_poly / ggrav ;
733
734 // Polytropic unit of density in terms of rho_unit :
735 double rho_poly = 1. / (ggrav * r_poly * r_poly) ;
736
737 ost.precision(10) ;
738 ost << endl << "Quantities in polytropic units : " << endl ;
739 ost << "==============================" << endl ;
740 ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
741 ost << " n_c : " << nbar()(0, 0, 0, 0) / rho_poly << endl ;
742 ost << " e_c : " << ener()(0, 0, 0, 0) / rho_poly << endl ;
743 ost << " Omega_c : " << get_omega_c() * t_poly << endl ;
744 ost << " P_c : " << 2.*M_PI / get_omega_c() / t_poly << endl ;
745 ost << " M_bar : " << mass_b() / m_poly << endl ;
746 ost << " M : " << mass_g() / m_poly << endl ;
747 ost << " J : " << angu_mom() / j_poly << endl ;
748 ost << " r_eq : " << ray_eq() / r_poly << endl ;
749 ost << " R_circ : " << r_circ() / r_poly << endl ;
750
751
752 }
753
754
755}
756
757
758
759
760
761
762 //-------------------------//
763 // Computational routines //
764 //-------------------------//
765
767
768 Tenseur d_khi = khi_shift.gradient() ;
769
770 if (d_khi.get_etat() == ETATQCQ) {
771 d_khi.dec2_dzpuis() ; // divide by r^2 in the external compactified
772 // domain
773 }
774
775 // x_k dW^k/dx_i
776
777 Tenseur x_d_w = skxk( w_shift.gradient() ) ;
778 x_d_w.dec_dzpuis() ;
779
780 double lambda = double(1) / double(3) ;
781
782 if ( mp.get_mg()->get_np(0) == 1 ) {
783 lambda = 0 ;
784 }
785
786 // The final computation is done component by component because
787 // d_khi and x_d_w are covariant comp. whereas w_shift is
788 // contravariant
789
790 shift.set_etat_qcq() ;
791
792 for (int i=0; i<3; i++) {
793 shift.set(i) = (lambda+2)/2./(lambda+1) * w_shift(i)
794 - (lambda/2./(lambda+1)) * (d_khi(i) + x_d_w(i)) ;
795 }
796
797 shift.set_triad( *(w_shift.get_triad()) ) ;
798
799}
800
801
802
804
805 if ( shift.get_etat() == ETATZERO ) {
806 tnphi = 0 ;
807 nphi = 0 ;
808 return ;
809 }
810
811 assert( shift.get_etat() == ETATQCQ ) ;
812
813 // Computation of tnphi
814 // --------------------
815 tnphi.set_etat_qcq() ;
816
817 mp.comp_p_from_cartesian( shift(0), shift(1), tnphi.set() ) ;
818
819 // Computation of nphi
820 // -------------------
821
822 nphi = tnphi ;
823 (nphi.set()).div_rsint() ;
824
825}
826}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Polytropic equation of state (relativistic case).
Definition eos.h:809
double get_gam() const
Returns the adiabatic index (cf. Eq. (3)).
Definition eos_poly.C:271
double get_kap() const
Returns the pressure coefficient (cf.
Definition eos_poly.C:275
Equation of state base class.
Definition eos.h:206
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Definition etoile.h:1628
double * p_r_circ_merid
Circumferential meridional radius.
Definition etoile.h:1639
Tenseur uuu
Norm of u_euler.
Definition etoile.h:1521
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile_rot.C:344
double omega
Rotation angular velocity ([f_unit] ).
Definition etoile.h:1504
double * p_z_pole
Redshift factor at North pole.
Definition etoile.h:1644
virtual double mom_quad() const
Quadrupole moment.
double * p_z_eqf
Forward redshift factor at equator.
Definition etoile.h:1642
void operator=(const Etoile_rot &)
Assignment to another Etoile_rot.
Definition etoile_rot.C:415
virtual double r_circ() const
Circumferential equatorial radius.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile_rot.C:374
Tenseur & logn
Metric potential = logn_auto.
Definition etoile.h:1524
double * p_mom_quad_old
Part of the quadrupole moment.
Definition etoile.h:1646
Cmp ssjm1_dzeta
Effective source at the previous step for the resolution of the Poisson equation for dzeta .
Definition etoile.h:1606
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition etoile.h:1534
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition etoile_rot.C:477
double * p_z_eqb
Backward redshift factor at equator.
Definition etoile.h:1643
virtual double mass_g() const
Gravitational mass.
virtual double f_isco() const
Orbital frequency at the innermost stable circular orbit (ISCO).
double * p_r_circ
Circumferential equatorial radius.
Definition etoile.h:1638
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1563
virtual double aplat() const
Flatening r_pole/r_eq.
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Etoile_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition etoile_rot.C:146
virtual double tsw() const
Ratio T/W.
Tenseur tggg
Metric potential .
Definition etoile.h:1540
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Definition etoile.h:1529
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Definition etoile.h:1611
Tenseur nphi
Metric coefficient .
Definition etoile.h:1513
virtual double mass_b() const
Baryon mass.
virtual double r_isco(ostream *ost=0x0) const
Circumferential radius of the innermost stable circular orbit (ISCO).
Definition et_rot_isco.C:87
Tenseur bbb
Metric factor B.
Definition etoile.h:1507
double * p_area
Surface area.
Definition etoile.h:1640
Tenseur ak_car
Scalar .
Definition etoile.h:1589
virtual ~Etoile_rot()
Destructor.
Definition etoile_rot.C:334
Tenseur & dzeta
Metric potential = beta_auto.
Definition etoile.h:1537
virtual double mean_radius() const
Mean radius.
double * p_grv3
Error on the virial identity GRV3.
Definition etoile.h:1637
double * p_grv2
Error on the virial identity GRV2.
Definition etoile.h:1636
double * p_aplat
Flatening r_pole/r_eq.
Definition etoile.h:1641
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition etoile.h:1595
double * p_mom_quad_Bo
Part of the quadrupole moment.
Definition etoile.h:1647
virtual double grv2() const
Error on the virial identity GRV2.
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Definition etoile_rot.C:803
virtual double r_circ_merid() const
Circumferential meridional radius.
double * p_mom_quad
Quadrupole moment.
Definition etoile.h:1645
double * p_f_eq
Orbital frequency at the equator.
Definition etoile.h:1654
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Definition etoile_rot.C:400
virtual double area() const
Surface area.
double * p_angu_mom
Angular momentum.
Definition etoile.h:1634
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:1619
virtual double angu_mom() const
Angular momentum.
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Definition etoile.h:1570
virtual double z_eqf() const
Forward redshift factor at equator.
virtual void display_poly(ostream &) const
Display in polytropic units.
Definition etoile_rot.C:708
double * p_f_isco
Orbital frequency of the ISCO.
Definition etoile.h:1649
double * p_lspec_isco
Specific angular momentum of a particle on the ISCO.
Definition etoile.h:1653
double * p_tsw
Ratio T/W.
Definition etoile.h:1635
Tenseur b_car
Square of the metric factor B.
Definition etoile.h:1510
virtual double z_eqb() const
Backward redshift factor at equator.
Tenseur tnphi
Component of the shift vector.
Definition etoile.h:1518
double * p_espec_isco
Specific energy of a particle on the ISCO.
Definition etoile.h:1651
double * p_r_isco
Circumferential radius of the ISCO.
Definition etoile.h:1648
virtual void sauve(FILE *) const
Save in a file.
Definition etoile_rot.C:452
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] ).
Definition etoile_rot.C:698
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition etoile.h:1601
virtual double z_pole() const
Redshift factor at North pole.
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1553
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Definition etoile_rot.C:645
void fait_shift()
Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog.
Definition etoile_rot.C:766
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile.C:399
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:435
void operator=(const Etoile &)
Assignment to another Etoile.
Definition etoile.C:433
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Definition etoile.h:512
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:487
Tenseur nbar
Baryon density in the fluid frame.
Definition etoile.h:462
const Eos & eos
Equation of state of the stellar matter.
Definition etoile.h:454
Map & mp
Mapping associated with the star.
Definition etoile.h:432
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile.C:381
Tenseur ener
Total energy density in the fluid frame.
Definition etoile.h:463
Etoile(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition etoile.C:146
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition etoile.C:514
Tenseur press
Fluid pressure.
Definition etoile.h:464
virtual void sauve(FILE *) const
Save in a file.
Definition etoile.C:486
Tenseur shift
Total shift vector.
Definition etoile.h:515
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Definition etoile.h:460
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Definition etoile.C:413
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition etoile.h:509
Tenseur a_car
Total conformal factor .
Definition etoile.h:518
Basic array class.
Definition tbl.h:161
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
void dec2_dzpuis()
dzpuis -= 2 ;
Definition tenseur.C:1136
void dec_dzpuis()
dzpuis -= 1 ;
Definition tenseur.C:1110
int get_etat() const
Returns the logical state.
Definition tenseur.h:710
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:507
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
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 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
Tenseur skxk(const Tenseur &)
Contraction of the last index of (*this) with or , depending on the type of S .
Lorene prototypes.
Definition app_hor.h:67
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
Definition map.h:803
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.