LORENE
Lorene::Mtbl_cf Class Reference

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. More...
 
 Mtbl_cf (const Mg3d *p_mgrid, const Base_val &basis)
 Constructor. More...
 
 Mtbl_cf (const Mg3d &, FILE *)
 Constructor from a file (see sauve(FILE*) ) More...
 
 Mtbl_cf (const Mtbl_cf &)
 Copy constructor. More...
 
 ~Mtbl_cf ()
 Destructor. More...
 
void operator= (const Mtbl_cf &)
 Assignement to another Mtbl_cf. More...
 
void operator= (double)
 Assignement to a double. More...
 
void operator= (int)
 Assignement to a int. 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_qcq ()
 Sets the logical state to ETATQCQ (ordinary state). More...
 
void annule_hard ()
 Sets the Mtbl_cf to zero in a hard way. More...
 
void annule (int l_min, int l_max)
 Sets the Mtbl_cf to zero in some domains. More...
 
Tblset (int l)
 Read/write of the Tbl containing the coefficients in a given domain. More...
 
const Tbloperator() (int l) const
 Read-only of the Tbl containing the coefficients in a given domain. More...
 
double & set (int l, int k, int j, int i)
 Read/write of a particular element. More...
 
double operator() (int l, int k, int j, int i) const
 Read-only of a particular element. 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_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. More...
 
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. 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...
 
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 $\xi$, but collocation point in $(\theta', \phi')$, by means of the spectral expansion. More...
 
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 $\xi$, but collocation point in $(\theta', \phi')$, by means of the spectral expansion. More...
 
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 $\xi = 1 $ by means of the spectral expansion. More...
 
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 $\xi = -1 $ by means of the spectral expansion. More...
 
const Mg3dget_mg () const
 Returns the Mg3d on which the Mtbl_cf is defined. More...
 
int get_etat () const
 Returns the logical state. More...
 
int get_nzone () const
 Returns the number of zones (domains) More...
 
void sauve (FILE *) const
 Save in a file. More...
 
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. More...
 
void affiche_seuil (ostream &ostr, int precision=4, double threshold=1.e-7) const
 Prints only the values greater than a given threshold. More...
 
void operator+= (const Mtbl_cf &)
 += Mtbl_cf More...
 
void operator-= (const Mtbl_cf &)
 -= Mtbl_cf More...
 
void operator*= (double)
 *= double More...
 
void operator/= (double)
 /= double More...
 
void dsdx ()
 ${\partial \over \partial \xi} $ More...
 
void d2sdx2 ()
 ${\partial^2\over \partial \xi^2} $ More...
 
void sx ()
 ${1 \over \xi} $ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over \xi-1} $ (r -sampling = UNSURR ) More...
 
void sx2 ()
 ${1 \over \xi^2}$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ ${1 \over (\xi-1)^2}$ (r -sampling = UNSURR ) More...
 
void mult_x ()
 $\xi \, Id$ (r -sampling = RARE ) \ Id (r sampling = FIN ) \ $(\xi-1) \, Id $ (r -sampling = UNSURR ) More...
 
void sxm1_zec ()
 Id (r sampling = RARE, FIN ) \ ${1 \over (\xi-1)}$ (r -sampling = UNSURR ) More...
 
void mult_xm1_shell (int i)
 Id (r sampling = RARE, UNSURR ) \ $(\xi-1) \, Id$ (r -sampling = FIN ) More...
 
void mult_xp1_shell (int i)
 Id (r sampling = RARE, UNSURR ) \ $(\xi+1) \, Id$ (r -sampling = FIN ) More...
 
void mult_x_shell (int i)
 Id (r sampling = RARE, UNSURR ) \ $\xi \, Id$ (r -sampling = FIN ) More...
 
void mult_xm1_zec ()
 Id (r sampling = RARE, FIN ) \ $(\xi-1) \, Id$ (r -sampling = UNSURR ) More...
 
void mult2_xm1_zec ()
 Id (r sampling = RARE, FIN ) \ $(\xi-1)^2 \, Id$ (r -sampling = UNSURR ) More...
 
void dsdt ()
 ${\partial \over \partial \theta}$ More...
 
void d2sdt2 ()
 ${\partial^2 \over \partial \theta^2}$ More...
 
void ssint ()
 $Id\over\sin\theta$ More...
 
void scost ()
 $Id\over\cos\theta$ More...
 
void mult_ct ()
 $\cos\theta \, Id$ More...
 
void mult_st ()
 $\sin\theta \, Id$ More...
 
void dsdp ()
 ${\partial \over \partial \phi}$ More...
 
void d2sdp2 ()
 ${\partial^2 \over \partial \phi^2}$ More...
 
void mult_cp ()
 $\cos\phi \, Id$ More...
 
void mult_sp ()
 $\sin\phi \, Id$ More...
 
void lapang ()
 Angular Laplacian. More...
 
void poisson_angu (double lambda=0)
 Resolution of the generalized angular Poisson equation. More...
 

Public Attributes

Base_val base
 Bases of the spectral expansions. More...
 
Tbl ** t
 Array (size nzone ) of pointers on the Tbl 's which contain the spectral coefficients in each domain. More...
 

Private Member Functions

void del_t ()
 Logical destructor: dellocates the memory occupied by the Tbl array t . More...
 

Private Attributes

const Mg3dmg
 Pointer on the multi-grid Mgd3 on which this is defined. More...
 
int nzone
 Number of domains (zones) More...
 
int etat
 Logical state (ETATNONDEF , ETATQCQ or ETATZERO ). More...
 

Friends

ostream & operator<< (ostream &, const Mtbl_cf &)
 Display. More...
 

Detailed Description

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 $\phi$-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 196 of file mtbl_cf.h.

Constructor & Destructor Documentation

◆ Mtbl_cf() [1/4]

Lorene::Mtbl_cf::Mtbl_cf ( const Mg3d mgrid,
const Base_val basis 
)

Constructor.

Definition at line 127 of file mtbl_cf.C.

◆ Mtbl_cf() [2/4]

Lorene::Mtbl_cf::Mtbl_cf ( const Mg3d p_mgrid,
const Base_val basis 
)

Constructor.

Definition at line 135 of file mtbl_cf.C.

◆ Mtbl_cf() [3/4]

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

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

Definition at line 176 of file mtbl_cf.C.

References base, etat, Lorene::fread_be(), Lorene::Mg3d::get_nzone(), Lorene::Base_val::get_nzone(), mg, nzone, and t.

◆ Mtbl_cf() [4/4]

Lorene::Mtbl_cf::Mtbl_cf ( const Mtbl_cf mtc)

Copy constructor.

Definition at line 153 of file mtbl_cf.C.

References get_etat().

◆ ~Mtbl_cf()

Lorene::Mtbl_cf::~Mtbl_cf ( )

Destructor.

Definition at line 146 of file mtbl_cf.C.

References del_t().

Member Function Documentation

◆ affiche_seuil()

void Lorene::Mtbl_cf::affiche_seuil ( ostream &  ostr,
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
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 400 of file mtbl_cf.C.

References base, etat, and nzone.

◆ annule()

void Lorene::Mtbl_cf::annule ( int  l_min,
int  l_max 
)

Sets the Mtbl_cf to zero in some domains.

Parameters
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 335 of file mtbl_cf.C.

References etat, nzone, and set_etat_zero().

◆ annule_hard()

void Lorene::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 315 of file mtbl_cf.C.

References Lorene::Tbl::annule_hard(), etat, Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), mg, nzone, and t.

◆ d2sdp2()

void Lorene::Mtbl_cf::d2sdp2 ( )

${\partial^2 \over \partial \phi^2}$

Definition at line 135 of file valeur_d2sdp2.C.

References etat, MAX_BASE, P_COSSIN, and TRA_P.

◆ d2sdt2()

void Lorene::Mtbl_cf::d2sdt2 ( )

◆ d2sdx2()

void Lorene::Mtbl_cf::d2sdx2 ( )

◆ del_t()

void Lorene::Mtbl_cf::del_t ( )
private

Logical destructor: dellocates the memory occupied by the Tbl array t .

Definition at line 279 of file mtbl_cf.C.

References etat, nzone, and t.

◆ display()

void Lorene::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.

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 62 of file mtbl_cf_display.C.

References base, and etat.

◆ dsdp()

void Lorene::Mtbl_cf::dsdp ( )

${\partial \over \partial \phi}$

Definition at line 142 of file valeur_dsdp.C.

References etat, MAX_BASE_2, P_COSSIN, P_COSSIN_I, P_COSSIN_P, and TRA_P.

◆ dsdt()

void Lorene::Mtbl_cf::dsdt ( )

◆ dsdx()

void Lorene::Mtbl_cf::dsdx ( )

◆ get_etat()

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

Returns the logical state.

Definition at line 466 of file mtbl_cf.h.

References etat.

◆ get_mg()

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

Returns the Mg3d on which the Mtbl_cf is defined.

Definition at line 463 of file mtbl_cf.h.

References mg.

◆ get_nzone()

int Lorene::Mtbl_cf::get_nzone ( ) const
inline

Returns the number of zones (domains)

Definition at line 469 of file mtbl_cf.h.

References nzone.

◆ lapang()

void Lorene::Mtbl_cf::lapang ( )

◆ mult2_xm1_zec()

void Lorene::Mtbl_cf::mult2_xm1_zec ( )

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

Definition at line 106 of file valeur_mult2_xm1.C.

References etat, MAX_BASE, R_CHEBU, and TRA_R.

◆ mult_cp()

void Lorene::Mtbl_cf::mult_cp ( )

$\cos\phi \, Id$

Definition at line 123 of file valeur_mult_cp.C.

References etat, MAX_BASE_2, P_COSSIN, P_COSSIN_I, P_COSSIN_P, and TRA_P.

◆ mult_ct()

void Lorene::Mtbl_cf::mult_ct ( )

◆ mult_sp()

void Lorene::Mtbl_cf::mult_sp ( )

$\sin\phi \, Id$

Definition at line 123 of file valeur_mult_sp.C.

References etat, MAX_BASE_2, P_COSSIN, P_COSSIN_I, P_COSSIN_P, and TRA_P.

◆ mult_st()

void Lorene::Mtbl_cf::mult_st ( )

◆ mult_x()

void Lorene::Mtbl_cf::mult_x ( )

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

Definition at line 151 of file valeur_mult_x.C.

References etat, MAX_BASE, 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, and TRA_R.

◆ mult_x_shell()

void Lorene::Mtbl_cf::mult_x_shell ( int  i)

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

Parameters
i[input] : the index of the domain where the computation is done.

Definition at line 218 of file valeur_mult_x.C.

References Lorene::Mg3d::get_type_r(), and mg.

◆ mult_xm1_shell()

void Lorene::Mtbl_cf::mult_xm1_shell ( int  i)

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

Parameters
i[input] : the index of the domain where the computation is done.

Definition at line 110 of file valeur_mult_xm1.C.

References Lorene::Mg3d::get_type_r(), and mg.

◆ mult_xm1_zec()

void Lorene::Mtbl_cf::mult_xm1_zec ( )

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

Definition at line 156 of file valeur_mult_xm1.C.

References etat, MAX_BASE, R_CHEBU, and TRA_R.

◆ mult_xp1_shell()

void Lorene::Mtbl_cf::mult_xp1_shell ( int  i)

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

Parameters
i[input] : the index of the domain where the computation is done.

Definition at line 72 of file valeur_mult_xp1.C.

References Lorene::Mg3d::get_type_r(), and mg.

◆ operator()() [1/2]

const Tbl& Lorene::Mtbl_cf::operator() ( int  l) const
inline

Read-only of the Tbl containing the coefficients in a given domain.

Parameters
l[input] domain index

Definition at line 316 of file mtbl_cf.h.

References etat, and nzone.

◆ operator()() [2/2]

double Lorene::Mtbl_cf::operator() ( int  l,
int  k,
int  j,
int  i 
) const
inline

Read-only of a particular element.

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

Definition at line 342 of file mtbl_cf.h.

References etat.

◆ operator*=()

void Lorene::Mtbl_cf::operator*= ( double  )

*= double

Definition at line 432 of file mtbl_cf.C.

References Lorene::c_est_pas_fait().

◆ operator+=()

void Lorene::Mtbl_cf::operator+= ( const Mtbl_cf mi)

+= Mtbl_cf

Definition at line 305 of file mtbl_cf_arithm.C.

References base, etat, get_mg(), and mg.

◆ operator-=()

void Lorene::Mtbl_cf::operator-= ( const Mtbl_cf mi)

-= Mtbl_cf

Definition at line 328 of file mtbl_cf_arithm.C.

References base, etat, get_mg(), and mg.

◆ operator/=()

void Lorene::Mtbl_cf::operator/= ( double  )

/= double

Definition at line 437 of file mtbl_cf.C.

References Lorene::c_est_pas_fait().

◆ operator=() [1/3]

void Lorene::Mtbl_cf::operator= ( const Mtbl_cf mtc)

Assignement to another Mtbl_cf.

Definition at line 223 of file mtbl_cf.C.

References get_etat(), and mg.

◆ operator=() [2/3]

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

Assignement to a double.

Definition at line 245 of file mtbl_cf.C.

References nzone, set_etat_qcq(), set_etat_zero(), and t.

◆ operator=() [3/3]

void Lorene::Mtbl_cf::operator= ( int  m)

Assignement to a int.

Definition at line 259 of file mtbl_cf.C.

References nzone, set_etat_qcq(), set_etat_zero(), and t.

◆ poisson_angu()

void Lorene::Mtbl_cf::poisson_angu ( double  lambda = 0)

Resolution of the generalized angular Poisson equation.

The generalized angular Poisson equation is $\Delta_{\theta\varphi} u + \lambda u = \sigma$, where $\Delta_{\theta\varphi} u := \frac{\partial^2 u} {\partial \theta^2} + \frac{1}{\tan \theta} \frac{\partial u} {\partial \theta} +\frac{1}{\sin^2 \theta}\frac{\partial^2 u} {\partial \varphi^2}$.

Before the call to poisson_angu() , *this contains the coefficients of the source $\sigma$; after the call, it contains the coefficients of the solution $u$.

Parameters
lambda[input] coefficient $\lambda$ in the above equation (default value = 0)

Definition at line 86 of file mtbl_cf_pde.C.

References Lorene::Base_val::b, base, get_mg(), Lorene::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.

◆ sauve()

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

Save in a file.

Definition at line 207 of file mtbl_cf.C.

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

◆ scost()

void Lorene::Mtbl_cf::scost ( )

◆ set() [1/2]

Tbl& Lorene::Mtbl_cf::set ( int  l)
inline

Read/write of the Tbl containing the coefficients in a given domain.

Parameters
l[input] domain index

Definition at line 304 of file mtbl_cf.h.

References etat, and nzone.

◆ set() [2/2]

double& Lorene::Mtbl_cf::set ( int  l,
int  k,
int  j,
int  i 
)
inline

Read/write of a particular element.

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

Definition at line 329 of file mtbl_cf.h.

References etat, and nzone.

◆ set_etat_nondef()

void Lorene::Mtbl_cf::set_etat_nondef ( )

Sets the logical state to ETATNONDEF (undefined).

Deallocates the memory occupied by the Tbl array t .

Definition at line 297 of file mtbl_cf.C.

References etat.

◆ set_etat_qcq()

void Lorene::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 303 of file mtbl_cf.C.

References etat.

◆ set_etat_zero()

void Lorene::Mtbl_cf::set_etat_zero ( )

Sets the logical state to ETATZERO (zero).

Deallocates the memory occupied by the Tbl array t .

Definition at line 291 of file mtbl_cf.C.

References etat.

◆ ssint()

void Lorene::Mtbl_cf::ssint ( )

◆ sx()

void Lorene::Mtbl_cf::sx ( )

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

Definition at line 156 of file valeur_sx.C.

References etat, MAX_BASE, 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, and TRA_R.

◆ sx2()

void Lorene::Mtbl_cf::sx2 ( )

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

Definition at line 160 of file valeur_sx2.C.

References etat, MAX_BASE, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_LEG, R_LEGI, R_LEGP, and TRA_R.

◆ sxm1_zec()

void Lorene::Mtbl_cf::sxm1_zec ( )

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

Definition at line 115 of file valeur_sxm1_zec.C.

References etat, MAX_BASE, R_CHEBU, and TRA_R.

◆ val_in_bound_jk()

double Lorene::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 $\xi = -1 $ by means of the spectral expansion.

Not defined in the nucleus.

Parameters
l[input] index of the domain
j[input] index in $\theta$-direction
k[input] index in $\phi$-direction
Returns
coefficient at the boundary $\xi = -1$ in the domain no. l of the field whose spectral coefficients are stored in *this .

Definition at line 455 of file mtbl_cf_val_point.C.

References Lorene::Base_val::b, base, Lorene::Mg3d::get_nr(), mg, MSQ_R, operator()(), R_CHEB, R_CHEBU, and TRA_R.

◆ val_out_bound_jk()

double Lorene::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 $\xi = 1 $ by means of the spectral expansion.

Parameters
l[input] index of the domain
j[input] index in $\theta$-direction
k[input] index in $\phi$-direction
Returns
coefficient at the boundary $\xi = 1$ in the domain no. l of the field whose spectral coefficients are stored in *this .

Definition at line 431 of file mtbl_cf_val_point.C.

References Lorene::Base_val::b, base, Lorene::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.

◆ val_point()

double Lorene::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.

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 whose spectral coefficients are stored in *this .

Definition at line 118 of file mtbl_cf_val_point.C.

References etat, MAX_BASE, MAX_BASE_2, 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_SIN, T_SIN_I, T_SIN_P, TRA_P, TRA_R, and TRA_T.

◆ val_point_asymy()

double Lorene::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.

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 whose spectral coefficients are stored in *this .

Definition at line 74 of file mtbl_cf_vp_asymy.C.

References etat, MAX_BASE, MAX_BASE_2, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, T_COS, T_COS_P, T_COSSIN_CI, T_COSSIN_CP, T_SIN, T_SIN_P, TRA_P, TRA_R, and TRA_T.

◆ val_point_jk()

double Lorene::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 $\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 whose spectral coefficients are stored in *this .

Definition at line 261 of file mtbl_cf_val_point.C.

References etat, MAX_BASE, 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.

◆ val_point_jk_asymy()

double Lorene::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 $\xi$, but collocation point in $(\theta', \phi')$, by means of the spectral expansion.

Case where the field is antisymmetric with respect to the y=0 plane.

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 whose spectral coefficients are stored in *this .

Definition at line 206 of file mtbl_cf_vp_asymy.C.

References etat, MAX_BASE, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, and TRA_R.

◆ val_point_jk_symy()

double Lorene::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 $\xi$, but collocation point in $(\theta', \phi')$, by means of the spectral expansion.

Case where the field is symmetric with respect to the y=0 plane.

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 whose spectral coefficients are stored in *this .

Definition at line 204 of file mtbl_cf_vp_symy.C.

References etat, MAX_BASE, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, and TRA_R.

◆ val_point_symy()

double Lorene::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.

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 whose spectral coefficients are stored in *this .

Definition at line 72 of file mtbl_cf_vp_symy.C.

References etat, MAX_BASE, MAX_BASE_2, P_COSSIN, P_COSSIN_I, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, T_COS, T_COS_P, T_COSSIN_CI, T_COSSIN_CP, T_SIN, T_SIN_P, TRA_P, TRA_R, and TRA_T.

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  o,
const Mtbl_cf mt 
)
friend

Display.

Definition at line 369 of file mtbl_cf.C.

Member Data Documentation

◆ base

Base_val Lorene::Mtbl_cf::base

Bases of the spectral expansions.

Definition at line 210 of file mtbl_cf.h.

◆ etat

int Lorene::Mtbl_cf::etat
private

Logical state (ETATNONDEF , ETATQCQ or ETATZERO ).

Definition at line 206 of file mtbl_cf.h.

◆ mg

const Mg3d* Lorene::Mtbl_cf::mg
private

Pointer on the multi-grid Mgd3 on which this is defined.

Definition at line 202 of file mtbl_cf.h.

◆ nzone

int Lorene::Mtbl_cf::nzone
private

Number of domains (zones)

Definition at line 204 of file mtbl_cf.h.

◆ t

Tbl** Lorene::Mtbl_cf::t

Array (size nzone ) of pointers on the Tbl 's which contain the spectral coefficients in each domain.

Definition at line 215 of file mtbl_cf.h.


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