25char scalar_manip_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Tensor/Scalar/scalar_manip.C,v 1.19 2014/10/13 08:53:46 j_novak Exp $" ;
104#include "utilitaires.h"
113 assert (
etat != ETATNONDEF) ;
114 assert (l_min <= l_max) ;
115 assert (l_min >= 0) ;
116 if (
etat == ETATZERO )
119 if (
etat == ETATUN) {
127 int l_q, m_q, base_r ;
128 for (
int lz=0; lz<
mp->get_mg()->get_nzone(); lz++)
129 if (m_coef(lz).
get_etat() != ETATZERO)
130 for (
int k=0; k<
mp->get_mg()->get_np(lz)+1; k++)
131 for (
int j=0; j<
mp->get_mg()->get_nt(lz); j++)
132 for (
int i=0; i<
mp->get_mg()->get_nr(lz); i++) {
134 if ((l_min <= l_q) && (l_q<= l_max))
135 m_coef.
set(lz, k, j, i) = 0 ;
153 assert (
etat != ETATNONDEF) ;
154 if ( (
etat == ETATZERO) || (
etat == ETATUN) )
157 int nz =
mp->get_mg()->get_nzone() ;
158 int np =
mp->get_mg()->get_np(nz-1) ;
159 int nt =
mp->get_mg()->get_nt(nz-1) ;
160 int nr =
mp->get_mg()->get_nr(nz-1) ;
165 va.set_etat_cf_qcq() ;
168 for (
int k=0 ; k<np+1 ; k++)
170 for (
int j=0 ; j<nt ; j++)
171 for (
int i=nr-1 ; i>nr-1-n ; i--)
172 va.c_cf->set(nz-1, k, j, i) = 0 ;
181 assert (
etat != ETATNONDEF) ;
182 if ( (
etat == ETATZERO) || (
etat == ETATUN) )
188 va.set_etat_cf_qcq() ;
189 int nz =
mp->get_mg()->get_nzone() ;
190 int* nr =
new int[nz];
191 int* nt =
new int[nz];
192 int* np =
new int[nz];
193 for (
int l=0; l<=nz-1; l++) {
194 nr[l] =
mp->get_mg()->get_nr(l) ;
195 nt[l] =
mp->get_mg()->get_nt(l) ;
196 np[l] =
mp->get_mg()->get_np(l) ;
199 for (
int l=0; l<=nz-1; l++)
200 for (
int k=0 ; k<np[l]+1 ; k++)
202 for (
int j=0 ; j<nt[l] ; j++)
203 for (
int i=nr[l]-1; i>nr[l]-1-nn[l] ; i--)
204 va.c_cf->set(l, k, j, i) = 0 ;
219 assert (
etat != ETATNONDEF) ;
220 if ( (
etat == ETATZERO) || (
etat == ETATUN) )
226 va.set_etat_cf_qcq() ;
227 int nr =
mp->get_mg()->get_nr(nz) ;
228 int nt =
mp->get_mg()->get_nt(nz) ;
229 int np =
mp->get_mg()->get_np(nz) ;
231 for (
int k=0 ; k<np+1 ; k++)
233 for (
int j=0 ; j<nt ; j++)
234 for (
int i=nr-1; i>nr-1-n ; i--)
235 va.c_cf->set(nz, k, j, i) = 0 ;
250 assert (
etat != ETATNONDEF) ;
251 if ( (
etat == ETATZERO) || (
etat == ETATUN) )
257 va.set_etat_cf_qcq() ;
258 int np =
mp->get_mg()->get_np(nz) ;
259 int nt =
mp->get_mg()->get_nt(nz) ;
260 int nr =
mp->get_mg()->get_nr(nz) ;
262 for (
int k=np+1-n ; k<np+1 ; k++)
263 for (
int j=0 ; j<nt ; j++)
264 for (
int i=0 ; i<nr ; i++)
265 va.c_cf->set(nz, k, j, i) = 0 ;
272 va.filtre_tp(nn, nz1, nz2) ;
285 assert (
etat != ETATNONDEF) ;
286 if (
etat == ETATZERO) {
287 if (x0 ==
double(0)) return ;
291 if (
etat == ETATUN) {
292 if (x0 ==
double(1)) return ;
293 else etat = ETATQCQ ;
298 int nt =
mp->get_mg()->get_nt(l0) ;
299 int np =
mp->get_mg()->get_np(l0) ;
302 va.set_etat_c_qcq() ;
304 for (
int k=0 ; k<np ; k++)
305 for (
int j=0 ; j<nt ; j++)
306 va.set(l0, k, j, 0) = x0 ;
317 assert (
etat != ETATNONDEF) ;
318 if (
etat == ETATZERO) {
319 if (x0 ==
double(0)) return ;
323 if (
etat == ETATUN) {
324 if (x0 ==
double(1)) return ;
325 else etat = ETATQCQ ;
330 int nrm1 =
mp->get_mg()->get_nr(l0) - 1 ;
331 int nt =
mp->get_mg()->get_nt(l0) ;
332 int np =
mp->get_mg()->get_np(l0) ;
335 va.set_etat_c_qcq() ;
337 for (
int k=0 ; k<np ; k++)
338 for (
int j=0 ; j<nt ; j++)
339 va.set(l0, k, j, nrm1) = x0 ;
355 int nz =
mp->get_mg()->get_nzone() ;
356 int np =
mp->get_mg()->get_np(nz-1) ;
357 int nt =
mp->get_mg()->get_nt(nz-1) ;
358 int nr =
mp->get_mg()->get_nr(nz-1) ;
362 cout <<
"Le mapping doit etre affine" << endl ;
371 va.set_etat_cf_qcq() ;
373 for (
int conte=0 ; conte<nbre ; conte++) {
380 double* coloc =
new double [nr] ;
381 int * deg =
new int[3] ;
386 for (
int i=0 ; i<nr ; i++)
387 coloc[i] =
pow(alpha,
double(conte))*
388 pow(-1-
cos(M_PI*i/(nr-1)),
double(conte)) ;
390 cfrcheb(deg, deg, coloc, deg, coloc) ;
392 for (
int k=0 ; k<np+1 ; k++)
394 for (
int j=0 ; j<nt ; j++) {
397 double* coef =
new double [nr] ;
398 double* auxi =
new double[1] ;
399 for (
int i=0 ; i<nr ; i++)
400 coef[i] = (*courant.
va.
c_cf)(nz-1, k, j, i) ;
403 som_r_chebu (coef, nr, 1, 1, 1, auxi) ;
406 som_r_pas_prevu (coef, nr, 1, 1, 1, auxi) ;
413 courant.
va.
c_cf->
set(nz-1, k, j, 0) -= *auxi ;
415 for (
int i=0 ; i<nr ; i++)
416 this->
va.c_cf->set(nz-1, k, j, i) -= *auxi * coloc[i] ;
432 if (output_ylm)
va.ylm() ;
434 int np =
mp->get_mg()->get_np(l_zone) ;
435 int nt =
mp->get_mg()->get_nt(l_zone) ;
439 assert(
etat == ETATQCQ) ;
441 for (
int k=0; k<np+2; k++)
442 for (
int j=0; j<nt; j++)
443 resu.
set(k, j) =
va.c_cf->val_out_bound_jk(l_zone, j, k) ;
449 assert(
mp->get_mg()->get_type_r(l_zone) != RARE) ;
451 if (output_ylm)
va.ylm() ;
453 int np =
mp->get_mg()->get_np(l_zone) ;
454 int nt =
mp->get_mg()->get_nt(l_zone) ;
458 assert(
etat == ETATQCQ) ;
460 for (
int k=0; k<np+2; k++)
461 for (
int j=0; j<nt; j++)
462 resu.
set(k, j) =
va.c_cf->val_in_bound_jk(l_zone, j, k) ;
469 if (output_ylm)
va.ylm() ;
478 assert(
etat == ETATQCQ) ;
480 int np =
mp->get_mg()->get_np(l_zone) ;
481 int nt =
mp->get_mg()->get_nt(l_zone) ;
482 for (
int k=0; k<np+2; k++)
483 for (
int j=0; j<nt; j++)
485 =
va.c_cf->val_out_bound_jk(l_zone, j, k) ;
Bases of the spectral expansions.
void set_base_t(int base_t)
Sets the expansion basis for functions in all domains.
int get_base_r(int l) const
Returns the expansion basis for r ( ) functions in the domain of index l (e.g.
void give_quant_numbers(int, int, int, int &, int &, int &) const
Computes the various quantum numbers and 1d radial base.
const double * get_alpha() const
Returns the pointer on the array alpha.
Coefficients storage for the multi-domain spectral method.
Tbl & set(int l)
Read/write of the Tbl containing the coefficients in a given domain.
void annule_l(int l_min, int l_max, bool ylm_output=false)
Sets all the multipolar components between l_min and l_max to zero.
virtual void del_deriv() const
Logical destructor of the derivatives.
Scalar scalar_out_bound(int n, bool leave_ylm=false)
Returns the Scalar containing the values of angular coefficients at the outer boundary.
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
void fixe_decroissance(int puis)
Substracts all the components behaving like in the external domain, with n strictly lower than puis ...
Tbl tbl_in_bound(int n, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the inner boundary.
friend Scalar cos(const Scalar &)
Cosine.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Tbl tbl_out_bound(int l_dom, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the outer boundary.
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Scalar(const Map &mpi)
Constructor from mapping.
Valeur & set_spectral_va()
Returns va (read/write version).
void annule_hard()
Sets the Scalar to zero in a hard way.
void mult_r_ced()
Multiplication by r in the compactified external domain (CED), the dzpuis flag is not changed.
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...
int etat
The logical state ETATNONDEF (undefined), ETATZERO (null), ETATUN (one), or ETATQCQ (ordinary).
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
const Base_val & get_spectral_base() const
Returns the spectral bases of the Valeur va.
void set_inner_boundary(int l, double x)
Sets the value of the Scalar at the inner boundary of a given domain.
void filtre_r(int *nn)
Sets the n lasts coefficients in r to 0 in all domains.
void filtre_phi(int n, int zone)
Sets the n lasts coefficients in to 0 in the domain zone .
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
Valeur va
The numerical value of the Scalar.
friend Scalar pow(const Scalar &, int)
Power .
void set_spectral_base(const Base_val &)
Sets the spectral bases of the Valeur va.
int dzpuis
Power of r by which the quantity represented by this must be divided in the compactified external d...
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
double & set(int i)
Read/write of a particular element (index i) (1D case).
void set_etat_cf_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_c...
Mtbl * c
Values of the function at the points of the multi-grid.
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
void coef() const
Computes the coeffcients of *this.
Base_val base
Bases on which the spectral expansion is performed.
#define R_CHEBU
base de Chebychev ordinaire (fin), dev. en 1/r
const Map *const mp
Mapping on which the numerical values at the grid points are defined.