LORENE
star_bin.C
1  /*
2  * Methods for the class Star_bin
3  *
4  * (see file star.h for documentation)
5  */
6 
7 /*
8  * Copyright (c) 2004 Francois Limousin
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: star_bin.C,v 1.20 2016/12/05 16:18:14 j_novak Exp $
33  * $Log: star_bin.C,v $
34  * Revision 1.20 2016/12/05 16:18:14 j_novak
35  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36  *
37  * Revision 1.19 2014/10/13 08:53:37 j_novak
38  * Lorene classes and functions now belong to the namespace Lorene.
39  *
40  * Revision 1.18 2006/04/11 14:24:44 f_limousin
41  * New version of the code : improvement of the computation of some
42  * critical sources, estimation of the dirac gauge, helical symmetry...
43  *
44  * Revision 1.17 2005/09/13 19:38:31 f_limousin
45  * Reintroduction of the resolution of the equations in cartesian coordinates.
46  *
47  * Revision 1.16 2005/04/08 12:36:44 f_limousin
48  * Just to avoid warnings...
49  *
50  * Revision 1.15 2005/02/24 16:03:01 f_limousin
51  * Change the name of some variables (for instance dcov_logn --> dlogn).
52  *
53  * Revision 1.14 2005/02/17 17:29:28 f_limousin
54  * Change the name of some quantities to be consistent with other classes
55  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
56  *
57  * Revision 1.13 2005/02/11 18:11:57 f_limousin
58  * Introduction of a new member Map_af.
59  *
60  * Revision 1.12 2004/11/11 16:29:49 j_novak
61  * set_der_0x0 is no longer virtual (to be coherent with Tensor/Scalar classes).
62  *
63  * Revision 1.11 2004/07/21 11:49:03 f_limousin
64  * Remove function sprod.
65  *
66  * Revision 1.10 2004/06/22 12:48:52 f_limousin
67  * Change qq, qq_auto and qq_comp to beta, beta_auto and beta_comp.
68  *
69  * Revision 1.9 2004/04/08 16:32:28 f_limousin
70  * The new variable is ln(Q) instead of Q=psi^2*N. It improves the
71  * convergence of the code.
72  *
73  * Revision 1.8 2004/03/25 10:29:26 j_novak
74  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
75  *
76  * Revision 1.7 2004/03/23 09:54:54 f_limousin
77  * Add comments
78  *
79  * Revision 1.6 2004/02/27 09:50:57 f_limousin
80  * Scalars ssjm1_logn, ssjm1_qq ... have been added for all metric
81  * quantities for the resolution of Poisson equations.
82  * Members bsn and d_psi are now constructed on a cartesian triad.
83  *
84  * Revision 1.5 2004/02/21 17:05:13 e_gourgoulhon
85  * Method Scalar::point renamed Scalar::val_grid_point.
86  * Method Scalar::set_point renamed Scalar::set_grid_point.
87  *
88  * Revision 1.4 2004/01/22 10:07:18 f_limousin
89  * Add methods set_logn_comp() and set_shift_auto().
90  *
91  * Revision 1.3 2004/01/20 15:17:34 f_limousin
92  * First version
93  *
94  *
95  * $Header: /cvsroot/Lorene/C++/Source/Star/star_bin.C,v 1.20 2016/12/05 16:18:14 j_novak Exp $
96  *
97  */
98 
99 // Headers C
100 #include "math.h"
101 
102 // Headers Lorene
103 #include "etoile.h"
104 #include "star.h"
105 #include "eos.h"
106 #include "unites.h"
107 
108 // Local prototype
109 namespace Lorene {
110 Cmp raccord_c1(const Cmp& uu, int l1) ;
111 
112  //--------------//
113  // Constructors //
114  //--------------//
115 
116 // Standard constructor
117 // --------------------
118 Star_bin::Star_bin(Map& mpi, int nzet_i, const Eos& eos_i,
119  bool irrot, bool conf_flat0)
120  : Star(mpi, nzet_i, eos_i),
121  irrotational(irrot),
122  psi0(mpi),
123  d_psi(mpi, COV, mpi.get_bvect_cart()),
124  wit_w(mpi, CON, mpi.get_bvect_cart()),
125  loggam(mpi),
126  bsn(mpi, CON, mpi.get_bvect_cart()),
127  pot_centri(mpi),
128  logn_auto(mpi),
129  logn_comp(mpi),
130  dcov_logn(mpi, COV, mpi.get_bvect_cart()),
131  dcon_logn(mpi, CON, mpi.get_bvect_cart()),
132  lnq_auto(mpi),
133  lnq_comp(mpi),
134  psi4(mpi),
135  dcov_phi(mpi, COV, mpi.get_bvect_cart()),
136  dcon_phi(mpi, CON, mpi.get_bvect_cart()),
137  flat(mpi, mpi.get_bvect_cart()),
138  gtilde(flat),
139  beta_auto(mpi, CON, mpi.get_bvect_cart()),
140  beta_comp(mpi, CON, mpi.get_bvect_cart()),
141  hij(mpi, CON, mpi.get_bvect_cart()),
142  hij_auto(mpi, CON, mpi.get_bvect_cart()),
143  hij_comp(mpi, CON, mpi.get_bvect_cart()),
144  tkij_auto(mpi, CON, mpi.get_bvect_cart()),
145  tkij_comp(mpi, CON, mpi.get_bvect_cart()),
146  kcar_auto(mpi),
147  kcar_comp(mpi),
148  ssjm1_logn(mpi),
149  ssjm1_lnq(mpi),
150  ssjm1_khi(mpi),
151  ssjm1_wbeta(mpi, CON, mpi.get_bvect_cart()),
152  ssjm1_h11(mpi),
153  ssjm1_h21(mpi),
154  ssjm1_h31(mpi),
155  ssjm1_h22(mpi),
156  ssjm1_h32(mpi),
157  ssjm1_h33(mpi),
158  decouple(mpi),
159  conf_flat(conf_flat0){
160 
161  // Pointers of derived quantities initialized to zero :
162  set_der_0x0() ;
163 
164  // Quantities defined on a spherical triad in star.C are put on
165  // a cartesian one
169  Metric temp_met (mp.flat_met_cart()) ;
170  gamma = temp_met ;
171 
172  // All quantities are initialized to zero :
173  psi0 = 0 ;
174  d_psi.set_etat_zero() ;
175  wit_w.set_etat_zero() ;
176  loggam = 0 ;
177  bsn.set_etat_zero() ;
178  pot_centri = 0 ;
179 
180  logn_auto = 0 ;
181  logn_comp = 0 ;
186  lnq_auto = 0 ;
187  lnq_comp = 0 ;
188  psi4 = 1 ;
191  hij.set_etat_zero() ;
194 
197  kcar_auto = 0 ;
198  kcar_comp = 0 ;
199  ssjm1_logn = 0 ;
200  ssjm1_lnq = 0 ;
201  ssjm1_khi = 0 ;
202  ssjm1_wbeta.set_etat_zero() ;
203  ssjm1_h11 = 0 ;
204  ssjm1_h21 = 0 ;
205  ssjm1_h31 = 0 ;
206  ssjm1_h22 = 0 ;
207  ssjm1_h32 = 0 ;
208  ssjm1_h33 = 0 ;
209 }
210 
211 // Copy constructor
212 // ----------------
214  : Star(star),
215  irrotational(star.irrotational),
216  psi0(star.psi0),
217  d_psi(star.d_psi),
218  wit_w(star.wit_w),
219  loggam(star.loggam),
220  bsn(star.bsn),
221  pot_centri(star.pot_centri),
222  logn_auto(star.logn_auto),
223  logn_comp(star.logn_comp),
224  dcov_logn(star.dcov_logn),
225  dcon_logn(star.dcon_logn),
226  lnq_auto(star.lnq_auto),
227  lnq_comp(star.lnq_comp),
228  psi4(star.psi4),
229  dcov_phi(star.dcov_phi),
230  dcon_phi(star.dcon_phi),
231  flat(star.flat),
232  gtilde(star.gtilde),
233  beta_auto(star.beta_auto),
234  beta_comp(star.beta_comp),
235  hij(star.hij),
236  hij_auto(star.hij_auto),
237  hij_comp(star.hij_comp),
238  tkij_auto(star.tkij_auto),
239  tkij_comp(star.tkij_comp),
240  kcar_auto(star.kcar_auto),
241  kcar_comp(star.kcar_comp),
242  ssjm1_logn(star.ssjm1_logn),
243  ssjm1_lnq(star.ssjm1_lnq),
244  ssjm1_khi(star.ssjm1_khi),
245  ssjm1_wbeta(star.ssjm1_wbeta),
246  ssjm1_h11(star.ssjm1_h11),
247  ssjm1_h21(star.ssjm1_h21),
248  ssjm1_h31(star.ssjm1_h31),
249  ssjm1_h22(star.ssjm1_h22),
250  ssjm1_h32(star.ssjm1_h32),
251  ssjm1_h33(star.ssjm1_h33),
252  decouple(star.decouple),
253  conf_flat(star.conf_flat)
254 {
255  set_der_0x0() ;
256 
257 }
258 
259 // Constructor from a file
260 // -----------------------
261 Star_bin::Star_bin(Map& mpi, const Eos& eos_i, FILE* fich)
262  : Star(mpi, eos_i, fich),
263  psi0(mpi),
264  d_psi(mpi, COV, mpi.get_bvect_cart()),
265  wit_w(mpi, CON, mpi.get_bvect_cart()),
266  loggam(mpi),
267  bsn(mpi, CON, mpi.get_bvect_cart()),
268  pot_centri(mpi),
269  logn_auto(mpi, *(mpi.get_mg()), fich),
270  logn_comp(mpi),
271  dcov_logn(mpi, COV, mpi.get_bvect_cart()),
272  dcon_logn(mpi, CON, mpi.get_bvect_cart()),
273  lnq_auto(mpi, *(mpi.get_mg()), fich),
274  lnq_comp(mpi),
275  psi4(mpi),
276  dcov_phi(mpi, COV, mpi.get_bvect_cart()),
277  dcon_phi(mpi, CON, mpi.get_bvect_cart()),
278  flat(mpi, mpi.get_bvect_cart()),
279  gtilde(flat),
280  beta_auto(mpi, mpi.get_bvect_cart(), fich),
281  beta_comp(mpi, CON, mpi.get_bvect_cart()),
282  hij(mpi, CON, mpi.get_bvect_cart()),
283  hij_auto(mpi, mpi.get_bvect_cart(), fich),
284  hij_comp(mpi, CON, mpi.get_bvect_cart()),
285  tkij_auto(mpi, CON, mpi.get_bvect_cart()),
286  tkij_comp(mpi, CON, mpi.get_bvect_cart()),
287  kcar_auto(mpi),
288  kcar_comp(mpi),
289  ssjm1_logn(mpi, *(mpi.get_mg()), fich),
290  ssjm1_lnq(mpi, *(mpi.get_mg()), fich),
291  ssjm1_khi(mpi, *(mpi.get_mg()), fich),
292  ssjm1_wbeta(mpi, mpi.get_bvect_cart(), fich),
293  ssjm1_h11(mpi, *(mpi.get_mg()), fich),
294  ssjm1_h21(mpi, *(mpi.get_mg()), fich),
295  ssjm1_h31(mpi, *(mpi.get_mg()), fich),
296  ssjm1_h22(mpi, *(mpi.get_mg()), fich),
297  ssjm1_h32(mpi, *(mpi.get_mg()), fich),
298  ssjm1_h33(mpi, *(mpi.get_mg()), fich),
299  decouple(mpi){
300 
301  // Etoile parameters
302  // -----------------
303 
304  // irrotational and conf_flat is read in the file:
305  fread(&irrotational, sizeof(bool), 1, fich) ;
306  fread(&conf_flat, sizeof(bool), 1, fich) ;
307 
308 
309  // Read of the saved fields:
310  // ------------------------
311 
312  if (irrotational) {
313  Scalar gam_euler_file(mp, *(mp.get_mg()), fich) ;
314  gam_euler = gam_euler_file ;
315 
316  Scalar psi0_file(mp, *(mp.get_mg()), fich) ;
317  psi0 = psi0_file ;
318  }
319 
320  // Quantities defined on a spherical triad in star.C are put on
321  // a cartesian one
325  Metric temp_met (mp.flat_met_cart()) ;
326  gamma = temp_met ;
327 
328  // All other fields are initialized to zero :
329  // ----------------------------------------
330 
331  d_psi.set_etat_zero() ;
332  wit_w.set_etat_zero() ;
333  loggam = 0 ;
334  bsn.set_etat_zero() ;
335  pot_centri = 0 ;
336  logn_comp = 0 ;
340  lnq_comp = 0 ;
341  psi4 = 1 ;
344  hij.set_etat_zero() ;
346 
349  kcar_auto = 0 ;
350  kcar_comp = 0 ;
351 
352  // Pointers of derived quantities initialized to zero
353  // --------------------------------------------------
354  set_der_0x0() ;
355 
356 }
357 
358  //------------//
359  // Destructor //
360  //------------//
361 
363 
364  del_deriv() ;
365 
366 }
367 
368  //----------------------------------//
369  // Management of derived quantities //
370  //----------------------------------//
371 
372 void Star_bin::del_deriv() const {
373 
374  Star::del_deriv() ;
375 
376  if (p_xa_barycenter != 0x0) delete p_xa_barycenter ;
377 
378  set_der_0x0() ;
379 }
380 
381 
382 
383 
384 void Star_bin::set_der_0x0() const {
385 
387 
388  p_xa_barycenter = 0x0 ;
389 
390 }
391 
393 
395 
396  del_deriv() ;
397 
398 }
399 
400 
401  //--------------//
402  // Assignment //
403  //--------------//
404 
405 // Assignment to another Star_bin
406 // --------------------------------
407 void Star_bin::operator=(const Star_bin& star) {
408 
409  // Assignment of quantities common to the derived classes of Star
410  Star::operator=(star) ;
411 
412  // Assignement of proper quantities of class Star_bin
413  irrotational = star.irrotational ;
414  psi0 = star.psi0 ;
415  d_psi = star.d_psi ;
416  wit_w = star.wit_w ;
417  loggam = star.loggam ;
418  bsn = star.bsn ;
419  pot_centri = star.pot_centri ;
420  logn_auto = star.logn_auto ;
421  logn_comp = star.logn_comp ;
422  dcov_logn = star.dcov_logn ;
423  dcon_logn = star.dcon_logn ;
424  lnq_auto = star.lnq_auto ;
425  lnq_comp = star.lnq_comp ;
426  psi4 = star.psi4 ;
427  dcov_phi = star.dcov_phi ;
428  dcon_phi = star.dcon_phi ;
429  flat = star.flat ;
430  gtilde = star.gtilde ;
431  beta_auto = star.beta_auto ;
432  beta_comp = star.beta_comp ;
433  hij = star.hij ;
434  hij_auto = star.hij_auto ;
435  hij_comp = star.hij_comp ;
436  tkij_auto = star.tkij_auto ;
437  tkij_comp = star.tkij_comp ;
438  kcar_auto = star.kcar_auto ;
439  kcar_comp = star.kcar_comp ;
440  ssjm1_logn = star.ssjm1_logn ;
441  ssjm1_lnq = star.ssjm1_lnq ;
442  ssjm1_khi = star.ssjm1_khi ;
443  ssjm1_wbeta = star.ssjm1_wbeta ;
444  ssjm1_h11 = star.ssjm1_h11 ;
445  ssjm1_h21 = star.ssjm1_h21 ;
446  ssjm1_h31 = star.ssjm1_h31 ;
447  ssjm1_h22 = star.ssjm1_h22 ;
448  ssjm1_h32 = star.ssjm1_h32 ;
449  ssjm1_h33 = star.ssjm1_h33 ;
450  decouple = star.decouple ;
451  conf_flat = star.conf_flat ;
452 
453  del_deriv() ; // Deletes all derived quantities
454 
455 }
456 
458 
459  del_deriv() ; // sets to 0x0 all the derived quantities
460  return pot_centri ;
461 
462 }
463 
465 
466  del_deriv() ; // sets to 0x0 all the derived quantities
467  return logn_comp ;
468 
469 }
470 
472 
473  del_deriv() ; // sets to 0x0 all the derived quantities
474  return beta_auto ;
475 
476 }
477 
479 
480  del_deriv() ; // sets to 0x0 all the derived quantities
481  return beta ;
482 
483 }
484 
485 
486  //--------------//
487  // Outputs //
488  //--------------//
489 
490 // Save in a file
491 // --------------
492 void Star_bin::sauve(FILE* fich) const {
493 
494  Star::sauve(fich) ;
495 
496  logn_auto.sauve(fich) ;
497  lnq_auto.sauve(fich) ;
498  beta_auto.sauve(fich) ;
499  hij_auto.sauve(fich) ;
500 
501  ssjm1_logn.sauve(fich) ;
502  ssjm1_lnq.sauve(fich) ;
503  ssjm1_khi.sauve(fich) ;
504  ssjm1_wbeta.sauve(fich) ;
505  ssjm1_h11.sauve(fich) ;
506  ssjm1_h21.sauve(fich) ;
507  ssjm1_h31.sauve(fich) ;
508  ssjm1_h22.sauve(fich) ;
509  ssjm1_h32.sauve(fich) ;
510  ssjm1_h33.sauve(fich) ;
511 
512  fwrite(&irrotational, sizeof(bool), 1, fich) ;
513  fwrite(&conf_flat, sizeof(bool), 1, fich) ;
514 
515  if (irrotational) {
516  gam_euler.sauve(fich) ; // required to construct d_psi from psi0
517  psi0.sauve(fich) ;
518  }
519 
520 }
521 
522 // Printing
523 // --------
524 
525 ostream& Star_bin::operator>>(ostream& ost) const {
526 
527  using namespace Unites ;
528 
529  Star::operator>>(ost) ;
530 
531  ost << endl ;
532  ost << "Star in a binary system" << endl ;
533  ost << "-----------------------" << endl ;
534 
535  if (irrotational) {
536  ost << "irrotational configuration" << endl ;
537  }
538  else {
539  ost << "corotating configuration" << endl ;
540  }
541 
542  ost << "Absolute abscidia of the stellar center: " <<
543  mp.get_ori_x() / km << " km" << endl ;
544 
545  ost << "Absolute abscidia of the barycenter of the baryon density : " <<
546  xa_barycenter() / km << " km" << endl ;
547 
548  double r_0 = 0.5 * ( ray_eq() + ray_eq_pi() ) ;
549  double d_ns = fabs( mp.get_ori_x() ) + ray_eq_pi() - r_0 ;
550  double d_tilde = 2 * d_ns / r_0 ;
551 
552  ost << "d_tilde : " << d_tilde << endl ;
553 
554  ost << "Central value of gam_euler : "
555  << gam_euler.val_grid_point(0, 0, 0, 0) << endl ;
556 
557  ost << "Central u_euler (U^r, U^t, U^p) [c] : "
558  << u_euler(1).val_grid_point(0, 0, 0, 0) << " "
559  << u_euler(2).val_grid_point(0, 0, 0, 0) << " "
560  << u_euler(3).val_grid_point(0, 0, 0, 0) << endl ;
561 
562  if (irrotational) {
563  ost << "Central d_psi (r, t, p) [c] : "
564  << d_psi(1).val_grid_point(0, 0, 0, 0) << " "
565  << d_psi(2).val_grid_point(0, 0, 0, 0) << " "
566  << d_psi(3).val_grid_point(0, 0, 0, 0) << endl ;
567 
568  ost << "Central vel. / co-orb. (W^r, W^t, W^p) [c] : "
569  << wit_w(1).val_grid_point(0, 0, 0, 0) << " "
570  << wit_w(2).val_grid_point(0, 0, 0, 0) << " "
571  << wit_w(3).val_grid_point(0, 0, 0, 0) << endl ;
572 
573  ost << "Max vel. / co-orb. (W^r, W^t, W^p) [c] : "
574  << max(max(wit_w(1))) << " "
575  << max(max(wit_w(2))) << " "
576  << max(max(wit_w(3))) << endl ;
577 
578  ost << "Min vel. / co-orb. (W^r, W^t, W^p) [c] : "
579  << min(min(wit_w(1))) << " "
580  << min(min(wit_w(2))) << " "
581  << min(min(wit_w(3))) << endl ;
582 
583  double r_surf = mp.val_r(0,1.,M_PI/4,M_PI/4) ;
584 
585  ost << "Velocity at (r_surf,pi/4,pi/4) / co-orb. [c] : "
586  << wit_w(1).val_point(r_surf,M_PI/4,M_PI/4) << " "
587  << wit_w(2).val_point(r_surf,M_PI/4,M_PI/4) << " "
588  << wit_w(3).val_point(r_surf,M_PI/4,M_PI/4) << endl ;
589 
590  ost << "Central value of loggam : "
591  << loggam.val_grid_point(0, 0, 0, 0) << endl ;
592  }
593 
594 
595  ost << "Central value of log(N) auto, comp : "
596  << logn_auto.val_grid_point(0, 0, 0, 0) << " "
597  << logn_comp.val_grid_point(0, 0, 0, 0) << endl ;
598 
599  ost << "Central value of beta (N^r, N^t, N^p) [c] : "
600  << beta(1).val_grid_point(0, 0, 0, 0) << " "
601  << beta(2).val_grid_point(0, 0, 0, 0) << " "
602  << beta(3).val_grid_point(0, 0, 0, 0) << endl ;
603 
604  ost << " ... beta_auto part of it [c] : "
605  << beta_auto(1).val_grid_point(0, 0, 0, 0) << " "
606  << beta_auto(2).val_grid_point(0, 0, 0, 0) << " "
607  << beta_auto(3).val_grid_point(0, 0, 0, 0) << endl ;
608 
609  ost << endl << "Central value of (B^r, B^t, B^p)/N [c] : "
610  << bsn(1).val_grid_point(0, 0, 0, 0) << " "
611  << bsn(2).val_grid_point(0, 0, 0, 0) << " "
612  << bsn(3).val_grid_point(0, 0, 0, 0) << endl ;
613 
614 
615  ost << endl << "Central A^{ij} [c/km] : " << endl ;
616  ost << " A^{xx} auto, comp : "
617  << tkij_auto(1, 1).val_grid_point(0, 0, 0, 0) * km << " "
618  << tkij_comp(1, 1).val_grid_point(0, 0, 0, 0) * km << endl ;
619  ost << " A^{xy} auto, comp : "
620  << tkij_auto(1, 2).val_grid_point(0, 0, 0, 0) * km << " "
621  << tkij_comp(1, 2).val_grid_point(0, 0, 0, 0) * km << endl ;
622  ost << " A^{xz} auto, comp : "
623  << tkij_auto(1, 3).val_grid_point(0, 0, 0, 0) * km << " "
624  << tkij_comp(1, 3).val_grid_point(0, 0, 0, 0) * km << endl ;
625  ost << " A^{yy} auto, comp : "
626  << tkij_auto(2, 2).val_grid_point(0, 0, 0, 0) * km << " "
627  << tkij_comp(2, 2).val_grid_point(0, 0, 0, 0) * km << endl ;
628  ost << " A^{yz} auto, comp : "
629  << tkij_auto(2, 3).val_grid_point(0, 0, 0, 0) * km << " "
630  << tkij_comp(2, 3).val_grid_point(0, 0, 0, 0) * km << endl ;
631  ost << " A^{zz} auto, comp : "
632  << tkij_auto(3, 3).val_grid_point(0, 0, 0, 0) * km << " "
633  << tkij_comp(3, 3).val_grid_point(0, 0, 0, 0) * km << endl ;
634 
635  ost << endl << "Central A_{ij} A^{ij} [c^2/km^2] : " << endl ;
636  ost << " A_{ij} A^{ij} auto, comp : "
637  << kcar_auto.val_grid_point(0, 0, 0, 0) * km*km << " "
638  << kcar_comp.val_grid_point(0, 0, 0, 0) * km*km << endl ;
639 
640 
641  return ost ;
642 }
643 
644  //-------------------------//
645  // Computational routines //
646  //-------------------------//
647 
649 
650  if (!irrotational) {
652  return ;
653  }
654 
655  // Specific relativistic enthalpy ---> hhh
656  //----------------------------------
657 
658  Scalar hhh = exp(ent) ; // = 1 at the Newtonian limit
659 
660  // Computation of W^i = - h Gamma_n B^i/N
661  //----------------------------------------------
662 
663  Vector www = hhh * gam_euler * bsn * psi4 ;
664 
665  // Constant value of W^i at the center of the star
666  //-------------------------------------------------
667 
668  Vector v_orb(mp, COV, mp.get_bvect_cart()) ;
669 
670  for (int i=1; i<=3; i++) {
671  v_orb.set(i) = www(i).val_grid_point(0, 0, 0, 0) ;
672  }
673 
674  // Gradient of psi
675  //----------------
676 
677  Vector d_psi0 = psi0.derive_cov(flat) ;
678 
679  d_psi0.change_triad( mp.get_bvect_cart() ) ;
680  d_psi0.std_spectral_base() ;
681 
682  d_psi = d_psi0 + v_orb ;
683  for (int i=1; i<=3; i++) {
684  if (d_psi(i).get_etat() == ETATZERO)
685  d_psi.set(i).annule_hard() ;
686  }
687 
688  // C^1 continuation of d_psi outside the star
689  // (to ensure a smooth enthalpy field accross the stellar surface)
690  // ----------------------------------------------------------------
691 
692  int nzm1 = mp.get_mg()->get_nzone() - 1 ;
693  d_psi.annule(nzet, nzm1) ;
694  for (int i=1; i<=3; i++) {
695  Cmp d_psi_i (d_psi(i)) ;
696  d_psi_i.va.set_base( d_psi0(i).get_spectral_va().base ) ;
697  d_psi_i = raccord_c1(d_psi_i, nzet) ;
698  d_psi.set(i) = d_psi_i ;
699  }
700 
701 }
702 
703 
704 void Star_bin::relaxation(const Star_bin& star_jm1, double relax_ent,
705  double relax_met, int mer, int fmer_met) {
706 
707  double relax_ent_jm1 = 1. - relax_ent ;
708  double relax_met_jm1 = 1. - relax_met ;
709 
710  ent = relax_ent * ent + relax_ent_jm1 * star_jm1.ent ;
711 
712  if ( (mer != 0) && (mer % fmer_met == 0)) {
713 
714  logn_auto = relax_met * logn_auto + relax_met_jm1 * star_jm1.logn_auto ;
715  lnq_auto = relax_met * lnq_auto + relax_met_jm1 * star_jm1.lnq_auto ;
716 
717  beta_auto = relax_met * beta_auto
718  + relax_met_jm1 * star_jm1.beta_auto ;
719 
720  hij_auto = relax_met * hij_auto + relax_met_jm1 * star_jm1.hij_auto ;
721 
722  }
723 
724  del_deriv() ;
725 
726  equation_of_state() ;
727 
728 }
729 
730 void Star_bin::test_K_Hi() const {
731 
732  int nr = mp.get_mg()->get_nr(0) ;
733  int nt = mp.get_mg()->get_nt(0) ;
734  int np = mp.get_mg()->get_np(0) ;
735 
736  cout << "La jauge de Dirac est elle bien satisfaite ??" << endl ;
737  cout << "Vector Hi" << endl ;
738  for (int i=1; i<=3; i++)
739  cout << " Comp. " << i << " : " << norme((gtilde.con()
740  .divergence(flat))(i)/(nr*nt*np)) << endl ;
741 
742 
743  cout << "Pour comparaison valeur de D_i(g^1i)" << endl ;
744  for (int i=1; i<=3; i++)
745  cout << " i = " << i << " : " << norme((gtilde.con().derive_cov(flat))
746  (1, i, i)/(nr*nt*np)) << endl ;
747 
748 
749 
750 }
751 }
Sym_tensor hij_comp
Deviation of the inverse conformal metric from the inverse flat metric generated principally by the ...
Definition: star.h:594
Vector dcov_phi
Covariant derivative of the logarithm of the conformal factor.
Definition: star.h:555
bool irrotational
true for an irrotational star, false for a corotating one
Definition: star.h:491
Metric for tensor calculation.
Definition: metric.h:90
virtual void sauve(FILE *) const
Save in a binary file.
Definition: tensor_sym.C:375
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition: metric.C:293
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: star.C:319
void operator=(const Star_bin &)
Assignment to another Star_bin.
Definition: star_bin.C:407
Scalar psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case) ...
Definition: star.h:496
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
Vector wit_w
Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Definition: star.h:507
Map & mp
Mapping associated with the star.
Definition: star.h:180
Scalar ssjm1_lnq
Effective source at the previous step for the resolution of the Poisson equation for lnq_auto...
Definition: star.h:628
Sym_tensor hij_auto
Deviation of the inverse conformal metric from the inverse flat metric generated principally by the ...
Definition: star.h:588
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
Metric gamma
3-metric
Definition: star.h:235
Vector d_psi
Gradient of (in the irrotational case) (Spherical components with respect to the mapping of the star...
Definition: star.h:501
Scalar ssjm1_h31
Effective source at the previous step for the resolution of the Poisson equation for h20_auto...
Definition: star.h:651
Lorene prototypes.
Definition: app_hor.h:67
Scalar ssjm1_h11
Effective source at the previous step for the resolution of the Poisson equation for h00_auto...
Definition: star.h:641
Standard units of space, time and mass.
Equation of state base class.
Definition: eos.h:209
Sym_tensor stress_euler
Spatial part of the stress-energy tensor with respect to the Eulerian observer.
Definition: star.h:212
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:777
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
virtual void set_etat_nondef()
Sets the logical state of all components to ETATNONDEF (undefined state).
Definition: tensor.C:498
Scalar & set_logn_comp()
Read/write of the logarithm of the lapse generated principally by the companion.
Definition: star_bin.C:464
Base class for coordinate mappings.
Definition: map.h:682
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:780
Scalar ssjm1_h33
Effective source at the previous step for the resolution of the Poisson equation for h22_auto...
Definition: star.h:666
Scalar logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by...
Definition: star.h:532
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition: star.h:512
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition: star_bin.C:384
Scalar kcar_auto
Part of the scalar generated by beta_auto, i.e.
Definition: star.h:612
virtual void sauve(FILE *) const
Save in a binary file.
Definition: tensor.C:915
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Scalar ent
Log-enthalpy.
Definition: star.h:190
Base class for stars.
Definition: star.h:175
double * p_xa_barycenter
Absolute coordinate X of the barycenter of the baryon density.
Definition: star.h:687
Tbl min(const Cmp &)
Minimum values of a Cmp in each domain.
Definition: cmp_math.C:461
Vector beta
Shift vector.
Definition: star.h:228
Tensor field of valence 1.
Definition: vector.h:188
Scalar logn_auto
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by...
Definition: star.h:527
const Tensor_sym & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base )
Definition: valeur.C:813
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:334
Scalar ssjm1_h21
Effective source at the previous step for the resolution of the Poisson equation for h10_auto...
Definition: star.h:646
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:643
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: vector.C:322
void annule_hard()
Sets the Scalar to zero in a hard way.
Definition: scalar.C:386
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
void test_K_Hi() const
Test if the gauge conditions we impose are well satisfied.
Definition: star_bin.C:730
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: star.h:204
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
Scalar decouple
Function used to construct the part generated by the star from the total .
Definition: star.h:676
const Vector & divergence(const Metric &) const
Returns the divergence of this with respect to a Metric .
Definition: sym_tensor.C:352
Star_bin(Map &mp_i, int nzet_i, const Eos &eos_i, bool irrot, bool conf_flat)
Standard constructor.
Definition: star_bin.C:118
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition: cmp_math.C:484
Sym_tensor hij
Total deviation of the inverse conformal metric from the inverse flat metric.
Definition: star.h:581
Scalar pot_centri
Centrifugal potential.
Definition: star.h:521
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: star_bin.C:392
virtual ~Star_bin()
Destructor.
Definition: star_bin.C:362
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star) ...
Definition: star.h:562
Scalar lnq_auto
Scalar field generated principally by the star.
Definition: star.h:543
Scalar lnq_comp
Scalar field generated principally by the companion star.
Definition: star.h:548
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star_bin.C:372
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: star.h:207
Sym_tensor tkij_comp
Part of the extrinsic curvature tensor generated by beta_comp.
Definition: star.h:606
double ray_eq_pi() const
Coordinate radius at , [r_unit].
Definition: star_global.C:189
virtual void sauve(FILE *) const
Save in a file.
Definition: star_bin.C:492
void fait_d_psi()
Computes the gradient of the total velocity potential .
Definition: star_bin.C:648
Vector dcov_logn
Covariant derivative of the total logarithm of the lapse.
Definition: star.h:535
Class for stars in binary system.
Definition: star.h:483
Metric gtilde
Conformal metric .
Definition: star.h:565
Vector dcon_logn
Contravariant derivative of the total logarithm of the lapse.
Definition: star.h:538
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star.C:301
virtual void del_hydro_euler()
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer...
Definition: star.C:333
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:803
Scalar ssjm1_h22
Effective source at the previous step for the resolution of the Poisson equation for h11_auto...
Definition: star.h:656
Vector beta_comp
Part of the shift vector generated principally by the star (Spherical components with respect to the ...
Definition: star.h:575
Sym_tensor tkij_auto
Part of the extrinsic curvature tensor generated by beta_auto.
Definition: star.h:600
bool conf_flat
true if the 3-metric is conformally flat, false for a more general metric.
Definition: star.h:681
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
Scalar ssjm1_h32
Effective source at the previous step for the resolution of the Poisson equation for h21_auto...
Definition: star.h:661
Vector beta_auto
Part of the shift vector generated principally by the star (Spherical components with respect to the ...
Definition: star.h:570
Scalar kcar_comp
Part of the scalar generated by beta_auto and beta_comp, i.e.
Definition: star.h:618
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star_bin.C:525
virtual void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: star.C:420
Scalar & set_pot_centri()
Read/write the centrifugal potential.
Definition: star_bin.C:457
virtual void sauve(FILE *) const
Save in a file.
Definition: star.C:400
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:506
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:474
void relaxation(const Star_bin &star_prev, double relax_ent, double relax_met, int mer, int fmer_met)
Performs a relaxation on ent, logn_auto, lnq_auto, beta_auto and hij_auto.
Definition: star_bin.C:704
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for khi...
Definition: star.h:634
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
Vector & set_beta()
Read/write of .
Definition: star_bin.C:478
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: star.C:465
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:390
virtual void annule(int l_min, int l_max)
Sets the Tensor to zero in several domains.
Definition: tensor.C:680
Vector bsn
3-vector shift, divided by N, of the rotating coordinates, .
Definition: star.h:518
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:464
Scalar psi4
Conformal factor .
Definition: star.h:552
Vector dcon_phi
Contravariant derivative of the logarithm of the conformal factor.
Definition: star.h:557
Vector & set_beta_auto()
Read/write of .
Definition: star_bin.C:471
virtual double xa_barycenter() const
Absolute coordinate X of the barycenter of the baryon density,.
Scalar ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto...
Definition: star.h:623