LORENE
star_bhns.C
1/*
2 * Methods of class Star_bhns
3 *
4 * (see file star_bhns.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2005-2007 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 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
29
30/*
31 * $Id: star_bhns.C,v 1.5 2016/12/05 16:18:16 j_novak Exp $
32 * $Log: star_bhns.C,v $
33 * Revision 1.5 2016/12/05 16:18:16 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.4 2014/10/13 08:53:40 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.3 2014/10/06 15:13:16 j_novak
40 * Modified #include directives to use c++ syntax.
41 *
42 * Revision 1.2 2008/05/15 19:12:38 k_taniguchi
43 * Change of some parameters.
44 *
45 * Revision 1.1 2007/06/22 01:30:10 k_taniguchi
46 * *** empty log message ***
47 *
48 *
49 * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns.C,v 1.5 2016/12/05 16:18:16 j_novak Exp $
50 *
51 */
52
53// C++ headers
54//#include <>
55
56// C headers
57#include <cmath>
58
59// Lorene headers
60#include "etoile.h"
61#include "star.h"
62#include "star_bhns.h"
63#include "eos.h"
64#include "unites.h"
65
66// Local prototype
67namespace Lorene {
68Cmp raccord_c1(const Cmp& uu, int l1) ;
69
70 //---------------------//
71 // Constructor //
72 //---------------------//
73
74// Standard constructor
75// --------------------
76Star_bhns::Star_bhns(Map& mp_i, int nzet_i, const Eos& eos_i, bool irrot_i)
77 : Star(mp_i, nzet_i, eos_i),
78 mp_aff(mp_i),
79 irrotational(irrot_i),
80 psi0(mp_i),
81 d_psi(mp_i, COV, mp_i.get_bvect_cart()),
82 wit_w(mp_i, CON, mp_i.get_bvect_cart()),
83 loggam(mp_i),
84 bsn(mp_i, CON, mp_i.get_bvect_cart()),
85 gam(mp_i),
86 gam0(mp_i),
87 pot_centri(mp_i),
88 lapconf_auto(mp_i),
89 lapconf_comp(mp_i),
90 lapconf_tot(mp_i),
91 lapse_auto(mp_i),
92 lapse_tot(mp_i),
93 d_lapconf_auto(mp_i, COV, mp_i.get_bvect_cart()),
94 d_lapconf_comp(mp_i, COV, mp_i.get_bvect_cart()),
95 shift_auto(mp_i, CON, mp_i.get_bvect_cart()),
96 shift_comp(mp_i, CON, mp_i.get_bvect_cart()),
97 shift_tot(mp_i, CON, mp_i.get_bvect_cart()),
98 d_shift_auto(mp_i, 2, CON, mp_i.get_bvect_cart()),
99 d_shift_comp(mp_i, 2, CON, mp_i.get_bvect_cart()),
100 confo_auto(mp_i),
101 confo_comp(mp_i),
102 confo_tot(mp_i),
103 d_confo_auto(mp_i, COV, mp_i.get_bvect_cart()),
104 d_confo_comp(mp_i, COV, mp_i.get_bvect_cart()),
105 psi4(mp_i),
106 taij_auto(mp_i, CON, mp_i.get_bvect_cart()),
107 taij_quad_auto(mp_i),
108 flat(mp_i, mp_i.get_bvect_cart()),
109 ssjm1_lapconf(mp_i),
110 ssjm1_confo(mp_i),
111 ssjm1_khi(mp_i),
112 ssjm1_wshift(mp_i, CON, mp_i.get_bvect_cart())
113{
114
115 // Pointers of derived quantities initialized to zero :
116 set_der_0x0() ;
117
118 // Quantities defined on a spherical triad in star.C are put on
119 // a cartesian one
120 u_euler.change_triad(mp_i.get_bvect_cart()) ;
121
122 // All quantities are initialized to zero
123 psi0 = 0. ;
124 psi0.std_spectral_base() ;
125 d_psi.set_etat_zero() ;
126 wit_w.set_etat_zero() ;
127 loggam = 0. ;
128 loggam.std_spectral_base() ;
129 bsn.set_etat_zero() ;
130 gam = 0. ;
131 gam.std_spectral_base() ;
132 gam0 = 0. ;
133 gam0.std_spectral_base() ;
134 pot_centri = 0. ;
135 pot_centri.std_spectral_base() ;
136
137 lapconf_auto = 1. ;
138 lapconf_auto.std_spectral_base() ;
139 lapconf_comp = 0. ;
140 lapconf_comp.std_spectral_base() ;
141 lapconf_tot = 1. ;
142 lapconf_tot.std_spectral_base() ;
143 lapse_auto = 1. ;
144 lapse_auto.std_spectral_base() ;
145 lapse_tot = 1. ;
146 lapse_tot.std_spectral_base() ;
147 d_lapconf_auto.set_etat_zero() ;
148 d_lapconf_comp.set_etat_zero() ;
149 shift_auto.set_etat_zero() ;
150 shift_comp.set_etat_zero() ;
151 shift_tot.set_etat_zero() ;
152 d_shift_auto.set_etat_zero() ;
153 d_shift_comp.set_etat_zero() ;
154 confo_auto = 1. ;
155 confo_auto.std_spectral_base() ;
156 confo_comp = 0. ;
157 confo_comp.std_spectral_base() ;
158 confo_tot = 1. ;
159 confo_tot.std_spectral_base() ;
160 d_confo_auto.set_etat_zero() ;
161 d_confo_comp.set_etat_zero() ;
162 psi4 = 1. ;
163 psi4.std_spectral_base() ;
164
165 taij_auto.set_etat_zero() ;
166 taij_quad_auto = 0. ;
167 taij_quad_auto.std_spectral_base() ;
168
169 ssjm1_lapconf = 0. ;
170 ssjm1_lapconf.std_spectral_base() ;
171 ssjm1_confo = 0. ;
172 ssjm1_confo.std_spectral_base() ;
173 ssjm1_khi = 0. ;
174 ssjm1_khi.std_spectral_base() ;
175 ssjm1_wshift.set_etat_zero() ;
176
177}
178
179// Copy constructor
180// ----------------
182 : Star(star),
183 mp_aff(star.mp_aff),
185 psi0(star.psi0),
186 d_psi(star.d_psi),
187 wit_w(star.wit_w),
188 loggam(star.loggam),
189 bsn(star.bsn),
190 gam(star.gam),
191 gam0(star.gam0),
197 lapse_tot(star.lapse_tot),
202 shift_tot(star.shift_tot),
207 confo_tot(star.confo_tot),
210 psi4(star.psi4),
211 taij_auto(star.taij_auto),
213 flat(star.flat),
216 ssjm1_khi(star.ssjm1_khi),
218{
219 set_der_0x0() ;
220}
221
222// Constructor from a file
223// -----------------------
224Star_bhns::Star_bhns(Map& mp_i, const Eos& eos_i, FILE* fich)
225 : Star(mp_i, eos_i, fich),
226 mp_aff(mp_i),
227 psi0(mp_i),
228 d_psi(mp_i, COV, mp_i.get_bvect_cart()),
229 wit_w(mp_i, CON, mp_i.get_bvect_cart()),
230 loggam(mp_i),
231 bsn(mp_i, CON, mp_i.get_bvect_cart()),
232 gam(mp_i),
233 gam0(mp_i),
234 pot_centri(mp_i),
235 lapconf_auto(mp_i, *(mp_i.get_mg()), fich),
236 lapconf_comp(mp_i),
237 lapconf_tot(mp_i),
238 lapse_auto(mp_i),
239 lapse_tot(mp_i),
240 d_lapconf_auto(mp_i, COV, mp_i.get_bvect_cart()),
241 d_lapconf_comp(mp_i, COV, mp_i.get_bvect_cart()),
242 shift_auto(mp_i, mp_i.get_bvect_cart(), fich),
243 shift_comp(mp_i, CON, mp_i.get_bvect_cart()),
244 shift_tot(mp_i, CON, mp_i.get_bvect_cart()),
245 d_shift_auto(mp_i, 2, CON, mp_i.get_bvect_cart()),
246 d_shift_comp(mp_i, 2, CON, mp_i.get_bvect_cart()),
247 confo_auto(mp_i, *(mp_i.get_mg()), fich),
248 confo_comp(mp_i),
249 confo_tot(mp_i),
250 d_confo_auto(mp_i, COV, mp_i.get_bvect_cart()),
251 d_confo_comp(mp_i, COV, mp_i.get_bvect_cart()),
252 psi4(mp_i),
253 taij_auto(mp_i, CON, mp_i.get_bvect_cart()),
254 taij_quad_auto(mp_i),
255 flat(mp_i, mp_i.get_bvect_cart()),
256 ssjm1_lapconf(mp_i, *(mp_i.get_mg()), fich),
257 ssjm1_confo(mp_i, *(mp_i.get_mg()), fich),
258 ssjm1_khi(mp_i, *(mp_i.get_mg()), fich),
259 ssjm1_wshift(mp_i, mp_i.get_bvect_cart(), fich)
260{
261
262 // Star parameter
263 fread(&irrotational, sizeof(bool), 1, fich) ;
264
265 // Read of the saved fields
266 // ------------------------
267
268 if (irrotational) {
269 Scalar gam_euler_file(mp, *(mp.get_mg()), fich) ;
270 gam_euler = gam_euler_file ;
271
272 Scalar psi0_file(mp, *(mp.get_mg()), fich) ;
273 psi0 = psi0_file ;
274 }
275
276 // Quantities defined on a spherical triad in star.C are put on
277 // a cartesian one
278 u_euler.change_triad(mp_i.get_bvect_cart()) ;
279
280 // All other quantities are initialized to zero
281 // --------------------------------------------
282
283 d_psi.set_etat_zero() ;
284 wit_w.set_etat_zero() ;
285 loggam = 0. ;
286 loggam.std_spectral_base() ;
287 bsn.set_etat_zero() ;
288 gam = 0. ;
289 gam.std_spectral_base() ;
290 gam0 = 0. ;
291 gam0.std_spectral_base() ;
292 pot_centri = 0. ;
293 pot_centri.std_spectral_base() ;
294
295 lapconf_comp = 0. ;
296 lapconf_comp.std_spectral_base() ;
297 lapconf_tot = 1. ;
298 lapconf_tot.std_spectral_base() ;
299 lapse_auto = 1. ;
300 lapse_auto.std_spectral_base() ;
301 lapse_tot = 1. ;
302 lapse_tot.std_spectral_base() ;
303 d_lapconf_auto.set_etat_zero() ;
304 d_lapconf_comp.set_etat_zero() ;
305 shift_comp.set_etat_zero() ;
306 shift_tot.set_etat_zero() ;
307 d_shift_auto.set_etat_zero() ;
308 d_shift_comp.set_etat_zero() ;
309 confo_comp = 0. ;
310 confo_comp.std_spectral_base() ;
311 confo_tot = 1. ;
312 confo_tot.std_spectral_base() ;
313 d_confo_auto.set_etat_zero() ;
314 d_confo_comp.set_etat_zero() ;
315 psi4 = 1. ;
316 psi4.std_spectral_base() ;
317 taij_auto.set_etat_zero() ;
318 taij_quad_auto = 0. ;
319 taij_quad_auto.std_spectral_base() ;
320
321 // Pointers of derived quantities initialized to zero
322 // --------------------------------------------------
323 set_der_0x0() ;
324
325}
326
327
328 //--------------------//
329 // Destructor //
330 //--------------------//
331
333{
334
335 del_deriv() ;
336
337}
338
339
340 //------------------------------------------//
341 // Management of derived quantities //
342 //------------------------------------------//
343
345
347
348 if (p_mass_b_bhns != 0x0) delete p_mass_b_bhns ;
349 if (p_mass_g_bhns != 0x0) delete p_mass_g_bhns ;
350
351 set_der_0x0() ;
352
353}
354
356
358
359 p_mass_b_bhns = 0x0 ;
360 p_mass_g_bhns = 0x0 ;
361
362}
363
364
365 //--------------------//
366 // Assignment //
367 //--------------------//
368
369// Assignment to another Star_bhns
370// -------------------------------
372
373 // Assignment of quantities common to the derived classes of Star
374 Star::operator=(star) ;
375
376 // Assignment of proper quantities of class Star_bhns
377 mp_aff = star.mp_aff ;
379 psi0 = star.psi0 ;
380 d_psi = star.d_psi ;
381 wit_w = star.wit_w ;
382 loggam = star.loggam ;
383 bsn = star.bsn ;
384 gam = star.gam ;
385 gam0 = star.gam0 ;
386 pot_centri = star.pot_centri ;
389 lapconf_tot = star.lapconf_tot ;
390 lapse_auto = star.lapse_auto ;
391 lapse_tot = star.lapse_tot ;
394 shift_auto = star.shift_auto ;
395 shift_comp = star.shift_comp ;
396 shift_tot = star.shift_tot ;
399 confo_auto = star.confo_auto ;
400 confo_comp = star.confo_comp ;
401 confo_tot = star.confo_tot ;
404 psi4 = star.psi4 ;
405 taij_auto = star.taij_auto ;
407 flat = star.flat ;
409 ssjm1_confo = star.ssjm1_confo ;
410 ssjm1_khi = star.ssjm1_khi ;
412
413 del_deriv() ;
414
415}
416
418
419 del_deriv() ;
420 return pot_centri ;
421
422}
423
425
426 del_deriv() ;
427 return lapconf_auto ;
428
429}
430
432
433 del_deriv() ;
434 return lapconf_comp ;
435
436}
437
439
440 del_deriv() ;
441 return shift_auto ;
442
443}
444
446
447 del_deriv() ;
448 return shift_comp ;
449
450}
451
453
454 del_deriv() ;
455 return confo_auto ;
456
457}
458
460
461 del_deriv() ;
462 return confo_comp ;
463
464}
465
466
467 //-----------------//
468 // Outputs //
469 //-----------------//
470
471// Save in a file
472// --------------
473void Star_bhns::sauve(FILE* fich) const {
474
475 Star::sauve(fich) ;
476
477 lapconf_auto.sauve(fich) ;
478 shift_auto.sauve(fich) ;
479 confo_auto.sauve(fich) ;
480
481 ssjm1_lapconf.sauve(fich) ;
482 ssjm1_confo.sauve(fich) ;
483 ssjm1_khi.sauve(fich) ;
484 ssjm1_wshift.sauve(fich) ;
485
486 fwrite(&irrotational, sizeof(bool), 1, fich) ;
487
488 if (irrotational) {
489 gam_euler.sauve(fich) ; // required to construct d_psi from psi0
490 psi0.sauve(fich) ;
491 }
492
493}
494
495// Printing
496// --------
497
498ostream& Star_bhns::operator>>(ostream& ost) const {
499
500 using namespace Unites ;
501
502 // Star::operator>>(ost) ;
503
504 ost << endl ;
505 ost << "Neutron star in a BHNS binary" << endl ;
506 ost << "-----------------------------" << endl ;
507
508 ost << "Coordinate radius R_eq_tow : "
509 << ray_eq_pi() / km << " [km]" << endl ;
510 ost << "Coordinate radius R_eq_opp : "
511 << ray_eq() / km << " [km]" << endl ;
512 ost << "Coordinate radius R_eq_orb : "
513 << ray_eq_pis2() / km << " [km]" << endl ;
514 ost << "Coordinate radius R_pole : "
515 << ray_pole() / km << " [km]" << endl ;
516 ost << "Central enthalph H_ent : "
517 << ent.val_grid_point(0,0,0,0) << endl ;
518 ost << "Lapse function at the center of NS : "
519 << lapse_tot.val_grid_point(0,0,0,0) << endl ;
520 ost << "Conformal factor at the center of NS : "
521 << confo_tot.val_grid_point(0,0,0,0) << endl ;
522 ost << "shift(1) at the center of NS : "
523 << shift_tot(1).val_grid_point(0,0,0,0) << endl ;
524 ost << "shift(2) at the center of NS : "
525 << shift_tot(2).val_grid_point(0,0,0,0) << endl ;
526 ost << "shift(3) at the center of NS : "
527 << shift_tot(3).val_grid_point(0,0,0,0) << endl ;
528
529 return ost ;
530
531}
532
533
534 //--------------------------------//
535 // Computational routines //
536 //--------------------------------//
537
539
540 if (!irrotational) {
541 d_psi.set_etat_nondef() ;
542 return ;
543 }
544
545 // Specific relativistic enthalpy ---> hhh
546 // ------------------------------
547
548 Scalar hhh = exp(ent) ; // = 1 at the Newtonian limit
549
550 // Computation of W^i = h Gamma_n B^i / N
551 // --------------------------------------
552
553 Vector www = hhh * gam_euler * bsn * psi4 ;
554
555 // Constant value of W^i at the center of the star
556 // -----------------------------------------------
557
558 Vector v_orb(mp, COV, mp.get_bvect_cart()) ;
559
560 for (int i=1; i<=3; i++) {
561 v_orb.set(i) = www(i).val_grid_point(0,0,0,0) ;
562 }
563
564 // Gradient of psi
565 // ---------------
566
567 Vector d_psi0(mp, COV, mp.get_bvect_cart()) ;
568 d_psi0.set_etat_qcq() ;
569 for (int i=1; i<=3; i++)
570 d_psi0.set(i) = psi0.deriv(i) ;
571
572 d_psi0.std_spectral_base() ;
573
574 d_psi = d_psi0 + v_orb ;
575
576 for (int i=1; i<=3; i++) {
577 if (d_psi(i).get_etat() == ETATZERO)
578 d_psi.set(i).annule_hard() ;
579 }
580
581 // C^1 continuation of d_psi outside the star
582 // (to ensure a smooth enthalpy field accross the stellar surface)
583 // ---------------------------------------------------------------
584
585 int nzm1 = mp.get_mg()->get_nzone() - 1 ;
586 d_psi.annule(nzet, nzm1) ;
587 for (int i=1; i<=3; i++) {
588 Cmp d_psi_i (d_psi(i)) ;
589 d_psi_i.va.set_base( d_psi0(i).get_spectral_va().base ) ;
590 d_psi_i = raccord_c1(d_psi_i, nzet) ;
591 d_psi.set(i) = d_psi_i ;
592 }
593
594}
595
596
597void Star_bhns::relax_bhns(const Star_bhns& star_prev, double relax_ent,
598 double relax_met, int mer, int fmer_met) {
599
600 double relax_ent_jm1 = 1. - relax_ent ;
601 double relax_met_jm1 = 1. - relax_met ;
602
603 ent = relax_ent * ent + relax_ent_jm1 * star_prev.ent ;
604
605 if ( (mer != 0) && (mer % fmer_met == 0)) {
606
607 lapconf_auto = relax_met * lapconf_auto
608 + relax_met_jm1 * star_prev.lapconf_auto ;
609
610 shift_auto = relax_met * shift_auto
611 + relax_met_jm1 * star_prev.shift_auto ;
612
613 confo_auto = relax_met * confo_auto
614 + relax_met_jm1 * star_prev.confo_auto ;
615
616 }
617
618 del_deriv() ;
619
621
622}
623}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Valeur va
The numerical value of the Cmp.
Definition cmp.h:464
Equation of state base class.
Definition eos.h:206
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star ).
Definition star_bhns.h:192
Scalar & set_confo_comp()
Read/write of the conformal factor generated by the companion black hole.
Definition star_bhns.C:459
virtual void sauve(FILE *) const
Save in a file.
Definition star_bhns.C:473
Vector shift_tot
Total shift vector.
Definition star_bhns.h:144
Sym_tensor taij_auto
Part of the extrinsic curvature tensor generated by shift_auto , lapse_auto , and confo_auto .
Definition star_bhns.h:182
Scalar lapconf_auto
Lapconf function generated by the star.
Definition star_bhns.h:113
Scalar confo_comp
Conformal factor generated by the companion black hole.
Definition star_bhns.h:160
Vector wit_w
Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Definition star_bhns.h:88
Tensor d_shift_auto
Derivative of the shift vector generated by the star .
Definition star_bhns.h:149
Vector bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition star_bhns.h:99
Scalar & set_confo_auto()
Read/write of the conformal factor generated by the neutron star.
Definition star_bhns.C:452
Vector & set_shift_auto()
Read/write of the shift vector generated by the neutron star.
Definition star_bhns.C:438
Scalar & set_lapconf_auto()
Read/write of the lapconf function generated by the neutron star.
Definition star_bhns.C:424
Scalar & set_pot_centri()
Read/write the centrifugal potential.
Definition star_bhns.C:417
void fait_d_psi_bhns()
Computes the gradient of the total velocity potential .
Definition star_bhns.C:538
Scalar lapse_tot
Total lapse function.
Definition star_bhns.h:125
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Definition star_bhns.h:220
virtual void del_deriv() const
Deletes all the derived quantities.
Definition star_bhns.C:344
Vector d_confo_comp
Derivative of the conformal factor generated by the companion black hole.
Definition star_bhns.h:173
void operator=(const Star_bhns &)
Assignment to another Star_bhns.
Definition star_bhns.C:371
Scalar confo_tot
Total conformal factor.
Definition star_bhns.h:163
Scalar psi4
Fourth power of the total conformal factor.
Definition star_bhns.h:176
Vector d_lapconf_auto
Derivative of the lapconf function generated by the star .
Definition star_bhns.h:130
Vector d_psi
Gradient of (in the irrotational case) (Spherical components with respect to the mapping of the star...
Definition star_bhns.h:82
Vector d_confo_auto
Derivative of the conformal factor generated by the star .
Definition star_bhns.h:168
Scalar & set_lapconf_comp()
Read/write of the lapconf function generated by the companion black hole.
Definition star_bhns.C:431
Scalar gam0
Lorentz factor between the co-orbiting observer and the Eulerian one.
Definition star_bhns.h:107
Vector shift_comp
Shift vector generated by the companion black hole.
Definition star_bhns.h:141
Scalar lapconf_tot
Total lapconf function.
Definition star_bhns.h:119
Star_bhns(Map &mp_i, int nzet_i, const Eos &eos_i, bool irrot_i)
Standard constructor.
Definition star_bhns.C:76
Vector d_lapconf_comp
Derivative of the lapconf function generated by the companion black hole.
Definition star_bhns.h:135
Vector & set_shift_comp()
Read/write of the shift vector generated by the companion black hole.
Definition star_bhns.C:445
Scalar gam
Lorentz factor between the fluid and the co-orbiting observer.
Definition star_bhns.h:102
Scalar ssjm1_lapconf
Effective source at the previous step for the resolution of the Poisson equation for lapconf_auto .
Definition star_bhns.h:197
bool irrotational
true for an irrotational star, false for a corotating one
Definition star_bhns.h:72
Scalar pot_centri
Centrifugal potential.
Definition star_bhns.h:110
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition star_bhns.C:355
double * p_mass_b_bhns
Baryon mass.
Definition star_bhns.h:225
Scalar lapse_auto
Lapse function generated by the "star".
Definition star_bhns.h:122
void relax_bhns(const Star_bhns &star_prev, double relax_ent, double relax_met, int mer, int fmer_met)
Performs a relaxation on ent , lapse_auto , shift_auto , confo_auto .
Definition star_bhns.C:597
Scalar psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case).
Definition star_bhns.h:77
Scalar lapconf_comp
Lapconf function generated by the companion black hole.
Definition star_bhns.h:116
Scalar taij_quad_auto
Part of the scalar generated by .
Definition star_bhns.h:187
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition star_bhns.C:498
Vector shift_auto
Shift vector generated by the star.
Definition star_bhns.h:138
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition star_bhns.h:210
Scalar confo_auto
Conformal factor generated by the star.
Definition star_bhns.h:157
Scalar ssjm1_confo
Effective source at the previous step for the resolution of the Poisson equation for confo_auto .
Definition star_bhns.h:202
Map_af mp_aff
Affine mapping for solving poisson's equations of metric quantities.
Definition star_bhns.h:67
Tensor d_shift_comp
Derivative of the shift vector generated by the companion black hole.
Definition star_bhns.h:154
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition star_bhns.h:93
virtual ~Star_bhns()
Destructor.
Definition star_bhns.C:332
double * p_mass_g_bhns
Gravitational mass.
Definition star_bhns.h:226
Star(Map &mp_i, int nzet_i, const Eos &eos_i)
Standard constructor.
Definition star.C:127
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition star.C:465
virtual void del_deriv() const
Deletes all the derived quantities.
Definition star.C:301
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition star.h:204
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition star.C:319
double ray_eq() const
Coordinate radius at , [r_unit].
virtual void sauve(FILE *) const
Save in a file.
Definition star.C:400
double ray_eq_pis2() const
Coordinate radius at , [r_unit].
double ray_eq_pi() const
Coordinate radius at , [r_unit].
Scalar ent
Log-enthalpy.
Definition star.h:190
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition star.h:207
Map & mp
Mapping associated with the star.
Definition star.h:180
void operator=(const Star &)
Assignment to another Star.
Definition star.C:354
int nzet
Number of domains of *mp occupied by the star.
Definition star.h:183
double ray_pole() const
Coordinate radius at [r_unit].
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base ).
Definition valeur.C:813
Tensor field of valence 1.
Definition vector.h:188
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition vector.C:319
Scalar & set(int)
Read/write access to a component.
Definition vector.C:299
Cmp exp(const Cmp &)
Exponential.
Definition cmp_math.C:273
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition tensor.C:490
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 Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:777
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.