29char et_magnetisation_C[] =
"$Header $" ;
73#include "et_rot_mag.h"
74#include "utilitaires.h"
88 const Eos& eos_i,
bool magnetisation,
91 use_B_in_eos(B_eos), include_magnetisation(magnetisation),
94 J_I(mp_i, COV, mp_i.get_bvect_spher()),
95 Sij_I(mp_i, COV, mp_i.get_bvect_spher())
98 if (tmp_p_eos == 0x0) {
99 cerr <<
"Et_magnetisation::Et_magnetisation : " << endl ;
100 cerr <<
"Only magnetised EoS is admitted for this class!" << endl ;
101 cerr <<
"Aborting ... " << endl ;
105 cerr <<
"Et_magnetisation::Et_magnetisation : " << endl ;
106 cerr <<
"Magnetisation terms can be included only if "
107 <<
"the magnetic field is used in the EoS!" << endl;
108 cerr <<
"Aborting ... " << endl ;
122 use_B_in_eos(true), include_magnetisation(true),
125 J_I(mp_i, COV, mp_i.get_bvect_spher()),
126 Sij_I(mp_i, COV, mp_i.get_bvect_spher())
156 use_B_in_eos(et.use_B_in_eos), include_magnetisation(et.include_magnetisation),
192 assert (mageos != 0x0) ;
197 double epsilon = 1.e-12 ;
204 for (
int l=0; l<nz; l++) {
206 for (
int k=0; k<mg->
get_np(l); k++) {
207 for (
int j=0; j<mg->
get_nt(l); j++) {
208 for (
int i=0; i<mg->
get_nr(l); i++) {
219 fact_ent.
set(0) = 1 + epsilon * xi(0) * xi(0) ;
220 fact_ent.
set(1) = 1 - 0.25 * epsilon * (xi(1) - 1) * (xi(1) - 1) ;
222 for (
int l=
nzet; l<nz; l++) {
223 fact_ent.
set(l) = 1 ;
229 fact_ent.
set(1) = 1 - 0.5 * epsilon * (xi(1) - 0.5) * (xi(1) - 0.5) ;
230 fact_ent.
set(2) = 1 - 0.25 * epsilon * (xi(2) - 1) * (xi(2) - 1) ;
234 cerr <<
"Et_magnetisation::equation_of_state: "
235 <<
"not ready yet for nzet > 3 !" << endl ;
238 ent_eos = fact_ent * ent_eos ;
252 norm_b =
sqrt(
a_car() * ( Bmag(0)*Bmag(0) + Bmag(1)*Bmag(1) ))
267 for (
int l=0; l< nz; l++) {
268 Tbl* tent = ent_eos.
va.
c->
t[l] ;
280 for (
int k=0; k < mg->
get_np(l); k++) {
281 for (
int j=0; j < mg->
get_nt(l); j++) {
282 for (
int i=0; i < mg->
get_nr(l); i++) {
283 magb0 = norm_b(l, k, j, i) ;
284 double ent0 = ent_eos(l, k, j, i) ;
348 ost <<
"Rotating magnetized neutron star"
351 ost <<
"Using magnetic field in the EoS" << endl ;
353 ost <<
"Including magnetisation terms in the equations" << endl ;
354 ost <<
"Maximal value of the magnetization scalar x : "
368 double fact_omega,
int nzadapt,
const Tbl& ent_limit,
369 const Itbl& icontrol,
const Tbl& control,
double mbar_wanted,
370 double magmom_wanted,
371 double aexp_mass,
Tbl& diff,
double Q0,
double a_j0,
372 Cmp (*f_j)(
const Cmp&,
const double),
373 Cmp (*M_j)(
const Cmp& x,
const double)) {
381 char display_bold[]=
"x[1m" ; display_bold[0] = 27 ;
382 char display_normal[] =
"x[0m" ; display_normal[0] = 27 ;
397 int i_b = mg->
get_nr(l_b) - 1 ;
398 int j_b = mg->
get_nt(l_b) - 1 ;
402 double ent_b = ent_limit(
nzet-1) ;
407 int mer_max = icontrol(0) ;
408 int mer_rot = icontrol(1) ;
409 int mer_change_omega = icontrol(2) ;
410 int mer_fix_omega = icontrol(3) ;
411 int mer_mass = icontrol(4) ;
412 int mermax_poisson = icontrol(5) ;
413 int delta_mer_kep = icontrol(6) ;
414 int mer_mag = icontrol(7) ;
415 int mer_change_mag = icontrol(8) ;
416 int mer_fix_mag = icontrol(9) ;
417 int mer_magmom = icontrol(10) ;
420 if (mer_change_omega < mer_rot) {
421 cerr <<
"Et_magnetisation::equilibrium_mag: mer_change_omega < mer_rot !"
423 cerr <<
" mer_change_omega = " << mer_change_omega << endl ;
424 cerr <<
" mer_rot = " << mer_rot << endl ;
427 if (mer_fix_omega < mer_change_omega) {
428 cerr <<
"Et_magnetisation::equilibrium_mag: mer_fix_omega < mer_change_omega !"
430 cerr <<
" mer_fix_omega = " << mer_fix_omega << endl ;
431 cerr <<
" mer_change_omega = " << mer_change_omega << endl ;
437 bool change_ent = true ;
440 mer_mass =
abs(mer_mass) ;
443 double precis = control(0) ;
444 double omega_ini = control(1) ;
445 double relax = control(2) ;
446 double relax_prev = double(1) - relax ;
447 double relax_poisson = control(3) ;
448 double thres_adapt = control(4) ;
449 double precis_adapt = control(5) ;
450 double Q_ini = control(6) ;
451 double a_j_ini = control (7) ;
457 double& diff_ent = diff.
set(0) ;
468 int nz_search =
nzet + 1 ;
471 double reg_map = 1. ;
473 par_adapt.
add_int(nitermax, 0) ;
475 par_adapt.
add_int(nzadapt, 1) ;
478 par_adapt.
add_int(nz_search, 2) ;
480 par_adapt.
add_int(adapt_flag, 3) ;
500 par_adapt.
add_tbl(ent_limit, 0) ;
507 double precis_poisson = 1.e-16 ;
509 Param par_poisson_nuf ;
510 par_poisson_nuf.
add_int(mermax_poisson, 0) ;
511 par_poisson_nuf.
add_double(relax_poisson, 0) ;
512 par_poisson_nuf.
add_double(precis_poisson, 1) ;
516 Param par_poisson_nuq ;
517 par_poisson_nuq.
add_int(mermax_poisson, 0) ;
518 par_poisson_nuq.
add_double(relax_poisson, 0) ;
519 par_poisson_nuq.
add_double(precis_poisson, 1) ;
523 Param par_poisson_tggg ;
524 par_poisson_tggg.
add_int(mermax_poisson, 0) ;
525 par_poisson_tggg.
add_double(relax_poisson, 0) ;
526 par_poisson_tggg.
add_double(precis_poisson, 1) ;
532 Param par_poisson_dzeta ;
539 Param par_poisson_vect ;
541 par_poisson_vect.
add_int(mermax_poisson, 0) ;
542 par_poisson_vect.
add_double(relax_poisson, 0) ;
543 par_poisson_vect.
add_double(precis_poisson, 1) ;
551 Param par_poisson_At ;
554 par_poisson_At.
add_int(mermax_poisson, 0) ;
556 par_poisson_At.
add_double(precis_poisson, 1) ;
560 Param par_poisson_Avect ;
565 par_poisson_Avect.
add_int(mermax_poisson, 0) ;
566 par_poisson_Avect.
add_double(relax_poisson, 0) ;
567 par_poisson_Avect.
add_double(precis_poisson, 1) ;
580 double accrois_omega = (omega0 - omega_ini)
581 /
double(mer_fix_omega - mer_change_omega) ;
582 double accrois_Q = (Q0 - Q_ini)
583 /
double(mer_fix_mag - mer_change_mag);
584 double accrois_a_j = (a_j0 - a_j_ini)
585 /
double(mer_fix_mag - mer_change_mag);
636 ofstream fichconv(
"convergence.d") ;
637 fichconv <<
"# diff_ent GRV2 " << endl ;
639 ofstream fichfreq(
"frequency.d") ;
640 fichfreq <<
"# f [Hz]" << endl ;
642 ofstream fichevol(
"evolution.d") ;
644 "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
648 double err_grv2 = 1 ;
654 for(
int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
656 cout <<
"-----------------------------------------------" << endl ;
657 cout <<
"step: " << mer << endl ;
658 cout <<
"diff_ent = " << display_bold << diff_ent << display_normal
660 cout <<
"err_grv2 = " << err_grv2 << endl ;
665 if (mer >= mer_rot) {
667 if (mer < mer_change_omega) {
671 if (mer <= mer_fix_omega) {
672 omega = omega_ini + accrois_omega *
673 (mer - mer_change_omega) ;
679 if (mer >= mer_mag) {
680 if (mer < mer_change_mag) {
685 if (mer <= mer_fix_mag) {
686 Q = Q_ini + accrois_Q * (mer - mer_change_mag) ;
687 a_j = a_j_ini + accrois_a_j * (mer - mer_change_mag) ;
695 magnet_comput(adapt_flag, f_j, par_poisson_At, par_poisson_Avect) ;
701 SrrplusStt = SrrplusStt /
a_car ;
728 source_nuf = qpig *
nbar ;
737 source_dzf = 2 * qpig *
a_car
750 * ( 2*
press + SrrplusStt ) ;
753 (source_tggg.
set()).mult_rsint() ;
788 for (
int i=0; i<3; i++) {
789 source_shift.
set(i) += squad(i) ;
793 if (source_shift.
get_etat() == ETATZERO) {
795 for (
int i=0; i<3; i++) {
796 source_shift.
set(i) = mtmp(i) ;
801 for (
int i=0; i<3; i++)
802 source_shift.
set(i) += mtmp(i) ;
811 source_nuf().poisson(par_poisson_nuf,
nuf.
set()) ;
819 source_nuq().poisson(par_poisson_nuq,
nuq.
set()) ;
826 if (source_shift.
get_etat() != ETATZERO) {
828 for (
int i=0; i<3; i++) {
829 if(source_shift(i).dz_nonzero()) {
830 assert( source_shift(i).get_dzpuis() == 4 ) ;
833 (source_shift.
set(i)).set_dzpuis(4) ;
841 double lambda_shift = double(1) / double(3) ;
843 if ( mg->
get_np(0) == 1 ) {
847 source_shift.
poisson_vect(lambda_shift, par_poisson_vect,
862 if (mer > mer_fix_omega + delta_mer_kep) {
864 omega *= fact_omega ;
867 bool omega_trop_grand = false ;
874 bool superlum = true ;
890 ((
uuu.
set()).va).set_base( (tmp.
va).base ) ;
897 for (
int l=0; l<
nzet; l++) {
898 for (
int i=0; i<mg->
get_nr(l); i++) {
900 double u1 =
uuu()(l, 0, j_b, i) ;
903 cout <<
"U > c for l, i : " << l <<
" " << i
904 <<
" U = " << u1 << endl ;
909 cout <<
"**** VELOCITY OF LIGHT REACHED ****" << endl ;
910 omega /= fact_omega ;
911 cout <<
"New rotation frequency : "
912 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
913 omega_trop_grand = true ;
933 mlngamma = - 0.5 *
uuu*
uuu ;
943 double nuf_b =
nuf()(l_b, k_b, j_b, i_b) ;
944 double nuq_b =
nuq()(l_b, k_b, j_b, i_b) ;
945 double mlngamma_b = mlngamma()(l_b, k_b, j_b, i_b) ;
946 double mag_b = mag()(l_b, k_b, j_b, i_b) ;
949 double nuf_c =
nuf()(0,0,0,0) ;
950 double nuq_c =
nuq()(0,0,0,0) ;
951 double mlngamma_c = 0 ;
952 double mag_c = mag()(0,0,0,0) ;
956 double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
957 + nuq_c - nuq_b + mag_c - mag_b)
958 / ( nuf_b - nuf_c ) ;
959 alpha_r =
sqrt(alpha_r2) ;
960 cout <<
"alpha_r = " << alpha_r << endl ;
966 double nu_c =
logn()(0,0,0,0) ;
970 ent = (ent_c + nu_c + mlngamma_c + mag_c) -
logn - mlngamma - mag ;
977 for (
int l=0; l<
nzet; l++) {
978 int imax = mg->
get_nr(l) - 1 ;
979 if (l == l_b) imax-- ;
980 for (
int i=0; i<imax; i++) {
981 if (
ent()(l, 0, j_b, i) < 0. ) {
983 cout <<
"ent < 0 for l, i : " << l <<
" " << i
984 <<
" ent = " <<
ent()(l, 0, j_b, i) << endl ;
990 cout <<
"**** KEPLERIAN VELOCITY REACHED ****" << endl ;
991 omega /= fact_omega ;
992 cout <<
"New rotation frequency : "
993 <<
omega/(2.*M_PI) * f_unit <<
" Hz" << endl ;
994 omega_trop_grand = true ;
999 if ( omega_trop_grand ) {
1001 fact_omega =
sqrt( fact_omega ) ;
1002 cout <<
"**** New fact_omega : " << fact_omega << endl ;
1012 double dent_eq =
ent().dsdr()(l_b, k_b, j_b, i_b) ;
1013 double dent_pole =
ent().dsdr()(l_b, k_b, 0, i_b) ;
1014 double rap_dent = fabs( dent_eq / dent_pole ) ;
1015 cout <<
"| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
1017 if ( rap_dent < thres_adapt ) {
1019 cout <<
"******* FROZEN MAPPING *********" << endl ;
1071 err_grv2 = lbda_grv2 - 1;
1072 cout <<
"GRV2: " << err_grv2 << endl ;
1086 logn = relax *
logn + relax_prev * logn_prev ;
1088 dzeta = relax *
dzeta + relax_prev * dzeta_prev ;
1100 fichfreq <<
" " <<
omega / (2*M_PI) * f_unit ;
1101 fichevol <<
" " << rap_dent ;
1103 fichevol <<
" " << ent_c ;
1109 if (mer > mer_mass) {
1112 if (mbar_wanted > 0.) {
1113 xx =
mass_b() / mbar_wanted - 1. ;
1114 cout <<
"Discrep. baryon mass <-> wanted bar. mass : " << xx
1118 xx =
mass_g() / fabs(mbar_wanted) - 1. ;
1119 cout <<
"Discrep. grav. mass <-> wanted grav. mass : " << xx
1122 double xprog = ( mer > 2*mer_mass) ? 1. :
1123 double(mer-mer_mass)/double(mer_mass) ;
1125 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
1126 double fact =
pow(ax, aexp_mass) ;
1127 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
1128 xx <<
" " << ax <<
" " << fact << endl ;
1134 if (mer%4 == 0)
omega *= fact ;
1142 if (mer > mer_magmom) {
1151 double xx =
MagMom() / magmom_wanted - 1. ;
1152 cout <<
"Discrep. mag. moment <-> wanted mag. moment : " << xx
1155 double xprog = ( mer > 2*mer_magmom) ? 1. :
1156 double(mer-mer_magmom)/double(mer_magmom) ;
1158 double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
1159 double fact =
pow(ax, aexp_mass) ;
1160 cout <<
" xprog, xx, ax, fact : " << xprog <<
" " <<
1161 xx <<
" " << ax <<
" " << fact << endl ;
1171 diff_ent = diff_ent_tbl(0) ;
1172 for (
int l=1; l<
nzet; l++) {
1173 diff_ent += diff_ent_tbl(l) ;
1177 fichconv <<
" " <<
log10( fabs(diff_ent) + 1.e-16 ) ;
1178 fichconv <<
" " <<
log10( fabs(err_grv2) + 1.e-16 ) ;
1186 dzeta_prev =
dzeta ;
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
void mult_rsint()
Multiplication by .
void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
int get_etat() const
Returns the logical state.
Valeur va
The numerical value of the Cmp
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
void annule(int l)
Sets the Cmp to zero in a given domain.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Tbl & set(int l)
Read/write of the value in a given domain.
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Class for a magnetized (tabulated) equation of state.
double mag_ent_p(double ent, const Param *par=0x0) const
Computes the magnetisation.
virtual double press_ent_p(double ent, const Param *par=0x0) const
Computes the pressure from the log-enthalpy.
virtual double nbar_ent_p(double ent, const Param *par=0x0) const
Computes the baryon density from the log-enthalpy.
virtual double ener_ent_p(double ent, const Param *par=0x0) const
Computes the total energy density from the log-enthalpy.
Equation of state base class.
virtual 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 operator=(const Et_magnetisation &)
Assignment to another Et_rot_mag.
Sym_tensor Sij_I
Interaction stress 3-tensor.
Scalar xmag
The magnetisation scalar.
Vector J_I
Interaction momentum density 3-vector.
virtual ~Et_magnetisation()
Destructor.
virtual double grv2() const
Error on the virial identity GRV2.
bool use_B_in_eos
Flag : true if the value of the magnetic field is used in the Eos.
virtual void MHD_comput()
Computes the electromagnetic part of the stress-energy tensor.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Et_magnetisation(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, bool include_mag=true, bool use_B=true)
Standard constructor.
virtual void sauve(FILE *) const
Save in a file.
bool include_magnetisation
Flag : true if magnetisation terms are included in the equations.
Scalar E_I
Interaction (magnetisation) energy density.
virtual double mass_g() const
Gravitational mass.
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 magmom_wanted, double aexp_mass, Tbl &diff, double Q0, double a_j0, Cmp(*f_j)(const Cmp &x, const double), Cmp(*M_j)(const Cmp &x, const double))
Computes an equilibrium configuration.
Class for magnetized (isolator or perfect conductor), rigidly rotating stars.
void operator=(const Et_rot_mag &)
Assignment to another Et_rot_mag.
Cmp A_phi
-component of the electromagnetic potential 1-form divided by .
double a_j
Amplitude of the curent/charge function.
double MagMom() const
Magnetic Momentum in SI units.
Cmp A_t
t-component of the elecctromagnetic potential 1-form, divided by .
virtual void del_deriv() const
Deletes all the derived quantities.
bool is_conduct() const
Tells if the star is made of conducting or isolating material.
virtual void sauve(FILE *) const
Save in a file.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Tenseur E_em
electromagnetic energy density in the Eulerian frame
Tenseur Jp_em
component of the electromagnetic momentum density 3-vector, as measured in the Eulerian frame.
double Q
In the case of a perfect conductor, the requated baryonic charge.
Tenseur Magn() const
Computes the magnetic field spherical components in Lorene's units.
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Tenseur uuu
Norm of u_euler.
double omega
Rotation angular velocity ([f_unit] )
Tenseur & logn
Metric potential = logn_auto.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Tenseur tggg
Metric potential .
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Tenseur nphi
Metric coefficient .
virtual double mass_b() const
Baryon mass.
Tenseur bbb
Metric factor B.
void update_metric()
Computes metric coefficients from known potentials.
Tenseur & dzeta
Metric potential = beta_auto.
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Tenseur b_car
Square of the metric factor B.
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
int nzet
Number of domains of *mp occupied by the star.
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Tenseur nbar
Baryon density in the fluid frame.
const Eos & eos
Equation of state of the stellar matter.
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.
Tenseur ener
Total energy density in the fluid frame.
Tenseur press
Fluid pressure.
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 a_car
Total conformal factor .
double ray_pole() const
Coordinate radius at [r_unit].
double * x
Array of values of at the nr collocation points.
Basic integer array class.
Radial mapping of rather general form.
virtual void homothetie(double lambda)
Sets a new radial scale.
Base class for coordinate mappings.
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
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.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
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.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
const Grille3d * get_grille3d(int l) const
Returns a pointer on the 3D mono-grid for domain no. l.
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.
Tbl & set(int l)
Read/write of the Tbl in a given domain.
Tbl ** t
Array (size nzone ) of pointers on the Tbl 's.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Tensor field of valence 0 (or component of a tensorial field).
virtual void sauve(FILE *) const
Save in a file.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Tbl & set_domain(int l)
Read/write of the value in a given domain.
virtual void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
double & set_grid_point(int l, int k, int j, int i)
Setting the value of the field at a given grid point.
virtual void exponential_filter_r(int lzmin, int lzmax, int p, double alpha=-16.)
Applies an exponential filter to the spectral coefficients in the radial direction.
Class intended to describe valence-2 symmetric tensors.
int get_etat() const
Gives the logical state.
void annule_hard()
Sets the Tbl to zero in a hard way.
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
double & set(int i)
Read/write of a particular element (index i) (1D case)
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 set_std_base()
Set the standard spectal basis of decomposition for each component.
void set_etat_zero()
Sets the logical state to ETATZERO (zero state).
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
int get_etat() const
Returns the logical state.
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c ).
Mtbl * c
Values of the function at the points of the multi-grid
void coef_i() const
Computes the physical value of *this.
Tensor field of valence 1.
Cmp sqrt(const Cmp &)
Square root.
Cmp log10(const Cmp &)
Basis 10 logarithm.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
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 sauve(FILE *) const
Save in a binary file.
virtual void sauve(FILE *) const
Save in a binary file.
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
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...
Tbl maxabs(const Tensor &aa, const char *comment=0x0, ostream &ost=cout, bool verb=true)
Maxima in each domain of the absolute values of the tensor components.
Standard electro-magnetic units.