ITK v6 Migration Guide¶
This guide documents the changes required to migrate a code base which uses ITK v5 to use ITK v6. The migration guide for transition from v4 to v5 can be found here.
Require modern C++ language feature use¶
Many backward compatible/ forward enabling compiler features are now required to be used.
Replace ITKv5_CONST with const
Minimum supported VTK version: 9.1¶
ITKVtkGlue (and any code that links it) now requires VTK 9.1 or
newer. The previous floor of VTK 8.1 has been raised because:
VTK 9.0 was the first stable release of the new module system, which replaces the legacy
vtk_module_config(NS …)macro with imported targets such asVTK::CommonCoreandVTK::IOImage.VTK 9.0 removed the
OpenGLrendering backend; onlyOpenGL2is shipped.VTK 9.1 stabilised the imported-target names and the
vtk_module_autoinithelper that ITKVtkGlue now relies on unconditionally, removing a significant amount of legacy compatibility code fromModules/Bridge/VtkGlue/.
VTK 9.4 or newer is recommended for current security/CVE patches and bug fixes, but not required — every imported-target name and CMake helper that ITKVtkGlue depends on was already stable in 9.1, so ITK 6 is happy to build against any 9.1.x release. Use 9.4+ if you have the choice.
What you need to do¶
If your project links ITKVtkGlue (directly or via find_package(ITK)):
Build / link against VTK 9.1 or newer. Distribution packages that ship a supported VTK as of this release: Ubuntu 22.04 LTS and 24.04 LTS (
libvtk9.1), Fedora 38+ (9.2+), Homebrew (brew install vtktracks current — currently 9.5). Anaconda’sconda-forge::vtkis at 9.4+. 3D Slicer and ParaView build VTK 9.4+ from source as part of their superbuild.Remove any local
vtk_module_config(...)calls in your CMake. The shim that ITKVtkGlue used to inject when consuming legacy VTK has been deleted. Replace such calls with a directtarget_link_libraries(<your-target> PUBLIC VTK::IOImage VTK::…), then callvtk_module_autoinitso VTK’s object factories and I/O readers/writers register at startup. Without this call, factory- dispatched features (e.g.,vtkImageReader2Factory) silently produce empty results:target_link_libraries(<your-target> PUBLIC VTK::IOImage VTK::ImagingSources # …other VTK modules you use ) vtk_module_autoinit( TARGETS <your-target> MODULES VTK::IOImage VTK::ImagingSources # …same list as above )
Remove any
VTK_RENDERING_BACKENDselection logic that fell back toOpenGL. OnlyOpenGL2is supported now. Code such asif(VTK_RENDERING_BACKEND STREQUAL "OpenGL") # ... endif()
can be deleted; the special cases for
OpenGLno longer exist upstream.Update
find_package(VTK …)calls to require 9.1 explicitly:find_package(VTK 9.1 NO_MODULE REQUIRED)
If a project also exports its own
<Project>Config.cmake, mirror the version pin so downstream consumers get a clear failure rather than a confusing missing-target error. Bump to9.4instead of9.1if you want to enforce the recommended-not-required floor in your own build.
Note for downstream projects that pin both ITK and VTK¶
ITK 6 + VTK 9.1+ is the supported combination, with VTK 9.4+
recommended. If you maintain a SuperBuild that pins both
(e.g., 3D Slicer, BRAINSTools, ANTs, ITK Remote Modules), bump the
pinned VTK tag to at least v9.1.0 (or, preferably, the latest
9.4.x) before bumping the ITK pin to 6. Building ITK 6 against VTK
8.x or VTK 9.0 will fail at find_package time with a clear
version-mismatch error.
Distributions still on VTK 8.x or VTK 9.0 should track ITK 5.4.x until their VTK is updated; ITK 5.4 patch releases continue to support the wider VTK 8.1+/9.0+ range.
Remove support for ITKv4 interfaces¶
ITKV4_COMPATIBILITY is no longer a supported option. Previously this option
was off by default. Previously when enabled, the ITK API was modified to
provide support for ITKV4 functions.
Remove support for ITKv5¶
All contents of the Deprecated module were removed. This includes TreeContainer and related classes; atomic primitives, mutexes and related classes which are now part of C++ standard; specialized Vector filters - specialized versions are no longer needed, as regular filters can work with vector images. For details, see ITKv5 Migration Guide.
Prefer standard CXX language features rather than ITK macros¶
Replace ITK aliases (left column) with CXX standard feature (right column)
ITK_FALLTHROUGH [[fallthrough]]
ITK_DELETE_FUNCTION = delete
ITK_CONSTEXPR_FUNC constexpr
ITK_CONSTEXPR_VAR constexpr
ITK_ALIGNAS(X) alignas(X)
ITK_ALIGNOF(X) alignof(X)
ITK_DEPRECATED [[deprecated]]
ITK_DEPRECATED_MSG(MSG) [[deprecated(MSG)]]
ITK_CONSTEXPR constexpr
ITK_DELETED_FUNCTION = delete
ITK_EXTERN_TEMPLATE extern
ITK_FINAL final
ITK_NOEXCEPT noexcept
ITK_NOEXCEPT_EXPR(X) noexcept(X)
ITK_NULLPTR nullptr
ITK_OVERRIDE override
ITK_STATIC_ASSERT(X) static_assert(X, #X)
ITK_STATIC_ASSERT_MSG(X, MSG) static_assert(X, MSG)
ITK_THREAD_LOCAL thread_local
Removed ITKv5 migration/maintenance scripts¶
The following scripts used for migrating to ITKv5 were removed from the ITKv6.
CheckForOutdatedDefines.sh
EnumPrintFunction.py
Move_DISALLOW_COPY_to_public_section.cpp
ReplaceITK_NULLPTRMacroNames.sh
ReplaceITK_OVERRIDEMacroNames.sh
ReplaceitkGetObjectMacro.sh
UpdateAllC++Headers.sh
UseNativeC++Syntax.sh
misc-unused-parameters.sh
modernize-loop-convert.sh
modernize-pass-by-value.sh
modernize-return-braced-init-list.sh
modernize-use-auto.sh
modernize-use-bool-literals.sh
modernize-use-default-member-init.sh
modernize-use-emplace.sh
modernize-use-equals-default.sh
modernize-use-equals-delete.sh
modernize-use-nullptr.sh
modernize-use-override.sh
performance-general.sh
prefer-type-alias-over-typedef.sh
prefer_constexpr_for_const_literals.sh
readability-container-size-empty.sh
replaceClassWithTypename.py
replace_itkStaticConstMacro.sh
replace_vnl_math_XXX.sh
Accessing outdated ITKv5 migration scripts¶
git worktree add .../ITKv5.4 v5.4.0
ls ../ITKv5/Utilities/ITKv5Preparation
Class changes¶
The Clone() member function of itk::PointSet now does a “deep copy” of its
data, creating a new instance that has a copy of the points, the point data and
the region information properties of the original PointSet object. With previous
ITK versions, PointSet::Clone() did not copy any data. (It previously just
created a default-constructed PointSet object, like PointSet::CreateAnother()
does.)
For the sake of code readability, a new CoordinateType alias is added for
each nested CoordRepType alias. The old CoordRepType aliases will still be
available with ITK 6.0, but it is recommended to use CoordinateType instead.
The CoordRepType aliases will be removed when ITK_FUTURE_LEGACY_REMOVE is
enabled. Similarly, InputCoordinateType, OutputCoordinateType, and
ImagePointCoordinateType replace InputCoordRepType, OutputCoordRepType,
and ImagePointCoordRepType, respectively.
FEM LinearSystemWrapperDenseVNL public type aliases removed¶
These two nested type aliases are removed from the public interface of
fem::LinearSystemWrapperDenseVNL:
using MatrixRepresentation = vnl_matrix<Float>;
using MatrixHolder = std::vector<MatrixRepresentation *>;
ITKVNLInstantiation library is removed¶
The usage of ITKVNLInstantiation library should be replaced directly with the ITKVNL module. The ITKVNLInstantiation library was an empty library used for compatibility and provided transitive linking to ITKVNL.
itk::Math::sqrteps changed by 1 ULP for accuracy¶
itk::Math::sqrteps is now 0x1p-26, the exact sqrt(eps);
vnl_math::sqrteps (1.490116119384766e-08) was 1 ULP higher.
Prefer itk::AnatomicalOrientation over itk::SpatialOrientation¶
The enumeration defined in itk::SpatialOrientation, including itk::SpatialOrientation::CoordinateTerms,
itk::SpatialOrientation::CoordinateMajornessTerms, and itk::SpatialOrientation::ValidCoordinateOrientations may be
deprecated in the future. These terms described an orientation of an axis in a coordinate system, by the “from” or
“negative” direction. For example CoordinateTerms::ITK_COORDINATE_Right describes an axis moving from the “Right” to
the “Left” side of the body.
This is the opposite of the DICOM Patient Orientation (0020,0020) tag.
The itk::AnatomicalOrientation class now represents the anatomical orientation. The class can provide a representation
of itself as an unambiguous enumeration, string and a matrix. It provides both a PositiveEnum and a NegativeEnum for
three letter descriptions of the anatomical orientation, which is ambiguous.
The recommended unambiguous way to define an anatomical orientation is the following:
itk::AnatomicalOrientation(itk::AnatomicalOrientation::CoordinateEnum::RightToLeft,
itk::AnatomicalOrientation::CoordinateEnum::AnteriorToPosterior,
itk::AnatomicalOrientation::CoordinateEnum::InferiorToSuperior);
The itk::SpatialOrientation::ValidCoordinateOrientations enumerations can be explicitly or implicitly converted to the
new AnatomicalOrientation object:
itk::AnatomicalOrientation orientation = itk::SpatialOrientation::ValidCoordinateOrientations::ITK_COORDINATE_ORIENTATION_RAI;
orientation = itk::AnatomicalOrientation(itk::SpatialOrientation::ValidCoordinateOrientations::ITK_COORDINATE_ORIENTATION_RIP;)
Implicit conversion is not available in ITK Python. An error messaging like the following can occur:
TypeError: in method 'itkOrientImageFilterID3ID3_SetDesiredCoordinateOrientation', argument 2 of type 'itkAnatomicalOrientation'
Explicitly converting to AnatomicalOrientation is required in ITK Python:
reoriented_image = itk.orient_image_filter(
image,
use_image_direction=True,
desired_coordinate_orientation=itk.AnatomicalOrientation(
itk.SpatialOrientationEnums.ValidCoordinateOrientations_ITK_COORDINATE_ORIENTATION_RAS
),
)
GetIndex() replaced with ComputeIndex(), for iterators without index¶
ImageConstIterator::GetIndex() is marked to be removed in the future. This iterator does not have an index,
internally, so its GetIndex() member function does a potentially expensive computation. For an iterator that is
derived from ImageConstIterator, please use its newly introduced ComputeIndex() member function instead. Or use an
iterator with index (like ImageIteratorWithIndex).
Remove support for Python wrapped long double types¶
The swig wrapping of long double types into python
resulted in implicit type conversions to double,
which results in silent loss of precision.
There has never been an option for “ITK_WRAP_long_double” configuration, and manually wrapped functions in vnl for for long double were never needed by ITK based functions.
Given the undefined behavior of wrapping long double types with swig, and given that there is no use case for long double support directly from the wrapped ITK API, these manual vnl wrappings could not be reached from ITK interfaces.
The handful of manually wrapped long double functions were removed from python wrapping.
IO modules write floating-point metadata with lossless precision¶
float and double metadata values are now serialized using
itk::ConvertNumberToString(), which produces the shortest decimal string
that round-trips exactly — replacing the previous 6-significant-digit default.
Affected modules: NIfTI (scl_slope, scl_inter, pixdim, spacing,
and offset fields), MetaImage, NRRD, VoxBoCUB, SpatialObject,
and Mesh IO.
Tests that compare raw header text against golden baselines may fail because
values previously written as e.g. "1" are now written as "1.0000001".
Replace exact string matches with numeric comparisons or regenerate baselines.
Files not written by ITK are unaffected.
Legacy GoogleTest Target Names Removed¶
ITK 6 now uses the standard CMake FindGTest target names for GoogleTest libraries, aligning with upstream project and CMake.
Target Name Changes¶
Before (ITK 5):
target_link_libraries(MyTest
GTest::GTest # legacy target
GTest::Main # legacy target
)
After (ITK 6):
target_link_libraries(MyTest
GTest::gtest # compatible target
GTest::gtest_main # compatible target
)
Rationale¶
These names were deprecated in CMake 3.20 and removed in CMake 4.1.0. Additionally, the GoogleTest project itself uses the lowercase target names (GTest::gtest and GTest::gtest_main), meaning the old ITK-specific aliases were not compatible when using GoogleTest directly from its upstream repository. ITK 6 adopts the standard lowercase target names to ensure compatibility with modern CMake versions, the GoogleTest project, and consistency with other projects.
Python Global Interpreter Lock (GIL) Release¶
ITK now releases the Python Global Interpreter Lock (GIL) during C++ operations by default,
allowing for true multi-threaded execution of ITK operations from Python. This enables
parallel filter invocation when using ITK in parallel computing frameworks like Dask, Ray, or
Python’s standard threading module.
Key Changes¶
New CMake Option:
ITK_PYTHON_RELEASE_GIL(default:ON) - Controls whether the GIL is released during ITK operationsWhen enabled, the
-threadsflag is passed to SWIG to generate thread-safe wrappers
Benefits:
Multiple Python threads can execute ITK operations concurrently
Improves performance in parallel computing scenarios
Prevents thread blocking when using frameworks like Dask
Example:
import itk
import threading
image_paths = ["image1.mha", "image2.mha"]
def process_image(image_path):
# GIL is released during ITK operations
image = itk.imread(image_path)
smoothed = itk.median_image_filter(image, radius=5)
return smoothed
# Multiple threads can now execute ITK operations concurrently
threads = [
threading.Thread(target=process_image, args=(path,))
for path in image_paths
]
for thread in threads:
thread.start()
for thread in threads:
thread.join()
Note: ITK callbacks and event monitoring may be affected by GIL release. If you encounter
issues with callbacks, you can disable GIL release by setting -DITK_PYTHON_RELEASE_GIL=OFF
when building ITK.
Python lazy loading is always on¶
The Python wrapping layer’s lazy-loading mechanism has been rewritten on top
of the standard-library PEP 562
module-level __getattr__ / __dir__ protocol. Lazy first-touch resolution
is now the only mode of operation: there is no eager-import alternative.
Removed¶
itkConfig.LazyLoading— previously a boolean attribute onitkConfigthat, when set toFalse, forced every wrapped submodule to be imported atimport itktime.ITK_PYTHON_LAZYLOADING— the environment variable that seeded the value ofitkConfig.LazyLoadingat startup.The custom
itk.support.lazy.LazyITKModuletypes.ModuleTypesubclass (and itsITKLazyLoadLock,_lazy_itk_module_reconstructor, andnot_loadedsentinel) that the previous implementation relied on.
What you need to do¶
If your code or test suite contains either of the following, remove them:
itkConfig.LazyLoading = False— assignment alone is silently ignored (module attributes are mutable, but no ITK code reads the value anymore). On a freshimport itkConfigwith no prior assignment, a read ofitkConfig.LazyLoadingraisesAttributeErrorbecause the attribute is no longer defined. Delete both the writes and any reads:import itkConfig # Remove any assignments — they are silently no-ops: itkConfig.LazyLoading = False # Remove any reads — without a prior assignment they raise # AttributeError, since itkConfig no longer defines this attribute: if itkConfig.LazyLoading: ...
ITK_PYTHON_LAZYLOADING— the environment variable is no longer consulted. Any export is silently ignored; remove it from launch scripts and CI configurations:ITK_PYTHON_LAZYLOADING=0 python my_script.py # drop the env var
The observable contract from a caller’s point of view is otherwise
unchanged: import itk returns the package, attribute access on itk
or any itk.<Submodule> resolves on first touch, dir(itk) lists the
full set of lazily-discoverable attributes, and pickle / cloudpickle
round-trips through sys.modules continue to work.
Why this changed¶
The previous LazyITKModule.__getattribute__ intercepted every attribute
read on the itk package just to check a sentinel; the PEP 562 hook only
fires on a miss, so steady-state attribute access now hits the standard
module fast path. The legacy ITKLazyLoadLock also combined
threading.RLock and multiprocessing.RLock, and constructing the
multiprocessing half pinned the multiprocessing start method at
import itk time — breaking downstream callers that wanted to call
multiprocessing.set_start_method(...) after import. The replacement uses
only threading.RLock.
Modern CMake Interface Libraries¶
ITK 6 introduces modern CMake interface libraries for all ITK modules, providing improved dependency management and
simplified linking. Each module now exports a namespaced interface library following the pattern ITK::{ModuleName}Module.
Key Changes¶
Interface Library Naming:
All ITK module libraries are now exported with the
ITK::namespaceModule interface libraries follow the pattern:
ITK::{ModuleName}ModuleExample:
ITKCommonlibrary is accessed asITK::ITKCommonModule
Variable Behavior Changes:
${ModuleName}_LIBRARIESnow contains only libraries produced by that specific moduleTransitive dependencies are automatically handled through the interface library
Previously,
${ModuleName}_LIBRARIESincluded some transitive dependencies
Use Target-Specific Properties:
CMake global and directory scoped properties such and include directories and libraries have been replaced with target-specific properties.
The interface libraries handle transitive dependencies, include paths, and compiler flags automatically through CMake’s
INTERFACEproperties.Eliminates the need for
UseITK.cmakewhich modified global CMake state.
Benefits:
Cleaner, more maintainable CMake code
Prevents conflicts between different dependency requirements
Explicit dependency relationships
Better IDE integration and IntelliSense support
Application Configuration¶
The usage of ${ITK_USE_FILE} or UseITK.cmake` is now deprecated, and provides the legacy compatibility of setting global properties for include directory, CXX flags etc.
For applications using ITK, use the interface library for proper dependency linking. See the Installation Example for a complete working example.
Before:
find_package(ITK REQUIRED COMPONENTS MyModule)
include(${ITK_USE_FILE})
add_executable(Example Example.cxx)
target_link_libraries(Example ${ITK_LIBRARIES})
After:
find_package(ITK REQUIRED COMPONENTS MyModule)
itk_generate_factory_registration()
add_executable(Example Example.cxx)
target_link_libraries(Example ITK::MyModuleModule)
Alternatively, the CMake variable ITK_INTERFACE_LIBRARIES can be used to link against all loaded modules.
Factory Registration with Meta-Modules:
ITK uses factory registration to load IO formats and pluggable components at runtime. ITK 6 uses meta-modules that simplify factory registration by grouping related modules together.
Factory meta-modules include:
ITKImageIO- All image IO modules (JPEG, PNG, NIFTI, DICOM, etc.)ITKMeshIO- All mesh IO modulesITKTransformIO- Transform IO modulesITKFFTImageFilterInit- FFT implementations
To register all factories:
find_package(ITK REQUIRED COMPONENTS ITKCommon ITKImageIO )
itk_generate_factory_registration( )
...
target_link_libraries(MyTarget ${ITK_INTERFACE_LIBRARIES} ITK::ITKImageIO)
The itk_generate_factory_registration() macro must be called before adding executables. It optionally takes the factory types to register. The CMake macro generates registration code for the IO modules and factory-enabled components from the modules loaded in find_package(). The meta-module must be linked to the target to include and enable the generated registration code.
In the above example, because ITKImageIO is provided as a required component, all available ImageIO modules are registered. Note that if no components are requested in find_package(), then all components are loaded and registered.
To explicitly register factories:
find_package(ITK REQUIRED COMPONENTS ITKCommon ITKIOGDCM ITKIONRRD )
itk_generate_factory_registration( ImageIO )
...
target_link_libraries(MyTarget ${ITK_INTERFACE_LIBRARIES} ITK::ITKImageIO)
In this example, the GDCM and NRRD ImageIO modules are explicitly loaded. Only the ImageIO factory type is registered and generated.
Determining Required Modules:
To identify which ITK modules your code depends on, use the WhatModulesITK.py utility:
python Utilities/Maintenance/WhatModulesITK.py /path/to/ITK/source file1.cxx file2.h
This script analyzes your source files and reports the required ITK modules based on the headers included.
The output lists the modules you should specify in find_package(ITK REQUIRED COMPONENTS ...).
Use the --link option to generate target_link_libraries() commands with the interface library names.
This outputs the namespaced interface libraries (e.g., ITK::ITKCommonModule) suitable for direct use in
target_link_libraries() commands.
Migration for ITK Modules¶
Remote modules should be updated to use modern CMake patterns for better integration with ITK 6.
Using itk_module_add_library¶
Replace direct add_library() calls with itk_module_add_library():
Before:
set(MyModule_SRCS
itkClass1.cxx
itkClass2.cxx
)
add_library(MyModule ${ITK_LIBRARY_BUILD_TYPE} ${MyModule_SRCS})
itk_module_link_dependencies(MyModule)
itk_module_target(MyModule)
After:
set(MyModule_SRCS
itkClass1.cxx
itkClass2.cxx
)
itk_module_add_library(MyModule ${MyModule_SRCS})
The itk_module_add_library() macro automatically:
Sets the appropriate library type (SHARED/STATIC)
Configures include directories using generator expressions
Links dependencies via the interface library
Sets up export targets
In some cases when the dependencies of the module are not needed for the library, only itk_module_target() is needed.
Using Interface Libraries for Executables¶
When linking executables or examples, use the interface library instead of module variables:
Before:
add_executable(MyExample MyExample.cxx)
target_link_libraries(MyExample ${MyModule_LIBRARIES})
After:
add_executable(MyExample MyExample.cxx)
target_link_libraries(MyExample ITK::MyModuleModule)
Backward Compatibility¶
For backward compatibility, non-namespaced aliases are created with deprecation warnings. However, new code should use the namespaced ITK:: targets exclusively.
NumPy bridge: image_from_array(arr.T) is no longer transpose-equivalent to image_from_array(arr)¶
Rationale¶
Previously, itk.image_from_array (and the underlying
itk.PyBuffer.GetImageViewFromArray / GetImageFromArray) special-cased
F-contiguous inputs, with the result that image_from_array(arr) and
image_from_array(arr.T) produced ITK images with identical sizes and
content – the transpose was effectively a no-op through this API. The
docstrings explicitly documented this behavior and pointed users at
np.reshape for any actual flipping.
This was buggy in practice: the C++ side performed its own shape reversal
on F-contiguous buffers, so size handling double-reversed and image.shape
disagreed with array.shape[::-1] for F-contiguous inputs. Round-tripping
through array_from_image(image_from_array(arr)) could silently transpose
the data.
What changed in ITKv6¶
GetImageViewFromArray / GetImageFromArray now deep-copy any
non-C-contiguous input via np.ascontiguousarray() before constructing
the image. The resulting image owns an independent buffer that does
not reference the user’s original array; mutations to the original
will no longer be visible through the image. The shape rule is now
uniform regardless of memory layout:
itk.size(image) == array.shape[::-1]
image.GetPixel((i, j, k)) == array[k, j, i]
As a consequence, arr and arr.T produce ITK images with different
(reversed) sizes – np.transpose / .T is now the correct operation
to flip the dimension order through this API.
What you need to do¶
Audit any code that relied on image_from_array(arr) and
image_from_array(arr.T) returning the same image. The most common idioms
to inspect:
Code that called
itk.image_from_array(some_view.T)expecting the transpose to be ignored. Drop the.Tto recover the prior behavior, or accept the new transposed-image semantics if that was the intent.Code that round-trips arrays through ITK: previously
array_from_image(image_from_array(arr_f))could differ fromarr_ffor F-contiguous inputs; now it equalsnp.ascontiguousarray(arr_f), preserving values exactly with the originalarr_f.shape.A new warning fires when
GetImageViewFromArrayis called on a non-C-contiguous input (since a deep-copy is performed). Setneed_contiguous=Falseto silence it if the deep-copy is intended.GetImageFromArrayalready silences the warning internally because it always deep-copies viaImageDuplicator; the kwarg is not exposed on that API.