LORENE
etoile.C
1/*
2 * Methods of class Etoile
3 *
4 * (see file etoile.h for documentation)
5 */
6
7/*
8 * Copyright (c) 2000-2001 Eric Gourgoulhon
9 * Copyright (c) 2000-2001 Keisuke Taniguchi
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: etoile.C,v 1.11 2016/12/05 16:17:54 j_novak Exp $
34 * $Log: etoile.C,v $
35 * Revision 1.11 2016/12/05 16:17:54 j_novak
36 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
37 *
38 * Revision 1.10 2014/10/13 08:52:58 j_novak
39 * Lorene classes and functions now belong to the namespace Lorene.
40 *
41 * Revision 1.9 2012/08/12 17:48:35 p_cerda
42 * Magnetstar: New classes for magnetstar. Allowing for non-equatorial symmetry in Etoile et al. Adding B_phi in Et_rot_mag.
43 *
44 * Revision 1.8 2005/01/18 22:36:50 k_taniguchi
45 * Delete a pointer for ray_eq(int kk).
46 *
47 * Revision 1.7 2005/01/18 20:35:05 k_taniguchi
48 * Addition of ray_eq(int kk).
49 *
50 * Revision 1.6 2005/01/17 20:40:25 k_taniguchi
51 * Addition of ray_eq_3pis2().
52 *
53 * Revision 1.5 2004/03/25 10:29:06 j_novak
54 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
55 *
56 * Revision 1.4 2003/10/13 15:23:56 f_limousin
57 * *** empty log message ***
58 *
59 * Revision 1.3 2002/04/09 14:32:15 e_gourgoulhon
60 * 1/ Added extra parameters in EOS computational functions (argument par)
61 * 2/ New class MEos for multi-domain EOS
62 *
63 * Revision 1.2 2001/12/04 21:27:53 e_gourgoulhon
64 *
65 * All writing/reading to a binary file are now performed according to
66 * the big endian convention, whatever the system is big endian or
67 * small endian, thanks to the functions fwrite_be and fread_be
68 *
69 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
70 * LORENE
71 *
72 * Revision 2.14 2000/11/24 13:27:44 eric
73 * Dans eqution_of_state(): changement leger de ent dans le cas ou l'on a
74 * deux domaine avant d'appeler l'EOS.
75 *
76 * Revision 2.13 2000/09/25 12:22:02 keisuke
77 * *** empty log message ***
78 *
79 * Revision 2.12 2000/09/22 15:50:58 keisuke
80 * Ajout du membre d_logn_auto_div.
81 *
82 * Revision 2.11 2000/09/07 14:34:09 keisuke
83 * Ajout du membre logn_auto_regu.
84 *
85 * Revision 2.10 2000/08/31 15:36:54 eric
86 * Bases spectrales standards pour nnn, a_car et gam_euler dans le
87 * constructeur (initialisation a la metrique plate).
88 *
89 * Revision 2.9 2000/08/29 11:37:49 eric
90 * Ajout des membres k_div et logn_auto_div.
91 *
92 * Revision 2.8 2000/07/21 12:01:11 eric
93 * Modif dans Etoile::del_deriv() :
94 * appel de Etoile::set_der_0x0() et non de la fonction virtuelle set_der_0x0().
95 *
96 * Revision 2.7 2000/03/21 12:39:34 eric
97 * Le constructeur standard teste la compatibilite de l'EOS avec le
98 * caractere relativiste de l'etoile.
99 *
100 * Revision 2.6 2000/02/21 14:32:40 eric
101 * gam_euler est initialise a 1 dans le constructeur standard.
102 * Suppression de l'appel a del_hydro_euler dans equation_of_state().
103 *
104 * Revision 2.5 2000/02/09 19:30:47 eric
105 * La triade de decomposition doit desormais figurer en argument des
106 * constructeurs de Tenseur.
107 *
108 * Revision 2.4 2000/02/02 09:23:34 eric
109 * Affichage de la masse.
110 *
111 * Revision 2.3 2000/01/28 17:18:10 eric
112 * Modifs noms des quantites globales.
113 * Affichage.
114 *
115 * Revision 2.2 2000/01/24 17:13:36 eric
116 * Le mapping mp n'est plus constant.
117 *
118 * Revision 2.1 2000/01/24 13:37:22 eric
119 * *** empty log message ***
120 *
121 * Revision 2.0 2000/01/20 17:04:45 eric
122 * *** empty log message ***
123 *
124 *
125 * $Header: /cvsroot/Lorene/C++/Source/Etoile/etoile.C,v 1.11 2016/12/05 16:17:54 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 "utilitaires.h"
136#include "param.h"
137#include "unites.h"
138
139 //--------------//
140 // Constructors //
141 //--------------//
142
143// Standard constructor
144// --------------------
145namespace Lorene {
146Etoile::Etoile(Map& mpi, int nzet_i, bool relat, const Eos& eos_i)
147 : mp(mpi),
148 nzet(nzet_i),
149 relativistic(relat),
150 k_div(0),
151 eos(eos_i),
152 ent(mpi),
153 nbar(mpi),
154 ener(mpi),
155 press(mpi),
156 ener_euler(mpi),
157 s_euler(mpi),
158 gam_euler(mpi),
159 u_euler(mpi, 1, CON, mp.get_bvect_cart()),
160 logn_auto(mpi),
161 logn_auto_regu(mpi),
162 logn_auto_div(mpi),
163 d_logn_auto_div(mpi, 1, COV, mp.get_bvect_spher()),
164 beta_auto(mpi),
165 nnn(mpi),
166 shift(mpi, 1, CON, mp.get_bvect_cart()),
167 a_car(mpi) {
168
169 // Check of the EOS
170 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;
171 const Eos_poly_newt* p_eos_poly_newt =
172 dynamic_cast<const Eos_poly_newt*>( &eos ) ;
173 const Eos_incomp* p_eos_incomp = dynamic_cast<const Eos_incomp*>( &eos ) ;
174 const Eos_incomp_newt* p_eos_incomp_newt =
175 dynamic_cast<const Eos_incomp_newt*>( &eos ) ;
176
177 if (relativistic) {
178
179 if (p_eos_poly_newt != 0x0) {
180 cout <<
181 "Etoile::Etoile : the EOS Eos_poly_newt must not be employed"
182 << " for a relativistic star ! " << endl ;
183 cout << "(Use Eos_poly instead)" << endl ;
184 abort() ;
185 }
186 if (p_eos_incomp_newt != 0x0) {
187 cout <<
188 "Etoile::Etoile : the EOS Eos_incomp_newt must not be employed"
189 << " for a relativistic star ! " << endl ;
190 cout << "(Use Eos_incomp instead)" << endl ;
191 abort() ;
192 }
193
194 }
195 else{
196
197 if ( (p_eos_poly != 0x0) && (p_eos_poly_newt == 0x0) ) {
198 cout <<
199 "Etoile::Etoile : the EOS Eos_poly must not be employed"
200 << " for a Newtonian star ! " << endl ;
201 cout << "(Use Eos_poly_newt instead)" << endl ;
202 abort() ;
203 }
204 if ( (p_eos_incomp != 0x0) && (p_eos_incomp_newt == 0x0) ) {
205 cout <<
206 "Etoile::Etoile : the EOS Eos_incomp must not be employed"
207 << " for a relativistic star ! " << endl ;
208 cout << "(Use Eos_incomp_newt instead)" << endl ;
209 abort() ;
210 }
211
212 }
213
214
215 // Parameter 1/c^2
216 unsurc2 = relativistic ? double(1) : double(0) ;
217
218 // Pointers of derived quantities initialized to zero :
219 set_der_0x0() ;
220
221 // All the matter quantities are initialized to zero :
222 nbar = 0 ;
223 ener = 0 ;
224 press = 0 ;
225 ent = 0 ;
226 ener_euler = 0 ;
227 s_euler = 0 ;
228 gam_euler = 1 ;
229 gam_euler.set_std_base() ;
230 u_euler = 0 ;
231
232 // The metric is initialized to the flat one :
233 logn_auto = 0 ;
234 logn_auto_regu = 0 ;
235 logn_auto_div = 0 ;
236 d_logn_auto_div = 0 ;
237 beta_auto = 0 ;
238 nnn = 1 ;
239 nnn.set_std_base() ;
240 shift = 0 ;
241 a_car = 1 ;
242 a_car.set_std_base() ;
243
244}
245
246// Copy constructor
247// ----------------
249 : mp(et.mp),
250 nzet(et.nzet),
252 unsurc2(et.unsurc2),
253 k_div(et.k_div),
254 eos(et.eos),
255 ent(et.ent),
256 nbar(et.nbar),
257 ener(et.ener),
258 press(et.press),
260 s_euler(et.s_euler),
261 gam_euler(et.gam_euler),
262 u_euler(et.u_euler),
263 logn_auto(et.logn_auto),
267 beta_auto(et.beta_auto),
268 nnn(et.nnn),
269 shift(et.shift),
270 a_car(et.a_car) {
271
272 set_der_0x0() ;
273
274}
275
276// Constructor from a file
277// -----------------------
278Etoile::Etoile(Map& mpi, const Eos& eos_i, FILE* fich)
279 : mp(mpi),
280 eos(eos_i),
281 ent(mpi),
282 nbar(mpi),
283 ener(mpi),
284 press(mpi),
285 ener_euler(mpi),
286 s_euler(mpi),
287 gam_euler(mpi),
288 u_euler(mpi, 1, CON, mp.get_bvect_cart()),
289 logn_auto(mpi),
290 logn_auto_regu(mpi),
291 logn_auto_div(mpi),
292 d_logn_auto_div(mpi, 1, COV, mp.get_bvect_spher()),
293 beta_auto(mpi),
294 nnn(mpi),
295 shift(mpi, 1, CON, mp.get_bvect_cart()),
296 a_car(mpi) {
297
298 // Etoile parameters
299 // -----------------
300
301 // nzet and relativistic are read in the file:
302 int xx ;
303 fread_be(&xx, sizeof(int), 1, fich) ;
304 k_div = xx / 1000 ; // integer part
305 nzet = xx - k_div * 1000 ;
306
307 fread(&relativistic, sizeof(bool), 1, fich) ;
308
309 // Parameter 1/c^2 is deduced from relativistic:
310 unsurc2 = relativistic ? double(1) : double(0) ;
311
312
313 // Equation of state
314 // -----------------
315
316 // Read of the saved EOS
317 Eos* p_eos_file = Eos::eos_from_file(fich) ;
318
319 // Comparison with the assigned EOS:
320 if (eos != *p_eos_file) {
321 cout <<
322 "Etoile::Etoile(const Map&, const Eos&, FILE*) : the EOS given in "
323 << endl <<
324 " argument and that read in the file are different !" << endl ;
325 abort() ;
326 }
327
328 // p_eos_file is no longer required (it was used only for checking the
329 // EOS compatibility)
330 delete p_eos_file ;
331
332 // Read of the saved fields:
333 // ------------------------
334 Tenseur ent_file(mp, fich) ;
335 ent = ent_file ;
336
337 Tenseur logn_auto_file(mp, fich) ;
338 logn_auto = logn_auto_file ;
339
340 Tenseur beta_auto_file(mp, fich) ;
341 beta_auto = beta_auto_file ;
342
343 if (k_div == 0) {
344 logn_auto_div = 0 ;
345 d_logn_auto_div = 0 ;
346 }
347 else {
348
349 Tenseur logn_auto_div_file(mp, fich) ;
350 logn_auto_div = logn_auto_div_file ;
351
352 Tenseur d_logn_auto_div_file(mp, mp.get_bvect_spher(), fich) ;
353 d_logn_auto_div = d_logn_auto_div_file ;
354 }
355
357
358 shift = 0 ;
359
360 // Pointers of derived quantities initialized to zero
361 // --------------------------------------------------
362 set_der_0x0() ;
363
364}
365
366 //------------//
367 // Destructor //
368 //------------//
369
371
372 del_deriv() ;
373
374}
375
376
377 //----------------------------------//
378 // Management of derived quantities //
379 //----------------------------------//
380
381void Etoile::del_deriv() const {
382
383 if (p_mass_b != 0x0) delete p_mass_b ;
384 if (p_mass_g != 0x0) delete p_mass_g ;
385 if (p_ray_eq != 0x0) delete p_ray_eq ;
386 if (p_ray_eq_pis2 != 0x0) delete p_ray_eq_pis2 ;
387 if (p_ray_eq_pi != 0x0) delete p_ray_eq_pi ;
388 if (p_ray_eq_3pis2 != 0x0) delete p_ray_eq_3pis2 ;
389 if (p_ray_pole != 0x0) delete p_ray_pole ;
390 if (p_l_surf != 0x0) delete p_l_surf ;
391 if (p_xi_surf != 0x0) delete p_xi_surf ;
392
394}
395
396
397
398
400
401 p_mass_b = 0x0 ;
402 p_mass_g = 0x0 ;
403 p_ray_eq = 0x0 ;
404 p_ray_eq_pis2 = 0x0 ;
405 p_ray_eq_pi = 0x0 ;
406 p_ray_eq_3pis2 = 0x0 ;
407 p_ray_pole = 0x0 ;
408 p_l_surf = 0x0 ;
409 p_xi_surf = 0x0 ;
410
411}
412
414
415 ener_euler.set_etat_nondef() ;
416 s_euler.set_etat_nondef() ;
417 gam_euler.set_etat_nondef() ;
418 u_euler.set_etat_nondef() ;
419
420 del_deriv() ;
421
422}
423
424
425
426
427 //--------------//
428 // Assignment //
429 //--------------//
430
431// Assignment to another Etoile
432// ----------------------------
433void Etoile::operator=(const Etoile& et) {
434
435 assert( &(et.mp) == &mp ) ; // Same mapping
436 assert( &(et.eos) == &eos ) ; // Same EOS
437
438 nzet = et.nzet ;
440 k_div = et.k_div ;
441 unsurc2 = et.unsurc2 ;
442
443 ent = et.ent ;
444 nbar = et.nbar ;
445 ener = et.ener ;
446 press = et.press ;
448 s_euler = et.s_euler ;
449 gam_euler = et.gam_euler ;
450 u_euler = et.u_euler ;
451 logn_auto = et.logn_auto ;
455 beta_auto = et.beta_auto ;
456 nnn = et.nnn ;
457 shift = et.shift ;
458 a_car = et.a_car ;
459
460
461 del_deriv() ; // Deletes all derived quantities
462
463}
464
465// Assignment of the enthalpy field
466// --------------------------------
467
468void Etoile::set_enthalpy(const Cmp& ent_i) {
469
470 ent = ent_i ;
471
472 // Update of (nbar, ener, press) :
474
475 // The derived quantities are obsolete:
476 del_deriv() ;
477
478}
479
480 //--------------//
481 // Outputs //
482 //--------------//
483
484// Save in a file
485// --------------
486void Etoile::sauve(FILE* fich) const {
487
488 int xx = nzet + k_div * 1000 ;
489 fwrite_be(&xx, sizeof(int), 1, fich) ;
490
491 fwrite(&relativistic, sizeof(bool), 1, fich) ;
492
493 eos.sauve(fich) ;
494
495 ent.sauve(fich) ;
496 logn_auto.sauve(fich) ;
497 beta_auto.sauve(fich) ;
498
499 if (k_div != 0) {
500 logn_auto_div.sauve(fich) ;
501 d_logn_auto_div.sauve(fich) ;
502 }
503
504}
505
506// Printing
507// --------
508
509ostream& operator<<(ostream& ost, const Etoile& et) {
510 et >> ost ;
511 return ost ;
512}
513
514ostream& Etoile::operator>>(ostream& ost) const {
515
516 using namespace Unites ;
517
518 ost << endl ;
519 if (relativistic) {
520 ost << "Relativistic star" << endl ;
521 ost << "-----------------" << endl ;
522 }
523 else {
524 ost << "Newtonian star" << endl ;
525 ost << "--------------" << endl ;
526 }
527
528 ost << "Number of domains occupied by the star : " << nzet << endl ;
529
530 ost << "Equation of state : " << endl ;
531 ost << eos << endl ;
532
533 ost << endl << "Central enthalpy : " << ent()(0,0,0,0) << " c^2" << endl ;
534 ost << "Central proper baryon density : " << nbar()(0,0,0,0)
535 << " x 0.1 fm^-3" << endl ;
536 ost << "Central proper energy density : " << ener()(0,0,0,0)
537 << " rho_nuc c^2" << endl ;
538 ost << "Central pressure : " << press()(0,0,0,0)
539 << " rho_nuc c^2" << endl ;
540
541 ost << endl
542 << "Regularization index of the gravitational potential : k_div = "
543 << k_div << endl ;
544 ost << "Central lapse N : " << nnn()(0,0,0,0) << endl ;
545 ost << "Central value of A^2 : " << a_car()(0,0,0,0) << endl ;
546
547 ost << endl
548 << "Coordinate equatorial radius (phi=0) a1 = "
549 << ray_eq()/km << " km" << endl ;
550 ost << "Coordinate equatorial radius (phi=pi/2) a2 = "
551 << ray_eq_pis2()/km << " km" << endl ;
552 ost << "Coordinate equatorial radius (phi=pi): "
553 << ray_eq_pi()/km << " km" << endl ;
554 ost << "Coordinate polar radius a3 = "
555 << ray_pole()/km << " km" << endl ;
556 ost << "Axis ratio a2/a1 = " << ray_eq_pis2() / ray_eq()
557 << " a3/a1 = " << ray_pole() / ray_eq() << endl ;
558
559 ost << endl << "Baryon mass : " << mass_b() / msol << " M_sol" << endl ;
560 ost << "Gravitational mass : " << mass_g() / msol << " M_sol" << endl ;
561
562 return ost ;
563}
564
565 //-----------------------------------------//
566 // Computation of hydro quantities //
567 //-----------------------------------------//
568
570
571 Cmp ent_eos = ent() ;
572
573
574 // Slight rescale of the enthalpy field in case of 2 domains inside the
575 // star
576
577
578 double epsilon = 1.e-12 ;
579
580 const Mg3d* mg = mp.get_mg() ;
581 int nz = mg->get_nzone() ;
582
583 Mtbl xi(mg) ;
584 xi.set_etat_qcq() ;
585 for (int l=0; l<nz; l++) {
586 xi.t[l]->set_etat_qcq() ;
587 for (int k=0; k<mg->get_np(l); k++) {
588 for (int j=0; j<mg->get_nt(l); j++) {
589 for (int i=0; i<mg->get_nr(l); i++) {
590 xi.set(l,k,j,i) =
591 mg->get_grille3d(l)->x[i] ;
592 }
593 }
594 }
595
596 }
597
598 Cmp fact_ent(mp) ;
599 fact_ent.allocate_all() ;
600
601 fact_ent.set(0) = 1 + epsilon * xi(0) * xi(0) ;
602 fact_ent.set(1) = 1 - 0.25 * epsilon * (xi(1) - 1) * (xi(1) - 1) ;
603
604 for (int l=nzet; l<nz; l++) {
605 fact_ent.set(l) = 1 ;
606 }
607
608 if (nzet > 1) {
609
610 if(nzet == 3) {
611 fact_ent.set(1) = 1 - 0.5 * epsilon * (xi(1) - 0.5) * (xi(1) - 0.5) ;
612 fact_ent.set(2) = 1 - 0.25 * epsilon * (xi(2) - 1) * (xi(2) - 1) ;
613 }
614
615 if (nzet > 3) {
616
617 cout << "Etoile::equation_of_state: not ready yet for nzet > 3 !"
618 << endl ;
619 }
620
621 ent_eos = fact_ent * ent_eos ;
622 ent_eos.std_base_scal() ;
623 }
624
625
626
627
628
629 // Call to the EOS (the EOS is called domain by domain in order to
630 // allow for the use of MEos)
631
632 Cmp tempo(mp) ;
633
634 nbar.set_etat_qcq() ;
635 nbar.set() = 0 ;
636 for (int l=0; l<nzet; l++) {
637
638 Param par ; // Paramater for multi-domain equation of state
639 par.add_int(l) ;
640
641 tempo = eos.nbar_ent(ent_eos, 1, l, &par) ;
642
643 nbar = nbar() + tempo ;
644
645 }
646
647 ener.set_etat_qcq() ;
648 ener.set() = 0 ;
649 for (int l=0; l<nzet; l++) {
650
651 Param par ; // Paramater for multi-domain equation of state
652 par.add_int(l) ;
653
654 tempo = eos.ener_ent(ent_eos, 1, l, &par) ;
655
656 ener = ener() + tempo ;
657
658 }
659
660 press.set_etat_qcq() ;
661 press.set() = 0 ;
662 for (int l=0; l<nzet; l++) {
663
664 Param par ; // Paramater for multi-domain equation of state
665 par.add_int(l) ;
666
667 tempo = eos.press_ent(ent_eos, 1, l, &par) ;
668
669 press = press() + tempo ;
670
671 }
672
673
674 // Set the bases for spectral expansion
675 nbar.set_std_base() ;
676 ener.set_std_base() ;
677 press.set_std_base() ;
678
679 // The Eulerian quantities are obsolete
680 //## del_hydro_euler() ;
681
682 // The derived quantities are obsolete
683 del_deriv() ;
684
685}
686
688
689 cout <<
690 "Etoile::hydro_euler : hydro_euler must be called via a derived class"
691 << endl << " of Etoile !" << endl ;
692
693 abort() ;
694
695}
696}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition cmp.C:326
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition cmp.C:647
Tbl & set(int l)
Read/write of the value in a given domain.
Definition cmp.h:724
Equation of state of incompressible matter (Newtonian case).
Definition eos.h:1816
Equation of state of incompressible matter (relativistic case).
Definition eos.h:1632
Polytropic equation of state (Newtonian case).
Definition eos.h:1107
Polytropic equation of state (relativistic case).
Definition eos.h:809
Equation of state base class.
Definition eos.h:206
static Eos * eos_from_file(FILE *)
Construction of an EOS from a binary file.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile.C:399
Tenseur logn_auto_div
Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by t...
Definition etoile.h:500
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:435
double * p_mass_b
Baryon mass.
Definition etoile.h:550
double * p_mass_g
Gravitational mass.
Definition etoile.h:551
Tbl * p_xi_surf
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the su...
Definition etoile.h:548
void operator=(const Etoile &)
Assignment to another Etoile.
Definition etoile.C:433
Tenseur logn_auto_regu
Regular part of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:494
double ray_eq_pi() const
Coordinate radius at , [r_unit].
int k_div
Index of regularity of the gravitational potential logn_auto .
Definition etoile.h:452
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
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition etoile.C:569
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition etoile.h:477
Itbl * p_l_surf
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surfa...
Definition etoile.h:542
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:474
friend ostream & operator<<(ostream &, const Etoile &)
Display.
Definition etoile.C:509
Map & mp
Mapping associated with the star.
Definition etoile.h:432
double * p_ray_eq_3pis2
Coordinate radius at , .
Definition etoile.h:533
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile.C:381
double * p_ray_eq_pi
Coordinate radius at , .
Definition etoile.h:530
virtual double mass_b() const
Baryon mass.
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
double ray_eq_pis2() const
Coordinate radius at , [r_unit].
virtual ~Etoile()
Destructor.
Definition etoile.C:370
double * p_ray_pole
Coordinate radius at .
Definition etoile.h:536
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 ).
Definition etoile.h:504
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 double mass_g() const
Gravitational mass.
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:440
double * p_ray_eq_pis2
Coordinate radius at , .
Definition etoile.h:527
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition etoile.h:468
virtual void sauve(FILE *) const
Save in a file.
Definition etoile.C:486
Tenseur shift
Total shift vector.
Definition etoile.h:515
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition etoile.C:687
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition etoile.h:471
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
void set_enthalpy(const Cmp &)
Assignment of the enthalpy field.
Definition etoile.C:468
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition etoile.h:445
double * p_ray_eq
Coordinate radius at , .
Definition etoile.h:524
double ray_pole() const
Coordinate radius at [r_unit].
double * x
Array of values of at the nr collocation points.
Definition grilles.h:215
Multi-domain grid.
Definition grilles.h:279
const Grille3d * get_grille3d(int l) const
Returns a pointer on the 3D mono-grid for domain no. l.
Definition grilles.h:517
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:479
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:474
int get_nzone() const
Returns the number of domains.
Definition grilles.h:465
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:469
Multi-domain array.
Definition mtbl.h:118
Tbl & set(int l)
Read/write of the Tbl in a given domain.
Definition mtbl.h:221
Tbl ** t
Array (size nzone ) of pointers on the Tbl 's.
Definition mtbl.h:132
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition mtbl.C:302
Parameter storage.
Definition param.h:125
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:249
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:364
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
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
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
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition map.h:795
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.