Skip to content

Advanced Registration

These utilities expose lower-level transform estimation and application functions. They are useful for custom registration experiments, but most users should use coregix.align_image_pair() instead.

coregix.preprocess.registration

Image registration utilities.

estimate_elastix_transform(fixed_image_path, moving_image_path, *, parameter_map='rigid', force_linear_resample=False, force_nearest_resample=False, fixed_mask_path=None, moving_mask_path=None, log_to_console=False)

Estimate transform parameters that map a source image to a reference image.

Parameters:

Name Type Description Default
fixed_image_path str

Path to the reference image.

required
moving_image_path str

Path to the source image to be aligned.

required
parameter_map Union[str, Sequence[str]]

Elastix parameter map name(s). Typical values include "translation", "rigid", "affine", and "bspline".

'rigid'
force_linear_resample bool

If True, enforce linear final resampling in all loaded parameter maps (mimics explicit wrapper settings).

False
force_nearest_resample bool

If True, enforce nearest-neighbor final resampling in all loaded parameter maps.

False
fixed_mask_path Optional[str]

Optional reference-image mask path.

None
moving_mask_path Optional[str]

Optional source-image mask path.

None
log_to_console bool

If True, emit registration backend logs to stdout.

False

Returns:

Type Description
Any

ITK transform parameter object produced by elastix.

Raises:

Type Description
RuntimeError

If itk-elastix is not installed.

Source code in coregix/preprocess/registration.py
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
def estimate_elastix_transform(
    fixed_image_path: str,
    moving_image_path: str,
    *,
    parameter_map: Union[str, Sequence[str]] = "rigid",
    force_linear_resample: bool = False,
    force_nearest_resample: bool = False,
    fixed_mask_path: Optional[str] = None,
    moving_mask_path: Optional[str] = None,
    log_to_console: bool = False,
) -> Any:
    """Estimate transform parameters that map a source image to a reference image.

    Args:
        fixed_image_path: Path to the reference image.
        moving_image_path: Path to the source image to be aligned.
        parameter_map: Elastix parameter map name(s). Typical values include
            ``"translation"``, ``"rigid"``, ``"affine"``, and ``"bspline"``.
        force_linear_resample: If ``True``, enforce linear final resampling in all
            loaded parameter maps (mimics explicit wrapper settings).
        force_nearest_resample: If ``True``, enforce nearest-neighbor final
            resampling in all loaded parameter maps.
        fixed_mask_path: Optional reference-image mask path.
        moving_mask_path: Optional source-image mask path.
        log_to_console: If ``True``, emit registration backend logs to stdout.

    Returns:
        ITK transform parameter object produced by elastix.

    Raises:
        RuntimeError: If ``itk-elastix`` is not installed.
    """
    itk = _require_itk()

    fixed = itk.imread(fixed_image_path, itk.F)
    moving = itk.imread(moving_image_path, itk.F)

    parameter_object = itk.ParameterObject.New()
    map_names = [parameter_map] if isinstance(parameter_map, str) else list(parameter_map)
    for map_name in map_names:
        parameter_object.AddParameterMap(parameter_object.GetDefaultParameterMap(map_name))
    if force_linear_resample:
        for idx in range(int(parameter_object.GetNumberOfParameterMaps())):
            parameter_object.SetParameter(idx, "ResampleInterpolator", "FinalLinearInterpolator")
            parameter_object.SetParameter(idx, "DefaultPixelValue", "0")
            parameter_object.SetParameter(idx, "FinalBSplineInterpolationOrder", "1")
    if force_nearest_resample:
        for idx in range(int(parameter_object.GetNumberOfParameterMaps())):
            parameter_object.SetParameter(idx, "ResampleInterpolator", "FinalNearestNeighborInterpolator")
            parameter_object.SetParameter(idx, "DefaultPixelValue", "0")
            parameter_object.SetParameter(idx, "FinalBSplineInterpolationOrder", "0")

    kwargs = {
        "fixed_image": fixed,
        "moving_image": moving,
        "parameter_object": parameter_object,
        "log_to_console": log_to_console,
    }

    if fixed_mask_path:
        kwargs["fixed_mask"] = itk.imread(fixed_mask_path, itk.UC)
    if moving_mask_path:
        kwargs["moving_mask"] = itk.imread(moving_mask_path, itk.UC)

    _, transform_parameter_object = itk.elastix_registration_method(**kwargs)
    return transform_parameter_object

apply_elastix_transform(moving_image_path, output_image_path, transform_parameter_object, *, reference_image_path=None, log_to_console=False)

Apply a precomputed registration transform to an image and write the result.

Parameters:

Name Type Description Default
moving_image_path str

Path to the source image to warp.

required
output_image_path str

Output path for the transformed image.

required
transform_parameter_object Any

Transform parameter object from :func:estimate_elastix_transform.

required
reference_image_path Optional[str]

Optional raster whose CRS/transform are copied to output_image_path after transformix writes the data.

None
log_to_console bool

If True, emit transformix logs to stdout.

False

Returns:

Type Description
str

The written output_image_path.

Raises:

Type Description
RuntimeError

If itk-elastix is not installed.

Source code in coregix/preprocess/registration.py
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
def apply_elastix_transform(
    moving_image_path: str,
    output_image_path: str,
    transform_parameter_object: Any,
    *,
    reference_image_path: Optional[str] = None,
    log_to_console: bool = False,
) -> str:
    """Apply a precomputed registration transform to an image and write the result.

    Args:
        moving_image_path: Path to the source image to warp.
        output_image_path: Output path for the transformed image.
        transform_parameter_object: Transform parameter object from
            :func:`estimate_elastix_transform`.
        reference_image_path: Optional raster whose CRS/transform are copied to
            ``output_image_path`` after transformix writes the data.
        log_to_console: If ``True``, emit transformix logs to stdout.

    Returns:
        The written ``output_image_path``.

    Raises:
        RuntimeError: If ``itk-elastix`` is not installed.
    """
    itk = _require_itk()
    moving = itk.imread(moving_image_path, itk.F)
    transformed = itk.transformix_filter(
        moving,
        transform_parameter_object=transform_parameter_object,
        log_to_console=log_to_console,
    )
    itk.imwrite(transformed, output_image_path)
    if reference_image_path:
        import rasterio

        with rasterio.open(reference_image_path) as ref, rasterio.open(output_image_path, "r+") as out:
            out.crs = ref.crs
            out.transform = ref.transform
            if ref.nodata is not None:
                out.nodata = ref.nodata
    return output_image_path

apply_elastix_transform_array(moving_image, transform_parameter_object, *, log_to_console=False)

Apply a precomputed registration transform to an in-memory image array.

Source code in coregix/preprocess/registration.py
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
def apply_elastix_transform_array(
    moving_image: np.ndarray,
    transform_parameter_object: Any,
    *,
    log_to_console: bool = False,
) -> np.ndarray:
    """Apply a precomputed registration transform to an in-memory image array."""
    itk = _require_itk()
    moving = itk.image_from_array(np.asarray(moving_image, dtype=np.float32))
    transformed = itk.transformix_filter(
        moving,
        transform_parameter_object=transform_parameter_object,
        log_to_console=log_to_console,
    )
    return itk.array_from_image(transformed)

deformation_field_from_transform(reference_image_path, transform_parameter_object, *, output_directory=None)

Return the transformix deformation field as a NumPy array.

Source code in coregix/preprocess/registration.py
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
def deformation_field_from_transform(
    reference_image_path: str,
    transform_parameter_object: Any,
    *,
    output_directory: Optional[str] = None,
) -> np.ndarray:
    """Return the transformix deformation field as a NumPy array."""
    itk = _require_itk()
    reference = itk.imread(reference_image_path, itk.F)
    kwargs = {
        "transform_parameter_object": transform_parameter_object,
    }
    if output_directory is not None:
        kwargs["output_directory"] = output_directory
    field = itk.transformix_deformation_field(reference, **kwargs)
    return itk.array_from_image(field)

deformation_field_from_transform_region(transform_parameter_object, *, row_off, col_off, height, width, output_directory=None)

Return a deformation field for a fixed-grid subregion.

Source code in coregix/preprocess/registration.py
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
def deformation_field_from_transform_region(
    transform_parameter_object: Any,
    *,
    row_off: int,
    col_off: int,
    height: int,
    width: int,
    output_directory: Optional[str] = None,
) -> np.ndarray:
    """Return a deformation field for a fixed-grid subregion."""
    itk = _require_itk()
    reference = itk.image_from_array(np.zeros((height, width), dtype=np.float32))
    reference.SetOrigin((float(col_off), float(row_off)))
    reference.SetSpacing((1.0, 1.0))
    kwargs = {
        "transform_parameter_object": transform_parameter_object,
    }
    if output_directory is not None:
        kwargs["output_directory"] = output_directory
    field = itk.transformix_deformation_field(reference, **kwargs)
    return itk.array_from_image(field)

run_elastix_registration(fixed_image_path, moving_image_path, output_image_path, *, parameter_map='rigid', fixed_mask_path=None, moving_mask_path=None, log_to_console=False)

Register a source image to a reference image using the registration backend.

Parameters:

Name Type Description Default
fixed_image_path str

Path to the reference image.

required
moving_image_path str

Path to the source image that will be warped.

required
output_image_path str

Path where the registered image will be written.

required
parameter_map Union[str, Sequence[str]]

Elastix parameter map name(s). Common values: "rigid", "affine", "bspline".

'rigid'
fixed_mask_path Optional[str]

Optional reference-image mask path.

None
moving_mask_path Optional[str]

Optional source-image mask path.

None
log_to_console bool

Whether the registration backend should print logs to stdout.

False

Returns:

Type Description
str

The written output_image_path.

Source code in coregix/preprocess/registration.py
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
def run_elastix_registration(
    fixed_image_path: str,
    moving_image_path: str,
    output_image_path: str,
    *,
    parameter_map: Union[str, Sequence[str]] = "rigid",
    fixed_mask_path: Optional[str] = None,
    moving_mask_path: Optional[str] = None,
    log_to_console: bool = False,
) -> str:
    """Register a source image to a reference image using the registration backend.

    Args:
        fixed_image_path: Path to the reference image.
        moving_image_path: Path to the source image that will be warped.
        output_image_path: Path where the registered image will be written.
        parameter_map: Elastix parameter map name(s). Common values:
            ``"rigid"``, ``"affine"``, ``"bspline"``.
        fixed_mask_path: Optional reference-image mask path.
        moving_mask_path: Optional source-image mask path.
        log_to_console: Whether the registration backend should print logs to stdout.

    Returns:
        The written ``output_image_path``.
    """
    transform_parameter_object = estimate_elastix_transform(
        fixed_image_path=fixed_image_path,
        moving_image_path=moving_image_path,
        parameter_map=parameter_map,
        fixed_mask_path=fixed_mask_path,
        moving_mask_path=moving_mask_path,
        log_to_console=log_to_console,
    )
    apply_elastix_transform(
        moving_image_path=moving_image_path,
        output_image_path=output_image_path,
        transform_parameter_object=transform_parameter_object,
        log_to_console=log_to_console,
    )
    return output_image_path