26char star_bin_equilibrium_xcts_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Star/star_bin_equilibrium_xcts.C,v 1.13 2014/10/13 08:53:38 j_novak Exp $" ;
84#include "utilitaires.h"
98 const Tbl& fact_resize,
283 cout <<
"-----------------------------------------------" <<
endl ;
338 cout <<
"k_b, j_b, alpha_r: " <<
k_b <<
" " <<
j_b <<
" "
388 cout <<
"******* FROZEN MAPPING *********" <<
endl ;
427 for (
int i=
nzet-1;
i<nz-4;
i++) {
506 cout <<
"Max. abs(enthalpy) at the boundary between domains nzet-1"
507 <<
" and nzet : " <<
endl ;
534 Psi3.std_spectral_base() ;
537 sPsi7.std_spectral_base() ;
553 cout <<
"Resolution of the Poisson equation for Psi_auto:" <<
endl ;
564 "Relative error in the resolution of the equation for Psi_auto : "
566 for (
int l=0;
l<nz;
l++) {
570 <<
"==========================================================="
589 cout <<
"Resolution of the Poisson equation for chi_auto:" <<
endl ;
600 "Relative error in the resolution of the equation for chi_auto : "
603 for (
int l=0;
l<nz;
l++) {
607 <<
"==========================================================="
630 for (
int i=0;
i<3;
i++) {
660 for (
int i=1;
i<=3;
i++) {
689 "Relative error in the resolution of the equation for beta_auto : "
691 cout <<
"x component : " ;
692 for (
int l=0;
l<nz;
l++) {
696 cout <<
"y component : " ;
697 for (
int l=0;
l<nz;
l++) {
701 cout <<
"z component : " ;
702 for (
int l=0;
l<nz;
l++) {
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Time evolution with partial storage (*** under development ***).
Radial mapping of rather general form.
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
virtual void resize(int l, double lambda)=0
Rescales the outer boundary of one domain.
virtual void reevaluate_symy(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.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
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_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.
Tensor field of valence 0 (or component of a tensorial field).
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
const Scalar & dsdr() const
Returns of *this .
Valeur & set_spectral_va()
Returns va (read/write version)
Scalar Psi
Total conformal factor .
Scalar pot_centri
Centrifugal potential.
Sym_tensor haij_auto
Part of the extrinsic curvature tensor generated by beta_auto.
Scalar psi4
Conformal factor .
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star) .
void equilibrium(double ent_c, int mermax, int mermax_potvit, int mermax_poisson, double relax_poisson, double relax_potvit, double thres_adapt, const Tbl &fact, const Tbl *pent_limit, Tbl &diff)
Computes an equilibrium configuration.
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for khi.
Scalar khi
Solution for the scalar part of the vector Poisson equation for .
Vector ssjm1_wbeta
Effective source at the previous step for wbeta of the vector Poisson equation for wbeta (needed for ...
Scalar Psi_comp
Scalar field generated principally by the companion star.
Scalar ssjm1_chi
Effective source at the previous step for the resolution of the Poisson equation for \chi_auto .
Vector dcov_Psi
Covariant derivative of the conformal factor .
Vector beta_auto
Part of the shift vector generated principally by the star (Spherical components with respect to the ...
Scalar Psi_auto
Scalar field generated principally by the star.
Scalar hacar_comp
Part of the scalar generated by beta_auto and beta_comp, i.e.
Scalar chi
Total function .
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Vector dcov_chi
Covariant derivative of the function .
Scalar chi_auto
Scalar field generated principally by the star.
Vector w_beta
Solution for the vector part of the vector Poisson equation for .
Scalar chi_comp
Scalar field generated principally by the companion star.
Scalar ssjm1_psi
Effective source at the previous step for the resolution of the Poisson equation for \Psi_auto .
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
bool irrotational
true for an irrotational star, false for a corotating one
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
Scalar hacar_auto
Part of the scalar generated by beta_auto, i.e.
Scalar nn
Lapse function N .
Scalar ener_euler
Total energy density in the Eulerian frame.
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
double ray_eq() const
Coordinate radius at , [r_unit].
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Scalar press
Fluid pressure.
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Map & mp
Mapping associated with the star.
int nzet
Number of domains of *mp occupied by the star.
double ray_pole() const
Coordinate radius at [r_unit].
Tensor handling *** DEPRECATED : use class Tensor instead ***.
void smooth(int nzet, Valeur &uuva) const
Changes the function *this as a smooth one when there exists a discontinuity between the nucleus and ...
Tensor field of valence 1.
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
Scalar & set(int)
Read/write access to a component.
Cmp sqrt(const Cmp &)
Square root.
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.
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Cmp pow(const Cmp &, int)
Power .
Cmp abs(const Cmp &)
Absolute value.
Cmp log(const Cmp &)
Neperian logarithm.
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Standard units of space, time and mass.