Slicer 5.11
Slicer is a multi-platform, free and open source software package for visualization and medical image computing
Loading...
Searching...
No Matches
vtkITKTransformConverter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: $RCSfile: vtkITKTransformConverter.h,v $
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or https://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
15
16#ifndef __vtkITKTransformConverter_h
17#define __vtkITKTransformConverter_h
18
20
21// VTK includes
22#include <vtkCollection.h>
23#include <vtkImageData.h>
24#include <vtkOrientedBSplineTransform.h>
25#include <vtkPoints.h>
26#include <vtkThinPlateSplineTransform.h>
27#include <vtkTransform.h>
28
29// ITK includes
30#include <itkAffineTransform.h>
31#include <itkBSplineDeformableTransform.h> // ITKv3 style
32#include <itkBSplineTransform.h> // ITKv4 style
33#include <itkCompositeTransform.h>
34#include <itkCompositeTransformIOHelper.h>
35#include <itkDisplacementFieldTransform.h>
36#include <itkIdentityTransform.h>
37#include <itkTransformFileWriter.h>
38#include <itkTransformFileReader.h>
39#include <itkImageFileReader.h>
40#include <itkImageFileWriter.h>
41#include <itkTranslationTransform.h>
42#include <itkScaleTransform.h>
43#include <itkTransformFactory.h>
44#include <itkThinPlateSplineKernelTransform.h>
45
46// Constants and typedefs
47
48static const unsigned int VTKDimension = 3;
49
50static const int BSPLINE_TRANSFORM_ORDER = 3;
51
52typedef itk::TransformFileWriter TransformWriterType;
53
54typedef itk::DisplacementFieldTransform<double, 3> DisplacementFieldTransformDoubleType;
55typedef DisplacementFieldTransformDoubleType::DisplacementFieldType GridImageDoubleType;
56
57typedef itk::ThinPlateSplineKernelTransform<double, 3> ThinPlateSplineTransformDoubleType;
58
60
62{
63public:
65
66 template <typename T>
67 static vtkAbstractTransform* CreateVTKTransformFromITK(vtkObject* loggerObject,
68 typename itk::TransformBaseTemplate<T>::Pointer transformItk,
69 double center_LocalRAS[3] = nullptr);
70
81 static itk::Object::Pointer CreateITKTransformFromVTK(vtkObject* loggerObject,
82 vtkAbstractTransform* transformVtk,
83 itk::Object::Pointer& secondaryTransformItk,
84 int preferITKv3CompatibleTransforms,
85 bool initialize = true,
86 double center_LocalRAS[3] = nullptr);
87
88 template <typename T>
89 static bool SetVTKBSplineFromITKv3Generic(vtkObject* loggerObject,
90 vtkOrientedBSplineTransform* bsplineVtk,
91 typename itk::TransformBaseTemplate<T>::Pointer warpTransformItk,
92 typename itk::TransformBaseTemplate<T>::Pointer bulkTransformItk);
93
94 template <typename T>
95 static bool SetVTKOrientedGridTransformFromITKImage(vtkObject* loggerObject,
96 vtkOrientedGridTransform* grid_Ras,
97 typename itk::DisplacementFieldTransform<T, 3>::DisplacementFieldType::Pointer gridImage_Lps);
98 static bool SetITKImageFromVTKOrientedGridTransform(vtkObject* loggerObject, GridImageDoubleType::Pointer& gridImage_Lps, vtkOrientedGridTransform* grid_Ras);
99
100 // Convert a 2D ITK matrix to a 3D ITK matrix.
101 // Spacing for the third dimension is geometric mean of the other two dimensions.
102 // This is chosen to minimize geometric distortion.
103 // Sign of the third dimension is chosen to create a right-handed coordinate system.
104 template <typename T>
105 static itk::Matrix<T, 3, 3> Matrix2Dto3D(itk::Matrix<T, 2, 2> m2D);
106
107protected:
108 template <typename T>
109 static typename itk::AffineTransform<T, 3>::Pointer ConvertITKLinearTransformFrom2Dto3D(typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS);
110 template <typename T>
111 static bool SetVTKLinearTransformFromITK(vtkObject* loggerObject,
112 vtkMatrix4x4* transformVtk_RAS,
113 typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS,
114 double center_LocalRAS[3] = nullptr);
115 static bool SetITKLinearTransformFromVTK(vtkObject* loggerObject, itk::Object::Pointer& transformItk_LPS, vtkMatrix4x4* transformVtk_RAS, double center_LocalRAS[3] = nullptr);
116
117 template <typename T>
118 static bool SetVTKOrientedGridTransformFrom2DITKImage(vtkObject* loggerObject,
119 vtkOrientedGridTransform* grid_Ras,
120 typename itk::DisplacementFieldTransform<T, 2>::DisplacementFieldType::Pointer gridImage_Lps);
121
122 template <typename T, unsigned Dimension>
123 static bool SetVTKOrientedGridTransformFromITK(vtkObject* loggerObject,
124 vtkOrientedGridTransform* transformVtk_RAS,
125 typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS);
126 static bool SetITKOrientedGridTransformFromVTK(vtkObject* loggerObject, itk::Object::Pointer& transformItk_LPS, vtkOrientedGridTransform* transformVtk_RAS);
127
128 static bool SetITKv3BSplineFromVTK(vtkObject* loggerObject,
129 itk::Object::Pointer& warpTransformItk,
130 itk::Object::Pointer& bulkTransformItk,
131 vtkOrientedBSplineTransform* bsplineVtk,
132 bool alwaysAddBulkTransform);
133 static bool SetITKv4BSplineFromVTK(vtkObject* loggerObject, itk::Object::Pointer& warpTransformItk, vtkOrientedBSplineTransform* bsplineVtk);
134
135 template <typename T, unsigned VDimension>
136 static bool SetVTKThinPlateSplineTransformFromITK(vtkObject* loggerObject,
137 vtkThinPlateSplineTransform* transformVtk_RAS,
138 typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS);
139 static bool SetITKThinPlateSplineTransformFromVTK(vtkObject* loggerObject,
140 itk::Object::Pointer& transformItk_LPS,
141 vtkThinPlateSplineTransform* transformVtk_RAS,
142 bool initialize = true);
143
144 static bool IsIdentityMatrix(vtkMatrix4x4* matrix);
145
146 template <typename BSplineTransformType>
147 static bool SetVTKBSplineParametersFromITKGeneric(vtkObject* loggerObject,
148 vtkOrientedBSplineTransform* bsplineVtk,
149 typename itk::TransformBaseTemplate<typename BSplineTransformType::ScalarType>::Pointer warpTransformItk);
150 template <typename T>
151 static bool SetVTKBSplineFromITKv4Generic(vtkObject* loggerObject, vtkOrientedBSplineTransform* bsplineVtk, typename itk::TransformBaseTemplate<T>::Pointer warpTransformItk);
152
153 template <typename BSplineTransformType>
154 static bool SetITKBSplineParametersFromVTKGeneric(vtkObject* loggerObject,
155 typename itk::Transform<typename BSplineTransformType::ScalarType, VTKDimension, VTKDimension>::Pointer& warpTransformItk,
156 vtkOrientedBSplineTransform* bsplineVtk);
157 template <typename T>
158 static bool SetITKv3BSplineFromVTKGeneric(vtkObject* loggerObject,
159 typename itk::Transform<T, VTKDimension, VTKDimension>::Pointer& warpTransformItk,
160 typename itk::Transform<T, VTKDimension, VTKDimension>::Pointer& bulkTransformItk,
161 vtkOrientedBSplineTransform* bsplineVtk,
162 bool alwaysAddBulkTransform);
163 template <typename T>
164 static bool SetITKv4BSplineFromVTKGeneric(vtkObject* loggerObject,
165 typename itk::Transform<T, VTKDimension, VTKDimension>::Pointer& warpTransformItk,
166 vtkOrientedBSplineTransform* bsplineVtk);
167};
168
169//----------------------------------------------------------------------------
171{
172 itk::TransformFactory<InverseDisplacementFieldTransformFloatType>::RegisterTransform();
173 itk::TransformFactory<InverseDisplacementFieldTransformDoubleType>::RegisterTransform();
174
175 itk::TransformFactory<InverseBSplineTransformFloatITKv3Type>::RegisterTransform();
176 itk::TransformFactory<InverseBSplineTransformFloatITKv4Type>::RegisterTransform();
177 itk::TransformFactory<InverseBSplineTransformDoubleITKv3Type>::RegisterTransform();
178 itk::TransformFactory<InverseBSplineTransformDoubleITKv4Type>::RegisterTransform();
179
180 itk::TransformFactory<InverseThinPlateSplineTransformFloatType>::RegisterTransform();
181 itk::TransformFactory<InverseThinPlateSplineTransformDoubleType>::RegisterTransform();
182
183 itk::TransformFactory<itk::InverseThinPlateSplineKernelTransform<float, 2>>::RegisterTransform();
184 itk::TransformFactory<itk::InverseThinPlateSplineKernelTransform<double, 2>>::RegisterTransform();
185
186 typedef itk::ThinPlateSplineKernelTransform<float, 3> ThinPlateSplineTransformFloatType;
187 typedef itk::ThinPlateSplineKernelTransform<double, 3> ThinPlateSplineTransformDoubleType;
188
189 // by default they are not registered
190 itk::TransformFactory<ThinPlateSplineTransformFloatType>::RegisterTransform();
191 itk::TransformFactory<ThinPlateSplineTransformDoubleType>::RegisterTransform();
192
193 itk::TransformFactory<itk::ThinPlateSplineKernelTransform<float, 2>>::RegisterTransform();
194 itk::TransformFactory<itk::ThinPlateSplineKernelTransform<double, 2>>::RegisterTransform();
195}
196
197//----------------------------------------------------------------------------
198template <typename T>
199itk::Matrix<T, 3, 3> vtkITKTransformConverter::Matrix2Dto3D(itk::Matrix<T, 2, 2> m2D)
200{
201 itk::Matrix<T, 3, 3> direction3d;
202 direction3d.SetIdentity();
203 for (unsigned int i = 0; i < 2; i++)
204 {
205 for (unsigned int j = 0; j < 2; j++)
206 {
207 direction3d[i][j] = m2D(i, j);
208 }
209 }
210
211 // compute scale for Z direction
212 itk::Vector<T, VTKDimension> directionCosine{};
213 directionCosine[0] = direction3d[0][0];
214 directionCosine[1] = direction3d[0][1];
215 auto scaleX = directionCosine.GetNorm();
216 directionCosine[0] = direction3d[1][0];
217 directionCosine[1] = direction3d[1][1];
218 auto scaleY = directionCosine.GetNorm();
219 T scaleZ = std::sqrt(scaleX * scaleY); // geometric mean
220 auto det = vnl_determinant(m2D.GetVnlMatrix());
221 direction3d(2, 2) = (det < 0.0) ? -scaleZ : scaleZ;
222
223 return direction3d;
224}
225
226//----------------------------------------------------------------------------
227template <typename T>
228typename itk::AffineTransform<T, 3>::Pointer vtkITKTransformConverter::ConvertITKLinearTransformFrom2Dto3D(typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS)
229{
230 using Affine3D = itk::AffineTransform<T, 3>;
231 using LinearTransformType = itk::MatrixOffsetTransformBase<T, 2, 2>;
232
233 typename LinearTransformType::ConstPointer t2d = dynamic_cast<const LinearTransformType*>(transformItk_LPS.GetPointer());
234 if (t2d == nullptr)
235 {
236 throw std::runtime_error("vtkITKTransformConverter::ConvertITKLinearTransformFrom2Dto3D: Unable to cast ITK transform to 2D linear transform");
237 }
238
239 // convert into 3D affine which Slicer can handle
240 typename Affine3D::MatrixType m = vtkITKTransformConverter::Matrix2Dto3D<T>(t2d->GetMatrix());
241 typename Affine3D::TranslationType t;
242 t.Fill(0);
243 typename Affine3D::InputPointType c;
244 c.Fill(0);
245 for (unsigned int i = 0; i < 2; i++)
246 {
247 t[i] = t2d->GetTranslation()[i];
248 c[i] = t2d->GetCenter()[i];
249 }
250
251 typename Affine3D::Pointer aTr = Affine3D::New();
252 aTr->SetCenter(c);
253 aTr->SetMatrix(m);
254 aTr->SetTranslation(t);
255 return aTr;
256}
257
258//----------------------------------------------------------------------------
259template <typename T>
261 vtkMatrix4x4* transformVtk_RAS,
262 typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS,
263 double center_LocalRAS[3] /*=nullptr*/)
264{
265 static const unsigned int D = VTKDimension;
266 typedef itk::MatrixOffsetTransformBase<T, D, D> LinearTransformType;
267 typedef itk::TranslationTransform<T, D> TranslateTransformType;
268
269 vtkSmartPointer<vtkMatrix4x4> transformVtk_LPS = vtkSmartPointer<vtkMatrix4x4>::New();
270 const unsigned itkDim = transformItk_LPS->GetOutputSpaceDimension();
271
272 bool convertedToVtkMatrix = false;
273
274 std::string itkTransformClassName = transformItk_LPS->GetNameOfClass();
275
276 // Linear transform of doubles or floats, dimension 2 or 3
277
278 // ITKIO transform libraries are built as shared and dynamic_cast
279 // can NOT be used with templated classes that are
280 // instantiated in a translation unit different than the one where they are
281 // defined. It will work only if the classes are explicitly instantiated
282 // and exported.
283 // To workaround the issue, instead of using dynamic_cast:
284 // (1) to ensure the objects are of the right type string comparison is done
285 // (2) static_cast is used instead of dynamic_cast.
286 // See InsightSoftwareConsortium/ITK@d1e9fe2
287 // and see https://stackoverflow.com/questions/8024010/why-do-template-class-functions-have-to-be-declared-in-the-same-translation-unit
288 //
289 // The disadvantage of this approach is that each supported class name has to be explicitly listed here and if the class hierarchy changes in ITK
290 // then the static cast may produce invalid results. Also, even if the transform class name is matching,
291 // we may cast the transform to a wrong type due to mismatch in dimensions (not 3) or data type (not double or float).
292
293 if (itkTransformClassName.find("AffineTransform") != std::string::npos || //
294 itkTransformClassName == "MatrixOffsetTransformBase" || //
295 itkTransformClassName == "Rigid3DTransform" || //
296 itkTransformClassName == "Euler3DTransform" || //
297 itkTransformClassName == "CenteredEuler3DTransform" || //
298 itkTransformClassName == "QuaternionRigidTransform" || //
299 itkTransformClassName == "VersorTransform" || //
300 itkTransformClassName == "VersorRigid3DTransform" || //
301 itkTransformClassName == "ScaleSkewVersor3DTransform" || //
302 itkTransformClassName == "ScaleVersor3DTransform" || //
303 itkTransformClassName == "Similarity3DTransform" || //
304 itkTransformClassName == "ScaleTransform" || //
305 itkTransformClassName == "ScaleLogarithmicTransform")
306 {
307 typename LinearTransformType::Pointer dlt;
308 if (itkDim == 2)
309 {
310 dlt = ConvertITKLinearTransformFrom2Dto3D<T>(transformItk_LPS);
311 }
312 else
313 {
314 dlt = static_cast<LinearTransformType*>(transformItk_LPS.GetPointer());
315 }
316
317 convertedToVtkMatrix = true;
318 for (unsigned int i = 0; i < D; i++)
319 {
320 for (unsigned int j = 0; j < D; j++)
321 {
322 transformVtk_LPS->SetElement(i, j, dlt->GetMatrix()[i][j]);
323 }
324 transformVtk_LPS->SetElement(i, D, dlt->GetOffset()[i]);
325 }
326
327 if (center_LocalRAS)
328 {
329 auto center_LocalLPS = dlt->GetCenter();
330 center_LocalRAS[0] = -center_LocalLPS[0];
331 center_LocalRAS[1] = -center_LocalLPS[1];
332 center_LocalRAS[2] = center_LocalLPS[2];
333 }
334 }
335
336 // Identity transform of doubles or floats, any dimension
337 if (itkTransformClassName == "IdentityTransform")
338 {
339 // nothing to do, matrix is already the identity
340 convertedToVtkMatrix = true;
341 }
342
343 // Translate transform of doubles or floats, dimension 2 or 3
344 // Separated out as it does not inherit from MatrixOffsetTransformBase
345 if (itkTransformClassName == "TranslationTransform")
346 {
347 if (itkDim == 2)
348 {
349 using Translate2DTransformType = itk::TranslationTransform<T, 2>;
350 typename Translate2DTransformType::Pointer dtt = static_cast<Translate2DTransformType*>(transformItk_LPS.GetPointer());
351 for (unsigned int i = 0; i < 2; i++)
352 {
353 transformVtk_LPS->SetElement(i, D, dtt->GetOffset()[i]);
354 }
355 // third element remains zero
356 }
357 else
358 {
359 typename TranslateTransformType::Pointer dtt = static_cast<TranslateTransformType*>(transformItk_LPS.GetPointer());
360 for (unsigned int i = 0; i < D; i++)
361 {
362 transformVtk_LPS->SetElement(i, D, dtt->GetOffset()[i]);
363 }
364 }
365 convertedToVtkMatrix = true;
366 }
367
368 // Convert from LPS (ITK) to RAS (Slicer)
369 //
370 // Tras = lps2ras * Tlps * ras2lps
371 //
372
373 vtkSmartPointer<vtkMatrix4x4> lps2ras = vtkSmartPointer<vtkMatrix4x4>::New();
374 lps2ras->SetElement(0, 0, -1);
375 lps2ras->SetElement(1, 1, -1);
376 vtkMatrix4x4* ras2lps = lps2ras; // lps2ras is diagonal therefore the inverse is identical
377
378 vtkMatrix4x4::Multiply4x4(lps2ras, transformVtk_LPS, transformVtk_LPS);
379 vtkMatrix4x4::Multiply4x4(transformVtk_LPS, ras2lps, transformVtk_RAS);
380
381 if (center_LocalRAS)
382 {
383 for (int i = 0; i < 3; ++i)
384 {
385 // Apply the offset to the center of transformation to account for the expected behavior of the ITK transform center.
386 // For image registration, the center of the transform is the center of the fixed image.
387 // See ITK software guide section 3.6.2.
388 double offset_RAS = transformVtk_RAS->GetElement(i, 3);
389 center_LocalRAS[i] = center_LocalRAS[i] + offset_RAS;
390 }
391 }
392
393 return convertedToVtkMatrix;
394}
395
396//----------------------------------------------------------------------------
398 itk::Object::Pointer& transformItk_LPS,
399 vtkMatrix4x4* transformVtk_RAS,
400 double center_LocalRAS[3] /*=nullptr*/)
401{
402 typedef itk::AffineTransform<double, VTKDimension> AffineTransformType;
403
404 if (transformVtk_RAS == nullptr)
405 {
406 vtkErrorWithObjectMacro(loggerObject, "vtkITKTransformConverter::SetITKLinearTransformFromVTK failed: invalid input transform");
407 return false;
408 }
409
410 vtkSmartPointer<vtkMatrix4x4> lps2ras = vtkSmartPointer<vtkMatrix4x4>::New();
411 lps2ras->SetElement(0, 0, -1);
412 lps2ras->SetElement(1, 1, -1);
413 vtkMatrix4x4* ras2lps = lps2ras; // lps2ras is diagonal therefore the inverse is identical
414
415 // Convert from RAS (Slicer) to LPS (ITK)
416 //
417 // Tlps = ras2lps * Tras * lps2ras
418 //
419 vtkSmartPointer<vtkMatrix4x4> vtkmat = vtkSmartPointer<vtkMatrix4x4>::New();
420
421 vtkMatrix4x4::Multiply4x4(ras2lps, transformVtk_RAS, vtkmat);
422 vtkMatrix4x4::Multiply4x4(vtkmat, lps2ras, vtkmat);
423
424 typedef AffineTransformType::MatrixType MatrixType;
425 typedef AffineTransformType::OutputVectorType OffsetType;
426
427 MatrixType itkmat;
428 OffsetType itkoffset;
429
430 for (unsigned int i = 0; i < VTKDimension; i++)
431 {
432 for (unsigned int j = 0; j < VTKDimension; j++)
433 {
434 itkmat[i][j] = vtkmat->GetElement(i, j);
435 }
436 itkoffset[i] = vtkmat->GetElement(i, VTKDimension);
437 }
438
439 AffineTransformType::Pointer affine = AffineTransformType::New();
440 // Matrix and offset are stored last as they represent desired transform of the VTK rotation and translation.
441 if (center_LocalRAS)
442 {
443 // For itk::AffineTransform, the center is not normally specified with matrix/offset, however it can be set explicitly before setting the matrix/offset.
444 // When using matrix/offset it is only safe to explicitly set the center of transformation before setting the matrix/offset. Setting the center after
445 // the offset would change the contents of the transform.
446 // Translation will be recomputed from the offset to represent the specified center of the transform.
447 // For more information, see the documentation of itk::AffineTransform::SetCenter.
448
449 // Convert the center from the local RAS to the parent RAS
450 double center_LocalLPS[3] = { -center_LocalRAS[0], -center_LocalRAS[1], center_LocalRAS[2] };
451
452 // Account for the offset to the center of transformation to account for the expected behavior of the ITK transform center.
453 // For image registration, the center of the transform is the center of the fixed image.
454 // See ITK software guide section 3.6.2.
455 vtkMath::Subtract(center_LocalLPS, itkoffset, center_LocalLPS);
456 affine->SetCenter(center_LocalLPS);
457 }
458 affine->SetMatrix(itkmat);
459 affine->SetOffset(itkoffset);
460
461 transformItk_LPS = affine;
462 return true;
463}
464
465//----------------------------------------------------------------------------
467{
468 static double identity[16] = { 1, 0, 0, 0, 0, 1, 0, 0, 0, 0, 1, 0, 0, 0, 0, 1 };
469 int i, j;
470
471 for (i = 0; i < 4; i++)
472 {
473 for (j = 0; j < 4; j++)
474 {
475 if (matrix->GetElement(i, j) != identity[4 * i + j])
476 {
477 return false;
478 }
479 }
480 }
481 return true;
482}
483
484//----------------------------------------------------------------------------
485template <typename BSplineTransformType>
487 vtkOrientedBSplineTransform* bsplineVtk,
488 typename itk::TransformBaseTemplate<typename BSplineTransformType::ScalarType>::Pointer warpTransformItk)
489{
490 //
491 // this version uses the itk::BSplineTransform not the itk::BSplineDeformableTransform
492 //
493 typedef typename BSplineTransformType::ScalarType T;
494 if (bsplineVtk == nullptr)
495 {
496 vtkErrorWithObjectMacro(loggerObject, "vtkMRMLTransformStorageNode::SetVTKBSplineFromITKv4 failed: bsplineVtk is invalid");
497 return false;
498 }
499 bool isDoublePrecisionInput = true;
500 if (sizeof(T) == sizeof(float))
501 {
502 isDoublePrecisionInput = false;
503 }
504 else if (sizeof(T) == sizeof(double))
505 {
506 isDoublePrecisionInput = true;
507 }
508 else
509 {
510 vtkErrorWithObjectMacro(loggerObject, "Unsupported scalar type in BSpline transform file (only float and double are supported)");
511 return false;
512 }
513
514 typename BSplineTransformType::Pointer bsplineItk = BSplineTransformType::New();
515 std::string warpTransformItkName = warpTransformItk->GetNameOfClass();
516 std::string requestedWarpTransformItkName = bsplineItk->GetNameOfClass();
517 if (warpTransformItkName != requestedWarpTransformItkName)
518 {
519 return false;
520 }
521 const unsigned itkDim = warpTransformItk->GetOutputSpaceDimension();
522 constexpr unsigned instDim = BSplineTransformType::SpaceDimension;
523 if (itkDim != instDim)
524 {
525 vtkErrorWithObjectMacro(loggerObject, "Actual and instantiated dimensionality differ: actual = " << itkDim << ", instantiated = " << instDim);
526 return false;
527 }
528 bsplineItk = static_cast<BSplineTransformType*>(warpTransformItk.GetPointer());
529
530 // now get the fixed parameters and map them to the vtk analogs
531
532 // it turns out that for a BSplineTransform, the TransformDomain information
533 // is not populated when reading from a file, so instead we access these
534 // fields from one of the coefficient images (they are assumed to be
535 // identical)
536 const typename BSplineTransformType::CoefficientImageArray coefficientImages = bsplineItk->GetCoefficientImages();
537
538 // * mesh size X, Y, Z (including the BSPLINE_TRANSFORM_ORDER=3 boundary nodes,
539 // 1 before and 2 after the grid)
540 typename BSplineTransformType::MeshSizeType meshSize = coefficientImages[0]->GetLargestPossibleRegion().GetSize();
541
542 // * mesh origin X, Y, Z (position of the boundary node before the grid)
543 typename BSplineTransformType::OriginType origin = coefficientImages[0]->GetOrigin();
544
545 // * mesh spacing X, Y, Z
546 typename BSplineTransformType::SpacingType spacing = coefficientImages[0]->GetSpacing();
547
548 // * mesh direction 3x3 matrix (first row, second row, third row)
549 typename BSplineTransformType::DirectionType direction = coefficientImages[0]->GetDirection();
550
551 itk::Matrix<T, 3, 3> direction3d;
552 if constexpr (instDim == 2)
553 {
554 direction3d = Matrix2Dto3D(direction);
555 }
556 else
557 {
558 direction3d = direction;
559 }
560
561 vtkNew<vtkMatrix4x4> gridDirectionMatrix_LPS;
562 for (unsigned int row = 0; row < VTKDimension; row++)
563 {
564 for (unsigned int column = 0; column < VTKDimension; column++)
565 {
566 gridDirectionMatrix_LPS->SetElement(row, column, direction3d[row][column]);
567 }
568 }
569
570 // Set bspline grid and coefficients
571 bsplineVtk->SetBorderModeToZero();
572
573 vtkNew<vtkImageData> bsplineCoefficients;
574
575 if constexpr (instDim == 2)
576 {
577 bsplineCoefficients->SetExtent(0, meshSize[0] - 1, 0, meshSize[1] - 1, 0, 0);
578 bsplineCoefficients->SetSpacing(spacing[0], spacing[1], std::sqrt(spacing[0] * spacing[1]));
579 }
580 else
581 {
582 bsplineCoefficients->SetExtent(0, meshSize[0] - 1, 0, meshSize[1] - 1, 0, meshSize[2] - 1);
583 bsplineCoefficients->SetSpacing(spacing[0], spacing[1], spacing[2]);
584 }
585
586 // convert the direction matrix from LPS (itk) to RAS (slicer)
587 vtkNew<vtkMatrix4x4> lpsToRas;
588 lpsToRas->SetElement(0, 0, -1);
589 lpsToRas->SetElement(1, 1, -1);
590
591 vtkNew<vtkMatrix4x4> rasToLps;
592 rasToLps->SetElement(0, 0, -1);
593 rasToLps->SetElement(1, 1, -1);
594
595 vtkNew<vtkMatrix4x4> gridDirectionMatrix_RAS;
596 vtkMatrix4x4::Multiply4x4(lpsToRas.GetPointer(), gridDirectionMatrix_LPS.GetPointer(), gridDirectionMatrix_RAS.GetPointer());
597 bsplineVtk->SetGridDirectionMatrix(gridDirectionMatrix_RAS.GetPointer());
598
599 // convert the origin from LPS (itk) to RAS (slicer)
600 double gridOrigin_RAS[3] = { -origin[0], -origin[1], 0 };
601 if constexpr (instDim == 3)
602 {
603 gridOrigin_RAS[2] = origin[2];
604 }
605 bsplineCoefficients->SetOrigin(gridOrigin_RAS);
606
607 int bsplineCoefficientsScalarType = VTK_FLOAT;
608 if (isDoublePrecisionInput)
609 {
610 bsplineCoefficientsScalarType = VTK_DOUBLE;
611 }
612
613 bsplineCoefficients->AllocateScalars(bsplineCoefficientsScalarType, 3);
614
615 const unsigned int expectedNumberOfVectors = meshSize.CalculateProductOfElements();
616 const unsigned int expectedNumberOfParameters = expectedNumberOfVectors * instDim;
617 const unsigned int actualNumberOfParameters = bsplineItk->GetNumberOfParameters();
618
619 if (actualNumberOfParameters != expectedNumberOfParameters)
620 {
621 vtkErrorWithObjectMacro(loggerObject, "Mismatch in number of BSpline parameters in the transform file and the MRML node");
622 return false;
623 }
624 const T* itkBSplineParams_LPS = static_cast<const T*>(bsplineItk->GetParameters().data_block());
625 T* vtkBSplineParams_RAS = static_cast<T*>(bsplineCoefficients->GetScalarPointer());
626 for (unsigned int i = 0; i < expectedNumberOfVectors; i++)
627 {
628 *(vtkBSplineParams_RAS++) = -(*(itkBSplineParams_LPS));
629 *(vtkBSplineParams_RAS++) = -(*(itkBSplineParams_LPS + expectedNumberOfVectors));
630 if constexpr (instDim == 3)
631 {
632 *(vtkBSplineParams_RAS++) = (*(itkBSplineParams_LPS + expectedNumberOfVectors * 2));
633 }
634 else // 2D
635 {
636 *(vtkBSplineParams_RAS++) = 0; // Z component
637 }
638 itkBSplineParams_LPS++;
639 }
640 bsplineVtk->SetCoefficientData(bsplineCoefficients.GetPointer());
641
642 // Success
643 return true;
644}
645
646//----------------------------------------------------------------------------
647template <typename T>
649 vtkOrientedBSplineTransform* bsplineVtk,
650 typename itk::TransformBaseTemplate<T>::Pointer warpTransformItk,
651 typename itk::TransformBaseTemplate<T>::Pointer bulkTransformItk)
652{
653 if (bsplineVtk == nullptr)
654 {
655 vtkErrorWithObjectMacro(loggerObject, "vtkMRMLTransformStorageNode::SetVTKBSplineFromITK failed: bsplineVtk is invalid");
656 return false;
657 }
658
659 const unsigned int itkDim = warpTransformItk->GetOutputSpaceDimension();
660 if (itkDim < 2 || itkDim > 3)
661 {
662 vtkErrorWithObjectMacro(loggerObject, "Unsupported number of dimensions in BSpline transform file (expected = 2 or 3, actual = " << itkDim << ")");
663 return false;
664 }
665 bool inverse = false;
666 if (itkDim == 2)
667 {
668 // inverse class is derived from forward class, so it has to be checked first
670 {
671 inverse = true;
672 }
673 else if (SetVTKBSplineParametersFromITKGeneric<itk::BSplineDeformableTransform<T, 2, BSPLINE_TRANSFORM_ORDER>>(loggerObject, bsplineVtk, warpTransformItk))
674 {
675 inverse = false;
676 }
677 else
678 {
679 vtkDebugWithObjectMacro(loggerObject, "Not an ITKv3 BSpline transform");
680 return false;
681 }
682 }
683 else
684 {
685 // inverse class is derived from forward class, so it has to be checked first
687 {
688 inverse = true;
689 }
690 else if (SetVTKBSplineParametersFromITKGeneric<itk::BSplineDeformableTransform<T, 3, BSPLINE_TRANSFORM_ORDER>>(loggerObject, bsplineVtk, warpTransformItk))
691 {
692 inverse = false;
693 }
694 else
695 {
696 vtkDebugWithObjectMacro(loggerObject, "Not an ITKv3 BSpline transform");
697 return false;
698 }
699 }
700
701 // Set the bulk transform
702 if (bulkTransformItk)
703 {
704 std::string bulkTransformItkTransformName = bulkTransformItk->GetNameOfClass();
705
706 if (bulkTransformItkTransformName == "AffineTransform")
707 {
708 vtkNew<vtkMatrix4x4> bulkMatrix_RAS;
709 if (!SetVTKLinearTransformFromITK<T>(loggerObject, bulkMatrix_RAS, bulkTransformItk))
710 {
711 vtkErrorWithObjectMacro(loggerObject, "Cannot read the bulk transform in BSplineTransform");
712 return false;
713 }
714 bsplineVtk->SetBulkTransformMatrix(bulkMatrix_RAS.GetPointer());
715 }
716 else if (bulkTransformItkTransformName == "IdentityTransform")
717 {
718 // bulk transform is identity, which is equivalent to no bulk transform
719 }
720 else
721 {
722 vtkErrorWithObjectMacro(loggerObject, "Cannot read the 2nd transform in BSplineTransform (expected AffineTransform_double_3_3 or IdentityTransform)");
723 return false;
724 }
725 }
726
727 if (inverse)
728 {
729 bsplineVtk->Inverse();
730 }
731
732 // Success
733 return true;
734}
735
736//----------------------------------------------------------------------------
737template <typename T>
739 vtkOrientedBSplineTransform* bsplineVtk,
740 typename itk::TransformBaseTemplate<T>::Pointer warpTransformItk)
741{
742 bool inverse = false;
743 const unsigned int itkDim = warpTransformItk->GetOutputSpaceDimension();
744
745 if (itkDim == 2)
746 {
747 // inverse class is derived from forward class, so it has to be checked first
749 {
750 inverse = true;
751 }
752 else if (SetVTKBSplineParametersFromITKGeneric<itk::BSplineTransform<T, 2, BSPLINE_TRANSFORM_ORDER>>(loggerObject, bsplineVtk, warpTransformItk))
753 {
754 inverse = false;
755 }
756 else
757 {
758 vtkDebugWithObjectMacro(loggerObject, "Not an ITKv4 BSpline transform");
759 return false;
760 }
761 }
762 else
763 {
764 // inverse class is derived from forward class, so it has to be checked first
766 {
767 inverse = true;
768 }
769 else if (SetVTKBSplineParametersFromITKGeneric<itk::BSplineTransform<T, 3, BSPLINE_TRANSFORM_ORDER>>(loggerObject, bsplineVtk, warpTransformItk))
770 {
771 inverse = false;
772 }
773 else
774 {
775 vtkDebugWithObjectMacro(loggerObject, "Not an ITKv4 BSpline transform");
776 return false;
777 }
778 }
779
780 if (inverse)
781 {
782 bsplineVtk->Inverse();
783 }
784 return true;
785}
786
787//----------------------------------------------------------------------------
788template <typename BSplineTransformType>
790 vtkObject* loggerObject,
791 typename itk::Transform<typename BSplineTransformType::ScalarType, VTKDimension, VTKDimension>::Pointer& warpTransformItk,
792 vtkOrientedBSplineTransform* bsplineVtk)
793{
794 typedef typename BSplineTransformType::ScalarType T;
795 if (bsplineVtk == nullptr)
796 {
797 vtkErrorWithObjectMacro(loggerObject, "vtkMRMLTransformStorageNode::SetITKBSplineFromVTK failed: bsplineVtk is invalid");
798 return false;
799 }
800 vtkImageData* bsplineCoefficients_RAS = bsplineVtk->GetCoefficientData();
801 if (bsplineCoefficients_RAS == nullptr)
802 {
803 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file: coefficients are not specified");
804 return false;
805 }
806
807 typename BSplineTransformType::Pointer bsplineItk = BSplineTransformType::New();
808 warpTransformItk = bsplineItk;
809
810 // Fixed parameters:
811 // * mesh size X, Y, Z (including the BSPLINE_TRANSFORM_ORDER=3 boundary nodes, 1 before and 2 after the grid)
812 // * mesh origin X, Y, Z (position of the boundary node before the grid)
813 // * mesh spacing X, Y, Z
814 // * mesh direction 3x3 matrix (first row, second row, third row)
815 typename BSplineTransformType::FixedParametersType transformFixedParamsItk;
816 const unsigned int numberOfFixedParameters = VTKDimension * (VTKDimension + 3);
817 transformFixedParamsItk.SetSize(numberOfFixedParameters);
818
819 int* gridExtent = bsplineCoefficients_RAS->GetExtent();
820 int gridSize[3] = { gridExtent[1] - gridExtent[0] + 1, gridExtent[3] - gridExtent[2] + 1, gridExtent[5] - gridExtent[4] + 1 };
821 transformFixedParamsItk[0] = gridSize[0];
822 transformFixedParamsItk[1] = gridSize[1];
823 transformFixedParamsItk[2] = gridSize[2];
824
825 double* gridOrigin_RAS = bsplineCoefficients_RAS->GetOrigin();
826 double gridOrigin_LPS[3] = { -gridOrigin_RAS[0], -gridOrigin_RAS[1], gridOrigin_RAS[2] };
827 transformFixedParamsItk[3] = gridOrigin_LPS[0];
828 transformFixedParamsItk[4] = gridOrigin_LPS[1];
829 transformFixedParamsItk[5] = gridOrigin_LPS[2];
830
831 double* gridSpacing = bsplineCoefficients_RAS->GetSpacing();
832 transformFixedParamsItk[6] = gridSpacing[0];
833 transformFixedParamsItk[7] = gridSpacing[1];
834 transformFixedParamsItk[8] = gridSpacing[2];
835
836 vtkNew<vtkMatrix4x4> gridDirectionMatrix_RAS;
837 if (bsplineVtk->GetGridDirectionMatrix() != nullptr)
838 {
839 gridDirectionMatrix_RAS->DeepCopy(bsplineVtk->GetGridDirectionMatrix());
840 }
841 vtkNew<vtkMatrix4x4> lpsToRas;
842 lpsToRas->SetElement(0, 0, -1);
843 lpsToRas->SetElement(1, 1, -1);
844 vtkNew<vtkMatrix4x4> rasToLps;
845 rasToLps->SetElement(0, 0, -1);
846 rasToLps->SetElement(1, 1, -1);
847 vtkNew<vtkMatrix4x4> gridDirectionMatrix_LPS;
848 vtkMatrix4x4::Multiply4x4(rasToLps.GetPointer(), gridDirectionMatrix_RAS.GetPointer(), gridDirectionMatrix_LPS.GetPointer());
849 int fpIndex = 9;
850 for (unsigned int row = 0; row < VTKDimension; row++)
851 {
852 for (unsigned int column = 0; column < VTKDimension; column++)
853 {
854 transformFixedParamsItk[fpIndex++] = gridDirectionMatrix_LPS->GetElement(row, column);
855 }
856 }
857
858 bsplineItk->SetFixedParameters(transformFixedParamsItk);
859
860 // BSpline coefficients
861
862 const unsigned int expectedNumberOfVectors = gridSize[0] * gridSize[1] * gridSize[2];
863 const unsigned int expectedNumberOfParameters = expectedNumberOfVectors * VTKDimension;
864 if (bsplineItk->GetNumberOfParameters() != expectedNumberOfParameters)
865 {
866 vtkErrorWithObjectMacro(loggerObject, "Mismatch in number of BSpline parameters in the ITK transform and the VTK transform");
867 return false;
868 }
869
870 typename BSplineTransformType::ParametersType transformParamsItk(expectedNumberOfParameters);
871 T* itkBSplineParams_LPS = static_cast<T*>(transformParamsItk.data_block());
872 T* vtkBSplineParams_RAS = static_cast<T*>(bsplineCoefficients_RAS->GetScalarPointer());
873 double coefficientScale = bsplineVtk->GetDisplacementScale();
874 for (unsigned int i = 0; i < expectedNumberOfVectors; i++)
875 {
876 *(itkBSplineParams_LPS) = -coefficientScale * (*(vtkBSplineParams_RAS++));
877 *(itkBSplineParams_LPS + expectedNumberOfVectors) = -coefficientScale * (*(vtkBSplineParams_RAS++));
878 *(itkBSplineParams_LPS + expectedNumberOfVectors * 2) = coefficientScale * (*(vtkBSplineParams_RAS++));
879 itkBSplineParams_LPS++;
880 }
881
882 bsplineItk->SetParameters(transformParamsItk);
883 return true;
884}
885
886//----------------------------------------------------------------------------
887template <typename T>
889 typename itk::Transform<T, VTKDimension, VTKDimension>::Pointer& warpTransformItk,
890 typename itk::Transform<T, VTKDimension, VTKDimension>::Pointer& bulkTransformItk,
891 vtkOrientedBSplineTransform* bsplineVtk,
892 bool alwaysAddBulkTransform)
893{
894 if (bsplineVtk == nullptr)
895 {
896 vtkErrorWithObjectMacro(loggerObject, "vtkMRMLTransformStorageNode::SetITKBSplineFromVTK failed: bsplineVtk is invalid");
897 return false;
898 }
899 // Update is needed because it refreshes the inverse flag (the flag may be out-of-date if the transform depends on its inverse)
900 bsplineVtk->Update();
901 bool itkTransformSetSuccessfully = false;
902 if (bsplineVtk->GetInverseFlag())
903 {
904 itkTransformSetSuccessfully =
906 }
907 else
908 {
909 itkTransformSetSuccessfully = SetITKBSplineParametersFromVTKGeneric<itk::BSplineDeformableTransform<T, VTKDimension, VTKDimension>>(loggerObject, warpTransformItk, bsplineVtk);
910 }
911 if (!itkTransformSetSuccessfully)
912 {
913 vtkErrorWithObjectMacro(loggerObject, "vtkMRMLTransformStorageNode::SetITKBSplineFromVTK failed: cannot determine BSpline parameters");
914 return false;
915 }
916
917 vtkMatrix4x4* bulkMatrix_RAS = bsplineVtk->GetBulkTransformMatrix();
918 if (bulkMatrix_RAS || alwaysAddBulkTransform)
919 {
920 vtkNew<vtkMatrix4x4> lpsToRas;
921 lpsToRas->SetElement(0, 0, -1);
922 lpsToRas->SetElement(1, 1, -1);
923 vtkNew<vtkMatrix4x4> rasToLps;
924 rasToLps->SetElement(0, 0, -1);
925 rasToLps->SetElement(1, 1, -1);
926 vtkNew<vtkMatrix4x4> bulkMatrix_LPS; // bulk_LPS = rasToLps * bulk_RAS * lpsToRas
927 // If bulk transform is available then use it, otherwise just write an identity matrix (we just write it because
928 // alwaysAddBulkTransform was requested, due to backward compatibility reasons)
929 if (bulkMatrix_RAS != nullptr)
930 {
931 vtkMatrix4x4::Multiply4x4(rasToLps.GetPointer(), bulkMatrix_RAS, bulkMatrix_LPS.GetPointer());
932 vtkMatrix4x4::Multiply4x4(bulkMatrix_LPS.GetPointer(), lpsToRas.GetPointer(), bulkMatrix_LPS.GetPointer());
933 }
934 typedef itk::AffineTransform<T, VTKDimension> BulkTransformType;
935 typename BulkTransformType::Pointer affineItk = BulkTransformType::New();
936 bulkTransformItk = affineItk;
937
938 typename BulkTransformType::MatrixType affineMatrix;
939 typename BulkTransformType::OffsetType affineOffset;
940 for (unsigned int i = 0; i < VTKDimension; i++)
941 {
942 for (unsigned int j = 0; j < VTKDimension; j++)
943 {
944 affineMatrix[i][j] = bulkMatrix_LPS->GetElement(i, j);
945 }
946 affineOffset[i] = bulkMatrix_LPS->GetElement(i, VTKDimension);
947 }
948
949 affineItk->SetMatrix(affineMatrix);
950 affineItk->SetOffset(affineOffset);
951 }
952 else
953 {
954 bulkTransformItk = nullptr;
955 }
956
957 return true;
958}
959
960//----------------------------------------------------------------------------
961template <typename T>
963 typename itk::Transform<T, VTKDimension, VTKDimension>::Pointer& warpTransformItk,
964 vtkOrientedBSplineTransform* bsplineVtk)
965{
966 // Update is needed because it refreshes the inverse flag (the flag may be out-of-date if the transform depends on its inverse)
967 bsplineVtk->Update();
968 bool itkTransformSetSuccessfully = false;
969 if (bsplineVtk->GetInverseFlag())
970 {
971 itkTransformSetSuccessfully = SetITKBSplineParametersFromVTKGeneric<itk::InverseBSplineTransform<T, VTKDimension, VTKDimension>>(loggerObject, warpTransformItk, bsplineVtk);
972 }
973 else
974 {
975 itkTransformSetSuccessfully = SetITKBSplineParametersFromVTKGeneric<itk::BSplineTransform<T, VTKDimension, VTKDimension>>(loggerObject, warpTransformItk, bsplineVtk);
976 }
977 if (!itkTransformSetSuccessfully)
978 {
979 vtkErrorWithObjectMacro(loggerObject, "vtkMRMLTransformStorageNode::SetITKv4BSplineFromVTKGeneric failed: cannot determine BSpline parameters");
980 return false;
981 }
982 return true;
983}
984
985//----------------------------------------------------------------------------
987 itk::Object::Pointer& warpTransformItk,
988 itk::Object::Pointer& bulkTransformItk,
989 vtkOrientedBSplineTransform* bsplineVtk,
990 bool alwaysAddBulkTransform)
991{
992 if (bsplineVtk == nullptr)
993 {
994 vtkErrorWithObjectMacro(loggerObject, "Cannot retrieve BSpline transform from node");
995 return false;
996 }
997
998 vtkImageData* bsplineCoefficients = bsplineVtk->GetCoefficientData();
999
1000 if (bsplineCoefficients == nullptr)
1001 {
1002 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file: coefficients are not specified");
1003 return false;
1004 }
1005
1006 if (bsplineCoefficients->GetScalarType() == VTK_FLOAT)
1007 {
1008 typedef itk::Transform<float, VTKDimension, VTKDimension> ITKTransformType;
1009 ITKTransformType::Pointer floatWarpTransformItk;
1010 ITKTransformType::Pointer floatBulkTransformItk;
1011 if (!SetITKv3BSplineFromVTKGeneric<float>(loggerObject, floatWarpTransformItk, floatBulkTransformItk, bsplineVtk, alwaysAddBulkTransform))
1012 {
1013 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file");
1014 return false;
1015 }
1016 warpTransformItk = floatWarpTransformItk.GetPointer();
1017 bulkTransformItk = floatBulkTransformItk.GetPointer();
1018 }
1019 else if (bsplineCoefficients->GetScalarType() == VTK_DOUBLE)
1020 {
1021 typedef itk::Transform<double, VTKDimension, VTKDimension> ITKTransformType;
1022 ITKTransformType::Pointer doubleWarpTransformItk;
1023 ITKTransformType::Pointer doubleBulkTransformItk;
1024 if (!SetITKv3BSplineFromVTKGeneric<double>(loggerObject, doubleWarpTransformItk, doubleBulkTransformItk, bsplineVtk, alwaysAddBulkTransform))
1025 {
1026 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file");
1027 return false;
1028 }
1029 warpTransformItk = doubleWarpTransformItk;
1030 bulkTransformItk = doubleBulkTransformItk;
1031 }
1032 else
1033 {
1034 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file: only float and double coefficient types are supported");
1035 return false;
1036 }
1037
1038 return true;
1039}
1040
1041//----------------------------------------------------------------------------
1042bool vtkITKTransformConverter::SetITKv4BSplineFromVTK(vtkObject* loggerObject, itk::Object::Pointer& warpTransformItk, vtkOrientedBSplineTransform* bsplineVtk)
1043{
1044 if (bsplineVtk == nullptr)
1045 {
1046 vtkErrorWithObjectMacro(loggerObject, "Cannot retrieve BSpline transform from node");
1047 return false;
1048 }
1049
1050 vtkImageData* bsplineCoefficients = bsplineVtk->GetCoefficientData();
1051
1052 if (bsplineCoefficients == nullptr)
1053 {
1054 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file: coefficients are not specified");
1055 return false;
1056 }
1057
1058 if (bsplineCoefficients->GetScalarType() == VTK_FLOAT)
1059 {
1060 typedef itk::Transform<float, VTKDimension, VTKDimension> ITKTransformType;
1061 ITKTransformType::Pointer floatWarpTransformItk;
1062 if (!SetITKv4BSplineFromVTKGeneric<float>(loggerObject, floatWarpTransformItk, bsplineVtk))
1063 {
1064 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file");
1065 return false;
1066 }
1067 warpTransformItk = floatWarpTransformItk.GetPointer();
1068 }
1069 else if (bsplineCoefficients->GetScalarType() == VTK_DOUBLE)
1070 {
1071 typedef itk::Transform<double, VTKDimension, VTKDimension> ITKTransformType;
1072 ITKTransformType::Pointer doubleWarpTransformItk;
1073 if (!SetITKv4BSplineFromVTKGeneric<double>(loggerObject, doubleWarpTransformItk, bsplineVtk))
1074 {
1075 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file");
1076 return false;
1077 }
1078 warpTransformItk = doubleWarpTransformItk;
1079 }
1080 else
1081 {
1082 vtkErrorWithObjectMacro(loggerObject, "Cannot write BSpline transform to file: only float and double coefficient types are supported");
1083 return false;
1084 }
1085
1086 return true;
1087}
1088
1089//----------------------------------------------------------------------------
1090template <typename T, unsigned Dimension>
1092 vtkOrientedGridTransform* transformVtk_RAS,
1093 typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS)
1094{
1095 typedef itk::DisplacementFieldTransform<T, Dimension> DisplacementFieldTransformType;
1096 typedef itk::InverseDisplacementFieldTransform<T, Dimension> InverseDisplacementFieldTransformType;
1097
1098 std::string transformItkClassName = transformItk_LPS->GetNameOfClass();
1099
1100 bool inverse = false;
1101 typename DisplacementFieldTransformType::DisplacementFieldType* gridImageItk_Lps = nullptr;
1102 if (transformItkClassName == "InverseDisplacementFieldTransform") // inverse class is derived from forward class, so it has to be checked first
1103 {
1104 DisplacementFieldTransformType* inverseDisplacementFieldTransform = static_cast<InverseDisplacementFieldTransformType*>(transformItk_LPS.GetPointer());
1105 inverse = true;
1106 gridImageItk_Lps = inverseDisplacementFieldTransform->GetDisplacementField();
1107 }
1108 else if (transformItkClassName == "DisplacementFieldTransform")
1109 {
1110 DisplacementFieldTransformType* displacementFieldTransform = static_cast<DisplacementFieldTransformType*>(transformItk_LPS.GetPointer());
1111 inverse = false;
1112 gridImageItk_Lps = displacementFieldTransform->GetDisplacementField();
1113 }
1114 else
1115 {
1116 vtkDebugWithObjectMacro(loggerObject, "Not a grid transform");
1117 return false;
1118 }
1119 if constexpr (Dimension != VTKDimension)
1120 {
1121 if (!SetVTKOrientedGridTransformFrom2DITKImage<T>(loggerObject, transformVtk_RAS, gridImageItk_Lps))
1122 {
1123 return false;
1124 }
1125 }
1126 else
1127 {
1128 if (!SetVTKOrientedGridTransformFromITKImage<T>(loggerObject, transformVtk_RAS, gridImageItk_Lps))
1129 {
1130 return false;
1131 }
1132 }
1133 if (inverse)
1134 {
1135 transformVtk_RAS->Inverse();
1136 }
1137 return true;
1138}
1139
1140//----------------------------------------------------------------------------
1141bool vtkITKTransformConverter::SetITKOrientedGridTransformFromVTK(vtkObject* loggerObject, itk::Object::Pointer& transformItk_LPS, vtkOrientedGridTransform* transformVtk_RAS)
1142{
1143 GridImageDoubleType::Pointer gridImageItk_Lps;
1144 if (!SetITKImageFromVTKOrientedGridTransform(loggerObject, gridImageItk_Lps, transformVtk_RAS))
1145 {
1146 vtkErrorWithObjectMacro(loggerObject, "vtkITKTransformConverter::SetITKOrientedGridTransformFromVTK failed: input transform is invalid");
1147 return false;
1148 }
1149 // Update is needed because it refreshes the inverse flag (the flag may be out-of-date if the transform depends on its inverse)
1150 transformVtk_RAS->Update();
1151 if (transformVtk_RAS->GetInverseFlag())
1152 {
1153 InverseDisplacementFieldTransformDoubleType::Pointer gridTransformItk = InverseDisplacementFieldTransformDoubleType::New();
1154 gridTransformItk->SetDisplacementField(gridImageItk_Lps);
1155 transformItk_LPS = gridTransformItk;
1156 }
1157 else
1158 {
1159 DisplacementFieldTransformDoubleType::Pointer gridTransformItk = DisplacementFieldTransformDoubleType::New();
1160 gridTransformItk->SetDisplacementField(gridImageItk_Lps);
1161 transformItk_LPS = gridTransformItk;
1162 }
1163 return true;
1164}
1165
1166//----------------------------------------------------------------------------
1167template <typename T>
1169 vtkOrientedGridTransform* grid_Ras,
1170 typename itk::DisplacementFieldTransform<T, 3>::DisplacementFieldType::Pointer gridImage_Lps)
1171{
1172 typedef itk::DisplacementFieldTransform<T, 3> DisplacementFieldTransformType;
1173 typedef typename DisplacementFieldTransformType::DisplacementFieldType GridImageType;
1174
1175 vtkNew<vtkImageData> gridImage_Ras;
1176
1177 // Origin
1178 gridImage_Ras->SetOrigin(-gridImage_Lps->GetOrigin()[0], -gridImage_Lps->GetOrigin()[1], gridImage_Lps->GetOrigin()[2]);
1179
1180 // Spacing
1181 gridImage_Ras->SetSpacing(gridImage_Lps->GetSpacing()[0], gridImage_Lps->GetSpacing()[1], gridImage_Lps->GetSpacing()[2]);
1182
1183 // Direction
1184 vtkNew<vtkMatrix4x4> gridDirectionMatrix_LPS;
1185 for (unsigned int row = 0; row < VTKDimension; row++)
1186 {
1187 for (unsigned int column = 0; column < VTKDimension; column++)
1188 {
1189 gridDirectionMatrix_LPS->SetElement(row, column, gridImage_Lps->GetDirection()(row, column));
1190 }
1191 }
1192 vtkNew<vtkMatrix4x4> lpsToRas;
1193 lpsToRas->SetElement(0, 0, -1);
1194 lpsToRas->SetElement(1, 1, -1);
1195 vtkNew<vtkMatrix4x4> gridDirectionMatrix_RAS;
1196 vtkMatrix4x4::Multiply4x4(lpsToRas.GetPointer(), gridDirectionMatrix_LPS.GetPointer(), gridDirectionMatrix_RAS.GetPointer());
1197 grid_Ras->SetGridDirectionMatrix(gridDirectionMatrix_RAS.GetPointer());
1198
1199 // Vectors
1200 typename GridImageType::SizeType size = gridImage_Lps->GetBufferedRegion().GetSize();
1201 gridImage_Ras->SetDimensions(size[0], size[1], size[2]);
1202 unsigned int numberOfScalarComponents = GridImageType::PixelType::Dimension;
1203 if (numberOfScalarComponents != VTKDimension)
1204 {
1205 vtkErrorWithObjectMacro(loggerObject,
1206 "Cannot load grid transform: the input displacement field expected to contain " << VTKDimension << " components but it actually contains "
1207 << numberOfScalarComponents);
1208 return false;
1209 }
1210 gridImage_Ras->AllocateScalars(VTK_DOUBLE, 3);
1211
1212 double* displacementVectors_Ras = reinterpret_cast<double*>(gridImage_Ras->GetScalarPointer());
1213 itk::ImageRegionConstIterator<GridImageType> inputIt(gridImage_Lps, gridImage_Lps->GetRequestedRegion());
1214 while (!inputIt.IsAtEnd())
1215 {
1216 typename GridImageType::PixelType displacementVectorLps = inputIt.Get();
1217 *(displacementVectors_Ras++) = -displacementVectorLps[0];
1218 *(displacementVectors_Ras++) = -displacementVectorLps[1];
1219 *(displacementVectors_Ras++) = displacementVectorLps[2];
1220 ++inputIt;
1221 }
1222
1223 grid_Ras->SetDisplacementGridData(gridImage_Ras.GetPointer());
1224
1225 // Set the interpolation to cubic to have smooth derivatives
1226 grid_Ras->SetInterpolationModeToCubic();
1227
1228 return true;
1229}
1230
1231//----------------------------------------------------------------------------
1232template <typename T>
1234 vtkOrientedGridTransform* grid_Ras,
1235 typename itk::DisplacementFieldTransform<T, 2>::DisplacementFieldType::Pointer gridImage_Lps)
1236{
1237 typedef itk::DisplacementFieldTransform<T, 2> DisplacementFieldTransformType;
1238 typedef typename DisplacementFieldTransformType::DisplacementFieldType GridImageType;
1239
1240 vtkNew<vtkImageData> gridImage_Ras;
1241
1242 // Origin (convert LPS -> RAS for X and Y, set Z to 0)
1243 const auto& itkOrigin = gridImage_Lps->GetOrigin();
1244 double originRasX = -itkOrigin[0];
1245 double originRasY = -itkOrigin[1];
1246 double originRasZ = 0.0;
1247 gridImage_Ras->SetOrigin(originRasX, originRasY, originRasZ);
1248
1249 // Spacing (X,Y from ITK, Z=geometric mean)
1250 const auto& itkSpacing = gridImage_Lps->GetSpacing();
1251 gridImage_Ras->SetSpacing(itkSpacing[0], itkSpacing[1], std::sqrt(itkSpacing[0] * itkSpacing[1]));
1252
1253 itk::Matrix<itk::SpacePrecisionType, 3, 3> direction3d = vtkITKTransformConverter::Matrix2Dto3D(gridImage_Lps->GetDirection());
1254
1255 vtkNew<vtkMatrix4x4> gridDirectionMatrix_LPS;
1256 for (unsigned int row = 0; row < 3; ++row)
1257 {
1258 for (unsigned int col = 0; col < 3; ++col)
1259 {
1260 gridDirectionMatrix_LPS->SetElement(row, col, direction3d(row, col));
1261 }
1262 }
1263
1264 // Convert direction from LPS to RAS
1265 vtkNew<vtkMatrix4x4> lpsToRas;
1266 lpsToRas->SetElement(0, 0, -1);
1267 lpsToRas->SetElement(1, 1, -1);
1268 vtkNew<vtkMatrix4x4> gridDirectionMatrix_RAS;
1269 vtkMatrix4x4::Multiply4x4(lpsToRas.GetPointer(), gridDirectionMatrix_LPS.GetPointer(), gridDirectionMatrix_RAS.GetPointer());
1270 grid_Ras->SetGridDirectionMatrix(gridDirectionMatrix_RAS.GetPointer());
1271
1272 // Expand from 2->3 components
1273 typename GridImageType::SizeType size = gridImage_Lps->GetBufferedRegion().GetSize();
1274 int dimX = static_cast<int>(size[0]);
1275 int dimY = static_cast<int>(size[1]);
1276 gridImage_Ras->SetDimensions(dimX, dimY, 1);
1277
1278 unsigned int numberOfScalarComponents = GridImageType::PixelType::Dimension;
1279 if (numberOfScalarComponents != 2)
1280 {
1281 vtkErrorWithObjectMacro(loggerObject,
1282 "Cannot load 2D grid transform: the input displacement field expected to contain 2 components but it actually contains " << numberOfScalarComponents);
1283 return false;
1284 }
1285
1286 if constexpr (std::is_same_v<T, float>)
1287 {
1288 gridImage_Ras->AllocateScalars(VTK_FLOAT, 3);
1289 }
1290 else if constexpr (std::is_same_v<T, double>)
1291 {
1292 gridImage_Ras->AllocateScalars(VTK_DOUBLE, 3);
1293 }
1294 else
1295 {
1296 vtkErrorWithObjectMacro(loggerObject, "Cannot load 2D grid transform: unsupported pixel type. Only float and double are supported.");
1297 return false;
1298 }
1299
1300 // Vectors: copy displacement vectors, converting LPS -> RAS and expand from 2->3 components
1301 T* displacementVectors_Ras = reinterpret_cast<T*>(gridImage_Ras->GetScalarPointer());
1302 itk::ImageRegionConstIterator<GridImageType> inputIt(gridImage_Lps, gridImage_Lps->GetRequestedRegion());
1303 while (!inputIt.IsAtEnd())
1304 {
1305 typename GridImageType::PixelType displacementVectorLps = inputIt.Get();
1306 // Convert LPS -> RAS for X and Y, Z is zero for 2D
1307 *(displacementVectors_Ras++) = -static_cast<T>(displacementVectorLps[0]);
1308 *(displacementVectors_Ras++) = -static_cast<T>(displacementVectorLps[1]);
1309 *(displacementVectors_Ras++) = 0.0;
1310 ++inputIt;
1311 }
1312
1313 grid_Ras->SetDisplacementGridData(gridImage_Ras.GetPointer());
1314
1315 // Use cubic interpolation for smooth derivatives
1316 grid_Ras->SetInterpolationModeToCubic();
1317
1318 return true;
1319}
1320
1321//----------------------------------------------------------------------------
1322bool vtkITKTransformConverter::SetITKImageFromVTKOrientedGridTransform(vtkObject* loggerObject, GridImageDoubleType::Pointer& gridImage_Lps, vtkOrientedGridTransform* grid_Ras)
1323{
1324 if (grid_Ras == nullptr)
1325 {
1326 vtkErrorWithObjectMacro(loggerObject, "Cannot save grid transform: the input vtkOrientedGridTransform is invalid");
1327 return false;
1328 }
1329
1330 // Update is needed because DisplacementGrid may be out-of-date if the transform depends on its inverse
1331 grid_Ras->Update();
1332
1333 vtkImageData* gridImage_Ras = grid_Ras->GetDisplacementGrid();
1334 if (gridImage_Ras == nullptr)
1335 {
1336 vtkErrorWithObjectMacro(loggerObject, "Cannot save grid transform: the input vtkOrientedGridTransform does not contain a valid displacement grid");
1337 return false;
1338 }
1339 if (gridImage_Ras->GetNumberOfScalarComponents() != static_cast<int>(VTKDimension))
1340 {
1341 vtkErrorWithObjectMacro(loggerObject,
1342 "Cannot save grid transform: the input vtkOrientedGridTransform expected to contain " << VTKDimension << " components but it actually contains "
1343 << gridImage_Ras->GetNumberOfScalarComponents());
1344 return false;
1345 }
1346
1347 gridImage_Lps = GridImageDoubleType::New();
1348
1349 // Origin
1350 double* origin_Ras = gridImage_Ras->GetOrigin();
1351 double origin_Lps[3] = { -origin_Ras[0], -origin_Ras[1], origin_Ras[2] };
1352 gridImage_Lps->SetOrigin(origin_Lps);
1353
1354 // Spacing
1355 double* spacing = gridImage_Ras->GetSpacing();
1356 // GridType::SpacingType spacing( spacing );
1357 gridImage_Lps->SetSpacing(spacing);
1358
1359 // Direction
1360 vtkNew<vtkMatrix4x4> gridDirectionMatrix_Ras;
1361 if (grid_Ras->GetGridDirectionMatrix() != nullptr)
1362 {
1363 gridDirectionMatrix_Ras->DeepCopy(grid_Ras->GetGridDirectionMatrix());
1364 }
1365 vtkNew<vtkMatrix4x4> rasToLps;
1366 rasToLps->SetElement(0, 0, -1);
1367 rasToLps->SetElement(1, 1, -1);
1368 vtkNew<vtkMatrix4x4> gridDirectionMatrix_Lps;
1369 vtkMatrix4x4::Multiply4x4(rasToLps.GetPointer(), gridDirectionMatrix_Ras.GetPointer(), gridDirectionMatrix_Lps.GetPointer());
1370 GridImageDoubleType::DirectionType gridDirectionMatrixItk_Lps;
1371 for (unsigned int row = 0; row < VTKDimension; row++)
1372 {
1373 for (unsigned int column = 0; column < VTKDimension; column++)
1374 {
1375 gridDirectionMatrixItk_Lps(row, column) = gridDirectionMatrix_Lps->GetElement(row, column);
1376 }
1377 }
1378 gridImage_Lps->SetDirection(gridDirectionMatrixItk_Lps);
1379
1380 // Vectors
1381 GridImageDoubleType::IndexType start;
1382 start[0] = start[1] = start[2] = 0;
1383 int* Nijk = gridImage_Ras->GetDimensions();
1384 GridImageDoubleType::SizeType size;
1385 size[0] = Nijk[0];
1386 size[1] = Nijk[1];
1387 size[2] = Nijk[2];
1388 GridImageDoubleType::RegionType region;
1389 region.SetSize(size);
1390 region.SetIndex(start);
1391 gridImage_Lps->SetRegions(region);
1392 gridImage_Lps->Allocate();
1393 itk::ImageRegionIterator<GridImageDoubleType> gridImageIt_Lps(gridImage_Lps, region);
1394 gridImageIt_Lps.GoToBegin();
1395 GridImageDoubleType::PixelType displacementVectorLps;
1396 double displacementScale = grid_Ras->GetDisplacementScale();
1397 double displacementShift = grid_Ras->GetDisplacementShift();
1398
1399 if (gridImage_Ras->GetScalarType() == VTK_DOUBLE)
1400 {
1401 double* displacementVectors_Ras = reinterpret_cast<double*>(gridImage_Ras->GetScalarPointer());
1402 while (!gridImageIt_Lps.IsAtEnd())
1403 {
1404 displacementVectorLps[0] = -(displacementScale * (*(displacementVectors_Ras++)) + displacementShift);
1405 displacementVectorLps[1] = -(displacementScale * (*(displacementVectors_Ras++)) + displacementShift);
1406 displacementVectorLps[2] = (displacementScale * (*(displacementVectors_Ras++)) + displacementShift);
1407 gridImageIt_Lps.Set(displacementVectorLps);
1408 ++gridImageIt_Lps;
1409 }
1410 }
1411 else if (gridImage_Ras->GetScalarType() == VTK_FLOAT)
1412 {
1413 float* displacementVectors_Ras = reinterpret_cast<float*>(gridImage_Ras->GetScalarPointer());
1414 while (!gridImageIt_Lps.IsAtEnd())
1415 {
1416 displacementVectorLps[0] = -(displacementScale * (*(displacementVectors_Ras++)) + displacementShift);
1417 displacementVectorLps[1] = -(displacementScale * (*(displacementVectors_Ras++)) + displacementShift);
1418 displacementVectorLps[2] = (displacementScale * (*(displacementVectors_Ras++)) + displacementShift);
1419 gridImageIt_Lps.Set(displacementVectorLps);
1420 ++gridImageIt_Lps;
1421 }
1422 }
1423 else
1424 {
1425 vtkErrorWithObjectMacro(loggerObject, "Cannot save grid transform: only float and double scalar types are supported");
1426 return false;
1427 }
1428 return true;
1429}
1430
1431//----------------------------------------------------------------------------
1432template <typename T, unsigned VDimension>
1434 vtkThinPlateSplineTransform* transformVtk_RAS,
1435 typename itk::TransformBaseTemplate<T>::Pointer transformItk_LPS)
1436{
1437 typedef itk::ThinPlateSplineKernelTransform<T, VDimension> ThinPlateSplineTransformType;
1438 typedef itk::InverseThinPlateSplineKernelTransform<T, VDimension> InverseThinPlateSplineTransformType;
1439
1440 std::string transformItkClassName = transformItk_LPS->GetNameOfClass();
1441
1442 bool inverse = false;
1443 typename ThinPlateSplineTransformType::PointSetType::Pointer sourceLandmarksItk_Lps;
1444 typename ThinPlateSplineTransformType::PointSetType::Pointer targetLandmarksItk_Lps;
1445 if (transformItkClassName == "InverseThinPlateSplineKernelTransform") // inverse class is derived from forward class, so it has to be checked first
1446 {
1447 ThinPlateSplineTransformType* inverseTpsTransform = static_cast<InverseThinPlateSplineTransformType*>(transformItk_LPS.GetPointer());
1448 inverse = true;
1449 sourceLandmarksItk_Lps = inverseTpsTransform->GetSourceLandmarks();
1450 targetLandmarksItk_Lps = inverseTpsTransform->GetTargetLandmarks();
1451 }
1452 else if (transformItkClassName == "ThinPlateSplineKernelTransform")
1453 {
1454 ThinPlateSplineTransformType* tpsTransform = static_cast<ThinPlateSplineTransformType*>(transformItk_LPS.GetPointer());
1455 inverse = false;
1456 sourceLandmarksItk_Lps = tpsTransform->GetSourceLandmarks();
1457 targetLandmarksItk_Lps = tpsTransform->GetTargetLandmarks();
1458 }
1459 else
1460 {
1461 vtkDebugWithObjectMacro(loggerObject, "Not a ThinPlateSpline transform");
1462 return false;
1463 }
1464
1465 vtkNew<vtkPoints> sourceLandmarksVtk_Ras;
1466 unsigned int numberOfSourceLandmarks = sourceLandmarksItk_Lps->GetNumberOfPoints();
1467 vtkNew<vtkPoints> targetLandmarksVtk_Ras;
1468 unsigned int numberOfTargetLandmarks = targetLandmarksItk_Lps->GetNumberOfPoints();
1469 if constexpr (std::is_same_v<T, float>)
1470 {
1471 sourceLandmarksVtk_Ras->SetDataTypeToFloat();
1472 sourceLandmarksVtk_Ras->Allocate(numberOfSourceLandmarks);
1473
1474 targetLandmarksVtk_Ras->SetDataTypeToFloat();
1475 targetLandmarksVtk_Ras->Allocate(numberOfTargetLandmarks);
1476 }
1477 else if constexpr (std::is_same_v<T, double>)
1478 {
1479 sourceLandmarksVtk_Ras->SetDataTypeToDouble();
1480 sourceLandmarksVtk_Ras->Allocate(numberOfSourceLandmarks);
1481
1482 targetLandmarksVtk_Ras->SetDataTypeToDouble();
1483 targetLandmarksVtk_Ras->Allocate(numberOfTargetLandmarks);
1484 }
1485 else
1486 {
1487 vtkErrorWithObjectMacro(loggerObject, "Cannot load 2D TPS transform: unsupported pixel type. Only float and double are supported.");
1488 return false;
1489 }
1490
1491 for (unsigned int i = 0; i < numberOfSourceLandmarks; i++)
1492 {
1493 typename ThinPlateSplineTransformType::InputPointType pointItk_Lps;
1494 bool pointExists = sourceLandmarksItk_Lps->GetPoint(i, &pointItk_Lps);
1495 if (!pointExists)
1496 {
1497 continue;
1498 }
1499 T pointVtk_Ras[3] = { 0 };
1500 pointVtk_Ras[0] = -pointItk_Lps[0];
1501 pointVtk_Ras[1] = -pointItk_Lps[1];
1502 if constexpr (VDimension == 2)
1503 {
1504 pointVtk_Ras[2] = 0.0;
1505 }
1506 else
1507 {
1508 pointVtk_Ras[2] = pointItk_Lps[2];
1509 }
1510 sourceLandmarksVtk_Ras->InsertNextPoint(pointVtk_Ras);
1511 }
1512
1513 for (unsigned int i = 0; i < numberOfTargetLandmarks; i++)
1514 {
1515 typename ThinPlateSplineTransformType::InputPointType pointItk_Lps;
1516 bool pointExists = targetLandmarksItk_Lps->GetPoint(i, &pointItk_Lps);
1517 if (!pointExists)
1518 {
1519 continue;
1520 }
1521 T pointVtk_Ras[3] = { 0 };
1522 pointVtk_Ras[0] = -pointItk_Lps[0];
1523 pointVtk_Ras[1] = -pointItk_Lps[1];
1524 if constexpr (VDimension == 2)
1525 {
1526 pointVtk_Ras[2] = 0.0;
1527 }
1528 else
1529 {
1530 pointVtk_Ras[2] = pointItk_Lps[2];
1531 }
1532 targetLandmarksVtk_Ras->InsertNextPoint(pointVtk_Ras);
1533 }
1534
1535 transformVtk_RAS->SetBasisToR(); // need to set this because data is 3D
1536 transformVtk_RAS->SetSourceLandmarks(sourceLandmarksVtk_Ras.GetPointer());
1537 transformVtk_RAS->SetTargetLandmarks(targetLandmarksVtk_Ras.GetPointer());
1538
1539 if (inverse)
1540 {
1541 transformVtk_RAS->Inverse();
1542 }
1543 return true;
1544}
1545
1546//----------------------------------------------------------------------------
1548 itk::Object::Pointer& transformItk_LPS,
1549 vtkThinPlateSplineTransform* transformVtk_RAS,
1550 bool initialize /*= true*/)
1551{
1552 if (transformVtk_RAS == nullptr)
1553 {
1554 vtkErrorWithObjectMacro(loggerObject, "Cannot set ITK thin-plate spline transform from VTK: the input vtkThinPlateSplineTransform is invalid");
1555 return false;
1556 }
1557
1558 // Update is needed because the Basis value and Inverse flag may be out-of-date if the transform depends on its inverse
1559 transformVtk_RAS->Update();
1560
1561 if (transformVtk_RAS->GetBasis() != VTK_RBF_R)
1562 {
1563 vtkErrorWithObjectMacro(loggerObject,
1564 "Cannot set ITK thin-plate spline transform from VTK: basis function must be R."
1565 " Call SetBasisToR() method of the vtkThinPlateSplineTransform object before attempting to write it to file.");
1566 return false;
1567 }
1568
1569 ThinPlateSplineTransformDoubleType::PointSetType::Pointer sourceLandmarksItk_Lps = ThinPlateSplineTransformDoubleType::PointSetType::New();
1570 vtkPoints* sourceLandmarksVtk_Ras = transformVtk_RAS->GetSourceLandmarks();
1571 if (sourceLandmarksVtk_Ras != nullptr)
1572 {
1573 for (int i = 0; i < sourceLandmarksVtk_Ras->GetNumberOfPoints(); i++)
1574 {
1575 double posVtk_Ras[3] = { 0 };
1576 sourceLandmarksVtk_Ras->GetPoint(i, posVtk_Ras);
1577 ThinPlateSplineTransformDoubleType::InputPointType posItk_Lps;
1578 posItk_Lps[0] = -posVtk_Ras[0];
1579 posItk_Lps[1] = -posVtk_Ras[1];
1580 posItk_Lps[2] = posVtk_Ras[2];
1581 sourceLandmarksItk_Lps->GetPoints()->InsertElement(i, posItk_Lps);
1582 }
1583 }
1584 ThinPlateSplineTransformDoubleType::PointSetType::Pointer targetLandmarksItk_Lps = ThinPlateSplineTransformDoubleType::PointSetType::New();
1585 vtkPoints* targetLandmarksVtk_Ras = transformVtk_RAS->GetTargetLandmarks();
1586 if (targetLandmarksVtk_Ras != nullptr)
1587 {
1588 for (int i = 0; i < targetLandmarksVtk_Ras->GetNumberOfPoints(); i++)
1589 {
1590 double posVtk_Ras[3] = { 0 };
1591 targetLandmarksVtk_Ras->GetPoint(i, posVtk_Ras);
1592 ThinPlateSplineTransformDoubleType::InputPointType posItk_Lps;
1593 posItk_Lps[0] = -posVtk_Ras[0];
1594 posItk_Lps[1] = -posVtk_Ras[1];
1595 posItk_Lps[2] = posVtk_Ras[2];
1596 targetLandmarksItk_Lps->GetPoints()->InsertElement(i, posItk_Lps);
1597 }
1598 }
1599
1600 if (transformVtk_RAS->GetInverseFlag())
1601 {
1602 InverseThinPlateSplineTransformDoubleType::Pointer tpsTransformItk = InverseThinPlateSplineTransformDoubleType::New();
1603 tpsTransformItk->SetSourceLandmarks(sourceLandmarksItk_Lps);
1604 tpsTransformItk->SetTargetLandmarks(targetLandmarksItk_Lps);
1605 if (initialize)
1606 {
1607 tpsTransformItk->ComputeWMatrix();
1608 }
1609 transformItk_LPS = tpsTransformItk;
1610 }
1611 else
1612 {
1613 ThinPlateSplineTransformDoubleType::Pointer tpsTransformItk = ThinPlateSplineTransformDoubleType::New();
1614 tpsTransformItk->SetSourceLandmarks(sourceLandmarksItk_Lps);
1615 tpsTransformItk->SetTargetLandmarks(targetLandmarksItk_Lps);
1616 if (initialize)
1617 {
1618 tpsTransformItk->ComputeWMatrix();
1619 }
1620 transformItk_LPS = tpsTransformItk;
1621 }
1622 return true;
1623}
1624
1625//----------------------------------------------------------------------------
1626template <typename T>
1627vtkAbstractTransform* vtkITKTransformConverter::CreateVTKTransformFromITK(vtkObject* loggerObject,
1628 typename itk::TransformBaseTemplate<T>::Pointer transformItk,
1629 double center_LocalRAS[3] /*=nullptr*/)
1630{
1631 bool conversionSuccess = false;
1632
1633 // Perform these sanity checks early
1634 if (!transformItk)
1635 {
1636 vtkErrorWithObjectMacro(loggerObject, "CreateVTKTransformFromITK failed: nullptr provided");
1637 return nullptr;
1638 }
1639 const unsigned itkDim = transformItk->GetOutputSpaceDimension();
1640 if (itkDim != transformItk->GetInputSpaceDimension())
1641 {
1642 vtkErrorWithObjectMacro(loggerObject,
1643 "We do not support transforms with different number of input and output dimensions, " << "input dimension = " << transformItk->GetInputSpaceDimension()
1644 << ", output dimension = " << itkDim << ")");
1645 return nullptr;
1646 }
1647 if (itkDim < 2 || itkDim > 3) // do we need to support 4D transforms?
1648 {
1649 vtkErrorWithObjectMacro(loggerObject, "Unsupported transform dimension: " << itkDim);
1650 return nullptr;
1651 }
1652
1653 std::string itkTransformClassName = transformItk->GetNameOfClass();
1654
1655 // Composite Transform
1656 if (itkTransformClassName == "CompositeTransform")
1657 {
1658 // Composition is handled by the caller (vtkMRMLTransformStorageNode)
1659 return nullptr;
1660 }
1661
1662 // Linear
1663 vtkNew<vtkMatrix4x4> transformMatrixVtk;
1664 conversionSuccess = SetVTKLinearTransformFromITK<T>(loggerObject, transformMatrixVtk.GetPointer(), transformItk, center_LocalRAS);
1665 if (conversionSuccess)
1666 {
1667 vtkNew<vtkTransform> linearTransformVtk;
1668 linearTransformVtk->SetMatrix(transformMatrixVtk.GetPointer());
1669 linearTransformVtk->Register(nullptr);
1670 return linearTransformVtk.GetPointer();
1671 }
1672
1673 // Grid
1674 vtkNew<vtkOrientedGridTransform> gridTransformVtk;
1675 if (itkDim == 3)
1676 {
1677 conversionSuccess = SetVTKOrientedGridTransformFromITK<T, 3>(loggerObject, gridTransformVtk.GetPointer(), transformItk);
1678 }
1679 else if (itkDim == 2)
1680 {
1681 conversionSuccess = SetVTKOrientedGridTransformFromITK<T, 2>(loggerObject, gridTransformVtk.GetPointer(), transformItk);
1682 }
1683 if (conversionSuccess)
1684 {
1685 gridTransformVtk->Register(nullptr);
1686 return gridTransformVtk.GetPointer();
1687 }
1688
1689 // BSpline
1690 vtkNew<vtkOrientedBSplineTransform> bsplineTransformVtk;
1691 conversionSuccess = SetVTKBSplineFromITKv4Generic<T>(loggerObject, bsplineTransformVtk.GetPointer(), transformItk);
1692 if (conversionSuccess)
1693 {
1694 bsplineTransformVtk->Register(nullptr);
1695 return bsplineTransformVtk.GetPointer();
1696 }
1697
1698 // ThinPlateSpline
1699 vtkNew<vtkThinPlateSplineTransform> tpsTransformVtk;
1700 if (itkDim == 3)
1701 {
1702 conversionSuccess = SetVTKThinPlateSplineTransformFromITK<T, 3>(loggerObject, tpsTransformVtk.GetPointer(), transformItk);
1703 }
1704 else if (itkDim == 2)
1705 {
1706 conversionSuccess = SetVTKThinPlateSplineTransformFromITK<T, 2>(loggerObject, tpsTransformVtk.GetPointer(), transformItk);
1707 }
1708 if (conversionSuccess)
1709 {
1710 tpsTransformVtk->Register(nullptr);
1711 return tpsTransformVtk.GetPointer();
1712 }
1713
1714 return nullptr;
1715}
1716
1717//----------------------------------------------------------------------------
1718itk::Object::Pointer vtkITKTransformConverter::CreateITKTransformFromVTK(vtkObject* loggerObject,
1719 vtkAbstractTransform* transformVtk,
1720 itk::Object::Pointer& secondaryTransformItk,
1721 int preferITKv3CompatibleTransforms,
1722 bool initialize /*= true*/,
1723 double center_LocalRAS[3] /*=nullptr*/)
1724{
1725 typedef itk::CompositeTransform<double> CompositeTransformType;
1726
1727 if (transformVtk == nullptr)
1728 {
1729 vtkErrorWithObjectMacro(loggerObject, "CreateITKTransformFromVTK failed: invalid VTK transform");
1730 return nullptr;
1731 }
1732 vtkNew<vtkCollection> transformList;
1733 vtkMRMLTransformNode::FlattenGeneralTransform(transformList.GetPointer(), transformVtk);
1734 if (transformList->GetNumberOfItems() == 0)
1735 {
1736 // no transformation means identity transform
1737 vtkNew<vtkTransform> identity;
1738 transformList->AddItem(identity.GetPointer());
1739 }
1740
1741 itk::Object::Pointer primaryTransformItk;
1742 if (transformList->GetNumberOfItems() == 1)
1743 {
1744 // Simple, non-composite transform (one item in the input transform)
1745 vtkObject* singleTransformVtk = transformList->GetItemAsObject(0);
1746 // Linear
1747 if (vtkHomogeneousTransform::SafeDownCast(singleTransformVtk))
1748 {
1749 vtkHomogeneousTransform* linearTransformVtk = vtkHomogeneousTransform::SafeDownCast(singleTransformVtk);
1750 vtkMatrix4x4* transformMatrix = linearTransformVtk->GetMatrix();
1751 if (!SetITKLinearTransformFromVTK(loggerObject, primaryTransformItk, transformMatrix, center_LocalRAS))
1752 {
1753 // conversion failed
1754 return nullptr;
1755 }
1756 return primaryTransformItk;
1757 }
1758 // BSpline
1759 else if (vtkOrientedBSplineTransform::SafeDownCast(singleTransformVtk))
1760 {
1761 vtkOrientedBSplineTransform* bsplineTransformVtk = vtkOrientedBSplineTransform::SafeDownCast(singleTransformVtk);
1762 vtkMatrix4x4* bulkMatrix = bsplineTransformVtk->GetBulkTransformMatrix(); // non-zero for ITKv3 bspline transform only
1763 if (preferITKv3CompatibleTransforms || (bulkMatrix != nullptr && !IsIdentityMatrix(bulkMatrix)))
1764 {
1765 if (!SetITKv3BSplineFromVTK(loggerObject, primaryTransformItk, secondaryTransformItk, bsplineTransformVtk, preferITKv3CompatibleTransforms))
1766 {
1767 // conversion failed
1768 return nullptr;
1769 }
1770 return primaryTransformItk;
1771 }
1772 else
1773 {
1774 if (!SetITKv4BSplineFromVTK(loggerObject, primaryTransformItk, bsplineTransformVtk))
1775 {
1776 // conversion failed
1777 return nullptr;
1778 }
1779 return primaryTransformItk;
1780 }
1781 }
1782 // Grid
1783 else if (vtkOrientedGridTransform::SafeDownCast(singleTransformVtk))
1784 {
1785 vtkOrientedGridTransform* gridTransformVtk = vtkOrientedGridTransform::SafeDownCast(singleTransformVtk);
1786 if (!SetITKOrientedGridTransformFromVTK(loggerObject, primaryTransformItk, gridTransformVtk))
1787 {
1788 // conversion failed
1789 return nullptr;
1790 }
1791 return primaryTransformItk;
1792 }
1793 // ThinPlateSpline
1794 else if (vtkThinPlateSplineTransform::SafeDownCast(singleTransformVtk))
1795 {
1796 vtkThinPlateSplineTransform* tpsTransformVtk = vtkThinPlateSplineTransform::SafeDownCast(singleTransformVtk);
1797 if (!SetITKThinPlateSplineTransformFromVTK(loggerObject, primaryTransformItk, tpsTransformVtk, initialize))
1798 {
1799 // conversion failed
1800 return nullptr;
1801 }
1802 return primaryTransformItk;
1803 }
1804 else
1805 {
1806 if (singleTransformVtk == nullptr)
1807 {
1808 vtkErrorWithObjectMacro(loggerObject, "vtkITKTransformConverter::CreateITKTransformFromVTK failed: invalid input transform");
1809 return nullptr;
1810 }
1811 vtkErrorWithObjectMacro(
1812 loggerObject, "vtkITKTransformConverter::CreateITKTransformFromVTK failed: conversion of transform type " << singleTransformVtk->GetClassName() << " is not supported");
1813 return nullptr;
1814 }
1815 }
1816 else
1817 {
1818 // Composite transform (more than one items in the input transform)
1819 // Create an ITK composite transform from the VTK general transform
1820 CompositeTransformType::Pointer compositeTransformItk = CompositeTransformType::New();
1821 primaryTransformItk = compositeTransformItk;
1822 // We need to iterate through the list in reverse order, as ITK CompositeTransform can only append transform on one side
1823 for (int transformIndex = transformList->GetNumberOfItems() - 1; transformIndex >= 0; --transformIndex)
1824 {
1825 vtkAbstractTransform* singleTransformVtk = vtkAbstractTransform::SafeDownCast(transformList->GetItemAsObject(transformIndex));
1826 itk::Object::Pointer secondaryTransformItkTmp;
1827 // We use ITKv4 format (PreferITKv3Transform format is set to false), because
1828 // legacy ITKv3 code does not know how to interpret composite transforms,
1829 // and also ITKv3 bspline transform with bulk component cannot be saved in a composite transform
1830 itk::Object::Pointer singleTransformItk = CreateITKTransformFromVTK(loggerObject, singleTransformVtk, secondaryTransformItkTmp, false);
1831 if (secondaryTransformItkTmp.IsNotNull())
1832 {
1833 vtkErrorWithObjectMacro(loggerObject,
1834 "vtkITKTransformConverter::CreateITKTransformFromVTK failed:"
1835 " composite transforms cannot contain legacy transforms (that contains secondary transforms)."
1836 " Do not harden transforms on legacy ITK transforms to avoid this error.");
1837 return nullptr;
1838 }
1839
1840 if (singleTransformItk.IsNull() //
1841 || std::string(singleTransformItk->GetNameOfClass()).find("Transform") == std::string::npos)
1842 {
1843 vtkErrorWithObjectMacro(loggerObject,
1844 "vtkITKTransformConverter::CreateITKTransformFromVTK failed:"
1845 " invalid element found while trying to create a composite transform");
1846 return nullptr;
1847 }
1848 CompositeTransformType::TransformType::Pointer singleTransformItkTypeChecked = static_cast<CompositeTransformType::TransformType*>(singleTransformItk.GetPointer());
1849 compositeTransformItk->AddTransform(singleTransformItkTypeChecked.GetPointer());
1850 }
1851 return primaryTransformItk;
1852 }
1853 return nullptr;
1854}
1855
1856#endif
static itk::AffineTransform< T, 3 >::Pointer ConvertITKLinearTransformFrom2Dto3D(typename itk::TransformBaseTemplate< T >::Pointer transformItk_LPS)
static bool SetVTKBSplineFromITKv3Generic(vtkObject *loggerObject, vtkOrientedBSplineTransform *bsplineVtk, typename itk::TransformBaseTemplate< T >::Pointer warpTransformItk, typename itk::TransformBaseTemplate< T >::Pointer bulkTransformItk)
static bool SetITKThinPlateSplineTransformFromVTK(vtkObject *loggerObject, itk::Object::Pointer &transformItk_LPS, vtkThinPlateSplineTransform *transformVtk_RAS, bool initialize=true)
static bool SetITKv4BSplineFromVTKGeneric(vtkObject *loggerObject, typename itk::Transform< T, VTKDimension, VTKDimension >::Pointer &warpTransformItk, vtkOrientedBSplineTransform *bsplineVtk)
static bool SetVTKLinearTransformFromITK(vtkObject *loggerObject, vtkMatrix4x4 *transformVtk_RAS, typename itk::TransformBaseTemplate< T >::Pointer transformItk_LPS, double center_LocalRAS[3]=nullptr)
static itk::Object::Pointer CreateITKTransformFromVTK(vtkObject *loggerObject, vtkAbstractTransform *transformVtk, itk::Object::Pointer &secondaryTransformItk, int preferITKv3CompatibleTransforms, bool initialize=true, double center_LocalRAS[3]=nullptr)
static bool SetITKv3BSplineFromVTKGeneric(vtkObject *loggerObject, typename itk::Transform< T, VTKDimension, VTKDimension >::Pointer &warpTransformItk, typename itk::Transform< T, VTKDimension, VTKDimension >::Pointer &bulkTransformItk, vtkOrientedBSplineTransform *bsplineVtk, bool alwaysAddBulkTransform)
static bool SetITKBSplineParametersFromVTKGeneric(vtkObject *loggerObject, typename itk::Transform< typename BSplineTransformType::ScalarType, VTKDimension, VTKDimension >::Pointer &warpTransformItk, vtkOrientedBSplineTransform *bsplineVtk)
static bool SetITKLinearTransformFromVTK(vtkObject *loggerObject, itk::Object::Pointer &transformItk_LPS, vtkMatrix4x4 *transformVtk_RAS, double center_LocalRAS[3]=nullptr)
static bool SetITKImageFromVTKOrientedGridTransform(vtkObject *loggerObject, GridImageDoubleType::Pointer &gridImage_Lps, vtkOrientedGridTransform *grid_Ras)
static bool SetVTKOrientedGridTransformFromITK(vtkObject *loggerObject, vtkOrientedGridTransform *transformVtk_RAS, typename itk::TransformBaseTemplate< T >::Pointer transformItk_LPS)
static bool SetITKOrientedGridTransformFromVTK(vtkObject *loggerObject, itk::Object::Pointer &transformItk_LPS, vtkOrientedGridTransform *transformVtk_RAS)
static bool SetVTKBSplineParametersFromITKGeneric(vtkObject *loggerObject, vtkOrientedBSplineTransform *bsplineVtk, typename itk::TransformBaseTemplate< typename BSplineTransformType::ScalarType >::Pointer warpTransformItk)
static bool SetITKv4BSplineFromVTK(vtkObject *loggerObject, itk::Object::Pointer &warpTransformItk, vtkOrientedBSplineTransform *bsplineVtk)
static bool SetVTKThinPlateSplineTransformFromITK(vtkObject *loggerObject, vtkThinPlateSplineTransform *transformVtk_RAS, typename itk::TransformBaseTemplate< T >::Pointer transformItk_LPS)
static vtkAbstractTransform * CreateVTKTransformFromITK(vtkObject *loggerObject, typename itk::TransformBaseTemplate< T >::Pointer transformItk, double center_LocalRAS[3]=nullptr)
static bool SetVTKOrientedGridTransformFromITKImage(vtkObject *loggerObject, vtkOrientedGridTransform *grid_Ras, typename itk::DisplacementFieldTransform< T, 3 >::DisplacementFieldType::Pointer gridImage_Lps)
static bool SetITKv3BSplineFromVTK(vtkObject *loggerObject, itk::Object::Pointer &warpTransformItk, itk::Object::Pointer &bulkTransformItk, vtkOrientedBSplineTransform *bsplineVtk, bool alwaysAddBulkTransform)
static itk::Matrix< T, 3, 3 > Matrix2Dto3D(itk::Matrix< T, 2, 2 > m2D)
static bool SetVTKOrientedGridTransformFrom2DITKImage(vtkObject *loggerObject, vtkOrientedGridTransform *grid_Ras, typename itk::DisplacementFieldTransform< T, 2 >::DisplacementFieldType::Pointer gridImage_Lps)
static bool SetVTKBSplineFromITKv4Generic(vtkObject *loggerObject, vtkOrientedBSplineTransform *bsplineVtk, typename itk::TransformBaseTemplate< T >::Pointer warpTransformItk)
static bool IsIdentityMatrix(vtkMatrix4x4 *matrix)
static void FlattenGeneralTransform(vtkCollection *outputTransformList, vtkAbstractTransform *inputTransform)
itk::AffineTransform< TransformRealType, 3 > AffineTransformType
Definition dtitypes.h:64
static const unsigned int VTKDimension
itk::DisplacementFieldTransform< double, 3 > DisplacementFieldTransformDoubleType
itk::TransformFileWriter TransformWriterType
DisplacementFieldTransformDoubleType::DisplacementFieldType GridImageDoubleType
itk::ThinPlateSplineKernelTransform< double, 3 > ThinPlateSplineTransformDoubleType
static const int BSPLINE_TRANSFORM_ORDER