ITK  6.0.0
Insight Toolkit
Examples/RegistrationITKv4/DisplacementFieldInitialization.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.
*
*=========================================================================*/
// Software Guide : BeginLatex
//
// In order to initialize deformable registration algorithm it is often
// convenient to generate a displacement field from a set of feature
// correspondences provided by the user. The following example illustrates
// how to use the \doxygen{itkLandmarkDisplacementFieldSource} class in order
// to generate a displacement field from the specification of two sets of
// landmarks. Landmarks from one set are associated one-to-one to the
// landmarks in the other set. Each landmark pair defines one displacement
// vector. This class interpolates the values of all other displacement
// vectors using \doxygen{KernelBasedTransform}
//
//
// \index{Registration!Finite Element-Based}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "itkImage.h"
#include "itkVector.h"
#include <fstream>
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Missing Parameters " << std::endl;
std::cerr << "Usage: " << argv[0];
std::cerr << " landmarksFile fixedImage outputDisplacementField"
<< std::endl;
return EXIT_FAILURE;
}
constexpr unsigned int Dimension = 2;
using VectorComponentType = float;
using DisplacementFieldType = itk::Image<VectorType, Dimension>;
using PixelType = unsigned char;
using FixedImageType = itk::Image<PixelType, Dimension>;
using FixedReaderType = itk::ImageFileReader<FixedImageType>;
auto fixedReader = FixedReaderType::New();
fixedReader->SetFileName(argv[2]);
try
{
fixedReader->Update();
}
catch (const itk::ExceptionObject & excp)
{
std::cerr << "Exception thrown " << std::endl;
std::cerr << excp << std::endl;
return EXIT_FAILURE;
}
FixedImageType::ConstPointer fixedImage = fixedReader->GetOutput();
using FilterType =
auto filter = FilterType::New();
filter->SetOutputSpacing(fixedImage->GetSpacing());
filter->SetOutputOrigin(fixedImage->GetOrigin());
filter->SetOutputRegion(fixedImage->GetLargestPossibleRegion());
filter->SetOutputDirection(fixedImage->GetDirection());
// Create source and target landmarks.
//
using LandmarkContainerType = FilterType::LandmarkContainer;
using LandmarkPointType = FilterType::LandmarkPointType;
auto sourceLandmarks = LandmarkContainerType::New();
auto targetLandmarks = LandmarkContainerType::New();
std::ifstream pointsFile;
pointsFile.open(argv[1]);
LandmarkPointType sourcePoint;
pointsFile >> sourcePoint;
LandmarkPointType targetPoint;
pointsFile >> targetPoint;
unsigned int pointId = 0;
while (!pointsFile.fail())
{
sourceLandmarks->InsertElement(pointId, sourcePoint);
targetLandmarks->InsertElement(pointId, targetPoint);
pointId++;
pointsFile >> sourcePoint;
pointsFile >> targetPoint;
}
pointsFile.close();
filter->SetSourceLandmarks(sourceLandmarks);
filter->SetTargetLandmarks(targetLandmarks);
try
{
filter->UpdateLargestPossibleRegion();
}
catch (const itk::ExceptionObject & excp)
{
std::cerr << "Exception thrown " << std::endl;
std::cerr << excp << std::endl;
return EXIT_FAILURE;
}
// Write an image for regression testing
auto writer = WriterType::New();
writer->SetInput(filter->GetOutput());
writer->SetFileName(argv[3]);
filter->Print(std::cout);
try
{
writer->Update();
}
catch (const itk::ExceptionObject & excp)
{
std::cerr << "Exception thrown by writer" << std::endl;
std::cerr << excp << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
// Software Guide : EndLatex
}
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::GTest::TypedefsAndConstructors::Dimension2::VectorType
ImageBaseType::SpacingType VectorType
Definition: itkGTestTypedefsAndConstructors.h:53
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itkImageFileReader.h
itkImage.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itkImageFileWriter.h
itkVector.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itkLandmarkDisplacementFieldSource.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::LandmarkDisplacementFieldSource
Computes a displacement field from two sets of landmarks.
Definition: itkLandmarkDisplacementFieldSource.h:49