go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkAdvancedBSplineDeformableTransformBase.h
Go to the documentation of this file.
1 /*======================================================================
2 
3 This file is part of the elastix software.
4 
5 Copyright (c) University Medical Center Utrecht. All rights reserved.
6 See src/CopyrightElastix.txt or http://elastix.isi.uu.nl/legal.php for
7 details.
8 
9 This software is distributed WITHOUT ANY WARRANTY; without even
10 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11 PURPOSE. See the above copyright notices for more information.
12 
13 ======================================================================*/
14 
15 #ifndef __itkAdvancedBSplineDeformableTransformBase_h
16 #define __itkAdvancedBSplineDeformableTransformBase_h
17 
18 #include "itkAdvancedTransform.h"
19 #include "itkImage.h"
20 #include "itkImageRegion.h"
21 
22 
23 namespace itk
24 {
25 
36 template <
37  class TScalarType = double, // Data type for scalars
38  unsigned int NDimensions = 3 > // Number of dimensions
40  : public AdvancedTransform< TScalarType, NDimensions, NDimensions >
41 {
42 public:
45  typedef AdvancedTransform<
46  TScalarType, NDimensions, NDimensions > Superclass;
47  typedef SmartPointer<Self> Pointer;
48  typedef SmartPointer<const Self> ConstPointer;
49 
52 
54  itkStaticConstMacro( SpaceDimension, unsigned int, NDimensions );
55 
70 
71  typedef typename Superclass
74  typedef typename Superclass
77  typedef typename Superclass
80 
99  void SetParameters( const ParametersType & parameters );
100 
116  void SetFixedParameters( const ParametersType & parameters );
117 
134  void SetParametersByValue( const ParametersType & parameters );
135 
144  void SetIdentity( void );
145 
147  virtual const ParametersType& GetParameters( void ) const;
148 
150  virtual const ParametersType& GetFixedParameters( void ) const;
151 
153  typedef typename ParametersType::ValueType PixelType;
154  typedef Image< PixelType,
155  itkGetStaticConstMacro( SpaceDimension )> ImageType;
156  typedef typename ImageType::Pointer ImagePointer;
157 
159  virtual const ImagePointer * GetCoefficientImage( void ) const
160  { return this->m_CoefficientImage; }
161 
173  virtual void SetCoefficientImage( ImagePointer images[] );
174 
176  typedef ImageRegion< itkGetStaticConstMacro( SpaceDimension ) > RegionType;
177 
178  typedef typename RegionType::IndexType IndexType;
179  typedef typename RegionType::SizeType SizeType;
180  typedef typename ImageType::SpacingType SpacingType;
181  typedef typename ImageType::DirectionType DirectionType;
182  typedef typename ImageType::PointType OriginType;
184 
186  virtual void SetGridRegion( const RegionType& region ) = 0;
187  //itkGetMacro( GridRegion, RegionType );
188  itkGetConstMacro( GridRegion, RegionType );
189 
191  virtual void SetGridSpacing( const SpacingType & spacing );
192  //itkGetMacro( GridSpacing, SpacingType );
193  itkGetConstMacro( GridSpacing, SpacingType );
194 
196  virtual void SetGridDirection( const DirectionType & direction );
197  //itkGetMacro( GridDirection, DirectionType );
198  itkGetConstMacro( GridDirection, DirectionType );
199 
201  virtual void SetGridOrigin( const OriginType& origin );
202  //itkGetMacro( GridOrigin, OriginType );
203  itkGetConstMacro( GridOrigin, OriginType );
204 
206  typedef Array<unsigned long> ParameterIndexArrayType;
207 
211  virtual OutputVectorType TransformVector( const InputVectorType & ) const
212  {
213  itkExceptionMacro( << "Method not applicable for deformable transform." );
214  return OutputVectorType();
215  }
216 
220  virtual OutputVnlVectorType TransformVector( const InputVnlVectorType & ) const
221  {
222  itkExceptionMacro( << "Method not applicable for deformable transform. ");
223  return OutputVnlVectorType();
224  }
225 
229  virtual OutputCovariantVectorType TransformCovariantVector(
230  const InputCovariantVectorType & ) const
231  {
232  itkExceptionMacro( << "Method not applicable for deformable transform. ");
233  return OutputCovariantVectorType();
234  }
235 
237  virtual unsigned int GetNumberOfParameters( void ) const;
238 
240  virtual unsigned int GetNumberOfParametersPerDimension( void ) const;
241 
243  itkGetConstReferenceMacro( ValidRegion, RegionType );
244 
250  virtual bool IsLinear( void ) const { return false; }
251 
252  virtual unsigned int GetNumberOfAffectedWeights( void ) const = 0;
253 
254  virtual unsigned long GetNumberOfNonZeroJacobianIndices( void ) const = 0;
255 
259  typedef ContinuousIndex<ScalarType, SpaceDimension> ContinuousIndexType;
260 
261 protected:
263  virtual void PrintSelf( std::ostream &os, Indent indent ) const;
264 
267 
269  void WrapAsImages( void );
270 
272  void TransformPointToContinuousGridIndex(
273  const InputPointType & point, ContinuousIndexType & index ) const;
274 
275  virtual void ComputeNonZeroJacobianIndices(
276  NonZeroJacobianIndicesType & nonZeroJacobianIndices,
277  const RegionType & supportRegion ) const = 0;
278 
280  virtual bool InsideValidRegion( const ContinuousIndexType& index ) const;
281 
285  ImagePointer m_CoefficientImage[ NDimensions ];
286 
293 
299 
301 
303  unsigned long m_Offset;
307 
310 
313 
315  typedef typename JacobianType::ValueType JacobianPixelType;
316  typedef Image< JacobianPixelType,
317  itkGetStaticConstMacro( SpaceDimension ) > JacobianImageType;
318 
319  typename JacobianImageType::Pointer m_JacobianImage[ NDimensions ];
320 
325 
327  ImagePointer m_WrappedImage[ NDimensions ];
328 
331 
332  void UpdateGridOffsetTable( void );
333 
334 private:
335  AdvancedBSplineDeformableTransformBase(const Self&); //purposely not implemented
336  void operator=(const Self&); //purposely not implemented
337 
338 }; //class AdvancedBSplineDeformableTransformBase
339 
340 
341 } // namespace itk
342 
343 // Define instantiation macro for this template.
344 #define ITK_TEMPLATE_AdvancedBSplineDeformableTransformBase(_, EXPORT, x, y) namespace itk { \
345  _(3(class EXPORT AdvancedBSplineDeformableTransformBase< ITK_TEMPLATE_3 x >)) \
346  namespace Templates { typedef AdvancedBSplineDeformableTransformBase< ITK_TEMPLATE_3 x > \
347  AdvancedBSplineDeformableTransformBase##y; } \
348  }
349 
350 #if ITK_TEMPLATE_EXPLICIT
351 # include "Templates/itkAdvancedBSplineDeformableTransformBase+-.h"
352 #endif
353 
354 #if ITK_TEMPLATE_TXX
355 # include "itkAdvancedBSplineDeformableTransformBase.txx"
356 #endif
357 
358 #endif /* __itkAdvancedBSplineDeformableTransformBase_h */


Generated on 21-03-2014 for elastix by doxygen 1.8.1.2 elastix logo