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