ITK  4.13.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:
149 
151  itkNewMacro(Self);
152 
155 
158  itkStaticConstMacro(ImageDimension, unsigned int,
159  TOutputImage::ImageDimension);
160 
162  typedef TInputImage InputImageType;
163  typedef typename InputImageType::RegionType InputRegionType;
164  typedef typename InputImageType::Pointer InputImagePointer;
165  typedef typename InputImageType::ConstPointer InputImageConstPointer;
168 
169  typedef TOutputImage OutputImageType;
170  typedef typename OutputImageType::Pointer OutputImagePointer;
171  typedef typename OutputImageType::PixelType OutputPixelType;
172 
180 
181  typedef TMaskImage MaskImageType;
182  typedef typename MaskImageType::Pointer MaskImagePointer;
183 
186 
188  itkSetInputMacro(FixedImage, InputImageType);
189  itkGetInputMacro(FixedImage, InputImageType);
191 
193  itkSetInputMacro(MovingImage, InputImageType);
194  itkGetInputMacro(MovingImage, InputImageType);
196 
198  itkSetInputMacro(FixedImageMask, MaskImageType);
199  itkGetInputMacro(FixedImageMask, MaskImageType);
201 
203  itkSetInputMacro(MovingImageMask, MaskImageType);
204  itkGetInputMacro(MovingImageMask, MaskImageType);
206 
208  itkSetMacro(RequiredNumberOfOverlappingPixels,SizeValueType);
209  itkGetMacro(RequiredNumberOfOverlappingPixels,SizeValueType);
211 
213  itkGetMacro(RequiredFractionOfOverlappingPixels,RealPixelType);
214  itkSetClampMacro(RequiredFractionOfOverlappingPixels, RealPixelType, 0.0f, 1.0f);
216 
218  itkGetMacro(MaximumNumberOfOverlappingPixels,SizeValueType);
219 
220 #ifdef ITK_USE_CONCEPT_CHECKING
221  // Begin concept checking
222  itkConceptMacro( OutputPixelTypeIsFloatingPointCheck,
224  // End concept checking
225 #endif
226 
227 protected:
228  MaskedFFTNormalizedCorrelationImageFilter():m_TotalForwardAndInverseFFTs(12)
229  {
230  // #0 "FixedImage" required
231  Self::SetPrimaryInputName("FixedImage");
232 
233  // #1 "MaskImage" required
234  Self::AddRequiredInputName("MovingImage", 1);
235 
236  // #2 "FixedImageMask" optional
237  Self::AddOptionalInputName("FixedImageMask", 2);
238 
239  // #3 "MaskImageMask" optional
240  Self::AddOptionalInputName("MovingImageMask", 3);
241 
242  m_RequiredNumberOfOverlappingPixels = 0;
243  m_RequiredFractionOfOverlappingPixels = 0;
244  m_MaximumNumberOfOverlappingPixels = 0;
245  m_AccumulatedProgress = 0.0;
246  }
248  void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
249 
251  void VerifyInputInformation() ITK_OVERRIDE;
252 
254  void GenerateData() ITK_OVERRIDE;
255 
261  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
262 
267  void GenerateOutputInformation() ITK_OVERRIDE;
268 
269  void EnlargeOutputRequestedRegion( DataObject *output ) ITK_OVERRIDE;
270 
271  typename TMaskImage::Pointer PreProcessMask( const InputImageType * inputImage, const MaskImageType * inputMask );
272 
273  typename TInputImage::Pointer PreProcessImage( const InputImageType * inputImage, const MaskImageType * inputMask );
274 
275  template< typename LocalInputImageType >
276  typename LocalInputImageType::Pointer RotateImage( LocalInputImageType * inputImage );
277 
278  template< typename LocalInputImageType, typename LocalOutputImageType >
279  typename LocalOutputImageType::Pointer CalculateForwardFFT( LocalInputImageType * inputImage, InputSizeType & FFTImageSize );
280 
281  template< typename LocalInputImageType, typename LocalOutputImageType >
282  typename LocalOutputImageType::Pointer CalculateInverseFFT( LocalInputImageType * inputImage, RealSizeType & combinedImageSize );
283 
284  // Helper math methods.
285  template< typename LocalInputImageType, typename LocalOutputImageType >
286  typename LocalOutputImageType::Pointer ElementProduct( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
287 
288  template< typename LocalInputImageType >
289  typename LocalInputImageType::Pointer ElementQuotient( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
290 
291  template< typename LocalInputImageType >
292  typename LocalInputImageType::Pointer ElementSubtraction( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
293 
294  template< typename LocalInputImageType >
295  typename LocalInputImageType::Pointer ElementPositive( LocalInputImageType * inputImage );
296 
297  template< typename LocalInputImageType, typename LocalOutputImageType >
298  typename LocalOutputImageType::Pointer ElementRound( LocalInputImageType * inputImage );
299 
300  // This function factorizes the image size uses factors of 2, 3, and
301  // 5. After this factorization, if there are any remaining values,
302  // the function returns this value.
303  int FactorizeNumber( int n );
304 
305  // Find the closest valid dimension above the desired dimension. This
306  // will be a combination of 2s, 3s, and 5s.
307  int FindClosestValidDimension( int n );
308 
309  template< typename LocalInputImageType >
310  double CalculatePrecisionTolerance( LocalInputImageType * inputImage );
311 
312 private:
313  ITK_DISALLOW_COPY_AND_ASSIGN(MaskedFFTNormalizedCorrelationImageFilter);
314 
318  SizeValueType m_RequiredNumberOfOverlappingPixels;
319 
323  RealPixelType m_RequiredFractionOfOverlappingPixels;
324 
326  SizeValueType m_MaximumNumberOfOverlappingPixels;
327 
329  const unsigned int m_TotalForwardAndInverseFFTs;
330 
332  float m_AccumulatedProgress;
333 };
334 } // end namespace itk
335 
336 #ifndef ITK_MANUAL_INSTANTIATION
337 #include "itkMaskedFFTNormalizedCorrelationImageFilter.hxx"
338 #endif
339 
340 #endif
Superclass::RegionType RegionType
Definition: itkImage.h:137
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Base class for all process objects that output image data.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
Image< std::complex< RealPixelType >, ImageDimension > FFTImageType
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)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75