LORENE
et_bin_nsbh.C
1/*
2 * Methods for the class Et_bin_nsbh
3 *
4 * (see file et_bin_nsbh.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2003 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: et_bin_nsbh.C,v 1.13 2016/12/05 16:17:53 j_novak Exp $
32 * $Log: et_bin_nsbh.C,v $
33 * Revision 1.13 2016/12/05 16:17:53 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.12 2014/10/13 08:52:56 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.11 2014/10/06 15:13:08 j_novak
40 * Modified #include directives to use c++ syntax.
41 *
42 * Revision 1.10 2006/09/25 10:01:49 p_grandclement
43 * Addition of N-dimensional Tbl
44 *
45 * Revision 1.9 2006/06/01 12:47:53 p_grandclement
46 * update of the Bin_ns_bh project
47 *
48 * Revision 1.8 2005/08/29 15:10:16 p_grandclement
49 * Addition of things needed :
50 * 1) For BBH with different masses
51 * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
52 * WORKING YET !!!
53 *
54 * Revision 1.7 2004/06/09 06:39:35 k_taniguchi
55 * Introduce set_n_auto() and set_confpsi_auto().
56 *
57 * Revision 1.6 2004/03/25 10:29:04 j_novak
58 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
59 *
60 * Revision 1.5 2003/12/16 05:30:18 k_taniguchi
61 * Some changes in the constructor from file and the saving to file.
62 *
63 * Revision 1.4 2003/11/28 04:30:22 k_taniguchi
64 * Add the output operator.
65 *
66 * Revision 1.3 2003/10/24 12:30:58 k_taniguchi
67 * Add some derivative terms
68 *
69 * Revision 1.2 2003/10/22 08:12:46 k_taniguchi
70 * Supress some terms in Et_bin_nsbh::sauve.
71 *
72 * Revision 1.1 2003/10/21 11:50:38 k_taniguchi
73 * Methods for class Et_bin_nsbh
74 *
75 *
76 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh.C,v 1.13 2016/12/05 16:17:53 j_novak Exp $
77 *
78 */
79
80// C headers
81#include <cmath>
82
83// Lorene headers
84#include "et_bin_nsbh.h"
85#include "etoile.h"
86#include "eos.h"
87#include "utilitaires.h"
88#include "param.h"
89#include "unites.h"
90
91 //--------------//
92 // Constructors //
93 //--------------//
94
95// Standard constructor
96// --------------------
97namespace Lorene {
98Et_bin_nsbh::Et_bin_nsbh(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i,
99 bool irrot, const Base_vect& ref_triad_i)
100 : Etoile_bin(mp_i, nzet_i, relat, eos_i, irrot, ref_triad_i),
101 n_auto(mp_i),
102 n_comp(mp_i),
103 d_n_auto(mp_i, 1, COV, ref_triad_i),
104 d_n_comp(mp_i, 1, COV, ref_triad_i),
105 confpsi(mp_i),
106 confpsi_auto(mp_i),
107 confpsi_comp(mp_i),
108 d_confpsi_auto(mp_i, 1, COV, ref_triad_i),
109 d_confpsi_comp(mp_i, 1, COV, ref_triad_i),
110 taij_auto(mp_i, 2, CON, ref_triad_i) ,
111 taij_comp(mp_i, 2, CON, ref_triad_i),
112 taij_tot(mp_i, 2, CON, ref_triad_i) ,
113 tkij_auto(mp_i, 2, CON, ref_triad_i),
114 tkij_tot(mp_i, 2, CON, ref_triad_i),
115 ssjm1_lapse(mp_i),
116 ssjm1_confpsi(mp_i) {
117
118 // Pointers of derived quantities initialized to zero :
119 set_der_0x0() ;
120
121 // The metric is initialized to the flat one :
122 n_auto = 0.5 ;
123 n_auto.set_std_base() ;
124 n_comp = 0.5 ;
125 n_comp.set_std_base() ;
126 d_n_auto = 0 ;
127 d_n_comp = 0 ;
128 confpsi = 1 ;
129 confpsi.set_std_base() ;
130 confpsi_auto = 0.5 ;
131 confpsi_auto.set_std_base() ;
132 confpsi_comp = 0.5 ;
133 confpsi_comp.set_std_base() ;
134 d_confpsi_auto = 0 ;
135 d_confpsi_comp = 0 ;
136
137 taij_auto.set_etat_zero() ;
138 taij_comp.set_etat_zero() ;
139 taij_tot.set_etat_zero() ;
140 tkij_auto.set_etat_zero() ;
141 tkij_tot.set_etat_zero() ;
142
143 ssjm1_lapse.set_etat_qcq() ;
144 ssjm1_lapse = 0. ;
145 ssjm1_confpsi.set_etat_qcq() ;
146 ssjm1_confpsi = 0. ;
147
148}
149
150// Copy constructor
151// ----------------
153 : Etoile_bin(et),
154 n_auto(et.n_auto),
155 n_comp(et.n_comp),
156 d_n_auto(et.d_n_auto),
157 d_n_comp(et.d_n_comp),
158 confpsi(et.confpsi),
159 confpsi_auto(et.confpsi_auto),
160 confpsi_comp(et.confpsi_comp),
161 d_confpsi_auto(et.d_confpsi_auto),
162 d_confpsi_comp(et.d_confpsi_comp),
163 taij_auto(et.taij_auto),
164 taij_comp(et.taij_comp),
165 taij_tot(et.taij_tot),
166 tkij_auto(et.tkij_auto),
167 tkij_tot(et.tkij_tot),
168 ssjm1_lapse(et.ssjm1_lapse),
169 ssjm1_confpsi(et.ssjm1_confpsi) {
170
171 set_der_0x0() ;
172}
173
174// Constructor from a file
175// -----------------------
176Et_bin_nsbh::Et_bin_nsbh(Map& mp_i, const Eos& eos_i,
177 const Base_vect& ref_triad_i, FILE* fich, bool old)
178 : Etoile_bin(mp_i, eos_i, ref_triad_i, fich),
179 n_auto(mp_i),
180 n_comp(mp_i),
181 d_n_auto(mp_i, 1, COV, ref_triad_i),
182 d_n_comp(mp_i, 1, COV, ref_triad_i ),
183 confpsi(mp_i),
184 confpsi_auto(mp_i),
185 confpsi_comp(mp_i),
186 d_confpsi_auto(mp_i, 1, COV, ref_triad_i),
187 d_confpsi_comp(mp_i, 1, COV, ref_triad_i),
188 taij_auto(mp_i, 2, CON, ref_triad_i),
189 taij_comp(mp_i, 2, CON, ref_triad_i),
190 taij_tot(mp_i, 2, CON, ref_triad_i),
191 tkij_auto(mp_i, 2, CON, ref_triad_i),
192 tkij_tot(mp_i, 2, CON, ref_triad_i),
193 ssjm1_lapse(mp_i),
194 ssjm1_confpsi(mp_i) {
195
196 // Construct from data in Etoile_bin
197 Cmp n_from_file (mp_i, *(mp_i.get_mg()), fich) ;
198 n_auto.set_etat_qcq() ;
199 n_auto.set() = n_from_file ;
200
201 Cmp psi_from_file (mp_i, *(mp_i.get_mg()), fich) ;
202 confpsi_auto.set_etat_qcq() ;
203 confpsi_auto.set() = psi_from_file ;
204
205 if (!old) {
206 Tenseur shift_auto_file (mp, ref_triad_i, fich) ;
207 shift_auto = shift_auto_file ;
208 }
209
210 // All other fields are initialized to zero or some constants :
211 // ----------------------------------------------------------
212 n_comp = 0.5 ;
213 n_comp.set_std_base() ;
214 d_n_auto = 0 ;
215 d_n_comp = 0 ;
216 confpsi = 1 ;
217 confpsi.set_std_base() ;
218 confpsi_comp = 0.5 ;
219 confpsi_comp.set_std_base() ;
220 d_confpsi_auto = 0 ;
221 d_confpsi_comp = 0 ;
222 taij_auto.set_etat_zero() ;
223 taij_comp.set_etat_zero() ;
224 taij_tot.set_etat_zero() ;
225 tkij_auto.set_etat_zero() ;
226 tkij_tot.set_etat_zero() ;
227
228 // Read of the saved fields :
229 Cmp ssjm1_lapse_file(mp_i, *(mp_i.get_mg()), fich) ;
230 ssjm1_lapse = ssjm1_lapse_file ;
231
232 Cmp ssjm1_confpsi_file(mp_i, *(mp_i.get_mg()), fich) ;
233 ssjm1_confpsi = ssjm1_confpsi_file ;
234
235 // Pointers of derived quantities initialized to zero
236 // --------------------------------------------------
237 set_der_0x0() ;
238}
239
240 //------------//
241 // Destructor //
242 //------------//
243
244Et_bin_nsbh::~Et_bin_nsbh(){
245
246 del_deriv() ;
247
248}
249
250 //--------------//
251 // Assignment //
252 //--------------//
253
254// Assignment to another Et_bin_nsbh
255// ---------------------------------
257
258 // Assignment of quantities common to the derived classes of Etoile_bin
260
261 n_auto = et.n_auto ;
262 n_comp = et.n_comp ;
263 d_n_auto = et.d_n_auto ;
264 d_n_comp = et.d_n_comp ;
265 confpsi = et.confpsi ;
270 taij_auto = et.taij_auto ;
271 taij_comp = et.taij_comp ;
272 taij_tot = et.taij_tot ;
273 tkij_auto = et.tkij_auto ;
274 tkij_tot = et.tkij_tot ;
277
278 del_deriv() ; // Deletes all derived quantities
279}
280
282
283 del_deriv() ; // sets to 0x0 all the derived quantities
284 return n_auto ;
285
286}
287
289
290 del_deriv() ; // sets to 0x0 all the derived quantities
291 return n_comp ;
292
293}
294
296
297 del_deriv() ; // sets to 0x0 all the derived quantities
298 return confpsi_auto ;
299
300}
301
303
304 del_deriv() ; // sets to 0x0 all the derived quantities
305 return confpsi_comp ;
306
307}
308
309 //--------------//
310 // Outputs //
311 //--------------//
312
313// Save in a file
314// --------------
315void Et_bin_nsbh::sauve(FILE* fich) const {
316
317 Etoile_bin::sauve(fich) ;
318
319 n_auto().sauve(fich) ;
320 confpsi_auto().sauve(fich) ;
321 shift_auto.sauve(fich) ;
322
323 ssjm1_lapse.sauve(fich) ;
324 ssjm1_confpsi.sauve(fich) ;
325}
326
327
328// Printing
329// --------
330
331ostream& Et_bin_nsbh::operator>>(ostream& ost) const {
332
333 using namespace Unites ;
334
335 Etoile::operator>>(ost) ;
336
337 ost << endl ;
338 ost << "Neutron star in a binary system" << endl ;
339 ost << "-------------------------------" << endl ;
340
341 if (irrotational) {
342 ost << "irrotational configuration" << endl ;
343 }
344 else {
345 ost << "corotating configuration" << endl ;
346 }
347
348 ost << "Absolute abscidia of the stellar center: " <<
349 mp.get_ori_x() / km << " km" << endl ;
350
351 ost << "Absolute abscidia of the barycenter of the baryon density : " <<
352 xa_barycenter() / km << " km" << endl ;
353
354 double r_0 = 0.5 * ( ray_eq() + ray_eq_pi() ) ;
355 double d_ns = fabs( mp.get_ori_x() ) + ray_eq_pi() - r_0 ;
356 double d_tilde = 2 * d_ns / r_0 ;
357
358 ost << "d_tilde : " << d_tilde << endl ;
359
360 ost << "Orientation with respect to the absolute frame : " <<
361 mp.get_rot_phi() << " rad" << endl ;
362
363 ost << "Central value of gam_euler : "
364 << gam_euler()(0, 0, 0, 0) << endl ;
365
366 ost << "Central u_euler (U^X, U^Y, U^Z) [c] : "
367 << u_euler(0)(0, 0, 0, 0) << " "
368 << u_euler(1)(0, 0, 0, 0) << " "
369 << u_euler(2)(0, 0, 0, 0) << endl ;
370
371 if (irrotational) {
372 ost << "Central d_psi (X, Y, Z) [c] : "
373 << d_psi(0)(0, 0, 0, 0) << " "
374 << d_psi(1)(0, 0, 0, 0) << " "
375 << d_psi(2)(0, 0, 0, 0) << endl ;
376
377 ost << "Central vel. / co-orb. (W^X, W^Y, W^Z) [c] : "
378 << wit_w(0)(0, 0, 0, 0) << " "
379 << wit_w(1)(0, 0, 0, 0) << " "
380 << wit_w(2)(0, 0, 0, 0) << endl ;
381
382 ost << "Max vel. / co-orb. (W^X, W^Y, W^Z) [c] : "
383 << max(max(wit_w(0))) << " "
384 << max(max(wit_w(1))) << " "
385 << max(max(wit_w(2))) << endl ;
386
387 ost << "Min vel. / co-orb. (W^X, W^Y, W^Z) [c] : "
388 << min(min(wit_w(0))) << " "
389 << min(min(wit_w(1))) << " "
390 << min(min(wit_w(2))) << endl ;
391
392 double r_surf = mp.val_r(0,1.,M_PI/4,M_PI/4) ;
393
394 ost << "Velocity at (r_surf,pi/4,pi/4) / co-orb. [c] : "
395 << wit_w(0).val_point(r_surf,M_PI/4,M_PI/4) << " "
396 << wit_w(1).val_point(r_surf,M_PI/4,M_PI/4) << " "
397 << wit_w(2).val_point(r_surf,M_PI/4,M_PI/4) << endl ;
398
399 ost << "Central value of loggam : "
400 << loggam()(0, 0, 0, 0) << endl ;
401 }
402
403 ost << "Central value of lapse(N) auto : "
404 << n_auto()(0, 0, 0, 0) << endl ;
405
406 ost << "Central value of confpsi auto : "
407 << confpsi_auto()(0, 0, 0, 0) << endl ;
408
409 ost << "Central value of shift (N^X, N^Y, N^Z) [c] : "
410 << shift(0)(0, 0, 0, 0) << " "
411 << shift(1)(0, 0, 0, 0) << " "
412 << shift(2)(0, 0, 0, 0) << endl ;
413
414 ost << " ... shift_auto part of it [c] : "
415 << shift_auto(0)(0, 0, 0, 0) << " "
416 << shift_auto(1)(0, 0, 0, 0) << " "
417 << shift_auto(2)(0, 0, 0, 0) << endl ;
418
419 ost << " ... shift_comp part of it [c] : "
420 << shift_comp(0)(0, 0, 0, 0) << " "
421 << shift_comp(1)(0, 0, 0, 0) << " "
422 << shift_comp(2)(0, 0, 0, 0) << endl ;
423
424 ost << " ... w_shift (NB: components in the star Cartesian frame) [c] : "
425 << endl
426 << w_shift(0)(0, 0, 0, 0) << " "
427 << w_shift(1)(0, 0, 0, 0) << " "
428 << w_shift(2)(0, 0, 0, 0) << endl ;
429
430 ost << "Central value of khi_shift [km c] : "
431 << khi_shift()(0, 0, 0, 0) / km << endl ;
432
433 ost << endl << "Central value of (B^X, B^Y, B^Z)/N [c] : "
434 << bsn(0)(0, 0, 0, 0) << " "
435 << bsn(1)(0, 0, 0, 0) << " "
436 << bsn(2)(0, 0, 0, 0) << endl ;
437
438 ost << endl <<
439 "Central (d/dX,d/dY,d/dZ)(logn_auto) [km^{-1}] : "
440 << d_n_auto(0)(0, 0, 0, 0) * km << " "
441 << d_n_auto(1)(0, 0, 0, 0) * km << " "
442 << d_n_auto(2)(0, 0, 0, 0) * km << endl ;
443
444 ost << endl <<
445 "Central (d/dX,d/dY,d/dZ)(beta_auto) [km^{-1}] : "
446 << d_confpsi_auto(0)(0, 0, 0, 0) * km << " "
447 << d_confpsi_auto(1)(0, 0, 0, 0) * km << " "
448 << d_confpsi_auto(2)(0, 0, 0, 0) * km << endl ;
449
450 return ost ;
451}
452}
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition base_vect.h:105
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Equation of state base class.
Definition eos.h:206
Class for a star in a NS-BH binary system.
Definition et_bin_nsbh.h:79
Tenseur confpsi_auto
Part of the conformal factor $\Psi$ generated principaly by the star.
Tenseur d_confpsi_comp
Gradient of {\tt confpsi_comp} (Cartesian components with respect to {\tt ref_triad}...
Cmp ssjm1_confpsi
Effective source at the previous step for the resolution of the Poisson equation for {\tt confpsi_aut...
virtual ostream & operator>>(ostream &) const
Save in a file.
Tenseur confpsi_comp
Part of the conformal factor $\Psi$ generated principaly by the companion star.
Tenseur_sym tkij_tot
Total extrinsic curvature tensor $K^{ij}$ generated by {\tt shift_auto} and {\tt shift_comp}...
Tenseur d_n_comp
Gradient of {\tt n_comp} (Cartesian components with respect to {\tt ref_triad}).
Definition et_bin_nsbh.h:98
Tenseur & set_n_comp()
Read/write the lapse {\it N} generated principaly by the companion star.
Tenseur_sym taij_tot
Total extrinsic curvature tensor $\tilde A^{ij} = 2 N K^{ij}$ generated by {\tt shift_auto}...
Tenseur & set_n_auto()
Read/write the lapse {\it N} generated principaly by the star.
Tenseur d_n_auto
Gradient of {\tt n_auto} (Cartesian components with respect to {\tt ref_triad}).
Definition et_bin_nsbh.h:93
void operator=(const Et_bin_nsbh &)
Destructor.
Cmp ssjm1_lapse
Effective source at the previous step for the resolution of the Poisson equation for {\tt n_auto}...
Et_bin_nsbh(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, bool irrot, const Base_vect &ref_triad_i)
Standard constructor.
Definition et_bin_nsbh.C:98
Tenseur_sym taij_auto
Part of the extrinsic curvature tensor $\tilde A^{ij} = 2 N K^{ij}$ generated by {\tt shift_auto}...
Tenseur & set_confpsi_comp()
Read/write the conformal factor $\Psi$ generated principaly by the companion star.
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor $K^{ij}$ generated by {\tt shift_auto}.
Tenseur & set_confpsi_auto()
Read/write the conformal factor $\Psi$ generated principaly by the star.
Tenseur_sym taij_comp
Part of the extrinsic curvature tensor $\tilde A^{ij} = 2 N K^{ij}$ generated by {\tt shift_comp}...
Tenseur n_auto
Part of the lapse {\it N} generated principaly by the star.
Definition et_bin_nsbh.h:85
Tenseur n_comp
Part of the lapse {\it N} generated principaly by the companion star.
Definition et_bin_nsbh.h:88
Tenseur confpsi
Total conformal factor $\Psi$.
virtual void sauve(FILE *) const
Save in a file.
Tenseur d_confpsi_auto
Gradient of {\tt confpsi_auto} (Cartesian components with respect to {\tt ref_triad}...
Class for stars in binary system.
Definition etoile.h:817
Tenseur shift_comp
Part of the shift vector generated principaly by the companion star.
Definition etoile.h:898
virtual void sauve(FILE *) const
Save in a file.
Definition etoile_bin.C:562
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition etoile.h:953
Tenseur d_psi
Gradient of (in the irrotational case) (Cartesian components with respect to ref_triad ).
Definition etoile.h:841
bool irrotational
true for an irrotational star, false for a corotating one
Definition etoile.h:825
Tenseur w_shift
Vector used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:911
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition etoile_bin.C:462
virtual void del_deriv() const
Deletes all the derived quantities.
Definition etoile_bin.C:450
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition etoile.h:852
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition etoile.h:892
virtual double xa_barycenter() const
Absolute coordinate X of the barycenter of the baryon density, defined according to the formula.
Tenseur wit_w
Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Definition etoile.h:847
Etoile_bin(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, bool irrot, const Base_vect &ref_triad_i)
Standard constructor.
Definition etoile_bin.C:210
void operator=(const Etoile_bin &)
Assignment to another Etoile_bin.
Definition etoile_bin.C:485
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:921
double ray_eq_pi() const
Coordinate radius at , [r_unit].
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition etoile.h:477
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:474
Map & mp
Mapping associated with the star.
Definition etoile.h:432
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition etoile.C:514
Tenseur shift
Total shift vector.
Definition etoile.h:515
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
void sauve(FILE *) const
Save in a file.
Definition tenseur.C:1331
Tbl min(const Cmp &)
Minimum values of a Cmp in each domain.
Definition cmp_math.C:461
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition cmp_math.C:438
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.