32 #include "hot_star_rot_diff_cfc.h" 34 #include "graphique.h" 35 #include "utilitaires.h" 41 int ,
const Tbl& ent_limit,
43 const Tbl& control,
double mbar_wanted,
44 double aexp_mass,
Tbl& diff,
Param*){
52 char display_bold[]=
"x[1m" ; display_bold[0] = 27 ;
53 char display_normal[] =
"x[0m" ; display_normal[0] = 27 ;
65 assert( ( type_t == SYM) || (type_t == NONSYM) ) ;
67 int i_b = mg->
get_nr(l_b) - 1 ;
68 int j_b = (type_t == SYM ? mg->
get_nt(l_b) - 1 : mg->
get_nt(l_b)/2 ) ;
72 double ent_b = ent_limit(
nzet-1) ;
77 int mer_max = icontrol(0) ;
78 int mer_rot = icontrol(1) ;
79 int mer_change_omega = icontrol(2) ;
80 int mer_fix_omega = icontrol(3) ;
81 int mer_mass = icontrol(4) ;
82 int mer_triax = icontrol(5) ;
83 int delta_mer_kep = icontrol(6) ;
86 if (mer_change_omega < mer_rot) {
87 cout <<
"Hot_star_rot_diff_CFC::equilibrium: mer_change_omega < mer_rot !" << endl ;
88 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
89 cout <<
" mer_rot = " << mer_rot << endl ;
92 if (mer_fix_omega < mer_change_omega) {
93 cout <<
"Hot_star_rot_diff_CFC::equilibrium: mer_fix_omega < mer_change_omega !" 95 cout <<
" mer_fix_omega = " << mer_fix_omega << endl ;
96 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
102 bool change_ent = true ;
105 mer_mass =
abs(mer_mass) ;
108 double precis = control(0) ;
109 double omega_ini = control(1) ;
110 double relax = control(2) ;
111 double relax_prev = double(1) - relax ;
112 double ampli_triax = control(3) ;
119 double& diff_ent = diff.
set(0) ;
120 double& vit_triax = diff.
set(7) ;
130 double accrois_omega = (omega_c0 - omega_ini) /
131 double(mer_fix_omega - mer_change_omega) ;
152 ofstream fichconv(
"convergence.d") ;
153 fichconv <<
"# diff_ent GRV2 max_triax vit_triax" << endl ;
155 ofstream fichfreq(
"frequency.d") ;
156 fichfreq <<
"# f [Hz]" << endl ;
158 ofstream fichevol(
"evolution.d") ;
160 "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c" 164 double err_grv2 = 1 ;
165 double max_triax_prev = 0 ;
172 for(
int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
174 cout <<
"-----------------------------------------------" << endl ;
175 cout <<
"step: " << mer << endl ;
176 cout <<
"diff_ent = " << display_bold << diff_ent << display_normal
178 cout <<
"err_grv2 = " << err_grv2 << endl ;
184 if (mer >= mer_rot) {
186 if (mer < mer_change_omega) {
187 omega_c = omega_ini ;
190 if (mer <= mer_fix_omega) {
191 omega_c = omega_ini + accrois_omega *
192 (mer - mer_change_omega) ;
229 if (mer == mer_triax) {
231 if ( mg->
get_np(0) == 1 ) {
233 "Hot_star_rot_diff::equilibrium: np must be stricly greater than 1" 234 << endl <<
" to set a triaxial perturbation !" << endl ;
241 perturb = 1 + ampli_triax * sint*sint *
cos(2*phi) ;
242 ln_f_new = ln_f_new * perturb ;
255 double max_triax = 0 ;
257 if ( mg->
get_np(0) > 1 ) {
259 for (
int l=0; l<nz; l++) {
260 for (
int j=0; j<mg->
get_nt(l); j++) {
261 for (
int i=0; i<mg->
get_nr(l); i++) {
264 double xcos2p = (*(va_nuf.
c_cf))(l, 2, j, i) ;
267 double xsin2p = (*(va_nuf.
c_cf))(l, 3, j, i) ;
269 double xx =
sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
271 max_triax = ( xx > max_triax ) ? xx : max_triax ;
278 cout <<
"Triaxial part of nuf : " << max_triax << endl ;
293 if (mer > mer_fix_omega + delta_mer_kep) {
295 omega_c *= fact_omega ;
299 bool omega_trop_grand = false ;
306 bool superlum = true ;
318 for (
int l=
nzet+1; l<nz; l++) {
340 +
beta(3).val_grid_point(l_b, k_b, j_b, i_b) ;
345 if (F_e !=
double(0)) {
350 for (l=0; l<
nzet+1; l++) {
351 for (k=0; k<mg->
get_np(l); k++) {
352 for (j=0; j<mg->
get_nt(l); j++) {
353 for (i=0; i<mg->
get_nr(l); i++) {
355 if (om_lkji > om_max) {
369 if (lambda2 == -1.) {
373 lambda1 = (1 +
pow(F_max/(B2*omega_c), p))/(1 +
pow(F_max/(A2*omega_c), q+p)) ;
376 cout <<
"Warning: F_e=" << F_e <<
" < 0. Aborting." << endl ;
399 / ((lambda1-1)*
pow(F_e, p) - (lambda2-1)*
pow(F_max, p)), 1./(p+q) ) ;
401 / (lambda2*(lambda1-1)*
pow(F_e, p+q) - lambda1*(lambda2-1)*
pow(F_max, p+q)), 1./p) ;
403 cout << F_e <<
" " << F_max << endl ;
410 cout <<
"Parameters of Omega(F) : " <<
par_frot << endl ;
412 double omeg_min = lambda2 * omega_c ;
413 double omeg_max = lambda1 * omega_c ;
414 double precis1 = 1.e-14 ;
415 int nitermax1 = 1000 ;
425 double omeg_min = 0 ;
426 double omeg_max = omega_c ;
427 double precis1 = 1.e-14 ;
428 int nitermax1 = 100 ;
458 for (
int l=0; l<
nzet; l++) {
459 for (
int i=0; i<mg->
get_nr(l); i++) {
464 cout <<
"U > c for l, i : " << l <<
" " << i
465 <<
" U = " <<
sqrt(u2) << endl ;
470 cout <<
"**** VELOCITY OF LIGHT REACHED ****" << endl ;
471 omega_c /= fact_omega ;
472 cout <<
"New central rotation frequency : " 473 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
474 omega_trop_grand = true ;
498 Scalar prim_entropy = expmHT ;
502 prim_entropy = prim_entropy.
primr(
false) ;
508 double prim_entr_b = prim_entropy.
val_grid_point(l_b, k_b, j_b, i_b) ;
515 double mlngamma_c = 0 ;
521 double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
522 - prim_entr_c + prim_entr_b + ln_q_c - ln_q_b + primf_c - primf_b)
523 / ( ln_f_b - ln_f_c ) ;
524 alpha_r =
sqrt(alpha_r2) ;
526 cout <<
"alpha_r = " << alpha_r << endl ;
535 logn = alpha_r2 * ln_f_new + ln_q_new ;
542 ent = (ent_c + logn_c + mlngamma_c - prim_entr_c)
557 for (
int l=0; l<
nzet; l++) {
558 int imax = mg->
get_nr(l) - 1 ;
559 if (l == l_b) imax-- ;
560 for (
int i=0; i<imax; i++) {
563 cout <<
"ent < 0 for l, i : " << l <<
" " << i
570 cout <<
"**** KEPLERIAN VELOCITY REACHED ****" << endl ;
571 omega_c /= fact_omega ;
572 cout <<
"New central rotation frequency : " 573 << omega_c/(2.*M_PI) * f_unit <<
" Hz" << endl ;
574 omega_trop_grand = true ;
579 if ( omega_trop_grand ) {
581 fact_omega =
sqrt( fact_omega ) ;
582 cout <<
"**** New fact_omega : " << fact_omega << endl ;
602 logn = relax *
logn + relax_prev * logn_prev ;
604 psi = relax * psi_new + relax_prev * psi_prev ;
606 beta = relax *
beta + relax_prev * beta_prev ;
620 fichfreq <<
" " << omega_c / (2*M_PI) * f_unit ;
622 fichevol <<
" " << ent_c ;
628 if (mer > mer_mass) {
631 if (mbar_wanted > 0.) {
632 xx =
mass_b() / mbar_wanted - 1. ;
633 cout <<
"Discrep. baryon mass <-> wanted bar. mass : " << xx
637 xx =
mass_g() / fabs(mbar_wanted) - 1. ;
638 cout <<
"Discrep. grav. mass <-> wanted grav. mass : " << xx
641 double xprog = ( mer > 2*mer_mass) ? 1. :
642 double(mer-mer_mass)/double(mer_mass) ;
644 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
645 double fact =
pow(ax, aexp_mass) ;
646 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
647 xx <<
" " << ax <<
" " << fact << endl ;
653 if (mer%4 == 0) omega_c *= fact ;
664 diff_ent = diff_ent_tbl(0) ;
665 for (
int l=1; l<
nzet; l++) {
666 diff_ent += diff_ent_tbl(l) ;
670 fichconv <<
" " <<
log10( fabs(diff_ent) + 1.e-16 ) ;
671 fichconv <<
" " <<
log10( fabs(err_grv2) + 1.e-16 ) ;
672 fichconv <<
" " <<
log10( fabs(max_triax) + 1.e-16 ) ;
675 if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
676 vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
679 fichconv <<
" " << vit_triax ;
688 max_triax_prev = max_triax ;
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
void solve_shift(Vector &shift_new) const
Solution of the shift equation for rotating stars in CFC.
Cmp log(const Cmp &)
Neperian logarithm.
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Scalar temp
Temperature field (in MeV)
Cmp exp(const Cmp &)
Exponential.
Scalar psi
Conformal factor .
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.
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Standard units of space, time and mass.
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)
virtual double grv2() const
Error on the virial identity GRV2.
Tensor field of valence 0 (or component of a tensorial field).
Tbl par_frot
Parameters of the function .
void solve_logn_q(Scalar &ln_q_new) const
Solution of the quadratic part of the Poisson equation for the lapse for rotating stars in CFC...
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Scalar primr(bool null_infty=true) const
Computes the radial primitive which vanishes for .
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.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Values and coefficients of a (real-value) function.
Scalar entropy
Entropy per baryon field (in $k_B$)
Scalar nn
Lapse function N .
Tensor field of valence 1.
Cmp cos(const Cmp &)
Cosine.
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
void solve_logn_f(Scalar &ln_f_new) const
Solution of the `‘matter’' part of the Poisson equation for the lapse for rotating stars in CFC...
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 logn
Logarithm of the lapse N .
Coord phi
coordinate centered on the grid
double ray_eq() const
Coordinate radius at , [r_unit].
int nzet
Number of domains of *mp occupied by the star.
virtual double mass_b() const
Baryonic mass.
virtual void homothetie(double lambda)=0
Sets a new radial scale.
double omega
Rotation angular velocity ([f_unit] )
int get_nzone() const
Returns the number of domains.
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
virtual double mass_g() const
Gravitational mass.
Map & mp
Mapping associated with the star.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
void solve_psi(Scalar &psi_new)
Solution of the equations for the conformal factor for rotating stars in CFC.
Scalar psi4
Conformal factor .
void mult_rsint()
Multiplication by everywhere; dzpuis is not changed.
Cmp pow(const Cmp &, int)
Power .
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
Active physical coordinates and mapping derivatives.
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.
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
void fait_omega_field(double omeg_min, double omeg_max, double precis, int nitermax)
Computes (member omega_field ).
double omega_max
Maximum value of .
Cmp log10(const Cmp &)
Basis 10 logarithm.
Cmp abs(const Cmp &)
Absolute value.
int get_taille() const
Gives the total size (ie dim.taille)
void update_metric()
Computes metric quantities from known potentials.
const Scalar & dsdr() const
Returns of *this .
Scalar omega_field
Field .
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Scalar & set(int)
Read/write access to a component.
double ray_pole() const
Coordinate radius at [r_unit].
void update_entropy()
Updates the part of the entropy from entropy_param.
void annule_hard()
Sets the Tbl to zero in a hard way.
double omega_min
Minimum value of .
const Valeur & get_spectral_va() const
Returns va (read only version)
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.