83#include "et_rot_diff.h"
87#include "utilitaires.h"
92 int nzadapt,
const Tbl& ent_limit,
94 const Tbl& control,
double mbar_wanted,
95 double aexp_mass,
Tbl& diff,
Param*) {
103 char display_bold[]=
"x[1m" ; display_bold[0] = 27 ;
104 char display_normal[] =
"x[0m" ; display_normal[0] = 27 ;
109 const Mg3d* mg =
mp.get_mg() ;
119 int i_b = mg->
get_nr(l_b) - 1 ;
120 int j_b = mg->
get_nt(l_b) - 1 ;
124 double ent_b = ent_limit(
nzet-1) ;
129 int mer_max = icontrol(0) ;
130 int mer_rot = icontrol(1) ;
131 int mer_change_omega = icontrol(2) ;
132 int mer_fix_omega = icontrol(3) ;
133 int mer_mass = icontrol(4) ;
134 int mermax_poisson = icontrol(5) ;
135 int mer_triax = icontrol(6) ;
136 int delta_mer_kep = icontrol(7) ;
139 if (mer_change_omega < mer_rot) {
140 cout <<
"Et_rot_diff::equilibrium: mer_change_omega < mer_rot !" << endl ;
141 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
142 cout <<
" mer_rot = " << mer_rot << endl ;
145 if (mer_fix_omega < mer_change_omega) {
146 cout <<
"Et_rot_diff::equilibrium: mer_fix_omega < mer_change_omega !"
148 cout <<
" mer_fix_omega = " << mer_fix_omega << endl ;
149 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
155 bool change_ent = true ;
158 mer_mass =
abs(mer_mass) ;
161 double precis = control(0) ;
162 double omega_ini = control(1) ;
163 double relax = control(2) ;
164 double relax_prev = double(1) - relax ;
165 double relax_poisson = control(3) ;
166 double thres_adapt = control(4) ;
167 double ampli_triax = control(5) ;
168 double precis_adapt = control(6) ;
175 double& diff_ent = diff.
set(0) ;
176 double& diff_nuf = diff.
set(1) ;
177 double& diff_nuq = diff.
set(2) ;
180 double& diff_shift_x = diff.
set(5) ;
181 double& diff_shift_y = diff.
set(6) ;
182 double& vit_triax = diff.
set(7) ;
193 int nz_search =
nzet + 1 ;
196 double reg_map = 1. ;
198 par_adapt.
add_int(nitermax, 0) ;
200 par_adapt.
add_int(nzadapt, 1) ;
203 par_adapt.
add_int(nz_search, 2) ;
205 par_adapt.
add_int(adapt_flag, 3) ;
224 par_adapt.
add_tbl(ent_limit, 0) ;
230 double precis_poisson = 1.e-16 ;
232 Param par_poisson_nuf ;
233 par_poisson_nuf.
add_int(mermax_poisson, 0) ;
234 par_poisson_nuf.
add_double(relax_poisson, 0) ;
235 par_poisson_nuf.
add_double(precis_poisson, 1) ;
239 Param par_poisson_nuq ;
240 par_poisson_nuq.
add_int(mermax_poisson, 0) ;
241 par_poisson_nuq.
add_double(relax_poisson, 0) ;
242 par_poisson_nuq.
add_double(precis_poisson, 1) ;
246 Param par_poisson_tggg ;
247 par_poisson_tggg.
add_int(mermax_poisson, 0) ;
248 par_poisson_tggg.
add_double(relax_poisson, 0) ;
249 par_poisson_tggg.
add_double(precis_poisson, 1) ;
255 Param par_poisson_dzeta ;
262 Param par_poisson_vect ;
264 par_poisson_vect.
add_int(mermax_poisson, 0) ;
265 par_poisson_vect.
add_double(relax_poisson, 0) ;
266 par_poisson_vect.
add_double(precis_poisson, 1) ;
278 double accrois_omega = (omega_c0 - omega_ini) /
279 double(mer_fix_omega - mer_change_omega) ;
301 Tenseur source_shift(
mp, 1, CON,
mp.get_bvect_cart()) ;
307 if (
nuf.get_etat() == ETATZERO) {
313 if (
nuq.get_etat() == ETATZERO) {
318 if (
tggg.get_etat() == ETATZERO) {
319 tggg.set_etat_qcq() ;
323 if (
dzeta.get_etat() == ETATZERO) {
324 dzeta.set_etat_qcq() ;
329 ofstream fichconv(
"convergence.d") ;
330 fichconv <<
"# diff_ent GRV2 max_triax vit_triax" << endl ;
332 ofstream fichfreq(
"frequency.d") ;
333 fichfreq <<
"# f [Hz]" << endl ;
335 ofstream fichevol(
"evolution.d") ;
337 "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
341 double err_grv2 = 1 ;
342 double max_triax_prev = 0 ;
348 for(
int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
350 cout <<
"-----------------------------------------------" << endl ;
351 cout <<
"step: " << mer << endl ;
352 cout <<
"diff_ent = " << display_bold << diff_ent << display_normal
354 cout <<
"err_grv2 = " << err_grv2 << endl ;
359 if (mer >= mer_rot) {
361 if (mer < mer_change_omega) {
362 omega_c = omega_ini ;
365 if (mer <= mer_fix_omega) {
366 omega_c = omega_ini + accrois_omega *
367 (mer - mer_change_omega) ;
389 source_nuf = qpig *
nbar ;
402 logn.gradient_spher() ) ;
411 (source_tggg.
set()).mult_rsint() ;
432 for (
int i=0; i<3; i++) {
433 source_shift.
set(i) += squad(i) ;
443 source_nuf().poisson(par_poisson_nuf,
nuf.set()) ;
445 cout <<
"Test of the Poisson equation for nuf :" << endl ;
446 Tbl err = source_nuf().test_poisson(
nuf(), cout,
true) ;
447 diff_nuf = err(0, 0) ;
453 if (mer == mer_triax) {
455 if ( mg->
get_np(0) == 1 ) {
457 "Et_rot_diff::equilibrium: np must be stricly greater than 1"
458 << endl <<
" to set a triaxial perturbation !" << endl ;
466 nuf.set() =
nuf() * perturb ;
479 double max_triax = 0 ;
481 if ( mg->
get_np(0) > 1 ) {
483 for (
int l=0; l<nz; l++) {
484 for (
int j=0; j<mg->
get_nt(l); j++) {
485 for (
int i=0; i<mg->
get_nr(l); i++) {
488 double xcos2p = (*(va_nuf.
c_cf))(l, 2, j, i) ;
491 double xsin2p = (*(va_nuf.
c_cf))(l, 3, j, i) ;
493 double xx =
sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
495 max_triax = ( xx > max_triax ) ? xx : max_triax ;
502 cout <<
"Triaxial part of nuf : " << max_triax << endl ;
510 source_nuq().poisson(par_poisson_nuq,
nuq.set()) ;
512 cout <<
"Test of the Poisson equation for nuq :" << endl ;
513 err = source_nuq().test_poisson(
nuq(), cout,
true) ;
514 diff_nuq = err(0, 0) ;
521 if (source_shift.
get_etat() != ETATZERO) {
523 for (
int i=0; i<3; i++) {
524 if(source_shift(i).dz_nonzero()) {
525 assert( source_shift(i).get_dzpuis() == 4 ) ;
528 (source_shift.
set(i)).set_dzpuis(4) ;
536 double lambda_shift = double(1) / double(3) ;
538 if ( mg->
get_np(0) == 1 ) {
542 source_shift.
poisson_vect(lambda_shift, par_poisson_vect,
545 cout <<
"Test of the Poisson equation for shift_x :" << endl ;
546 err = source_shift(0).test_poisson(
shift(0), cout,
true) ;
547 diff_shift_x = err(0, 0) ;
549 cout <<
"Test of the Poisson equation for shift_y :" << endl ;
550 err = source_shift(1).test_poisson(
shift(1), cout,
true) ;
551 diff_shift_y = err(0, 0) ;
565 if (mer > mer_fix_omega + delta_mer_kep) {
567 omega_c *= fact_omega ;
571 bool omega_trop_grand = false ;
578 bool superlum = true ;
592 double omeg_min = 0 ;
593 double omeg_max = omega_c ;
594 double precis1 = 1.e-14 ;
595 int nitermax1 = 100 ;
610 if (
uuu.get_etat() == ETATQCQ) {
612 ((
uuu.set()).va).set_base( (tmp.
va).base ) ;
620 for (
int l=0; l<
nzet; l++) {
621 for (
int i=0; i<mg->
get_nr(l); i++) {
623 double u1 =
uuu()(l, 0, j_b, i) ;
626 cout <<
"U > c for l, i : " << l <<
" " << i
627 <<
" U = " << u1 << endl ;
632 cout <<
"**** VELOCITY OF LIGHT REACHED ****" << endl ;
633 omega_c /= fact_omega ;
634 cout <<
"New central rotation frequency : "
635 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
636 omega_trop_grand = true ;
657 mlngamma = - 0.5 *
uuu*
uuu ;
661 double nuf_b =
nuf()(l_b, k_b, j_b, i_b) ;
662 double nuq_b =
nuq()(l_b, k_b, j_b, i_b) ;
663 double mlngamma_b = mlngamma()(l_b, k_b, j_b, i_b) ;
664 double primf_b =
prim_field()(l_b, k_b, j_b, i_b) ;
668 double nuf_c =
nuf()(0,0,0,0) ;
669 double nuq_c =
nuq()(0,0,0,0) ;
670 double mlngamma_c = 0 ;
675 double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
676 + nuq_c - nuq_b + primf_c - primf_b)
677 / ( nuf_b - nuf_c ) ;
678 alpha_r =
sqrt(alpha_r2) ;
679 cout <<
"alpha_r = " << alpha_r << endl ;
685 double nu_c =
logn()(0,0,0,0) ;
697 for (
int l=0; l<
nzet; l++) {
698 int imax = mg->
get_nr(l) - 1 ;
699 if (l == l_b) imax-- ;
700 for (
int i=0; i<imax; i++) {
701 if (
ent()(l, 0, j_b, i) < 0. ) {
703 cout <<
"ent < 0 for l, i : " << l <<
" " << i
704 <<
" ent = " <<
ent()(l, 0, j_b, i) << endl ;
710 cout <<
"**** KEPLERIAN VELOCITY REACHED ****" << endl ;
711 omega_c /= fact_omega ;
712 cout <<
"New central rotation frequency : "
713 << omega_c/(2.*M_PI) * f_unit <<
" Hz" << endl ;
714 omega_trop_grand = true ;
719 if ( omega_trop_grand ) {
721 fact_omega =
sqrt( fact_omega ) ;
722 cout <<
"**** New fact_omega : " << fact_omega << endl ;
733 double dent_eq =
ent().dsdr()(l_b, k_b, j_b, i_b) ;
734 double dent_pole =
ent().dsdr()(l_b, k_b, 0, i_b) ;
735 double rap_dent = fabs( dent_eq / dent_pole ) ;
736 cout <<
"| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
738 if ( rap_dent < thres_adapt ) {
740 cout <<
"******* FROZEN MAPPING *********" << endl ;
749 mp.adapt(
ent(), par_adapt) ;
757 mp.reevaluate(&mp_prev,
nzet+1,
ent.set()) ;
785 mp.poisson2d(source_tggg(),
mp.cmp_zero(), par_poisson_tggg,
792 mp.poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta,
795 err_grv2 = lbda_grv2 - 1;
796 cout <<
"GRV2: " << err_grv2 << endl ;
811 logn = relax *
logn + relax_prev * logn_prev ;
813 dzeta = relax *
dzeta + relax_prev * dzeta_prev ;
825 fichfreq <<
" " << omega_c / (2*M_PI) * f_unit ;
826 fichevol <<
" " << rap_dent ;
828 fichevol <<
" " << ent_c ;
834 if (mer > mer_mass) {
837 if (mbar_wanted > 0.) {
838 xx =
mass_b() / mbar_wanted - 1. ;
839 cout <<
"Discrep. baryon mass <-> wanted bar. mass : " << xx
843 xx =
mass_g() / fabs(mbar_wanted) - 1. ;
844 cout <<
"Discrep. grav. mass <-> wanted grav. mass : " << xx
847 double xprog = ( mer > 2*mer_mass) ? 1. :
848 double(mer-mer_mass)/double(mer_mass) ;
850 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
851 double fact =
pow(ax, aexp_mass) ;
852 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
853 xx <<
" " << ax <<
" " << fact << endl ;
859 if (mer%4 == 0) omega_c *= fact ;
869 diff_ent = diff_ent_tbl(0) ;
870 for (
int l=1; l<
nzet; l++) {
871 diff_ent += diff_ent_tbl(l) ;
875 fichconv <<
" " <<
log10( fabs(diff_ent) + 1.e-16 ) ;
876 fichconv <<
" " <<
log10( fabs(err_grv2) + 1.e-16 ) ;
877 fichconv <<
" " <<
log10( fabs(max_triax) + 1.e-16 ) ;
880 if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
881 vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
884 fichconv <<
" " << vit_triax ;
893 max_triax_prev = max_triax ;
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
void mult_rsint()
Multiplication by .
Valeur va
The numerical value of the Cmp.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
void annule(int l)
Sets the Cmp to zero in a given domain.
Active physical coordinates and mapping derivatives.
void fait_omega_field(double omeg_min, double omeg_max, double precis, int nitermax)
Computes (member omega_field ).
Tenseur prim_field
Field .
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Tenseur omega_field
Field .
Tbl par_frot
Parameters of the function .
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.
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Tenseur uuu
Norm of u_euler.
double omega
Rotation angular velocity ([f_unit] ).
Tenseur & logn
Metric potential = logn_auto.
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
virtual double mass_g() const
Gravitational mass.
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Tenseur tggg
Metric potential .
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Tenseur nphi
Metric coefficient .
virtual double mass_b() const
Baryon mass.
Tenseur bbb
Metric factor B.
void update_metric()
Computes metric coefficients from known potentials.
Tenseur & dzeta
Metric potential = beta_auto.
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
virtual double grv2() const
Error on the virial identity GRV2.
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
int nzet
Number of domains of *mp occupied by the star.
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Tenseur nbar
Baryon density in the fluid frame.
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Map & mp
Mapping associated with the star.
Tenseur press
Fluid pressure.
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Tenseur ener_euler
Total energy density in the Eulerian frame.
Tenseur shift
Total shift vector.
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Tenseur a_car
Total conformal factor .
double ray_pole() const
Coordinate radius at [r_unit].
Basic integer array class.
Radial mapping of rather general form.
virtual void homothetie(double lambda)
Sets a new radial scale.
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
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_double_mod(double &x, int position=0)
Adds the address of a new modifiable double 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.
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_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 poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
int get_etat() const
Returns the logical state.
Values and coefficients of a (real-value) function.
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
void coef() const
Computes the coeffcients of *this.
Cmp sqrt(const Cmp &)
Square root.
Cmp log10(const Cmp &)
Basis 10 logarithm.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Cmp pow(const Cmp &, int)
Power .
Cmp cos(const Cmp &)
Cosine.
Cmp abs(const Cmp &)
Absolute value.
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...
Coord phi
coordinate centered on the grid
Standard units of space, time and mass.