|
LORENE
|
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... | |
| Tbl & | set (int l) |
| Read/write of the value in a given domain (configuration space). More... | |
| const Tbl & | operator() (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 , but collocation point in , 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 of *this. More... | |
| void | ylm_i () |
Inverse of ylm() More... | |
| void | val_propre_1d () |
Set the basis to the eigenvalues of . More... | |
| void | val_propre_1d_i () |
Inverse transformation of val_propre_1d. More... | |
| const Base_val & | get_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 ( ) functions in a given domain. More... | |
| void | set_base_t (int base_t) |
Sets the expansion basis for functions in all domains. More... | |
| void | set_base_p (int base_p) |
Sets the expansion basis for functions in all domains. More... | |
| 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. More... | |
| const Valeur & | dsdx () const |
Returns of *this. More... | |
| const Valeur & | d2sdx2 () const |
Returns of *this. More... | |
| const Valeur & | dsdt () const |
Returns of *this. More... | |
| const Valeur & | d2sdt2 () const |
Returns of *this. More... | |
| const Valeur & | ssint () const |
Returns of *this. More... | |
| const Valeur & | scost () const |
Returns of *this. More... | |
| const Valeur & | mult_ct () const |
Returns applied to *this. More... | |
| const Valeur & | mult_st () const |
Returns applied to *this. More... | |
| const Valeur & | dsdp () const |
Returns of *this. More... | |
| const Valeur & | stdsdp () const |
Returns of *this. More... | |
| const Valeur & | d2sdp2 () const |
Returns of *this. More... | |
| const Valeur & | mult_cp () const |
Returns applied to *this. More... | |
| const Valeur & | mult_sp () const |
Returns applied to *this. More... | |
| const Valeur & | lapang () const |
Returns the angular Laplacian of *this. More... | |
| const Valeur & | sx () const |
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) More... | |
| const Valeur & | sx2 () const |
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) More... | |
| const Valeur & | mult_x () const |
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) More... | |
| void | sxm1_zec () |
Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ (r -sampling = UNSURR ) More... | |
| void | mult_xm1_shell (int index) |
Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ (r -sampling = FIN ) \. More... | |
| void | mult_xp1_shell (int index) |
Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ (r -sampling = FIN ) \. More... | |
| void | mult_x_shell (int index) |
Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \ (r -sampling = FIN ) \. More... | |
| void | mult_xm1_zec () |
Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ (r -sampling = UNSURR ) More... | |
| void | mult2_xm1_zec () |
Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \ (r -sampling = UNSURR ) More... | |
| void | va_x () |
Returns (r -sampling = RARE ) \ (r sampling = FIN ) \ (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 Mg3d * | get_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 | |
| Mtbl * | c |
| Values of the function at the points of the multi-grid. More... | |
| Mtbl_cf * | c_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 Mg3d * | mg |
Multi-grid Mgd3 on which this is defined. More... | |
| int | etat |
Logical state (ETATNONDEF , ETATQCQ or ETATZERO ). More... | |
| Valeur * | p_dsdx |
Pointer on . More... | |
| Valeur * | p_d2sdx2 |
Pointer on . More... | |
| Valeur * | p_sx |
Pointer on . More... | |
| Valeur * | p_sx2 |
Pointer on . More... | |
| Valeur * | p_mult_x |
Pointer on . More... | |
| Valeur * | p_dsdt |
Pointer on . More... | |
| Valeur * | p_d2sdt2 |
Pointer on . More... | |
| Valeur * | p_ssint |
Pointer on . More... | |
| Valeur * | p_scost |
Pointer on . More... | |
| Valeur * | p_mult_ct |
Pointer on . More... | |
| Valeur * | p_mult_st |
Pointer on . More... | |
| Valeur * | p_dsdp |
Pointer on . More... | |
| Valeur * | p_stdsdp |
Pointer on . More... | |
| Valeur * | p_d2sdp2 |
Pointer on . More... | |
| Valeur * | p_mult_cp |
Pointer on . More... | |
| Valeur * | p_mult_sp |
Pointer on . More... | |
| Valeur * | p_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... | |
|
explicit |
|
explicit |
| 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.
| Lorene::Valeur::Valeur | ( | const Valeur & | vc | ) |
| Lorene::Valeur::~Valeur | ( | ) |
| 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.
| 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.
| void Lorene::Valeur::annule | ( | int | l | ) |
| void Lorene::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 757 of file valeur.C.
References etat, Lorene::Mg3d::get_nzone(), mg, and set_etat_zero().
| 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.
| void Lorene::Valeur::coef | ( | ) | const |
Computes the coeffcients of *this.
Definition at line 151 of file valeur_coef.C.
References etat, MAX_BASE, MAX_BASE_2, 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, 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 Lorene::Valeur::coef_i | ( | ) | const |
Computes the physical value of *this.
Definition at line 143 of file valeur_coef_i.C.
References etat, MAX_BASE, MAX_BASE_2, 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, 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 & Lorene::Valeur::d2sdp2 | ( | ) | const |
| const Valeur & Lorene::Valeur::d2sdt2 | ( | ) | const |
| const Valeur & Lorene::Valeur::d2sdx2 | ( | ) | const |
|
private |
|
private |
| 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.
| 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.
| const Valeur & Lorene::Valeur::dsdp | ( | ) | const |
| const Valeur & Lorene::Valeur::dsdt | ( | ) | const |
| const Valeur & Lorene::Valeur::dsdx | ( | ) | const |
| 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 \
\
\ 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 84 of file valeur_equipot.C.
References etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), and mg.
| 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 \
\
\ 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 86 of file valeur_equipot_out.C.
References etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), and mg.
| void Lorene::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 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.
|
inline |
|
inline |
|
inline |
| const Valeur & Lorene::Valeur::lapang | ( | ) | const |
Returns the angular Laplacian of *this.
Definition at line 75 of file valeur_lapang.C.
References etat.
| void Lorene::Valeur::mult2_xm1_zec | ( | ) |
Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \
(r -sampling = UNSURR )
Definition at line 79 of file valeur_mult2_xm1.C.
References etat.
| const Valeur & Lorene::Valeur::mult_cp | ( | ) | const |
| const Valeur & Lorene::Valeur::mult_ct | ( | ) | const |
| const Valeur & Lorene::Valeur::mult_sp | ( | ) | const |
| const Valeur & Lorene::Valeur::mult_st | ( | ) | const |
| const Valeur & Lorene::Valeur::mult_x | ( | ) | const |
Returns
(r -sampling = RARE ) \ Id (r sampling = FIN ) \
(r -sampling = UNSURR )
Definition at line 108 of file valeur_mult_x.C.
References etat.
| void Lorene::Valeur::mult_x_shell | ( | int | index | ) |
Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \
(r -sampling = FIN ) \.
| index | [input] Index of the domain where computation is performed. |
Definition at line 192 of file valeur_mult_x.C.
References etat.
| void Lorene::Valeur::mult_xm1_shell | ( | int | index | ) |
Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \
(r -sampling = FIN ) \.
| index | [input] Index of the domain where computation is performed. |
Definition at line 83 of file valeur_mult_xm1.C.
References etat.
| void Lorene::Valeur::mult_xm1_zec | ( | ) |
Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \
(r -sampling = UNSURR )
Definition at line 129 of file valeur_mult_xm1.C.
References etat.
| void Lorene::Valeur::mult_xp1_shell | ( | int | index | ) |
Applies the following operator to *this : \ Id (r sampling = RARE, UNSURR ) \
(r -sampling = FIN ) \.
| index | [input] Index of the domain where computation is performed. |
Definition at line 45 of file valeur_mult_xp1.C.
References etat.
|
private |
|
inline |
|
inline |
| void Lorene::Valeur::operator*= | ( | const Valeur & | vi | ) |
| void Lorene::Valeur::operator+= | ( | const Valeur & | vi | ) |
| void Lorene::Valeur::operator-= | ( | const Valeur & | vi | ) |
| void Lorene::Valeur::operator= | ( | const Valeur & | a | ) |
| void Lorene::Valeur::operator= | ( | const Mtbl & | mt | ) |
Assignement to a Mtbl.
Definition at line 390 of file valeur.C.
References Lorene::Mtbl::get_etat().
| 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().
| 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().
| 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().
| const Valeur & Lorene::Valeur::scost | ( | ) | const |
|
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 373 of file valeur.h.
References etat.
|
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 411 of file valeur.h.
References etat.
| 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.
| void Lorene::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 865 of file valeur.C.
References Lorene::Mtbl_cf::base, base, c_cf, and Lorene::Base_val::set_base_p().
| void Lorene::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 839 of file valeur.C.
References Lorene::Mtbl_cf::base, base, c_cf, and Lorene::Base_val::set_base_r().
| void Lorene::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 852 of file valeur.C.
References Lorene::Mtbl_cf::base, base, c_cf, and Lorene::Base_val::set_base_t().
|
private |
| 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.
| 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.
| void Lorene::Valeur::set_etat_nondef | ( | ) |
| void Lorene::Valeur::set_etat_zero | ( | ) |
| 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.
| 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.
| const Valeur & Lorene::Valeur::ssint | ( | ) | const |
| 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().
| 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().
| const Valeur & Lorene::Valeur::stdsdp | ( | ) | const |
| const Valeur & Lorene::Valeur::sx | ( | ) | const |
Returns
(r -sampling = RARE ) \ Id (r sampling = FIN ) \
(r -sampling = UNSURR )
Definition at line 113 of file valeur_sx.C.
References etat.
| const Valeur & Lorene::Valeur::sx2 | ( | ) | const |
Returns
(r -sampling = RARE ) \ Id (r sampling = FIN ) \
(r -sampling = UNSURR )
Definition at line 117 of file valeur_sx2.C.
References etat.
| void Lorene::Valeur::sxm1_zec | ( | ) |
Applies the following operator to *this : \ Id (r sampling = RARE, FIN ) \
(r -sampling = UNSURR )
Definition at line 84 of file valeur_sxm1_zec.C.
References etat.
| void Lorene::Valeur::va_x | ( | ) |
Returns
(r -sampling = RARE ) \ (r sampling = FIN ) \
(r -sampling = UNSURR )
Definition at line 43 of file valeur_x.C.
References etat.
| 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.
| l | [input] index of the domain |
| x | [input] value of the coordinate |
| theta | [input] value of the coordinate |
| phi | [input] value of the coordinate |
in the domain no. l of the field represented by *this . Definition at line 885 of file valeur.C.
References etat.
| 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
, 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 |
in the domain no. l of the field represented by *this . Definition at line 903 of file valeur.C.
References etat.
| void Lorene::Valeur::val_propre_1d | ( | ) |
Set the basis to the eigenvalues of
.
Definition at line 325 of file valeur_val_propre_1d.C.
References base, Lorene::Base_val::get_base_t(), Lorene::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 Lorene::Valeur::val_propre_1d_i | ( | ) |
Inverse transformation of val_propre_1d.
Definition at line 361 of file valeur_val_propre_1d.C.
References base, Lorene::Base_val::get_base_t(), Lorene::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 Lorene::Valeur::ylm | ( | ) |
Computes the coefficients
of *this.
Definition at line 141 of file valeur_ylm.C.
References etat, get_mg(), Lorene::Mg3d::get_nzone(), MAX_BASE, NONDEF, 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 Lorene::Valeur::ylm_i | ( | ) |
Inverse of ylm()
Definition at line 134 of file valeur_ylm_i.C.
References etat, get_mg(), Lorene::Mg3d::get_nzone(), MAX_BASE, NONDEF, 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.
|
friend |
|
friend |
Friend fonction.
Definition at line 199 of file valeur_val_propre_1d.C.
|
friend |
Friend fonction.
Definition at line 69 of file valeur_val_propre_1d.C.
| Base_val Lorene::Valeur::base |
|
mutable |
|
mutable |
|
private |
|
private |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |