LORENE
Lorene::Valeur Class Reference

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

#include <valeur.h>

Public Member Functions

 Valeur (const Mg3d &mgrid)
 Constructor. More...
 
 Valeur (const Mg3d *p_mgrid)
 Constructor. More...
 
 Valeur (const Mg3d &, FILE *)
 Constructor from a file (see sauve(FILE*) ) More...
 
 Valeur (const Valeur &)
 Copy constructor. More...
 
 ~Valeur ()
 Destructor. More...
 
void operator= (const Valeur &a)
 Assignement to another Valeur. More...
 
void operator= (const Mtbl &mt)
 Assignement to a Mtbl. More...
 
void operator= (const Mtbl_cf &mtcf)
 Assignement to a Mtbl_cf. More...
 
void operator= (double)
 Assignement to a double. More...
 
Tblset (int l)
 Read/write of the value in a given domain (configuration space). More...
 
const Tbloperator() (int l) const
 Read-only of the value in a given domain (configuration space). More...
 
double & set (int l, int k, int j, int i)
 Read/write of a particular element (configuration space). More...
 
double operator() (int l, int k, int j, int i) const
 Read-only of a particular element (configuration space). More...
 
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. More...
 
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. More...
 
void coef () const
 Computes the coeffcients of *this. More...
 
void coef_i () const
 Computes the physical value of *this. More...
 
void ylm ()
 Computes the coefficients $Y_l^m$ of *this. More...
 
void ylm_i ()
 Inverse of ylm() More...
 
void val_propre_1d ()
 Set the basis to the eigenvalues of $ \partial^2_theta - \cos\theta/\sin\theta \partial_\theta$. More...
 
void val_propre_1d_i ()
 Inverse transformation of val_propre_1d. More...
 
const Base_valget_base () const
 Return the bases for spectral expansions (member base ) More...
 
void set_base (const Base_val &)
 Sets the bases for spectral expansions (member base ) More...
 
void std_base_scal ()
 Sets the bases for spectral expansions (member base ) to the standard ones for a scalar. More...
 
void std_base_scal_odd ()
 Sets the bases for spectral expansions (member base ) to the standard odd ones for a scalar. More...
 
void set_base_r (int l, int base_r)
 Sets the expansion basis for r ( $\xi$) functions in a given domain. More...
 
void set_base_t (int base_t)
 Sets the expansion basis for $\theta$ functions in all domains. More...
 
void set_base_p (int base_p)
 Sets the expansion basis for $\phi$ functions in all domains. More...
 
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. More...
 
const Valeurdsdx () const
 Returns $\partial / \partial \xi$ of *this. More...
 
const Valeurd2sdx2 () const
 Returns $\partial^2 / \partial \xi^2$ of *this. More...
 
const Valeurdsdt () const
 Returns $\partial / \partial \theta$ of *this. More...
 
const Valeurd2sdt2 () const
 Returns $\partial^2 / \partial \theta^2$ of *this. More...
 
const Valeurssint () const
 Returns $1 / \sin(\theta)$ of *this. More...
 
const Valeurscost () const
 Returns $1 / \cos(\theta)$ of *this. More...
 
const Valeurmult_ct () const
 Returns $\cos(\theta) \, Id$ applied to *this. More...
 
const Valeurmult_st () const
 Returns $\sin(\theta) \, Id$ applied to *this. More...
 
const Valeurdsdp () const
 Returns $\partial / \partial \phi$ of *this. More...
 
const Valeurstdsdp () const
 Returns $1/\sin(\theta) \, \partial / \partial \phi$ of *this. More...
 
const Valeurd2sdp2 () const
 Returns $\partial^2 / \partial \phi^2$ of *this. More...
 
const Valeurmult_cp () const
 Returns $\cos(\phi) \, Id$ applied to *this. More...
 
const Valeurmult_sp () const
 Returns $\sin(\phi) \, Id$ applied to *this. More...
 
const Valeurlapang () const
 Returns the angular Laplacian of *this. More...
 
const Valeursx () const
 Returns ${1 \over \xi}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over \xi-1}$ (r -sampling = UNSURR ) More...
 
const Valeursx2 () const
 Returns ${1 \over \xi^2}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over (\xi-1)^2}$ (r -sampling = UNSURR ) More...
 
const Valeurmult_x () const
 Returns $\xi \, Id$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ $(\xi-1) \, Id $ (r -sampling = UNSURR ) More...
 
void sxm1_zec ()
 Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ ${1 \over (\xi-1)}$ (r -sampling = UNSURR ) More...
 
void mult_xm1_shell (int index)
 Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ $(\xi-1) \, Id$ (r -sampling = FIN ) \. More...
 
void mult_xp1_shell (int index)
 Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ $(\xi+1) \, Id$ (r -sampling = FIN ) \. More...
 
void mult_x_shell (int index)
 Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ $\xi \, Id$ (r -sampling = FIN ) \. More...
 
void mult_xm1_zec ()
 Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ $(\xi-1) \, Id$ (r -sampling = UNSURR ) More...
 
void mult2_xm1_zec ()
 Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ $(\xi-1)^2 \, Id$ (r -sampling = UNSURR ) More...
 
void va_x ()
 Returns ${\xi}$ (r -sampling = RARE ) \ (r sampling = FIN ) \ ${1 \over \xi-1}$ (r -sampling = UNSURR ) More...
 
void sauve (FILE *) const
 Save in a file. More...
 
void display_coef (double threshold=1.e-7, int precision=4, ostream &ostr=cout) const
 Displays the spectral coefficients and the associated basis functions. More...
 
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. More...
 
void set_etat_nondef ()
 Sets the logical state to ETATNONDEF (undefined). More...
 
void set_etat_zero ()
 Sets the logical state to ETATZERO (zero). More...
 
void set_etat_c_qcq ()
 Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c ). More...
 
void set_etat_cf_qcq ()
 Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_cf ). More...
 
void annule_hard ()
 Sets the Valeur to zero in a hard way. More...
 
void annule (int l)
 Sets the Valeur to zero in a given domain. More...
 
void annule (int l_min, int l_max)
 Sets the Valeur to zero in several domains. More...
 
int get_etat () const
 Returns the logical state. More...
 
const Mg3dget_mg () const
 Returns the Mg3d on which the this is defined. More...
 
void operator+= (const Valeur &)
 += Valeur More...
 
void operator-= (const Valeur &)
 -= Valeur More...
 
void operator*= (const Valeur &)
 *= Valeur More...
 
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). More...
 
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). More...
 
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. More...
 

Public Attributes

Mtblc
 Values of the function at the points of the multi-grid. More...
 
Mtbl_cfc_cf
 Coefficients of the spectral expansion of the function. More...
 
Base_val base
 Bases on which the spectral expansion is performed. More...
 

Private Member Functions

void nouveau ()
 Memory allocation. More...
 
void del_t ()
 Logical destructor. More...
 
void del_deriv ()
 Logical destructor of the derivatives. More...
 
void set_der_0x0 ()
 Sets the pointers for derivatives to 0x0. More...
 

Private Attributes

const Mg3dmg
 Multi-grid Mgd3 on which this is defined. More...
 
int etat
 Logical state (ETATNONDEF , ETATQCQ or ETATZERO ). More...
 
Valeurp_dsdx
 Pointer on $\partial / \partial \xi$. More...
 
Valeurp_d2sdx2
 Pointer on $\partial^2 / \partial \xi^2$. More...
 
Valeurp_sx
 Pointer on $1 / \xi$. More...
 
Valeurp_sx2
 Pointer on $1 / \xi^2$. More...
 
Valeurp_mult_x
 Pointer on $\xi \, Id$. More...
 
Valeurp_dsdt
 Pointer on $\partial / \partial \theta$. More...
 
Valeurp_d2sdt2
 Pointer on $\partial^2 / \partial \theta^2$. More...
 
Valeurp_ssint
 Pointer on $1 / \sin(\theta)$. More...
 
Valeurp_scost
 Pointer on $1 / \cos(\theta)$. More...
 
Valeurp_mult_ct
 Pointer on $\cos(\theta) \, Id$. More...
 
Valeurp_mult_st
 Pointer on $\sin(\theta) \, Id$. More...
 
Valeurp_dsdp
 Pointer on $\partial / \partial \phi$. More...
 
Valeurp_stdsdp
 Pointer on $1/\sin\theta \partial / \partial \phi$. More...
 
Valeurp_d2sdp2
 Pointer on $\partial^2 / \partial \phi^2$. More...
 
Valeurp_mult_cp
 Pointer on $\cos(\phi) \, Id$. More...
 
Valeurp_mult_sp
 Pointer on $\sin(\phi) \, Id$. More...
 
Valeurp_lapang
 Pointer on the angular Laplacian. More...
 

Friends

class Cmp
 Friend class. More...
 
class Scalar
 Friend class. More...
 
ostream & operator<< (ostream &, const Valeur &)
 Display. More...
 
void rotate_propre_pair (Valeur &, bool)
 Friend fonction. More...
 
void rotate_propre_impair (Valeur &, bool)
 Friend fonction. More...
 

Detailed Description

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

()

Definition at line 297 of file valeur.h.

Constructor & Destructor Documentation

◆ Valeur() [1/4]

Lorene::Valeur::Valeur ( const Mg3d mgrid)
explicit

Constructor.

Definition at line 203 of file valeur.C.

References nouveau().

◆ Valeur() [2/4]

Lorene::Valeur::Valeur ( const Mg3d p_mgrid)
explicit

Constructor.

Definition at line 209 of file valeur.C.

References nouveau().

◆ Valeur() [3/4]

Lorene::Valeur::Valeur ( const Mg3d g,
FILE *  fd 
)

Constructor from a file (see sauve(FILE*) )

Definition at line 266 of file valeur.C.

References etat, Lorene::fread_be(), and mg.

◆ Valeur() [4/4]

Lorene::Valeur::Valeur ( const Valeur vc)

Copy constructor.

Definition at line 216 of file valeur.C.

References get_etat(), get_mg(), mg, and nouveau().

◆ ~Valeur()

Lorene::Valeur::~Valeur ( )

Destructor.

Definition at line 300 of file valeur.C.

References del_t().

Member Function Documentation

◆ affiche_seuil()

void Lorene::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 561 of file valeur.C.

References etat.

◆ annule() [1/2]

void Lorene::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 747 of file valeur.C.

◆ annule() [2/2]

void Lorene::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 757 of file valeur.C.

References etat, Lorene::Mg3d::get_nzone(), mg, and set_etat_zero().

◆ annule_hard()

void Lorene::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 726 of file valeur.C.

References Lorene::Mtbl::annule_hard(), Lorene::Mtbl_cf::annule_hard(), base, c, c_cf, del_deriv(), etat, and mg.

◆ coef()

◆ coef_i()

◆ d2sdp2()

const Valeur & Lorene::Valeur::d2sdp2 ( ) const

Returns $\partial^2 / \partial \phi^2$ of *this.

Definition at line 95 of file valeur_d2sdp2.C.

References etat.

◆ d2sdt2()

const Valeur & Lorene::Valeur::d2sdt2 ( ) const

Returns $\partial^2 / \partial \theta^2$ of *this.

Definition at line 110 of file valeur_d2sdt2.C.

References etat.

◆ d2sdx2()

const Valeur & Lorene::Valeur::d2sdx2 ( ) const

Returns $\partial^2 / \partial \xi^2$ of *this.

Definition at line 117 of file valeur_d2sdx2.C.

References etat.

◆ del_deriv()

void Lorene::Valeur::del_deriv ( )
private

Logical destructor of the derivatives.

Definition at line 641 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().

◆ del_t()

void Lorene::Valeur::del_t ( )
private

Logical destructor.

Definition at line 629 of file valeur.C.

References c, c_cf, and del_deriv().

◆ display_coef()

void Lorene::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 539 of file valeur.C.

References etat.

◆ dsdp()

const Valeur & Lorene::Valeur::dsdp ( ) const

Returns $\partial / \partial \phi$ of *this.

Definition at line 101 of file valeur_dsdp.C.

References etat.

◆ dsdt()

const Valeur & Lorene::Valeur::dsdt ( ) const

Returns $\partial / \partial \theta$ of *this.

Definition at line 115 of file valeur_dsdt.C.

References etat.

◆ dsdx()

const Valeur & Lorene::Valeur::dsdx ( ) const

Returns $\partial / \partial \xi$ of *this.

Definition at line 114 of file valeur_dsdx.C.

References etat.

◆ equipot()

void Lorene::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 84 of file valeur_equipot.C.

References etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), and mg.

◆ equipot_outward()

void Lorene::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 86 of file valeur_equipot_out.C.

References etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), and mg.

◆ filtre_tp()

void Lorene::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 927 of file valeur.C.

References etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), and mg.

◆ get_base()

const Base_val& Lorene::Valeur::get_base ( ) const
inline

Return the bases for spectral expansions (member base )

Definition at line 490 of file valeur.h.

References base.

◆ get_etat()

int Lorene::Valeur::get_etat ( ) const
inline

Returns the logical state.

Definition at line 760 of file valeur.h.

References etat.

◆ get_mg()

const Mg3d* Lorene::Valeur::get_mg ( ) const
inline

Returns the Mg3d on which the this is defined.

Definition at line 763 of file valeur.h.

References mg.

◆ lapang()

const Valeur & Lorene::Valeur::lapang ( ) const

Returns the angular Laplacian of *this.

Definition at line 75 of file valeur_lapang.C.

References etat.

◆ mult2_xm1_zec()

void Lorene::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 79 of file valeur_mult2_xm1.C.

References etat.

◆ mult_cp()

const Valeur & Lorene::Valeur::mult_cp ( ) const

Returns $\cos(\phi) \, Id$ applied to *this.

Definition at line 80 of file valeur_mult_cp.C.

References etat.

◆ mult_ct()

const Valeur & Lorene::Valeur::mult_ct ( ) const

Returns $\cos(\theta) \, Id$ applied to *this.

Definition at line 101 of file valeur_mult_ct.C.

References etat.

◆ mult_sp()

const Valeur & Lorene::Valeur::mult_sp ( ) const

Returns $\sin(\phi) \, Id$ applied to *this.

Definition at line 80 of file valeur_mult_sp.C.

References etat.

◆ mult_st()

const Valeur & Lorene::Valeur::mult_st ( ) const

Returns $\sin(\theta) \, Id$ applied to *this.

Definition at line 100 of file valeur_mult_st.C.

References etat.

◆ mult_x()

const Valeur & Lorene::Valeur::mult_x ( ) const

Returns $\xi \, Id$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ $(\xi-1) \, Id $ (r -sampling = UNSURR )

Definition at line 108 of file valeur_mult_x.C.

References etat.

◆ mult_x_shell()

void Lorene::Valeur::mult_x_shell ( int  index)

Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ $\xi \, Id$ (r -sampling = FIN ) \.

Parameters
index[input] Index of the domain where computation is performed.

Definition at line 192 of file valeur_mult_x.C.

References etat.

◆ mult_xm1_shell()

void Lorene::Valeur::mult_xm1_shell ( int  index)

Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ $(\xi-1) \, Id$ (r -sampling = FIN ) \.

Parameters
index[input] Index of the domain where computation is performed.

Definition at line 83 of file valeur_mult_xm1.C.

References etat.

◆ mult_xm1_zec()

void Lorene::Valeur::mult_xm1_zec ( )

Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ $(\xi-1) \, Id$ (r -sampling = UNSURR )

Definition at line 129 of file valeur_mult_xm1.C.

References etat.

◆ mult_xp1_shell()

void Lorene::Valeur::mult_xp1_shell ( int  index)

Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ $(\xi+1) \, Id$ (r -sampling = FIN ) \.

Parameters
index[input] Index of the domain where computation is performed.

Definition at line 45 of file valeur_mult_xp1.C.

References etat.

◆ nouveau()

void Lorene::Valeur::nouveau ( )
private

Memory allocation.

Definition at line 620 of file valeur.C.

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

◆ operator()() [1/2]

const Tbl& Lorene::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 391 of file valeur.h.

References etat.

◆ operator()() [2/2]

double Lorene::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 431 of file valeur.h.

References etat.

◆ operator*=()

void Lorene::Valeur::operator*= ( const Valeur vi)

*= Valeur

Definition at line 964 of file valeur_arithm.C.

References etat, get_mg(), and mg.

◆ operator+=()

void Lorene::Valeur::operator+= ( const Valeur vi)

+= Valeur

Definition at line 840 of file valeur_arithm.C.

References etat, get_mg(), and mg.

◆ operator-=()

void Lorene::Valeur::operator-= ( const Valeur vi)

-= Valeur

Definition at line 903 of file valeur_arithm.C.

References etat, get_mg(), and mg.

◆ operator=() [1/4]

void Lorene::Valeur::operator= ( const Valeur a)

Assignement to another Valeur.

Definition at line 330 of file valeur.C.

References base, c, c_cf, etat, and mg.

◆ operator=() [2/4]

void Lorene::Valeur::operator= ( const Mtbl mt)

Assignement to a Mtbl.

Definition at line 390 of file valeur.C.

References Lorene::Mtbl::get_etat().

◆ operator=() [3/4]

void Lorene::Valeur::operator= ( const Mtbl_cf mtcf)

Assignement to a Mtbl_cf.

Definition at line 433 of file valeur.C.

References Lorene::Mtbl_cf::get_etat().

◆ operator=() [4/4]

void Lorene::Valeur::operator= ( double  x)

Assignement to a double.

Definition at line 310 of file valeur.C.

References base, c, Lorene::Base_val::set_base_nondef(), set_etat_c_qcq(), and set_etat_zero().

◆ sauve()

void Lorene::Valeur::sauve ( FILE *  fd) const

Save in a file.

Definition at line 478 of file valeur.C.

References base, etat, Lorene::fwrite_be(), mg, Lorene::Base_val::sauve(), and Lorene::Mg3d::sauve().

◆ scost()

const Valeur & Lorene::Valeur::scost ( ) const

Returns $1 / \cos(\theta)$ of *this.

Definition at line 102 of file valeur_scost.C.

References etat.

◆ set() [1/2]

Tbl& Lorene::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 373 of file valeur.h.

References etat.

◆ set() [2/2]

double& Lorene::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 411 of file valeur.h.

References etat.

◆ set_base()

void Lorene::Valeur::set_base ( const Base_val base_i)

Sets the bases for spectral expansions (member base )

Definition at line 813 of file valeur.C.

References Lorene::Mtbl_cf::base, base, and c_cf.

◆ set_base_p()

void Lorene::Valeur::set_base_p ( int  base_p)

Sets the expansion basis for $\phi$ functions in all domains.

Parameters
base_ptype 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 865 of file valeur.C.

References Lorene::Mtbl_cf::base, base, c_cf, and Lorene::Base_val::set_base_p().

◆ set_base_r()

void Lorene::Valeur::set_base_r ( int  l,
int  base_r 
)

Sets the expansion basis for r ( $\xi$) functions in a given domain.

Parameters
lDomain index
base_rtype 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 839 of file valeur.C.

References Lorene::Mtbl_cf::base, base, c_cf, and Lorene::Base_val::set_base_r().

◆ set_base_t()

void Lorene::Valeur::set_base_t ( int  base_t)

Sets the expansion basis for $\theta$ functions in all domains.

Parameters
base_ttype 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 852 of file valeur.C.

References Lorene::Mtbl_cf::base, base, c_cf, and Lorene::Base_val::set_base_t().

◆ set_der_0x0()

void Lorene::Valeur::set_der_0x0 ( )
private

Sets the pointers for derivatives to 0x0.

Definition at line 668 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.

◆ set_etat_c_qcq()

void Lorene::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 704 of file valeur.C.

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

◆ set_etat_cf_qcq()

void Lorene::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 715 of file valeur.C.

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

◆ set_etat_nondef()

void Lorene::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 698 of file valeur.C.

References etat.

◆ set_etat_zero()

void Lorene::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 692 of file valeur.C.

References etat.

◆ smooth()

void Lorene::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 80 of file valeur_smooth.C.

References etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), and mg.

◆ ssint()

const Valeur & Lorene::Valeur::ssint ( ) const

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

Definition at line 115 of file valeur_ssint.C.

References etat.

◆ std_base_scal()

void Lorene::Valeur::std_base_scal ( )

Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.

Definition at line 827 of file valeur.C.

References mg, set_base(), and Lorene::Mg3d::std_base_scal().

◆ std_base_scal_odd()

void Lorene::Valeur::std_base_scal_odd ( )

Sets the bases for spectral expansions (member base ) to the standard odd ones for a scalar.

Definition at line 833 of file valeur.C.

References mg, set_base(), and Lorene::Mg3d::std_base_scal_odd().

◆ stdsdp()

const Valeur & Lorene::Valeur::stdsdp ( ) const

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

Definition at line 63 of file valeur_stdsdp.C.

References etat.

◆ sx()

const Valeur & Lorene::Valeur::sx ( ) const

Returns ${1 \over \xi}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over \xi-1}$ (r -sampling = UNSURR )

Definition at line 113 of file valeur_sx.C.

References etat.

◆ sx2()

const Valeur & Lorene::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 117 of file valeur_sx2.C.

References etat.

◆ sxm1_zec()

void Lorene::Valeur::sxm1_zec ( )

Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ ${1 \over (\xi-1)}$ (r -sampling = UNSURR )

Definition at line 84 of file valeur_sxm1_zec.C.

References etat.

◆ va_x()

void Lorene::Valeur::va_x ( )

Returns ${\xi}$ (r -sampling = RARE ) \ (r sampling = FIN ) \ ${1 \over \xi-1}$ (r -sampling = UNSURR )

Definition at line 43 of file valeur_x.C.

References etat.

◆ val_point()

double Lorene::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 885 of file valeur.C.

References etat.

◆ val_point_jk()

double Lorene::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 903 of file valeur.C.

References etat.

◆ val_propre_1d()

void Lorene::Valeur::val_propre_1d ( )

◆ val_propre_1d_i()

void Lorene::Valeur::val_propre_1d_i ( )

◆ ylm()

void Lorene::Valeur::ylm ( )

◆ ylm_i()

Friends And Related Function Documentation

◆ Cmp

friend class Cmp
friend

Friend class.

Definition at line 859 of file valeur.h.

◆ operator<<

ostream& operator<< ( ostream &  o,
const Valeur vi 
)
friend

Display.

Definition at line 499 of file valeur.C.

◆ rotate_propre_impair

void rotate_propre_impair ( Valeur so,
bool  sens 
)
friend

Friend fonction.

Definition at line 199 of file valeur_val_propre_1d.C.

◆ rotate_propre_pair

void rotate_propre_pair ( Valeur so,
bool  sens 
)
friend

Friend fonction.

Definition at line 69 of file valeur_val_propre_1d.C.

◆ Scalar

friend class Scalar
friend

Friend class.

Definition at line 860 of file valeur.h.

Member Data Documentation

◆ base

Base_val Lorene::Valeur::base

Bases on which the spectral expansion is performed.

Definition at line 315 of file valeur.h.

◆ c

Mtbl* Lorene::Valeur::c
mutable

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

Definition at line 309 of file valeur.h.

◆ c_cf

Mtbl_cf* Lorene::Valeur::c_cf
mutable

Coefficients of the spectral expansion of the function.

Definition at line 312 of file valeur.h.

◆ etat

int Lorene::Valeur::etat
private

Logical state (ETATNONDEF , ETATQCQ or ETATZERO ).

Definition at line 305 of file valeur.h.

◆ mg

const Mg3d* Lorene::Valeur::mg
private

Multi-grid Mgd3 on which this is defined.

Definition at line 302 of file valeur.h.

◆ p_d2sdp2

Valeur* Lorene::Valeur::p_d2sdp2
mutableprivate

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

Definition at line 335 of file valeur.h.

◆ p_d2sdt2

Valeur* Lorene::Valeur::p_d2sdt2
mutableprivate

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

Definition at line 327 of file valeur.h.

◆ p_d2sdx2

Valeur* Lorene::Valeur::p_d2sdx2
mutableprivate

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

Definition at line 321 of file valeur.h.

◆ p_dsdp

Valeur* Lorene::Valeur::p_dsdp
mutableprivate

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

Definition at line 333 of file valeur.h.

◆ p_dsdt

Valeur* Lorene::Valeur::p_dsdt
mutableprivate

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

Definition at line 326 of file valeur.h.

◆ p_dsdx

Valeur* Lorene::Valeur::p_dsdx
mutableprivate

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

Definition at line 320 of file valeur.h.

◆ p_lapang

Valeur* Lorene::Valeur::p_lapang
mutableprivate

Pointer on the angular Laplacian.

Definition at line 339 of file valeur.h.

◆ p_mult_cp

Valeur* Lorene::Valeur::p_mult_cp
mutableprivate

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

Definition at line 336 of file valeur.h.

◆ p_mult_ct

Valeur* Lorene::Valeur::p_mult_ct
mutableprivate

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

Definition at line 330 of file valeur.h.

◆ p_mult_sp

Valeur* Lorene::Valeur::p_mult_sp
mutableprivate

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

Definition at line 337 of file valeur.h.

◆ p_mult_st

Valeur* Lorene::Valeur::p_mult_st
mutableprivate

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

Definition at line 331 of file valeur.h.

◆ p_mult_x

Valeur* Lorene::Valeur::p_mult_x
mutableprivate

Pointer on $\xi \, Id$.

Definition at line 324 of file valeur.h.

◆ p_scost

Valeur* Lorene::Valeur::p_scost
mutableprivate

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

Definition at line 329 of file valeur.h.

◆ p_ssint

Valeur* Lorene::Valeur::p_ssint
mutableprivate

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

Definition at line 328 of file valeur.h.

◆ p_stdsdp

Valeur* Lorene::Valeur::p_stdsdp
mutableprivate

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

Definition at line 334 of file valeur.h.

◆ p_sx

Valeur* Lorene::Valeur::p_sx
mutableprivate

Pointer on $1 / \xi$.

Definition at line 322 of file valeur.h.

◆ p_sx2

Valeur* Lorene::Valeur::p_sx2
mutableprivate

Pointer on $1 / \xi^2$.

Definition at line 323 of file valeur.h.


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