ITK  5.2.0
Insight Toolkit
Public Types | Public Member Functions | List of all members
itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage > Class Template Reference

#include <itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h>

+ Inheritance diagram for itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >:
+ Collaboration diagram for itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >:

Public Types

using AccessorType = typename Superclass::AccessorType
 
using FrequencyType = typename ImageType::SpacingType
 
using FrequencyValueType = typename ImageType::SpacingValueType
 
using ImageType = typename Superclass::ImageType
 
using IndexType = typename Superclass::IndexType
 
using InternalPixelType = typename Superclass::InternalPixelType
 
using OffsetType = typename Superclass::OffsetType
 
using PixelContainer = typename Superclass::PixelContainer
 
using PixelContainerPointer = typename Superclass::PixelContainerPointer
 
using PixelType = typename Superclass::PixelType
 
using RegionType = typename Superclass::RegionType
 
using Self = FrequencyFFTLayoutImageRegionConstIteratorWithIndex
 
using SizeType = typename Superclass::SizeType
 
using Superclass = ImageRegionConstIteratorWithIndex< TImage >
 
- Public Types inherited from itk::ImageRegionConstIteratorWithIndex< TImage >
using AccessorType = typename Superclass::AccessorType
 
using ImageType = typename Superclass::ImageType
 
using IndexType = typename Superclass::IndexType
 
using InternalPixelType = typename Superclass::InternalPixelType
 
using OffsetType = typename Superclass::OffsetType
 
using PixelContainer = typename Superclass::PixelContainer
 
using PixelContainerPointer = typename Superclass::PixelContainerPointer
 
using PixelType = typename Superclass::PixelType
 
using RegionType = typename Superclass::RegionType
 
using Self = ImageRegionConstIteratorWithIndex
 
using SizeType = typename Superclass::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
 

Public Member Functions

 FrequencyFFTLayoutImageRegionConstIteratorWithIndex ()
 
 FrequencyFFTLayoutImageRegionConstIteratorWithIndex (const Superclass &it)
 
 FrequencyFFTLayoutImageRegionConstIteratorWithIndex (const TImage *ptr, const RegionType &region)
 
IndexType GetFrequencyBin () const
 
- Public Member Functions inherited from itk::ImageRegionConstIteratorWithIndex< TImage >
 ImageRegionConstIteratorWithIndex ()
 
 ImageRegionConstIteratorWithIndex (const ImageConstIteratorWithIndex< TImage > &it)
 
 ImageRegionConstIteratorWithIndex (const TImage *ptr, const RegionType &region)
 
Selfoperator++ ()
 
Selfoperator-- ()
 
- Public Member Functions inherited from itk::ImageConstIteratorWithIndex< TImage >
const IndexTypeGetIndex () const
 
const RegionTypeGetRegion () const
 
 ImageConstIteratorWithIndex ()
 
 ImageConstIteratorWithIndex (const Self &it)
 
 ImageConstIteratorWithIndex (const TImage *ptr, const RegionType &region)
 
bool operator!= (const Self &it) const
 
bool operator< (const Self &it) const
 
bool operator<= (const Self &it) const
 
Selfoperator= (const Self &it)
 
bool operator== (const Self &it) const
 
bool operator> (const Self &it) const
 
bool operator>= (const Self &it) const
 
virtual ~ImageConstIteratorWithIndex ()=default
 
void SetIndex (const IndexType &ind)
 
PixelType Get () const
 
const PixelTypeValue () const
 
void GoToBegin ()
 
void GoToReverseBegin ()
 
bool IsAtReverseEnd () const
 
bool IsAtEnd () const
 
bool Remaining ()
 
IndexType m_LargestPositiveFrequencyIndex
 
IndexType m_MinIndex
 
IndexType m_MaxIndex
 
FrequencyType m_FrequencyOrigin
 
FrequencyType m_FrequencySpacing
 
bool m_ActualXDimensionIsOdd
 
FrequencyType GetFrequency () const
 
FrequencyValueType GetFrequencyModuloSquare () const
 
virtual const IndexTypeGetLargestPositiveFrequencyIndex () const
 
virtual const IndexTypeGetMinIndex () const
 
virtual const IndexTypeGetMaxIndex () const
 
virtual const FrequencyTypeGetFrequencyOrigin () const
 
virtual const FrequencyTypeGetFrequencySpacing () const
 
void SetActualXDimensionIsOdd (bool value)
 
virtual bool GetActualXDimensionIsOdd ()
 
virtual void ActualXDimensionIsOddOn ()
 
virtual void ActualXDimensionIsOddOff ()
 
void Init ()
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::ImageConstIteratorWithIndex< TImage >
static 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 >
TImage::ConstWeakPointer m_Image
 
IndexType m_PositionIndex
 
IndexType m_BeginIndex
 
IndexType m_EndIndex
 
RegionType m_Region
 
OffsetValueType m_OffsetTable [ImageDimension+1]
 
const InternalPixelTypem_Position
 
const InternalPixelTypem_Begin
 
const InternalPixelTypem_End
 
bool m_Remaining
 
AccessorType m_PixelAccessor
 
AccessorFunctorType m_PixelAccessorFunctor
 

Detailed Description

template<typename TImage>
class itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >

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.

See also
ForwardFFTImageFilter
MORE INFORMATION
For a complete description of the ITK Image Iterators and their API, please see the Iterators chapter in the ITK Software Guide. The ITK Software Guide is available in print and as a free .pdf download from https://www.itk.org.
See also
ImageRegionIteratorWithIndex
ImageConstIterator
ConditionalConstIterator
ConstNeighborhoodIterator
ConstShapedNeighborhoodIterator
ConstSliceIterator
CorrespondenceDataStructureIterator
FloodFilledFunctionConditionalConstIterator
FloodFilledImageFunctionConditionalConstIterator
FloodFilledImageFunctionConditionalIterator
FloodFilledSpatialFunctionConditionalConstIterator
FloodFilledSpatialFunctionConditionalIterator
ImageConstIterator
ImageConstIteratorWithIndex
ImageIterator
ImageIteratorWithIndex
ImageLinearConstIteratorWithIndex
ImageLinearIteratorWithIndex
ImageRandomConstIteratorWithIndex
ImageRandomIteratorWithIndex
ImageRegionConstIterator
ImageRegionConstIteratorWithIndex
ImageRegionExclusionConstIteratorWithIndex
ImageRegionExclusionIteratorWithIndex
ImageRegionIterator
ImageRegionIteratorWithIndex
ImageRegionReverseConstIterator
ImageRegionReverseIterator
ImageReverseConstIterator
ImageReverseIterator
ImageSliceConstIteratorWithIndex
ImageSliceIteratorWithIndex
NeighborhoodIterator
PathConstIterator
PathIterator
ShapedNeighborhoodIterator
SliceIterator
ImageConstIteratorWithIndex

Definition at line 115 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Member Typedef Documentation

◆ AccessorType

◆ FrequencyType

template<typename TImage>
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::FrequencyType = typename ImageType::SpacingType

◆ FrequencyValueType

template<typename TImage>
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::FrequencyValueType = typename ImageType::SpacingValueType

◆ ImageType

template<typename TImage>
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::ImageType = typename Superclass::ImageType

◆ IndexType

template<typename TImage>
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::IndexType = typename Superclass::IndexType

Types inherited from the Superclass

Definition at line 124 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

◆ InternalPixelType

◆ OffsetType

◆ PixelContainer

◆ PixelContainerPointer

◆ PixelType

template<typename TImage>
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::PixelType = typename Superclass::PixelType

◆ RegionType

◆ Self

Standard class type alias.

Definition at line 120 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

◆ SizeType

template<typename TImage>
using itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::SizeType = typename Superclass::SizeType

◆ Superclass

Constructor & Destructor Documentation

◆ FrequencyFFTLayoutImageRegionConstIteratorWithIndex() [1/3]

Default constructor. Needed since we provide a cast constructor.

Definition at line 138 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

References itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

◆ FrequencyFFTLayoutImageRegionConstIteratorWithIndex() [2/3]

template<typename TImage>
itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::FrequencyFFTLayoutImageRegionConstIteratorWithIndex ( const TImage *  ptr,
const RegionType region 
)
inline

Constructor establishes an iterator to walk a particular image and a particular region of that image.

Definition at line 146 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

References itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

◆ FrequencyFFTLayoutImageRegionConstIteratorWithIndex() [3/3]

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 158 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

References itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

Member Function Documentation

◆ ActualXDimensionIsOddOff()

template<typename TImage>
virtual void itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::ActualXDimensionIsOddOff ( )
virtual

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;

◆ ActualXDimensionIsOddOn()

template<typename TImage>
virtual void itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::ActualXDimensionIsOddOn ( )
virtual

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;

◆ GetActualXDimensionIsOdd()

template<typename TImage>
virtual bool itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetActualXDimensionIsOdd ( )
virtual

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;

◆ GetFrequency()

template<typename TImage>
FrequencyType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequency ( ) const
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 209 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

References itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyBin(), itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_FrequencyOrigin, and itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_FrequencySpacing.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyModuloSquare().

◆ GetFrequencyBin()

template<typename TImage>
IndexType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyBin ( ) const
inline

◆ GetFrequencyModuloSquare()

template<typename TImage>
FrequencyValueType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyModuloSquare ( ) const
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 223 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

References itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequency().

◆ GetFrequencyOrigin()

template<typename TImage>
virtual const FrequencyType& itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyOrigin ( ) const
virtual

Origin of frequencies is zero for FFT output.

◆ GetFrequencySpacing()

template<typename TImage>
virtual const FrequencyType& itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencySpacing ( ) const
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.

◆ GetLargestPositiveFrequencyIndex()

template<typename TImage>
virtual const IndexType& itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetLargestPositiveFrequencyIndex ( ) const
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

◆ GetMaxIndex()

template<typename TImage>
virtual const IndexType& itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetMaxIndex ( ) const
virtual

Default to UpperIndex of the largest possible Region.

◆ GetMinIndex()

template<typename TImage>
virtual const IndexType& itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetMinIndex ( ) const
virtual

Default to first index of the largest possible Region.

◆ Init()

template<typename TImage>
void itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init ( )
inlineprivate

◆ SetActualXDimensionIsOdd()

template<typename TImage>
void itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::SetActualXDimensionIsOdd ( bool  value)
inline

Does nothing. This member only affects HalfHermitianFrequencyIterator. Provided for homogeneous interface between iterators.

Definition at line 266 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

References itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_ActualXDimensionIsOdd.

Member Data Documentation

◆ m_ActualXDimensionIsOdd

template<typename TImage>
bool itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_ActualXDimensionIsOdd
private

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 303 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::SetActualXDimensionIsOdd().

◆ m_FrequencyOrigin

template<typename TImage>
FrequencyType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_FrequencyOrigin
private

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 301 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequency(), and itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

◆ m_FrequencySpacing

template<typename TImage>
FrequencyType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_FrequencySpacing
private

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 302 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequency(), and itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

◆ m_LargestPositiveFrequencyIndex

template<typename TImage>
IndexType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_LargestPositiveFrequencyIndex
private

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 298 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyBin(), and itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

◆ m_MaxIndex

template<typename TImage>
IndexType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_MaxIndex
private

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 300 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyBin(), and itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().

◆ m_MinIndex

template<typename TImage>
IndexType itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::m_MinIndex
private

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 299 of file itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h.

Referenced by itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::GetFrequencyBin(), and itk::FrequencyFFTLayoutImageRegionConstIteratorWithIndex< TImage >::Init().


The documentation for this class was generated from the following file: