Registration and Transformation
Overview
This example demonstrates image registration and transformation using SimpleElastix. SimpleElastix provides access to the powerful elastix registration toolbox through SimpleITK’s ElastixImageFilter and TransformixImageFilter classes.
This example includes two programs:
elx - Performs image registration using ElastixImageFilter
tfx - Applies transformations using TransformixImageFilter
The registration workflow typically consists of:
Reading fixed and moving images
Loading or creating a parameter map that defines the registration strategy
Executing the registration to obtain a registered image and transform parameters
Optionally applying the transform to other images using Transformix
Parameter maps can be loaded from text files, created using default presets, or customized programmatically. The elastix parameter files use a simple text format that allows fine-grained control over the registration process.
Code
Registration (elx)
using System;
using itk.simple;
namespace itk.simple.examples
{
class Elastix
{
static void Main(string[] args)
{
if (args.Length < 5)
{
Console.WriteLine("Usage: elx <fixedImage> <movingImage> <parameterFile> <outputImage> <outputParameterFile>");
return;
}
// Instantiate SimpleElastix
ElastixImageFilter elastixImageFilter = new ElastixImageFilter();
// Read input
elastixImageFilter.SetFixedImage(SimpleITK.ReadImage(args[0]));
elastixImageFilter.SetMovingImage(SimpleITK.ReadImage(args[1]));
elastixImageFilter.SetParameterMap(SimpleITK.ReadParameterFile(args[2]));
elastixImageFilter.LogToConsoleOn();
// Run registration
elastixImageFilter.Execute();
// Write result image
SimpleITK.WriteImage(elastixImageFilter.GetResultImage(), args[3]);
// Write parameter file. This example only supports one parameter map and one transform parameter map.
SimpleITK.WriteParameterFile(elastixImageFilter.GetTransformParameterMaps()[0], args[4]);
}
}
}
// This one header will include all SimpleITK filters and external
// objects.
#include <SimpleITK.h>
#include <iostream>
#include <stdlib.h>
namespace sitk = itk::simple;
int
main(int argc, char * argv[])
{
if (argc < 5)
{
std::cerr << "Usage: " << argv[0]
<< " <fixedImage> <movingImage> <parameterFile> <outputImage> <outputParameterFile> \n";
return 1;
}
// Instantiate SimpleElastix
sitk::ElastixImageFilter elastixImageFilter;
// Read input
sitk::ImageFileReader reader;
reader.SetFileName(std::string(argv[1]));
elastixImageFilter.SetFixedImage(reader.Execute());
reader.SetFileName(std::string(argv[2]));
elastixImageFilter.SetMovingImage(reader.Execute());
elastixImageFilter.SetParameterMap(sitk::ReadParameterFile(std::string(argv[3])));
elastixImageFilter.LogToConsoleOn();
// Run registration
elastixImageFilter.Execute();
// Write result image
sitk::ImageFileWriter writer;
writer.SetFileName(std::string(argv[4]));
writer.Execute(elastixImageFilter.GetResultImage());
// Write parameter file. This test executable only supports one parameter map and one transform parameter map.
sitk::WriteParameterFile(elastixImageFilter.GetTransformParameterMaps()[0], std::string(argv[5]));
return 0;
}
import org.itk.simple.*;
public class elx {
public static void main(String[] args) {
if (args.length < 5) {
System.out.println("Usage: elx <fixedImage> <movingImage> <parameterFile> <outputImage> <outputParameterFile>");
System.exit(1);
}
// Instantiate SimpleElastix
ElastixImageFilter elastixImageFilter = new ElastixImageFilter();
// Read input
elastixImageFilter.setFixedImage(SimpleITK.readImage(args[0]));
elastixImageFilter.setMovingImage(SimpleITK.readImage(args[1]));
elastixImageFilter.setParameterMap(SimpleITK.readParameterFile(args[2]));
elastixImageFilter.logToConsoleOn();
// Run registration
elastixImageFilter.execute();
// Write result image
SimpleITK.writeImage(elastixImageFilter.getResultImage(), args[3]);
// Write parameter file. This example only supports one parameter map and one transform parameter map.
SimpleITK.writeParameterFile(elastixImageFilter.getTransformParameterMaps().get(0), args[4]);
}
}
require "SimpleITK"
if #arg < 5 then
print("Usage: elx <fixedImage> <movingImage> <parameterFile> <outputImage> <outputParameterFile>")
os.exit(1)
end
-- Instantiate SimpleElastix
elastixImageFilter = SimpleITK.ElastixImageFilter()
-- Read input
elastixImageFilter:SetFixedImage(SimpleITK.ReadImage(arg[1]))
elastixImageFilter:SetMovingImage(SimpleITK.ReadImage(arg[2]))
elastixImageFilter:SetParameterMap(SimpleITK.ReadParameterFile(arg[3]))
elastixImageFilter:LogToConsoleOn()
-- Run registration
elastixImageFilter:Execute()
-- Write result image
SimpleITK.WriteImage(elastixImageFilter:GetResultImage(), arg[4])
-- Write parameter file. This example only supports one parameter map and one transform parameter map.
parameterMaps = elastixImageFilter:GetTransformParameterMaps()
SimpleITK.WriteParameterFile(parameterMaps[0], arg[5])
#!/usr/bin/env python
# This example demonstrates the usage of the SimpleElastix framework
# for image registration using the ElastixImageFilter.
import SimpleITK as sitk
import sys
if len(sys.argv) < 6:
print(
f"Usage: {sys.argv[0]} <fixedImage> <movingImage> <parameterFile> <outputImage> <outputParameterFile>"
)
sys.exit(1)
# Instantiate SimpleElastix
elastix_image_filter = sitk.ElastixImageFilter()
# Read input images
fixed_image = sitk.ReadImage(sys.argv[1])
moving_image = sitk.ReadImage(sys.argv[2])
# Set images
elastix_image_filter.SetFixedImage(fixed_image)
elastix_image_filter.SetMovingImage(moving_image)
# Read and set parameter map
parameter_map = sitk.ReadParameterFile(sys.argv[3])
elastix_image_filter.SetParameterMap(parameter_map)
# Enable logging to console
elastix_image_filter.LogToConsoleOn()
# Run registration
elastix_image_filter.Execute()
# Write result image
sitk.WriteImage(elastix_image_filter.GetResultImage(), sys.argv[4])
# Write parameter file
# This example only supports one parameter map and one transform parameter map
sitk.WriteParameterFile(elastix_image_filter.GetTransformParameterMaps()[0], sys.argv[5])
library(SimpleITK)
args <- commandArgs(TRUE)
if (length(args) < 5) {
stop("Usage: elx <fixedImage> <movingImage> <parameterFile> <outputImage> <outputParameterFile>")
}
# Instantiate SimpleElastix
elastixImageFilter <- ElastixImageFilter()
# Read input
elastixImageFilter$SetFixedImage(ReadImage(args[1]))
elastixImageFilter$SetMovingImage(ReadImage(args[2]))
elastixImageFilter$SetParameterMap(elastixImageFilter$ReadParameterFile(args[3]))
elastixImageFilter$LogToConsoleOn()
# Run registration
output_image <- elastixImageFilter$Execute()
# Write result image
WriteImage(output_image, args[4])
# Write parameter file. This example only supports one parameter map and one transform parameter map.
elastixImageFilter$WriteParameterFile(elastixImageFilter$GetTransformParameterMap(0L), args[5])
Transformation (tfx)
using System;
using itk.simple;
namespace itk.simple.examples
{
class Transformix
{
static void Main(string[] args)
{
if (args.Length < 3)
{
Console.WriteLine("Usage: tfx <inputImage> <parameterFile> <outputImage>");
return;
}
// Instantiate transformix
TransformixImageFilter transformixImageFilter = new TransformixImageFilter();
transformixImageFilter.LogToConsoleOn();
// Read input
transformixImageFilter.SetMovingImage(SimpleITK.ReadImage(args[0]));
transformixImageFilter.SetTransformParameterMap(SimpleITK.ReadParameterFile(args[1]));
// Run warp
transformixImageFilter.Execute();
// Write result image
SimpleITK.WriteImage(transformixImageFilter.GetResultImage(), args[2]);
}
}
}
// This one header will include all SimpleITK filters and external
// objects.
#include <SimpleITK.h>
#include <iostream>
#include <stdlib.h>
namespace sitk = itk::simple;
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Usage: " << argv[0] << " <inputImage> <parameterFile> <outputImage>\n";
return 1;
}
// Instantiate transformix
sitk::TransformixImageFilter transformixImageFilter;
transformixImageFilter.LogToConsoleOn();
// Read input
sitk::ImageFileReader reader;
reader.SetFileName(std::string(argv[1]));
transformixImageFilter.SetMovingImage(reader.Execute());
transformixImageFilter.SetTransformParameterMap(sitk::ReadParameterFile(argv[2]));
// Run warp
transformixImageFilter.Execute();
// Write result image
sitk::ImageFileWriter writer;
writer.SetFileName(std::string(argv[3]));
writer.Execute(transformixImageFilter.GetResultImage());
return 0;
}
import org.itk.simple.*;
public class tfx {
public static void main(String[] args) {
if (args.length < 3) {
System.out.println("Usage: tfx <inputImage> <parameterFile> <outputImage>");
System.exit(1);
}
// Instantiate transformix
TransformixImageFilter transformixImageFilter = new TransformixImageFilter();
transformixImageFilter.logToConsoleOn();
// Read input
transformixImageFilter.setMovingImage(SimpleITK.readImage(args[0]));
transformixImageFilter.setTransformParameterMap(SimpleITK.readParameterFile(args[1]));
// Run warp
transformixImageFilter.execute();
// Write result image
SimpleITK.writeImage(transformixImageFilter.getResultImage(), args[2]);
}
}
require "SimpleITK"
if #arg < 3 then
print("Usage: tfx <inputImage> <parameterFile> <outputImage>")
os.exit(1)
end
-- Instantiate transformix
transformixImageFilter = SimpleITK.TransformixImageFilter()
transformixImageFilter:LogToConsoleOn()
-- Read input
transformixImageFilter:SetMovingImage(SimpleITK.ReadImage(arg[1]))
transformixImageFilter:SetTransformParameterMap(SimpleITK.ReadParameterFile(arg[2]))
-- Run warp
transformixImageFilter:Execute()
-- Write result image
SimpleITK.WriteImage(transformixImageFilter:GetResultImage(), arg[3])
#!/usr/bin/env python
import SimpleITK as sitk
import sys
if len(sys.argv) < 4:
print(f"Usage: {sys.argv[0]} <inputImage> <parameterFile> <outputImage>")
sys.exit(1)
# Instantiate transformix
transformix_image_filter = sitk.TransformixImageFilter()
# Enable logging to console
transformix_image_filter.LogToConsoleOn()
# Read input image
moving_image = sitk.ReadImage(sys.argv[1])
transformix_image_filter.SetMovingImage(moving_image)
# Read and set transform parameter map
transform_parameter_map = sitk.ReadParameterFile(sys.argv[2])
transformix_image_filter.SetTransformParameterMap(transform_parameter_map)
# Run warp
transformix_image_filter.Execute()
# Write result image
sitk.WriteImage(transformix_image_filter.GetResultImage(), sys.argv[3])
library(SimpleITK)
args <- commandArgs(TRUE)
if (length(args) < 3) {
stop("Usage: tfx <inputImage> <parameterFile> <outputImage>")
}
# Instantiate transformix
transformixImageFilter <- TransformixImageFilter()
transformixImageFilter$LogToConsoleOn()
# Read input
transformixImageFilter$SetMovingImage(ReadImage(args[1]))
transformixImageFilter$SetTransformParameterMap(transformixImageFilter$ReadParameterFile(args[2]))
# Run warp
resampled_image <- transformixImageFilter$Execute()
# Write result image
WriteImage(resampled_image, args[3])