ITK  6.0.0
Insight Toolkit
itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * https://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18#ifndef itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex_h
19#define itkFrequencyFFTLayoutImageRegionConstIteratorWithIndex_h
20
22
23namespace itk
24{
114template <typename TImage>
116 : public ImageRegionConstIteratorWithIndex<TImage>
117{
118public:
119
123
125 using typename Superclass::IndexType;
126 using typename Superclass::SizeType;
127 using typename Superclass::OffsetType;
128 using typename Superclass::RegionType;
129 using typename Superclass::ImageType;
130 using typename Superclass::PixelContainer;
132 using typename Superclass::InternalPixelType;
133 using typename Superclass::PixelType;
134 using typename Superclass::AccessorType;
135
136 using FrequencyType = typename ImageType::SpacingType;
137 using FrequencyValueType = typename ImageType::SpacingValueType;
141 {
142 this->Init();
143 }
144
148 : ImageRegionConstIteratorWithIndex<TImage>(ptr, region)
149 {
150 this->Init();
151 }
152
161 {
162 this->Init();
163 }
164
165 /*
166 * Image Index [0, N - 1] returns [0 to N/2] (positive) union [-N/2 + 1, -1] (negative). So index N/2 + 1 returns the
167 * bin -N/2 + 1. If first index of the image is not zero, it stills returns values in the same range. f = [0, 1, ...,
168 * N/2-1, -N/2, ..., -1] if N is even f = [0, 1, ..., (N-1)/2, -(N-1)/2, ..., -1] if N is odd
169 */
172 {
173 IndexType freqInd;
174
175 freqInd.Fill(0);
176 for (unsigned int dim = 0; dim < TImage::ImageDimension; ++dim)
177 {
178 if (this->m_PositionIndex[dim] <= m_LargestPositiveFrequencyIndex[dim])
179 {
180 freqInd[dim] = this->m_PositionIndex[dim] - this->m_MinIndex[dim];
181 }
182 else // -. From -N/2 + 1 (Nyquist if even) to -1 (-df in frequency)
183 {
184 freqInd[dim] = this->m_PositionIndex[dim] - (this->m_MaxIndex[dim] + 1);
185 }
186 }
187 return freqInd;
188 }
189
209 FrequencyType
211 {
212 FrequencyType freq;
213 IndexType freqInd = this->GetFrequencyBin();
216 for (unsigned int dim = 0; dim < TImage::ImageDimension; ++dim)
217 {
218 freq[dim] = this->m_FrequencyOrigin[dim] + this->m_FrequencySpacing[dim] * freqInd[dim];
219 }
220 return freq;
221 }
222
223 FrequencyValueType
225 {
226 FrequencyValueType w2(0);
227 FrequencyType w(this->GetFrequency());
228
229 for (unsigned int dim = 0; dim < TImage::ImageDimension; ++dim)
230 {
231 w2 += w[dim] * w[dim];
232 }
233 return w2;
234 }
235
245 itkGetConstReferenceMacro(LargestPositiveFrequencyIndex, IndexType);
246
248 itkGetConstReferenceMacro(MinIndex, IndexType);
249
251 itkGetConstReferenceMacro(MaxIndex, IndexType);
252
254 itkGetConstReferenceMacro(FrequencyOrigin, FrequencyType);
255
262 itkGetConstReferenceMacro(FrequencySpacing, FrequencyType);
263
266 void
268 {
269 this->m_ActualXDimensionIsOdd = value;
270 }
271 itkGetMacro(ActualXDimensionIsOdd, bool);
272 itkBooleanMacro(ActualXDimensionIsOdd);
275private:
279 void
281 {
282 SizeType sizeImage = this->m_Image->GetLargestPossibleRegion().GetSize();
283 this->m_MinIndex = this->m_Image->GetLargestPossibleRegion().GetIndex();
284 this->m_MaxIndex = this->m_Image->GetLargestPossibleRegion().GetUpperIndex();
285 for (unsigned int dim = 0; dim < ImageType::ImageDimension; ++dim)
286 {
287 this->m_LargestPositiveFrequencyIndex[dim] =
288 static_cast<FrequencyValueType>(this->m_MinIndex[dim] + std::floor(sizeImage[dim] / 2.0));
289 // Set frequency metadata.
290 // Origin of frequencies is zero after a FFT.
291 this->m_FrequencyOrigin[dim] = 0.0;
292 // SamplingFrequency = 1.0 / SpatialImageSpacing
293 // Freq_BinSize = SamplingFrequency / Size
294 this->m_FrequencySpacing[dim] = 1.0 / (this->m_Image->GetSpacing()[dim] * sizeImage[dim]);
295 }
296 }
305};
306} // end namespace itk
307#endif
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
A base class for multi-dimensional iterators templated over image type that are designed to efficient...
typename PixelContainer::Pointer PixelContainerPointer
typename TImage::InternalPixelType InternalPixelType
typename TImage::PixelContainer PixelContainer
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....