Coefficients storage for the multi-domain spectral method. More...
#include <mtbl_cf.h>
Public Member Functions | |
Mtbl_cf (const Mg3d &mgrid, const Base_val &basis) | |
Constructor. | |
Mtbl_cf (const Mg3d *p_mgrid, const Base_val &basis) | |
Constructor. | |
Mtbl_cf (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE*) ). | |
Mtbl_cf (const Mtbl_cf &) | |
Copy constructor. | |
~Mtbl_cf () | |
Destructor. | |
void | operator= (const Mtbl_cf &) |
Assignement to another Mtbl_cf . | |
void | operator= (double) |
Assignement to a double . | |
void | operator= (int) |
Assignement to a int . | |
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_qcq () |
Sets the logical state to ETATQCQ (ordinary state). | |
void | annule_hard () |
Sets the Mtbl_cf to zero in a hard way. | |
void | annule (int l_min, int l_max) |
Sets the Mtbl_cf to zero in some domains. | |
Tbl & | set (int l) |
Read/write of the Tbl containing the coefficients in a given domain. | |
const Tbl & | operator() (int l) const |
Read-only of the Tbl containing the coefficients in a given domain. | |
double & | set (int l, int k, int j, int i) |
Read/write of a particular element. | |
double | operator() (int l, int k, int j, int i) const |
Read-only of a particular element. | |
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_symy (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_asymy (int l, double x, double theta, double phi) const |
Computes the value of the field represented by *this at an arbitrary point, by means of the spectral expansion. | |
double | val_point_jk (int l, double x, int j, int k) const |
Computes the value of the field represented by *this at an arbitrary point in , but collocation point in , by means of the spectral expansion. | |
double | val_point_jk_symy (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. | |
double | val_point_jk_asymy (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. | |
double | val_out_bound_jk (int l, int j, int k) const |
Computes the angular coefficient of index j ,k of the field represented by *this at by means of the spectral expansion. | |
double | val_in_bound_jk (int l, int j, int k) const |
Computes the angular coefficient of index j ,k of the field represented by *this at by means of the spectral expansion. | |
const Mg3d * | get_mg () const |
Returns the Mg3d on which the Mtbl_cf is defined. | |
int | get_etat () const |
Returns the logical state. | |
int | get_nzone () const |
Returns the number of zones (domains). | |
void | sauve (FILE *) const |
Save in a file. | |
void | display (double threshold=1.e-7, int precision=4, ostream &ostr=cout) const |
Prints the coefficients whose values are greater than a given threshold, as well as the corresponding basis. | |
void | affiche_seuil (ostream &ostr, int precision=4, double threshold=1.e-7) const |
Prints only the values greater than a given threshold. | |
void | operator+= (const Mtbl_cf &) |
+= Mtbl_cf | |
void | operator-= (const Mtbl_cf &) |
-= Mtbl_cf | |
void | operator*= (double) |
*= double | |
void | operator/= (double) |
/= double | |
void | dsdx () |
| |
void | d2sdx2 () |
| |
void | sx () |
(r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) | |
void | sx2 () |
(r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) | |
void | mult_x () |
(r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) | |
void | sxm1_zec () |
Id (r sampling = RARE , FIN ) \ (r -sampling = UNSURR ). | |
void | mult_xm1_zec () |
Id (r sampling = RARE , FIN ) \ (r -sampling = UNSURR ). | |
void | mult2_xm1_zec () |
Id (r sampling = RARE , FIN ) \ (r -sampling = UNSURR ). | |
void | dsdt () |
| |
void | d2sdt2 () |
| |
void | ssint () |
| |
void | scost () |
| |
void | mult_ct () |
| |
void | mult_st () |
| |
void | dsdp () |
| |
void | d2sdp2 () |
| |
void | mult_cp () |
| |
void | mult_sp () |
| |
void | lapang () |
Angular Laplacian. | |
void | poisson_angu (double lambda=0) |
Resolution of the generalized angular Poisson equation. | |
Public Attributes | |
Base_val | base |
Bases of the spectral expansions. | |
Tbl ** | t |
Array (size nzone ) of pointers on the Tbl 's which contain the spectral coefficients in each domain. | |
Private Member Functions | |
void | del_t () |
Logical destructor: dellocates the memory occupied by the Tbl array t . | |
Private Attributes | |
const Mg3d * | mg |
Pointer on the multi-grid Mgd3 on which this is defined. | |
int | nzone |
Number of domains (zones). | |
int | etat |
Logical state (ETATNONDEF , ETATQCQ or ETATZERO ). | |
Friends | |
ostream & | operator<< (ostream &, const Mtbl_cf &) |
Display. |
Coefficients storage for the multi-domain spectral method.
This class is essentially an array (on the various physical domains) of Tbl
specially designed for storage of the coefficients of the spectral expansions in each domain. It is intended to be used in conjunction with the class Mtbl
(see class Valeur
). A difference between a Mtbl
and a Mtbl_cf
, both defined one the same grid Mg3d
, is that each Tbl
of the Mtbl_cf
has 2 more elements in the -dimension (Dim_tbl::dim[2]) than the corresponding Tbl
of the Mtbl
. A Mbl_cf
is initialy created with a logical state ETATZERO
. Arithmetic operations are provided with the usual meaning (see below).
()
Definition at line 182 of file mtbl_cf.h.
Constructor.
Definition at line 120 of file mtbl_cf.C.
References base, Base_val::get_nzone(), Mg3d::get_nzone(), and nzone.
Constructor.
Definition at line 128 of file mtbl_cf.C.
References base, Base_val::get_nzone(), Mg3d::get_nzone(), and nzone.
Mtbl_cf::Mtbl_cf | ( | const Mg3d & | g, | |
FILE * | fd | |||
) |
Constructor from a file (see sauve(FILE*)
).
Definition at line 169 of file mtbl_cf.C.
References base, etat, fread_be(), Base_val::get_nzone(), Mg3d::get_nzone(), mg, nzone, and t.
Mtbl_cf::Mtbl_cf | ( | const Mtbl_cf & | mtc | ) |
Copy constructor.
Definition at line 146 of file mtbl_cf.C.
References etat, get_etat(), nzone, set_etat_qcq(), and t.
void Mtbl_cf::affiche_seuil | ( | ostream & | ostr, | |
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 | |
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 393 of file mtbl_cf.C.
References Tbl::affiche_seuil(), base, etat, and t.
void Mtbl_cf::annule | ( | int | l_min, | |
int | l_max | |||
) |
Sets the Mtbl_cf
to zero in some domains.
l_min | [input] The Mtbl_cf 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, nzone-1) is equivalent to set_etat_zero()
.
Definition at line 328 of file mtbl_cf.C.
References etat, nzone, Tbl::set_etat_zero(), set_etat_zero(), and t.
void Mtbl_cf::annule_hard | ( | ) |
Sets the Mtbl_cf
to zero in a hard way.
1/ Sets the logical state to ETATQCQ
, i.e. to an ordinary state. 2/ Allocates the memory of the Tbl
array t
, 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 308 of file mtbl_cf.C.
References Tbl::annule_hard(), etat, Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), mg, nzone, and t.
void Mtbl_cf::d2sdp2 | ( | ) |
void Mtbl_cf::d2sdt2 | ( | ) |
Definition at line 144 of file valeur_d2sdt2.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_T, nzone, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
void Mtbl_cf::d2sdx2 | ( | ) |
Definition at line 151 of file valeur_d2sdx2.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_R, nzone, 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, and TRA_R.
void Mtbl_cf::del_t | ( | ) | [private] |
void Mtbl_cf::display | ( | double | threshold = 1.e-7 , |
|
int | precision = 4 , |
|||
ostream & | ostr = cout | |||
) | const |
Prints the coefficients whose values are greater than a given threshold, as well as the corresponding basis.
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 52 of file mtbl_cf_display.C.
References base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), mg, Base_val::name_phi(), Base_val::name_r(), Base_val::name_theta(), nzone, and t.
void Mtbl_cf::dsdp | ( | ) |
Definition at line 135 of file valeur_dsdp.C.
References Base_val::b, base, etat, MAX_BASE_2, MSQ_P, nzone, P_COSSIN, P_COSSIN_I, P_COSSIN_P, t, and TRA_P.
void Mtbl_cf::dsdt | ( | ) |
Definition at line 148 of file valeur_dsdt.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_T, nzone, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
void Mtbl_cf::dsdx | ( | ) |
Definition at line 149 of file valeur_dsdx.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_R, nzone, 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, and TRA_R.
int Mtbl_cf::get_etat | ( | ) | const [inline] |
const Mg3d* Mtbl_cf::get_mg | ( | ) | const [inline] |
int Mtbl_cf::get_nzone | ( | ) | const [inline] |
void Mtbl_cf::lapang | ( | ) |
void Mtbl_cf::mult2_xm1_zec | ( | ) |
void Mtbl_cf::mult_cp | ( | ) |
Definition at line 116 of file valeur_mult_cp.C.
References Base_val::b, base, etat, MAX_BASE_2, MSQ_P, nzone, P_COSSIN, P_COSSIN_I, P_COSSIN_P, t, and TRA_P.
void Mtbl_cf::mult_ct | ( | ) |
Definition at line 137 of file valeur_mult_ct.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_T, nzone, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
void Mtbl_cf::mult_sp | ( | ) |
Definition at line 116 of file valeur_mult_sp.C.
References Base_val::b, base, etat, MAX_BASE_2, MSQ_P, nzone, P_COSSIN, P_COSSIN_I, P_COSSIN_P, t, and TRA_P.
void Mtbl_cf::mult_st | ( | ) |
Definition at line 136 of file valeur_mult_st.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_T, nzone, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
void Mtbl_cf::mult_x | ( | ) |
(r -sampling = RARE
) \ Id (r sampling = FIN
) \ (r -sampling = UNSURR
)
Definition at line 135 of file valeur_mult_x.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_R, nzone, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, t, and TRA_R.
void Mtbl_cf::mult_xm1_zec | ( | ) |
double Mtbl_cf::operator() | ( | int | l, | |
int | k, | |||
int | j, | |||
int | i | |||
) | const [inline] |
const Tbl& Mtbl_cf::operator() | ( | int | l | ) | const [inline] |
void Mtbl_cf::operator*= | ( | double | ) |
void Mtbl_cf::operator+= | ( | const Mtbl_cf & | mi | ) |
+= Mtbl_cf
Definition at line 298 of file mtbl_cf_arithm.C.
References annule_hard(), base, etat, get_etat(), get_mg(), mg, nzone, and t.
void Mtbl_cf::operator-= | ( | const Mtbl_cf & | mi | ) |
-= Mtbl_cf
Definition at line 321 of file mtbl_cf_arithm.C.
References annule_hard(), base, etat, get_etat(), get_mg(), mg, nzone, and t.
void Mtbl_cf::operator/= | ( | double | ) |
void Mtbl_cf::operator= | ( | int | m | ) |
Assignement to a int
.
Definition at line 252 of file mtbl_cf.C.
References nzone, set_etat_qcq(), set_etat_zero(), and t.
void Mtbl_cf::operator= | ( | double | x | ) |
Assignement to a double
.
Definition at line 238 of file mtbl_cf.C.
References nzone, set_etat_qcq(), set_etat_zero(), and t.
void Mtbl_cf::operator= | ( | const Mtbl_cf & | mtc | ) |
Assignement to another Mtbl_cf
.
Definition at line 216 of file mtbl_cf.C.
References base, get_etat(), mg, nzone, set_etat_qcq(), set_etat_zero(), and t.
void Mtbl_cf::poisson_angu | ( | double | lambda = 0 |
) |
Resolution of the generalized angular Poisson equation.
The generalized angular Poisson equation is , where .
Before the call to poisson_angu()
, *this
contains the coefficients of the source ; after the call, it contains the coefficients of the solution .
lambda | [input] coefficient in the above equation (default value = 0) |
Definition at line 79 of file mtbl_cf_pde.C.
References Base_val::b, base, get_mg(), Mg3d::get_nzone(), MAX_BASE, MSQ_T, T_LEG, T_LEG_I, T_LEG_IP, T_LEG_MI, T_LEG_MP, T_LEG_P, T_LEG_PI, T_LEG_PP, and TRA_T.
void Mtbl_cf::sauve | ( | FILE * | fd | ) | const |
Save in a file.
Definition at line 200 of file mtbl_cf.C.
References base, etat, fwrite_be(), mg, nzone, Tbl::sauve(), Mg3d::sauve(), Base_val::sauve(), and t.
void Mtbl_cf::scost | ( | ) |
Definition at line 138 of file valeur_scost.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_T, nzone, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
double& Mtbl_cf::set | ( | int | l, | |
int | k, | |||
int | j, | |||
int | i | |||
) | [inline] |
Tbl& Mtbl_cf::set | ( | int | l | ) | [inline] |
void Mtbl_cf::set_etat_nondef | ( | ) |
void Mtbl_cf::set_etat_qcq | ( | ) |
Sets the logical state to ETATQCQ
(ordinary state).
If the state (member etat
) is already ETATQCQ
, this function does nothing. Otherwise, it performs the memory allocation for the Tbl
array t
.
Definition at line 296 of file mtbl_cf.C.
References etat, Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), mg, nzone, and t.
void Mtbl_cf::set_etat_zero | ( | ) |
void Mtbl_cf::ssint | ( | ) |
Definition at line 150 of file valeur_ssint.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_T, nzone, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, and TRA_T.
void Mtbl_cf::sx | ( | ) |
(r -sampling = RARE
) \ Id (r sampling = FIN
) \ (r -sampling = UNSURR
)
Definition at line 139 of file valeur_sx.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_R, nzone, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, t, and TRA_R.
void Mtbl_cf::sx2 | ( | ) |
(r -sampling = RARE
) \ Id (r sampling = FIN
) \ (r -sampling = UNSURR
)
Definition at line 148 of file valeur_sx2.C.
References Base_val::b, base, etat, MAX_BASE, MSQ_R, nzone, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, t, and TRA_R.
void Mtbl_cf::sxm1_zec | ( | ) |
double Mtbl_cf::val_in_bound_jk | ( | int | l, | |
int | j, | |||
int | k | |||
) | const |
Computes the angular coefficient of index j
,k of the field represented by *this
at by means of the spectral expansion.
Not defined in the nucleus.
l | [input] index of the domain | |
j | [input] index in -direction | |
k | [input] index in -direction |
*this
. Definition at line 448 of file mtbl_cf_val_point.C.
References Base_val::b, base, Mg3d::get_nr(), mg, MSQ_R, operator()(), R_CHEB, R_CHEBU, and TRA_R.
double Mtbl_cf::val_out_bound_jk | ( | int | l, | |
int | j, | |||
int | k | |||
) | const |
Computes the angular coefficient of index j
,k of the field represented by *this
at by means of the spectral expansion.
l | [input] index of the domain | |
j | [input] index in -direction | |
k | [input] index in -direction |
*this
. Definition at line 424 of file mtbl_cf_val_point.C.
References Base_val::b, base, Mg3d::get_nr(), mg, MSQ_R, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, and TRA_R.
double Mtbl_cf::val_point | ( | int | l, | |
double | x, | |||
double | theta, | |||
double | phi | |||
) | const |
Computes the value of the field represented by *this
at an arbitrary point, by means of the spectral expansion.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
theta | [input] value of the coordinate | |
phi | [input] value of the coordinate |
*this
. Definition at line 111 of file mtbl_cf_val_point.C.
References Base_val::b, base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), MAX_BASE, MAX_BASE_2, mg, MSQ_P, MSQ_R, MSQ_T, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, R_LEG, R_LEGI, R_LEGP, Tbl::t, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, TRA_P, TRA_R, and TRA_T.
double Mtbl_cf::val_point_asymy | ( | 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.
Case where the field is antisymmetric with respect to the y=0 plane.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
theta | [input] value of the coordinate | |
phi | [input] value of the coordinate |
*this
. Definition at line 67 of file mtbl_cf_vp_asymy.C.
References Base_val::b, base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), MAX_BASE, MAX_BASE_2, mg, MSQ_P, MSQ_R, MSQ_T, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, Tbl::t, t, T_COS, T_COS_P, T_COSSIN_CI, T_COSSIN_CP, T_SIN, T_SIN_P, TRA_P, TRA_R, and TRA_T.
double Mtbl_cf::val_point_jk | ( | int | l, | |
double | x, | |||
int | j, | |||
int | k | |||
) | const |
Computes the value of the field represented by *this
at an arbitrary point in , but collocation point in , by means of the spectral expansion.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
j | [input] index of the collocation point in | |
k | [input] index of the collocation point in |
*this
. Definition at line 254 of file mtbl_cf_val_point.C.
References Base_val::b, base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), MAX_BASE, mg, MSQ_R, MSQ_T, Base_val::phi_functions(), R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, Tbl::t, t, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, T_COSSIN_S, T_COSSIN_SI, T_COSSIN_SP, T_SIN, T_SIN_I, T_SIN_P, Base_val::theta_functions(), and TRA_R.
double Mtbl_cf::val_point_jk_asymy | ( | 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.
Case where the field is antisymmetric with respect to the y=0 plane.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
j | [input] index of the collocation point in | |
k | [input] index of the collocation point in |
*this
. Definition at line 199 of file mtbl_cf_vp_asymy.C.
References Base_val::b, base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), MAX_BASE, mg, MSQ_R, MSQ_T, Base_val::phi_functions(), R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, Tbl::t, t, T_COSSIN_CP, Base_val::theta_functions(), and TRA_R.
double Mtbl_cf::val_point_jk_symy | ( | 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.
Case where the field is symmetric with respect to the y=0 plane.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
j | [input] index of the collocation point in | |
k | [input] index of the collocation point in |
*this
. Definition at line 197 of file mtbl_cf_vp_symy.C.
References Base_val::b, base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), MAX_BASE, mg, MSQ_R, MSQ_T, Base_val::phi_functions(), R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, Tbl::t, t, T_COSSIN_CP, Base_val::theta_functions(), and TRA_R.
double Mtbl_cf::val_point_symy | ( | 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.
Case where the field is symmetric with respect to the y=0 plane.
l | [input] index of the domain | |
x | [input] value of the coordinate | |
theta | [input] value of the coordinate | |
phi | [input] value of the coordinate |
*this
. Definition at line 65 of file mtbl_cf_vp_symy.C.
References Base_val::b, base, etat, Tbl::get_etat(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), MAX_BASE, MAX_BASE_2, mg, MSQ_P, MSQ_R, MSQ_T, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, Tbl::t, t, T_COS, T_COS_P, T_COSSIN_CI, T_COSSIN_CP, T_SIN, T_SIN_P, TRA_P, TRA_R, and TRA_T.
ostream& operator<< | ( | ostream & | , | |
const Mtbl_cf & | ||||
) | [friend] |
Display.
int Mtbl_cf::etat [private] |
const Mg3d* Mtbl_cf::mg [private] |
int Mtbl_cf::nzone [private] |