87#include "utilitaires.h"
101 const Tbl& fact_resize,
102 const Tbl* pent_limit,
113 const Mg3d* mg =
mp.get_mg() ;
121 int i_b = mg->
get_nr(l_b) - 1 ;
131 double& diff_ent = diff.
set(0) ;
132 double& diff_vel_pot = diff.
set(1) ;
133 double& diff_psi = diff.
set(2) ;
134 double& diff_chi = diff.
set(3) ;
135 double& diff_beta_x = diff.
set(4) ;
136 double& diff_beta_y = diff.
set(5) ;
137 double& diff_beta_z = diff.
set(6) ;
151 int nz_search =
nzet ;
155 double precis_secant = 1.e-14 ;
157 double reg_map = 1. ;
160 par_adapt.
add_int(nitermax, 0) ;
166 par_adapt.
add_int(nz_search, 2) ;
169 par_adapt.
add_int(adapt_flag, 3) ;
192 if (pent_limit != 0x0) ent_limit = *pent_limit ;
194 par_adapt.
add_tbl(ent_limit, 0) ;
197 double precis_poisson = 1.e-16 ;
207 par_psi.
add_int(mermax_poisson, 0) ;
218 par_chi.
add_int(mermax_poisson, 0) ;
229 par_beta.
add_int(mermax_poisson, 0) ;
239 Tenseur ssjm1wbeta(
mp, 1, CON,
mp.get_bvect_cart()) ;
255 Scalar logn_auto =
log(chi_auto_p1) -
log(Psi_auto_p1) ;
276 Vector source_beta(
mp, CON,
mp.get_bvect_cart()) ;
284 for(
int mer=0 ; mer<mermax ; mer++ ) {
286 cout <<
"-----------------------------------------------" << endl ;
287 cout <<
"step: " << mer << endl ;
288 cout <<
"diff_ent = " << diff_ent << endl ;
318 double alpha_r2 = 0 ;
319 for (
int k=0; k<mg->
get_np(l_b); k++) {
320 for (
int j=0; j<mg->
get_nt(l_b); j++) {
325 double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b) /
326 ( logn_auto_b - logn_auto_c ) ;
330 if (alpha_r2_jk > alpha_r2) {
331 alpha_r2 = alpha_r2_jk ;
339 alpha_r =
sqrt(alpha_r2) ;
341 cout <<
"k_b, j_b, alpha_r: " << k_b <<
" " << j_b <<
" "
371 ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
373 cout <<
"pot" <<
norme(pot_ext) << endl ;
375 (
ent.set_spectral_va()).smooth(
nzet,
ent.set_spectral_va()) ;
384 double dent_eq =
ent.dsdr().val_point(
ray_eq(),M_PI/2.,0.) ;
385 double dent_pole =
ent.dsdr().val_point(
ray_pole(),0.,0.) ;
386 double rap_dent = fabs( dent_eq / dent_pole ) ;
387 cout <<
"| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
389 if ( rap_dent < thres_adapt ) {
391 cout <<
"******* FROZEN MAPPING *********" << endl ;
399 for (
int l=0; l<
nzet; l++) {
400 ent_limit.
set(l) =
ent.val_grid_point(l, k_b, j_b, i_b) ;
402 } ent_limit.
set(
nzet-1) = ent_b ;
407 mp.adapt(ent_cmp, par_adapt) ;
413 double rr_in_1 =
mp.val_r(
nzet, -1., M_PI/2., 0.) ;
416 double rr_out_nm2 =
mp.val_r(nz-2, 1., M_PI/2., 0.) ;
418 mp.resize(nz-2, rr_in_1/rr_out_nm2 * fact_resize(1)) ;
421 double rr_out_nm3 =
mp.val_r(nz-3, 1., M_PI/2., 0.) ;
423 mp.resize(nz-3, rr_in_1/rr_out_nm3 * fact_resize(0)) ;
428 double rr_out_nm3_new =
mp.val_r(nz-3, 1., M_PI/2., 0.) ;
430 for (
int i=
nzet-1; i<nz-4; i++) {
432 double rr_out_i =
mp.val_r(i, 1., M_PI/2., 0.) ;
434 double rr_mid = rr_out_i
435 + (rr_out_nm3_new - rr_out_i) /
double(nz - 3 - i) ;
437 double rr_2timesi = 2. * rr_out_i ;
439 if (rr_2timesi < rr_mid) {
441 double rr_out_ip1 =
mp.val_r(i+1, 1., M_PI/2., 0.) ;
443 mp.resize(i+1, rr_2timesi / rr_out_ip1) ;
447 double rr_out_ip1 =
mp.val_r(i+1, 1., M_PI/2., 0.) ;
449 mp.resize(i+1, rr_mid / rr_out_ip1) ;
494 mp.reevaluate_symy(&mp_prev,
nzet+1,
ent) ;
496 double ent_s_max = -1 ;
499 for (
int k=0; k<mg->
get_np(l_b); k++) {
500 for (
int j=0; j<mg->
get_nt(l_b); j++) {
501 double xx = fabs(
ent.val_grid_point(l_b, k, j, i_b) ) ;
502 if (xx > ent_s_max) {
509 cout <<
"Max. abs(enthalpy) at the boundary between domains nzet-1"
510 <<
" and nzet : " << endl ;
511 cout <<
" " << ent_s_max <<
" reached for k = " << k_s_max <<
512 " and j = " << j_s_max << endl ;
532 int nr =
mp.get_mg()->get_nr(0) ;
533 int nt =
mp.get_mg()->get_nt(0) ;
534 int np =
mp.get_mg()->get_np(0) ;
556 cout <<
"Resolution of the Poisson equation for Psi_auto:" << endl ;
560 cout <<
"Psi_auto: " << endl <<
norme(
Psi_auto/(nr*nt*np)) << endl ;
567 "Relative error in the resolution of the equation for Psi_auto : "
569 for (
int l=0; l<nz; l++) {
570 cout << tdiff_psi(l) <<
" " ;
573 <<
"==========================================================="
576 diff_psi =
max(
abs(tdiff_psi)) ;
592 cout <<
"Resolution of the Poisson equation for chi_auto:" << endl ;
596 cout <<
"chi_auto: " << endl <<
norme(
chi_auto/(nr*nt*np)) << endl ;
603 "Relative error in the resolution of the equation for chi_auto : "
606 for (
int l=0; l<nz; l++) {
607 cout << tdiff_chi(l) <<
" " ;
610 <<
"==========================================================="
613 diff_chi =
max(
abs(tdiff_chi)) ;
631 Tenseur source_p(
mp, 1, CON,
mp.get_bvect_cart() ) ;
633 for (
int i=0; i<3; i++) {
635 source_p.
set(i) =
Cmp(source_beta(i+1)) ;
640 Tenseur vect_auxi (
mp, 1, CON,
mp.get_bvect_cart()) ;
642 for (
int i=0; i<3 ;i++) vect_auxi.
set(i) =
Cmp(
w_beta(i+1)) ;
659 double lambda = double(1) / double(3) ;
660 source_p.
poisson_vect(lambda, par_beta, resu_p, vect_auxi, scal_auxi) ;
663 for (
int i=1; i<=3; i++) {
665 w_beta.set(i) = vect_auxi(i-1) ;
671 for (
int i=0; i<3; i++)
ssjm1_wbeta.set(i+1) = ssjm1wbeta(i) ;
687 Tbl tdiff_beta_x =
diffrel(lap_beta(1), source_beta(1)) ;
688 Tbl tdiff_beta_y =
diffrel(lap_beta(2), source_beta(2)) ;
689 Tbl tdiff_beta_z =
diffrel(lap_beta(3), source_beta(3)) ;
692 "Relative error in the resolution of the equation for beta_auto : "
694 cout <<
"x component : " ;
695 for (
int l=0; l<nz; l++) {
696 cout << tdiff_beta_x(l) <<
" " ;
699 cout <<
"y component : " ;
700 for (
int l=0; l<nz; l++) {
701 cout << tdiff_beta_y(l) <<
" " ;
704 cout <<
"z component : " ;
705 for (
int l=0; l<nz; l++) {
706 cout << tdiff_beta_z(l) <<
" " ;
710 diff_beta_x =
max(
abs(tdiff_beta_x)) ;
711 diff_beta_y =
max(
abs(tdiff_beta_y)) ;
712 diff_beta_z =
max(
abs(tdiff_beta_z)) ;
721 diff_ent = diff_ent_tbl(0) ;
722 for (
int l=1; l<
nzet; l++) diff_ent += diff_ent_tbl(l) ;
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Radial mapping of rather general form.
virtual void homothetie(double lambda)
Sets a new radial scale.
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
int get_nzone() const
Returns the number of domains.
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
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_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Tensor field of valence 0 (or component of a tensorial field).
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Scalar Psi
Total conformal factor .
Scalar pot_centri
Centrifugal potential.
Sym_tensor haij_auto
Part of the extrinsic curvature tensor generated by beta_auto.
Scalar psi4
Conformal factor .
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star) .
void equilibrium(double ent_c, int mermax, int mermax_potvit, int mermax_poisson, double relax_poisson, double relax_potvit, double thres_adapt, const Tbl &fact, const Tbl *pent_limit, Tbl &diff)
Computes an equilibrium configuration.
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for khi.
Scalar khi
Solution for the scalar part of the vector Poisson equation for .
Vector ssjm1_wbeta
Effective source at the previous step for wbeta of the vector Poisson equation for wbeta (needed for ...
Scalar Psi_comp
Scalar field generated principally by the companion star.
Scalar ssjm1_chi
Effective source at the previous step for the resolution of the Poisson equation for \chi_auto .
Vector dcov_Psi
Covariant derivative of the conformal factor .
Vector beta_auto
Part of the shift vector generated principally by the star (Spherical components with respect to the ...
Scalar Psi_auto
Scalar field generated principally by the star.
Scalar hacar_comp
Part of the scalar generated by beta_auto and beta_comp, i.e.
Scalar chi
Total function .
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Vector dcov_chi
Covariant derivative of the function .
Scalar chi_auto
Scalar field generated principally by the star.
Vector w_beta
Solution for the vector part of the vector Poisson equation for .
Scalar chi_comp
Scalar field generated principally by the companion star.
Scalar ssjm1_psi
Effective source at the previous step for the resolution of the Poisson equation for \Psi_auto .
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
bool irrotational
true for an irrotational star, false for a corotating one
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
Scalar hacar_auto
Part of the scalar generated by beta_auto, i.e.
Scalar nn
Lapse function N .
Scalar ener_euler
Total energy density in the Eulerian frame.
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
double ray_eq() const
Coordinate radius at , [r_unit].
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Scalar press
Fluid pressure.
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Map & mp
Mapping associated with the star.
int nzet
Number of domains of *mp occupied by the star.
double ray_pole() const
Coordinate radius at [r_unit].
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
double & set(int i)
Read/write of a particular element (index i) (1D case).
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void set_std_base()
Set the standard spectal basis of decomposition for each component.
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
Tensor field of valence 1.
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
Cmp sqrt(const Cmp &)
Square root.
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.
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Cmp pow(const Cmp &, int)
Power .
Cmp abs(const Cmp &)
Absolute value.
Cmp log(const Cmp &)
Neperian logarithm.
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Standard units of space, time and mass.