LORENE
FFTW3/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
23char cipcossin_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cipcossin.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $" ;
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
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.3 2014/10/13 08:53:19 j_novak Exp $
84 * $Log: cipcossin.C,v $
85 * Revision 1.3 2014/10/13 08:53:19 j_novak
86 * Lorene classes and functions now belong to the namespace Lorene.
87 *
88 * Revision 1.2 2014/10/06 15:18:49 j_novak
89 * Modified #include directives to use c++ syntax.
90 *
91 * Revision 1.1 2004/12/21 17:06:02 j_novak
92 * Added all files for using fftw3.
93 *
94 * Revision 1.4 2003/01/31 10:31:23 e_gourgoulhon
95 * Suppressed the directive #include <malloc.h> for malloc is defined
96 * in <stdlib.h>
97 *
98 * Revision 1.3 2002/10/16 14:36:53 j_novak
99 * Reorganization of #include instructions of standard C++, in order to
100 * use experimental version 3 of gcc.
101 *
102 * Revision 1.2 2002/09/09 13:00:40 e_gourgoulhon
103 * Modification of declaration of Fortran 77 prototypes for
104 * a better portability (in particular on IBM AIX systems):
105 * All Fortran subroutine names are now written F77_* and are
106 * defined in the new file C++/Include/proto_f77.h.
107 *
108 * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
109 * LORENE
110 *
111 * Revision 2.0 1999/02/22 15:43:58 hyc
112 * *** empty log message ***
113 *
114 *
115 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/cipcossin.C,v 1.3 2014/10/13 08:53:19 j_novak Exp $
116 *
117 */
118
119
120// headers du C
121#include <cstdlib>
122#include <fftw3.h>
123
124//Lorene prototypes
125#include "tbl.h"
126
127// Prototypage des sous-routines utilisees:
128namespace Lorene {
129fftw_plan back_fft(int, Tbl*&) ;
130//*****************************************************************************
131
132void cipcossin(const int* deg, const int* dimc, const int* dimf,
133 double* cf, double* ff)
134{
135
136// Dimensions des tableaux ff et cf :
137 int n1f = dimf[0] ;
138 int n2f = dimf[1] ;
139 int n3f = dimf[2] ;
140 int n1c = dimc[0] ;
141 int n2c = dimc[1] ;
142 int n3c = dimc[2] ;
143
144// Nombres de degres de liberte en phi:
145 int np = deg[0] ;
146
147// Tests de dimension:
148 if (np+2 > n1c) {
149 cout << "cipcossin: np+2 > n1c : np = " << np << " , n1c = "
150 << n1c << endl ;
151 abort () ;
152 exit(-1) ;
153 }
154 if (np > n1f) {
155 cout << "cipcossin: np > n1f : np = " << np << " , n1f = "
156 << n1f << endl ;
157 abort () ;
158 exit(-1) ;
159 }
160 if (n3f > n3c) {
161 cout << "cipcossin: n3f > n3c : n3f = " << n3f << " , n3c = "
162 << n3c << endl ;
163 abort () ;
164 exit(-1) ;
165 }
166 if (n2f > n2c) {
167 cout << "cipcossin: n2f > n2c : n2f = " << n2f << " , n2c = "
168 << n2c << endl ;
169 abort () ;
170 exit(-1) ;
171 }
172
173 // Recherche des tables
174 Tbl* pg = 0x0 ;
175 fftw_plan p = back_fft(np, pg) ;
176 Tbl& g = *pg ;
177 int index ;
178 int n2n3c = n2c*n3c ;
179 int n2n3f = n2f*n3f ;
180
181// boucle sur r et theta
182
183 for (int j=0; j<n2c; j++) {
184 for (int k=0; k<n3c; k++) {
185 index = n3c * j + k ;
186 double* debut = cf + index ;
187 g.set(0) = *debut ;
188 debut += 2*n2n3c ;
189 for (int i=1; i<np/2; i++) {
190 int isym = np - i ;
191 g.set(i) = 0.5 * (*debut) ; debut += n2n3c ;
192 g.set(isym) = -0.5 * (*debut) ; debut += n2n3c ;
193 }
194 g.set(np/2) = *debut ;
195
196// FFT inverse
197
198 fftw_execute(p) ;
199
200// On range dans ff
201 if ((j<n2f) && (k<n3f)) {
202 debut = ff + n3f * j + k ;
203 for (int i=0; i<np; i++) {
204 *debut = g(i) ;
205 debut += n2n3f ;
206 }
207 }
208 } // fin de la boucle sur r
209 } // fin de la boucle sur theta
210
211}
212}
Lorene prototypes.
Definition app_hor.h:64