LORENE
FFT991/cipcossin.C
1 /*
2  * Copyright (c) 1999-2002 Eric Gourgoulhon
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 /*
27  * Transformation de Fourier inverse sur le premier indice d'un tableau 3-D
28  * par le biais de la routine FFT Fortran FFT991
29  *
30  * Entree:
31  * -------
32  * int* deg : tableau du nombre effectif de degres de liberte dans chacune
33  * des 3 dimensions: le nombre de points de collocation
34  * en phi est np = deg[0] et doit etre de la forme
35  * np = 2^p 3^q 5^r
36  * int* dimc : tableau du nombre d'elements dans chacune des trois
37  * dimensions du tableau cf
38  * On doit avoir dimc[0] >= deg[0] + 2 = np + 2.
39  *
40  * int* dimf : tableau du nombre d'elements dans chacune des trois
41  * dimensions du tableau ff
42  * On doit avoir dimf[0] >= deg[0] = np .
43  *
44  *
45  * Entree / sortie :
46  * -----------------
47  * double* cf : entree: tableau des coefficients de la fonction f;
48  * L'espace memoire correspondant a ce
49  * pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit
50  * avoir ete alloue avant l'appel a la routine.
51  * La convention de stokage est la suivante:
52  * cf[ dimc[2]*dimc[1]*k + dimc[2]*j + i ] = c_k 0<= k <= np,
53  * ou les indices j et i correspondent respectivement
54  * a theta et a r et ou les c_k sont les coefficients
55  * du developpement de f en series de Fourier:
56  *
57  * f(phi) = som_{l=0}^{np/2} c_{2l} cos( 2 pi/T l phi )
58  * + c_{2l+1} sin( 2 pi/T l phi ),
59  *
60  * ou T est la periode de f.
61  * En particulier cf[1] = 0.
62  * De plus, cf[np+1] n'est pas egal a c_{np+1}
63  * mais a zero.
64  * !!!! CE TABLEAU EST DETRUIT EN SORTIE !!!!!
65  *
66  * Sortie:
67  * -------
68  * double* ff : tableau des valeurs de la fonction aux points de
69  * de collocation
70  *
71  * phi_k = T k/np 0 <= k <= np-1
72  *
73  * suivant la convention de stokage:
74  * ff[ dimf[2]*dimf[1]*k + dimf[2]*j + i ] = f(phi_k) 0 <= k <= np-1,
75  * les indices j et i correspondant respectivement
76  * a theta et a r.
77  * L'espace memoire correspondant a ce
78  * pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit
79  * avoir ete alloue avant l'appel a la routine.
80  */
81 
82 /*
83  * $Id: cipcossin.C,v 1.5 2016/12/05 16:18:03 j_novak Exp $
84  * $Log: cipcossin.C,v $
85  * Revision 1.5 2016/12/05 16:18:03 j_novak
86  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
87  *
88  * Revision 1.4 2014/10/15 12:48:21 j_novak
89  * Corrected namespace declaration.
90  *
91  * Revision 1.3 2014/10/13 08:53:16 j_novak
92  * Lorene classes and functions now belong to the namespace Lorene.
93  *
94  * Revision 1.2 2014/10/06 15:18:45 j_novak
95  * Modified #include directives to use c++ syntax.
96  *
97  * Revision 1.1 2004/12/21 17:06:01 j_novak
98  * Added all files for using fftw3.
99  *
100  * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
101  * Suppressed the directive #include <malloc.h> for malloc is defined
102  * in <stdlib.h>
103  *
104  * Revision 1.3 2002/10/16 14:36:53 j_novak
105  * Reorganization of #include instructions of standard C++, in order to
106  * use experimental version 3 of gcc.
107  *
108  * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
109  * Modification of declaration of Fortran 77 prototypes for
110  * a better portability (in particular on IBM AIX systems):
111  * All Fortran subroutine names are now written F77_* and are
112  * defined in the new file C++/Include/proto_f77.h.
113  *
114  * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
115  * LORENE
116  *
117  * Revision 2.0 1999/02/22 15:43:58 hyc
118  * *** empty log message ***
119  *
120  *
121  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFT991/cipcossin.C,v 1.5 2016/12/05 16:18:03 j_novak Exp $
122  *
123  */
124 
125 
126 // headers du C
127 #include <cassert>
128 #include <cstdlib>
129 
130 // Prototypes of F77 subroutines
131 #include "headcpp.h"
132 #include "proto_f77.h"
133 
134 // Prototypage des sous-routines utilisees:
135 namespace Lorene {
136 int* facto_ini(int ) ;
137 double* trigo_ini(int ) ;
138 //*****************************************************************************
139 
140 void cipcossin(const int* deg, const int* dimc, const int* dimf,
141  double* cf, double* ff)
142 {
143 
144 int i, j, k, index ;
145 
146 // Dimensions des tableaux ff et cf :
147  int n1f = dimf[0] ;
148  int n2f = dimf[1] ;
149  int n3f = dimf[2] ;
150  int n1c = dimc[0] ;
151  int n2c = dimc[1] ;
152  int n3c = dimc[2] ;
153 
154 // Nombres de degres de liberte en phi:
155  int np = deg[0] ;
156 
157 // Tests de dimension:
158  if (np+2 > n1c) {
159  cout << "cipcossin: np+2 > n1c : np = " << np << " , n1c = "
160  << n1c << endl ;
161  abort () ;
162  exit(-1) ;
163  }
164  if (np > n1f) {
165  cout << "cipcossin: np > n1f : np = " << np << " , n1f = "
166  << n1f << endl ;
167  abort () ;
168  exit(-1) ;
169  }
170  if (n3f > n3c) {
171  cout << "cipcossin: n3f > n3c : n3f = " << n3f << " , n3c = "
172  << n3c << endl ;
173  abort () ;
174  exit(-1) ;
175  }
176  if (n2f > n2c) {
177  cout << "cipcossin: n2f > n2c : n2f = " << n2f << " , n2c = "
178  << n2c << endl ;
179  abort () ;
180  exit(-1) ;
181  }
182 
183  // Recherche des tables
184  int* facto = facto_ini(np) ;
185  double* trigo = trigo_ini(np) ;
186 
187  // Tableau de travail
188  double* t1 = (double*)( malloc( (np+2)*sizeof(double) ) ) ;
189 
190 // Denormalisation des cosinus
191  int n2n3c = n2c * n3c ;
192  for (i=2; i<np; i += 2 ) {
193  for (j=0; j<n2c; j++) {
194  for (k=0; k<n3c; k++) {
195  index = n2n3c * i + n3c * j + k ;
196  cf[index] *= .5 ;
197  }
198  }
199  }
200 
201 // Normalisation des sinus (les termes sin(0*phi) et sin(np/2 *phi) ne sont pas
202 // traites)
203  for (i=3; i<np+1; i += 2 ) {
204  for (j=0; j<n2c; j++) {
205  for (k=0; k<n3c; k++) {
206  index = n2n3c * i + n3c * j + k ;
207  cf[index] *= -.5 ;
208  }
209  }
210  }
211 
212 // Parametres pour la routine FFT991
213  int jump = 1 ;
214  int inc = n2n3c ;
215  int lot = 1 ;
216  int isign = 1 ;
217 
218 // boucle sur r et theta
219 
220  for (j=0; j<n2c; j++) {
221  for (k=0; k<n3c; k++) {
222 
223  index = n3c * j + k ;
224 
225 // FFT inverse
226  double* debut = cf + index ;
227 
228  F77_fft991( debut, t1, trigo, facto, &inc, &jump, &np,
229  &lot, &isign) ;
230  } // fin de la boucle sur r
231  } // fin de la boucle sur theta
232 
233 // On recopie le resultat dans ff si besoin est (c'est-a-dire si le
234 // pointeur ff est different de cf) :
235 
236  if ( ff != cf ) {
237  int n2n3f = n2f * n3f ;
238  for (i=0; i<np; i++) {
239  for (j=0; j<n2f; j++) {
240  for (k=0; k<n3f; k++) {
241  int indexc = n2n3c * i + n3c * j + k ;
242  int indexf = n2n3f * i + n3f * j + k ;
243  ff[indexf] = cf[indexc] ;
244  }
245  }
246  }
247 
248  }
249 
250  // Menage
251  free (t1) ;
252 
253 }
254 }
Lorene prototypes.
Definition: app_hor.h:67