LORENE
et_bin_nsbh_equilibrium.C
1 /*
2  * Method of class Etoile to compute a static spherical configuration
3  * of a neutron star in a NS-BH binary system.
4  *
5  * (see file etoile.h for documentation).
6  *
7  */
8 
9 /*
10  * Copyright (c) 2003 Keisuke Taniguchi
11  *
12  * This file is part of LORENE.
13  *
14  * LORENE is free software; you can redistribute it and/or modify
15  * it under the terms of the GNU General Public License version 2
16  * as published by the Free Software Foundation.
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 char et_bin_nsbh_equilibrium_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_equilibrium.C,v 1.13 2014/10/13 08:52:56 j_novak Exp $" ;
30 
31 /*
32  * $Id: et_bin_nsbh_equilibrium.C,v 1.13 2014/10/13 08:52:56 j_novak Exp $
33  * $Log: et_bin_nsbh_equilibrium.C,v $
34  * Revision 1.13 2014/10/13 08:52:56 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.12 2014/10/06 15:13:08 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.11 2008/09/26 08:38:45 p_grandclement
41  * get rid of desaliasing
42  *
43  * Revision 1.10 2006/09/05 13:39:45 p_grandclement
44  * update of the bin_ns_bh project
45  *
46  * Revision 1.9 2006/06/01 12:47:53 p_grandclement
47  * update of the Bin_ns_bh project
48  *
49  * Revision 1.8 2006/04/25 07:21:58 p_grandclement
50  * Various changes for the NS_BH project
51  *
52  * Revision 1.7 2006/03/30 07:33:47 p_grandclement
53  * *** empty log message ***
54  *
55  * Revision 1.6 2005/10/18 13:12:33 p_grandclement
56  * update of the mixted binary codes
57  *
58  * Revision 1.5 2005/08/29 15:10:17 p_grandclement
59  * Addition of things needed :
60  * 1) For BBH with different masses
61  * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
62  * WORKING YET !!!
63  *
64  * Revision 1.4 2004/03/25 10:29:04 j_novak
65  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
66  *
67  * Revision 1.3 2003/10/24 12:34:06 k_taniguchi
68  * Change the notation as it should be
69  *
70  * Revision 1.2 2003/10/21 11:49:33 k_taniguchi
71  * Change the class from Etoile_bin to sub-class Et_bin_nsbh.
72  *
73  * Revision 1.1 2003/10/20 15:01:55 k_taniguchi
74  * Computation of an equilibrium configuration of a neutron star
75  * in a NS-BH binary system.
76  *
77  *
78  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_equilibrium.C,v 1.13 2014/10/13 08:52:56 j_novak Exp $
79  *
80  */
81 
82 // Headers C
83 #include <cmath>
84 
85 // Headers Lorene
86 #include "etoile.h"
87 #include "map.h"
88 #include "nbr_spx.h"
89 #include "et_bin_nsbh.h"
90 #include "param.h"
91 
92 #include "graphique.h"
93 #include "utilitaires.h"
94 #include "unites.h"
95 
96 namespace Lorene {
97 void Et_bin_nsbh::equilibrium_nsbh(bool adapt, double ent_c, int& niter, int mermax,
98  int mermax_poisson, double relax_poisson,
99  int mermax_potvit, double relax_potvit,
100  Tbl& diff) {
101 
102  // Fundamental constants and units
103  // -------------------------------
104  using namespace Unites ;
105 
106  // Initializations
107  // --------------
108 
109  const Mg3d* mg = mp.get_mg() ;
110  int nz = mg->get_nzone() ; // total number of domains
111 
112  // The following is required to initialize mp_prev as a Map_et:
113  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
114 
115  // Error indicators
116  // ----------------
117  double& diff_ent = diff.set(0) ;
118  double& diff_vel_pot = diff.set(1) ;
119  double& diff_lapse = diff.set(2) ;
120  double& diff_confpsi = diff.set(3) ;
121  double& diff_shift_x = diff.set(4) ;
122  double& diff_shift_y = diff.set(5) ;
123  double& diff_shift_z = diff.set(6) ;
124 
125 
126  // Parameters for the function Map_et::poisson for n_auto
127  // ------------------------------------------------------
128  double precis_poisson = 1.e-16 ;
129 
130  Param par_poisson1 ;
131  par_poisson1.add_int(mermax_poisson, 0) ; // maximum number of iterations
132  par_poisson1.add_double(relax_poisson, 0) ; // relaxation parameter
133  par_poisson1.add_double(precis_poisson, 1) ; // required precision
134  par_poisson1.add_int_mod(niter, 0) ; // number of iterations actually used
135  par_poisson1.add_cmp_mod( ssjm1_lapse ) ;
136 
137  // Parameters for the function Map_et::poisson for confpsi_auto
138  // ------------------------------------------------------------
139 
140  Param par_poisson2 ;
141  par_poisson2.add_int(mermax_poisson, 0) ; // maximum number of iterations
142  par_poisson2.add_double(relax_poisson, 0) ; // relaxation parameter
143  par_poisson2.add_double(precis_poisson, 1) ; // required precision
144  par_poisson2.add_int_mod(niter, 0) ; // number of iterations actually used
145  par_poisson2.add_cmp_mod( ssjm1_confpsi ) ;
146 
147  // Parameters for the function Tenseur::poisson_vect
148  // -------------------------------------------------
149 
150  Param par_poisson_vect ;
151  par_poisson_vect.add_int(mermax_poisson, 0) ;
152  // maximum number of iterations
153  par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
154  par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
155  par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
156  par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
157  par_poisson_vect.add_int_mod(niter, 0) ;
158 
159  // Parameters for the adaptation
160  Param par_adapt ;
161  int nitermax = 100 ;
162  int niter_adapt ;
163  int adapt_flag = (adapt) ? 1 : 0 ;
164  int nz_search = nzet + 1 ;
165  double precis_secant = 1.e-14 ;
166  double alpha_r ;
167  double reg_map = 1. ;
168  int k_b ;
169  int j_b ;
170  Tbl ent_limit(nzet) ;
171 
172  par_adapt.add_int(nitermax, 0) ;
173  par_adapt.add_int(nzet, 1) ;
174  par_adapt.add_int(nz_search, 2) ;
175  par_adapt.add_int(adapt_flag, 3) ;
176  par_adapt.add_int(j_b, 4) ;
177  par_adapt.add_int(k_b, 5) ;
178  par_adapt.add_int_mod(niter_adapt, 0) ;
179  par_adapt.add_double(precis_secant, 0) ;
180  par_adapt.add_double(reg_map, 1) ;
181  par_adapt.add_double(alpha_r, 2) ;
182  par_adapt.add_tbl(ent_limit, 0) ;
183 
184  // External potential
185  // See Eq (99) from Gourgoulhon et al. (2001)
186  // -----------------------------------------
187 
188 
189  Tenseur ent_jm1 = ent ; // Enthalpy at previous step
190  Tenseur source(mp) ; // source term in the equation for logn_auto
191  // and beta_auto
192  Tenseur source_shift(mp, 1, CON, ref_triad) ; // source term in the
193  // equation for shift_auto
194 
195  //=========================================================================
196  // Start of iteration
197  //=========================================================================
198 
199  for(int mer=0 ; mer<mermax ; mer++ ) {
200 
201  cout << "-----------------------------------------------" << endl ;
202  cout << "step: " << mer << endl ;
203  cout << "diff_ent = " << diff_ent << endl ;
204  //-----------------------------------------------------
205  // Resolution of the elliptic equation for the velocity
206  // scalar potential
207  //-----------------------------------------------------
208 
209  if (irrotational) {
210  diff_vel_pot = velocity_potential(mermax_potvit, precis_poisson,
211  relax_potvit) ;
212  }
213 
214  // Equation de la surface
215  //if (adapt) {
216 
217  // Rescaling of the radius : (Be carefull !)
218  int nt = mg->get_nt(nzet-1) ;
219  int np = mg->get_np(nzet-1) ;
220  int nr = mg->get_nr(nzet-1) ;
221 
222  // valeurs au centre
223  double hc = exp(ent_c) ;
224  double gamma_c = exp(loggam())(0,0,0,0) ;
225  double gamma_0_c = exp(-pot_centri())(0,0,0,0) ;
226  double n_auto_c = n_auto()(0,0,0,0) ;
227  double n_comp_c = n_comp()(0,0,0,0) ;
228 
229  double alpha_square = 0 ;
230  double constante = 0;
231  for (int k=0; k<np; k++) {
232  for (int j=0; j<nt; j++) {
233 
234  // valeurs au bord
235  double gamma_b = exp(loggam())(nzet-1,k,j,nr-1) ;
236  double gamma_0_b = exp(-pot_centri())(nzet-1,k,j,nr-1) ;
237  double n_auto_b = n_auto()(nzet-1,k,j,nr-1) ;
238  double n_comp_b = n_comp()(nzet-1,k,j,nr-1) ;
239 
240  // Les solutions :
241  double alpha_square_courant = (gamma_0_c*gamma_b*n_comp_b - hc*gamma_c*gamma_0_b*n_comp_c) /
242  (hc*gamma_c*gamma_0_b*n_auto_c-gamma_0_c*gamma_b*n_auto_b) ;
243  double constante_courant = gamma_b*(n_comp_b+alpha_square_courant*n_auto_b)/gamma_0_b ;
244 
245  if (alpha_square_courant > alpha_square) {
246  alpha_square = alpha_square_courant ;
247  k_b = k ;
248  j_b = j ;
249  constante = constante_courant ;
250  }
251  }
252  }
253 
254  alpha_r = sqrt(alpha_square) ;
255  cout << "Adaptation : " << k_b << " " << j_b << " " << alpha_r << endl ;
256 
257  // Le potentiel :
258  Tenseur potentiel (constante*exp(-loggam-pot_centri)/(n_auto*alpha_square+n_comp)) ;
259  potentiel.set_std_base() ;
260  for (int l=nzet+1 ; l<nz ; l++)
261  potentiel.set().va.set(l) = 1 ;
262 
263  Map_et mp_prev = mp_et ;
264  ent = log(potentiel) ;
265  ent.set_std_base() ;
266  ent().va.smooth(nzet, (ent.set().va)) ;
267 
268  ent_limit.set_etat_qcq() ;
269  for (int l=0; l<nzet; l++) { // loop on domains inside the star
270  ent_limit.set(l) = ent()(l, k_b, j_b, nr-1) ;
271  }
272 
273  // On adapte :
274  mp.adapt(ent(), par_adapt, 4) ;
275  mp_prev.homothetie(alpha_r) ;
276 
277  for (int l=nzet ; l<nz-1 ; l++)
278  mp.resize(l, 1./alpha_r) ;
279  mp.reevaluate_symy (&mp_prev, nzet, ent.set()) ;
280 
281 
282 
283  // Equation of state
284  //----------------------------------------------------
285  equation_of_state() ; // computes new values for nbar (n), ener (e)
286  // and press (p) from the new ent (H)
287 
288  //---------------------------------------------------------
289  // Matter source terms in the gravitational field equations
290  //---------------------------------------------------------
291  hydro_euler() ; // computes new values for ener_euler (E),
292  // s_euler (S) and u_euler (U^i)
293 
294 
295  //-------------------------------------------------
296  // Relative change in enthalpy
297  //-------------------------------------------------
298 
299  Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
300  diff_ent = diff_ent_tbl(0) ;
301  for (int l=1; l<nzet; l++) {
302  diff_ent += diff_ent_tbl(l) ;
303  }
304  diff_ent /= nzet ;
305 
306  ent_jm1 = ent ;
307 
308 
309  //--------------------------------------------------------
310  // Poisson equation for n_auto
311  //--------------------------------------------------------
312 
313  // Source
314  // See Eq (50) from Gourgoulhon et al. (2001)
315  // ------------------------------------------
316 
317  Tenseur confpsi_q = pow(confpsi, 4.) ;
318  Tenseur confpsi_c = pow(confpsi, 5.) ;
319 
320  if (relativistic) {
321  Tenseur tmp = flat_scalar_prod(tkij_tot, tkij_auto) ;
322  Tenseur kk (mp) ;
323  kk = 0 ;
324  Tenseur tmp2(mp) ;
325  tmp2.set_etat_qcq() ;
326  for (int i=0 ; i<3 ; i++) {
327  tmp2.set() = tmp(i, i) ;
328  kk = kk + tmp2 ;
329  }
330 
331  source = qpig * nnn * confpsi_q * (ener_euler + s_euler)
332  + nnn * confpsi_q * kk
334  confpsi ;
335  }
336  else {
337  cout <<
338  "WARNING : Et_bin_nsbh is for the relativistic calculation"
339  << endl ;
340  abort() ;
341  }
342 
343  source.set_std_base() ;
344 
345  // Resolution of the Poisson equation
346  // ----------------------------------
347  Cmp n_auto_old (n_auto()) ;
348  source().poisson(par_poisson1, n_auto.set()) ;
349  n_auto.set() = n_auto() + 0.5 ;
350 
351  // Difference pas précédent
352  // -----------------------------------------------------
353 
354  Tbl tdiff_lapse = diffrel(n_auto(), n_auto_old) ;
355  cout <<
356  "Relative difference on n_auto : "
357  << endl ;
358  for (int l=0; l<nz; l++) {
359  cout << tdiff_lapse(l) << " " ;
360  }
361  cout << endl ;
362  diff_lapse = max(abs(tdiff_lapse)) ;
363 
364  if (relativistic) {
365 
366 
367  //--------------------------------------------------------
368  // Poisson equation for confpsi_auto
369  //--------------------------------------------------------
370 
371  // Source
372  // See Eq (51) from Gourgoulhon et al. (2001)
373  // ------------------------------------------
374 
375  Tenseur tmp = flat_scalar_prod(tkij_tot, tkij_auto) ;
376  Tenseur kk (mp) ;
377  kk = 0 ;
378  Tenseur tmp2(mp) ;
379  tmp2.set_etat_qcq() ;
380  for (int i=0 ; i<3 ; i++) {
381  tmp2.set() = tmp(i, i) ;
382  kk = kk + tmp2 ;
383  }
384 
385  source = -0.5 * qpig * confpsi_c * ener_euler
386  - 0.125 * confpsi_c * kk ;
387 
388  source.set_std_base() ;
389 
390  // Resolution of the Poisson equation
391  // ----------------------------------
392  Cmp psi_old (confpsi_auto()) ;
393  source().poisson(par_poisson2, confpsi_auto.set()) ;
394  confpsi_auto.set() = confpsi_auto() + 0.5 ;
395 
396 
397  // Check: has the Poisson equation been correctly solved ?
398  // -----------------------------------------------------
399 
400  Tbl tdiff_confpsi = diffrel(confpsi_auto(), psi_old) ;
401  cout <<
402  "Relative difference on confpsi_auto : "
403  << endl ;
404  for (int l=0; l<nz; l++) {
405  cout << tdiff_confpsi(l) << " " ;
406  }
407  cout << endl ;
408  diff_confpsi = max(abs(tdiff_confpsi)) ;
409 
410  //--------------------------------------------------------
411  // Vector Poisson equation for shift_auto
412  //--------------------------------------------------------
413 
414  // Source
415  // See Eq (52) from Gourgoulhon et al. (2001)
416  // ------
417  Tenseur vtmp = d_n_auto -6. * nnn * d_confpsi_auto / confpsi ;
418  source_shift = 4.*qpig * nnn *confpsi_q *(ener_euler + press)
419  * u_euler ;
420  if (tkij_tot.get_etat() != ETATZERO)
421  source_shift = source_shift + 2.* flat_scalar_prod(tkij_tot, vtmp) ;
422  source_shift.set_std_base() ;
423  // Resolution of the Poisson equation
424  // ----------------------------------
425  // Filter for the source of shift vector
426  for (int i=0 ; i<3 ; i++)
427  if (source_shift(i).get_etat() != ETATZERO)
428  source_shift.set(i).va.coef_i() ;
429 
430 for (int i=0; i<3; i++)
431  if ((source_shift(i).get_etat() != ETATZERO) && (source_shift(i).va.c->t[nz-1]->get_etat() != ETATZERO))
432  source_shift.set(i).filtre(4) ;
433  for (int i=0; i<3; i++) {
434  if(source_shift(i).dz_nonzero()) {
435  assert( source_shift(i).get_dzpuis() == 4 ) ;
436  }
437  else{
438  (source_shift.set(i)).set_dzpuis(4) ;
439  }
440  }
441  //##
442  // source_shift.dec2_dzpuis() ; // dzpuis 4 -> 2
443 
444  double lambda_shift = double(1) / double(3) ;
445  // ON DOIT CHANGER DE TRIADE
446  source_shift.change_triad(mp.get_bvect_cart()) ;
447  Tenseur shift_old (shift_auto) ;
448  source_shift.poisson_vect_oohara(lambda_shift, par_poisson_vect,
451 
452  // Check: has the equation for shift_auto been correctly solved ?
453  // --------------------------------------------------------------
454 
455 
456 
457  Tbl tdiff_shift_x = diffrel(shift_auto(0), shift_old(0)) ;
458  Tbl tdiff_shift_y = diffrel(shift_auto(1), shift_old(1)) ;
459  Tbl tdiff_shift_z = diffrel(shift_auto(2), shift_old(2)) ;
460 
461  cout <<
462  "Relative difference on shift_auto : "
463  << endl ;
464  cout << "x component : " ;
465  for (int l=0; l<nz; l++) {
466  cout << tdiff_shift_x(l) << " " ;
467  }
468  cout << endl ;
469  cout << "y component : " ;
470  for (int l=0; l<nz; l++) {
471  cout << tdiff_shift_y(l) << " " ;
472  }
473  cout << endl ;
474  cout << "z component : " ;
475  for (int l=0; l<nz; l++) {
476  cout << tdiff_shift_z(l) << " " ;
477  }
478  cout << endl ;
479 
480  diff_shift_x = max(abs(tdiff_shift_x)) ;
481  diff_shift_y = max(abs(tdiff_shift_y)) ;
482  diff_shift_z = max(abs(tdiff_shift_z)) ;
483  } // End of relativistic equations
484 
485 
486  } // End of main loop
487 
488  //=========================================================================
489  // End of iteration
490  //=========================================================================
491 }
492 
493 // Truc pourri
494 void Et_bin_nsbh::equilibrium_nsbh (double, int, int, double,
495  int, double, double, const Tbl&, Tbl&) {
496 
497  cout << "Not implemented !" << endl ;
498  abort() ;
499 }
500 }
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:296
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur &#39;s...
Definition: etoile.h:828
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:270
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:220
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1170
Tenseur pot_centri
Centrifugal potential.
Definition: etoile.h:953
Lorene prototypes.
Definition: app_hor.h:64
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:765
Tenseur nnn
Total lapse function.
Definition: etoile.h:509
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:468
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
Tenseur press
Fluid pressure.
Definition: etoile.h:461
Tenseur n_auto
Part of the lapse { N} generated principaly by the star.
Definition: et_bin_nsbh.h:85
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:504
Tenseur confpsi
Total conformal factor $$.
Definition: et_bin_nsbh.h:101
Tenseur d_confpsi_comp
Gradient of { confpsi_comp} (Cartesian components with respect to { ref_triad})
Definition: et_bin_nsbh.h:119
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition: etoile.h:889
bool irrotational
true for an irrotational star, false for a corotating one
Definition: etoile.h:822
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:824
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:474
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:668
Map & mp
Mapping associated with the star.
Definition: etoile.h:429
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:448
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: etoile.C:566
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:435
int get_etat() const
Returns the logical state.
Definition: tenseur.h:704
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:918
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:348
virtual void resize(int l, double lambda)=0
Rescales the outer boundary of one domain.
Tenseur d_confpsi_auto
Gradient of { confpsi_auto} (Cartesian components with respect to { ref_triad})
Definition: et_bin_nsbh.h:114
Tenseur d_n_auto
Gradient of { n_auto} (Cartesian components with respect to { ref_triad})
Definition: et_bin_nsbh.h:93
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:432
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:437
Tenseur confpsi_auto
Part of the conformal factor $$ generated principaly by the star.
Definition: et_bin_nsbh.h:104
virtual void reevaluate_symy(const Map *mp_prev, int nzet, Cmp &uu) const =0
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
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:791
Cmp ssjm1_confpsi
Effective source at the previous step for the resolution of the Poisson equation for { confpsi_auto} ...
Definition: et_bin_nsbh.h:164
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition: etoile.h:849
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:457
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:410
Cmp ssjm1_lapse
Effective source at the previous step for the resolution of the Poisson equation for { n_auto} by mea...
Definition: et_bin_nsbh.h:158
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
Tenseur_sym tkij_tot
Total extrinsic curvature tensor $K^{ij}$ generated by { shift_auto} and { shift_comp}.
Definition: et_bin_nsbh.h:152
Basic array class.
Definition: tbl.h:161
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor $K^{ij}$ generated by { shift_auto}.
Definition: et_bin_nsbh.h:146
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:636
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:465
Tenseur n_comp
Part of the lapse { N} generated principaly by the companion star.
Definition: et_bin_nsbh.h:88
virtual void equilibrium_nsbh(double ent_c, int mermax, int mermax_poisson, double relax_poisson, int mermax_potvit, double relax_potvit, double thres_adapt, const Tbl &fact, Tbl &diff)
Computes an equilibrium configuration in a NS-BH binary system.
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:461
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: etoile.h:973
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: et_bin_hydro.C:106
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Definition: etoile.h:983