File: itkBSplineBaseTransform.h

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (413 lines) | stat: -rw-r--r-- 16,143 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
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  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
 *
 *         https://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 itkBSplineBaseTransform_h
#define itkBSplineBaseTransform_h

#include <iostream>
#include "itkTransform.h"
#include "itkImage.h"
#include "itkBSplineInterpolationWeightFunction.h"

namespace itk
{
/** \class BSplineBaseTransform
 * \brief A base class with common elements of BSplineTransform and BSplineDeformableTransform
 *
 * \ingroup ITKTransform
 */
template <typename TParametersValueType = double, unsigned int VDimension = 3, unsigned int VSplineOrder = 3>
class ITK_TEMPLATE_EXPORT BSplineBaseTransform : public Transform<TParametersValueType, VDimension, VDimension>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(BSplineBaseTransform);

  /** Standard class type aliases. */
  using Self = BSplineBaseTransform;
  using Superclass = Transform<TParametersValueType, VDimension, VDimension>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(BSplineBaseTransform);

  /** Dimension of the domain space. */
  static constexpr unsigned int SpaceDimension = VDimension;

  /** The BSpline order. */
  static constexpr unsigned int SplineOrder = VSplineOrder;

  /** implement type-specific clone method*/
  itkCloneMacro(Self);

  /** Standard scalar type for this class. */
  using typename Superclass::ScalarType;

  /** Standard parameters container. */
  using typename Superclass::FixedParametersType;
  using typename Superclass::ParametersType;

  /** Standard Jacobian container. */
  using typename Superclass::JacobianType;
  using typename Superclass::JacobianPositionType;
  using typename Superclass::InverseJacobianPositionType;

  /** Transform category type. */
  using typename Superclass::TransformCategoryEnum;

  /** The number of parameters defining this transform. */
  using typename Superclass::NumberOfParametersType;

  /** Standard vector type for this class. */
  using InputVectorType = Vector<TParametersValueType, Self::SpaceDimension>;
  using OutputVectorType = Vector<TParametersValueType, Self::SpaceDimension>;

  /** Standard covariant vector type for this class. */
  using InputCovariantVectorType = CovariantVector<TParametersValueType, Self::SpaceDimension>;
  using OutputCovariantVectorType = CovariantVector<TParametersValueType, Self::SpaceDimension>;

  /** Standard vnl_vector type for this class. */
  using InputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
  using OutputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;

  /** Standard coordinate point type for this class. */
  using InputPointType = Point<TParametersValueType, Self::SpaceDimension>;
  using OutputPointType = Point<TParametersValueType, Self::SpaceDimension>;

  /** This method sets the parameters of the transform.
   * For a BSpline 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 scalars (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
   * fixed parameters.
   * NOTE: The transform domain must be set first.
   *
   */
  void
  SetParameters(const ParametersType & parameters) override;

  /** This method sets the fixed parameters of the transform.
   * For a BSpline deformation transform, the fixed parameters are the
   * following: grid size, grid origin, grid spacing, and grid direction.
   * However, all of these are set via the much more intuitive
   * SetTransformDomainXXX() functions
   *
   * The fixed parameters are the three times the size of the templated
   * dimensions.  This function has the effect of make the following non-
   * existing functional calls:
   *   transform->SetGridSpacing( spacing );
   *   transform->SetGridOrigin( origin );
   *   transform->SetGridDirection( direction );
   *   transform->SetGridRegion( bsplineRegion );
   *
   * With recent updates to this transform, however, all these parameters
   * are set indirectly by setting the transform domain parameters unless
   * the user sets them with SetFixedParameters().
   *
   * This function was added to allow the transform to work with the
   * itkTransformReader/Writer I/O filters.
   *
   */
  void
  SetFixedParameters(const FixedParametersType & parameters) override = 0;

  /** This method sets the parameters of the transform.
   * For a BSpline 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
   * fixed parameters.
   * NOTE: The fixed parameters 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 BSplineBaseTransform
   *  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 ParametersValueType = typename ParametersType::ValueType;
  using ImageType = Image<ParametersValueType, Self::SpaceDimension>;
  using ImagePointer = typename ImageType::Pointer;
  using CoefficientImageArray = FixedArray<ImagePointer, VDimension>;

  /** Set the array of coefficient images.
   *
   * This is an alternative API for setting the BSpline coefficients
   * as an array of SpaceDimension images. The fixed parameters are
   * taken from the first image. It is assumed 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(const CoefficientImageArray & images) = 0;

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

  using typename Superclass::DerivativeType;

  /** Update the transform's parameters by the adding values in \c update
   * to current parameter values.
   * We assume \c update is of the same length as Parameters. Throw
   * exception otherwise.
   * \c factor is a scalar multiplier for each value in update.
   * SetParameters is called at the end of this method, to allow transforms
   * to perform any required operations on the update parameters, typically
   * a conversion to member variables for use in TransformPoint.
   * Derived classes should override to provide specialized behavior.
   */
  void
  UpdateTransformParameters(const DerivativeType & update, TParametersValueType factor = 1.0) override;

  /** Typedefs for specifying the extent of 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;

  /** Transform points by a BSpline deformable transformation. */
  OutputPointType
  TransformPoint(const InputPointType & point) const override;

  /** Interpolation weights function type. */
  using WeightsFunctionType = BSplineInterpolationWeightFunction<ScalarType, Self::SpaceDimension, Self::SplineOrder>;

  using WeightsType = typename WeightsFunctionType::WeightsType;
  using ContinuousIndexType = typename WeightsFunctionType::ContinuousIndexType;

  /** Number of weights. */
  static constexpr unsigned int NumberOfWeights{ WeightsFunctionType::NumberOfWeights };

  /** Parameter index array type. */
  using ParameterIndexArrayType = FixedArray<unsigned long, NumberOfWeights>;

  /**
   * Transform points by a BSpline deformable transformation.
   * On return, weights contains the interpolation weights used to compute the
   * deformation and indices of the x (zeroth) dimension coefficient parameters
   * in the support region used to compute the deformation.
   * Parameter indices for the i-th dimension can be obtained by adding
   * ( i * this->GetNumberOfParametersPerDimension() ) to the indices array.
   */
  virtual void
  TransformPoint(const InputPointType &    inputPoint,
                 OutputPointType &         outputPoint,
                 WeightsType &             weights,
                 ParameterIndexArrayType & indices,
                 bool &                    inside) const = 0;

#if !defined(ITK_LEGACY_REMOVE)
  /** Get number of weights. */
  itkLegacyMacro(unsigned long GetNumberOfWeights() const) { return m_WeightsFunction->GetNumberOfWeights(); }
#endif

  /** Method to transform a vector -
   *  not applicable for this type of transform. */
  using Superclass::TransformVector;
  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 */
  using Superclass::TransformCovariantVector;
  OutputCovariantVectorType
  TransformCovariantVector(const InputCovariantVectorType &) const override
  {
    itkExceptionMacro("Method not applicable for deformable transform. ");
  }

  /** Get Jacobian at a point. A very specialized function just for BSplines */
  void
  ComputeJacobianFromBSplineWeightsWithRespectToPosition(const InputPointType &,
                                                         WeightsType &,
                                                         ParameterIndexArrayType &) const;

  void
  ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &) const override = 0;

  void
  ComputeJacobianWithRespectToPosition(const InputPointType &, JacobianPositionType &) const override
  {
    itkExceptionMacro("ComputeJacobianWithRespectToPosition not yet implemented "
                      "for "
                      << this->GetNameOfClass());
  }
  using Superclass::ComputeJacobianWithRespectToPosition;

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

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

  TransformCategoryEnum
  GetTransformCategory() const override
  {
    return Self::TransformCategoryEnum::BSpline;
  }

  unsigned int
  GetNumberOfAffectedWeights() const;

  using PhysicalDimensionsType = typename ImageType::SpacingType;
  using PixelType = typename ImageType::PixelType;

  using MeshSizeType = SizeType;

  /** Return the number of local parameters */
  NumberOfParametersType
  GetNumberOfLocalParameters() const override
  {
    return this->GetNumberOfParameters();
  }

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

  BSplineBaseTransform() = default;
  ~BSplineBaseTransform() override = default;

  /** Get/Set to allow subclasses to access and manipulate the weights function. */
  itkSetObjectMacro(WeightsFunction, WeightsFunctionType);
  itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType);

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

protected:
  /** Construct control point grid from transform domain information */
  void
  SetFixedParametersFromTransformDomainInformation() const;

  /** Construct control point grid size from transform domain information */
  virtual void
  SetFixedParametersGridSizeFromTransformDomainInformation() const = 0;

  /** Construct control point grid origin from transform domain information */
  virtual void
  SetFixedParametersGridOriginFromTransformDomainInformation() const = 0;

  /** Construct control point grid spacing from transform domain information */
  virtual void
  SetFixedParametersGridSpacingFromTransformDomainInformation() const = 0;

  /** Construct control point grid direction from transform domain information */
  virtual void
  SetFixedParametersGridDirectionFromTransformDomainInformation() const = 0;

  /** Construct control point grid size from transform domain information */
  virtual void
  SetCoefficientImageInformationFromFixedParameters() = 0;

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

  // NOTE:  There is a natural duality between the
  //       two representations of of the coefficients
  //       whereby the m_InternalParametersBuffer is
  //       needed to fit into the optimization framework
  //       and the m_CoefficientImages is needed for
  //       the Jacobian computations.  This implementation
  //       is an attempt to remove as much redundancy as possible
  //       and share as much information between the two
  //       instances as possible.
  //
  /** Array of images representing the B-spline coefficients
   *  in each dimension wrapped from the flat parameters in
   *  m_InternalParametersBuffer
   */
  CoefficientImageArray m_CoefficientImages{ Self::ArrayOfImagePointerGeneratorHelper() };

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

  /** Pointer to function used to compute Bspline interpolation weights. */
  typename WeightsFunctionType::Pointer m_WeightsFunction{ WeightsFunctionType::New() };

private:
  static CoefficientImageArray
  ArrayOfImagePointerGeneratorHelper();
}; // class BSplineBaseTransform
} // namespace itk

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

#endif /* itkBSplineBaseTransform_h */