source: Sophya/trunk/SophyaExt/IFFTW/fftw3server.cc@ 3009

Last change on this file since 3009 was 3000, checked in by ansari, 19 years ago

Passage a FFTW3 + suppression const des arguments FFTForward/Backward cmv+Reza 03/07/2006

File size: 6.7 KB
Line 
1#include "FFTW/fftw3.h"
2
3
4static void FillSizes4FFTW(const BaseArray & in, int * sz)
5{
6 int k1 = 0;
7 int k2 = 0;
8 for(k1=in.NbDimensions()-1; k1>=0; k1--) {
9 sz[k2] = in.Size(k1); k2++;
10 }
11}
12
13/* --Methode-- */
14//! Constructor - If preserve_input==true, input arrays will not be overwritten.
15FFTWServer::FFTWServer(bool preserve_input)
16 : FFTServerInterface("FFTServer using FFTW3 package")
17 , ckR8("FFTWServer - ", true, false),
18 _preserve_input(preserve_input)
19{
20}
21
22
23/* --Methode-- */
24FFTWServer::~FFTWServer()
25{
26}
27
28/* --Methode-- */
29FFTServerInterface * FFTWServer::Clone()
30{
31 return (new FFTWServer) ;
32}
33
34/* --Methode-- */
35void
36FFTWServer::FFTForward(TArray< complex<r_8> > & in, TArray< complex<r_8> > & out)
37{
38 int rank = ckR8.CheckResize(in, out);
39 if (rank == 1) { // One dimensional transform
40 fftw_plan p = fftw_plan_dft_1d(in.Size(), (fftw_complex *)in.Data(),
41 (fftw_complex *)out.Data(),
42 FFTW_FORWARD, FFTW_ESTIMATE);
43 if (p == NULL)
44 throw ParmError("FFTWServer::FFTForward( complex<r_8>, complex<r_8> ) Error creating fftw_plan");
45 fftw_execute(p);
46 fftw_destroy_plan(p);
47 }
48 else { // Multi dimensional
49 if (in.NbDimensions() > MAXND_FFTW)
50 throw ParmError("FFTWServer::FFTForward( complex<r_8>, complex<r_8> ) rank > MAXND_FFTW !");
51 int sz[MAXND_FFTW];
52 FillSizes4FFTW(in, sz);
53 fftw_plan p = fftw_plan_dft(rank, sz,
54 (fftw_complex *)in.Data(), (fftw_complex *)out.Data(),
55 FFTW_FORWARD, FFTW_ESTIMATE);
56 if (p == NULL)
57 throw ParmError("FFTWServer::FFTForward( complex<r_8>, complex<r_8> ) Error creating fftw_plan");
58 fftw_execute(p);
59 fftw_destroy_plan(p);
60 }
61 if(this->getNormalize()) out=out/complex<r_8>((double)in.Size(),0.);
62 return;
63}
64
65/* --Methode-- */
66void FFTWServer::FFTBackward(TArray< complex<r_8> > & in, TArray< complex<r_8> > & out)
67{
68 int rank = ckR8.CheckResize(in, out);
69 if (rank == 1) { // One dimensional transform
70 fftw_plan p = fftw_plan_dft_1d(in.Size(), (fftw_complex *)in.Data(),
71 (fftw_complex *)out.Data(),
72 FFTW_BACKWARD, FFTW_ESTIMATE);
73 if (p == NULL)
74 throw ParmError("FFTWServer::FFTBackward( complex<r_8>, complex<r_8> ) Error creating fftw_plan");
75 fftw_execute(p);
76 fftw_destroy_plan(p);
77 }
78 else { // Multi dimensional
79 if (in.NbDimensions() > MAXND_FFTW)
80 throw ParmError("FFTWServer::FFTBackward( complex<r_8>, complex<r_8> ) rank > MAXND_FFTW !");
81 int sz[MAXND_FFTW];
82 FillSizes4FFTW(in, sz);
83 fftw_plan p = fftw_plan_dft(rank, sz,
84 (fftw_complex *)in.Data(), (fftw_complex *)out.Data(),
85 FFTW_BACKWARD, FFTW_ESTIMATE);
86 if (p == NULL)
87 throw ParmError("FFTWServer::FFTBackward( complex<r_8>, complex<r_8> ) Error creating fftw_plan");
88 fftw_execute(p);
89 fftw_destroy_plan(p);
90 }
91
92 return;
93}
94
95
96void FFTWServer::FFTForward(TArray< r_8 > & in, TArray< complex<r_8> > & out)
97{
98 int rank = ckR8.CheckResize(in, out);
99 if (rank == 1) { // One dimensional transform
100 fftw_plan p = fftw_plan_dft_r2c_1d(in.Size(), in.Data(),
101 (fftw_complex *)out.Data(),
102 FFTW_ESTIMATE);
103 if (p == NULL)
104 throw ParmError("FFTWServer::FFTForward(r_8, complex<r_8> ) Error creating fftw_plan");
105 fftw_execute(p);
106 fftw_destroy_plan(p);
107
108
109 }
110 else { // Multi dimensional
111 if (in.NbDimensions() > MAXND_FFTW)
112 throw ParmError("FFTWServer::FFTForward(r_8, complex<r_8> ) rank > MAXND_FFTW !");
113 int sz[MAXND_FFTW];
114 FillSizes4FFTW(in, sz);
115 fftw_plan p = fftw_plan_dft_r2c(rank, sz, in.Data(),
116 (fftw_complex *)out.Data(),
117 FFTW_ESTIMATE);
118 if (p == NULL)
119 throw ParmError("FFTWServer::FFTForward(r_8, complex<r_8> ) Error creating fftw_plan");
120 fftw_execute(p);
121 fftw_destroy_plan(p);
122 }
123 if(this->getNormalize()) out=out/complex<r_8>((double)in.Size(),0.);
124 return;
125
126}
127
128
129
130void FFTWServer::FFTBackward(TArray< complex<r_8> > & in, TArray< r_8 > & out,
131 bool usoutsz)
132{
133 // ATTENTION dans les TF (Complex->Reel), c'est la taille logique de la TF
134 // qu'il faut indiquer lors de la creation des plans, cad taille tableau reel
135 int rank = ckR8.CheckResize(in, out, usoutsz);
136 bool share = (_preserve_input) ? false : true;
137 TArray< complex<r_8> > inp(in, share);
138 if (rank == 1) { // One dimensional transform
139 fftw_plan p = fftw_plan_dft_c2r_1d(out.Size(), (fftw_complex *)inp.Data(),
140 out.Data(),
141 FFTW_ESTIMATE);
142 if (p == NULL)
143 throw ParmError("FFTWServer::FFTBackward(r_8, complex<r_8> ) Error creating fftw_plan");
144 fftw_execute(p);
145 fftw_destroy_plan(p);
146 }
147 else { // Multi dimensional
148 if (in.NbDimensions() > MAXND_FFTW)
149 throw ParmError("FFTWServer::FFTBackward(r_8, complex<r_8> ) rank > MAXND_FFTW !");
150 int sz[MAXND_FFTW];
151 FillSizes4FFTW(out, sz);
152 fftw_plan p = fftw_plan_dft_c2r(rank, sz, (fftw_complex *)inp.Data(),
153 out.Data(),
154 FFTW_ESTIMATE);
155 if (p == NULL)
156 throw ParmError("FFTWServer::FFTBackward(r_8, complex<r_8> ) Error creating fftw_plan");
157 fftw_execute(p);
158 fftw_destroy_plan(p);
159 }
160 return;
161}
162
163
164
165/* --Methode-- */
166void FFTWServer::ReShapetoReal(TArray< complex<r_8> > const & ina, TArray< r_8 > & outa,
167 bool usz)
168{
169 TVector< complex<r_8> > in(ina);
170 TVector< r_8> out(outa);
171 int n = in.NElts();
172 r_8 thr = FFTArrayChecker<r_8>::ZeroThreshold();
173 sa_size_t ncs;
174 if (usz) {
175 if ( (out.NElts() != 2*n-1) && (out.NElts() != 2*n-2) )
176 throw SzMismatchError("FFTWServer::ReShapetoReal(..., true) - Wrong output array size ");
177 ncs = out.NElts();
178 }
179 else {
180 ncs = ( (in(n-1).imag() < -thr) || (in(n-1).imag() > thr) ) ?
181 ncs = 2*n-1 : ncs = 2*n-2;
182
183 if (out.NElts() != ncs) {
184 cerr << " DEBUG-FFTWServer::ReShapetoReal() ncs= " << ncs
185 << " out.NElts()= " << out.NElts() << endl;
186 throw SzMismatchError("FFTWServer::ReShapetoReal() - Wrong output array size !");
187 }
188 }
189 // if (ncs == 2*n-2) {
190 out(0) = in(0).real();
191 for(int k=1; k<n; k++) {
192 out(k) = in(k).real();
193 out(ncs-k) = in(k).imag();
194 }
195 if (ncs == 2*n-2)
196 out(n-1) = in(n-1).real();
197
198 // for(int k=0; k<out.NElts(); k++) cout << "ReShapetoReal out " << k << " " << out(k) << endl;
199}
200
201
202/* --Methode-- */
203void FFTWServer::ReShapetoCompl(TArray< r_8 > const & ina, TArray< complex<r_8> > & outa)
204{
205 TVector< r_8> in(ina);
206 TVector< complex<r_8> > out(outa);
207
208 sa_size_t n = in.NElts();
209 sa_size_t ncs = n/2+1;
210 sa_size_t nc = (n%2 != 0) ? n/2+1 : n/2;
211 if (out.NElts() != ncs)
212 throw SzMismatchError("FFTWServer::ReShapetoCompl() - Wrong output array size !");
213
214 out(0) = complex<r_8> (in(0),0.);
215 for(int k=1; k<ncs; k++)
216 out(k) = complex<r_8>(in(k), in(n-k));
217 if (n%2 == 0) out(ncs-1) = complex<r_8>(in(n/2), 0.);
218
219}
220
Note: See TracBrowser for help on using the repository browser.