ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Filtering/Convolution/NormalizedCorrelationUsingFFTWithMaskImages/Code.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#include "itkImage.h"
#include <iostream>
#include <string>
namespace
{
using ImageType = itk::Image<unsigned char, 2>;
using FloatImageType = itk::Image<float, 2>;
} // namespace
static void
CreateMask(MaskType * const mask);
static void
CreateImage(ImageType::Pointer image, const itk::Index<2> & cornerOfSquare);
int
main()
{
itk::Index<2> offset;
offset[0] = 5;
offset[1] = 6;
// Setup mask
auto mask = MaskType::New();
CreateMask(mask);
auto fixedImage = ImageType::New();
itk::Index<2> cornerOfFixedSquare;
cornerOfFixedSquare[0] = 3;
cornerOfFixedSquare[1] = 8;
CreateImage(fixedImage, cornerOfFixedSquare);
itk::WriteImage(fixedImage.GetPointer(), "fixedImage.png");
auto movingImage = ImageType::New();
itk::Index<2> cornerOfMovingSquare;
cornerOfMovingSquare[0] = cornerOfFixedSquare[0] + offset[0];
cornerOfMovingSquare[1] = cornerOfFixedSquare[1] + offset[1];
CreateImage(movingImage, cornerOfMovingSquare);
itk::WriteImage(movingImage.GetPointer(), "movingImage.png");
// Perform normalized correlation
auto correlationFilter = CorrelationFilterType::New();
correlationFilter->SetFixedImage(fixedImage);
correlationFilter->SetMovingImage(movingImage);
correlationFilter->SetMovingImageMask(mask);
// correlationFilter->SetFixedImageMask(mask);
correlationFilter->Update();
itk::WriteImage(correlationFilter->GetOutput(), "correlation.mha");
auto rescaleFilter = RescaleFilterType::New();
rescaleFilter->SetInput(correlationFilter->GetOutput());
rescaleFilter->SetOutputMinimum(0);
rescaleFilter->SetOutputMaximum(255);
rescaleFilter->Update();
itk::WriteImage(rescaleFilter->GetOutput(), "correlation.png");
using MinimumMaximumImageCalculatorType = itk::MinimumMaximumImageCalculator<FloatImageType>;
MinimumMaximumImageCalculatorType::Pointer minimumMaximumImageCalculatorFilter =
minimumMaximumImageCalculatorFilter->SetImage(correlationFilter->GetOutput());
minimumMaximumImageCalculatorFilter->Compute();
itk::Index<2> maximumCorrelationPatchCenter = minimumMaximumImageCalculatorFilter->GetIndexOfMaximum();
itk::Size<2> outputSize = correlationFilter->GetOutput()->GetLargestPossibleRegion().GetSize();
itk::Index<2> maximumCorrelationPatchCenterFixed;
maximumCorrelationPatchCenterFixed[0] = outputSize[0] / 2 - maximumCorrelationPatchCenter[0];
maximumCorrelationPatchCenterFixed[1] = outputSize[1] / 2 - maximumCorrelationPatchCenter[1];
std::cout << "Maximum location: " << maximumCorrelationPatchCenter << std::endl;
std::cout << "Maximum location fixed: " << maximumCorrelationPatchCenterFixed
<< std::endl; // This is the value we expect!
std::cout << "Maximum value: " << minimumMaximumImageCalculatorFilter->GetMaximum()
<< std::endl; // If the images can be perfectly aligned, the value is 1
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image, const itk::Index<2> & cornerOfSquare)
{
start.Fill(0);
size.Fill(51);
ImageType::RegionType region(start, size);
image->SetRegions(region);
image->Allocate();
image->FillBuffer(0);
itk::ImageRegionIterator<ImageType> imageIterator(image, region);
ImageType::IndexValueType squareSize = 8;
while (!imageIterator.IsAtEnd())
{
if (imageIterator.GetIndex()[0] > cornerOfSquare[0] &&
imageIterator.GetIndex()[0] < cornerOfSquare[0] + squareSize &&
imageIterator.GetIndex()[1] > cornerOfSquare[1] && imageIterator.GetIndex()[1] < cornerOfSquare[1] + squareSize)
{
imageIterator.Set(255);
}
++imageIterator;
}
}
void
CreateMask(MaskType * const mask)
{
start.Fill(0);
size.Fill(51);
ImageType::RegionType region(start, size);
mask->SetRegions(region);
mask->Allocate();
mask->FillBuffer(255); // Make the whole mask "valid"
// unsigned int squareSize = 8;
ImageType::IndexValueType squareSize = 3;
itk::Index<2> cornerOfSquare = { { 3, 8 } };
// Remove pixels from the mask in a small square. The correlationw will not be computed at these pixels.
itk::ImageRegionIterator<MaskType> maskIterator(mask, region);
while (!maskIterator.IsAtEnd())
{
if (maskIterator.GetIndex()[0] > cornerOfSquare[0] && maskIterator.GetIndex()[0] < cornerOfSquare[0] + squareSize &&
maskIterator.GetIndex()[1] > cornerOfSquare[1] && maskIterator.GetIndex()[1] < cornerOfSquare[1] + squareSize)
{
maskIterator.Set(0);
}
++maskIterator;
}
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Index
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:70
itk::Size< 2 >
itk::MinimumMaximumImageCalculator
Computes the minimum and the maximum intensity values of an image.
Definition: itkMinimumMaximumImageCalculator.h:45
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkMinimumMaximumImageCalculator.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:274
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:213
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:90
itk::FFTNormalizedCorrelationImageFilter
Calculate normalized cross correlation using FFTs.
Definition: itkFFTNormalizedCorrelationImageFilter.h:103
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::Index::GetIndex
const IndexValueType * GetIndex() const
Definition: itkIndex.h:232
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itkRegionOfInterestImageFilter.h
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition: itkRescaleIntensityImageFilter.h:133
itkImageKernelOperator.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:254
itkFFTNormalizedCorrelationImageFilter.h
itk::Size::GetSize
const SizeValueType * GetSize() const
Definition: itkSize.h:171