SimpleITK  
Elastix/ParameterMaps/ParameterMaps.py
1#!/usr/bin/env python
2
3import sys
4import SimpleITK as sitk
5
6if len(sys.argv) < 3:
7 print("Usage: ParameterMaps <fixedImage> <movingImage> [outputImage]")
8 sys.exit(1)
9
10# Read input images
11fixed_image = sitk.ReadImage(sys.argv[1])
12moving_image = sitk.ReadImage(sys.argv[2])
13
14print("SimpleITK Elastix Parameter Maps Example")
15print("=" * 60)
16
17#START_DEFAULT_PARAMETER_MAP
18# Get a default parameter map for a common registration type
19parameter_map = sitk.GetDefaultParameterMap("translation")
20
21# Print parameter map to see its contents
22sitk.PrintParameterMap(parameter_map)
23#END_DEFAULT_PARAMETER_MAP
24
25#START_MODIFY_PARAMETER_MAP
26# Modify parameter map entries (values are lists of strings)
27parameter_map["Transform"] = ["AffineTransform"]
28parameter_map["MaximumNumberOfIterations"] = ["512"]
29parameter_map["NumberOfSpatialSamples"] = ["8192"]
30#END_MODIFY_PARAMETER_MAP
31
32#START_MULTI_STAGE_REGISTRATION
33# Multi-stage registration: translation -> affine -> bspline
34elastix_image_filter = sitk.ElastixImageFilter()
35elastix_image_filter.SetFixedImage(fixed_image)
36elastix_image_filter.SetMovingImage(moving_image)
37
38elastix_image_filter.SetParameterMap(sitk.GetDefaultParameterMap("translation"))
39elastix_image_filter.AddParameterMap(sitk.GetDefaultParameterMap("affine"))
40
41bspline_map = sitk.GetDefaultParameterMap("bspline")
42bspline_map["FinalGridSpacingInPhysicalUnits"] = ["8.0"]
43elastix_image_filter.AddParameterMap(bspline_map)
44
45elastix_image_filter.LogToConsoleOff()
46result_image = elastix_image_filter.Execute()
47
48# Retrieve transform parameter maps from the result
49transform_parameter_maps = elastix_image_filter.GetTransformParameterMaps()
50#END_MULTI_STAGE_REGISTRATION
51
52#START_READ_WRITE_PARAMETER_MAP
53# Write a parameter map to a file
54sitk.WriteParameterFile(parameter_map, "parameters.txt")
55
56# Read a parameter map from a file
57parameter_map = sitk.ReadParameterFile("parameters.txt")
58#END_READ_WRITE_PARAMETER_MAP
59
60# Write results if output path provided
61if len(sys.argv) > 3:
62 output_dir = sys.argv[3]
63 for i, tpm in enumerate(transform_parameter_maps):
64 filename = f"{output_dir}_TransformParameters_{i}.txt"
65 sitk.WriteParameterFile(tpm, filename)
66 sitk.WriteImage(result_image, f"{output_dir}.nii")
67print("Example completed successfully!")
68print("=" * 60)
The class that wraps the elastix registration library.
SITKIO_EXPORT Image ReadImage(const PathType &filename, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageFileReader class which is convenient for most image r...
SITKElastix_EXPORT std::map< std::string, std::vector< std::string > > GetDefaultParameterMap(const std::string transform, const unsigned int numberOfResolutions=4, const double finalGridSpacingInPhysicalUnits=8.0)
SITKIO_EXPORT void WriteImage(const Image &image, const PathType &fileName, bool useCompression=false, int compressionLevel=-1)
WriteImage is a procedural interface to the ImageFileWriter. class which is convenient for many image...
SITKElastix_EXPORT void PrintParameterMap(const std::map< std::string, std::vector< std::string > > parameterMap)
SITKElastix_EXPORT std::map< std::string, std::vector< std::string > > ReadParameterFile(const std::string filename)
SITKElastix_EXPORT void WriteParameterFile(const std::map< std::string, std::vector< std::string > > parameterMap, const std::string filename)