ITK  5.4.0
Insight Toolkit
itkMaskedFFTNormalizedCorrelationImageFilter.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 itkMaskedFFTNormalizedCorrelationImageFilter_h
19 #define itkMaskedFFTNormalizedCorrelationImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkImage.h"
23 
24 namespace itk
25 {
144 template <typename TInputImage, typename TOutputImage, typename TMaskImage = TInputImage>
146  : public ImageToImageFilter<TInputImage, TOutputImage>
147 {
148 public:
149  ITK_DISALLOW_COPY_AND_MOVE(MaskedFFTNormalizedCorrelationImageFilter);
150 
156 
158  itkNewMacro(Self);
159 
161  itkOverrideGetNameOfClassMacro(MaskedFFTNormalizedCorrelationImageFilter);
162 
165  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
166 
168  using InputImageType = TInputImage;
174 
175  using OutputImageType = TOutputImage;
177  using OutputPixelType = typename OutputImageType::PixelType;
178 
186 
187  using MaskImageType = TMaskImage;
189 
192 
194  itkSetInputMacro(FixedImage, InputImageType);
195  itkGetInputMacro(FixedImage, InputImageType);
199  itkSetInputMacro(MovingImage, InputImageType);
200  itkGetInputMacro(MovingImage, InputImageType);
204  itkSetInputMacro(FixedImageMask, MaskImageType);
205  itkGetInputMacro(FixedImageMask, MaskImageType);
209  itkSetInputMacro(MovingImageMask, MaskImageType);
210  itkGetInputMacro(MovingImageMask, MaskImageType);
214  itkSetMacro(RequiredNumberOfOverlappingPixels, SizeValueType);
215  itkGetMacro(RequiredNumberOfOverlappingPixels, SizeValueType);
219  itkGetMacro(RequiredFractionOfOverlappingPixels, RealPixelType);
220  itkSetClampMacro(RequiredFractionOfOverlappingPixels, RealPixelType, 0.0f, 1.0f);
224  itkGetMacro(MaximumNumberOfOverlappingPixels, SizeValueType);
225 
226 #ifdef ITK_USE_CONCEPT_CHECKING
227  // Begin concept checking
228  itkConceptMacro(OutputPixelTypeIsFloatingPointCheck, (Concept::IsFloatingPoint<OutputPixelType>));
229  // End concept checking
230 #endif
231 
232 protected:
234  {
235  // #0 "FixedImage" required
236  Self::SetPrimaryInputName("FixedImage");
237 
238  // #1 "MaskImage" required
239  Self::AddRequiredInputName("MovingImage", 1);
240 
241  // #2 "FixedImageMask" optional
242  Self::AddOptionalInputName("FixedImageMask", 2);
243 
244  // #3 "MaskImageMask" optional
245  Self::AddOptionalInputName("MovingImageMask", 3);
246 
247  m_RequiredNumberOfOverlappingPixels = 0;
248  m_RequiredFractionOfOverlappingPixels = 0;
249  m_MaximumNumberOfOverlappingPixels = 0;
250  m_AccumulatedProgress = 0.0;
251  }
252  ~MaskedFFTNormalizedCorrelationImageFilter() override = default;
253  void
254  PrintSelf(std::ostream & os, Indent indent) const override;
255 
257  void
258  VerifyInputInformation() ITKv5_CONST override;
259 
261  void
262  GenerateData() override;
263 
269  void
270  GenerateInputRequestedRegion() override;
271 
276  void
277  GenerateOutputInformation() override;
278 
279  void
280  EnlargeOutputRequestedRegion(DataObject * output) override;
281 
282  typename TMaskImage::Pointer
283  PreProcessMask(const InputImageType * inputImage, const MaskImageType * inputMask);
284 
285  typename TInputImage::Pointer
286  PreProcessImage(const InputImageType * inputImage, const MaskImageType * inputMask);
287 
288  template <typename LocalInputImageType>
289  typename LocalInputImageType::Pointer
290  RotateImage(LocalInputImageType * inputImage);
291 
292  template <typename LocalInputImageType, typename LocalOutputImageType>
293  typename LocalOutputImageType::Pointer
294  CalculateForwardFFT(LocalInputImageType * inputImage, InputSizeType & FFTImageSize);
295 
296  template <typename LocalInputImageType, typename LocalOutputImageType>
297  typename LocalOutputImageType::Pointer
298  CalculateInverseFFT(LocalInputImageType * inputImage, RealSizeType & combinedImageSize);
299 
300  // Helper math methods.
301  template <typename LocalInputImageType, typename LocalOutputImageType>
302  typename LocalOutputImageType::Pointer
303  ElementProduct(LocalInputImageType * inputImage1, LocalInputImageType * inputImage2);
304 
305  template <typename LocalInputImageType>
306  typename LocalInputImageType::Pointer
307  ElementQuotient(LocalInputImageType * inputImage1, LocalInputImageType * inputImage2);
308 
309  template <typename LocalInputImageType>
310  typename LocalInputImageType::Pointer
311  ElementSubtraction(LocalInputImageType * inputImage1, LocalInputImageType * inputImage2);
312 
313  template <typename LocalInputImageType>
314  typename LocalInputImageType::Pointer
315  ElementPositive(LocalInputImageType * inputImage);
316 
317  template <typename LocalInputImageType, typename LocalOutputImageType>
318  typename LocalOutputImageType::Pointer
319  ElementRound(LocalInputImageType * inputImage);
320 
321  // This function factorizes the image size uses factors of 2, 3, and
322  // 5. After this factorization, if there are any remaining values,
323  // the function returns this value.
324  int
325  FactorizeNumber(int n);
326 
327  // Find the closest valid dimension above the desired dimension. This
328  // will be a combination of 2s, 3s, and 5s.
329  int
330  FindClosestValidDimension(int n);
331 
332  template <typename LocalInputImageType>
333  double
334  CalculatePrecisionTolerance(LocalInputImageType * inputImage);
335 
336 private:
340  SizeValueType m_RequiredNumberOfOverlappingPixels{};
341 
344  RealPixelType m_RequiredFractionOfOverlappingPixels{};
345 
347  SizeValueType m_MaximumNumberOfOverlappingPixels{};
348 
350  const unsigned int m_TotalForwardAndInverseFFTs{ 12 };
351 
353  float m_AccumulatedProgress{};
354 };
355 } // end namespace itk
356 
357 #ifndef ITK_MANUAL_INSTANTIATION
358 # include "itkMaskedFFTNormalizedCorrelationImageFilter.hxx"
359 #endif
360 
361 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::MaskImageType
TInputImage MaskImageType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:187
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::RealSizeType
typename RealImageType::SizeType RealSizeType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:183
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::SizeValueType
itk::SizeValueType SizeValueType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:173
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::MaskImagePointer
typename MaskImageType::Pointer MaskImagePointer
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:188
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:181
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::RealIndexType
typename RealImageType::IndexType RealIndexType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:182
itk::MaskedFFTNormalizedCorrelationImageFilter::MaskedFFTNormalizedCorrelationImageFilter
MaskedFFTNormalizedCorrelationImageFilter()
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:233
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::Concept::IsFloatingPoint
Definition: itkConceptChecking.h:946
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::MaskedFFTNormalizedCorrelationImageFilter
Calculate masked normalized cross correlation using FFTs.
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:145
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::RealPixelType
OutputPixelType RealPixelType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:179
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::OutputPixelType
typename OutputImageType::PixelType OutputPixelType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:177
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::RealRegionType
typename RealImageType::RegionType RealRegionType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:184
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::InputRegionType
typename InputImageType::RegionType InputRegionType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:169
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::InputSizeType
typename InputImageType::SizeType InputSizeType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:172
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::FFTImagePointer
typename FFTImageType::Pointer FFTImagePointer
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:191
itk::MaskedFFTNormalizedCorrelationImageFilter< TInputImage, TOutputImage >::RealPointType
typename RealImageType::PointType RealPointType
Definition: itkMaskedFFTNormalizedCorrelationImageFilter.h:185
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293