145 #include "graphique.h" 146 #include "utilitaires.h" 151 int nzadapt,
const Tbl& ent_limit,
const Itbl& icontrol,
152 const Tbl& control,
double mbar_wanted,
153 double aexp_mass,
Tbl& diff,
Param*) {
162 char display_bold[]=
"x[1m" ; display_bold[0] = 27 ;
163 char display_normal[] =
"x[0m" ; display_normal[0] = 27 ;
178 int i_b = mg->
get_nr(l_b) - 1 ;
179 int j_b = mg->
get_nt(l_b) - 1 ;
183 double ent_b = ent_limit(
nzet-1) ;
188 int mer_max = icontrol(0) ;
189 int mer_rot = icontrol(1) ;
190 int mer_change_omega = icontrol(2) ;
191 int mer_fix_omega = icontrol(3) ;
192 int mer_mass = icontrol(4) ;
193 int mermax_poisson = icontrol(5) ;
194 int mer_triax = icontrol(6) ;
195 int delta_mer_kep = icontrol(7) ;
198 if (mer_change_omega < mer_rot) {
199 cout <<
"Etoile_rot::equilibrium: mer_change_omega < mer_rot !" << endl ;
200 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
201 cout <<
" mer_rot = " << mer_rot << endl ;
204 if (mer_fix_omega < mer_change_omega) {
205 cout <<
"Etoile_rot::equilibrium: mer_fix_omega < mer_change_omega !" 207 cout <<
" mer_fix_omega = " << mer_fix_omega << endl ;
208 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
214 bool change_ent = true ;
217 mer_mass =
abs(mer_mass) ;
220 double precis = control(0) ;
221 double omega_ini = control(1) ;
222 double relax = control(2) ;
223 double relax_prev = double(1) - relax ;
224 double relax_poisson = control(3) ;
225 double thres_adapt = control(4) ;
226 double ampli_triax = control(5) ;
227 double precis_adapt = control(6) ;
234 double& diff_ent = diff.
set(0) ;
235 double& diff_nuf = diff.
set(1) ;
236 double& diff_nuq = diff.
set(2) ;
239 double& diff_shift_x = diff.
set(5) ;
240 double& diff_shift_y = diff.
set(6) ;
241 double& vit_triax = diff.
set(7) ;
252 int nz_search =
nzet + 1 ;
255 double reg_map = 1. ;
257 par_adapt.
add_int(nitermax, 0) ;
259 par_adapt.
add_int(nzadapt, 1) ;
262 par_adapt.
add_int(nz_search, 2) ;
264 par_adapt.
add_int(adapt_flag, 3) ;
283 par_adapt.
add_tbl(ent_limit, 0) ;
289 double precis_poisson = 1.e-16 ;
291 Param par_poisson_nuf ;
292 par_poisson_nuf.
add_int(mermax_poisson, 0) ;
293 par_poisson_nuf.
add_double(relax_poisson, 0) ;
294 par_poisson_nuf.
add_double(precis_poisson, 1) ;
298 Param par_poisson_nuq ;
299 par_poisson_nuq.
add_int(mermax_poisson, 0) ;
300 par_poisson_nuq.
add_double(relax_poisson, 0) ;
301 par_poisson_nuq.
add_double(precis_poisson, 1) ;
305 Param par_poisson_tggg ;
306 par_poisson_tggg.
add_int(mermax_poisson, 0) ;
307 par_poisson_tggg.
add_double(relax_poisson, 0) ;
308 par_poisson_tggg.
add_double(precis_poisson, 1) ;
314 Param par_poisson_dzeta ;
321 Param par_poisson_vect ;
323 par_poisson_vect.
add_int(mermax_poisson, 0) ;
324 par_poisson_vect.
add_double(relax_poisson, 0) ;
325 par_poisson_vect.
add_double(precis_poisson, 1) ;
337 double accrois_omega = (omega0 - omega_ini) /
338 double(mer_fix_omega - mer_change_omega) ;
388 ofstream fichconv(
"convergence.d") ;
389 fichconv <<
"# diff_ent GRV2 max_triax vit_triax" << endl ;
391 ofstream fichfreq(
"frequency.d") ;
392 fichfreq <<
"# f [Hz]" << endl ;
394 ofstream fichevol(
"evolution.d") ;
396 "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c" 400 double err_grv2 = 1 ;
401 double max_triax_prev = 0 ;
407 for(
int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
409 cout <<
"-----------------------------------------------" << endl ;
410 cout <<
"step: " << mer << endl ;
411 cout <<
"diff_ent = " << display_bold << diff_ent << display_normal
413 cout <<
"err_grv2 = " << err_grv2 << endl ;
418 if (mer >= mer_rot) {
420 if (mer < mer_change_omega) {
424 if (mer <= mer_fix_omega) {
425 omega = omega_ini + accrois_omega *
426 (mer - mer_change_omega) ;
448 source_nuf = qpig *
nbar ;
470 (source_tggg.
set()).mult_rsint() ;
491 for (
int i=0; i<3; i++) {
492 source_shift.set(i) += squad(i) ;
496 source_shift.set_std_base() ;
502 source_nuf().poisson(par_poisson_nuf,
nuf.
set()) ;
504 cout <<
"Test of the Poisson equation for nuf :" << endl ;
505 Tbl err = source_nuf().test_poisson(
nuf(), cout,
true) ;
506 diff_nuf = err(0, 0) ;
512 if (mer == mer_triax) {
514 if ( mg->
get_np(0) == 1 ) {
516 "Etoile_rot::equilibrium: np must be stricly greater than 1" 517 << endl <<
" to set a triaxial perturbation !" << endl ;
524 perturb = 1 + ampli_triax * sint*sint *
cos(2*phi) ;
538 double max_triax = 0 ;
540 if ( mg->
get_np(0) > 1 ) {
542 for (
int l=0; l<nz; l++) {
543 for (
int j=0; j<mg->
get_nt(l); j++) {
544 for (
int i=0; i<mg->
get_nr(l); i++) {
547 double xcos2p = (*(va_nuf.
c_cf))(l, 2, j, i) ;
550 double xsin2p = (*(va_nuf.
c_cf))(l, 3, j, i) ;
552 double xx =
sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
554 max_triax = ( xx > max_triax ) ? xx : max_triax ;
561 cout <<
"Triaxial part of nuf : " << max_triax << endl ;
569 source_nuq().poisson(par_poisson_nuq,
nuq.
set()) ;
571 cout <<
"Test of the Poisson equation for nuq :" << endl ;
572 err = source_nuq().test_poisson(
nuq(), cout,
true) ;
573 diff_nuq = err(0, 0) ;
580 if (source_shift.get_etat() != ETATZERO) {
582 for (
int i=0; i<3; i++) {
583 if(source_shift(i).dz_nonzero()) {
584 assert( source_shift(i).get_dzpuis() == 4 ) ;
587 (source_shift.set(i)).set_dzpuis(4) ;
595 double lambda_shift = double(1) / double(3) ;
597 if ( mg->
get_np(0) == 1 ) {
601 source_shift.poisson_vect(lambda_shift, par_poisson_vect,
604 cout <<
"Test of the Poisson equation for shift_x :" << endl ;
605 err = source_shift(0).test_poisson(
shift(0), cout,
true) ;
606 diff_shift_x = err(0, 0) ;
608 cout <<
"Test of the Poisson equation for shift_y :" << endl ;
609 err = source_shift(1).test_poisson(
shift(1), cout,
true) ;
610 diff_shift_y = err(0, 0) ;
629 if (mer > mer_fix_omega + delta_mer_kep) {
631 omega *= fact_omega ;
634 bool omega_trop_grand = false ;
641 bool superlum = true ;
657 ((
uuu.
set()).va).set_base( (tmp.
va).base ) ;
665 for (
int l=0; l<
nzet; l++) {
666 for (
int i=0; i<mg->
get_nr(l); i++) {
668 double u1 =
uuu()(l, 0, j_b, i) ;
671 cout <<
"U > c for l, i : " << l <<
" " << i
672 <<
" U = " << u1 << endl ;
677 cout <<
"**** VELOCITY OF LIGHT REACHED ****" << endl ;
678 omega /= fact_omega ;
679 cout <<
"New rotation frequency : " 680 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
681 omega_trop_grand = true ;
702 mlngamma = - 0.5 *
uuu*
uuu ;
706 double nuf_b =
nuf()(l_b, k_b, j_b, i_b) ;
707 double nuq_b =
nuq()(l_b, k_b, j_b, i_b) ;
708 double mlngamma_b = mlngamma()(l_b, k_b, j_b, i_b) ;
711 double nuf_c =
nuf()(0,0,0,0) ;
712 double nuq_c =
nuq()(0,0,0,0) ;
713 double mlngamma_c = 0 ;
717 double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
718 + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
719 alpha_r =
sqrt(alpha_r2) ;
720 cout <<
"alpha_r = " << alpha_r << endl ;
726 double nu_c =
logn()(0,0,0,0) ;
731 ent = (ent_c + nu_c + mlngamma_c) -
logn - mlngamma ;
738 for (
int l=0; l<
nzet; l++) {
739 int imax = mg->
get_nr(l) - 1 ;
740 if (l == l_b) imax-- ;
741 for (
int i=0; i<imax; i++) {
742 if (
ent()(l, 0, j_b, i) < 0. ) {
744 cout <<
"ent < 0 for l, i : " << l <<
" " << i
745 <<
" ent = " <<
ent()(l, 0, j_b, i) << endl ;
751 cout <<
"**** KEPLERIAN VELOCITY REACHED ****" << endl ;
752 omega /= fact_omega ;
753 cout <<
"New rotation frequency : " 754 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
755 omega_trop_grand = true ;
760 if ( omega_trop_grand ) {
762 fact_omega =
sqrt( fact_omega ) ;
763 cout <<
"**** New fact_omega : " << fact_omega << endl ;
778 double dent_eq =
ent().dsdr()(l_b, k_b, j_b, i_b) ;
779 double dent_pole =
ent().dsdr()(l_b, k_b, 0, i_b) ;
780 double rap_dent = fabs( dent_eq / dent_pole ) ;
781 cout <<
"| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
783 if ( rap_dent < thres_adapt ) {
785 cout <<
"******* FROZEN MAPPING *********" << endl ;
847 mp.
poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta,
850 err_grv2 = lbda_grv2 - 1;
851 cout <<
"GRV2: " << err_grv2 << endl ;
866 logn = relax *
logn + relax_prev * logn_prev ;
868 dzeta = relax *
dzeta + relax_prev * dzeta_prev ;
880 fichfreq <<
" " <<
omega / (2*M_PI) * f_unit ;
881 fichevol <<
" " << rap_dent ;
883 fichevol <<
" " << ent_c ;
889 if (mer > mer_mass) {
892 if (mbar_wanted > 0.) {
893 xx =
mass_b() / mbar_wanted - 1. ;
894 cout <<
"Discrep. baryon mass <-> wanted bar. mass : " << xx
898 xx =
mass_g() / fabs(mbar_wanted) - 1. ;
899 cout <<
"Discrep. grav. mass <-> wanted grav. mass : " << xx
902 double xprog = ( mer > 2*mer_mass) ? 1. :
903 double(mer-mer_mass)/double(mer_mass) ;
905 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
906 double fact =
pow(ax, aexp_mass) ;
907 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
908 xx <<
" " << ax <<
" " << fact << endl ;
914 if (mer%4 == 0)
omega *= fact ;
924 diff_ent = diff_ent_tbl(0) ;
925 for (
int l=1; l<
nzet; l++) {
926 diff_ent += diff_ent_tbl(l) ;
930 fichconv <<
" " <<
log10( fabs(diff_ent) + 1.e-16 ) ;
931 fichconv <<
" " <<
log10( fabs(err_grv2) + 1.e-16 ) ;
932 fichconv <<
" " <<
log10( fabs(max_triax) + 1.e-16 ) ;
935 if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
936 vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
939 fichconv <<
" " << vit_triax ;
948 max_triax_prev = max_triax ;
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Cmp log(const Cmp &)
Neperian logarithm.
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Radial mapping of rather general form.
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
void coef() const
Computes the coeffcients of *this.
Cmp sqrt(const Cmp &)
Square root.
void annule(int l)
Sets the Cmp to zero in a given domain.
void set_std_base()
Set the standard spectal basis of decomposition for each component.
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Standard units of space, time and mass.
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
double & set(int i)
Read/write of a particular element (index i) (1D case)
Tenseur nnn
Total lapse function.
Tenseur nphi
Metric coefficient .
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
double ray_eq() const
Coordinate radius at , [r_unit].
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.
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Basic integer array class.
Values and coefficients of a (real-value) function.
Tenseur press
Fluid pressure.
Cmp cos(const Cmp &)
Cosine.
Tenseur shift
Total shift vector.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
Coord phi
coordinate centered on the grid
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void update_metric()
Computes metric coefficients from known potentials.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
virtual void poisson2d(const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const =0
Computes the solution of a 2-D Poisson equation.
Tenseur nbar
Baryon density in the fluid frame.
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
virtual double grv2() const
Error on the virial identity GRV2.
virtual double mass_b() const
Baryon mass.
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
void mult_rsint()
Multiplication by .
Map & mp
Mapping associated with the star.
int get_nzone() const
Returns the number of domains.
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
virtual void homothetie(double lambda)
Sets a new radial scale.
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
int get_etat() const
Returns the logical state.
Tenseur bbb
Metric factor B.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Cmp pow(const Cmp &, int)
Power .
Active physical coordinates and mapping derivatives.
Tenseur uuu
Norm of u_euler.
Tenseur tggg
Metric potential .
double omega
Rotation angular velocity ([f_unit] )
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
int nzet
Number of domains of *mp occupied by the star.
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Tenseur a_car
Total conformal factor .
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
double ray_pole() const
Coordinate radius at [r_unit].
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Cmp log10(const Cmp &)
Basis 10 logarithm.
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
virtual void equilibrium(double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff, Param *=0x0)
Computes an equilibrium configuration.
Cmp abs(const Cmp &)
Absolute value.
virtual void reevaluate(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.
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Tenseur & logn
Metric potential = logn_auto.
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
virtual double mass_g() const
Gravitational mass.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Tenseur ener_euler
Total energy density in the Eulerian frame.
Tenseur & dzeta
Metric potential = beta_auto.
Valeur va
The numerical value of the Cmp.
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Tensor handling *** DEPRECATED : use class Tensor instead ***.
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.