72 #include "graphique.h" 73 #include "utilitaires.h" 78 int nzadapt,
const Tbl& ent_limit,
const Itbl& icontrol,
79 const Tbl& control,
double mbar_wanted,
80 double aexp_mass,
Tbl& diff,
Param*) {
89 char display_bold[]=
"x[1m" ; display_bold[0] = 27 ;
90 char display_normal[] =
"x[0m" ; display_normal[0] = 27 ;
105 int i_b = mg->
get_nr(l_b) - 1 ;
106 int j_b = mg->
get_nt(l_b) - 1 ;
110 double ent_b = ent_limit(
nzet-1) ;
115 int mer_max = icontrol(0) ;
116 int mer_rot = icontrol(1) ;
117 int mer_change_omega = icontrol(2) ;
118 int mer_fix_omega = icontrol(3) ;
119 int mer_mass = icontrol(4) ;
120 int mermax_poisson = icontrol(5) ;
121 int mer_triax = icontrol(6) ;
122 int delta_mer_kep = icontrol(7) ;
125 if (mer_change_omega < mer_rot) {
126 cout <<
"Star_rot::equilibrium: mer_change_omega < mer_rot !" << endl ;
127 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
128 cout <<
" mer_rot = " << mer_rot << endl ;
131 if (mer_fix_omega < mer_change_omega) {
132 cout <<
"Star_rot::equilibrium: mer_fix_omega < mer_change_omega !" 134 cout <<
" mer_fix_omega = " << mer_fix_omega << endl ;
135 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
141 bool change_ent = true ;
144 mer_mass =
abs(mer_mass) ;
147 double precis = control(0) ;
148 double omega_ini = control(1) ;
149 double relax = control(2) ;
150 double relax_prev = double(1) - relax ;
151 double relax_poisson = control(3) ;
152 double thres_adapt = control(4) ;
153 double ampli_triax = control(5) ;
154 double precis_adapt = control(6) ;
161 double& diff_ent = diff.
set(0) ;
162 double& diff_nuf = diff.
set(1) ;
163 double& diff_nuq = diff.
set(2) ;
166 double& diff_shift_x = diff.
set(5) ;
167 double& diff_shift_y = diff.
set(6) ;
168 double& vit_triax = diff.
set(7) ;
179 int nz_search =
nzet + 1 ;
182 double reg_map = 1. ;
184 par_adapt.
add_int(nitermax, 0) ;
186 par_adapt.
add_int(nzadapt, 1) ;
189 par_adapt.
add_int(nz_search, 2) ;
191 par_adapt.
add_int(adapt_flag, 3) ;
210 par_adapt.
add_tbl(ent_limit, 0) ;
216 double precis_poisson = 1.e-16 ;
224 for (
int i=1; i<=3; i++) {
230 for (
int i=1; i<=3; i++) {
231 cshift.set(i-1) = -
beta(i) ;
236 for (
int i=1; i<=3; i++) {
237 cw_shift.set(i-1) =
w_shift(i) ;
244 Param par_poisson_nuf ;
245 par_poisson_nuf.
add_int(mermax_poisson, 0) ;
246 par_poisson_nuf.
add_double(relax_poisson, 0) ;
247 par_poisson_nuf.
add_double(precis_poisson, 1) ;
251 Param par_poisson_nuq ;
252 par_poisson_nuq.
add_int(mermax_poisson, 0) ;
253 par_poisson_nuq.
add_double(relax_poisson, 0) ;
254 par_poisson_nuq.
add_double(precis_poisson, 1) ;
258 Param par_poisson_tggg ;
259 par_poisson_tggg.
add_int(mermax_poisson, 0) ;
260 par_poisson_tggg.
add_double(relax_poisson, 0) ;
261 par_poisson_tggg.
add_double(precis_poisson, 1) ;
267 Param par_poisson_dzeta ;
274 Param par_poisson_vect ;
276 par_poisson_vect.
add_int(mermax_poisson, 0) ;
277 par_poisson_vect.
add_double(relax_poisson, 0) ;
278 par_poisson_vect.
add_double(precis_poisson, 1) ;
290 double accrois_omega = (omega0 - omega_ini) /
291 double(mer_fix_omega - mer_change_omega) ;
318 ofstream fichconv(
"convergence.d") ;
319 fichconv <<
"# diff_ent GRV2 max_triax vit_triax" << endl ;
321 ofstream fichfreq(
"frequency.d") ;
322 fichfreq <<
"# f [Hz]" << endl ;
324 ofstream fichevol(
"evolution.d") ;
326 "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c" 330 double err_grv2 = 1 ;
331 double max_triax_prev = 0 ;
337 for(
int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
339 cout <<
"-----------------------------------------------" << endl ;
340 cout <<
"step: " << mer << endl ;
341 cout <<
"diff_ent = " << display_bold << diff_ent << display_normal
343 cout <<
"err_grv2 = " << err_grv2 << endl ;
348 if (mer >= mer_rot) {
350 if (mer < mer_change_omega) {
354 if (mer <= mer_fix_omega) {
355 omega = omega_ini + accrois_omega *
356 (mer - mer_change_omega) ;
377 source_nuq =
ak_car - d_logn(1)*(d_logn(1)+d_bet(1))
378 - d_logn(2)*(d_logn(2)+d_bet(2))
379 - d_logn(3)*(d_logn(3)+d_bet(3)) ;
382 source_nuf = qpig *
nbar ;
395 - d_logn(1)*d_logn(1) - d_logn(2)*d_logn(2) - d_logn(3)*d_logn(3) ;
429 cout <<
"Test of the Poisson equation for nuf :" << endl ;
431 diff_nuf = err(0, 0) ;
437 if (mer == mer_triax) {
439 if ( mg->
get_np(0) == 1 ) {
441 "Star_rot::equilibrium: np must be stricly greater than 1" 442 << endl <<
" to set a triaxial perturbation !" << endl ;
449 perturb = 1 + ampli_triax * sint*sint *
cos(2*phi) ;
463 double max_triax = 0 ;
465 if ( mg->
get_np(0) > 1 ) {
467 for (
int l=0; l<nz; l++) {
468 for (
int j=0; j<mg->
get_nt(l); j++) {
469 for (
int i=0; i<mg->
get_nr(l); i++) {
472 double xcos2p = (*(va_nuf.
c_cf))(l, 2, j, i) ;
475 double xsin2p = (*(va_nuf.
c_cf))(l, 3, j, i) ;
477 double xx =
sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
479 max_triax = ( xx > max_triax ) ? xx : max_triax ;
486 cout <<
"Triaxial part of nuf : " << max_triax << endl ;
496 cout <<
"Test of the Poisson equation for nuq :" << endl ;
498 diff_nuq = err(0, 0) ;
505 for (
int i=1; i<=3; i++) {
506 if(source_shift(i).get_etat() != ETATZERO) {
507 if(source_shift(i).dz_nonzero()) {
508 assert( source_shift(i).get_dzpuis() == 4 ) ;
511 (source_shift.set(i)).set_dzpuis(4) ;
516 double lambda_shift = double(1) / double(3) ;
518 if ( mg->
get_np(0) == 1 ) {
524 for (
int i=1; i<=3; i++) {
525 csource_shift.set(i-1) = source_shift(i) ;
527 csource_shift.set(2).set_etat_zero() ;
529 csource_shift.poisson_vect(lambda_shift, par_poisson_vect,
530 cshift, cw_shift, ckhi_shift) ;
532 for (
int i=1; i<=3; i++) {
539 cout <<
"Test of the Poisson equation for shift_x :" << endl ;
540 err = source_shift(1).test_poisson(-
beta(1), cout,
true) ;
541 diff_shift_x = err(0, 0) ;
543 cout <<
"Test of the Poisson equation for shift_y :" << endl ;
544 err = source_shift(2).test_poisson(-
beta(2), cout,
true) ;
545 diff_shift_y = err(0, 0) ;
559 if (mer > mer_fix_omega + delta_mer_kep) {
561 omega *= fact_omega ;
564 bool omega_trop_grand = false ;
571 bool superlum = true ;
594 for (
int l=0; l<
nzet; l++) {
595 for (
int i=0; i<mg->
get_nr(l); i++) {
600 cout <<
"U > c for l, i : " << l <<
" " << i
601 <<
" U = " << u1 << endl ;
606 cout <<
"**** VELOCITY OF LIGHT REACHED ****" << endl ;
607 omega /= fact_omega ;
608 cout <<
"New rotation frequency : " 609 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
610 omega_trop_grand = true ;
631 mlngamma = - 0.5 *
uuu*
uuu ;
637 double mlngamma_b = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ;
642 double mlngamma_c = 0 ;
646 double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
647 + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
648 alpha_r =
sqrt(alpha_r2) ;
649 cout <<
"alpha_r = " << alpha_r << endl ;
660 ent = (ent_c + nu_c + mlngamma_c) -
logn - mlngamma ;
667 for (
int l=0; l<
nzet; l++) {
668 int imax = mg->
get_nr(l) - 1 ;
669 if (l == l_b) imax-- ;
670 for (
int i=0; i<imax; i++) {
673 cout <<
"ent < 0 for l, i : " << l <<
" " << i
680 cout <<
"**** KEPLERIAN VELOCITY REACHED ****" << endl ;
681 omega /= fact_omega ;
682 cout <<
"New rotation frequency : " 683 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
684 omega_trop_grand = true ;
689 if ( omega_trop_grand ) {
691 fact_omega =
sqrt( fact_omega ) ;
692 cout <<
"**** New fact_omega : " << fact_omega << endl ;
704 double rap_dent = fabs( dent_eq / dent_pole ) ;
705 cout <<
"| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
707 if ( rap_dent < thres_adapt ) {
709 cout <<
"******* FROZEN MAPPING *********" << endl ;
755 Cmp csource_tggg(source_tggg) ;
766 Cmp csource_dzf(source_dzf) ;
767 Cmp csource_dzq(source_dzq) ;
769 mp.
poisson2d(csource_dzf, csource_dzq, par_poisson_dzeta,
773 err_grv2 = lbda_grv2 - 1;
774 cout <<
"GRV2: " << err_grv2 << endl ;
789 logn = relax *
logn + relax_prev * logn_prev ;
791 dzeta = relax *
dzeta + relax_prev * dzeta_prev ;
803 fichfreq <<
" " << setprecision(16) <<
omega / (2*M_PI) * f_unit ;
804 fichevol <<
" " << setprecision(16) << rap_dent ;
806 fichevol <<
" " << setprecision(16) << ent_c ;
812 if (mer > mer_mass) {
815 if (mbar_wanted > 0.) {
816 xx =
mass_b() / mbar_wanted - 1. ;
817 cout <<
"Discrep. baryon mass <-> wanted bar. mass : " << xx
821 xx =
mass_g() / fabs(mbar_wanted) - 1. ;
822 cout <<
"Discrep. grav. mass <-> wanted grav. mass : " << xx
825 double xprog = ( mer > 2*mer_mass) ? 1. :
826 double(mer-mer_mass)/double(mer_mass) ;
828 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
829 double fact =
pow(ax, aexp_mass) ;
830 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
831 xx <<
" " << ax <<
" " << fact << endl ;
837 if (mer%4 == 0)
omega *= fact ;
847 diff_ent = diff_ent_tbl(0) ;
848 for (
int l=1; l<
nzet; l++) {
849 diff_ent += diff_ent_tbl(l) ;
853 fichconv <<
" " <<
log10( fabs(diff_ent) + 1.e-16 ) ;
854 fichconv <<
" " <<
log10( fabs(err_grv2) + 1.e-16 ) ;
855 fichconv <<
" " <<
log10( fabs(max_triax) + 1.e-16 ) ;
858 if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
859 vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
862 fichconv <<
" " << vit_triax ;
871 max_triax_prev = max_triax ;
890 for (
int i=1; i<=3; i++) {
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.
virtual double mass_g() const
Gravitational mass.
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
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.
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Scalar dzeta
Metric potential .
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Scalar a_car
Square of the metric factor A.
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.
Map & mp
Mapping associated with the star.
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
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.
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Scalar bbb
Metric factor B.
Standard units of space, time and mass.
Tensor up(int ind, const Metric &gam) const
Computes a new tensor by raising an index of *this.
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)
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 double grv2() const
Error on the virial identity GRV2.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
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.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
void update_metric()
Computes metric coefficients from known potentials.
Values and coefficients of a (real-value) function.
Scalar nuq
Part of the Metric potential = logn generated by the quadratic terms.
double omega
Rotation angular velocity ([f_unit] )
Scalar nbar
Baryon density in the fluid frame.
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Tbl test_poisson(const Scalar &uu, ostream &ostr, bool detail=false) const
Checks if a Poisson equation with *this as a source has been correctly solved.
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Tensor field of valence 1.
Cmp cos(const Cmp &)
Cosine.
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
void set_dzpuis(int)
Modifies the dzpuis flag.
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Coord phi
coordinate centered on the grid
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
int nzet
Number of domains of *mp occupied by the star.
Scalar nphi
Metric coefficient .
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
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.
Scalar press
Fluid pressure.
virtual double mass_b() const
Baryon mass.
Scalar ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
int get_nzone() const
Returns the number of domains.
virtual void homothetie(double lambda)
Sets a new radial scale.
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Sym_tensor tkij
Tensor related to the extrinsic curvature tensor by .
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 .
Scalar logn
Logarithm of the lapse N .
Active physical coordinates and mapping derivatives.
double ray_pole() const
Coordinate radius at [r_unit].
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
virtual void partial_display(ostream &) const
Printing of some informations, excluding all global quantities.
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Scalar nn
Lapse function N .
Cmp log10(const Cmp &)
Basis 10 logarithm.
Cmp abs(const Cmp &)
Absolute value.
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Scalar uuu
Norm of u_euler.
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.
double ray_eq() const
Coordinate radius at , [r_unit].
const Scalar & dsdr() const
Returns of *this .
Scalar ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Valeur & set_spectral_va()
Returns va (read/write version)
Scalar & set(int)
Read/write access to a component.
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.
Scalar tggg
Metric potential .
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Vector w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Scalar ener_euler
Total energy density in the Eulerian frame.
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
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.
const Valeur & get_spectral_va() const
Returns va (read only version)