ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMaskedFFTNormalizedCorrelationImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 {
138 
139 template <typename TInputImage, typename TOutputImage, typename TMaskImage=TInputImage >
140 class ITK_TEMPLATE_EXPORT MaskedFFTNormalizedCorrelationImageFilter :
141  public ImageToImageFilter< TInputImage, TOutputImage >
142 {
143 public:
144  ITK_DISALLOW_COPY_AND_ASSIGN(MaskedFFTNormalizedCorrelationImageFilter);
145 
151 
153  itkNewMacro(Self);
154 
157 
160  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
161 
163  using InputImageType = TInputImage;
165  using InputImagePointer = typename InputImageType::Pointer;
166  using InputImageConstPointer = typename InputImageType::ConstPointer;
169 
170  using OutputImageType = TOutputImage;
171  using OutputImagePointer = typename OutputImageType::Pointer;
172  using OutputPixelType = typename OutputImageType::PixelType;
173 
181 
182  using MaskImageType = TMaskImage;
183  using MaskImagePointer = typename MaskImageType::Pointer;
184 
187 
189  itkSetInputMacro(FixedImage, InputImageType);
190  itkGetInputMacro(FixedImage, InputImageType);
192 
194  itkSetInputMacro(MovingImage, InputImageType);
195  itkGetInputMacro(MovingImage, InputImageType);
197 
199  itkSetInputMacro(FixedImageMask, MaskImageType);
200  itkGetInputMacro(FixedImageMask, MaskImageType);
202 
204  itkSetInputMacro(MovingImageMask, MaskImageType);
205  itkGetInputMacro(MovingImageMask, MaskImageType);
207 
209  itkSetMacro(RequiredNumberOfOverlappingPixels,SizeValueType);
210  itkGetMacro(RequiredNumberOfOverlappingPixels,SizeValueType);
212 
214  itkGetMacro(RequiredFractionOfOverlappingPixels,RealPixelType);
215  itkSetClampMacro(RequiredFractionOfOverlappingPixels, RealPixelType, 0.0f, 1.0f);
217 
219  itkGetMacro(MaximumNumberOfOverlappingPixels,SizeValueType);
220 
221 #ifdef ITK_USE_CONCEPT_CHECKING
222  // Begin concept checking
223  itkConceptMacro( OutputPixelTypeIsFloatingPointCheck,
225  // End concept checking
226 #endif
227 
228 protected:
230  {
231  // #0 "FixedImage" required
232  Self::SetPrimaryInputName("FixedImage");
233 
234  // #1 "MaskImage" required
235  Self::AddRequiredInputName("MovingImage", 1);
236 
237  // #2 "FixedImageMask" optional
238  Self::AddOptionalInputName("FixedImageMask", 2);
239 
240  // #3 "MaskImageMask" optional
241  Self::AddOptionalInputName("MovingImageMask", 3);
242 
243  m_RequiredNumberOfOverlappingPixels = 0;
244  m_RequiredFractionOfOverlappingPixels = 0;
245  m_MaximumNumberOfOverlappingPixels = 0;
246  m_AccumulatedProgress = 0.0;
247  }
248  ~MaskedFFTNormalizedCorrelationImageFilter() override = default;
249  void PrintSelf(std::ostream& os, Indent indent) const override;
250 
252  void VerifyInputInformation() ITKv5_CONST override;
253 
255  void GenerateData() override;
256 
262  void GenerateInputRequestedRegion() override;
263 
268  void GenerateOutputInformation() override;
269 
270  void EnlargeOutputRequestedRegion( DataObject *output ) override;
271 
272  typename TMaskImage::Pointer PreProcessMask( const InputImageType * inputImage, const MaskImageType * inputMask );
273 
274  typename TInputImage::Pointer PreProcessImage( const InputImageType * inputImage, const MaskImageType * inputMask );
275 
276  template< typename LocalInputImageType >
277  typename LocalInputImageType::Pointer RotateImage( LocalInputImageType * inputImage );
278 
279  template< typename LocalInputImageType, typename LocalOutputImageType >
280  typename LocalOutputImageType::Pointer CalculateForwardFFT( LocalInputImageType * inputImage, InputSizeType & FFTImageSize );
281 
282  template< typename LocalInputImageType, typename LocalOutputImageType >
283  typename LocalOutputImageType::Pointer CalculateInverseFFT( LocalInputImageType * inputImage, RealSizeType & combinedImageSize );
284 
285  // Helper math methods.
286  template< typename LocalInputImageType, typename LocalOutputImageType >
287  typename LocalOutputImageType::Pointer ElementProduct( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
288 
289  template< typename LocalInputImageType >
290  typename LocalInputImageType::Pointer ElementQuotient( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
291 
292  template< typename LocalInputImageType >
293  typename LocalInputImageType::Pointer ElementSubtraction( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
294 
295  template< typename LocalInputImageType >
296  typename LocalInputImageType::Pointer ElementPositive( LocalInputImageType * inputImage );
297 
298  template< typename LocalInputImageType, typename LocalOutputImageType >
299  typename LocalOutputImageType::Pointer ElementRound( LocalInputImageType * inputImage );
300 
301  // This function factorizes the image size uses factors of 2, 3, and
302  // 5. After this factorization, if there are any remaining values,
303  // the function returns this value.
304  int FactorizeNumber( int n );
305 
306  // Find the closest valid dimension above the desired dimension. This
307  // will be a combination of 2s, 3s, and 5s.
308  int FindClosestValidDimension( int n );
309 
310  template< typename LocalInputImageType >
311  double CalculatePrecisionTolerance( LocalInputImageType * inputImage );
312 
313 private:
317  SizeValueType m_RequiredNumberOfOverlappingPixels;
318 
322  RealPixelType m_RequiredFractionOfOverlappingPixels;
323 
325  SizeValueType m_MaximumNumberOfOverlappingPixels;
326 
328  const unsigned int m_TotalForwardAndInverseFFTs{12};
329 
332 };
333 } // end namespace itk
334 
335 #ifndef ITK_MANUAL_INSTANTIATION
336 #include "itkMaskedFFTNormalizedCorrelationImageFilter.hxx"
337 #endif
338 
339 #endif
typename OutputImageType::Pointer OutputImagePointer
unsigned long SizeValueType
Definition: itkIntTypes.h:83
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename InputImageType::Pointer InputImagePointer
TOutputImage OutputImageType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Calculate masked normalized cross correlation using FFTs.
#define itkConceptMacro(name, concept)
typename InputImageType::ConstPointer InputImageConstPointer
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75