LORENE
scalar_manip.C
1 /*
2  * Copyright (c) 2003 Eric Gourgoulhon & Jerome Novak
3  *
4  * Copyright (c) 2000-2001 Philippe Grandclement (Cmp version)
5  *
6  * This file is part of LORENE.
7  *
8  * LORENE is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * LORENE is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with LORENE; if not, write to the Free Software
20  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
21  *
22  */
23 
24 
25 
26 
27 /*
28  * $Id: scalar_manip.C,v 1.23 2025/03/31 08:29:50 j_novak Exp $
29  * $Log: scalar_manip.C,v $
30  * Revision 1.23 2025/03/31 08:29:50 j_novak
31  * Methods to copyt Tensors and Scalars to smaller grids...
32  *
33  * Revision 1.22 2025/03/28 09:12:50 j_novak
34  * News methods to copy Tensors and Scalars to larger grids. To be used for de-aliasing.
35  *
36  * Revision 1.21 2018/11/16 14:34:37 j_novak
37  * Changed minor points to avoid some compilation warnings.
38  *
39  * Revision 1.20 2016/12/05 16:18:18 j_novak
40  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
41  *
42  * Revision 1.19 2014/10/13 08:53:46 j_novak
43  * Lorene classes and functions now belong to the namespace Lorene.
44  *
45  * Revision 1.18 2014/10/06 15:16:15 j_novak
46  * Modified #include directives to use c++ syntax.
47  *
48  * Revision 1.17 2008/10/03 09:03:52 j_novak
49  * Correction of yet another mistake (the array values in physical space was not
50  * destroyed).
51  *
52  * Revision 1.16 2008/09/30 08:35:18 j_novak
53  * Correction of forgotten call to coef()
54  *
55  * Revision 1.15 2008/09/29 13:23:51 j_novak
56  * Implementation of the angular mapping associated with an affine
57  * mapping. Things must be improved to take into account the domain index.
58  *
59  * Revision 1.14 2008/09/22 19:08:01 j_novak
60  * New methods to deal with boundary conditions
61  *
62  * Revision 1.13 2007/06/21 20:00:00 k_taniguchi
63  * Addition of the method filtre_r (int n, int nz)
64  *
65  * Revision 1.12 2006/06/28 07:46:41 j_novak
66  * Better treatment in the case of a domain set to zero.
67  *
68  * Revision 1.11 2005/09/07 13:39:10 j_novak
69  * *** empty log message ***
70  *
71  * Revision 1.10 2005/09/07 13:10:48 j_novak
72  * Added a filter setting to zero all mulitpoles in a given range.
73  *
74  * Revision 1.9 2004/11/23 12:47:44 f_limousin
75  * Add function filtre_tp(int nn, int nz1, int nz2).
76  *
77  * Revision 1.8 2004/05/07 11:24:54 f_limousin
78  * Implement new method filtre_r (int* nn)
79  *
80  * Revision 1.7 2004/02/27 09:47:26 f_limousin
81  * New methods filtre_phi(int) and filtre_theta(int).
82  *
83  * Revision 1.6 2004/01/23 13:26:28 e_gourgoulhon
84  * Added methods set_inner_boundary and set_outer_boundary.
85  * Methods set_val_inf and set_val_hor, which are particular cases of
86  * the above, have been suppressed.
87  *
88  * Revision 1.5 2003/11/13 13:43:55 p_grandclement
89  * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
90  *
91  * Revision 1.4 2003/10/11 14:47:17 e_gourgoulhon
92  * Lines 112 and 145 : replaced "0" by "double(0)" in the comparison
93  * statement.
94  *
95  * Revision 1.3 2003/10/10 15:57:29 j_novak
96  * Added the state one (ETATUN) to the class Scalar
97  *
98  * Revision 1.2 2003/10/08 14:24:09 j_novak
99  * replaced mult_r_zec with mult_r_ced
100  *
101  * Revision 1.1 2003/09/25 09:33:36 j_novak
102  * Added methods for integral calculation and various manipulations
103  *
104  *
105  * $Header: /cvsroot/Lorene/C++/Source/Tensor/Scalar/scalar_manip.C,v 1.23 2025/03/31 08:29:50 j_novak Exp $
106  *
107  */
108 
109 //standard
110 #include <cstdlib>
111 #include <cmath>
112 
113 // Lorene
114 #include "tensor.h"
115 #include "proto.h"
116 #include "utilitaires.h"
117 
118 /*
119  * Annule tous les l entre l_min et l_max (compris)
120  */
121 
122 namespace Lorene {
123 void Scalar::annule_l (int l_min, int l_max, bool ylm_output) {
124 
125  assert (etat != ETATNONDEF) ;
126  assert (l_min <= l_max) ;
127  assert (l_min >= 0) ;
128  if (etat == ETATZERO )
129  return ;
130 
131  if (etat == ETATUN) {
132  if (l_min == 0) set_etat_zero() ;
133  else return ;
134  }
135 
136  va.ylm() ;
137  Mtbl_cf& m_coef = *va.c_cf ;
138  const Base_val& base = va.base ;
139  int l_q, m_q, base_r ;
140  for (int lz=0; lz<mp->get_mg()->get_nzone(); lz++)
141  if (m_coef(lz).get_etat() != ETATZERO)
142  for (int k=0; k<mp->get_mg()->get_np(lz)+1; k++)
143  for (int j=0; j<mp->get_mg()->get_nt(lz); j++)
144  for (int i=0; i<mp->get_mg()->get_nr(lz); i++) {
145  base.give_quant_numbers(lz, k, j, m_q, l_q, base_r) ;
146  if ((l_min <= l_q) && (l_q<= l_max))
147  m_coef.set(lz, k, j, i) = 0 ;
148  }
149  if (va.c != 0x0) {
150  delete va.c ;
151  va.c = 0x0 ;
152  }
153  if (!ylm_output)
154  va.ylm_i() ;
155 
156  return ;
157 }
158 
159 /*
160  * Annule les n derniers coefficients en r dans la derniere zone
161  */
162 
163 void Scalar::filtre (int n) {
164 
165  assert (etat != ETATNONDEF) ;
166  if ( (etat == ETATZERO) || (etat == ETATUN) )
167  return ;
168 
169  int nz = mp->get_mg()->get_nzone() ;
170  int np = mp->get_mg()->get_np(nz-1) ;
171  int nt = mp->get_mg()->get_nt(nz-1) ;
172  int nr = mp->get_mg()->get_nr(nz-1) ;
173 
174  del_deriv() ;
175 
176  va.coef() ;
177  va.set_etat_cf_qcq() ;
178 
179 
180  for (int k=0 ; k<np+1 ; k++)
181  if (k!=1)
182  for (int j=0 ; j<nt ; j++)
183  for (int i=nr-1 ; i>nr-1-n ; i--)
184  va.c_cf->set(nz-1, k, j, i) = 0 ;
185 }
186 
187 
188 /*
189  * Annule les n derniers coefficients en r dans toutes les zones
190  */
191 
192 void Scalar::filtre_r (int* nn) {
193  assert (etat != ETATNONDEF) ;
194  if ( (etat == ETATZERO) || (etat == ETATUN) )
195  return ;
196 
197  del_deriv() ;
198 
199  va.coef() ;
200  va.set_etat_cf_qcq() ;
201  int nz = mp->get_mg()->get_nzone() ;
202  int* nr = new int[nz];
203  int* nt = new int[nz];
204  int* np = new int[nz];
205  for (int l=0; l<=nz-1; l++) {
206  nr[l] = mp->get_mg()->get_nr(l) ;
207  nt[l] = mp->get_mg()->get_nt(l) ;
208  np[l] = mp->get_mg()->get_np(l) ;
209  }
210 
211  for (int l=0; l<=nz-1; l++)
212  for (int k=0 ; k<np[l]+1 ; k++)
213  if (k!=1)
214  for (int j=0 ; j<nt[l] ; j++)
215  for (int i=nr[l]-1; i>nr[l]-1-nn[l] ; i--)
216  va.c_cf->set(l, k, j, i) = 0 ;
217 
218  if (va.c != 0x0) {
219  delete va.c ;
220  va.c = 0x0 ;
221  }
222 
223 }
224 
225 
226 /*
227  * Annule les n derniers coefficients en r dans zone nz
228  */
229 
230 void Scalar::filtre_r (int n, int nz) {
231  assert (etat != ETATNONDEF) ;
232  if ( (etat == ETATZERO) || (etat == ETATUN) )
233  return ;
234 
235  del_deriv() ;
236 
237  va.coef() ;
238  va.set_etat_cf_qcq() ;
239  int nr = mp->get_mg()->get_nr(nz) ;
240  int nt = mp->get_mg()->get_nt(nz) ;
241  int np = mp->get_mg()->get_np(nz) ;
242 
243  for (int k=0 ; k<np+1 ; k++)
244  if (k!=1)
245  for (int j=0 ; j<nt ; j++)
246  for (int i=nr-1; i>nr-1-n ; i--)
247  va.c_cf->set(nz, k, j, i) = 0 ;
248 
249  if (va.c != 0x0) {
250  delete va.c ;
251  va.c = 0x0 ;
252  }
253 
254 }
255 
256 
257 /*
258  * Annule les n derniers coefficients en phi dans zone nz
259  */
260 
261 void Scalar::filtre_phi (int n, int nz) {
262  assert (etat != ETATNONDEF) ;
263  if ( (etat == ETATZERO) || (etat == ETATUN) )
264  return ;
265 
266  del_deriv() ;
267 
268  va.coef() ;
269  va.set_etat_cf_qcq() ;
270  int np = mp->get_mg()->get_np(nz) ;
271  int nt = mp->get_mg()->get_nt(nz) ;
272  int nr = mp->get_mg()->get_nr(nz) ;
273 
274  for (int k=np+1-n ; k<np+1 ; k++)
275  for (int j=0 ; j<nt ; j++)
276  for (int i=0 ; i<nr ; i++)
277  va.c_cf->set(nz, k, j, i) = 0 ;
278 
279 }
280 
281 
282 void Scalar::filtre_tp(int nn, int nz1, int nz2) {
283 
284  va.filtre_tp(nn, nz1, nz2) ;
285 
286 }
287 
288 
289  /* Sets the value of the {\tt Scalar} at the inner boundary of a given
290  * domain.
291  * @param l [input] domain index
292  * @param x [input] (constant) value at the inner boundary of domain no. {\tt l}
293  */
294 
295 void Scalar::set_inner_boundary(int l0, double x0) {
296 
297  assert (etat != ETATNONDEF) ;
298  if (etat == ETATZERO) {
299  if (x0 == double(0)) return ;
300  else annule_hard() ;
301  }
302 
303  if (etat == ETATUN) {
304  if (x0 == double(1)) return ;
305  else etat = ETATQCQ ;
306  }
307 
308  del_deriv() ;
309 
310  int nt = mp->get_mg()->get_nt(l0) ;
311  int np = mp->get_mg()->get_np(l0) ;
312 
313  va.coef_i() ;
314  va.set_etat_c_qcq() ;
315 
316  for (int k=0 ; k<np ; k++)
317  for (int j=0 ; j<nt ; j++)
318  va.set(l0, k, j, 0) = x0 ;
319 }
320 
321  /* Sets the value of the {\tt Scalar} at the outer boundary of a given
322  * domain.
323  * @param l [input] domain index
324  * @param x [input] (constant) value at the outer boundary of domain no. {\tt l}
325  */
326 
327 void Scalar::set_outer_boundary(int l0, double x0) {
328 
329  assert (etat != ETATNONDEF) ;
330  if (etat == ETATZERO) {
331  if (x0 == double(0)) return ;
332  else annule_hard() ;
333  }
334 
335  if (etat == ETATUN) {
336  if (x0 == double(1)) return ;
337  else etat = ETATQCQ ;
338  }
339 
340  del_deriv() ;
341 
342  int nrm1 = mp->get_mg()->get_nr(l0) - 1 ;
343  int nt = mp->get_mg()->get_nt(l0) ;
344  int np = mp->get_mg()->get_np(l0) ;
345 
346  va.coef_i() ;
347  va.set_etat_c_qcq() ;
348 
349  for (int k=0 ; k<np ; k++)
350  for (int j=0 ; j<nt ; j++)
351  va.set(l0, k, j, nrm1) = x0 ;
352 }
353 
354 /*
355  * Permet de fixer la decroissance du cmp a l infini en viurant les
356  * termes en 1/r^n
357  */
358 void Scalar::fixe_decroissance (int puis) {
359 
360  if (puis<dzpuis)
361  return ;
362  else {
363 
364  int nbre = puis-dzpuis ;
365 
366  // le confort avant tout ! (c'est bien le confort ...)
367  int nz = mp->get_mg()->get_nzone() ;
368  int np = mp->get_mg()->get_np(nz-1) ;
369  int nt = mp->get_mg()->get_nt(nz-1) ;
370  int nr = mp->get_mg()->get_nr(nz-1) ;
371 
372  const Map_af* map = dynamic_cast<const Map_af*>(mp) ;
373  if (map == 0x0) {
374  cout << "Le mapping doit etre affine" << endl ;
375  abort() ;
376  }
377 
378  double alpha = map->get_alpha()[nz-1] ;
379 
380  Scalar courant (*this) ;
381 
382  va.coef() ;
383  va.set_etat_cf_qcq() ;
384 
385  for (int conte=0 ; conte<nbre ; conte++) {
386 
387  int base_r = courant.va.base.get_base_r(nz-1) ;
388 
389  courant.va.coef() ;
390 
391  // On calcul les coefficients de 1/r^conte
392  double* coloc = new double [nr] ;
393  int * deg = new int[3] ;
394  deg[0] = 1 ;
395  deg[1] = 1 ;
396  deg[2] = nr ;
397 
398  for (int i=0 ; i<nr ; i++)
399  coloc[i] =pow(alpha, double(conte))*
400  pow(-1-cos(M_PI*i/(nr-1)), double(conte)) ;
401 
402  cfrcheb(deg, deg, coloc, deg, coloc) ;
403 
404  for (int k=0 ; k<np+1 ; k++)
405  if (k != 1)
406  for (int j=0 ; j<nt ; j++) {
407 
408  // On doit determiner le coefficient du truc courant :
409  double* coef = new double [nr] ;
410  double* auxi = new double[1] ;
411  for (int i=0 ; i<nr ; i++)
412  coef[i] = (*courant.va.c_cf)(nz-1, k, j, i) ;
413  switch (base_r) {
414  case R_CHEBU :
415  som_r_chebu (coef, nr, 1, 1, 1, auxi) ;
416  break ;
417  default :
418  som_r_pas_prevu (coef, nr, 1, 1, 1, auxi) ;
419  break ;
420  }
421 
422  // On modifie le cmp courant :
423  courant.va.coef() ;
424  courant.va.set_etat_cf_qcq() ;
425  courant.va.c_cf->set(nz-1, k, j, 0) -= *auxi ;
426 
427  for (int i=0 ; i<nr ; i++)
428  this->va.c_cf->set(nz-1, k, j, i) -= *auxi * coloc[i] ;
429 
430 
431  delete [] coef ;
432  delete [] auxi ;
433  }
434  delete [] coloc ;
435  delete [] deg ;
436 
437  courant.mult_r_ced() ;
438  }
439  }
440 }
441 
442 Tbl Scalar::tbl_out_bound(int l_zone, bool output_ylm) {
443  va.coef() ;
444  if (output_ylm) va.ylm() ;
445 
446  int np = mp->get_mg()->get_np(l_zone) ;
447  int nt = mp->get_mg()->get_nt(l_zone) ;
448  Tbl resu(np+2, nt) ;
449  if (etat == ETATZERO) resu.set_etat_zero() ;
450  else {
451  assert(etat == ETATQCQ) ;
452  resu.set_etat_qcq() ;
453  for (int k=0; k<np+2; k++)
454  for (int j=0; j<nt; j++)
455  resu.set(k, j) = va.c_cf->val_out_bound_jk(l_zone, j, k) ;
456  }
457  return resu ;
458 }
459 
460 Tbl Scalar::tbl_in_bound(int l_zone, bool output_ylm) {
461  assert(mp->get_mg()->get_type_r(l_zone) != RARE) ;
462  va.coef() ;
463  if (output_ylm) va.ylm() ;
464 
465  int np = mp->get_mg()->get_np(l_zone) ;
466  int nt = mp->get_mg()->get_nt(l_zone) ;
467  Tbl resu(np+2, nt) ;
468  if (etat == ETATZERO) resu.set_etat_zero() ;
469  else {
470  assert(etat == ETATQCQ) ;
471  resu.set_etat_qcq() ;
472  for (int k=0; k<np+2; k++)
473  for (int j=0; j<nt; j++)
474  resu.set(k, j) = va.c_cf->val_in_bound_jk(l_zone, j, k) ;
475  }
476  return resu ;
477 }
478 
479 Scalar Scalar::scalar_out_bound(int l_zone, bool output_ylm) {
480  va.coef() ;
481  if (output_ylm) va.ylm() ;
482 
483  Scalar resu(mp->mp_angu(l_zone)) ;
484  resu.std_spectral_base() ;
485  Base_val base = resu.get_spectral_base() ;
486  base.set_base_t(va.base.get_base_t(l_zone)) ;
487  resu.set_spectral_base(base) ;
488  if (etat == ETATZERO) resu.set_etat_zero() ;
489  else {
490  assert(etat == ETATQCQ) ;
491  resu.annule_hard() ;
492  int np = mp->get_mg()->get_np(l_zone) ;
493  int nt = mp->get_mg()->get_nt(l_zone) ;
494  for (int k=0; k<np+2; k++)
495  for (int j=0; j<nt; j++)
496  resu.set_spectral_va().c_cf->set(0, k, j, 0)
497  = va.c_cf->val_out_bound_jk(l_zone, j, k) ;
498  delete resu.set_spectral_va().c ;
499  resu.set_spectral_va().c = 0x0 ;
500  }
501  return resu ;
502 }
503 
504 
505  void Scalar::copy_coefs_from_small_grid(const Scalar& scal_small_grid) {
506 
507  assert(&scal_small_grid != this) ;
508 
509  const Base_val& b_i = scal_small_grid.get_spectral_base() ;
510 
511  const Mg3d& mgl = *mp->get_mg() ;
512  int nz = mgl.get_nzone( );
513 
514  const Mg3d& mgs = *scal_small_grid.get_mp().get_mg() ;
515 
516  // Some grid compatibility checks ...
517  //-------------------------------------
518  assert(scal_small_grid.etat != ETATNONDEF) ;
519  assert(mgs.get_nzone() == nz) ;
520 #ifndef NDEBUG
521  for (int l=0; l<nz; l++) {
522  assert(mgs.get_nr(l) <= mgl.get_nr(l)) ;
523  assert(mgs.get_nt(l) <= mgl.get_nt(l)) ;
524  assert(mgs.get_np(l) <= mgl.get_np(l)) ;
525  assert(mgs.get_type_r(l) == mgl.get_type_r(l)) ;
526  }
527  assert(mgs.get_type_t() == mgl.get_type_t()) ;
528  assert(mgs.get_type_p() == mgl.get_type_p()) ;
529 #endif
530 
531  if (scal_small_grid.get_etat() == ETATZERO) { // Nothing to do...
532  set_etat_zero() ;
533  return ;
534  }
535  else { // Copy must be done ...
536  set_etat_qcq() ;
537  va.set_etat_cf_qcq() ;
538  scal_small_grid.get_spectral_va().coef() ;
539  const Mtbl_cf& c_s = *(scal_small_grid.get_spectral_va().c_cf) ;
540  assert( c_s.get_etat() == ETATQCQ ) ;
541  Mtbl_cf& c_l = *(va.c_cf) ;
542  c_l.set_etat_qcq() ;
543  c_l.base = b_i ;
544 
545  for (int l=0; l<nz; l++) {
546  if (c_s.t[l]->get_etat() == ETATZERO) {
547  c_l.t[l]->set_etat_zero() ;
548  }
549  else {
550  int nrl = mgl.get_nr(l) ;
551  int nrs = mgs.get_nr(l) ;
552  int ntl = mgl.get_nt(l) ;
553  int nts = mgs.get_nt(l) ;
554  int npl = mgl.get_np(l)+2 ;
555  int nps = mgs.get_np(l)+2 ;
556  c_l.t[l]->set_etat_qcq() ;
557  double* p_out = c_l.t[l]->t ;
558  double* p_in = c_s.t[l]->t ;
559  for (int k=0; k<nps; k++) {
560  for (int j=0; j<nts; j++) {
561  for (int i=0; i<nrs; i++) {
562  *p_out = *p_in ;
563  p_out++ ; p_in++ ;
564  }
565  for (int i=nrs; i<nrl; i++) {
566  *p_out = 0. ; p_out++ ;
567  } // loop on i (r)
568  } // loop on j (theta)
569  for (int j=nts; j<ntl; j++) {
570  for (int i=0; i<nrl; i++) {
571  *p_out = 0. ; p_out++ ;
572  }
573  } // loop on j (theta)
574  } // loop on k (phi)
575  for (int k=nps; k<npl; k++) {
576  for (int j=0; j<ntl; j++) {
577  for (int i=0; i<nrl; i++) {
578  *p_out = 0. ; p_out++ ;
579  }
580  } // loop on j (theta)
581  } // loop on k (phi)
582  } // case etat != ETATZERO
583  } // loop on l (domains)
584  if (scal_small_grid.get_etat() == ETATUN) etat = ETATUN ;
585  } // case when copy must be done
586  dzpuis = scal_small_grid.get_dzpuis() ;
587  set_spectral_base(b_i ) ;
588  }
589 
590  void Scalar::copy_coefs_from_large_grid(const Scalar& scal_large_grid) {
591 
592  assert(&scal_large_grid != this) ;
593 
594  const Base_val& b_i = scal_large_grid.get_spectral_base() ;
595 
596  const Mg3d& mgs = *mp->get_mg() ;
597  int nz = mgs.get_nzone( );
598 
599  const Mg3d& mgl = *scal_large_grid.get_mp().get_mg() ;
600 
601  // Some grid compatibility checks ...
602  //-------------------------------------
603 #ifndef NDEBUG
604  assert(scal_large_grid.etat != ETATNONDEF) ;
605  assert(mgl.get_nzone() == nz) ;
606  for (int l=0; l<nz; l++) {
607  assert(mgs.get_nr(l) <= mgl.get_nr(l)) ;
608  assert(mgs.get_nt(l) <= mgl.get_nt(l)) ;
609  assert(mgs.get_np(l) <= mgl.get_np(l)) ;
610  assert(mgs.get_type_r(l) == mgl.get_type_r(l)) ;
611  }
612  assert(mgs.get_type_t() == mgl.get_type_t()) ;
613  assert(mgs.get_type_p() == mgl.get_type_p()) ;
614 #endif
615 
616  if (scal_large_grid.get_etat() == ETATZERO) { // Nothing to do...
617  set_etat_zero() ;
618  return ;
619  }
620  else { // Copy must be done ...
621  set_etat_qcq() ;
622  va.set_etat_cf_qcq() ;
623  scal_large_grid.get_spectral_va().coef() ;
624  const Mtbl_cf& c_l = *(scal_large_grid.get_spectral_va().c_cf) ;
625  assert( c_l.get_etat() == ETATQCQ ) ;
626  Mtbl_cf& c_s = *(va.c_cf) ;
627  c_s.set_etat_qcq() ;
628  c_s.base = b_i ;
629 
630  for (int l=0; l<nz; l++) {
631  if (c_l.t[l]->get_etat() == ETATZERO) {
632  c_s.t[l]->set_etat_zero() ;
633  }
634  else {
635  int nrs = mgs.get_nr(l) ;
636  int nts = mgs.get_nt(l) ;
637  int nps = mgs.get_np(l)+1 ;
638  int nrl = mgl.get_nr(l) ;
639  int dif_r = nrl - nrs ;
640  int dif_t = mgl.get_nt(l) - nts ;
641  c_s.t[l]->set_etat_qcq() ;
642  double* p_out = c_s.t[l]->t ;
643  double* p_in = c_l.t[l]->t ;
644  for (int k=0; k<nps; k++) {
645  for (int j=0; j<nts; j++) {
646  for (int i=0; i<nrs; i++) {
647  *p_out = *p_in ;
648  p_out++ ; p_in++ ;
649  }
650  p_in += dif_r ;
651  } // loop on j (theta)
652  p_in += dif_t*nrl ;
653  } // loop on k (phi)
654  // The last coef. in phi is set to 0
655  for (int j=0; j<nts; j++) {
656  for (int i=0; i<nrs; i++) {
657  *p_out = 0. ; p_out++ ;
658  }
659  } // loop on j (theta)
660  } // case etat != ETATZERO
661  } // loop on l (domains)
662  if (scal_large_grid.get_etat() == ETATUN) etat = ETATUN ;
663  } // case when copy must be done
664  dzpuis = scal_large_grid.get_dzpuis() ;
665  set_spectral_base(b_i ) ;
666  }
667 
668 
669 } // End of namespace Lorene
int get_type_p() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the transformatio...
Definition: grilles.h:514
const Base_val & get_spectral_base() const
Returns the spectral bases of the Valeur va.
Definition: scalar.h:1356
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.
Definition: scalar_manip.C:123
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
void filtre_r(int *nn)
Sets the n lasts coefficients in r to 0 in all domains.
Definition: scalar_manip.C:192
const double * get_alpha() const
Returns the pointer on the array alpha.
Definition: map_af.C:607
virtual const Map_af & mp_angu(int) const =0
Returns the "angular" mapping for the outside of domain l_zone.
void ylm_i()
Inverse of ylm()
Definition: valeur_ylm_i.C:134
void copy_coefs_from_small_grid(const Scalar &scal_small_grid)
Copies the content of the argument Scalar to this defined on a larger grid, with similar mappings...
Definition: scalar_manip.C:505
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
void set_etat_cf_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_c...
Definition: valeur.C:715
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:151
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
void give_quant_numbers(int, int, int, int &, int &, int &) const
Computes the various quantum numbers and 1d radial base.
Lorene prototypes.
Definition: app_hor.h:67
Tbl & set(int l)
Read/write of the Tbl containing the coefficients in a given domain.
Definition: mtbl_cf.h:304
int get_base_t(int l) const
Returns the expansion basis for functions in the domain of index l (e.g.
Definition: base_val.h:414
void ylm()
Computes the coefficients of *this.
Definition: valeur_ylm.C:141
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
void coef_i() const
Computes the physical value of *this.
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition: grilles.h:504
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
int get_etat() const
Gives the logical state.
Definition: tbl.h:414
Scalar scalar_out_bound(int n, bool leave_ylm=false)
Returns the Scalar containing the values of angular coefficients at the outer boundary.
Definition: scalar_manip.C:479
int dzpuis
Power of r by which the quantity represented by this must be divided in the compactified external do...
Definition: scalar.h:415
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:566
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
void set_base_t(int base_t)
Sets the expansion basis for functions in all domains.
Definition: base_val.C:198
Tbl tbl_in_bound(int n, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the inner boundary.
Definition: scalar_manip.C:460
int get_etat() const
Returns the logical state.
Definition: mtbl_cf.h:466
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 s...
void annule_hard()
Sets the Scalar to zero in a hard way.
Definition: scalar.C:386
int get_base_r(int l) const
Returns the expansion basis for r ( ) functions in the domain of index l (e.g.
Definition: base_val.h:403
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
friend Scalar pow(const Scalar &, int)
Power .
Definition: scalar_math.C:454
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 s...
friend Scalar cos(const Scalar &)
Cosine.
Definition: scalar_math.C:107
virtual void del_deriv() const
Logical destructor of the derivatives.
Definition: scalar.C:293
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:569
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
Definition: scalar_manip.C:327
double * t
The array of double.
Definition: tbl.h:176
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...
Definition: valeur.C:927
void mult_r_ced()
Multiplication by r in the compactified external domain (CED), the dzpuis flag is not changed...
Mtbl * c
Values of the function at the points of the multi-grid.
Definition: valeur.h:309
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
Valeur va
The numerical value of the Scalar.
Definition: scalar.h:417
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: mtbl_cf.C:303
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition: scalar_manip.C:163
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: tbl.C:350
void set_inner_boundary(int l, double x)
Sets the value of the Scalar at the inner boundary of a given domain.
Definition: scalar_manip.C:295
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...
Definition: scalar_manip.C:282
void fixe_decroissance(int puis)
Substracts all the components behaving like in the external domain, with n strictly lower than puis ...
Definition: scalar_manip.C:358
void set_spectral_base(const Base_val &)
Sets the spectral bases of the Valeur va
Definition: scalar.C:803
void copy_coefs_from_large_grid(const Scalar &scal_large_grid)
Copies the content of the argument Scalar to this defined on a smaller grid, with similar mappings...
Definition: scalar_manip.C:590
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
Multi-domain grid.
Definition: grilles.h:276
Bases of the spectral expansions.
Definition: base_val.h:325
void filtre_phi(int n, int zone)
Sets the n lasts coefficients in to 0 in the domain zone .
Definition: scalar_manip.C:261
Affine radial mapping.
Definition: map.h:2097
int etat
The logical state ETATNONDEF (undefined), ETATZERO (null), ETATUN (one), or ETATQCQ (ordinary)...
Definition: scalar.h:408
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c )...
Definition: valeur.C:704
Tbl tbl_out_bound(int l_dom, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the outer boundary.
Definition: scalar_manip.C:442
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:476
#define R_CHEBU
base de Chebychev ordinaire (fin), dev. en 1/r
Definition: type_parite.h:180
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:493
const Map *const mp
Mapping on which the numerical values at the grid points are defined.
Definition: tensor.h:307
const Map & get_mp() const
Returns the mapping.
Definition: tensor.h:902
Tbl ** t
Array (size nzone ) of pointers on the Tbl &#39;s which contain the spectral coefficients in each domain...
Definition: mtbl_cf.h:215
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:613
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition: valeur.h:373