LORENE
scalar_raccord_zec.C
1 /*
2  * Copyright (c) 2003 Eric Gourgoulhon & Jerome Novak
3  *
4  * Copyright (c) 2000-2001 Philippe Grandclement (for preceding Cmp version)
5  *
6  *
7  * This file is part of LORENE.
8  *
9  * LORENE is free software; you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation; either version 2 of the License, or
12  * (at your option) any later version.
13  *
14  * LORENE is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with LORENE; if not, write to the Free Software
21  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
22  *
23  */
24 
25 
26 
27 
28 /*
29  * $Id: scalar_raccord_zec.C,v 1.9 2016/12/05 16:18:19 j_novak Exp $
30  * $Log: scalar_raccord_zec.C,v $
31  * Revision 1.9 2016/12/05 16:18:19 j_novak
32  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
33  *
34  * Revision 1.8 2014/10/13 08:53:47 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.7 2014/10/06 15:16:16 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.6 2004/06/04 16:14:18 j_novak
41  * In smooth_decay, the configuration space was not up-to-date.
42  *
43  * Revision 1.5 2004/03/05 15:09:59 e_gourgoulhon
44  * Added method smooth_decay.
45  *
46  * Revision 1.4 2003/10/29 11:04:34 e_gourgoulhon
47  * dec2_dpzuis() replaced by dec_dzpuis(2).
48  * inc2_dpzuis() replaced by inc_dzpuis(2).
49  *
50  * Revision 1.3 2003/10/10 15:57:29 j_novak
51  * Added the state one (ETATUN) to the class Scalar
52  *
53  * Revision 1.2 2003/09/25 09:22:33 j_novak
54  * Added a #include<math.h>
55  *
56  * Revision 1.1 2003/09/25 08:58:10 e_gourgoulhon
57  * First version.
58  *
59  *
60  * $Header: /cvsroot/Lorene/C++/Source/Tensor/Scalar/scalar_raccord_zec.C,v 1.9 2016/12/05 16:18:19 j_novak Exp $
61  *
62  */
63 
64 //standard
65 #include <cstdlib>
66 #include <cmath>
67 
68 // LORENE
69 #include "matrice.h"
70 #include "tensor.h"
71 #include "proto.h"
72 #include "nbr_spx.h"
73 #include "utilitaires.h"
74 
75 // Fait le raccord C1 dans la zec ...
76 namespace Lorene {
77 // Suppose (pour le moment, le meme nbre de points sur les angles ...)
78 // et que la zone precedente est une coquille
79 
80 void Scalar::raccord_c1_zec(int puis, int nbre, int lmax) {
81 
82  assert (nbre>0) ;
83  assert (etat != ETATNONDEF) ;
84  if ((etat == ETATZERO) || (etat == ETATUN))
85  return ;
86 
87  // Le mapping doit etre affine :
88  const Map_af* map = dynamic_cast<const Map_af*>(mp) ;
89  if (map == 0x0) {
90  cout << "Le mapping doit etre affine" << endl ;
91  abort() ;
92  }
93 
94  int nz = map->get_mg()->get_nzone() ;
95  int nr = map->get_mg()->get_nr (nz-1) ;
96  int nt = map->get_mg()->get_nt (nz-1) ;
97  int np = map->get_mg()->get_np (nz-1) ;
98 
99  double alpha = map->get_alpha()[nz-1] ;
100  double r_cont = -1./2./alpha ; //Rayon de debut de la zec.
101 
102  // On calcul les coefficients des puissances de 1./r
103  Tbl coef (nbre+2*lmax, nr) ;
104  coef.set_etat_qcq() ;
105 
106  int* deg = new int[3] ;
107  deg[0] = 1 ; deg[1] = 1 ; deg[2] = nr ;
108  double* auxi = new double[nr] ;
109 
110  for (int conte=0 ; conte<nbre+2*lmax ; conte++) {
111  for (int i=0 ; i<nr ; i++)
112  auxi[i] = pow(-1-cos(M_PI*i/(nr-1)), (conte+puis)) ;
113 
114  cfrcheb(deg, deg, auxi, deg, auxi) ;
115  for (int i=0 ; i<nr ; i++)
116  coef.set(conte, i) = auxi[i]*pow (alpha, conte+puis) ;
117  }
118 
119  delete[] deg ;
120  // Maintenant on va calculer les valeurs de la ieme derivee :
121  Tbl valeurs (nbre, nt, np+1) ;
122  valeurs.set_etat_qcq() ;
123 
124  Scalar courant (*this) ;
125  double* res_val = new double[1] ;
126 
127  for (int conte=0 ; conte<nbre ; conte++) {
128 
129  courant.va.coef() ;
130  courant.va.ylm() ;
131  courant.va.c_cf->t[nz-1]->annule_hard() ;
132 
133  int base_r = courant.va.base.get_base_r(nz-2) ;
134  for (int k=0 ; k<np+1 ; k++)
135  for (int j=0 ; j<nt ; j++)
136  if (nullite_plm(j, nt, k, np, courant.va.base) == 1) {
137 
138  for (int i=0 ; i<nr ; i++)
139  auxi[i] = (*courant.va.c_cf)(nz-2, k, j, i) ;
140 
141  switch (base_r) {
142  case R_CHEB :
143  som_r_cheb (auxi, nr, 1, 1, 1, res_val) ;
144  break ;
145  default :
146  cout << "Cas non prevu dans raccord_zec" << endl ;
147  abort() ;
148  break ;
149  }
150  valeurs.set(conte, k, j) = res_val[0] ;
151  }
152  Scalar copie (courant) ;
153  copie.dec_dzpuis(2) ;
154  courant = copie.dsdr() ;
155  }
156 
157  delete [] auxi ;
158  delete [] res_val ;
159 
160  // On boucle sur les harmoniques : construction de la matrice
161  // et du second membre
162  va.coef() ;
163  va.ylm() ;
164  va.c_cf->t[nz-1]->annule_hard() ;
165  va.set_etat_cf_qcq() ;
166 
167  const Base_val& base = va.base ;
168  int base_r, l_quant, m_quant ;
169  for (int k=0 ; k<np+1 ; k++)
170  for (int j=0 ; j<nt ; j++)
171  if (nullite_plm(j, nt, k, np, va.base) == 1) {
172 
173  donne_lm (nz, nz-1, j, k, base, m_quant, l_quant, base_r) ;
174 
175  if (l_quant<=lmax) {
176 
177  Matrice systeme (nbre, nbre) ;
178  systeme.set_etat_qcq() ;
179 
180  for (int col=0 ; col<nbre ; col++)
181  for (int lig=0 ; lig<nbre ; lig++) {
182 
183  int facteur = (lig%2==0) ? 1 : -1 ;
184  for (int conte=0 ; conte<lig ; conte++)
185  facteur *= puis+col+conte+2*l_quant ;
186  systeme.set(lig, col) = facteur/pow(r_cont, puis+col+lig+2*l_quant) ;
187  }
188 
189  systeme.set_band(nbre, nbre) ;
190  systeme.set_lu() ;
191 
192  Tbl sec_membre (nbre) ;
193  sec_membre.set_etat_qcq() ;
194  for (int conte=0 ; conte<nbre ; conte++)
195  sec_membre.set(conte) = valeurs(conte, k, j) ;
196 
197  Tbl inv (systeme.inverse(sec_membre)) ;
198 
199  for (int conte=0 ; conte<nbre ; conte++)
200  for (int i=0 ; i<nr ; i++)
201  va.c_cf->set(nz-1, k, j, i)+=
202  inv(conte)*coef(conte+2*l_quant, i) ;
203  }
204  else for (int i=0 ; i<nr ; i++)
205  va.c_cf->set(nz-1, k, j, i)
206  = 0 ;
207  }
208 
209  va.ylm_i() ;
210  set_dzpuis (0) ;
211 }
212 
213 
214 //***************************************************************************
215 
216 void Scalar::smooth_decay(int kk, int nn) {
217 
218  assert(kk >= 0) ;
219 
220  if (etat == ETATZERO) return ;
221  if (va.get_etat() == ETATZERO) return ;
222 
223  const Mg3d& mgrid = *(mp->get_mg()) ;
224 
225  int nz = mgrid.get_nzone() ;
226  int nzm1 = nz - 1 ;
227  int np = mgrid.get_np(nzm1) ;
228  int nt = mgrid.get_nt(nzm1) ;
229  int nr_ced = mgrid.get_nr(nzm1) ;
230  int nr_shell = mgrid.get_nr(nzm1-1) ;
231 
232  assert(mgrid.get_np(nzm1-1) == np) ;
233  assert(mgrid.get_nt(nzm1-1) == nt) ;
234 
235  assert(mgrid.get_type_r(nzm1) == UNSURR) ;
236 
237  // In the present version, the mapping must be affine :
238  const Map_af* mapaf = dynamic_cast<const Map_af*>(mp) ;
239  if (mapaf == 0x0) {
240  cout << "Scalar::smooth_decay: present version supports only \n"
241  << " affine mappings !" << endl ;
242  abort() ;
243  }
244 
245 
246  assert(va.get_etat() == ETATQCQ) ;
247 
248 
249  // Spherical harmonics are required
250  va.ylm() ;
251 
252  // Array for the spectral coefficients in the CED:
253  assert( va.c_cf->t[nzm1] != 0x0) ;
254  Tbl& coefresu = *( va.c_cf->t[nzm1] ) ;
255  coefresu.set_etat_qcq() ;
256 
257  // 1-D grid
258  //----------
259  int nbr1[] = {nr_shell, nr_ced} ;
260  int nbt1[] = {1, 1} ;
261  int nbp1[] = {1, 1} ;
262  int typr1[] = {mgrid.get_type_r(nzm1-1), UNSURR} ;
263  Mg3d grid1d(2, nbr1, typr1, nbt1, SYM, nbp1, SYM) ;
264 
265  // 1-D mapping
266  //------------
267  double rbound = mapaf->get_alpha()[nzm1-1]
268  + mapaf->get_beta()[nzm1-1] ;
269  double rmin = - mapaf->get_alpha()[nzm1-1]
270  + mapaf->get_beta()[nzm1-1] ;
271  double bound1[] = {rmin, rbound, __infinity} ;
272 
273  Map_af map1d(grid1d, bound1) ;
274 
275  // 1-D scalars
276  // -----------
277  Scalar uu(map1d) ;
278  uu.std_spectral_base() ;
279  Scalar duu(map1d) ;
280  Scalar vv(map1d) ;
281  Scalar tmp(map1d) ;
282 
283  const Base_val& base = va.get_base() ;
284 
285  // Loop on the spherical harmonics
286  // -------------------------------
287  for (int k=0; k<=np; k++) {
288  for (int j=0; j<nt; j++) {
289 
290  if (nullite_plm(j, nt, k, np, base) != 1) {
291  for (int i=0 ; i<nr_ced ; i++) {
292  coefresu.set(k, j, i) = 0 ;
293  }
294  }
295  else {
296  int base_r_ced, base_r_shell , l_quant, m_quant ;
297  donne_lm(nz, nzm1, j, k, base, m_quant, l_quant, base_r_ced) ;
298  donne_lm(nz, nzm1-1, j, k, base, m_quant, l_quant, base_r_shell) ;
299 
300  int nn0 = l_quant + nn ; // lowerst power of decay
301 
302  uu.set_etat_qcq() ;
303  uu.va.set_etat_cf_qcq() ;
304  uu.va.c_cf->set_etat_qcq() ;
305  uu.va.c_cf->t[0]->set_etat_qcq() ;
306  uu.va.c_cf->t[1]->set_etat_qcq() ;
307  for (int i=0; i<nr_shell; i++) {
308  uu.va.c_cf->t[0]->set(0, 0, i) =
309  va.c_cf->operator()(nzm1-1, k, j, i) ;
310  uu.va.c_cf->t[0]->set(1, 0, i) = 0. ;
311  uu.va.c_cf->t[0]->set(2, 0, i) = 0. ;
312  }
313 
314  uu.va.c_cf->t[1]->set_etat_zero() ;
315 
316  uu.va.set_base_r(0, base_r_shell) ;
317  uu.va.set_base_r(1, base_r_ced) ;
318 
319  // Computation of the derivatives up to order k at the outer
320  // of the last shell:
321  // ----------------------------------------------------------
322  Tbl derivb(kk+1) ;
323  derivb.set_etat_qcq() ;
324  duu = uu ;
325  derivb.set(0) = uu.val_grid_point(0, 0, 0, nr_shell-1) ;
326  for (int p=1; p<=kk; p++) {
327  tmp = duu.dsdr() ;
328  duu = tmp ;
329  derivb.set(p) = duu.val_grid_point(0, 0, 0, nr_shell-1) ;
330  }
331 
332  // Matrix of the linear system to be solved
333  // ----------------------------------------
334 
335  Matrice mat(kk+1,kk+1) ;
336  mat.set_etat_qcq() ;
337 
338  for (int im=0; im<=kk; im++) {
339 
340  double fact = (im%2 == 0) ? 1. : -1. ;
341  fact /= pow(rbound, nn0 + im) ;
342 
343  for (int jm=0; jm<=kk; jm++) {
344 
345  double prod = 1 ;
346  for (int ir=0; ir<im; ir++) {
347  prod *= nn0 + jm + ir ;
348  }
349 
350  mat.set(im, jm) = fact * prod / pow(rbound, jm) ;
351 
352  }
353  }
354 
355  // Resolution of the linear system to get the coefficients
356  // of the 1/r^i functions
357  //---------------------------------------------------------
358  mat.set_band(kk+1, kk+1) ;
359  mat.set_lu() ;
360  Tbl aa = mat.inverse( derivb ) ;
361 
362  // Decaying function
363  // -----------------
364  vv = 0 ;
365  const Coord& r = map1d.r ;
366  for (int p=0; p<=kk; p++) {
367  tmp = aa(p) / pow(r, nn0 + p) ;
368  vv += tmp ;
369  }
370  vv.va.set_base( uu.va.get_base() ) ;
371 
372  // The coefficients of the decaying function
373  // are set to the result
374  //-------------------------------------------
375  vv.va.coef() ;
376 
377  if (vv.get_etat() == ETATZERO) {
378  for (int i=0; i<nr_ced; i++) {
379  coefresu.set(k,j,i) = 0 ;
380  }
381  }
382  else {
383  assert( vv.va.c_cf != 0x0 ) ;
384  for (int i=0; i<nr_ced; i++) {
385  coefresu.set(k,j,i) = vv.va.c_cf->operator()(1,0,0,i) ;
386  }
387  }
388 
389 
390  }
391 
392 
393  }
394  }
395 
396  if (va.c != 0x0) {
397  delete va.c ;
398  va.c = 0x0 ;
399  }
400  va.ylm_i() ;
401 
402  // Test of the computation
403  // -----------------------
404 
405  Scalar tmp1(*this) ;
406  Scalar tmp2(*mp) ;
407  for (int p=0; p<=kk; p++) {
408  double rd = pow(rbound, tmp1.get_dzpuis()) ;
409  double err = 0 ;
410  for (int k=0; k<np; k++) {
411  for (int j=0; j<nt; j++) {
412  double diff = fabs( tmp1.val_grid_point(nzm1, k, j, 0) / rd -
413  tmp1.val_grid_point(nzm1-1, k, j, nr_shell-1) ) ;
414  if (diff > err) err = diff ;
415  }
416  }
417  cout <<
418  "Scalar::smooth_decay: Max error matching of " << p
419  << "-th derivative: " << err << endl ;
420  tmp2 = tmp1.dsdr() ;
421  tmp1 = tmp2 ;
422  }
423 
424 
425 }
426 }
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
const double * get_alpha() const
Returns the pointer on the array alpha.
Definition: map_af.C:607
void ylm_i()
Inverse of ylm()
Definition: valeur_ylm_i.C:134
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
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
void set_lu() const
Calculate the LU-representation, assuming the band-storage has been done.
Definition: matrice.C:395
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
Tbl inverse(const Tbl &sec_membre) const
Solves the linear system represented by the matrix.
Definition: matrice.C:427
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:783
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:393
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
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:560
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base )
Definition: valeur.C:813
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:643
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
int get_etat() const
Returns the logical state.
Definition: valeur.h:760
friend Scalar pow(const Scalar &, int)
Power .
Definition: scalar_math.C:454
void smooth_decay(int k, int n)
Performs a matching of the last non-compactified shell with a decaying function where is the spher...
friend Scalar cos(const Scalar &)
Cosine.
Definition: scalar_math.C:107
Matrix handling.
Definition: matrice.h:152
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:563
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:465
Valeur va
The numerical value of the Scalar.
Definition: scalar.h:411
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: mtbl_cf.C:303
double & set(int j, int i)
Read/write of a particuliar element.
Definition: matrice.h:277
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: tbl.C:350
void set_band(int up, int low) const
Calculate the band storage of *std.
Definition: matrice.C:367
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
Multi-domain grid.
Definition: grilles.h:279
Bases of the spectral expansions.
Definition: base_val.h:325
void set_base_r(int l, int base_r)
Sets the expansion basis for r ( ) functions in a given domain.
Definition: valeur.C:839
Affine radial mapping.
Definition: map.h:2048
void raccord_c1_zec(int puis, int nbre, int lmax)
Performs the matching of the external domain with respect to the last shell using function like wit...
int etat
The logical state ETATNONDEF (undefined), ETATZERO (null), ETATUN (one), or ETATQCQ (ordinary)...
Definition: scalar.h:402
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:113
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: matrice.C:178
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:474
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:491
const Map *const mp
Mapping on which the numerical values at the grid points are defined.
Definition: tensor.h:301
void annule_hard()
Sets the Tbl to zero in a hard way.
Definition: tbl.C:375
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
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
Coord r
r coordinate centered on the grid
Definition: map.h:736
const Base_val & get_base() const
Return the bases for spectral expansions (member base )
Definition: valeur.h:490
#define R_CHEB
base de Chebychev ordinaire (fin)
Definition: type_parite.h:166