File: itkFFTTest.h

package info (click to toggle)
insighttoolkit4 4.13.3withdata-dfsg2-4
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 491,256 kB
  • sloc: cpp: 557,600; ansic: 180,546; fortran: 34,788; python: 16,572; sh: 2,187; lisp: 2,070; tcl: 993; java: 362; perl: 200; makefile: 133; csh: 81; pascal: 69; xml: 19; ruby: 10
file content (451 lines) | stat: -rw-r--r-- 15,371 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
/*=========================================================================
 *
 *  Copyright Insight Software Consortium
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         http://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkFFTTest_h
#define itkFFTTest_h

/* This test is build for testing forward and inverse Fast Fourier Transforms
 * using VNL and FFTW FFT libraries. */
#include "itkConfigure.h"
#include "itkImage.h"
#include "itkImageRegionIterator.h"
#include "itkVnlForwardFFTImageFilter.h"
#include "itkVnlInverseFFTImageFilter.h"
#if defined(ITK_USE_FFTWF) || defined(ITK_USE_FFTWD)
#include "itkFFTWInverseFFTImageFilter.h"
#include "itkFFTWForwardFFTImageFilter.h"
#endif
#include "itksys/SystemTools.hxx"
#include "vnl/vnl_sample.h"
#include <cmath>

/* test_fft is the test function and it is templated over the pixel, Image
 * dimensions and the FFT library to be used. */
template< typename TPixel, unsigned int VImageDimensions,
          typename R2CType, typename C2RType >
int
test_fft(unsigned int *SizeOfDimensions)
{
  typedef itk::Image< TPixel, VImageDimensions >                  RealImageType;
  typedef itk::Image< std::complex< TPixel >, VImageDimensions >  ComplexImageType;
  unsigned int counter = 0;
  typename RealImageType::SizeType  imageSize;
  typename RealImageType::IndexType imageIndex;

  typedef typename RealImageType::IndexType::IndexValueType indexValueType;

  // We are testing the FFT for 1D, 2D, and 3D images. An array
  // (SizeOfDimensions) containing the sizes of each dimension is
  // passed as an argument to this function. Based on the template
  // argument VImageDimensions, we create a 1D, 2D, or 3D image by
  // selecting the sizes of image dimensions from this array.
  for (unsigned int i = 0; i < VImageDimensions; i++)
    {
    imageSize.SetElement( i, SizeOfDimensions[i] );
    imageIndex.SetElement( i, -2*static_cast<indexValueType>(i) ); // Test for handling non-zero
                                      // image indices correctly
    }

  typename RealImageType::RegionType region;
  region.SetSize( imageSize );
  region.SetIndex( imageIndex );

  typename RealImageType::Pointer realImage = RealImageType::New();
  // Create the Real Image.
  realImage->SetLargestPossibleRegion( region );
  realImage->SetBufferedRegion( region );
  realImage->SetRequestedRegion( region );
  realImage->Allocate();

  // Set up spacing and origin to test passing of metadata
  typename RealImageType::PointType     origin;
  typename RealImageType::SpacingType   spacing;
  typename RealImageType::DirectionType direction;
  direction.Fill( 0.0 );
  for ( unsigned int i = 0; i < VImageDimensions; i++ )
    {
    origin[i]  = static_cast< typename RealImageType::PointValueType >( i ) + 1.0;
    spacing[i] = static_cast< typename RealImageType::SpacingValueType >( i ) + 2.0;
    direction[i][i] = static_cast< typename RealImageType::DirectionType::ValueType > ( i ) + 3.0;
    }
  realImage->SetOrigin( origin );
  realImage->SetSpacing( spacing );
  realImage->SetDirection( direction );
  vnl_sample_reseed( static_cast< int >( 123456 ) );

  // We use 2 region iterators for this test: the original image
  // iterator and another iterator for the resultant image after
  // performing FFT and IFFT.
  itk::ImageRegionIterator< RealImageType > originalImageIterator( realImage, region );

  // Allocate random pixel values to the image by iterating through it
  // and print out the image data.
  try
    {
    std::cerr << "---- Original image ----" << std::endl;
    while( !originalImageIterator.IsAtEnd() )
      {
      TPixel val = vnl_sample_uniform( 0.0, 16384.0 );
      if ( (counter + 1 ) % SizeOfDimensions[0] == 0 )
        {
        std::cout << val << std::endl;
        }
      else
        {
        std::cout << val << " ";
        }
      counter++;
      originalImageIterator.Set( val );
      ++originalImageIterator;
      }
    std::cout << std::endl << std::endl;
    }
  catch ( itk::ExceptionObject & ex )
    {
    ex.Print( std::cerr );
    return -1;
    }

  // Real to complex pointer. This computes the forward FFT.
  typename R2CType::Pointer R2C = R2CType::New();
  // Complex to Real pointer. This computes the Inverse FFT.
  typename C2RType::Pointer C2R = C2RType::New();

  // Set the real image created as the input to the forward FFT
  // filter.
  R2C->SetInput( realImage );
  R2C->Print( std::cout );

  try
    {
    R2C->Update();
    }
  catch ( itk::ExceptionObject & ex )
    {
    ex.Print( std::cerr );
    return -1;
    }

  // Get the size and the pointer to the complex image.
  typename ComplexImageType::Pointer complexImage = R2C->GetOutput();
  std::complex< TPixel > *fftbuf = complexImage->GetBufferPointer();
  const typename ComplexImageType::SizeType & complexImageSize =
    complexImage->GetLargestPossibleRegion().GetSize();

  unsigned int sizes[4] = { 1,1,1,1 };
  for( unsigned int i = 0; i < VImageDimensions; i++)
    {
    sizes[i] = complexImageSize[i];
    }
  /* Print out the the frequency domain data obtained after performing
   * the forward transform. */
  std::cout << "Frequency domain data after forward transform:" << std::endl;
  for( unsigned int i = 0; i < sizes[2]; i++)
    {
    unsigned int zStride = i * sizes[1] * sizes[0];
    for (unsigned int j = 0; j < sizes[1]; j++)
      {
      unsigned int yStride = j * sizes[0];
      for (unsigned int k = 0; k < sizes[0]; k++)
        {
        std::cout << fftbuf[zStride+yStride+k] << " ";
        }
      std::cout << std::endl;
      }
    }

  std::cout << std::endl << std::endl;

  // Check to see that the metadata has been copied
  if ( complexImage->GetOrigin() != origin )
    {
    std::cerr << "Origin in R2C output does not match!" << std::endl;
    return -1;
    }

  if ( complexImage->GetSpacing() != spacing )
    {
    std::cerr << "Spacing in R2C output does not match!" << std::endl;
    return -1;
    }

  if ( complexImage->GetDirection() != direction )
    {
    std::cerr << "Direction in R2C output does not match!" << std::endl;
    return -1;
    }

  // Perform the Inverse FFT to get back the Real Image. C2R is the
  // complex conjugate to real image filter and we give the resulting
  // complex image as input to this filter. This is the Inverse FFT of
  // the image.
  C2R->SetInput( complexImage );

  // Inform the filter that there's an odd # of pixels in the x
  // dimension.
  C2R->Print( std::cout );
  C2R->Update();
  std::cerr << "C2R region: " << C2R->GetOutput()->GetLargestPossibleRegion() << std::endl;
  typename RealImageType::Pointer imageAfterInverseFFT = C2R->GetOutput();

  // Check to see that the metadata has been copied
  if ( imageAfterInverseFFT->GetOrigin() != origin )
    {
    std::cerr << "Origin in C2R output does not match!" << std::endl;
    return -1;
    }

  if ( imageAfterInverseFFT->GetSpacing() != spacing )
    {
    std::cerr << "Spacing in C2R output does not match!" << std::endl;
    return -1;
    }

  if ( imageAfterInverseFFT->GetDirection() != direction )
    {
    std::cerr << "Direction in C2R output does not match!" << std::endl;
    return -1;
    }

  // The Inverse FFT image iterator is the resultant iterator after we
  // perform the FFT and Inverse FFT on the Original Image. */
  itk::ImageRegionIterator< RealImageType > inverseFFTImageIterator( imageAfterInverseFFT,
                                                                     region );
  counter = 0;
  inverseFFTImageIterator.GoToBegin();

  // Print the Image data obtained by performing the Inverse FFT.
  std::cerr << "---- Inverse FFT image ----" << std::endl;
  while ( !inverseFFTImageIterator.IsAtEnd() )
    {
    TPixel val = inverseFFTImageIterator.Value();
    if ( (counter + 1 ) % SizeOfDimensions[0] == 0 )
      {
      std::cerr << val << std::endl;
      }
    else
      {
      std::cerr << val << " ";
      }
    counter++;
    ++inverseFFTImageIterator;
    }

  std::cerr << std::endl << std::endl;

  // Subtract the Original image Pixel Values from the resultant image
  // values and test whether they are greater than 0.01 for the test
  // to pass.
  originalImageIterator.GoToBegin();
  inverseFFTImageIterator.GoToBegin();
  while ( !originalImageIterator.IsAtEnd() )
    {
    TPixel val = originalImageIterator.Value();
    TPixel val2 = inverseFFTImageIterator.Value();
    TPixel diff = itk::Math::abs( val - val2 );
    if ( itk::Math::NotAlmostEquals(val, 0.0) )
      {
      diff /= itk::Math::abs( val );
      }
    if ( diff > 0.01 )
      {
      std::cerr << "Diff found in test_fft: " << val << " " << val2 << " diff " << diff
                << std::endl;
      return -1;
      }
    ++originalImageIterator;
    ++inverseFFTImageIterator;
    }
  std::cout << std::endl << std::endl;

  return 0;
}


/* test_fft_rtc is the test function to compare two implementations
 * (Direct FFT only).  It is templated over the pixel, Image
 * dimensions and the FFT libraries to be used. */
template< typename TPixel, unsigned int VImageDimensions,
          typename R2CAType, typename R2CBType >
int
test_fft_rtc(unsigned int *SizeOfDimensions)
{
  typedef itk::Image< TPixel, VImageDimensions >                  RealImageType;
  typedef itk::Image< std::complex< TPixel >, VImageDimensions >  ComplexImageType;
  unsigned int counter = 0;
  typename RealImageType::SizeType  imageSize;
  typename RealImageType::IndexType imageIndex;

  // We are testing the FFT for 1D, 2D, and 3D images. An array
  // (SizeOfDimensions) containing the sizes of each dimension is
  // passed as an argument to this function. Based on the template
  // argument VImageDimensions, we create a 1D, 2D, or 3D image by
  // selecting the sizes of image dimensions from this array.
  for (unsigned int i = 0; i < VImageDimensions; i++)
    {
    imageSize.SetElement( i, SizeOfDimensions[i] );
    imageIndex.SetElement( i, 0 );
    }

  typename RealImageType::RegionType region;
  region.SetSize( imageSize );
  region.SetIndex( imageIndex );
  typename RealImageType::Pointer realImage = RealImageType::New();

  // Create the Real Image.
  realImage->SetLargestPossibleRegion( region );
  realImage->SetBufferedRegion( region );
  realImage->SetRequestedRegion( region );
  realImage->Allocate();
  vnl_sample_reseed( static_cast<int >( itksys::SystemTools::GetTime() / 10000.0 ) );

  // We use 2 region iterators for this test the original image
  // iterator and another iterator for the resultant image after
  // performing FFT and IFFT.
  itk::ImageRegionIterator< RealImageType > originalImageIterator( realImage, region );

  // Allocate random pixel values to the image by iterating through it
  // and Print out the image data.
  try
    {
    while( !originalImageIterator.IsAtEnd() )
      {
      TPixel val = vnl_sample_uniform( 0.0, 16384.0 );
      if ( (counter + 1 ) % SizeOfDimensions[0] == 0 )
        {
        std::cout << val << std::endl;
        }
      else
        {
        std::cout << val << " ";
        }
      counter++;
      originalImageIterator.Set( val );
      ++originalImageIterator;
      }
    std::cout << std::endl << std::endl;
    }
  catch( itk::ExceptionObject & ex )
    {
    ex.Print( std::cerr );
    return -1;
    }

  // Real to complex pointers. This computes the forward FFT.
  typename R2CAType::Pointer R2Ca = R2CAType::New();

  // Real to complex pointers. This computes the forward FFT.
  typename R2CBType::Pointer R2Cb = R2CBType::New();

  // Set the real image created as the input to the forward FFT
  // filter.
  R2Ca->SetInput( realImage );
  R2Ca->Update();

  R2Cb->SetInput( realImage );
  R2Cb->Update();

  // Get the size and the pointer to the complex image.
  typename ComplexImageType::Pointer complexImageA = R2Ca->GetOutput();
  std::complex< TPixel > *fftbufA = complexImageA->GetBufferPointer();
  const typename ComplexImageType::SizeType & complexImageSizeA =
    complexImageA->GetLargestPossibleRegion().GetSize();

  typename ComplexImageType::Pointer complexImageB = R2Cb->GetOutput();
  std::complex< TPixel > *fftbufB = complexImageB->GetBufferPointer();
  const typename ComplexImageType::SizeType & complexImageSizeB =
    complexImageB->GetLargestPossibleRegion().GetSize();


  unsigned int sizesA[4] = { 1,1,1,1 };
  unsigned int sizesB[4] = { 1,1,1,1 };
  for(unsigned int i = 0; i < VImageDimensions; i++)
    {
    // The size may be different if one implementation returns a
    // full matrix but not the other.
    sizesA[i] = complexImageSizeA[i];
    sizesB[i] = complexImageSizeB[i];
    }

  // Print out the the frequency domain data obtained after performing
  // the forward transform.
  std::cout << "Frequency domain data after forward transform:" << std::endl;
  for (unsigned int i = 0; i < sizesA[2]; i++)
    {
    unsigned int zStride = i * sizesA[1] * sizesA[0];
    for (unsigned int j = 0; j < sizesA[1]; j++)
      {
      unsigned int yStride = j * sizesA[0];
      for (unsigned int k = 0; k < sizesA[0]; k++)
        {
        std::cout << fftbufA[zStride+yStride+k] << " ";
        }
      std::cout << std::endl;
      }
    }
  std::cout << std::endl << std::endl;

  for (unsigned int i = 0; i < sizesB[2]; i++)
    {
    unsigned int zStride = i * sizesB[1] * sizesB[0];
    for (unsigned int j = 0; j < sizesB[1]; j++)
      {
      unsigned int yStride = j * sizesB[0];
      for (unsigned int k = 0; k < sizesB[0]; k++)
        {
        std::cout << fftbufB[zStride+yStride+k] << " ";
        }
      std::cout << std::endl;
      }
    }
  std::cout << std::endl << std::endl;


  // Subtract the pixel values from the two images. If one pixel
  // difference is greater than 0.01, the test is considered to have
  // failed.
  for (unsigned int i = 0; i < std::min( sizesA[2], sizesB[2] ); i++)
    {
    unsigned int zStrideA = i * sizesA[1] * sizesA[0];
    unsigned int zStrideB = i * sizesB[1] * sizesB[0];
    for (unsigned int j = 0; j < std::min( sizesA[1], sizesB[1] ); j++)
      {
      unsigned int yStrideA = j * sizesA[0];
      unsigned int yStrideB = j * sizesB[0];
      for (unsigned int k = 0; k < std::min( sizesA[0], sizesB[0] ); k++)
        {
        double val = std::abs(fftbufA[zStrideA+yStrideA+k]);
        double diff = std::abs(fftbufA[zStrideA+yStrideA+k] - fftbufB[zStrideB+yStrideB+k]);
        if ( itk::Math::NotAlmostEquals(val, 0.0) )
          {
          diff /= itk::Math::abs( val );
          }
        if ( diff > 0.01 )
          {
          std::cerr << "Diff found in test_fft_r2c: " << fftbufA[zStrideA+yStrideA+k]
                    << " " << fftbufB[zStrideB+yStrideB+k] << " diff " << diff << std::endl;
          return -1;
          }
        }
      }
    }

  std::cout << std::endl << std::endl;

  return 0;
}
#endif