32char et_equil_spher_regu_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Etoile/et_equil_spher_regu.C,v 1.4 2014/10/13 08:52:56 j_novak Exp $" ;
126 cout <<
"Input the regularization degree (k_div) : " ;
134 int i_b = mg->
get_nr(l_b) - 1 ;
135 int j_b = mg->
get_nt(l_b) - 1 ;
173 Cmp source_regu(
mp) ;
184 double diff_ent = 1 ;
193 for(
int mer=0 ; (diff_ent > precis) && (mer<mermax) ; mer++ ) {
195 cout <<
"-----------------------------------------------" << endl ;
196 cout <<
"step: " << mer << endl ;
197 cout <<
"alpha_r: " << alpha_r << endl ;
198 cout <<
"diff_ent = " << diff_ent << endl ;
225 (source.
set()).set_dzpuis(4) ;
254 source = - dlogn * dbeta ;
258 mpaff.
poisson(source(), par_nul, logn_quad.
set()) ;
262 mpaff.
dsdr(logn_quad(), dlogn_quad.
set()) ;
273 double nu_mat0_b =
logn_auto()(l_b, k_b, j_b, i_b) ;
274 double nu_mat0_c =
logn_auto()(0, 0, 0, 0) ;
276 double nu_quad0_b = logn_quad()(l_b, k_b, j_b, i_b) ;
277 double nu_quad0_c = logn_quad()(0, 0, 0, 0) ;
279 double alpha_r2 = ( ent_c - ent_b - nu_quad0_b + nu_quad0_c )
280 / ( qpig*(nu_mat0_b - nu_mat0_c) ) ;
282 alpha_r =
sqrt(alpha_r2) ;
303 double logn_c = logn()(0, 0, 0, 0) ;
304 ent = ent_c - logn() + logn_c ;
313 dlogn_auto = alpha_r*qpig * dlogn_auto ;
314 dlogn = dlogn_auto(0) + dlogn_quad() ;
327 - 0.5 * ( dlogn * dlogn + dbeta * dbeta ) ;
378 <<
"Characteristics of the star obtained by Etoile::equil_spher_regular : "
380 <<
"-----------------------------------------------------------------"
383 double ray =
mp.
val_r(l_b, 1., M_PI/2., 0) ;
384 cout <<
"Coordinate radius : " << ray / km <<
" km" << endl ;
386 double rcirc = ray *
sqrt(
a_car()(l_b, k_b, j_b, i_b) ) ;
388 double compact = qpig/(4.*M_PI) *
mass_g() / rcirc ;
390 cout <<
"Circumferential radius R : " << rcirc/km <<
" km" << endl ;
391 cout <<
"Baryon mass : " <<
mass_b()/msol <<
" Mo" << endl ;
392 cout <<
"Gravitational mass M : " <<
mass_g()/msol <<
" Mo" << endl ;
393 cout <<
"Compacity parameter GM/(c^2 R) : " << compact << endl ;
404 double vir_mat = source().integrale() ;
410 source = - ( dlogn * dlogn - 0.5 * tmp * tmp ) *
sqrt(
a_car()) ;
413 double vir_grav = source().integrale() ;
417 double grv3 = ( vir_mat + vir_grav ) / vir_mat ;
419 cout <<
"Virial theorem GRV3 : " << endl ;
420 cout <<
" 3P term : " << vir_mat << endl ;
421 cout <<
" grav. term : " << vir_grav << endl ;
422 cout <<
" relative error : " << grv3 << endl ;
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
virtual double der_ener_ent_p(double ent, const Param *par=0x0) const =0
Computes the logarithmic derivative from the log-enthalpy with extra parameters (virtual function im...
virtual double der_nbar_ent_p(double ent, const Param *par=0x0) const =0
Computes the logarithmic derivative from the log-enthalpy and extra parameters (virtual function imp...
Tenseur logn_auto_div
Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by t...
int nzet
Number of domains of *mp occupied by the star.
Tenseur logn_auto_regu
Regular part of the logarithm of the part of the lapse N generated principaly by the star.
int k_div
Index of regularity of the gravitational potential logn_auto .
Tenseur nnn
Total lapse function.
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Tenseur nbar
Baryon density in the fluid frame.
const Eos & eos
Equation of state of the stellar matter.
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.
virtual double mass_b() const
Baryon mass.
Tenseur ener
Total energy density in the fluid frame.
void equil_spher_regular(double ent_c, double precis=1.e-14)
Computes a spherical static configuration.
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 )
Tenseur press
Fluid pressure.
virtual double mass_g() const
Gravitational mass.
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 beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Tenseur a_car
Total conformal factor .
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
virtual void poisson_regular(const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const
Computes the solution of a scalar Poisson equation.
virtual void homothetie(double lambda)
Sets a new radial scale.
virtual void poisson(const Cmp &source, Param &par, Cmp &uu) const
Computes the solution of a scalar Poisson equation.
virtual void dsdr(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
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.
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void annule(int l)
Sets the Tenseur to zero in a given domain.
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 set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Cmp sqrt(const Cmp &)
Square root.
Cmp exp(const Cmp &)
Exponential.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Standard units of space, time and mass.