Values and coefficients of a (real-value) function. More...
#include <valeur.h>
Public Member Functions | |
Valeur (const Mg3d &mgrid) | |
Constructor. | |
Valeur (const Mg3d *p_mgrid) | |
Constructor. | |
Valeur (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE*) ). | |
Valeur (const Valeur &) | |
Copy constructor. | |
~Valeur () | |
Destructor. | |
void | operator= (const Valeur &a) |
Assignement to another Valeur . | |
void | operator= (const Mtbl &mt) |
Assignement to a Mtbl . | |
void | operator= (const Mtbl_cf &mtcf) |
Assignement to a Mtbl_cf . | |
void | operator= (double) |
Assignement to a double . | |
Tbl & | set (int l) |
Read/write of the value in a given domain (configuration space). | |
const Tbl & | operator() (int l) const |
Read-only of the value in a given domain (configuration space). | |
double & | set (int l, int k, int j, int i) |
Read/write of a particular element (configuration space). | |
double | operator() (int l, int k, int j, int i) const |
Read-only of a particular element (configuration space). | |
double | val_point (int l, double x, double theta, double phi) const |
Computes the value of the field represented by *this at an arbitrary point, by means of the spectral expansion. | |
double | val_point_jk (int l, double x, int j, int k) const |
Computes the value of the field represented by *this at an arbitrary point in , but collocation point in , by means of the spectral expansion. | |
void | coef () const |
Computes the coeffcients of *this . | |
void | coef_i () const |
Computes the physical value of *this . | |
void | ylm () |
Computes the coefficients of *this . | |
void | ylm_i () |
Inverse of ylm() . | |
void | val_propre_1d () |
Set the basis to the eigenvalues of . | |
void | val_propre_1d_i () |
Inverse transformation of val_propre_1d . | |
const Base_val & | get_base () const |
Return the bases for spectral expansions (member base ). | |
void | set_base (const Base_val &) |
Sets the bases for spectral expansions (member base ). | |
void | std_base_scal () |
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar. | |
void | std_base_scal_odd () |
Sets the bases for spectral expansions (member base ) to the standard odd ones for a scalar. | |
void | set_base_r (int l, int base_r) |
Sets the expansion basis for r () functions in a given domain. | |
void | set_base_t (int base_t) |
Sets the expansion basis for functions in all domains. | |
void | set_base_p (int base_p) |
Sets the expansion basis for functions in all domains. | |
void | filtre_tp (int nn, int nz1, int nz2) |
Sets the n lasts coefficients in to 0 in the domain nz1 to nz2 when expressed in spherical harmonics. | |
const Valeur & | dsdx () const |
Returns of *this . | |
const Valeur & | d2sdx2 () const |
Returns of *this . | |
const Valeur & | dsdt () const |
Returns of *this . | |
const Valeur & | d2sdt2 () const |
Returns of *this . | |
const Valeur & | ssint () const |
Returns of *this . | |
const Valeur & | scost () const |
Returns of *this . | |
const Valeur & | mult_ct () const |
Returns applied to *this . | |
const Valeur & | mult_st () const |
Returns applied to *this . | |
const Valeur & | dsdp () const |
Returns of *this . | |
const Valeur & | stdsdp () const |
Returns of *this . | |
const Valeur & | d2sdp2 () const |
Returns of *this . | |
const Valeur & | mult_cp () const |
Returns applied to *this . | |
const Valeur & | mult_sp () const |
Returns applied to *this . | |
const Valeur & | lapang () const |
Returns the angular Laplacian of *this . | |
const Valeur & | sx () const |
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ). | |
const Valeur & | sx2 () const |
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ). | |
const Valeur & | mult_x () const |
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ). | |
void | sxm1_zec () |
Applies the following operator to *this : \ Id (r sampling = RARE , FIN ) \ (r -sampling = UNSURR ). | |
void | mult_xm1_zec () |
Applies the following operator to *this : \ Id (r sampling = RARE , FIN ) \ (r -sampling = UNSURR ). | |
void | mult2_xm1_zec () |
Applies the following operator to *this : \ Id (r sampling = RARE , FIN ) \ (r -sampling = UNSURR ). | |
void | va_x () |
Returns (r -sampling = RARE ) \ (r sampling = FIN ) \ (r -sampling = UNSURR ). | |
void | sauve (FILE *) const |
Save in a file. | |
void | display_coef (double threshold=1.e-7, int precision=4, ostream &ostr=cout) const |
Displays the spectral coefficients and the associated basis functions. | |
void | affiche_seuil (ostream &ostr, int type=0, int precision=4, double threshold=1.e-7) const |
Prints only the values greater than a given threshold. | |
void | set_etat_nondef () |
Sets the logical state to ETATNONDEF (undefined). | |
void | set_etat_zero () |
Sets the logical state to ETATZERO (zero). | |
void | set_etat_c_qcq () |
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c ). | |
void | set_etat_cf_qcq () |
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_cf ). | |
void | annule_hard () |
Sets the Valeur to zero in a hard way. | |
void | annule (int l) |
Sets the Valeur to zero in a given domain. | |
void | annule (int l_min, int l_max) |
Sets the Valeur to zero in several domains. | |
int | get_etat () const |
Returns the logical state. | |
const Mg3d * | get_mg () const |
Returns the Mg3d on which the this is defined. | |
void | operator+= (const Valeur &) |
+= Valeur | |
void | operator-= (const Valeur &) |
-= Valeur | |
void | operator*= (const Valeur &) |
*= Valeur | |
void | equipot (double uu0, int nz_search, double precis, int nitermax, int &niter, Itbl &l_iso, Tbl &xi_iso) const |
Determines an equipotential surface of the field represented by *this (inward search). | |
void | equipot_outward (double uu0, int nz_search, double precis, int nitermax, int &niter, Itbl &l_iso, Tbl &xi_iso) const |
Determines an equipotential surface of the field represented by *this (outward search). | |
void | smooth (int nzet, Valeur &uuva) const |
Changes the function *this as a smooth one when there exists a discontinuity between the nucleus and the first shell. | |
Public Attributes | |
Mtbl * | c |
Values of the function at the points of the multi-grid. | |
Mtbl_cf * | c_cf |
Coefficients of the spectral expansion of the function. | |
Base_val | base |
Bases on which the spectral expansion is performed. | |
Private Member Functions | |
void | nouveau () |
Memory allocation. | |
void | del_t () |
Logical destructor. | |
void | del_deriv () |
Logical destructor of the derivatives. | |
void | set_der_0x0 () |
Sets the pointers for derivatives to 0x0. | |
Private Attributes | |
const Mg3d * | mg |
Multi-grid Mgd3 on which this is defined. | |
int | etat |
Logical state (ETATNONDEF , ETATQCQ or ETATZERO ). | |
Valeur * | p_dsdx |
Pointer on . | |
Valeur * | p_d2sdx2 |
Pointer on . | |
Valeur * | p_sx |
Pointer on . | |
Valeur * | p_sx2 |
Pointer on . | |
Valeur * | p_mult_x |
Pointer on . | |
Valeur * | p_dsdt |
Pointer on . | |
Valeur * | p_d2sdt2 |
Pointer on . | |
Valeur * | p_ssint |
Pointer on . | |
Valeur * | p_scost |
Pointer on . | |
Valeur * | p_mult_ct |
Pointer on . | |
Valeur * | p_mult_st |
Pointer on . | |
Valeur * | p_dsdp |
Pointer on . | |
Valeur * | p_stdsdp |
Pointer on . | |
Valeur * | p_d2sdp2 |
Pointer on . | |
Valeur * | p_mult_cp |
Pointer on . | |
Valeur * | p_mult_sp |
Pointer on . | |
Valeur * | p_lapang |
Pointer on the angular Laplacian. | |
Friends | |
class | Cmp |
Friend class. | |
class | Scalar |
Friend class. | |
ostream & | operator<< (ostream &, const Valeur &) |
Display. | |
void | rotate_propre_pair (Valeur &, bool) |
Friend fonction. | |
void | rotate_propre_impair (Valeur &, bool) |
Friend fonction. |
Values and coefficients of a (real-value) function.
()
Definition at line 283 of file valeur.h.
Valeur::Valeur | ( | const Mg3d & | mgrid | ) | [explicit] |
Valeur::Valeur | ( | const Mg3d * | p_mgrid | ) | [explicit] |
Valeur::Valeur | ( | const Mg3d & | g, | |
FILE * | fd | |||
) |
Constructor from a file (see sauve(FILE*)
).
Definition at line 259 of file valeur.C.
References c, c_cf, etat, fread_be(), mg, and set_der_0x0().
Valeur::Valeur | ( | const Valeur & | vc | ) |
Copy constructor.
Definition at line 209 of file valeur.C.
References c, c_cf, del_deriv(), etat, get_etat(), get_mg(), mg, nouveau(), set_etat_nondef(), and set_etat_zero().
void Valeur::affiche_seuil | ( | ostream & | ostr, | |
int | type = 0 , |
|||
int | precision = 4 , |
|||
double | threshold = 1.e-7 | |||
) | const |
Prints only the values greater than a given threshold.
ostr | [input] Output stream used for the printing | |
type | [input] Type of display : 0 = prints only the coefficients, 1 = prints only the values in configuration space, 2 = prints both | |
precision | [input] Number of printed digits (default: 4) | |
threshold | [input] Value above which an array element is printed (default: 1.e-7) |
Definition at line 554 of file valeur.C.
References Mtbl::affiche_seuil(), Mtbl_cf::affiche_seuil(), c, c_cf, coef(), coef_i(), and etat.
void Valeur::annule | ( | int | l_min, | |
int | l_max | |||
) |
Sets the Valeur
to zero in several domains.
l_min | [input] The Valeur will be set (logically) to zero in the domains whose indices are in the range [l_min,l_max] . | |
l_max | [input] see the comments for l_min . |
Note that annule
(0,mg->get_nzone()-1) is equivalent to set_etat_zero()
.
Definition at line 750 of file valeur.C.
References annule(), Mtbl_cf::annule(), Mtbl::annule(), c, c_cf, etat, Mg3d::get_nzone(), mg, p_d2sdp2, p_d2sdt2, p_d2sdx2, p_dsdp, p_dsdt, p_dsdx, p_lapang, p_mult_cp, p_mult_ct, p_mult_sp, p_mult_st, p_mult_x, p_scost, p_ssint, p_stdsdp, p_sx, p_sx2, and set_etat_zero().
void Valeur::annule | ( | int | l | ) |
void Valeur::annule_hard | ( | ) |
Sets the Valeur
to zero in a hard way.
1/ Sets the logical state to ETATQCQ
, i.e. to an ordinary state. 2/ Allocates the memory for c
and c_cf
, and fills it with zeros. NB: this function must be used for debugging purposes only. For other operations, the functions set_etat_zero()
or annule(int,int)
must be perferred.
Definition at line 719 of file valeur.C.
References Mtbl_cf::annule_hard(), Mtbl::annule_hard(), base, c, c_cf, del_deriv(), etat, and mg.
void Valeur::coef | ( | ) | const |
Computes the coeffcients of *this
.
Definition at line 144 of file valeur_coef.C.
References Base_val::b, base, c, c_cf, etat, Mg3d::get_colloc_r(), Tbl::get_dim(), Tbl::get_etat(), Mg3d::get_nzone(), MAX_BASE, MAX_BASE_2, mg, MSQ_P, MSQ_R, MSQ_T, NONDEF, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, R_LEG, R_LEGI, R_LEGP, Tbl::set_etat_qcq(), Mtbl_cf::set_etat_qcq(), Tbl::set_etat_zero(), sqrt(), Tbl::t, Mtbl_cf::t, Mtbl::t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_LEG, T_LEG_I, T_LEG_II, T_LEG_IP, T_LEG_MI, T_LEG_MP, T_LEG_P, T_LEG_PI, T_LEG_PP, T_SIN, T_SIN_I, T_SIN_P, TRA_P, TRA_R, and TRA_T.
void Valeur::coef_i | ( | ) | const |
Computes the physical value of *this
.
Definition at line 136 of file valeur_coef_i.C.
References Base_val::b, base, Mtbl_cf::base, c, c_cf, etat, Tbl::get_dim(), Tbl::get_etat(), Mg3d::get_nzone(), Tbl::get_taille(), MAX_BASE, MAX_BASE_2, mg, MSQ_P, MSQ_R, MSQ_T, NONDEF, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, R_LEG, R_LEGI, R_LEGP, Tbl::set_etat_qcq(), Mtbl::set_etat_qcq(), Tbl::set_etat_zero(), sqrt(), Tbl::t, Mtbl_cf::t, Mtbl::t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_LEG, T_LEG_I, T_LEG_II, T_LEG_IP, T_LEG_MI, T_LEG_MP, T_LEG_P, T_LEG_PI, T_LEG_PP, T_SIN, T_SIN_I, T_SIN_P, TRA_P, TRA_R, and TRA_T.
const Valeur & Valeur::d2sdp2 | ( | ) | const |
Returns of *this
.
Definition at line 88 of file valeur_d2sdp2.C.
References Mtbl_cf::base, base, c_cf, coef(), Mtbl_cf::d2sdp2(), etat, mg, p_d2sdp2, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::d2sdt2 | ( | ) | const |
Returns of *this
.
Definition at line 103 of file valeur_d2sdt2.C.
References Mtbl_cf::base, base, c_cf, coef(), Mtbl_cf::d2sdt2(), etat, mg, p_d2sdt2, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::d2sdx2 | ( | ) | const |
Returns of *this
.
Definition at line 110 of file valeur_d2sdx2.C.
References Mtbl_cf::base, base, c_cf, coef(), Mtbl_cf::d2sdx2(), etat, mg, p_d2sdx2, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
void Valeur::del_deriv | ( | ) | [private] |
void Valeur::del_t | ( | ) | [private] |
void Valeur::display_coef | ( | double | threshold = 1.e-7 , |
|
int | precision = 4 , |
|||
ostream & | ostr = cout | |||
) | const |
Displays the spectral coefficients and the associated basis functions.
This function shows only the values greater than a given threshold.
threshold | [input] Value above which a coefficient is printed (default: 1.e-7) | |
precision | [input] Number of printed digits (default: 4) | |
ostr | [input] Output stream used for the printing (default: cout) |
Definition at line 532 of file valeur.C.
References c_cf, coef(), Mtbl_cf::display(), and etat.
const Valeur & Valeur::dsdp | ( | ) | const |
Returns of *this
.
Definition at line 94 of file valeur_dsdp.C.
References Mtbl_cf::base, base, c_cf, coef(), Mtbl_cf::dsdp(), etat, mg, p_dsdp, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::dsdt | ( | ) | const |
Returns of *this
.
Definition at line 108 of file valeur_dsdt.C.
References Mtbl_cf::base, base, c_cf, coef(), Mtbl_cf::dsdt(), etat, mg, p_dsdt, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::dsdx | ( | ) | const |
Returns of *this
.
Definition at line 107 of file valeur_dsdx.C.
References Mtbl_cf::base, base, c_cf, coef(), Mtbl_cf::dsdx(), etat, mg, p_dsdx, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
void Valeur::equipot | ( | double | uu0, | |
int | nz_search, | |||
double | precis, | |||
int | nitermax, | |||
int & | niter, | |||
Itbl & | l_iso, | |||
Tbl & | xi_iso | |||
) | const |
Determines an equipotential surface of the field represented by *this
(inward search).
The equipotential is supposed to have the form \ \ \ where l is the domain index and the radial variable in each domain.
uu0 | [input] Value defining the equipotential by u = const = uu0 where u is the field represented by *this . | |
nz_search | [input] Number of domains where the equipotential is searched : the routine scans inward the nz_search innermost domains, starting from the domain of index nz_search-1 | |
precis | [input] Required absolute precision in the determination of zeros by the secant method (standard value: 1.e-14). | |
nitermax | [input] Maximum number of iterations in the secant method (standard value: 100). | |
niter | [output] Number of iterations effectively used in the secant method | |
l_iso | [output] 2-D Itbl containing the values of l on the equipotential surface at the collocation points in [Eq. (1)], with the following storage convention \ l_iso(k,j) = | |
xi_iso | [output] 2-D Tbl containing the values of on the equipotential surface at the collocation points in [Eq. (2)], with the following storage convention \ xi_iso(k,j) = |
Definition at line 77 of file valeur_equipot.C.
References Param::add_double(), Param::add_int_mod(), Param::add_mtbl_cf(), c_cf, coef(), etat, Mg3d::get_grille3d(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), mg, Tbl::set(), Itbl::set(), Tbl::set_etat_qcq(), Itbl::set_etat_qcq(), and zerosec().
void Valeur::equipot_outward | ( | double | uu0, | |
int | nz_search, | |||
double | precis, | |||
int | nitermax, | |||
int & | niter, | |||
Itbl & | l_iso, | |||
Tbl & | xi_iso | |||
) | const |
Determines an equipotential surface of the field represented by *this
(outward search).
The equipotential is supposed to have the form \ \ \ where l is the domain index and the radial variable in each domain.
uu0 | [input] Value defining the equipotential by u = const = uu0 where u is the field represented by *this . | |
nz_search | [input] Number of domains where the equipotential is searched : the routine scans outward the nz_search innermost domains, starting from the domain of index 0 | |
precis | [input] Required absolute precision in the determination of zeros by the secant method (standard value: 1.e-14). | |
nitermax | [input] Maximum number of iterations in the secant method (standard value: 100). | |
niter | [output] Number of iterations effectively used in the secant method | |
l_iso | [output] 2-D Itbl containing the values of l on the equipotential surface at the collocation points in [Eq. (1)], with the following storage convention \ l_iso(k,j) = | |
xi_iso | [output] 2-D Tbl containing the values of on the equipotential surface at the collocation points in [Eq. (2)], with the following storage convention \ xi_iso(k,j) = |
Definition at line 79 of file valeur_equipot_out.C.
References Param::add_double(), Param::add_int_mod(), Param::add_mtbl_cf(), c_cf, coef(), etat, Mg3d::get_grille3d(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), mg, Tbl::set(), Itbl::set(), Tbl::set_etat_qcq(), Itbl::set_etat_qcq(), and zerosec().
void Valeur::filtre_tp | ( | int | nn, | |
int | nz1, | |||
int | nz2 | |||
) |
Sets the n
lasts coefficients in to 0 in the domain nz1
to nz2
when expressed in spherical harmonics.
Definition at line 920 of file valeur.C.
References base, c, c_cf, etat, get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), mg, Mtbl_cf::set(), and ylm().
const Base_val& Valeur::get_base | ( | ) | const [inline] |
int Valeur::get_etat | ( | ) | const [inline] |
const Mg3d* Valeur::get_mg | ( | ) | const [inline] |
const Valeur & Valeur::lapang | ( | ) | const |
Returns the angular Laplacian of *this
.
Definition at line 68 of file valeur_lapang.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, Mtbl_cf::lapang(), mg, p_lapang, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
void Valeur::mult2_xm1_zec | ( | ) |
Applies the following operator to *this
: \ Id (r sampling = RARE
, FIN ) \ (r -sampling = UNSURR
).
Definition at line 72 of file valeur_mult2_xm1.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, Mtbl_cf::mult2_xm1_zec(), and set_etat_cf_qcq().
const Valeur & Valeur::mult_cp | ( | ) | const |
Returns applied to *this
.
Definition at line 73 of file valeur_mult_cp.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, Mtbl_cf::mult_cp(), p_mult_cp, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::mult_ct | ( | ) | const |
Returns applied to *this
.
Definition at line 94 of file valeur_mult_ct.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, Mtbl_cf::mult_ct(), p_mult_ct, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::mult_sp | ( | ) | const |
Returns applied to *this
.
Definition at line 73 of file valeur_mult_sp.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, Mtbl_cf::mult_sp(), p_mult_sp, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::mult_st | ( | ) | const |
Returns applied to *this
.
Definition at line 93 of file valeur_mult_st.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, Mtbl_cf::mult_st(), p_mult_st, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
const Valeur & Valeur::mult_x | ( | ) | const |
Returns (r -sampling = RARE
) \ Id (r sampling = FIN
) \ (r -sampling = UNSURR
).
Definition at line 92 of file valeur_mult_x.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, Mtbl_cf::mult_x(), p_mult_x, set_etat_cf_qcq(), set_etat_zero(), and Valeur().
void Valeur::mult_xm1_zec | ( | ) |
Applies the following operator to *this
: \ Id (r sampling = RARE
, FIN ) \ (r -sampling = UNSURR
).
Definition at line 72 of file valeur_mult_xm1.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, Mtbl_cf::mult_xm1_zec(), and set_etat_cf_qcq().
void Valeur::nouveau | ( | ) | [private] |
double Valeur::operator() | ( | int | l, | |
int | k, | |||
int | j, | |||
int | i | |||
) | const [inline] |
const Tbl& Valeur::operator() | ( | int | l | ) | const [inline] |
void Valeur::operator*= | ( | const Valeur & | vi | ) |
*= Valeur
Definition at line 954 of file valeur_arithm.C.
References base, c, c_cf, coef_i(), del_deriv(), etat, get_etat(), get_mg(), mg, and set_etat_zero().
void Valeur::operator+= | ( | const Valeur & | vi | ) |
+= Valeur
Definition at line 830 of file valeur_arithm.C.
References annule_hard(), c, c_cf, coef_i(), del_deriv(), etat, get_etat(), get_mg(), and mg.
void Valeur::operator-= | ( | const Valeur & | vi | ) |
-= Valeur
Definition at line 893 of file valeur_arithm.C.
References annule_hard(), c, c_cf, coef_i(), del_deriv(), etat, get_etat(), get_mg(), and mg.
void Valeur::operator= | ( | double | x | ) |
Assignement to a double
.
Definition at line 303 of file valeur.C.
References base, c, Base_val::set_base_nondef(), set_etat_c_qcq(), and set_etat_zero().
void Valeur::operator= | ( | const Mtbl_cf & | mtcf | ) |
Assignement to a Mtbl_cf
.
Definition at line 426 of file valeur.C.
References Mtbl_cf::base, base, c, c_cf, del_deriv(), etat, Mtbl_cf::get_etat(), Mtbl_cf::get_mg(), mg, and set_etat_zero().
void Valeur::operator= | ( | const Mtbl & | mt | ) |
Assignement to a Mtbl
.
Definition at line 383 of file valeur.C.
References base, c, c_cf, del_deriv(), etat, Mtbl::get_etat(), Mtbl::get_mg(), mg, Base_val::set_base_nondef(), and set_etat_zero().
void Valeur::operator= | ( | const Valeur & | a | ) |
Assignement to another Valeur
.
Definition at line 323 of file valeur.C.
References base, c, c_cf, del_deriv(), etat, mg, set_etat_nondef(), and set_etat_zero().
void Valeur::sauve | ( | FILE * | fd | ) | const |
Save in a file.
Definition at line 471 of file valeur.C.
References base, c, coef_i(), etat, fwrite_be(), mg, Mtbl::sauve(), Mg3d::sauve(), and Base_val::sauve().
const Valeur & Valeur::scost | ( | ) | const |
Returns of *this
.
Definition at line 95 of file valeur_scost.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, p_scost, Mtbl_cf::scost(), set_etat_cf_qcq(), set_etat_zero(), and Valeur().
double& Valeur::set | ( | int | l, | |
int | k, | |||
int | j, | |||
int | i | |||
) | [inline] |
Read/write of a particular element (configuration space).
NB: to gain in efficiency, the method del_deriv()
(to delete the derived members) is not called by this function. It must thus be invoqued by the user.
l | [input] domain index | |
k | [input] index | |
j | [input] index | |
i | [input] r () index |
Definition at line 397 of file valeur.h.
References c, c_cf, coef_i(), etat, and Mtbl::set().
Tbl& Valeur::set | ( | int | l | ) | [inline] |
Read/write of the value in a given domain (configuration space).
NB: to gain in efficiency, the method del_deriv()
(to delete the derived members) is not called by this function. It must thus be invoqued by the user.
l | [input] domain index |
l
. Definition at line 359 of file valeur.h.
References c, c_cf, coef_i(), etat, and Mtbl::set().
void Valeur::set_base | ( | const Base_val & | base_i | ) |
Sets the bases for spectral expansions (member base
).
Definition at line 806 of file valeur.C.
References Mtbl_cf::base, base, and c_cf.
void Valeur::set_base_p | ( | int | base_p | ) |
Sets the expansion basis for functions in all domains.
base_p | type of basis functions in (e.g. P_COSSIN , etc..., see documentation of class Base_val for denomination of the various bases). |
Definition at line 858 of file valeur.C.
References Mtbl_cf::base, base, c_cf, and Base_val::set_base_p().
void Valeur::set_base_r | ( | int | l, | |
int | base_r | |||
) |
Sets the expansion basis for r () functions in a given domain.
l | Domain index | |
base_r | type of basis functions in r () (e.g. R_CHEBP , etc..., see documentation of class Base_val for denomination of the various bases). |
Definition at line 832 of file valeur.C.
References Mtbl_cf::base, base, c_cf, and Base_val::set_base_r().
void Valeur::set_base_t | ( | int | base_t | ) |
Sets the expansion basis for functions in all domains.
base_t | type of basis functions in (e.g. T_COS_P , etc..., see documentation of class Base_val for denomination of the various bases). |
Definition at line 845 of file valeur.C.
References Mtbl_cf::base, base, c_cf, and Base_val::set_base_t().
void Valeur::set_der_0x0 | ( | ) | [private] |
void Valeur::set_etat_c_qcq | ( | ) |
Sets the logical state to ETATQCQ
(ordinary state) for values in the configuration space (Mtbl
c
).
If c
is not 0x0, this function does nothing on c
. Otherwise, it performs the memory allocation for c
. In all cases, this function deallocates the memory occupied by the Mtbl_cf
c_cf
, as well as by all the derivatives.
Definition at line 697 of file valeur.C.
References c, c_cf, del_deriv(), etat, and mg.
void Valeur::set_etat_cf_qcq | ( | ) |
Sets the logical state to ETATQCQ
(ordinary state) for values in the configuration space (Mtbl_cf
c_cf
).
If c_cf
is not 0x0, this function does nothing on c_cf
. Otherwise, it performs the memory allocation for c_cf
. In all cases, this function deallocates the memory occupied by the Mtbl
c
, as well as by all the derivatives.
void Valeur::set_etat_nondef | ( | ) |
void Valeur::set_etat_zero | ( | ) |
void Valeur::smooth | ( | int | nzet, | |
Valeur & | uuva | |||
) | const |
Changes the function *this
as a smooth one when there exists a discontinuity between the nucleus and the first shell.
nzet | [input] Number of domains covering a star. | |
uuva | [output] Smoothed function. |
Definition at line 72 of file valeur_smooth.C.
References c_cf, coef_i(), etat, Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), mg, pow(), Mtbl_cf::set(), Tbl::set(), set_etat_cf_qcq(), Mtbl_cf::set_etat_qcq(), Tbl::set_etat_qcq(), and Mtbl_cf::t.
const Valeur & Valeur::ssint | ( | ) | const |
Returns of *this
.
Definition at line 108 of file valeur_ssint.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, p_ssint, set_etat_cf_qcq(), set_etat_zero(), Mtbl_cf::ssint(), and Valeur().
void Valeur::std_base_scal | ( | ) |
Sets the bases for spectral expansions (member base
) to the standard ones for a scalar.
Definition at line 820 of file valeur.C.
References mg, set_base(), and Mg3d::std_base_scal().
void Valeur::std_base_scal_odd | ( | ) |
Sets the bases for spectral expansions (member base
) to the standard odd ones for a scalar.
Definition at line 826 of file valeur.C.
References mg, set_base(), and Mg3d::std_base_scal_odd().
const Valeur & Valeur::stdsdp | ( | ) | const |
const Valeur & Valeur::sx | ( | ) | const |
Returns (r -sampling = RARE
) \ Id (r sampling = FIN
) \ (r -sampling = UNSURR
).
Definition at line 96 of file valeur_sx.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, p_sx, set_etat_cf_qcq(), set_etat_zero(), Mtbl_cf::sx(), and Valeur().
const Valeur & Valeur::sx2 | ( | ) | const |
Returns (r -sampling = RARE
) \ Id (r sampling = FIN
) \ (r -sampling = UNSURR
).
Definition at line 105 of file valeur_sx2.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, mg, p_sx2, set_etat_cf_qcq(), set_etat_zero(), Mtbl_cf::sx2(), and Valeur().
void Valeur::sxm1_zec | ( | ) |
Applies the following operator to *this
: \ Id (r sampling = RARE
, FIN ) \ (r -sampling = UNSURR
).
Definition at line 77 of file valeur_sxm1_zec.C.
References Mtbl_cf::base, base, c_cf, coef(), etat, set_etat_cf_qcq(), and Mtbl_cf::sxm1_zec().
void Valeur::va_x | ( | ) |
Returns (r -sampling = RARE
) \ (r sampling = FIN
) \ (r -sampling = UNSURR
).
Definition at line 42 of file valeur_x.C.
References etat, Mg3d::get_grille3d(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), mg, and Grille3d::x.
double Valeur::val_point | ( | int | l, | |
double | x, | |||
double | theta, | |||
double | phi | |||
) | const |
Computes the value of the field represented by *this
at an arbitrary point, by means of the spectral expansion.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
theta | [input] value of the coordinate | |
phi | [input] value of the coordinate |
*this
. Definition at line 878 of file valeur.C.
References c_cf, coef(), etat, and Mtbl_cf::val_point().
double Valeur::val_point_jk | ( | int | l, | |
double | x, | |||
int | j, | |||
int | k | |||
) | const |
Computes the value of the field represented by *this
at an arbitrary point in , but collocation point in , by means of the spectral expansion.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
j | [input] index of the collocation point in | |
k | [input] index of the collocation point in |
*this
. Definition at line 896 of file valeur.C.
References c_cf, coef(), etat, and Mtbl_cf::val_point_jk().
void Valeur::val_propre_1d | ( | ) |
Set the basis to the eigenvalues of .
Definition at line 318 of file valeur_val_propre_1d.C.
References base, Base_val::get_base_t(), Mg3d::get_nzone(), mg, rotate_propre_impair, rotate_propre_pair, T_CL_COS_I, T_CL_COS_P, T_CL_SIN_I, T_CL_SIN_P, T_COS_I, T_COS_P, T_SIN_I, and T_SIN_P.
void Valeur::val_propre_1d_i | ( | ) |
Inverse transformation of val_propre_1d
.
Definition at line 354 of file valeur_val_propre_1d.C.
References base, Base_val::get_base_t(), Mg3d::get_nzone(), mg, rotate_propre_impair, rotate_propre_pair, T_CL_COS_I, T_CL_COS_P, T_CL_SIN_I, T_CL_SIN_P, T_COS_I, T_COS_P, T_SIN_I, and T_SIN_P.
void Valeur::ylm | ( | ) |
Computes the coefficients of *this
.
Definition at line 134 of file valeur_ylm.C.
References Base_val::b, Mtbl_cf::base, base, c_cf, coef(), etat, Tbl::get_etat(), get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Tbl::get_taille(), MAX_BASE, MSQ_P, MSQ_R, MSQ_T, NONDEF, P_COSSIN_I, P_COSSIN_P, sqrt(), Tbl::t, Mtbl_cf::t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_LEG, T_LEG_I, T_LEG_II, T_LEG_IP, T_LEG_MI, T_LEG_MP, T_LEG_P, T_LEG_PI, T_LEG_PP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
void Valeur::ylm_i | ( | ) |
Inverse of ylm()
.
Definition at line 127 of file valeur_ylm_i.C.
References Base_val::b, Mtbl_cf::base, base, c_cf, coef(), etat, Tbl::get_etat(), get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Tbl::get_taille(), MAX_BASE, MSQ_P, MSQ_R, MSQ_T, NONDEF, P_COSSIN_I, P_COSSIN_P, sqrt(), Tbl::t, Mtbl_cf::t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_LEG, T_LEG_I, T_LEG_II, T_LEG_IP, T_LEG_MI, T_LEG_MP, T_LEG_P, T_LEG_PI, T_LEG_PP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
ostream& operator<< | ( | ostream & | , | |
const Valeur & | ||||
) | [friend] |
Display.
void rotate_propre_impair | ( | Valeur & | , | |
bool | ||||
) | [friend] |
Friend fonction.
void rotate_propre_pair | ( | Valeur & | , | |
bool | ||||
) | [friend] |
Friend fonction.
Mtbl_cf* Valeur::c_cf [mutable] |
int Valeur::etat [private] |
const Mg3d* Valeur::mg [private] |
Valeur* Valeur::p_d2sdp2 [mutable, private] |
Valeur* Valeur::p_d2sdt2 [mutable, private] |
Valeur* Valeur::p_d2sdx2 [mutable, private] |
Valeur* Valeur::p_dsdp [mutable, private] |
Valeur* Valeur::p_dsdt [mutable, private] |
Valeur* Valeur::p_dsdx [mutable, private] |
Valeur* Valeur::p_lapang [mutable, private] |
Valeur* Valeur::p_mult_cp [mutable, private] |
Valeur* Valeur::p_mult_ct [mutable, private] |
Valeur* Valeur::p_mult_sp [mutable, private] |
Valeur* Valeur::p_mult_st [mutable, private] |
Valeur* Valeur::p_mult_x [mutable, private] |
Valeur* Valeur::p_scost [mutable, private] |
Valeur* Valeur::p_ssint [mutable, private] |
Valeur* Valeur::p_stdsdp [mutable, private] |
Valeur* Valeur::p_sx [mutable, private] |
Valeur* Valeur::p_sx2 [mutable, private] |