115 #include "et_rot_mag.h" 121 double fact_omega,
int nzadapt,
const Tbl& ent_limit,
122 const Itbl& icontrol,
const Tbl& control,
double mbar_wanted,
123 double aexp_mass,
Tbl& diff,
const double Q0,
const double a_j0,
124 Cmp (*f_j)(
const Cmp&,
const double),
125 Cmp (*M_j)(
const Cmp& x,
const double)) {
133 char display_bold[]=
"x[1m" ; display_bold[0] = 27 ;
134 char display_normal[] =
"x[0m" ; display_normal[0] = 27 ;
149 int i_b = mg->
get_nr(l_b) - 1 ;
150 int j_b = mg->
get_nt(l_b) - 1 ;
154 double ent_b = ent_limit(
nzet-1) ;
159 int mer_max = icontrol(0) ;
160 int mer_rot = icontrol(1) ;
161 int mer_change_omega = icontrol(2) ;
162 int mer_fix_omega = icontrol(3) ;
163 int mer_mass = icontrol(4) ;
164 int mermax_poisson = icontrol(5) ;
165 int delta_mer_kep = icontrol(6) ;
166 int mer_mag = icontrol(7) ;
167 int mer_change_mag = icontrol(8) ;
168 int mer_fix_mag = icontrol(9) ;
169 int mag_filter = icontrol(10) ;
172 if (mer_change_omega < mer_rot) {
173 cout <<
"Etoile_rot::equilibrium: mer_change_omega < mer_rot !" << endl ;
174 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
175 cout <<
" mer_rot = " << mer_rot << endl ;
178 if (mer_fix_omega < mer_change_omega) {
179 cout <<
"Etoile_rot::equilibrium: mer_fix_omega < mer_change_omega !" 181 cout <<
" mer_fix_omega = " << mer_fix_omega << endl ;
182 cout <<
" mer_change_omega = " << mer_change_omega << endl ;
188 bool change_ent = true ;
191 mer_mass =
abs(mer_mass) ;
194 double precis = control(0) ;
195 double omega_ini = control(1) ;
196 double relax = control(2) ;
197 double relax_prev = double(1) - relax ;
198 double relax_poisson = control(3) ;
199 double thres_adapt = control(4) ;
200 double precis_adapt = control(5) ;
201 double Q_ini = control(6) ;
202 double a_j_ini = control (7) ;
208 double& diff_ent = diff.
set(0) ;
219 int nz_search =
nzet + 1 ;
222 double reg_map = 1. ;
224 par_adapt.
add_int(nitermax, 0) ;
226 par_adapt.
add_int(nzadapt, 1) ;
229 par_adapt.
add_int(nz_search, 2) ;
231 par_adapt.
add_int(adapt_flag, 3) ;
250 par_adapt.
add_tbl(ent_limit, 0) ;
256 double precis_poisson = 1.e-16 ;
258 Param par_poisson_nuf ;
259 par_poisson_nuf.
add_int(mermax_poisson, 0) ;
260 par_poisson_nuf.
add_double(relax_poisson, 0) ;
261 par_poisson_nuf.
add_double(precis_poisson, 1) ;
265 Param par_poisson_nuq ;
266 par_poisson_nuq.
add_int(mermax_poisson, 0) ;
267 par_poisson_nuq.
add_double(relax_poisson, 0) ;
268 par_poisson_nuq.
add_double(precis_poisson, 1) ;
272 Param par_poisson_tggg ;
273 par_poisson_tggg.
add_int(mermax_poisson, 0) ;
274 par_poisson_tggg.
add_double(relax_poisson, 0) ;
275 par_poisson_tggg.
add_double(precis_poisson, 1) ;
281 Param par_poisson_dzeta ;
288 Param par_poisson_vect ;
290 par_poisson_vect.
add_int(mermax_poisson, 0) ;
291 par_poisson_vect.
add_double(relax_poisson, 0) ;
292 par_poisson_vect.
add_double(precis_poisson, 1) ;
301 Param par_poisson_At ;
304 par_poisson_At.
add_int(mermax_poisson, 0) ;
306 par_poisson_At.
add_double(precis_poisson, 1) ;
309 par_poisson_At.
add_int(mag_filter, 1) ;
311 Param par_poisson_Avect ;
316 par_poisson_Avect.
add_int(mermax_poisson, 0) ;
317 par_poisson_Avect.
add_double(relax_poisson, 0) ;
318 par_poisson_Avect.
add_double(precis_poisson, 1) ;
322 par_poisson_Avect.
add_int(mag_filter, 1) ;
333 double accrois_omega = (omega0 - omega_ini) /
334 double(mer_fix_omega - mer_change_omega) ;
335 double accrois_Q = (Q0 - Q_ini) /
336 double(mer_fix_mag - mer_change_mag);
337 double accrois_a_j = (a_j0 - a_j_ini) /
338 double(mer_fix_mag - mer_change_mag);
390 ofstream fichconv(
"convergence.d") ;
391 fichconv <<
"# diff_ent GRV2 " << endl ;
393 ofstream fichfreq(
"frequency.d") ;
394 fichfreq <<
"# f [Hz]" << endl ;
396 ofstream fichevol(
"evolution.d") ;
398 "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c" 402 double err_grv2 = 1 ;
408 for(
int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
410 cout <<
"-----------------------------------------------" << endl ;
411 cout <<
"step: " << mer << endl ;
412 cout <<
"diff_ent = " << display_bold << diff_ent << display_normal
414 cout <<
"err_grv2 = " << err_grv2 << endl ;
419 if (mer >= mer_rot) {
421 if (mer < mer_change_omega) {
425 if (mer <= mer_fix_omega) {
426 omega = omega_ini + accrois_omega *
427 (mer - mer_change_omega) ;
433 if (mer >= mer_mag) {
434 if (mer < mer_change_mag) {
439 if (mer <= mer_fix_mag) {
440 Q = Q_ini + accrois_Q * (mer - mer_change_mag) ;
441 a_j = a_j_ini + accrois_a_j * (mer - mer_change_mag) ;
451 f_j, par_poisson_At, par_poisson_Avect) ;
472 source_nuf = qpig *
nbar ;
494 (source_tggg.
set()).mult_rsint() ;
513 mtmp.set_etat_qcq() ;
529 for (
int i=0; i<3; i++) {
530 source_shift.set(i) += squad(i) ;
533 if (mtmp.get_etat() == ETATQCQ) {
534 if (source_shift.get_etat() == ETATZERO) {
535 source_shift.set_etat_qcq() ;
536 for (
int i=0; i<3; i++) {
537 source_shift.set(i) = mtmp(i) ;
538 source_shift.set(i).va.coef_i() ;
542 for (
int i=0; i<3; i++)
543 source_shift.set(i) += mtmp(i) ;
546 source_shift.set_std_base() ;
552 source_nuf().poisson(par_poisson_nuf,
nuf.
set()) ;
560 source_nuq().poisson(par_poisson_nuq,
nuq.
set()) ;
567 if (source_shift.get_etat() != ETATZERO) {
569 for (
int i=0; i<3; i++) {
570 if(source_shift(i).dz_nonzero()) {
571 assert( source_shift(i).get_dzpuis() == 4 ) ;
574 (source_shift.set(i)).set_dzpuis(4) ;
582 double lambda_shift = double(1) / double(3) ;
584 if ( mg->
get_np(0) == 1 ) {
588 source_shift.poisson_vect(lambda_shift, par_poisson_vect,
608 if (mer > mer_fix_omega + delta_mer_kep) {
610 omega *= fact_omega ;
613 bool omega_trop_grand = false ;
620 bool superlum = true ;
636 ((
uuu.
set()).va).set_base( (tmp.
va).base ) ;
644 for (
int l=0; l<
nzet; l++) {
645 for (
int i=0; i<mg->
get_nr(l); i++) {
647 double u1 =
uuu()(l, 0, j_b, i) ;
650 cout <<
"U > c for l, i : " << l <<
" " << i
651 <<
" U = " << u1 << endl ;
656 cout <<
"**** VELOCITY OF LIGHT REACHED ****" << endl ;
657 omega /= fact_omega ;
658 cout <<
"New rotation frequency : " 659 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
660 omega_trop_grand = true ;
681 mlngamma = - 0.5 *
uuu*
uuu ;
691 double nuf_b =
nuf()(l_b, k_b, j_b, i_b) ;
692 double nuq_b =
nuq()(l_b, k_b, j_b, i_b) ;
693 double mlngamma_b = mlngamma()(l_b, k_b, j_b, i_b) ;
694 double mag_b = mag()(l_b, k_b, j_b, i_b) ;
697 double nuf_c =
nuf()(0,0,0,0) ;
698 double nuq_c =
nuq()(0,0,0,0) ;
699 double mlngamma_c = 0 ;
700 double mag_c = mag()(0,0,0,0) ;
704 double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
705 + nuq_c - nuq_b + mag_c - mag_b)
706 / ( nuf_b - nuf_c ) ;
707 alpha_r =
sqrt(alpha_r2) ;
708 cout <<
"alpha_r = " << alpha_r << endl ;
714 double nu_c =
logn()(0,0,0,0) ;
719 ent = (ent_c + nu_c + mlngamma_c + mag_c) -
logn - mlngamma
727 for (
int l=0; l<
nzet; l++) {
728 int imax = mg->
get_nr(l) - 1 ;
729 if (l == l_b) imax-- ;
730 for (
int i=0; i<imax; i++) {
731 if (
ent()(l, 0, j_b, i) < 0. ) {
733 cout <<
"ent < 0 for l, i : " << l <<
" " << i
734 <<
" ent = " <<
ent()(l, 0, j_b, i) << endl ;
740 cout <<
"**** KEPLERIAN VELOCITY REACHED ****" << endl ;
741 omega /= fact_omega ;
742 cout <<
"New rotation frequency : " 743 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
744 omega_trop_grand = true ;
749 if ( omega_trop_grand ) {
751 fact_omega =
sqrt( fact_omega ) ;
752 cout <<
"**** New fact_omega : " << fact_omega << endl ;
762 double dent_eq =
ent().dsdr()(l_b, k_b, j_b, i_b) ;
763 double dent_pole =
ent().dsdr()(l_b, k_b, 0, i_b) ;
764 double rap_dent = fabs( dent_eq / dent_pole ) ;
765 cout <<
"| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
767 if ( rap_dent < thres_adapt ) {
769 cout <<
"******* FROZEN MAPPING *********" << endl ;
820 mp.
poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta,
823 err_grv2 = lbda_grv2 - 1;
824 cout <<
"GRV2: " << err_grv2 << endl ;
839 logn = relax *
logn + relax_prev * logn_prev ;
841 dzeta = relax *
dzeta + relax_prev * dzeta_prev ;
853 fichfreq <<
" " << setprecision(16) <<
omega / (2*M_PI) * f_unit ;
854 fichevol <<
" " << setprecision(16) << rap_dent ;
856 fichevol <<
" " << ent_c ;
862 if (mer > mer_mass) {
865 if (mbar_wanted > 0.) {
866 xx =
mass_b() / mbar_wanted - 1. ;
867 cout <<
"Discrep. baryon mass <-> wanted bar. mass : " << xx
871 xx =
mass_g() / fabs(mbar_wanted) - 1. ;
872 cout <<
"Discrep. grav. mass <-> wanted grav. mass : " << xx
875 double xprog = ( mer > 2*mer_mass) ? 1. :
876 double(mer-mer_mass)/double(mer_mass) ;
878 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
879 double fact =
pow(ax, aexp_mass) ;
880 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
881 xx <<
" " << ax <<
" " << fact << endl ;
887 if (mer%4 == 0)
omega *= fact ;
897 diff_ent = diff_ent_tbl(0) ;
898 for (
int l=1; l<
nzet; l++) {
899 diff_ent += diff_ent_tbl(l) ;
903 fichconv <<
" " << setprecision(16) <<
log10( fabs(diff_ent) + 1.e-16 ) ;
904 fichconv <<
" " <<
log10( fabs(err_grv2) + 1.e-16 ) ;
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Cmp log(const Cmp &)
Neperian logarithm.
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Tenseur E_em
electromagnetic energy density in the Eulerian frame
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.
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 ...
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
double Q
In the case of a perfect conductor, the requated baryonic charge.
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
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.
int get_etat() const
Returns the logical state.
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...
Tenseur b_car
Square of the metric factor B.
virtual double mass_g() const
Gravitational mass.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
Basic integer array class.
Tenseur press
Fluid pressure.
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.
Cmp A_phi
-component of the electromagnetic potential 1-form divided by .
void equilibrium_mag(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, const double Q0, const double a_j0, Cmp(*f_j)(const Cmp &x, const double), Cmp(*M_j)(const Cmp &x, const double))
Computes an equilibrium configuration.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void update_metric()
Computes metric coefficients from known potentials.
virtual double grv2() const
Error on the virial identity GRV2.
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
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 magnet_comput(const int adapt_flag, Cmp(*f_j)(const Cmp &x, const double), Param &par_poisson_At, Param &par_poisson_Avect)
Computes the electromagnetic quantities solving the Maxwell equations (6) and (7) of [Bocquet...
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
double a_j
Amplitude of the curent/charge function.
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.
bool is_conduct() const
Tells if the star is made of conducting or isolating material.
Tenseur nbar
Baryon density in the fluid frame.
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
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 .
Tenseur uuu
Norm of u_euler.
Tenseur tggg
Metric potential .
double omega
Rotation angular velocity ([f_unit] )
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)
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.
Cmp A_t
t-component of the elecctromagnetic potential 1-form, divided by .
virtual void MHD_comput()
Computes the electromagnetic part of the stress-energy tensor.
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.
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.
void set_etat_zero()
Sets the logical state to ETATZERO (zero state).
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 .
Standard electro-magnetic units.
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.
Tenseur Jp_em
component of the electromagnetic momentum density 3-vector, as measured in the Eulerian frame...