SimpleITK
|
Calculate normalized cross correlation using FFTs. More...
#include <sitkFFTNormalizedCorrelationImageFilter.h>
Public Types | |
using | PixelIDTypeList = BasicPixelIDTypeList |
using | Self = FFTNormalizedCorrelationImageFilter |
![]() | |
using | Self = ImageFilter |
![]() | |
using | Self = ProcessObject |
Public Member Functions | |
Image | Execute (const Image &fixedImage, const Image &movingImage) |
FFTNormalizedCorrelationImageFilter () | |
std::string | GetName () const |
double | GetRequiredFractionOfOverlappingPixels () const |
uint64_t | GetRequiredNumberOfOverlappingPixels () const |
Self & | SetRequiredFractionOfOverlappingPixels (double RequiredFractionOfOverlappingPixels) |
Self & | SetRequiredNumberOfOverlappingPixels (uint64_t RequiredNumberOfOverlappingPixels) |
std::string | ToString () const |
virtual | ~FFTNormalizedCorrelationImageFilter () |
![]() | |
ImageFilter () | |
virtual | ~ImageFilter ()=0 |
![]() | |
virtual void | Abort () |
virtual int | AddCommand (itk::simple::EventEnum event, const std::function< void()> &func) |
Directly add a callback to observe an event. More... | |
virtual int | AddCommand (itk::simple::EventEnum event, itk::simple::Command &cmd) |
Add a Command Object to observer the event. More... | |
virtual float | GetProgress () const |
An Active Measurement of the progress of execution. More... | |
virtual bool | HasCommand (itk::simple::EventEnum event) const |
Query of this object has any registered commands for event. More... | |
ProcessObject () | |
virtual void | RemoveAllCommands () |
Remove all registered commands. More... | |
virtual | ~ProcessObject () |
virtual void | DebugOn () |
virtual void | DebugOff () |
virtual bool | GetDebug () const |
virtual void | SetDebug (bool debugFlag) |
virtual void | SetNumberOfThreads (unsigned int n) |
virtual unsigned int | GetNumberOfThreads () const |
virtual void | SetNumberOfWorkUnits (unsigned int n) |
virtual unsigned int | GetNumberOfWorkUnits () const |
Private Types | |
using | MemberFunctionType = Image(Self::*)(const Image *fixedImage, const Image *movingImage) |
Private Member Functions | |
template<class TImageType > | |
Image | ExecuteInternal (const Image *fixedImage, const Image *movingImage) |
Private Attributes | |
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > | m_MemberFactory |
double | m_RequiredFractionOfOverlappingPixels {0.0} |
uint64_t | m_RequiredNumberOfOverlappingPixels {0u} |
Friends | |
struct | detail::MemberFunctionAddressor< MemberFunctionType > |
Additional Inherited Members | |
![]() | |
static bool | GetGlobalDefaultDebug () |
static void | GlobalDefaultDebugOff () |
static void | GlobalDefaultDebugOn () |
static void | SetGlobalDefaultDebug (bool debugFlag) |
static void | GlobalWarningDisplayOn () |
static void | GlobalWarningDisplayOff () |
static void | SetGlobalWarningDisplay (bool flag) |
static bool | GetGlobalWarningDisplay () |
static double | GetGlobalDefaultCoordinateTolerance () |
Access the global tolerance to determine congruent spaces. More... | |
static void | SetGlobalDefaultCoordinateTolerance (double) |
Access the global tolerance to determine congruent spaces. More... | |
static double | GetGlobalDefaultDirectionTolerance () |
Access the global tolerance to determine congruent spaces. More... | |
static void | SetGlobalDefaultDirectionTolerance (double) |
Access the global tolerance to determine congruent spaces. More... | |
static bool | SetGlobalDefaultThreader (const std::string &threader) |
Set/Get the default threader used for process objects. More... | |
static std::string | GetGlobalDefaultThreader () |
Set/Get the default threader used for process objects. More... | |
static void | SetGlobalDefaultNumberOfThreads (unsigned int n) |
static unsigned int | GetGlobalDefaultNumberOfThreads () |
Set/Get the default threader used for process objects. More... | |
![]() | |
void | CheckImageMatchingDimension (const Image &image1, const Image &image2, const std::string &image2Name) |
void | CheckImageMatchingPixelType (const Image &image1, const Image &image2, const std::string &image2Name) |
void | CheckImageMatchingSize (const Image &image1, const Image &image2, const std::string &image2Name) |
![]() | |
virtual unsigned long | AddITKObserver (const itk::EventObject &, itk::Command *) |
virtual itk::ProcessObject * | GetActiveProcess () |
virtual void | OnActiveProcessDelete () |
virtual void | onCommandDelete (const itk::simple::Command *cmd) noexcept |
virtual void | PreUpdate (itk::ProcessObject *p) |
virtual void | RemoveITKObserver (EventCommand &e) |
![]() | |
NonCopyable ()=default | |
NonCopyable (const NonCopyable &)=delete | |
NonCopyable & | operator= (const NonCopyable &)=delete |
![]() | |
template<class TImageType > | |
static void | FixNonZeroIndex (TImageType *img) |
![]() | |
template<class TImageType > | |
static TImageType::ConstPointer | CastImageToITK (const Image &img) |
template<class TPixelType , unsigned int VImageDimension, unsigned int VLength, template< typename, unsigned int > class TVector> | |
static Image | CastITKToImage (itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img) |
template<unsigned int VImageDimension, unsigned int VLength, template< unsigned int > class TVector> | |
static Image | CastITKToImage (itk::Image< TVector< VLength >, VImageDimension > *img) |
template<class TImageType > | |
static Image | CastITKToImage (TImageType *img) |
static const itk::EventObject & | GetITKEventObject (EventEnum e) |
template<typename T > | |
static std::ostream & | ToStringHelper (std::ostream &os, const T &v) |
static std::ostream & | ToStringHelper (std::ostream &os, const char &v) |
static std::ostream & | ToStringHelper (std::ostream &os, const signed char &v) |
static std::ostream & | ToStringHelper (std::ostream &os, const unsigned char &v) |
Calculate normalized cross correlation using FFTs.
This filter calculates the normalized cross correlation (NCC) of two images using FFTs instead of spatial correlation. It is much faster than spatial correlation for reasonably large structuring elements. This filter is a subclass of the more general MaskedFFTNormalizedCorrelationImageFilter and operates by essentially setting the masks in that algorithm to images of ones. As described in detail in the references below, there is no computational overhead to utilizing the more general masked algorithm because the FFTs of the images of ones are still necessary for the computations.
Inputs: Two images are required as inputs, fixedImage and movingImage. In the context of correlation, inputs are often defined as: "image" and "template". In this filter, the fixedImage plays the role of the image, and the movingImage plays the role of the template. However, this filter is capable of correlating any two images and is not restricted to small movingImages (templates).
Optional parameters: The RequiredNumberOfOverlappingPixels enables the user to specify how many voxels of the two images must overlap; any location in the correlation map that results from fewer than this number of voxels will be set to zero. Larger values zero-out pixels on a larger border around the correlation image. Thus, larger values remove less stable computations but also limit the capture range. If RequiredNumberOfOverlappingPixels is set to 0, the default, no zeroing will take place.
Image size: fixedImage and movingImage need not be the same size. Furthermore, whereas some algorithms require that the "template" be smaller than the "image" because of errors in the regions where the two are not fully overlapping, this filter has no such restriction.
Image spacing: Since the computations are done in the pixel domain, all input images must have the same spacing.
Outputs; The output is an image of RealPixelType that is the NCC of the two images and its values range from -1.0 to 1.0. The size of this NCC image is, by definition, size(fixedImage) + size(movingImage) - 1.
Example filter usage:
References: 1) D. Padfield. "Masked object registration in the Fourier domain." Transactions on Image Processing. 2) D. Padfield. "Masked FFT registration". In Proc. Computer Vision and Pattern Recognition, 2010.
Definition at line 74 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
private |
Setup for member function dispatching
Definition at line 124 of file sitkFFTNormalizedCorrelationImageFilter.h.
Define the pixels types supported by this filter
Definition at line 86 of file sitkFFTNormalizedCorrelationImageFilter.h.
Definition at line 76 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
virtual |
Destructor
itk::simple::FFTNormalizedCorrelationImageFilter::FFTNormalizedCorrelationImageFilter | ( | ) |
Default Constructor that takes no arguments and initializes default parameters
Image itk::simple::FFTNormalizedCorrelationImageFilter::Execute | ( | const Image & | fixedImage, |
const Image & | movingImage | ||
) |
Execute the filter on the input image
|
private |
|
inlinevirtual |
Name of this class
Implements itk::simple::ProcessObject.
Definition at line 110 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
inline |
Set and get the required fraction of overlapping pixels
Definition at line 107 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
inline |
Set and get the required number of overlapping pixels
Definition at line 97 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
inline |
Set and get the required fraction of overlapping pixels
Definition at line 102 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
inline |
Set and get the required number of overlapping pixels
Definition at line 92 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
virtual |
Print ourselves out
Reimplemented from itk::simple::ProcessObject.
|
friend |
Definition at line 128 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
private |
Definition at line 130 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
private |
Definition at line 135 of file sitkFFTNormalizedCorrelationImageFilter.h.
|
private |
Definition at line 133 of file sitkFFTNormalizedCorrelationImageFilter.h.