go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
elxTransformBase.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 __elxTransformBase_h
16 #define __elxTransformBase_h
17 
19 #include "elxMacro.h"
20 
21 #include "elxBaseComponentSE.h"
22 #include "itkAdvancedTransform.h"
24 #include "elxComponentDatabase.h"
25 #include "elxProgressCommand.h"
26 
27 #include <fstream>
28 #include <iomanip>
29 
30 namespace elastix
31 {
32  //using namespace itk; //Not here, because a TransformBase class was added to ITK...
33 
127 template <class TElastix>
129  : public BaseComponentSE<TElastix>
130 {
131 public:
132 
136 
138  itkTypeMacro( TransformBase, BaseComponentSE );
139 
141  typedef typename Superclass::ElastixType ElastixType;
145  typedef typename ConfigurationType
146  ::CommandLineArgumentMapType CommandLineArgumentMapType;
147  typedef typename ConfigurationType
148  ::CommandLineEntryType CommandLineEntryType;
151 
153  typedef typename ElastixType::CoordRepType CoordRepType;
154  typedef typename ElastixType::FixedImageType FixedImageType;
155  typedef typename ElastixType::MovingImageType MovingImageType;
156 
161 
164 
166  itkStaticConstMacro( FixedImageDimension,
167  unsigned int, FixedImageType::ImageDimension );
168 
170  itkStaticConstMacro( MovingImageDimension,
171  unsigned int, MovingImageType::ImageDimension );
172 
174  typedef itk::Object ObjectType;
175  typedef itk::AdvancedTransform<
176  CoordRepType,
177  itkGetStaticConstMacro( FixedImageDimension ),
178  itkGetStaticConstMacro( MovingImageDimension ) > ITKBaseType;
180  itkGetStaticConstMacro( FixedImageDimension ) > CombinationTransformType;
181  typedef typename
183 
186  typedef typename ParametersType::ValueType ValueType;
187 
191 
193  typedef typename RegistrationType::ITKBaseType ITKRegistrationType;
194  typedef typename ITKRegistrationType::OptimizerType OptimizerType;
195  typedef typename OptimizerType::ScalesType ScalesType;
196 
198  virtual ITKBaseType * GetAsITKBaseType( void )
199  {
200  return dynamic_cast<ITKBaseType *>( this );
201  }
202 
204  virtual const ITKBaseType * GetAsITKBaseType( void ) const
205  {
206  return dynamic_cast<const ITKBaseType *>( this );
207  }
208 
210  {
211  return dynamic_cast<const CombinationTransformType * >( this );
212  }
214  {
215  return dynamic_cast<CombinationTransformType * >( this );
216  }
217 
221  virtual int BeforeAllBase( void );
222 
226  virtual int BeforeAllTransformix( void );
227 
231  virtual void BeforeRegistrationBase( void );
232 
236  virtual void AfterRegistrationBase( void );
237 
239  virtual const InitialTransformType * GetInitialTransform( void ) const;
240 
242  virtual void SetInitialTransform( InitialTransformType * _arg );
243 
245  virtual void SetTransformParametersFileName( const char * filename );
246 
248  itkGetStringMacro( TransformParametersFileName );
249 
251  virtual void ReadFromFile( void );
252 
254  virtual void WriteToFile( const ParametersType & param ) const;
255 
257  virtual void WriteToFile( void ) const;
258 
260  virtual void SetReadWriteTransformParameters( const bool _arg );
261 
263  virtual void ReadInitialTransformFromFile(
264  const char * transformParameterFileName );
265 
267  virtual void TransformPoints( void ) const;
268 
270  virtual void TransformPointsSomePoints( const std::string filename ) const;
271 
273  virtual void TransformPointsSomePointsVTK( const std::string filename ) const;
274 
276  virtual void TransformPointsAllPoints( void ) const;
277 
279  virtual void ComputeDeterminantOfSpatialJacobian( void ) const;
280 
282  virtual void ComputeSpatialJacobian( void ) const;
283 
287  virtual void SetFinalParameters( void );
288 
289 protected:
290 
292  TransformBase();
294  virtual ~TransformBase();
295 
302  void AutomaticScalesEstimation( ScalesType & scales ) const;
303 
308 
309 private:
310 
312  TransformBase( const Self& ); // purposely not implemented
314  void operator=( const Self& ); // purposely not implemented
315 
318 
319 }; // end class TransformBase
320 
321 
322 } // end namespace elastix
323 
324 
325 #ifndef ITK_MANUAL_INSTANTIATION
326 #include "elxTransformBase.hxx"
327 #endif
328 
329 #endif // end #ifndef __elxTransformBase_h
This class combines two transforms: an &#39;initial transform&#39; with a &#39;current transform&#39;.
The BaseComponentSE class is a base class for elastix components that provides some basic functionali...
Superclass::RegistrationPointer RegistrationPointer
ElastixType::RegistrationBaseType RegistrationType
ElastixType::ConfigurationPointer ConfigurationPointer
ObjectPointer(* PtrToCreator)(void)
ComponentDatabaseType::ComponentDescriptionType ComponentDescriptionType
ParametersType * m_TransformParametersPointer
ElastixType::Pointer ElastixPointer
ElastixType::ConfigurationType ConfigurationType
virtual int BeforeAllTransformix(void)
std::string m_TransformParametersFileName
BaseComponentSE< TElastix > Superclass
ConfigurationType::CommandLineArgumentMapType CommandLineArgumentMapType
void AutomaticScalesEstimation(ScalesType &scales) const
OptimizerType::ScalesType ScalesType
ParametersType m_FinalParameters
virtual const CombinationTransformType * GetAsCombinationTransform(void) const
ParametersType::ValueType ValueType
virtual const ITKBaseType * GetAsITKBaseType(void) const
RegistrationType * RegistrationPointer
CombinationTransformType::InitialTransformType InitialTransformType
ElastixType::FixedImageType FixedImageType
virtual void SetReadWriteTransformParameters(const bool _arg)
virtual void WriteToFile(void) const
virtual int BeforeAllBase(void)
Superclass::ElastixType ElastixType
Superclass::ParametersType ParametersType
Transform maps points, vectors and covariant vectors from an input space to an output space...
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
virtual ITKBaseType * GetAsITKBaseType(void)
ComponentDatabase ComponentDatabaseType
ITKBaseType::ParametersType ParametersType
virtual void ReadFromFile(void)
virtual const InitialTransformType * GetInitialTransform(void) const
ElastixType::MovingImageType MovingImageType
virtual void TransformPointsSomePoints(const std::string filename) const
virtual void SetInitialTransform(InitialTransformType *_arg)
virtual void TransformPointsSomePointsVTK(const std::string filename) const
Superclass::ElastixPointer ElastixPointer
ElastixType::CoordRepType CoordRepType
elx::ProgressCommand ProgressCommandType
Superclass::RegistrationType RegistrationType
virtual void AfterRegistrationBase(void)
virtual void BeforeRegistrationBase(void)
Superclass::InputPointType InputPointType
virtual void ComputeDeterminantOfSpatialJacobian(void) const
itk::AdvancedCombinationTransform< CoordRepType, itkGetStaticConstMacro(FixedImageDimension) > CombinationTransformType
virtual void SetFinalParameters(void)
virtual void SetTransformParametersFileName(const char *filename)
ITKBaseType::InputPointType InputPointType
virtual void TransformPointsAllPoints(void) const
The ComponentDatabase class is a class that stores the New() functions of all components.
This class is the elastix base class for all Transforms.
ConfigurationType::CommandLineEntryType CommandLineEntryType
Superclass::ConfigurationType ConfigurationType
virtual CombinationTransformType * GetAsCombinationTransform(void)
Superclass::OutputPointType OutputPointType
A specialized Command object for updating the progress of a filter.
ITKRegistrationType::OptimizerType OptimizerType
Superclass::ConfigurationPointer ConfigurationPointer
virtual void ComputeSpatialJacobian(void) const
ComponentDatabase::PtrToCreator PtrToCreator
virtual void ReadInitialTransformFromFile(const char *transformParameterFileName)
virtual void TransformPoints(void) const
void operator=(const Self &)
ITKBaseType::OutputPointType OutputPointType
itk::AdvancedTransform< CoordRepType, itkGetStaticConstMacro(FixedImageDimension), itkGetStaticConstMacro(MovingImageDimension) > ITKBaseType
RegistrationType::ITKBaseType ITKRegistrationType


Generated on 05-12-2013 for elastix by doxygen 1.8.5 elastix logo