00001 /* 00002 * Copyright (c) 1999-2002 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 cfpcossin_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfpcossin.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $" ; 00024 00025 /* 00026 * Transformation de Fourier sur le premier indice d'un tableau 3-D 00027 * par le biais de la bibliotheque fftw. 00028 * 00029 * Entree: 00030 * ------- 00031 * int* deg : tableau du nombre effectif de degres de liberte dans chacune 00032 * des 3 dimensions: le nombre de points de collocation 00033 * en phi est np = deg[0] et doit etre de la forme 00034 * np = 2*p 00035 * int* dim : tableau du nombre d'elements de ff dans chacune des trois 00036 * dimensions. 00037 * On doit avoir dim[0] >= deg[0] + 2 = np + 2. 00038 * 00039 * Entree/Sortie : 00040 * --------------- 00041 * double* cf : entree: tableau des valeurs de la fonction f aux np points de 00042 * de collocation 00043 * phi_k = T k/np 0 <= k <= np-1 00044 * T etant la periode de f. La convention de stokage 00045 * utilisee est la suivante 00046 * cf[ dim[2]*dim[1]*k + dim[2]*j + i ] = f(phi_k) 0 <= k <= np-1, 00047 * les indices j et i correspondant respectivement 00048 * a theta et a r. 00049 * L'espace memoire correspondant au 00050 * pointeur cf doit etre dim[0]*dim[1]*dim[2] et doit 00051 * etre alloue avant l'appel a la routine. 00052 * 00053 * sortie: tableau des coefficients de la fonction suivant 00054 * la convention de stokage 00055 * cf[ dim[2]*dim[1]*k + dim[2]*j + i ] = c_k 0<= k <= np, 00056 * ou les indices j et i correspondent respectivement 00057 * a theta et a r et ou les c_k sont les coefficients 00058 * du developpement de f en series de Fourier: 00059 * 00060 * f(phi) = som_{l=0}^{np/2} c_{2l} cos( 2 pi/T l phi ) 00061 * + c_{2l+1} sin( 2 pi/T l phi ), 00062 * 00063 * ou T est la periode de f. 00064 * En particulier cf[1] = 0. 00065 * De plus, cf[np+1] n'est pas egal a c_{np+1} 00066 * mais a zero. 00067 * 00068 */ 00069 00070 /* 00071 * $Id: cfpcossin.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00072 * $Log: cfpcossin.C,v $ 00073 * Revision 1.1 2004/12/21 17:06:02 j_novak 00074 * Added all files for using fftw3. 00075 * 00076 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon 00077 * Suppressed the directive #include <malloc.h> for malloc is defined 00078 * in <stdlib.h> 00079 * 00080 * Revision 1.3 2002/10/16 14:36:43 j_novak 00081 * Reorganization of #include instructions of standard C++, in order to 00082 * use experimental version 3 of gcc. 00083 * 00084 * Revision 1.2 2002/09/09 13:00:39 e_gourgoulhon 00085 * Modification of declaration of Fortran 77 prototypes for 00086 * a better portability (in particular on IBM AIX systems): 00087 * All Fortran subroutine names are now written F77_* and are 00088 * defined in the new file C++/Include/proto_f77.h. 00089 * 00090 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon 00091 * LORENE 00092 * 00093 * Revision 2.0 1999/02/22 15:48:58 hyc 00094 * *** empty log message *** 00095 * 00096 * 00097 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cfpcossin.C,v 1.1 2004/12/21 17:06:02 j_novak Exp $ 00098 * 00099 */ 00100 // headers du C 00101 #include <stdlib.h> 00102 #include <fftw3.h> 00103 00104 //Lorene prototypes 00105 #include "tbl.h" 00106 00107 // Prototypage des sous-routines utilisees: 00108 fftw_plan prepare_fft(int, Tbl*&) ; 00109 //***************************************************************************** 00110 00111 void cfpcossin(const int* deg, const int* dim, double* cf) 00112 { 00113 // Dimensions du tableau cf : 00114 int n1 = dim[0] ; 00115 int n2 = dim[1] ; 00116 int n3 = dim[2] ; 00117 00118 // Nombres de degres de liberte en phi : 00119 int np = deg[0] ; 00120 00121 // Tests de dimension: 00122 if (np+2 > n1) { 00123 cout << "cfpcossin: np+2 > n1 : np+2 = " << np+2 << " , n1 = " 00124 << n1 << endl ; 00125 abort () ; 00126 exit(-1) ; 00127 } 00128 00129 // Recherche des tables 00130 Tbl* pg = 0x0 ; 00131 fftw_plan p = prepare_fft(np, pg) ; 00132 Tbl& g = *pg ; 00133 00134 int index = 0 ; 00135 int n2n3 = n2 * n3 ; 00136 double fac = 2./double(np) ; 00137 00138 // boucle sur theta et r 00139 for (int j=0; j<n2; j++) { 00140 for (int k=0; k<n3; k++) { 00141 index = n3 * j + k ; 00142 // FFT 00143 double* debut = cf + index ; 00144 double* tab = g.t ; 00145 for (int i=0; i<np; i++) { 00146 *tab = *debut ; 00147 tab++; 00148 debut += n2n3 ; 00149 } 00150 00151 fftw_execute(p) ; 00152 00153 debut = cf+index ; 00154 double* pcos = g.t ; 00155 double* psin = g.t + np - 1 ; 00156 (*debut) = (*pcos)/double(np) ; 00157 debut += n2n3 ; pcos++ ; 00158 (*debut) = 0. ; 00159 debut += n2n3 ; 00160 for (int i=1; i<np/2; i++){ 00161 *debut = (*pcos)*fac ; 00162 debut += n2n3 ; pcos++ ; 00163 *debut = -(*psin)*fac ; 00164 debut += n2n3 ; psin-- ; 00165 } 00166 (*debut) = (*pcos)/double(np) ; 00167 debut += n2n3 ; 00168 (*debut) = 0. ; 00169 00170 } // fin de la boucle sur r 00171 } // fin de la boucle sur theta 00172 00173 00174 }
1.4.6