LORENE
mat_cossinc_leg.C
1 /*
2  * Copyright (c) 2004 Michael Forot
3  *
4  * This file is part of LORENE.
5  *
6  * LORENE is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * LORENE is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with LORENE; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  *
20  */
21 
22 
23 
24 
25 /*
26  * Fournit la matrice de passage pour la transformation des coefficients du
27  * developpement en cos(j*theta) [m pair] / sin( j* theta) [m impair]
28  * dans les coefficients du developpement en fonctions associees de Legendre
29  * P_l^m(cos(theta)).
30  *
31  * Cette routine n'effectue le calcul de la matrice que si celui-ci n'a pas
32  * deja ete fait, sinon elle renvoie le pointeur sur une valeur precedemment
33  * calculee.
34  *
35  * Entree:
36  * -------
37  * int np : Nombre de degres de liberte en phi
38  * int nt : Nombre de degres de liberte en theta
39  *
40  * Sortie (valeur de retour) :
41  * ---------------------------
42  * double* mat_cossinc_leg : pointeur sur le tableau contenant l'ensemble
43  * (pour les np/2+1 valeurs de m) des
44  * matrices de passage.
45  * La dimension du tableau est (np/2+1)*nt^2
46  * Le stokage est le suivant:
47  *
48  * mat_cossinc_leg[ nt*nt* m + nt*l + j] = A_{mlj}
49  *
50  * ou A_{mlj} est defini par
51  *
52  * pour m pair :
53  * cos(j*theta) = som_{l=m}^{nt-1} A_{mlj} P_{l}^m( cos(theta) )
54  *
55  * pour m impair :
56  * sin(j*theta) = som_{l=m}^{nt-1} A_{mlj} P_{l}^m( cos(theta) )
57  *
58  * ou P_n^m(x) represente la fonction de Legendre associee de degre n et
59  * d'ordre m normalisee de facon a ce que
60  *
61  * int_0^pi [ P_n^m(cos(theta)) ]^2 sin(theta) dtheta = 1
62  *
63  *
64  */
65 
66 /*
67  * $Id: mat_cossinc_leg.C,v 1.8 2016/12/05 16:18:02 j_novak Exp $
68  * $Log: mat_cossinc_leg.C,v $
69  * Revision 1.8 2016/12/05 16:18:02 j_novak
70  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
71  *
72  * Revision 1.7 2014/10/13 08:53:13 j_novak
73  * Lorene classes and functions now belong to the namespace Lorene.
74  *
75  * Revision 1.6 2014/10/06 15:16:02 j_novak
76  * Modified #include directives to use c++ syntax.
77  *
78  * Revision 1.5 2009/10/13 13:49:36 j_novak
79  * New base T_LEG_MP.
80  *
81  * Revision 1.4 2005/02/18 13:14:13 j_novak
82  * Changing of malloc/free to new/delete + suppression of some unused variables
83  * (trying to avoid compilation warnings).
84  *
85  * Revision 1.3 2005/02/16 15:24:15 m_forot
86  * Replace int1d_chebp by int1d_cheb
87  *
88  * Revision 1.2 2004/12/17 15:42:02 e_gourgoulhon
89  * l_max = nt instead of nt2.
90  *
91  * Revision 1.1 2004/11/23 15:13:50 m_forot
92  * Added the bases for the cases without any equatorial symmetry
93  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
94  *
95  *
96  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/mat_cossinc_leg.C,v 1.8 2016/12/05 16:18:02 j_novak Exp $
97  *
98  */
99 
100 // headers du C
101 #include <cstdlib>
102 #include <cmath>
103 
104 // Prototypage
105 #include "headcpp.h"
106 #include "proto.h"
107 
108 namespace Lorene {
109 //******************************************************************************
110 
111 double* mat_cossinc_leg(int np, int nt) {
112 
113 #define NMAX 30 // Nombre maximun de couples(np,nt) differents
114 static double* tab[NMAX] ; // Tableau des pointeurs sur les tableaux
115 static int nb_dejafait = 0 ; // Nombre de tableaux deja initialises
116 static int np_dejafait[NMAX] ; // Valeurs de np pour lesquelles le
117  // calcul a deja ete fait
118 static int nt_dejafait[NMAX] ; // Valeurs de np pour lesquelles le
119  // calcul a deja ete fait
120 
121 int i, indice, j, j2, m, l ;
122 
123  {
124 
125  // Les matrices A_{mlj} pour ce couple (np,nt) ont-elles deja ete calculees ?
126  indice = -1 ;
127  for ( i=0 ; i < nb_dejafait ; i++ ) {
128  if ( (np_dejafait[i] == np) && (nt_dejafait[i] == nt) ) indice = i ;
129  }
130 
131 
132  // Si le calcul n'a pas deja ete fait, il faut le faire :
133  if (indice == -1) {
134  if ( nb_dejafait >= NMAX ) {
135  cout << "mat_cossinc_leg: nb_dejafait >= NMAX : "
136  << nb_dejafait << " <-> " << NMAX << endl ;
137  abort () ;
138  exit(-1) ;
139  }
140  indice = nb_dejafait ;
141  nb_dejafait++ ;
142  np_dejafait[indice] = np ;
143  nt_dejafait[indice] = nt ;
144 
145  tab[indice] = new double[(np/2+1)*nt*nt] ;
146 
147 //-----------------------
148 // Preparation du calcul
149 //-----------------------
150 
151 // Sur-echantillonnage pour calculer les produits scalaires sans aliasing:
152  int nt2 = 2*nt - 1 ;
153  int nt2m1 = nt2 - 1 ;
154 
155  int deg[3] ;
156  deg[0] = 1 ;
157  deg[1] = 1 ;
158  deg[2] = nt2 ;
159 
160 // Tableaux de travail
161  double* yy = new double[nt2] ;
162  double* cost = new double[nt*nt2] ;
163  double* sint = new double[nt*nt2] ;
164 
165 // Calcul des cos(j*theta) / sin( j*theta ) aux points de collocation
166 // de l'echantillonnage double :
167 
168  double dt = M_PI / double((nt2-1)) ;
169  for (j=0; j<nt; j++) {
170  for (j2=0; j2<nt2; j2++) {
171  double theta = j2*dt ;
172  cost[nt2*j + j2] = cos( j * theta ) ;
173  sint[nt2*j + j2] = sin( j * theta ) ;
174  }
175  }
176 
177 
178 //-------------------
179 // Boucle sur m
180 //-------------------
181 
182  for (m=0; m < np/2+1 ; m++) {
183 
184 // Recherche des fonctions de Legendre associees d'ordre m :
185 
186  double* leg = legendre_norm(m, nt) ;
187 
188  if (m%2==0) {
189 // Cas m pair
190 //-----------
191  for (l=m; l<nt; l++) { // boucle sur les P_{l}^m
192 
193  int parite = 1 - 2*((l-m)%2) ; // parite du P_l^m par rapport au plan theta=pi/2
194 
195  for (j=0; j<nt; j++) { // boucle sur les cos(j theta)
196 
197 //... produit scalaire de cos(j theta) par P_{l}^m(cos(theta))
198 
199  for (j2=0; j2<nt; j2++) {
200  yy[nt2m1-j2] = cost[nt2*j + j2] *
201  leg[nt2* (l-m) + 2*j2] ;
202  }
203 
204  for (j2=nt; j2<nt2; j2++) {
205  yy[nt2m1-j2] = cost[nt2*j + j2] *
206  parite * leg[nt2* (l-m) + 2*nt2 -2 -2*j2] ;
207  }
208 
209 //....... on passe en Tchebyshev vis-a-vis de x=cos(theta) pour calculer
210 // l'integrale (routine int1d_cheb) :
211  cfrcheb(deg, deg, yy, deg, yy) ;
212  tab[indice][ nt*nt* m + nt*l + j] =
213  int1d_cheb(nt2, yy) ;
214 
215  } // fin de la boucle sur j (indice de cos(j theta) )
216 
217  } // fin de la boucle sur l (indice de P_l^m)
218 
219 
220  } // fin du cas m pair
221  else {
222 
223 // Cas m impair
224 //-------------
225 
226  for (l=m; l<nt; l++) { // boucle sur les P_{l}^m
227 
228  int parite = 1 - 2*((l-m)%2) ; // parite du P_l^m par rapport au plan theta=pi/2
229 
230  for (j=0; j<nt; j++) { // boucle sur les sin(j)theta)
231 
232 //... produit scalaire de sin(j) theta) par P_{l}^m(cos(theta))
233 
234  for (j2=0; j2<nt; j2++) {
235  yy[nt2m1-j2] = sint[nt2*j + j2] *
236  leg[nt2* (l-m) + 2*j2] ;
237 
238  }
239 
240  for (j2=nt; j2<nt2; j2++) {
241  yy[nt2m1-j2] = sint[nt2*j + j2] *
242  parite * leg[nt2* (l-m) + 2*nt2 -2 - 2*j2] ;
243 
244  }
245 
246 //....... on passe en Tchebyshev vis-a-vis de x=cos(theta) pour calculer
247 // l'integrale (routine int1d_chebp) :
248  cfrcheb(deg, deg, yy, deg, yy) ;
249 
250  tab[indice][ nt*nt* m + nt*l + j] =
251  int1d_cheb(nt2, yy) ;
252 
253  } // fin de la boucle sur j (indice de sin(j*theta) )
254 
255  } // fin de la boucle sur l (indice de P_{l}^m)
256 
257 
258  } // fin du cas m impair
259 
260  delete [] leg ;
261 
262  } // fin de la boucle sur m
263 
264 // Liberation espace memoire
265 // -------------------------
266 
267  delete [] yy ;
268  delete [] cost ;
269  delete [] sint ;
270 
271  } // fin du cas ou le calcul etait necessaire
272 
273  } // Fin de zone critique
274 
275  return tab[indice] ;
276 
277 }
278 
279 
280 }
Lorene prototypes.
Definition: app_hor.h:67
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:97
Cmp sin(const Cmp &)
Sine.
Definition: cmp_math.C:72