ITK  5.2.0
Insight Toolkit
itkBSplineBaseTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkBSplineBaseTransform_h
19 #define itkBSplineBaseTransform_h
20 
21 #include <iostream>
22 #include "itkTransform.h"
23 #include "itkImage.h"
25 
26 namespace itk
27 {
33 template <typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
34 class ITK_TEMPLATE_EXPORT BSplineBaseTransform : public Transform<TParametersValueType, NDimensions, NDimensions>
35 {
36 public:
37  ITK_DISALLOW_COPY_AND_MOVE(BSplineBaseTransform);
38 
44 
46  itkTypeMacro(BSplineBaseTransform, Transform);
47 
49  static constexpr unsigned int SpaceDimension = NDimensions;
50 
52  static constexpr unsigned int SplineOrder = VSplineOrder;
53 
55  itkCloneMacro(Self);
56 
58  using ScalarType = typename Superclass::ScalarType;
59 
61  using FixedParametersType = typename Superclass::FixedParametersType;
62  using ParametersType = typename Superclass::ParametersType;
63 
65  using JacobianType = typename Superclass::JacobianType;
66  using JacobianPositionType = typename Superclass::JacobianPositionType;
67  using InverseJacobianPositionType = typename Superclass::InverseJacobianPositionType;
68 
70  using TransformCategoryEnum = typename Superclass::TransformCategoryEnum;
71 
73  using NumberOfParametersType = typename Superclass::NumberOfParametersType;
74 
78 
82 
84  using InputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
85  using OutputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
86 
90 
110  void
111  SetParameters(const ParametersType & parameters) override;
112 
135  void
136  SetFixedParameters(const FixedParametersType & parameters) override = 0;
138 
155  void
156  SetParametersByValue(const ParametersType & parameters) override;
157 
166  void
167  SetIdentity();
168 
170  const ParametersType &
171  GetParameters() const override;
172 
174  const FixedParametersType &
175  GetFixedParameters() const override;
176 
178  using ParametersValueType = typename ParametersType::ValueType;
182 
194  virtual void
195  SetCoefficientImages(const CoefficientImageArray & images) = 0;
196 
200  {
201  return this->m_CoefficientImages;
202  }
203 
204  using DerivativeType = typename Superclass::DerivativeType;
205 
216  void
217  UpdateTransformParameters(const DerivativeType & update, TParametersValueType factor = 1.0) override;
218 
221 
223  using SizeType = typename RegionType::SizeType;
227 
230  TransformPoint(const InputPointType & point) const override;
231 
234 
237 
240 
249  virtual void
250  TransformPoint(const InputPointType & inputPoint,
251  OutputPointType & outputPoint,
252  WeightsType & weights,
253  ParameterIndexArrayType & indices,
254  bool & inside) const = 0;
255 
257  unsigned long
259  {
260  return m_WeightsFunction->GetNumberOfWeights();
261  }
262 
265  using Superclass::TransformVector;
266  OutputVectorType
267  TransformVector(const InputVectorType &) const override
268  {
269  itkExceptionMacro(<< "Method not applicable for deformable transform.");
270  }
272 
275  OutputVnlVectorType
276  TransformVector(const InputVnlVectorType &) const override
277  {
278  itkExceptionMacro(<< "Method not applicable for deformable transform. ");
279  }
280 
283  using Superclass::TransformCovariantVector;
284  OutputCovariantVectorType
286  {
287  itkExceptionMacro(<< "Method not applicable for deformable transform. ");
288  }
290 
292  void
293  ComputeJacobianFromBSplineWeightsWithRespectToPosition(const InputPointType &,
294  WeightsType &,
295  ParameterIndexArrayType &) const;
296 
297  void
298  ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &) const override = 0;
299 
300  void
302  {
303  itkExceptionMacro(<< "ComputeJacobianWithRespectToPosition not yet implemented "
304  "for "
305  << this->GetNameOfClass());
306  }
307  using Superclass::ComputeJacobianWithRespectToPosition;
308 
310  NumberOfParametersType
311  GetNumberOfParameters() const override = 0;
312 
314  virtual NumberOfParametersType
315  GetNumberOfParametersPerDimension() const = 0;
316 
317  TransformCategoryEnum
318  GetTransformCategory() const override
319  {
320  return Self::TransformCategoryEnum::BSpline;
321  }
322 
323  unsigned int
324  GetNumberOfAffectedWeights() const;
325 
327  using PixelType = typename ImageType::PixelType;
328 
330 
333  GetNumberOfLocalParameters() const override
334  {
335  return this->GetNumberOfParameters();
336  }
337 
338 protected:
340  void
341  PrintSelf(std::ostream & os, Indent indent) const override;
342 
344  ~BSplineBaseTransform() override = default;
345 
347  itkSetObjectMacro(WeightsFunction, WeightsFunctionType);
348  itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType);
350 
352  void
353  WrapAsImages();
354 
355 protected:
357  void
358  SetFixedParametersFromTransformDomainInformation() const;
359 
361  virtual void
362  SetFixedParametersGridSizeFromTransformDomainInformation() const = 0;
363 
365  virtual void
366  SetFixedParametersGridOriginFromTransformDomainInformation() const = 0;
367 
369  virtual void
370  SetFixedParametersGridSpacingFromTransformDomainInformation() const = 0;
371 
373  virtual void
374  SetFixedParametersGridDirectionFromTransformDomainInformation() const = 0;
375 
377  virtual void
378  SetCoefficientImageInformationFromFixedParameters() = 0;
379 
381  virtual bool
382  InsideValidRegion(ContinuousIndexType &) const = 0;
383 
384  // NOTE: There is a natural duality between the
385  // two representations of of the coefficients
386  // whereby the m_InternalParametersBuffer is
387  // needed to fit into the optimization framework
388  // and the m_CoefficientImages is needed for
389  // the Jacobian computations. This implementation
390  // is an attempt to remove as much redundancy as possible
391  // and share as much information between the two
392  // instances as possible.
393  //
399 
402 
405 
406 private:
407  static CoefficientImageArray
408  ArrayOfImagePointerGeneratorHelper();
409 }; // class BSplineBaseTransform
410 } // namespace itk
411 
412 #ifndef ITK_MANUAL_INSTANTIATION
413 # include "itkBSplineBaseTransform.hxx"
414 #endif
415 
416 #endif /* itkBSplineBaseTransform_h */
itk::BSplineInterpolationWeightFunction
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Definition: itkBSplineInterpolationWeightFunction.h:47
itk::BSplineBaseTransform::ImagePointer
typename ImageType::Pointer ImagePointer
Definition: itkBSplineBaseTransform.h:180
itk::BSplineBaseTransform::TransformCovariantVector
OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
Definition: itkBSplineBaseTransform.h:285
itk::BSplineBaseTransform::m_WeightsFunction
WeightsFunctionType::Pointer m_WeightsFunction
Definition: itkBSplineBaseTransform.h:404
itk::BSplineBaseTransform::OutputVnlVectorType
vnl_vector_fixed< TParametersValueType, SpaceDimension > OutputVnlVectorType
Definition: itkBSplineBaseTransform.h:85
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::BSplineBaseTransform::NumberOfParametersType
typename Superclass::NumberOfParametersType NumberOfParametersType
Definition: itkBSplineBaseTransform.h:73
itk::BSplineBaseTransform::InverseJacobianPositionType
typename Superclass::InverseJacobianPositionType InverseJacobianPositionType
Definition: itkBSplineBaseTransform.h:67
itk::BSplineBaseTransform::IndexType
typename RegionType::IndexType IndexType
Definition: itkBSplineBaseTransform.h:222
itk::BSplineBaseTransform::OriginType
typename ImageType::PointType OriginType
Definition: itkBSplineBaseTransform.h:226
itk::BSplineBaseTransform::FixedParametersType
typename Superclass::FixedParametersType FixedParametersType
Definition: itkBSplineBaseTransform.h:61
itk::BSplineBaseTransform::GetNumberOfLocalParameters
NumberOfParametersType GetNumberOfLocalParameters() const override
Definition: itkBSplineBaseTransform.h:333
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:69
itk::BSplineBaseTransform::TransformCategoryEnum
typename Superclass::TransformCategoryEnum TransformCategoryEnum
Definition: itkBSplineBaseTransform.h:70
itk::BSplineBaseTransform::WeightsType
typename WeightsFunctionType::WeightsType WeightsType
Definition: itkBSplineBaseTransform.h:235
itk::BSplineBaseTransform::JacobianPositionType
typename Superclass::JacobianPositionType JacobianPositionType
Definition: itkBSplineBaseTransform.h:66
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::BSplineBaseTransform::DirectionType
typename ImageType::DirectionType DirectionType
Definition: itkBSplineBaseTransform.h:225
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itkBSplineInterpolationWeightFunction.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BSplineBaseTransform::InputVnlVectorType
vnl_vector_fixed< TParametersValueType, SpaceDimension > InputVnlVectorType
Definition: itkBSplineBaseTransform.h:84
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::BSplineBaseTransform::PixelType
typename ImageType::PixelType PixelType
Definition: itkBSplineBaseTransform.h:327
itk::BSplineBaseTransform::ScalarType
typename Superclass::ScalarType ScalarType
Definition: itkBSplineBaseTransform.h:58
itk::BSplineBaseTransform::MeshSizeType
SizeType MeshSizeType
Definition: itkBSplineBaseTransform.h:329
itk::BSplineBaseTransform::TransformVector
OutputVectorType TransformVector(const InputVectorType &) const override
Definition: itkBSplineBaseTransform.h:267
itk::BSplineBaseTransform::ComputeJacobianWithRespectToPosition
void ComputeJacobianWithRespectToPosition(const InputPointType &, JacobianPositionType &) const override
Definition: itkBSplineBaseTransform.h:301
itk::BSplineBaseTransform::SizeType
typename RegionType::SizeType SizeType
Definition: itkBSplineBaseTransform.h:223
itk::BSplineBaseTransform::ContinuousIndexType
typename WeightsFunctionType::ContinuousIndexType ContinuousIndexType
Definition: itkBSplineBaseTransform.h:236
itk::BSplineBaseTransform::DerivativeType
typename Superclass::DerivativeType DerivativeType
Definition: itkBSplineBaseTransform.h:204
itk::BSplineBaseTransform::SpacingType
typename ImageType::SpacingType SpacingType
Definition: itkBSplineBaseTransform.h:224
itk::FixedArray< ImagePointer, NDimensions >
itk::BSplineBaseTransform::TransformVector
OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
Definition: itkBSplineBaseTransform.h:276
itk::BSplineBaseTransform::ParametersValueType
typename ParametersType::ValueType ParametersValueType
Definition: itkBSplineBaseTransform.h:178
itk::Image::PixelType
TPixel PixelType
Definition: itkImage.h:106
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::Image::SpacingType
typename Superclass::SpacingType SpacingType
Definition: itkImage.h:154
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::BSplineBaseTransform::m_CoefficientImages
CoefficientImageArray m_CoefficientImages
Definition: itkBSplineBaseTransform.h:398
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::BSplineBaseTransform
A base class with common elements of BSplineTransform and BSplineDeformableTransform.
Definition: itkBSplineBaseTransform.h:34
itk::Array< double >
itk::BSplineBaseTransform::ParametersType
typename Superclass::ParametersType ParametersType
Definition: itkBSplineBaseTransform.h:62
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::BSplineBaseTransform::PhysicalDimensionsType
typename ImageType::SpacingType PhysicalDimensionsType
Definition: itkBSplineBaseTransform.h:326
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::BSplineBaseTransform::GetCoefficientImages
const CoefficientImageArray GetCoefficientImages() const
Definition: itkBSplineBaseTransform.h:199
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
itkTransform.h
itk::BSplineBaseTransform::JacobianType
typename Superclass::JacobianType JacobianType
Definition: itkBSplineBaseTransform.h:65
itk::BSplineBaseTransform::GetTransformCategory
TransformCategoryEnum GetTransformCategory() const override
Definition: itkBSplineBaseTransform.h:318
itk::BSplineBaseTransform::GetNumberOfWeights
unsigned long GetNumberOfWeights() const
Definition: itkBSplineBaseTransform.h:258
itk::BSplineBaseTransform::m_InternalParametersBuffer
ParametersType m_InternalParametersBuffer
Definition: itkBSplineBaseTransform.h:401