ITK  4.2.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 {
115 
116 template <class TInputImage, class TOutputImage >
118  public ImageToImageFilter< TInputImage, TOutputImage >
119 {
120 public:
126 
128  itkNewMacro(Self);
129 
131  itkTypeMacro(MaskedFFTNormalizedCorrelationImageFilter, MaskedFFTNormalizedCrossCorrelationImageFilter);
132 
135  itkStaticConstMacro(ImageDimension, unsigned int,
136  TOutputImage::ImageDimension);
137 
139  typedef TInputImage InputImageType;
140  typedef TOutputImage OutputImageType;
141  typedef typename InputImageType::RegionType InputRegionType;
142  typedef typename InputImageType::Pointer InputImagePointer;
143  typedef typename InputImageType::ConstPointer InputImageConstPointer;
144  typedef typename InputImageType::SizeType InputSizeType;
145  typedef typename OutputImageType::Pointer OutputImagePointer;
146  typedef typename OutputImageType::PixelType OutputPixelType;
147 
154 
157 
159  void SetFixedImage(const InputImageType *input)
160  {
161  this->SetNthInput(0, const_cast<InputImageType *>(input) );
162  }
163  InputImageType * GetFixedImage()
164  {
165  return static_cast<InputImageType*>(const_cast<DataObject *>(this->ProcessObject::GetInput(0)));
166  }
168 
170  void SetMovingImage(const InputImageType *input)
171  {
172  this->SetNthInput(1, const_cast<InputImageType *>(input) );
173  }
174  InputImageType * GetMovingImage()
175  {
176  return static_cast<InputImageType*>(const_cast<DataObject *>(this->ProcessObject::GetInput(1)));
177  }
179 
181  void SetFixedImageMask(const InputImageType *input)
182  {
183  this->SetNthInput(2, const_cast<InputImageType *>(input) );
184  }
185  InputImageType * GetFixedImageMask()
186  {
187  return static_cast<InputImageType*>(const_cast<DataObject *>(this->ProcessObject::GetInput(2)));
188  }
190 
192  void SetMovingImageMask(const InputImageType *input)
193  {
194  this->SetNthInput(3, const_cast<InputImageType *>(input) );
195  }
196  InputImageType * GetMovingImageMask()
197  {
198  return static_cast<InputImageType*>(const_cast<DataObject *>(this->ProcessObject::GetInput(3)));
199  }
201 
203  itkSetMacro(RequiredNumberOfOverlappingPixels,SizeValueType);
204  itkGetMacro(RequiredNumberOfOverlappingPixels,SizeValueType);
206 
207 #ifdef ITK_USE_CONCEPT_CHECKING
208 
209  itkConceptMacro( OutputPixelTypeIsFloatingPointCheck,
211 
213 #endif
214 
215 protected:
217  {
218  this->SetNumberOfRequiredInputs(2);
219  m_RequiredNumberOfOverlappingPixels = 0;
220  }
222  void PrintSelf(std::ostream& os, Indent indent) const;
223 
225  void VerifyInputInformation();
226 
228  void GenerateData();
229 
235  virtual void GenerateInputRequestedRegion();
236 
241  void GenerateOutputInformation();
242 
243  template< class LocalInputImageType >
244  typename LocalInputImageType::Pointer PreProcessMask( const LocalInputImageType * inputImage, const LocalInputImageType * inputMask );
245 
246  template< class LocalInputImageType >
247  typename LocalInputImageType::Pointer PreProcessImage( const LocalInputImageType * inputImage, LocalInputImageType * inputMask );
248 
249  typename TInputImage::Pointer RotateImage( InputImageType * inputImage );
250 
251  template< class LocalInputImageType, class LocalOutputImageType >
252  typename LocalOutputImageType::Pointer CalculateForwardFFT( LocalInputImageType * inputImage, InputSizeType & FFTImageSize );
253 
254  template< class LocalInputImageType, class LocalOutputImageType >
255  typename LocalOutputImageType::Pointer CalculateInverseFFT( LocalInputImageType * inputImage, RealSizeType & combinedImageSize );
256 
257  // Helper math methods.
258  template< class LocalInputImageType, class LocalOutputImageType >
259  typename LocalOutputImageType::Pointer ElementProduct( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
260 
261  template< class LocalInputImageType >
262  typename LocalInputImageType::Pointer ElementQuotient( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
263 
264  template< class LocalInputImageType >
265  typename LocalInputImageType::Pointer ElementSubtraction( LocalInputImageType * inputImage1, LocalInputImageType * inputImage2 );
266 
267  template< class LocalInputImageType >
268  typename LocalInputImageType::Pointer ElementPositive( LocalInputImageType * inputImage );
269 
270  template< class LocalInputImageType, class LocalOutputImageType >
271  typename LocalOutputImageType::Pointer ElementRound( LocalInputImageType * inputImage );
272 
273  // This function factorizes the image size uses factors of 2, 3, and
274  // 5. After this factorization, if there are any remaining values,
275  // the function returns this value.
276  int FactorizeNumber( int n );
277 
278  // Find the closest valid dimension above the desired dimension. This
279  // will be a combination of 2s, 3s, and 5s.
280  int FindClosestValidDimension( int n );
281 
282  template< class LocalInputImageType >
283  double CalculatePrecisionTolerance( LocalInputImageType * inputImage );
284 
285 private:
286  MaskedFFTNormalizedCorrelationImageFilter(const Self&); //purposely not implemented
287  void operator=(const Self&); //purposely not implemented
288 
293 };
294 } // end namespace itk
295 
296 #ifndef ITK_MANUAL_INSTANTIATION
297 #include "itkMaskedFFTNormalizedCorrelationImageFilter.hxx"
298 #endif
299 
300 #endif
301