00001 /* 00002 * Copyright (c) 1999-2001 Eric Gourgoulhon 00003 * 00004 * This file is part of LORENE. 00005 * 00006 * LORENE is free software; you can redistribute it and/or modify 00007 * it under the terms of the GNU General Public License as published by 00008 * the Free Software Foundation; either version 2 of the License, or 00009 * (at your option) any later version. 00010 * 00011 * LORENE is distributed in the hope that it will be useful, 00012 * but WITHOUT ANY WARRANTY; without even the implied warranty of 00013 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00014 * GNU General Public License for more details. 00015 * 00016 * You should have received a copy of the GNU General Public License 00017 * along with LORENE; if not, write to the Free Software 00018 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00019 * 00020 */ 00021 00022 00023 char mat_cossinci_legi_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/mat_cossinci_legi.C,v 1.3 2005/02/18 13:14:14 j_novak Exp $" ; 00024 00025 /* 00026 * Fournit la matrice de passage pour la transformation des coefficients du 00027 * developpement en cos((2*j+1)*theta) [m pair] / sin( 2*j * theta) [m impair] 00028 * dans les coefficients du developpement en fonctions associees de Legendre 00029 * P_l^m(cos(theta)) impaires (i.e. telles que l-m est impair). 00030 * 00031 * Cette routine n'effectue le calcul de la matrice que si celui-ci n'a pas 00032 * deja ete fait, sinon elle renvoie le pointeur sur une valeur precedemment 00033 * calculee. 00034 * 00035 * Entree: 00036 * ------- 00037 * int np : Nombre de degres de liberte en phi 00038 * int nt : Nombre de degres de liberte en theta 00039 * 00040 * Sortie (valeur de retour) : 00041 * --------------------------- 00042 * double* mat_cossinci_legi : pointeur sur le tableau contenant l'ensemble 00043 * (pour les np/2+1 valeurs de m) des 00044 * matrices de passage. 00045 * La dimension du tableau est (np/2+1)*nt^2 00046 * Le stokage est le suivant: 00047 * 00048 * mat_cossinci_legi[ nt*nt* m + nt*l + j] = A_{mlj} 00049 * 00050 * ou A_{mlj} est defini par 00051 * 00052 * pour m pair : 00053 * cos((2*j+1)*theta) = som_{l=m/2}^{nt-2} A_{mlj} P_{2l+1}^m( cos(theta) ) 00054 * pour 0 <= j <= nt-2 00055 * 00056 * pour m impair : 00057 * sin(2*j*theta) = som_{l=(m+1)/2}^{nt-2} A_{mlj} P_{2l}^m( cos(theta) ) 00058 * pour 1 <= j <= nt-2 00059 * 00060 * ou P_n^m(x) represente la fonction de Legendre associee de degre n et 00061 * d'ordre m normalisee de facon a ce que 00062 * 00063 * int_0^pi [ P_n^m(cos(theta)) ]^2 sin(theta) dtheta = 1 00064 * 00065 * 00066 */ 00067 00068 /* 00069 * $Id: mat_cossinci_legi.C,v 1.3 2005/02/18 13:14:14 j_novak Exp $ 00070 * $Log: mat_cossinci_legi.C,v $ 00071 * Revision 1.3 2005/02/18 13:14:14 j_novak 00072 * Changing of malloc/free to new/delete + suppression of some unused variables 00073 * (trying to avoid compilation warnings). 00074 * 00075 * Revision 1.2 2002/10/16 14:36:54 j_novak 00076 * Reorganization of #include instructions of standard C++, in order to 00077 * use experimental version 3 of gcc. 00078 * 00079 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon 00080 * LORENE 00081 * 00082 * Revision 2.0 1999/02/22 15:35:09 hyc 00083 * *** empty log message *** 00084 * 00085 * 00086 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/mat_cossinci_legi.C,v 1.3 2005/02/18 13:14:14 j_novak Exp $ 00087 * 00088 */ 00089 00090 // headers du C 00091 #include <stdlib.h> 00092 #include <math.h> 00093 00094 // Prototypage 00095 #include "headcpp.h" 00096 #include "proto.h" 00097 00098 // Variable de loch 00099 int loch_mat_cossinci_legi = 0 ; 00100 00101 //****************************************************************************** 00102 00103 double* mat_cossinci_legi(int np, int nt) { 00104 00105 #define NMAX 30 // Nombre maximun de couples(np,nt) differents 00106 static double* tab[NMAX] ; // Tableau des pointeurs sur les tableaux 00107 static int nb_dejafait = 0 ; // Nombre de tableaux deja initialises 00108 static int np_dejafait[NMAX] ; // Valeurs de np pour lesquelles le 00109 // calcul a deja ete fait 00110 static int nt_dejafait[NMAX] ; // Valeurs de np pour lesquelles le 00111 // calcul a deja ete fait 00112 00113 int i, indice, j, j2, m, l ; 00114 00115 // #pragma critical (loch_mat_cossinci_legi) 00116 { 00117 00118 // Les matrices A_{mlj} pour ce couple (np,nt) ont-elles deja ete calculees ? 00119 indice = -1 ; 00120 for ( i=0 ; i < nb_dejafait ; i++ ) { 00121 if ( (np_dejafait[i] == np) && (nt_dejafait[i] == nt) ) indice = i ; 00122 } 00123 00124 00125 // Si le calcul n'a pas deja ete fait, il faut le faire : 00126 if (indice == -1) { 00127 if ( nb_dejafait >= NMAX ) { 00128 cout << "mat_cossinci_legi: nb_dejafait >= NMAX : " 00129 << nb_dejafait << " <-> " << NMAX << endl ; 00130 abort () ; 00131 exit(-1) ; 00132 } 00133 indice = nb_dejafait ; 00134 nb_dejafait++ ; 00135 np_dejafait[indice] = np ; 00136 nt_dejafait[indice] = nt ; 00137 00138 tab[indice] = new double[(np/2+1)*nt*nt] ; 00139 00140 //----------------------- 00141 // Preparation du calcul 00142 //----------------------- 00143 00144 // Sur-echantillonnage pour calculer les produits scalaires sans aliasing: 00145 int nt2 = 2*nt - 1 ; 00146 int nt2m1 = nt2 - 1 ; 00147 00148 int deg[3] ; 00149 deg[0] = 1 ; 00150 deg[1] = 1 ; 00151 deg[2] = nt2 ; 00152 00153 // Tableaux de travail 00154 double* yy = new double[nt2] ; 00155 double* cost = new double[nt*nt2] ; 00156 double* sint = new double[nt*nt2] ; 00157 00158 // Calcul des cos(2*j*theta) / sin( (2*j+1)*theta ) aux points de collocation 00159 // de l'echantillonnage double : 00160 00161 double dt = M_PI / double(2*(nt2-1)) ; 00162 for (j=0; j<nt-1; j++) { 00163 for (j2=0; j2<nt2; j2++) { 00164 double theta = j2*dt ; 00165 cost[nt2*j + j2] = cos( (2*j+1) * theta ) ; 00166 sint[nt2*j + j2] = sin( 2*j * theta ) ; 00167 } 00168 } 00169 00170 00171 //------------------- 00172 // Boucle sur m 00173 //------------------- 00174 00175 for (m=0; m < np/2+1 ; m++) { 00176 00177 // Recherche des fonctions de Legendre associees d'ordre m : 00178 00179 double* leg = legendre_norm(m, nt) ; 00180 00181 if (m%2==0) { 00182 // Cas m pair 00183 //----------- 00184 for (l=m/2; l<nt-1; l++) { // boucle sur les P_{2l+1}^m 00185 00186 int ll = 2*l+1 ; // degre des fonctions de Legendre 00187 00188 for (j=0; j<nt-1; j++) { // boucle sur les cos((2j+1) theta) 00189 00190 //... produit scalaire de cos((2j+1) theta) par P_{2l+1}^m(cos(theta)) 00191 00192 for (j2=0; j2<nt2; j2++) { 00193 yy[nt2m1-j2] = cost[nt2*j + j2] * 00194 leg[nt2* (ll-m) + j2] ; 00195 } 00196 00197 //....... on passe en Tchebyshev vis-a-vis de x=cos(theta) pour calculer 00198 // l'integrale (routine int1d_chebp) : 00199 cfrchebp(deg, deg, yy, deg, yy) ; 00200 tab[indice][ nt*nt* m + nt*l + j] = 00201 2.*int1d_chebp(nt2, yy) ; 00202 00203 } // fin de la boucle sur j (indice de cos((2j+1) theta) ) 00204 00205 } // fin de la boucle sur l (indice de P_{2l+1}^m) 00206 00207 00208 } // fin du cas m pair 00209 else { 00210 00211 // Cas m impair 00212 //------------- 00213 00214 for (l=(m+1)/2; l<nt-1; l++) { // boucle sur les P_{2l}^m 00215 00216 int ll = 2*l ; // degre des fonctions de Legendre 00217 00218 for (j=0; j<nt-1; j++) { // boucle sur les sin(2j theta) 00219 00220 //... produit scalaire de sin((2j+1) theta) par P_{2l+1}^m(cos(theta)) 00221 00222 for (j2=0; j2<nt2; j2++) { 00223 yy[nt2m1-j2] = sint[nt2*j + j2] * 00224 leg[nt2* (ll-m) + j2] ; 00225 } 00226 00227 //....... on passe en Tchebyshev vis-a-vis de x=cos(theta) pour calculer 00228 // l'integrale (routine int1d_chebp) : 00229 cfrchebp(deg, deg, yy, deg, yy) ; 00230 tab[indice][ nt*nt* m + nt*l + j] = 00231 2.*int1d_chebp(nt2, yy) ; 00232 00233 } // fin de la boucle sur j (indice de sin(2j theta) ) 00234 00235 } // fin de la boucle sur l (indice de P_{2l}^m) 00236 00237 00238 } // fin du cas m impair 00239 00240 delete [] leg ; 00241 00242 } // fin de la boucle sur m 00243 00244 // Liberation espace memoire 00245 // ------------------------- 00246 00247 delete [] yy ; 00248 delete [] cost ; 00249 delete [] sint ; 00250 00251 } // fin du cas ou le calcul etait necessaire 00252 00253 } // Fin de zone critique 00254 00255 return tab[indice] ; 00256 00257 } 00258 00259
1.4.6