File: itkAdvancedBSplineDeformableTransformBase.h

package info (click to toggle)
elastix 5.2.0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 42,480 kB
  • sloc: cpp: 68,403; lisp: 4,118; python: 1,013; xml: 182; sh: 177; makefile: 33
file content (427 lines) | stat: -rw-r--r-- 15,438 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
/*=========================================================================
 *
 *  Copyright UMC Utrecht and contributors
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *        http://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkAdvancedBSplineDeformableTransformBase_h
#define itkAdvancedBSplineDeformableTransformBase_h

#include "itkAdvancedTransform.h"
#include "itkImage.h"
#include "itkImageRegion.h"

namespace itk
{

/** \class AdvancedBSplineDeformableTransformBase
 * \brief Base class for deformable transform using a B-spline representation
 *
 * This class is the base for the encapsulation of a deformable transform
 * of points from one N-dimensional one space to another N-dimensional space.
 *
 * This class is not templated over the spline order, which makes the use of
 * different spline orders more convenient in subsequent code.
 *
 */
template <class TScalarType = double, // Data type for scalars
          unsigned int NDimensions = 3>
// Number of dimensions
class ITK_TEMPLATE_EXPORT AdvancedBSplineDeformableTransformBase
  : public AdvancedTransform<TScalarType, NDimensions, NDimensions>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(AdvancedBSplineDeformableTransformBase);

  /** Standard class typedefs. */
  using Self = AdvancedBSplineDeformableTransformBase;
  using Superclass = AdvancedTransform<TScalarType, NDimensions, NDimensions>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

  /** Run-time type information (and related methods). */
  itkTypeMacro(AdvancedBSplineDeformableTransformBase, AdvancedTransform);

  /** Dimension of the domain space. */
  itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions);

  /** The number of fixed parameters. For Grid size, origin, spacing, and direction. */
  static constexpr unsigned int NumberOfFixedParameters = NDimensions * (NDimensions + 3);

  /** Typedefs from Superclass. */
  using typename Superclass::ParametersType;
  using typename Superclass::FixedParametersType;
  using typename Superclass::ParametersValueType;
  using typename Superclass::NumberOfParametersType;
  using typename Superclass::DerivativeType;
  using typename Superclass::JacobianType;
  using typename Superclass::ScalarType;
  using typename Superclass::InputPointType;
  using typename Superclass::OutputPointType;
  using typename Superclass::InputVectorType;
  using typename Superclass::OutputVectorType;
  using typename Superclass::InputVnlVectorType;
  using typename Superclass::OutputVnlVectorType;
  using typename Superclass::InputCovariantVectorType;
  using typename Superclass::OutputCovariantVectorType;
  using typename Superclass::TransformCategoryEnum;

  using typename Superclass::NonZeroJacobianIndicesType;
  using typename Superclass::SpatialJacobianType;
  using typename Superclass::JacobianOfSpatialJacobianType;
  using typename Superclass::SpatialHessianType;
  using typename Superclass::JacobianOfSpatialHessianType;
  using typename Superclass::InternalMatrixType;
  using typename Superclass::MovingImageGradientType;
  using typename Superclass::MovingImageGradientValueType;

  /* Creates a `BSplineDeformableTransform` of the specified derived type and spline order. */
  template <template <class, unsigned, unsigned> class TBSplineDeformableTransform>
  static Pointer
  Create(const unsigned splineOrder)
  {
    switch (splineOrder)
    {
      case 1:
      {
        return TBSplineDeformableTransform<TScalarType, NDimensions, 1>::New();
      }
      case 2:
      {
        return TBSplineDeformableTransform<TScalarType, NDimensions, 2>::New();
      }
      case 3:
      {
        return TBSplineDeformableTransform<TScalarType, NDimensions, 3>::New();
      }
    }
    itkGenericExceptionMacro("ERROR: The provided spline order (" << splineOrder << ") is not supported.");
  }


  unsigned
  GetSplineOrder() const
  {
    return m_SplineOrder;
  }


  /** This method sets the parameters of the transform.
   * For a B-spline deformation transform, the parameters are the BSpline
   * coefficients on a sparse grid.
   *
   * The parameters are N number of N-D grid of coefficients. Each N-D grid
   * is represented as a flat array of doubles
   * (in the same configuration as an itk::Image).
   * The N arrays are then concatenated to form one parameter array.
   *
   * For efficiency, this transform does not make a copy of the parameters.
   * It only keeps a pointer to the input parameters. It assumes that the memory
   * is managed by the caller. Use SetParametersByValue to force the transform
   * to call copy the parameters.
   *
   * This method wraps each grid as itk::Image's using the user specified
   * grid region, spacing and origin.
   * NOTE: The grid region, spacing and origin must be set first.
   */
  void
  SetParameters(const ParametersType & parameters) override;

  /** This method sets the fixed parameters of the transform.
   * For a B-spline deformation transform, the parameters are the following:
   *    Grid Size, Grid Origin, and Grid Spacing
   *
   * The fixed parameters are the three times the size of the templated
   * dimensions.
   * This function has the effect of make the following calls:
   *       transform->SetGridSpacing( spacing );
   *       transform->SetGridOrigin( origin );
   *       transform->SetGridDirection( direction );
   *       transform->SetGridRegion( bsplineRegion );
   *
   * This function was added to allow the transform to work with the
   * itkTransformReader/Writer I/O filters.
   */
  void
  SetFixedParameters(const FixedParametersType & parameters) override;

  /** This method sets the parameters of the transform.
   * For a B-spline deformation transform, the parameters are the BSpline
   * coefficients on a sparse grid.
   *
   * The parameters are N number of N-D grid of coefficients. Each N-D grid
   * is represented as a flat array of doubles
   * (in the same configuration as an itk::Image).
   * The N arrays are then concatenated to form one parameter array.
   *
   * This methods makes a copy of the parameters while for
   * efficiency the SetParameters method does not.
   *
   * This method wraps each grid as itk::Image's using the user specified
   * grid region, spacing and origin.
   * NOTE: The grid region, spacing and origin must be set first.
   */
  void
  SetParametersByValue(const ParametersType & parameters) override;

  /** This method can ONLY be invoked AFTER calling SetParameters().
   *  This restriction is due to the fact that the AdvancedBSplineDeformableTransform
   *  does not copy the array of parameters internally, instead it keeps a
   *  pointer to the user-provided array of parameters. This method is also
   *  in violation of the const-correctness of the parameters since the
   *  parameter array has been passed to the transform on a 'const' basis but
   *  the values get modified when the user invokes SetIdentity().
   */
  void
  SetIdentity();

  /** Get the Transformation Parameters. */
  const ParametersType &
  GetParameters() const override;

  /** Get the Transformation Fixed Parameters. */
  const FixedParametersType &
  GetFixedParameters() const override;

  /** Parameters as SpaceDimension number of images. */
  using PixelType = typename ParametersType::ValueType;
  using ImageType = Image<PixelType, Self::SpaceDimension>;
  using ImagePointer = typename ImageType::Pointer;

  /** Get the array of coefficient images. */
  virtual const ImagePointer *
  GetCoefficientImages() const
  {
    return this->m_CoefficientImages;
  }

  /** Set the array of coefficient images.
   *
   * This is an alternative API for setting the B-spline coefficients
   * as an array of SpaceDimension images. The grid region spacing
   * and origin is taken from the first image. It is assume that
   * the buffered region of all the subsequent images are the same
   * as the first image. Note that no error checking is done.
   *
   * Warning: use either the SetParameters() or SetCoefficientImages()
   * API. Mixing the two modes may results in unexpected results.
   */
  virtual void
  SetCoefficientImages(ImagePointer images[]);

  /** Typedefs for specifying the extend to the grid. */
  using RegionType = ImageRegion<Self::SpaceDimension>;

  using IndexType = typename RegionType::IndexType;
  using SizeType = typename RegionType::SizeType;
  using SpacingType = typename ImageType::SpacingType;
  using DirectionType = typename ImageType::DirectionType;
  using OriginType = typename ImageType::PointType;
  using GridOffsetType = IndexType;

  /** This method specifies the region over which the grid resides. */
  virtual void
  SetGridRegion(const RegionType & region) = 0;

  itkGetConstMacro(GridRegion, RegionType);

  /** This method specifies the grid spacing or resolution. */
  virtual void
  SetGridSpacing(const SpacingType & spacing);

  itkGetConstMacro(GridSpacing, SpacingType);

  /** This method specifies the grid directions . */
  virtual void
  SetGridDirection(const DirectionType & direction);

  itkGetConstMacro(GridDirection, DirectionType);

  /** This method specifies the grid origin. */
  virtual void
  SetGridOrigin(const OriginType & origin);

  itkGetConstMacro(GridOrigin, OriginType);

  /** Parameter index array type. */
  using ParameterIndexArrayType = Array<unsigned long>;

  /** Method to transform a vector -
   *  not applicable for this type of transform.
   */
  OutputVectorType
  TransformVector(const InputVectorType &) const override
  {
    itkExceptionMacro("Method not applicable for deformable transform.");
  }


  /** Method to transform a vnl_vector -
   *  not applicable for this type of transform.
   */
  OutputVnlVectorType
  TransformVector(const InputVnlVectorType &) const override
  {
    itkExceptionMacro("Method not applicable for deformable transform. ");
  }


  /** Method to transform a CovariantVector -
   *  not applicable for this type of transform.
   */
  OutputCovariantVectorType
  TransformCovariantVector(const InputCovariantVectorType &) const override
  {
    itkExceptionMacro("Method not applicable for deformable transform. ");
  }


  /** Return the number of parameters that completely define the Transform. */
  NumberOfParametersType
  GetNumberOfParameters() const override;

  /** Return the number of parameters per dimension */
  virtual NumberOfParametersType
  GetNumberOfParametersPerDimension() const;

  /** Return the region of the grid wholly within the support region */
  itkGetConstReferenceMacro(ValidRegion, RegionType);

  /** Indicates that this transform is linear. That is, given two
   * points P and Q, and scalar coefficients a and b, then
   *
   *           T( a*P + b*Q ) = a * T(P) + b * T(Q)
   */
  bool
  IsLinear() const override
  {
    return false;
  }

  /** Indicates the category transform.
   *  e.g. an affine transform, or a local one, e.g. a deformation field.
   */
  TransformCategoryEnum
  GetTransformCategory() const override
  {
    return TransformCategoryEnum::BSpline;
  }


  virtual unsigned int
  GetNumberOfAffectedWeights() const = 0;

  NumberOfParametersType
  GetNumberOfNonZeroJacobianIndices() const override = 0;

  /** This typedef should be equal to the typedef used
   * in derived classes based on the weights function.
   */
  using ContinuousIndexType = ContinuousIndex<ScalarType, SpaceDimension>;

protected:
  /** Print contents of an AdvancedBSplineDeformableTransformBase. */
  void
  PrintSelf(std::ostream & os, Indent indent) const override;

  AdvancedBSplineDeformableTransformBase() = delete;
  explicit AdvancedBSplineDeformableTransformBase(const unsigned splineOrder);
  ~AdvancedBSplineDeformableTransformBase() override = default;

  /** Wrap flat array into images of coefficients. */
  void
  WrapAsImages();

  /** Convert an input point to a continuous index inside the B-spline grid. */
  ContinuousIndexType
  TransformPointToContinuousGridIndex(const InputPointType & point) const;

  void
  UpdatePointIndexConversions();

  virtual void
  ComputeNonZeroJacobianIndices(NonZeroJacobianIndicesType & nonZeroJacobianIndices,
                                const RegionType &           supportRegion) const = 0;

  /** Check if a continuous index is inside the valid region. */
  virtual bool
  InsideValidRegion(const ContinuousIndexType & index) const;

private:
  const unsigned m_SplineOrder{};

  // Private using-declarations, to avoid `-Woverloaded-virtual` warnings from GCC (GCC 11.4) or clang (macos-12).
  using Superclass::TransformVector;
  using Superclass::TransformCovariantVector;

protected:
  /** Array of images representing the B-spline coefficients
   *  in each dimension.
   */
  ImagePointer m_CoefficientImages[NDimensions]{};

  /** Variables defining the coefficient grid extend. */
  RegionType     m_GridRegion{};
  SpacingType    m_GridSpacing{ 1.0 }; // default spacing is all ones
  DirectionType  m_GridDirection{ DirectionType::GetIdentity() };
  OriginType     m_GridOrigin{};
  GridOffsetType m_GridOffsetTable{};

  DirectionType                                     m_PointToIndexMatrix{};
  SpatialJacobianType                               m_PointToIndexMatrix2{};
  DirectionType                                     m_PointToIndexMatrixTransposed{};
  SpatialJacobianType                               m_PointToIndexMatrixTransposed2{};
  FixedArray<ScalarType, NDimensions>               m_PointToIndexMatrixDiagonal{};
  FixedArray<ScalarType, NDimensions * NDimensions> m_PointToIndexMatrixDiagonalProducts{};
  DirectionType                                     m_IndexToPoint{};
  bool                                              m_PointToIndexMatrixIsDiagonal{};

  RegionType m_ValidRegion{};

  /** Variables defining the interpolation support region. */
  unsigned long       m_Offset{};
  ContinuousIndexType m_ValidRegionBegin{};
  ContinuousIndexType m_ValidRegionEnd{};

  /** Keep a pointer to the input parameters. */
  const ParametersType * m_InputParametersPointer{};

  /** Jacobian as SpaceDimension number of images. */
  using JacobianPixelType = typename JacobianType::ValueType;
  using JacobianImageType = Image<JacobianPixelType, Self::SpaceDimension>;

  typename JacobianImageType::Pointer m_JacobianImage[NDimensions]{};

  /** Keep track of last support region used in computing the Jacobian
   * for fast resetting of Jacobian to zero.
   */
  mutable IndexType m_LastJacobianIndex{};

  /** Array holding images wrapped from the flat parameters. */
  ImagePointer m_WrappedImage[NDimensions]{};

  /** Internal parameters buffer. */
  ParametersType m_InternalParametersBuffer{};

  void
  UpdateGridOffsetTable();
};

} // namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkAdvancedBSplineDeformableTransformBase.hxx"
#endif

#endif /* itkAdvancedBSplineDeformableTransformBase_h */