ITK
6.0.0
Insight Toolkit
|
#include <itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h>
A multi-dimensional iterator templated over image type that walks pixels within a region and is specialized to keep track of its image index location.
This class is a specialization of ImageRegionConstIteratorWithIndex, adding method GetFrequencyBins to give the frequency bins corresponding to image indices, and GetFrequency to get the frequency of the bin. The frequency bins depends on the image size. The default assumes that the image to iterate over is the output of a forward FFT filter, where the first index corresponds to 0 frequency, and Nyquist Frequencies are in the middle, between positive and negative frequencies.
This class can be specialized further to iterate over other frequency layouts, for example shifted images (where 0 frequency is in the middle of the image, and Nyquist are in the border). For different layout, use other frequency iterator.
This iterator is for the frequency layout that results from applying a FFT (vnl or fftw) to an image. The default ImageInformation is: Origin = {{0}}, Spacing = {{1}}. In this case the frequency values will be in the range: [-1/2, 1/2] Hz Or [-pi, pi] rad/s To modify those ranges: a) Avoid modifying the origin. The origin index always corresponds to zero frequency after a FFT. The range should be always centered around zero. b) The spacing control the range of frequencies (always around zero). If the spacing is = {{0.5}} we get a frequency range of [-1/4, 1/4] or [-pi/2, pi/2].
The frequency layout is assumed to be: where fs is frequency sampling, or frequency spacing (1.0 by default). If N is even: Nyquist frequency at index=N/2 is shared between + and - regions. <---------—positive f ------------—><---------—negative f-----------—> 0(DC) fs/N 2*fs/N ... (N/2 -1)*fs/N fs/2 -(N/2-1)*fs/N ... -2*fs/N -fs/N
Example: Size 6:
0 <– DC Component (0 freq) 1 | 5 2 | 4 3 <– Shared between regions, unique Nyquist.
If N is odd: Nyquist frequency is not represented but there are symmetric largest frequencies at index=N/2, N/2 +1 <-------—positive f ------------—><---------—negative f--------------—> 0(DC) fs/N 2*fs/N ...... fs/2*(N-1)/N -fs/2*(N-1)/N ... -2*fs/N -fs/N
Example: Size 5:
0 <– DC Component (0 freq) 1 | 4 2 | 3 <– Absolute Largest Frequency bins (+, -)
Please see ImageRegionConstIteratorWithIndex for more information.
Definition at line 115 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
Public Types | |
using | AccessorType = typename TImage::AccessorType |
using | FrequencyType = typename ImageType::SpacingType |
using | FrequencyValueType = typename ImageType::SpacingValueType |
using | ImageType = TImage |
using | IndexType = typename TImage::IndexType |
using | InternalPixelType = typename TImage::InternalPixelType |
using | OffsetType = typename TImage::OffsetType |
using | PixelContainer = typename TImage::PixelContainer |
using | PixelContainerPointer = typename PixelContainer::Pointer |
using | PixelType = typename TImage::PixelType |
using | RegionType = typename TImage::RegionType |
using | Self = FrequencyFFTLayoutImageRegionConstIteratorWithIndex |
using | SizeType = typename TImage::SizeType |
using | Superclass = ImageRegionConstIteratorWithIndex< TImage > |
Public Types inherited from itk::ImageRegionConstIteratorWithIndex< TImage > | |
using | AccessorType = typename TImage::AccessorType |
using | ImageType = TImage |
using | IndexType = typename TImage::IndexType |
using | InternalPixelType = typename TImage::InternalPixelType |
using | OffsetType = typename TImage::OffsetType |
using | PixelContainer = typename TImage::PixelContainer |
using | PixelContainerPointer = typename PixelContainer::Pointer |
using | PixelType = typename TImage::PixelType |
using | RegionType = typename TImage::RegionType |
using | Self = ImageRegionConstIteratorWithIndex |
using | SizeType = typename TImage::SizeType |
using | Superclass = ImageConstIteratorWithIndex< TImage > |
Public Types inherited from itk::ImageConstIteratorWithIndex< TImage > | |
using | AccessorFunctorType = typename TImage::AccessorFunctorType |
using | AccessorType = typename TImage::AccessorType |
using | ImageType = TImage |
using | IndexType = typename TImage::IndexType |
using | IndexValueType = typename IndexType::IndexValueType |
using | InternalPixelType = typename TImage::InternalPixelType |
using | OffsetType = typename TImage::OffsetType |
using | OffsetValueType = typename OffsetType::OffsetValueType |
using | PixelContainer = typename TImage::PixelContainer |
using | PixelContainerPointer = typename PixelContainer::Pointer |
using | PixelType = typename TImage::PixelType |
using | RegionType = typename TImage::RegionType |
using | Self = ImageConstIteratorWithIndex |
using | SizeType = typename TImage::SizeType |
using | SizeValueType = typename SizeType::SizeValueType |
Private Member Functions | |
void | Init () |
Private Attributes | |
bool | m_ActualXDimensionIsOdd |
FrequencyType | m_FrequencyOrigin |
FrequencyType | m_FrequencySpacing |
IndexType | m_LargestPositiveFrequencyIndex |
IndexType | m_MaxIndex |
IndexType | m_MinIndex |
Additional Inherited Members | |
Static Public Member Functions inherited from itk::ImageConstIteratorWithIndex< TImage > | |
static constexpr unsigned int | GetImageDimension () |
Static Public Attributes inherited from itk::ImageConstIteratorWithIndex< TImage > | |
static constexpr unsigned int | ImageDimension = TImage::ImageDimension |
Protected Attributes inherited from itk::ImageConstIteratorWithIndex< TImage > | |
const InternalPixelType * | m_Begin { nullptr } |
IndexType | m_BeginIndex { { 0 } } |
const InternalPixelType * | m_End { nullptr } |
IndexType | m_EndIndex { { 0 } } |
TImage::ConstWeakPointer | m_Image {} |
OffsetValueType | m_OffsetTable [ImageDimension+1] {} |
AccessorType | m_PixelAccessor {} |
AccessorFunctorType | m_PixelAccessorFunctor {} |
const InternalPixelType * | m_Position { nullptr } |
IndexType | m_PositionIndex { { 0 } } |
RegionType | m_Region {} |
bool | m_Remaining { false } |
using itk::ImageConstIteratorWithIndex< TImage >::AccessorType = typename TImage::AccessorType |
Accessor type that converts data between internal and external representations.
Definition at line 132 of file itkImageConstIteratorWithIndex.h.
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::FrequencyType = typename ImageType::SpacingType |
Definition at line 136 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::FrequencyValueType = typename ImageType::SpacingValueType |
Definition at line 137 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::ImageType = TImage |
Image type alias support
Definition at line 116 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::IndexType = typename TImage::IndexType |
Index type alias support
Definition at line 105 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::InternalPixelType = typename TImage::InternalPixelType |
Internal Pixel Type
Definition at line 125 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::OffsetType = typename TImage::OffsetType |
Type of the Offset taken from the image
Definition at line 136 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::PixelContainer = typename TImage::PixelContainer |
PixelContainer type alias support Used to refer to the container for the pixel data. While this was already typedef'ed in the superclass, it needs to be redone here for this subclass to compile properly with gcc.
Definition at line 121 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::PixelContainerPointer = typename PixelContainer::Pointer |
Definition at line 122 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::PixelType = typename TImage::PixelType |
External Pixel Type
Definition at line 128 of file itkImageConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::RegionType = typename TImage::RegionType |
Region type alias support
Definition at line 113 of file itkImageConstIteratorWithIndex.h.
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Self = FrequencyFFTLayoutImageRegionConstIteratorWithIndex |
Standard class type alias.
Definition at line 121 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
using itk::ImageConstIteratorWithIndex< TImage >::SizeType = typename TImage::SizeType |
Size type alias support
Definition at line 109 of file itkImageConstIteratorWithIndex.h.
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Superclass = ImageRegionConstIteratorWithIndex<TImage> |
Definition at line 122 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
inline |
Default constructor. Needed since we provide a cast constructor.
Definition at line 139 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
inline |
Constructor establishes an iterator to walk a particular image and a particular region of that image. Initializes the iterator at the begin of the region.
Definition at line 147 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
inlineexplicit |
Constructor that can be used to cast from an ImageIterator to an ImageRegionIteratorWithIndex. Many routines return an ImageIterator, but for a particular task, you may want an ImageRegionIteratorWithIndex. Rather than provide overloaded APIs that return different types of Iterators, itk returns ImageIterators and uses constructors to cast from an ImageIterator to a ImageRegionIteratorWithIndex.
Definition at line 159 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
virtual |
Does nothing. This member only affects HalfHermitianFrequencyIterator. Provided for homogeneous interface between iterators.
|
virtual |
Does nothing. This member only affects HalfHermitianFrequencyIterator. Provided for homogeneous interface between iterators.
|
inline |
Note that this method is independent of the region in the constructor. It takes into account the ImageInformation of the Image in the frequency domain. This iterator is for the frequency layout that results from applying a FFT (vnl or fftw) to an image. If your image has a different layout, use other frequency iterator. The default ImageInformation is: Origin = {{0}}, Spacing = {{1}}. In this case the frequency values will be in the range: [-1/2, 1/2] Hz Or [-pi, pi] rad/s To modify those ranges: a) Avoid modifying the origin. The origin index always corresponds to zero frequency after a FFT. The range should be always centered around zero. b) The spacing control the range of frequencies (always around zero). If the spacing is = {{0.5}} we get a frequency range of [-1/4, 1/4] or [-pi/2, pi/2].
f = [0, 1, ..., N/2-1, -N/2, ..., -1] * FrequencySpacing if N is even f = [0, 1, ..., (N-1)/2, -(N-1)/2, ..., -1] * FrequencySpacing if N is odd
Where FrequencySpacing = samplingFrequency / N; and samplingFrequency = 1.0 / inputImageSpatialDomainSpacing;
Definition at line 208 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
inline |
Definition at line 171 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
inline |
Definition at line 222 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
virtual |
Origin of frequencies is zero for FFT output.
|
virtual |
This is the pixel width, or the bin size of the frequency in physical or world coordinates. SamplingFrequency = 1.0 / SpatialImageSpacing FrequencySpacing = SamplingFrequency / ImageSize FrequencySpacing = 1.0 / (SpatialImageSpacing * ImageSize) FrequencySpacing is computed at construction from the spacing of the input image, and cannot be modified.
|
virtual |
Index corresponding to the first highest frequency (Nyquist) after a FFT transform. If the size of the image is even, the Nyquist frequency = fs/2 is unique and shared between positive and negative frequencies. (Even Size) LargestPositiveFrequencyIndex = originIndex + N / 2 If odd, Nyquist frequency is not represented, but there is still a largest frequency at this index = fs/2 * (N-1)/N. (Odd Size) LargestPositiveFrequencyIndex = originIndex + (N + 1) / 2
|
virtual |
Default to UpperIndex of the largest possible Region.
|
virtual |
Default to first index of the largest possible Region.
|
inlineprivate |
Calculate Nyquist frequency index (m_LargestPositiveFrequencyIndex), Min/Max indices from LargestPossibleRegion. Also sets FrequencySpacing and FrequencyOrigin. Called by constructors.
Definition at line 278 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
inline |
Does nothing. This member only affects HalfHermitianFrequencyIterator. Provided for homogeneous interface between iterators.
Definition at line 265 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
private |
Definition at line 302 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
private |
Definition at line 300 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
private |
Definition at line 301 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
private |
Definition at line 297 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
private |
Definition at line 299 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.
|
private |
Definition at line 298 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.