ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkBSplineBaseTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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,
34  unsigned int NDimensions = 3,
35  unsigned int VSplineOrder = 3>
36 class ITK_TEMPLATE_EXPORT BSplineBaseTransform :
37  public Transform<TParametersValueType, NDimensions, NDimensions>
38 {
39 public:
40  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineBaseTransform);
41 
47 
49  itkTypeMacro( BSplineBaseTransform, Transform );
50 
52  static constexpr unsigned int SpaceDimension = NDimensions;
53 
55  static constexpr unsigned int SplineOrder = VSplineOrder;
56 
58  itkCloneMacro(Self);
59 
61  using ScalarType = typename Superclass::ScalarType;
62 
64  using FixedParametersType = typename Superclass::FixedParametersType;
66 
68  using JacobianType = typename Superclass::JacobianType;
69  using JacobianPositionType = typename Superclass::JacobianPositionType;
70  using InverseJacobianPositionType = typename Superclass::InverseJacobianPositionType;
71 
73  using TransformCategoryType = typename Superclass::TransformCategoryType;
74 
76  using NumberOfParametersType = typename Superclass::NumberOfParametersType;
77 
81 
85 
87  using InputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
88  using OutputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
89 
93 
113  void SetParameters( const ParametersType & parameters ) override;
114 
137  void SetFixedParameters( const FixedParametersType & parameters ) override = 0;
139 
156  void SetParametersByValue( const ParametersType & parameters ) override;
157 
166  void SetIdentity();
167 
169  const ParametersType & GetParameters() const override;
170 
172  const FixedParametersType & GetFixedParameters() const override;
173 
175  using ParametersValueType = typename ParametersType::ValueType;
179 
191  virtual void SetCoefficientImages( const CoefficientImageArray & images ) = 0;
192 
195  {
196  return this->m_CoefficientImages;
197  }
198 
199  using DerivativeType = typename Superclass::DerivativeType;
200 
211  void UpdateTransformParameters( const DerivativeType & update, TParametersValueType factor = 1.0 ) override;
212 
215 
217  using SizeType = typename RegionType::SizeType;
221 
223  OutputPointType TransformPoint( const InputPointType & point ) const override;
224 
227  Self::SpaceDimension ,
228  Self::SplineOrder >;
229 
232 
235 
244  virtual void TransformPoint( const InputPointType & inputPoint, OutputPointType & outputPoint,
245  WeightsType & weights, ParameterIndexArrayType & indices, bool & inside ) const = 0;
246 
248  unsigned long GetNumberOfWeights() const
249  {
250  return m_WeightsFunction->GetNumberOfWeights();
251  }
252 
255  using Superclass::TransformVector;
257  {
258  itkExceptionMacro( << "Method not applicable for deformable transform." );
259  }
261 
265  {
266  itkExceptionMacro( << "Method not applicable for deformable transform. " );
267  }
268 
271  using Superclass::TransformCovariantVector;
273  const InputCovariantVectorType & ) const override
274  {
275  itkExceptionMacro( << "Method not applicable for deformable transform. " );
276  }
278 
280  void ComputeJacobianFromBSplineWeightsWithRespectToPosition(
281  const InputPointType &, WeightsType &, ParameterIndexArrayType & ) const;
282 
283  void ComputeJacobianWithRespectToParameters( const InputPointType &, JacobianType & ) const override = 0;
284 
286  {
287  itkExceptionMacro( << "ComputeJacobianWithRespectToPosition not yet implemented "
288  "for " << this->GetNameOfClass() );
289  }
290  using Superclass::ComputeJacobianWithRespectToPosition;
291 
293  NumberOfParametersType GetNumberOfParameters() const override = 0;
294 
296  virtual NumberOfParametersType GetNumberOfParametersPerDimension() const = 0;
297 
299  {
300  return Self::BSpline;
301  }
302 
303  unsigned int GetNumberOfAffectedWeights() const;
304 
306  using PixelType = typename ImageType::PixelType;
307 
309 
312  {
313  return this->GetNumberOfParameters();
314  }
315 
316 protected:
318  void PrintSelf( std::ostream & os, Indent indent ) const override;
319 
321  ~BSplineBaseTransform() override = default;
322 
324  itkSetObjectMacro( WeightsFunction, WeightsFunctionType );
325  itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType );
327 
329  void WrapAsImages();
330 
331 protected:
333  void SetFixedParametersFromTransformDomainInformation() const;
334 
336  virtual void SetFixedParametersGridSizeFromTransformDomainInformation() const = 0;
337 
339  virtual void SetFixedParametersGridOriginFromTransformDomainInformation() const = 0;
340 
342  virtual void SetFixedParametersGridSpacingFromTransformDomainInformation() const = 0;
343 
345  virtual void SetFixedParametersGridDirectionFromTransformDomainInformation() const = 0;
346 
348  virtual void SetCoefficientImageInformationFromFixedParameters() =0;
349 
351  virtual bool InsideValidRegion( ContinuousIndexType & ) const = 0;
352 
353  // NOTE: There is a natural duality between the
354  // two representations of of the coefficients
355  // whereby the m_InternalParametersBuffer is
356  // needed to fit into the optimization framework
357  // and the m_CoefficientImages is needed for
358  // the Jacobian computations. This implementation
359  // is an attempt to remove as much redundancy as possible
360  // and share as much information between the two
361  // instances as possible.
362  //
368 
371 
374 
375 private:
376  static CoefficientImageArray ArrayOfImagePointerGeneratorHelper();
377 }; // class BSplineBaseTransform
378 } // namespace itk
379 
380 #ifndef ITK_MANUAL_INSTANTIATION
381 #include "itkBSplineBaseTransform.hxx"
382 #endif
383 
384 #endif /* itkBSplineBaseTransform_h */
typename Superclass::DerivativeType DerivativeType
vnl_vector_fixed< TParametersValueType, SpaceDimension > OutputVnlVectorType
typename Superclass::ParametersType ParametersType
Definition: itkTransform.h:119
TPixel PixelType
Definition: itkImage.h:95
Light weight base class for most itk classes.
typename ImageType::PixelType PixelType
unsigned long GetNumberOfWeights() const
typename ImageType::SpacingType SpacingType
CoefficientImageArray m_CoefficientImages
typename ParametersType::ValueType ParametersValueType
typename Superclass::NumberOfParametersType NumberOfParametersType
typename Superclass::ParametersType ParametersType
An image region represents a structured region of data.
typename Superclass::FixedParametersType FixedParametersType
typename ImageType::DirectionType DirectionType
typename ImageType::SpacingType PhysicalDimensionsType
TransformCategoryType GetTransformCategory() const override
OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
WeightsFunctionType::Pointer m_WeightsFunction
typename WeightsFunctionType::ContinuousIndexType ContinuousIndexType
typename ImageType::PointType OriginType
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
typename Superclass::JacobianPositionType JacobianPositionType
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
typename RegionType::IndexType IndexType
typename Superclass::SpacingType SpacingType
Definition: itkImage.h:143
A base class with common elements of BSplineTransform and BSplineDeformableTransform.
typename WeightsFunctionType::WeightsType WeightsType
typename Superclass::InverseJacobianPositionType InverseJacobianPositionType
typename RegionType::SizeType SizeType
typename Superclass::JacobianType JacobianType
A templated class holding a point in n-Dimensional image space.
const CoefficientImageArray GetCoefficientImages() const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void ComputeJacobianWithRespectToPosition(const InputPointType &, JacobianPositionType &) const override
NumberOfParametersType GetNumberOfLocalParameters() const override
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:52
typename Superclass::ScalarType ScalarType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
A templated class holding a n-Dimensional covariant vector.
vnl_vector_fixed< TParametersValueType, SpaceDimension > InputVnlVectorType
Templated n-dimensional image class.
Definition: itkImage.h:75
typename ImageType::Pointer ImagePointer
OutputVectorType TransformVector(const InputVectorType &) const override