ITK  4.13.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:
45 
47  itkTypeMacro( BSplineBaseTransform, Transform );
48 
50  itkStaticConstMacro( SpaceDimension, unsigned int, NDimensions );
51 
53  itkStaticConstMacro( SplineOrder, unsigned int, VSplineOrder );
54 
56  itkCloneMacro(Self);
57 
59  typedef typename Superclass::ScalarType ScalarType;
60 
62  typedef typename Superclass::FixedParametersType FixedParametersType;
63  typedef typename Superclass::ParametersType ParametersType;
64 
66  typedef typename Superclass::JacobianType JacobianType;
67 
69  typedef typename Superclass::TransformCategoryType TransformCategoryType;
70 
72  typedef typename Superclass::NumberOfParametersType NumberOfParametersType;
73 
78 
83 
85  typedef vnl_vector_fixed<TParametersValueType, SpaceDimension> InputVnlVectorType;
86  typedef vnl_vector_fixed<TParametersValueType, SpaceDimension> OutputVnlVectorType;
87 
92 
112  void SetParameters( const ParametersType & parameters ) ITK_OVERRIDE;
113 
136  virtual void SetFixedParameters( const FixedParametersType & parameters ) ITK_OVERRIDE = 0;
138 
155  void SetParametersByValue( const ParametersType & parameters ) ITK_OVERRIDE;
156 
165  void SetIdentity();
166 
168  virtual const ParametersType & GetParameters() const ITK_OVERRIDE;
169 
171  virtual const FixedParametersType & GetFixedParameters() const ITK_OVERRIDE;
172 
174  typedef typename ParametersType::ValueType ParametersValueType;
175  typedef Image<ParametersValueType, itkGetStaticConstMacro( SpaceDimension )> ImageType;
176  typedef typename ImageType::Pointer ImagePointer;
177  typedef FixedArray<ImagePointer, NDimensions> CoefficientImageArray;
178 
190  virtual void SetCoefficientImages( const CoefficientImageArray & images ) = 0;
191 
193  const CoefficientImageArray GetCoefficientImages() const
194  {
195  return this->m_CoefficientImages;
196  }
197 
198  typedef typename Superclass::DerivativeType DerivativeType;
199 
210  virtual void UpdateTransformParameters( const DerivativeType & update, TParametersValueType factor = 1.0 ) ITK_OVERRIDE;
211 
213  typedef ImageRegion<itkGetStaticConstMacro( SpaceDimension )> RegionType;
214 
215  typedef typename RegionType::IndexType IndexType;
216  typedef typename RegionType::SizeType SizeType;
217  typedef typename ImageType::SpacingType SpacingType;
219  typedef typename ImageType::PointType OriginType;
220 
222  OutputPointType TransformPoint( const InputPointType & point ) const ITK_OVERRIDE;
223 
226  itkGetStaticConstMacro( SpaceDimension ),
227  itkGetStaticConstMacro( SplineOrder )> WeightsFunctionType;
228 
229  typedef typename WeightsFunctionType::WeightsType WeightsType;
230  typedef typename WeightsFunctionType::ContinuousIndexType ContinuousIndexType;
231 
233  typedef Array<unsigned long> ParameterIndexArrayType;
234 
243  virtual void TransformPoint( const InputPointType & inputPoint, OutputPointType & outputPoint,
244  WeightsType & weights, ParameterIndexArrayType & indices, bool & inside ) const = 0;
245 
247  unsigned long GetNumberOfWeights() const
248  {
249  return m_WeightsFunction->GetNumberOfWeights();
250  }
251 
254  using Superclass::TransformVector;
255  virtual OutputVectorType TransformVector( const InputVectorType & ) const ITK_OVERRIDE
256  {
257  itkExceptionMacro( << "Method not applicable for deformable transform." );
258  }
260 
263  virtual OutputVnlVectorType TransformVector( const InputVnlVectorType & ) const ITK_OVERRIDE
264  {
265  itkExceptionMacro( << "Method not applicable for deformable transform. " );
266  }
267 
270  using Superclass::TransformCovariantVector;
272  const InputCovariantVectorType & ) const ITK_OVERRIDE
273  {
274  itkExceptionMacro( << "Method not applicable for deformable transform. " );
275  }
277 
279  void ComputeJacobianFromBSplineWeightsWithRespectToPosition(
280  const InputPointType &, WeightsType &, ParameterIndexArrayType & ) const;
281 
282  virtual void ComputeJacobianWithRespectToParameters( const InputPointType &, JacobianType & ) const ITK_OVERRIDE = 0;
283 
284  virtual void ComputeJacobianWithRespectToPosition( const InputPointType &, JacobianType & ) const ITK_OVERRIDE
285  {
286  itkExceptionMacro( << "ComputeJacobianWithRespectToPosition not yet implemented "
287  "for " << this->GetNameOfClass() );
288  }
289 
291  virtual NumberOfParametersType GetNumberOfParameters() const ITK_OVERRIDE = 0;
292 
294  virtual NumberOfParametersType GetNumberOfParametersPerDimension() const = 0;
295 
296  virtual TransformCategoryType GetTransformCategory() const ITK_OVERRIDE
297  {
298  return Self::BSpline;
299  }
300 
301  unsigned int GetNumberOfAffectedWeights() const;
302 
304  typedef typename ImageType::PixelType PixelType;
305 
307 
310  {
311  return this->GetNumberOfParameters();
312  }
313 
314 protected:
316  void PrintSelf( std::ostream & os, Indent indent ) const ITK_OVERRIDE;
317 
319  virtual ~BSplineBaseTransform() ITK_OVERRIDE;
320 
322  itkSetObjectMacro( WeightsFunction, WeightsFunctionType );
323  itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType );
325 
327  void WrapAsImages();
328 
329 protected:
331  void SetFixedParametersFromTransformDomainInformation() const;
332 
334  virtual void SetFixedParametersGridSizeFromTransformDomainInformation() const = 0;
335 
337  virtual void SetFixedParametersGridOriginFromTransformDomainInformation() const = 0;
338 
340  virtual void SetFixedParametersGridSpacingFromTransformDomainInformation() const = 0;
341 
343  virtual void SetFixedParametersGridDirectionFromTransformDomainInformation() const = 0;
344 
346  virtual void SetCoefficientImageInformationFromFixedParameters() =0;
347 
349  virtual bool InsideValidRegion( ContinuousIndexType & ) const = 0;
350 
351  // NOTE: There is a natural duality between the
352  // two representations of of the coefficients
353  // whereby the m_InternalParametersBuffer is
354  // needed to fit into the optimization framework
355  // and the m_CoefficientImages is needed for
356  // the Jacobian computations. This implementation
357  // is an attempt to remove as much redundancy as possible
358  // and share as much information between the two
359  // instances as possible.
360  //
365  CoefficientImageArray m_CoefficientImages;
366 
368  ParametersType m_InternalParametersBuffer;
369 
371  typename WeightsFunctionType::Pointer m_WeightsFunction;
372 
373 private:
374  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineBaseTransform);
375 
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 */
CovariantVector< TParametersValueType, itkGetStaticConstMacro(SpaceDimension)> InputCovariantVectorType
virtual NumberOfParametersType GetNumberOfLocalParameters() const override
Light weight base class for most itk classes.
Superclass::ScalarType ScalarType
Vector< TParametersValueType, itkGetStaticConstMacro(SpaceDimension)> InputVectorType
virtual OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
virtual OutputVectorType TransformVector(const InputVectorType &) const override
An image region represents a structured region of data.
ImageType::SpacingType PhysicalDimensionsType
Superclass::NumberOfParametersType NumberOfParametersType
virtual OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
TPixel PixelType
Definition: itkImage.h:89
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
ParametersType::ValueType ParametersValueType
vnl_vector_fixed< TParametersValueType, SpaceDimension > OutputVnlVectorType
vnl_vector_fixed< TParametersValueType, SpaceDimension > InputVnlVectorType
CovariantVector< TParametersValueType, itkGetStaticConstMacro(SpaceDimension)> OutputCovariantVectorType
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
Point< TParametersValueType, itkGetStaticConstMacro(SpaceDimension)> OutputPointType
Transform< TParametersValueType, NDimensions, NDimensions > Superclass
A base class with common elements of BSplineTransform and BSplineDeformableTransform.
SmartPointer< const Self > ConstPointer
Superclass::FixedParametersType FixedParametersType
Superclass::DerivativeType DerivativeType
A templated class holding a point in n-Dimensional image space.
Point< TParametersValueType, itkGetStaticConstMacro(SpaceDimension)> InputPointType
Superclass::TransformCategoryType TransformCategoryType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Vector< TParametersValueType, itkGetStaticConstMacro(SpaceDimension)> OutputVectorType
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:52
Returns the weights over the support region used for B-spline interpolation/reconstruction.
A templated class holding a n-Dimensional covariant vector.
Superclass::ParametersType ParametersType
Templated n-dimensional image class.
Definition: itkImage.h:75
Superclass::JacobianType JacobianType