Valeur Class Reference
[Spectral representation.]

Values and coefficients of a (real-value) function. More...

#include <valeur.h>

List of all members.

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.
Tblset (int l)
 Read/write of the value in a given domain (configuration space).
const Tbloperator() (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 $\xi$, but collocation point in $(\theta', \phi')$, 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 $Y_l^m$ of *this.
void ylm_i ()
 Inverse of ylm().
void val_propre_1d ()
 Set the basis to the eigenvalues of $ \partial^2_theta - \cos\theta/\sin\theta \partial_\theta$.
void val_propre_1d_i ()
 Inverse transformation of val_propre_1d.
const Base_valget_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 ($\xi$) functions in a given domain.
void set_base_t (int base_t)
 Sets the expansion basis for $\theta$ functions in all domains.
void set_base_p (int base_p)
 Sets the expansion basis for $\phi$ functions in all domains.
void filtre_tp (int nn, int nz1, int nz2)
 Sets the n lasts coefficients in $\theta$ to 0 in the domain nz1 to nz2 when expressed in spherical harmonics.
const Valeurdsdx () const
 Returns $\partial / \partial \xi$ of *this.
const Valeurd2sdx2 () const
 Returns $\partial^2 / \partial \xi^2$ of *this.
const Valeurdsdt () const
 Returns $\partial / \partial \theta$ of *this.
const Valeurd2sdt2 () const
 Returns $\partial^2 / \partial \theta^2$ of *this.
const Valeurssint () const
 Returns $1 / \sin(\theta)$ of *this.
const Valeurscost () const
 Returns $1 / \cos(\theta)$ of *this.
const Valeurmult_ct () const
 Returns $\cos(\theta) \, Id$ applied to *this.
const Valeurmult_st () const
 Returns $\sin(\theta) \, Id$ applied to *this.
const Valeurdsdp () const
 Returns $\partial / \partial \phi$ of *this.
const Valeurstdsdp () const
 Returns $1/\sin(\theta) \, \partial / \partial \phi$ of *this.
const Valeurd2sdp2 () const
 Returns $\partial^2 / \partial \phi^2$ of *this.
const Valeurmult_cp () const
 Returns $\cos(\phi) \, Id$ applied to *this.
const Valeurmult_sp () const
 Returns $\sin(\phi) \, Id$ applied to *this.
const Valeurlapang () const
 Returns the angular Laplacian of *this.
const Valeursx () const
 Returns ${1 \over \xi}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over \xi-1}$ (r -sampling = UNSURR ).
const Valeursx2 () const
 Returns ${1 \over \xi^2}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over (\xi-1)^2}$ (r -sampling = UNSURR ).
const Valeurmult_x () const
 Returns $\xi \, Id$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ $(\xi-1) \, Id $ (r -sampling = UNSURR ).
void sxm1_zec ()
 Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ ${1 \over (\xi-1)}$ (r -sampling = UNSURR ).
void mult_xm1_zec ()
 Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ $(\xi-1) \, Id$ (r -sampling = UNSURR ).
void mult2_xm1_zec ()
 Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ $(\xi-1)^2 \, Id$ (r -sampling = UNSURR ).
void va_x ()
 Returns ${\xi}$ (r -sampling = RARE ) \ (r sampling = FIN ) \ ${1 \over \xi-1}$ (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 Mg3dget_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

Mtblc
 Values of the function at the points of the multi-grid.
Mtbl_cfc_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 Mg3dmg
 Multi-grid Mgd3 on which this is defined.
int etat
 Logical state (ETATNONDEF , ETATQCQ or ETATZERO ).
Valeurp_dsdx
 Pointer on $\partial / \partial \xi$.
Valeurp_d2sdx2
 Pointer on $\partial^2 / \partial \xi^2$.
Valeurp_sx
 Pointer on $1 / \xi$.
Valeurp_sx2
 Pointer on $1 / \xi^2$.
Valeurp_mult_x
 Pointer on $\xi \, Id$.
Valeurp_dsdt
 Pointer on $\partial / \partial \theta$.
Valeurp_d2sdt2
 Pointer on $\partial^2 / \partial \theta^2$.
Valeurp_ssint
 Pointer on $1 / \sin(\theta)$.
Valeurp_scost
 Pointer on $1 / \cos(\theta)$.
Valeurp_mult_ct
 Pointer on $\cos(\theta) \, Id$.
Valeurp_mult_st
 Pointer on $\sin(\theta) \, Id$.
Valeurp_dsdp
 Pointer on $\partial / \partial \phi$.
Valeurp_stdsdp
 Pointer on $1/\sin\theta \partial / \partial \phi$.
Valeurp_d2sdp2
 Pointer on $\partial^2 / \partial \phi^2$.
Valeurp_mult_cp
 Pointer on $\cos(\phi) \, Id$.
Valeurp_mult_sp
 Pointer on $\sin(\phi) \, Id$.
Valeurp_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.

Detailed Description

Values and coefficients of a (real-value) function.

()

Definition at line 283 of file valeur.h.


Constructor & Destructor Documentation

Valeur::Valeur ( const Mg3d mgrid  )  [explicit]

Constructor.

Definition at line 196 of file valeur.C.

References nouveau().

Valeur::Valeur ( const Mg3d p_mgrid  )  [explicit]

Constructor.

Definition at line 202 of file valeur.C.

References nouveau().

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().

Valeur::~Valeur (  ) 

Destructor.

Definition at line 293 of file valeur.C.

References del_t().


Member Function Documentation

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.

Parameters:
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.

Parameters:
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  ) 

Sets the Valeur to zero in a given domain.

Parameters:
l [input] Index of the domain in which the Valeur will be set (logically) to zero.

Definition at line 740 of file valeur.C.

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
void Valeur::coef_i (  )  const
const Valeur & Valeur::d2sdp2 (  )  const

Returns $\partial^2 / \partial \phi^2$ 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 $\partial^2 / \partial \theta^2$ 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 $\partial^2 / \partial \xi^2$ 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]

Logical destructor of the derivatives.

Definition at line 634 of file valeur.C.

References 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_der_0x0().

void Valeur::del_t (  )  [private]

Logical destructor.

Definition at line 622 of file valeur.C.

References c, c_cf, and del_deriv().

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.

Parameters:
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 $\partial / \partial \phi$ 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 $\partial / \partial \theta$ 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 $\partial / \partial \xi$ 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 \ $l=L(\theta', \phi') \qquad (1) $ \ $\xi = X(\theta', \phi') \qquad (2)$ \ where l is the domain index and $\xi$ the radial variable in each domain.

Parameters:
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 $(\theta', \phi')$ [Eq. (1)], with the following storage convention \ l_iso(k,j) = $L({\theta'}_j, {\phi'}_k)$
xi_iso [output] 2-D Tbl containing the values of $\xi$ on the equipotential surface at the collocation points in $(\theta', \phi')$ [Eq. (2)], with the following storage convention \ xi_iso(k,j) = $X({\theta'}_j, {\phi'}_k)$

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 \ $l=L(\theta', \phi') \qquad (1) $ \ $\xi = X(\theta', \phi') \qquad (2)$ \ where l is the domain index and $\xi$ the radial variable in each domain.

Parameters:
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 $(\theta', \phi')$ [Eq. (1)], with the following storage convention \ l_iso(k,j) = $L({\theta'}_j, {\phi'}_k)$
xi_iso [output] 2-D Tbl containing the values of $\xi$ on the equipotential surface at the collocation points in $(\theta', \phi')$ [Eq. (2)], with the following storage convention \ xi_iso(k,j) = $X({\theta'}_j, {\phi'}_k)$

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 $\theta$ 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]

Return the bases for spectral expansions (member base ).

Definition at line 476 of file valeur.h.

References base.

int Valeur::get_etat (  )  const [inline]

Returns the logical state.

Definition at line 722 of file valeur.h.

References etat.

const Mg3d* Valeur::get_mg (  )  const [inline]

Returns the Mg3d on which the this is defined.

Definition at line 725 of file valeur.h.

References mg.

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 ) \ $(\xi-1)^2 \, Id$ (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 $\cos(\phi) \, Id$ 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 $\cos(\theta) \, Id$ 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 $\sin(\phi) \, Id$ 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 $\sin(\theta) \, Id$ 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 $\xi \, Id$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ $(\xi-1) \, Id $ (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 ) \ $(\xi-1) \, Id$ (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]

Memory allocation.

Definition at line 613 of file valeur.C.

References c, c_cf, etat, and set_der_0x0().

double Valeur::operator() ( int  l,
int  k,
int  j,
int  i 
) const [inline]

Read-only of a particular element (configuration space).

Parameters:
l [input] domain index
k [input] $\phi$ index
j [input] $\theta$ index
i [input] r ($\xi$) index

Definition at line 417 of file valeur.h.

References c, coef_i(), and etat.

const Tbl& Valeur::operator() ( int  l  )  const [inline]

Read-only of the value in a given domain (configuration space).

Parameters:
l [input] domain index
Returns:
Tbl containing the value of the field in domain l .

Definition at line 377 of file valeur.h.

References c, coef_i(), and etat.

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 $1 / \cos(\theta)$ 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.

Parameters:
l [input] domain index
k [input] $\phi$ index
j [input] $\theta$ index
i [input] r ($\xi$) 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.

Parameters:
l [input] domain index
Returns:
Tbl containing the value of the field in domain 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 $\phi$ functions in all domains.

Parameters:
base_p type of basis functions in $\phi$ (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 ($\xi$) functions in a given domain.

Parameters:
l Domain index
base_r type of basis functions in r ($\xi$) (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 $\theta$ functions in all domains.

Parameters:
base_t type of basis functions in $\theta$ (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]

Sets the pointers for derivatives to 0x0.

Definition at line 661 of file valeur.C.

References 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, and p_sx2.

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.

Definition at line 708 of file valeur.C.

References base, c, c_cf, del_deriv(), etat, and mg.

void Valeur::set_etat_nondef (  ) 

Sets the logical state to ETATNONDEF (undefined).

Deallocates the memory occupied by the Mtbl c and the Mtbl_cf c_cf , as well as by all the derivatives.

Definition at line 691 of file valeur.C.

References del_t(), and etat.

void Valeur::set_etat_zero (  ) 

Sets the logical state to ETATZERO (zero).

Deallocates the memory occupied by the Mtbl c and the Mtbl_cf c_cf , as well as by all the derivatives.

Definition at line 685 of file valeur.C.

References del_t(), and etat.

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.

Parameters:
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 $1 / \sin(\theta)$ 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

Returns $1/\sin(\theta) \, \partial / \partial \phi$ of *this.

Definition at line 56 of file valeur_stdsdp.C.

References dsdp(), etat, p_stdsdp, ssint(), and Valeur().

const Valeur & Valeur::sx (  )  const

Returns ${1 \over \xi}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over \xi-1}$ (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 ${1 \over \xi^2}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over (\xi-1)^2}$ (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 ) \ ${1 \over (\xi-1)}$ (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 ${\xi}$ (r -sampling = RARE ) \ (r sampling = FIN ) \ ${1 \over \xi-1}$ (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.

Parameters:
l [input] index of the domain
x [input] value of the coordinate $\xi$
theta [input] value of the coordinate $\theta'$
phi [input] value of the coordinate $\phi'$
Returns:
value at the point $(\xi, \theta', \phi')$ in the domain no. l of the field represented by *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 $\xi$, but collocation point in $(\theta', \phi')$, by means of the spectral expansion.

Parameters:
l [input] index of the domain
x [input] value of the coordinate $\xi$
j [input] index of the collocation point in $\theta'$
k [input] index of the collocation point in $\phi'$
Returns:
value at the point $(\xi, {\theta'}_j, {\phi'}_k)$ in the domain no. l of the field represented by *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 (  ) 
void Valeur::val_propre_1d_i (  ) 
void Valeur::ylm (  ) 
void Valeur::ylm_i (  ) 

Friends And Related Function Documentation

friend class Cmp [friend]

Friend class.

Definition at line 821 of file valeur.h.

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.

friend class Scalar [friend]

Friend class.

Definition at line 822 of file valeur.h.


Member Data Documentation

Bases on which the spectral expansion is performed.

Definition at line 301 of file valeur.h.

Mtbl* Valeur::c [mutable]

Values of the function at the points of the multi-grid.

Definition at line 295 of file valeur.h.

Mtbl_cf* Valeur::c_cf [mutable]

Coefficients of the spectral expansion of the function.

Definition at line 298 of file valeur.h.

int Valeur::etat [private]

Logical state (ETATNONDEF , ETATQCQ or ETATZERO ).

Definition at line 291 of file valeur.h.

const Mg3d* Valeur::mg [private]

Multi-grid Mgd3 on which this is defined.

Definition at line 288 of file valeur.h.

Valeur* Valeur::p_d2sdp2 [mutable, private]

Pointer on $\partial^2 / \partial \phi^2$.

Definition at line 321 of file valeur.h.

Valeur* Valeur::p_d2sdt2 [mutable, private]

Pointer on $\partial^2 / \partial \theta^2$.

Definition at line 313 of file valeur.h.

Valeur* Valeur::p_d2sdx2 [mutable, private]

Pointer on $\partial^2 / \partial \xi^2$.

Definition at line 307 of file valeur.h.

Valeur* Valeur::p_dsdp [mutable, private]

Pointer on $\partial / \partial \phi$.

Definition at line 319 of file valeur.h.

Valeur* Valeur::p_dsdt [mutable, private]

Pointer on $\partial / \partial \theta$.

Definition at line 312 of file valeur.h.

Valeur* Valeur::p_dsdx [mutable, private]

Pointer on $\partial / \partial \xi$.

Definition at line 306 of file valeur.h.

Valeur* Valeur::p_lapang [mutable, private]

Pointer on the angular Laplacian.

Definition at line 325 of file valeur.h.

Valeur* Valeur::p_mult_cp [mutable, private]

Pointer on $\cos(\phi) \, Id$.

Definition at line 322 of file valeur.h.

Valeur* Valeur::p_mult_ct [mutable, private]

Pointer on $\cos(\theta) \, Id$.

Definition at line 316 of file valeur.h.

Valeur* Valeur::p_mult_sp [mutable, private]

Pointer on $\sin(\phi) \, Id$.

Definition at line 323 of file valeur.h.

Valeur* Valeur::p_mult_st [mutable, private]

Pointer on $\sin(\theta) \, Id$.

Definition at line 317 of file valeur.h.

Valeur* Valeur::p_mult_x [mutable, private]

Pointer on $\xi \, Id$.

Definition at line 310 of file valeur.h.

Valeur* Valeur::p_scost [mutable, private]

Pointer on $1 / \cos(\theta)$.

Definition at line 315 of file valeur.h.

Valeur* Valeur::p_ssint [mutable, private]

Pointer on $1 / \sin(\theta)$.

Definition at line 314 of file valeur.h.

Valeur* Valeur::p_stdsdp [mutable, private]

Pointer on $1/\sin\theta \partial / \partial \phi$.

Definition at line 320 of file valeur.h.

Valeur* Valeur::p_sx [mutable, private]

Pointer on $1 / \xi$.

Definition at line 308 of file valeur.h.

Valeur* Valeur::p_sx2 [mutable, private]

Pointer on $1 / \xi^2$.

Definition at line 309 of file valeur.h.


The documentation for this class was generated from the following files:

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1