ITK  5.4.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  * https://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 VDimension = 3, unsigned int VSplineOrder = 3>
34 class ITK_TEMPLATE_EXPORT BSplineBaseTransform : public Transform<TParametersValueType, VDimension, VDimension>
35 {
36 public:
37  ITK_DISALLOW_COPY_AND_MOVE(BSplineBaseTransform);
38 
44 
46  itkOverrideGetNameOfClassMacro(BSplineBaseTransform);
47 
49  static constexpr unsigned int SpaceDimension = VDimension;
50 
52  static constexpr unsigned int SplineOrder = VSplineOrder;
53 
55  itkCloneMacro(Self);
56 
58  using typename Superclass::ScalarType;
59 
61  using typename Superclass::FixedParametersType;
62  using typename Superclass::ParametersType;
63 
65  using typename Superclass::JacobianType;
66  using typename Superclass::JacobianPositionType;
67  using typename Superclass::InverseJacobianPositionType;
68 
70  using typename Superclass::TransformCategoryEnum;
71 
73  using 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;
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 
182 
194  virtual void
195  SetCoefficientImages(const CoefficientImageArray & images) = 0;
196 
200  {
201  return this->m_CoefficientImages;
202  }
203 
204  using 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 
239  static constexpr unsigned int NumberOfWeights{ WeightsFunctionType::NumberOfWeights };
240 
243 
252  virtual void
253  TransformPoint(const InputPointType & inputPoint,
254  OutputPointType & outputPoint,
255  WeightsType & weights,
256  ParameterIndexArrayType & indices,
257  bool & inside) const = 0;
258 
259 #if !defined(ITK_LEGACY_REMOVE)
260 
261  itkLegacyMacro(unsigned long GetNumberOfWeights() const) { return m_WeightsFunction->GetNumberOfWeights(); }
262 #endif
263 
266  using Superclass::TransformVector;
267  OutputVectorType
268  TransformVector(const InputVectorType &) const override
269  {
270  itkExceptionMacro("Method not applicable for deformable transform.");
271  }
276  OutputVnlVectorType
277  TransformVector(const InputVnlVectorType &) const override
278  {
279  itkExceptionMacro("Method not applicable for deformable transform. ");
280  }
281 
284  using Superclass::TransformCovariantVector;
285  OutputCovariantVectorType
287  {
288  itkExceptionMacro("Method not applicable for deformable transform. ");
289  }
293  void
294  ComputeJacobianFromBSplineWeightsWithRespectToPosition(const InputPointType &,
295  WeightsType &,
296  ParameterIndexArrayType &) const;
297 
298  void
299  ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &) const override = 0;
300 
301  void
303  {
304  itkExceptionMacro("ComputeJacobianWithRespectToPosition not yet implemented "
305  "for "
306  << this->GetNameOfClass());
307  }
308  using Superclass::ComputeJacobianWithRespectToPosition;
309 
311  NumberOfParametersType
312  GetNumberOfParameters() const override = 0;
313 
315  virtual NumberOfParametersType
316  GetNumberOfParametersPerDimension() const = 0;
317 
318  TransformCategoryEnum
319  GetTransformCategory() const override
320  {
321  return Self::TransformCategoryEnum::BSpline;
322  }
323 
324  unsigned int
325  GetNumberOfAffectedWeights() const;
326 
328  using PixelType = typename ImageType::PixelType;
329 
331 
334  GetNumberOfLocalParameters() const override
335  {
336  return this->GetNumberOfParameters();
337  }
338 
339 protected:
341  void
342  PrintSelf(std::ostream & os, Indent indent) const override;
343 
344  BSplineBaseTransform() = default;
345  ~BSplineBaseTransform() override = default;
346 
348  itkSetObjectMacro(WeightsFunction, WeightsFunctionType);
349  itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType);
353  void
354  WrapAsImages();
355 
356 protected:
358  void
359  SetFixedParametersFromTransformDomainInformation() const;
360 
362  virtual void
363  SetFixedParametersGridSizeFromTransformDomainInformation() const = 0;
364 
366  virtual void
367  SetFixedParametersGridOriginFromTransformDomainInformation() const = 0;
368 
370  virtual void
371  SetFixedParametersGridSpacingFromTransformDomainInformation() const = 0;
372 
374  virtual void
375  SetFixedParametersGridDirectionFromTransformDomainInformation() const = 0;
376 
378  virtual void
379  SetCoefficientImageInformationFromFixedParameters() = 0;
380 
382  virtual bool
383  InsideValidRegion(ContinuousIndexType &) const = 0;
384 
385  // NOTE: There is a natural duality between the
386  // two representations of of the coefficients
387  // whereby the m_InternalParametersBuffer is
388  // needed to fit into the optimization framework
389  // and the m_CoefficientImages is needed for
390  // the Jacobian computations. This implementation
391  // is an attempt to remove as much redundancy as possible
392  // and share as much information between the two
393  // instances as possible.
394  //
399  CoefficientImageArray m_CoefficientImages{ Self::ArrayOfImagePointerGeneratorHelper() };
400 
402  ParametersType m_InternalParametersBuffer{};
403 
406 
407 private:
408  static CoefficientImageArray
409  ArrayOfImagePointerGeneratorHelper();
410 }; // class BSplineBaseTransform
411 } // namespace itk
412 
413 #ifndef ITK_MANUAL_INSTANTIATION
414 # include "itkBSplineBaseTransform.hxx"
415 #endif
416 
417 #endif /* itkBSplineBaseTransform_h */
itk::BSplineBaseTransform::GetTransformCategory
TransformCategoryEnum GetTransformCategory() const override
Definition: itkBSplineBaseTransform.h:319
itk::BSplineInterpolationWeightFunction
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Definition: itkBSplineInterpolationWeightFunction.h:48
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::OptimizerParameters
Class to hold and manage different parameter types used during optimization.
Definition: itkOptimizerParameters.h:36
itk::BSplineBaseTransform::MeshSizeType
SizeType MeshSizeType
Definition: itkBSplineBaseTransform.h:330
GetNameOfClass
const char * GetNameOfClass() const override
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::BSplineBaseTransform::PhysicalDimensionsType
typename ImageType::SpacingType PhysicalDimensionsType
Definition: itkBSplineBaseTransform.h:327
itk::BSplineBaseTransform::SizeType
typename RegionType::SizeType SizeType
Definition: itkBSplineBaseTransform.h:223
itk::BSplineBaseTransform::WeightsType
typename WeightsFunctionType::WeightsType WeightsType
Definition: itkBSplineBaseTransform.h:235
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:80
itk::BSplineBaseTransform::ComputeJacobianWithRespectToPosition
void ComputeJacobianWithRespectToPosition(const InputPointType &, JacobianPositionType &) const override
Definition: itkBSplineBaseTransform.h:302
itk::BSplineBaseTransform::DirectionType
typename ImageType::DirectionType DirectionType
Definition: itkBSplineBaseTransform.h:225
itk::BSplineBaseTransform::SpacingType
typename ImageType::SpacingType SpacingType
Definition: itkBSplineBaseTransform.h:224
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
itk::BSplineBaseTransform::ParametersValueType
typename ParametersType::ValueType ParametersValueType
Definition: itkBSplineBaseTransform.h:178
itkImage.h
itk::OptimizerParameters::ValueType
TParametersValueType ValueType
Definition: itkOptimizerParameters.h:40
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itkBSplineInterpolationWeightFunction.h
itk::BSplineBaseTransform::TransformVector
OutputVectorType TransformVector(const InputVectorType &) const override
Definition: itkBSplineBaseTransform.h:268
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BSplineBaseTransform::IndexType
typename RegionType::IndexType IndexType
Definition: itkBSplineBaseTransform.h:222
itk::BSplineBaseTransform::OriginType
typename ImageType::PointType OriginType
Definition: itkBSplineBaseTransform.h:226
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::BSplineInterpolationWeightFunction::WeightsType
typename Superclass::OutputType WeightsType
Definition: itkBSplineInterpolationWeightFunction.h:76
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::FixedArray< ImagePointer, VDimension >
itk::BSplineBaseTransform::PixelType
typename ImageType::PixelType PixelType
Definition: itkBSplineBaseTransform.h:328
itk::BSplineBaseTransform::GetNumberOfLocalParameters
NumberOfParametersType GetNumberOfLocalParameters() const override
Definition: itkBSplineBaseTransform.h:334
itk::Image::PixelType
TPixel PixelType
Definition: itkImage.h:108
itk::BSplineBaseTransform::OutputVnlVectorType
vnl_vector_fixed< TParametersValueType, SpaceDimension > OutputVnlVectorType
Definition: itkBSplineBaseTransform.h:85
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::BSplineBaseTransform::ImagePointer
typename ImageType::Pointer ImagePointer
Definition: itkBSplineBaseTransform.h:180
itk::Transform< TParametersValueType, VDimension, VDimension >::JacobianPositionType
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
Definition: itkTransform.h:145
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::BSplineBaseTransform::TransformVector
OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
Definition: itkBSplineBaseTransform.h:277
itk::BSplineBaseTransform
A base class with common elements of BSplineTransform and BSplineDeformableTransform.
Definition: itkBSplineBaseTransform.h:34
itk::BSplineBaseTransform::GetCoefficientImages
const CoefficientImageArray GetCoefficientImages() const
Definition: itkBSplineBaseTransform.h:199
itk::BSplineBaseTransform::TransformCovariantVector
OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
Definition: itkBSplineBaseTransform.h:286
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::TransformBaseTemplate::NumberOfParametersType
IdentifierType NumberOfParametersType
Definition: itkTransformBase.h:91
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
New
static Pointer New()
itkTransform.h
itk::BSplineBaseTransform::ContinuousIndexType
typename WeightsFunctionType::ContinuousIndexType ContinuousIndexType
Definition: itkBSplineBaseTransform.h:236
itk::BSplineBaseTransform::InputVnlVectorType
vnl_vector_fixed< TParametersValueType, SpaceDimension > InputVnlVectorType
Definition: itkBSplineBaseTransform.h:84
itk::images
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar images
Definition: itkVectorGradientMagnitudeImageFilter.h:107