138#include "et_bin_nsbh.h"
139#include "utilitaires.h"
143Cmp raccord_c1(
const Cmp& uu,
int l1) ;
150 if (
eos.identify() == 5 ||
eos.identify() == 4 ||
151 eos.identify() == 3) {
155 int nzm1 =
mp.get_mg()->get_nzone() - 1 ;
181 for (
int i=0; i<3; i++) {
182 v_orb.
set(i) = www(i)(0, 0, 0, 0) ;
185 v_orb.
annule(nzm1, nzm1) ;
200 for (
int l=
nzet; l <= nzm1; l++) {
201 dndh_log.
set(l) = 1 ;
209 Scalar zeta_h_scalar (zeta_h()) ;
211 for (
int l=1; l<=nzm1; l++)
214 Cmp zeta_h_cmp (zeta_h_scalar) ;
215 zeta_h.
set() = zeta_h_cmp ;
249 for (
int l=1; l<=nzm1; l++)
265 if (
psi0.get_etat() == ETATZERO) {
266 psi0.set_etat_qcq() ;
270 int nr =
mp.get_mg()->get_nr(0);
271 int nt =
mp.get_mg()->get_nt(0);
272 int np =
mp.get_mg()->get_np(0);
274 cout <<
"nr = " << nr <<
" nt = " << nt <<
" np = " << np << endl ;
276 cout <<
"psi0" << endl <<
norme(
psi0()/(nr*nt*np)) << endl ;
277 cout <<
"d(psi)/dr" << endl <<
norme(
psi0.set().dsdr()/(nr*nt*np)) << endl ;
279 Valeur lim(
mp.get_mg()->get_angu()) ;
283 Tenseur normal2 (
mp, 1, COV,
mp.get_bvect_cart()) ;
295 for (
int i=0; i<3; i++) {
296 for (
int j=0; j<i; j++) {
303 Metrique flat(plat,
true) ;
307 for (
int i=0; i<3; i++) {
308 normal.
set(i) = dcov_r(i) ;
309 normal2.
set(i) = dcov_r(i) ;
328 limite = ( - dcov_psi(1) * normal(1) - dcov_psi(2) * normal(2)
332 for (
int j=0; j<nt; j++)
333 for (
int k=0; k<np; k++)
334 lim.
set(0, k, j, 0) = limite(0, k, j, nr-1) ;
340 source().poisson_neumann_interne(lim, par, resu) ;
354 for (
int l=1; l<=nzm1; l++)
355 psi0.set().annule(l) ;
362 Cmp laplacien_psi0 =
psi0().laplacien() ;
364 erreur =
diffrel(laplacien_psi0, source())(0) ;
366 cout <<
"Check of the resolution of the continuity equation for strange stars: "
368 cout <<
"norme(source) : " <<
norme(source())(0) << endl
369 <<
"Error in the solution : " << erreur << endl ;
378 d_psi.set_etat_qcq() ;
380 for (
int i=0; i<3; i++) {
381 d_psi.set(i) = (
psi0.gradient())(i) + v_orb(i) ;
391 for (
int i=0; i<3; i++) {
395 assert(
d_psi.get_triad() == &(
mp.get_bvect_cart()) ) ;
408 int nzm1 =
mp.get_mg()->get_nzone() - 1 ;
434 for (
int i=0; i<3; i++) {
435 v_orb.
set(i) = www(i)(0, 0, 0, 0) ;
449 for (
int l=0; l<
nzet; l++) {
454 dndh_log = dndh_log +
eos.der_nbar_ent(
ent(), 1, l, &par) ;
462 for (
int l=
nzet; l <= nzm1; l++) {
463 dndh_log.
set(l) = 1 ;
509 if (
psi0.get_etat() == ETATZERO) {
510 psi0.set_etat_qcq() ;
516 mp.poisson_compact(
nzet, source(), zeta_h(), bb, par,
psi0.set() ) ;
524 Cmp oper = zeta_h() *
psi0().laplacien() + bb_dpsi0() ;
528 cout <<
"Check of the resolution of the continuity equation : " << endl ;
531 for (
int l=0; l<
nzet; l++) {
532 double err = terr(l) ;
533 cout <<
" domain " << l <<
" : norme(source) : " <<
norme(source())(l)
534 <<
" relative error : " << err << endl ;
535 if (err > erreur) erreur = err ;
546 d_psi.set_etat_qcq() ;
548 for (
int i=0; i<3; i++) {
549 d_psi.set(i) = (
psi0.gradient())(i) + v_orb(i) ;
559 for (
int i=0; i<3; i++) {
563 assert(
d_psi.get_triad() == &(
mp.get_bvect_cart()) ) ;
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Valeur va
The numerical value of the Cmp.
void annule(int l)
Sets the Cmp to zero in a given domain.
Tbl & set(int l)
Read/write of the value in a given domain.
Active physical coordinates and mapping derivatives.
Class for a star in a NS-BH binary system.
Tenseur psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case).
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Tenseur d_psi
Gradient of (in the irrotational case) (Cartesian components with respect to ref_triad ).
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur 's...
Tenseur beta_comp
Part of the logarithm of AN generated principaly by the companion star.
Cmp ssjm1_psi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
int nzet
Number of domains of *mp occupied by the star.
Tenseur nnn
Total lapse function.
const Eos & eos
Equation of state of the stellar matter.
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Map & mp
Mapping associated with the star.
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Tenseur a_car
Total conformal factor .
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Tensor field of valence 0 (or component of a tensorial field).
Tbl & set_domain(int l)
Read/write of the value in a given domain.
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
Class intended to describe tensors with a symmetry on the two last indices *** DEPRECATED : use class...
Tensor handling *** DEPRECATED : use class Tensor instead ***.
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
const Tenseur & derive_cov(const Metrique &met) const
Returns the covariant derivative of *this , with respect to met .
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void annule(int l)
Sets the Tenseur to zero in a given domain.
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates).
void set_std_base()
Set the standard spectal basis of decomposition for each component.
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Values and coefficients of a (real-value) function.
void ylm()
Computes the coefficients of *this.
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
void ylm_i()
Inverse of ylm().
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
void annule_hard()
Sets the Valeur to zero in a hard way.
Cmp sqrt(const Cmp &)
Square root.
Cmp exp(const Cmp &)
Exponential.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Cmp log(const Cmp &)
Neperian logarithm.
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...
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .