go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkTransformToSpatialJacobianSource.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright UMC Utrecht and contributors
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 __itkTransformToSpatialJacobianSource_h
19 #define __itkTransformToSpatialJacobianSource_h
20 
21 #include "itkAdvancedTransform.h"
22 #include "itkImageSource.h"
23 
24 namespace itk
25 {
26 
65 template< class TOutputImage,
66 class TTransformPrecisionType = double >
68  public ImageSource< TOutputImage >
69 {
70 public:
71 
74  typedef ImageSource< TOutputImage > Superclass;
76  typedef SmartPointer< const Self > ConstPointer;
77 
78  typedef TOutputImage OutputImageType;
79  typedef typename OutputImageType::Pointer OutputImagePointer;
80  typedef typename OutputImageType::ConstPointer OutputImageConstPointer;
81  typedef typename OutputImageType::RegionType OutputImageRegionType;
82 
84  itkNewMacro( Self );
85 
87  itkTypeMacro( TransformToSpatialJacobianSource, ImageSource );
88 
90  itkStaticConstMacro( ImageDimension, unsigned int,
91  TOutputImage::ImageDimension );
92 
94  typedef AdvancedTransform< TTransformPrecisionType,
95  itkGetStaticConstMacro( ImageDimension ),
96  itkGetStaticConstMacro( ImageDimension ) > TransformType;
99 
101  typedef typename OutputImageType::PixelType PixelType;
102  //typedef typename PixelType::ValueType PixelValueType;
103  typedef typename OutputImageType::RegionType RegionType;
104  typedef typename RegionType::SizeType SizeType;
105  typedef typename OutputImageType::IndexType IndexType;
106  typedef typename OutputImageType::PointType PointType;
107  typedef typename OutputImageType::SpacingType SpacingType;
108  typedef typename OutputImageType::PointType OriginType;
109  typedef typename OutputImageType::DirectionType DirectionType;
110 
112  typedef ImageBase< itkGetStaticConstMacro( ImageDimension ) > ImageBaseType;
113 
121  itkSetConstObjectMacro( Transform, TransformType );
122 
124  itkGetConstObjectMacro( Transform, TransformType );
125 
127  virtual void SetOutputSize( const SizeType & size );
128 
130  virtual const SizeType & GetOutputSize();
131 
134  virtual void SetOutputIndex( const IndexType & index );
135 
137  virtual const IndexType & GetOutputIndex();
138 
140  itkSetMacro( OutputRegion, OutputImageRegionType );
141 
143  itkGetConstReferenceMacro( OutputRegion, OutputImageRegionType );
144 
146  itkSetMacro( OutputSpacing, SpacingType );
147  virtual void SetOutputSpacing( const double * values );
148 
150  itkGetConstReferenceMacro( OutputSpacing, SpacingType );
151 
153  itkSetMacro( OutputOrigin, OriginType );
154  virtual void SetOutputOrigin( const double * values );
155 
157  itkGetConstReferenceMacro( OutputOrigin, OriginType );
158 
160  itkSetMacro( OutputDirection, DirectionType );
161  itkGetConstReferenceMacro( OutputDirection, DirectionType );
162 
165 
167  virtual void GenerateOutputInformation( void );
168 
171  virtual void BeforeThreadedGenerateData( void );
172 
174  unsigned long GetMTime( void ) const;
175 
176 protected:
177 
180 
181  void PrintSelf( std::ostream & os, Indent indent ) const;
182 
187  const OutputImageRegionType & outputRegionForThread,
188  ThreadIdType threadId );
189 
194  const OutputImageRegionType & outputRegionForThread,
195  ThreadIdType threadId );
196 
200  void LinearGenerateData( void );
201 
202 private:
203 
204  TransformToSpatialJacobianSource( const Self & ); // purposely not implemented
205  void operator=( const Self & ); // purposely not implemented
206 
208  RegionType m_OutputRegion; // region of the output image
209  TransformPointerType m_Transform; // Coordinate transform to use
210  SpacingType m_OutputSpacing; // output image spacing
211  OriginType m_OutputOrigin; // output image origin
212  DirectionType m_OutputDirection; // output image direction cosines
213 
214 };
215 
216 } // end namespace itk
217 
218 #ifndef ITK_MANUAL_INSTANTIATION
219 #include "itkTransformToSpatialJacobianSource.hxx"
220 #endif
221 
222 #endif // end #ifndef __itkTransformToSpatialJacobianSource_h
itk::TransformToSpatialJacobianSource::Self
TransformToSpatialJacobianSource Self
Definition: itkTransformToSpatialJacobianSource.h:73
itk::TransformToSpatialJacobianSource::m_Transform
TransformPointerType m_Transform
Definition: itkTransformToSpatialJacobianSource.h:209
itk::TransformToSpatialJacobianSource::TransformToSpatialJacobianSource
TransformToSpatialJacobianSource()
itk::TransformToSpatialJacobianSource::Pointer
SmartPointer< Self > Pointer
Definition: itkTransformToSpatialJacobianSource.h:75
itk::TransformToSpatialJacobianSource::itkStaticConstMacro
itkStaticConstMacro(ImageDimension, unsigned int, TOutputImage::ImageDimension)
itk::TransformToSpatialJacobianSource::Superclass
ImageSource< TOutputImage > Superclass
Definition: itkTransformToSpatialJacobianSource.h:74
itk::AdvancedTransform::SpatialJacobianType
Matrix< ScalarType, OutputSpaceDimension, InputSpaceDimension > SpatialJacobianType
Definition: itkAdvancedTransform.h:143
itkAdvancedTransform.h
SmartPointer< Self >
itk::TransformToSpatialJacobianSource::TransformPointerType
TransformType::ConstPointer TransformPointerType
Definition: itkTransformToSpatialJacobianSource.h:97
itk::TransformToSpatialJacobianSource::m_OutputSpacing
SpacingType m_OutputSpacing
Definition: itkTransformToSpatialJacobianSource.h:210
itk::TransformToSpatialJacobianSource::TransformType
AdvancedTransform< TTransformPrecisionType, itkGetStaticConstMacro(ImageDimension), itkGetStaticConstMacro(ImageDimension) > TransformType
Definition: itkTransformToSpatialJacobianSource.h:96
itk::TransformToSpatialJacobianSource::SetOutputOrigin
virtual void SetOutputOrigin(const double *values)
itk::TransformToSpatialJacobianSource::m_OutputOrigin
OriginType m_OutputOrigin
Definition: itkTransformToSpatialJacobianSource.h:211
itk::TransformToSpatialJacobianSource::SpacingType
OutputImageType::SpacingType SpacingType
Definition: itkTransformToSpatialJacobianSource.h:107
itk::TransformToSpatialJacobianSource::DirectionType
OutputImageType::DirectionType DirectionType
Definition: itkTransformToSpatialJacobianSource.h:109
itk::TransformToSpatialJacobianSource::SetOutputParametersFromImage
void SetOutputParametersFromImage(const ImageBaseType *image)
itk::TransformToSpatialJacobianSource::~TransformToSpatialJacobianSource
~TransformToSpatialJacobianSource()
Definition: itkTransformToSpatialJacobianSource.h:179
itk::TransformToSpatialJacobianSource::OutputImageType
TOutputImage OutputImageType
Definition: itkTransformToSpatialJacobianSource.h:78
itk::TransformToSpatialJacobianSource::OutputImageConstPointer
OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkTransformToSpatialJacobianSource.h:80
itk::TransformToSpatialJacobianSource::GetOutputIndex
virtual const IndexType & GetOutputIndex()
itk::TransformToSpatialJacobianSource::OriginType
OutputImageType::PointType OriginType
Definition: itkTransformToSpatialJacobianSource.h:108
itk::TransformToSpatialJacobianSource::RegionType
OutputImageType::RegionType RegionType
Definition: itkTransformToSpatialJacobianSource.h:103
itk::TransformToSpatialJacobianSource
Generate the spatial Jacobian matrix from a coordinate transform.
Definition: itkTransformToSpatialJacobianSource.h:69
itk::TransformToSpatialJacobianSource::SetOutputIndex
virtual void SetOutputIndex(const IndexType &index)
itk::TransformToSpatialJacobianSource::SetOutputSpacing
virtual void SetOutputSpacing(const double *values)
itk::TransformToSpatialJacobianSource::SizeType
RegionType::SizeType SizeType
Definition: itkTransformToSpatialJacobianSource.h:104
itk::TransformToSpatialJacobianSource::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
itk::TransformToSpatialJacobianSource::m_OutputRegion
RegionType m_OutputRegion
Definition: itkTransformToSpatialJacobianSource.h:208
itk::AdvancedTransform
Transform maps points, vectors and covariant vectors from an input space to an output space.
Definition: itkAdvancedTransform.h:87
itk::TransformToSpatialJacobianSource::SetOutputSize
virtual void SetOutputSize(const SizeType &size)
ThreadIdType
itk::TransformToSpatialJacobianSource::TransformToSpatialJacobianSource
TransformToSpatialJacobianSource(const Self &)
itk::TransformToSpatialJacobianSource::ImageBaseType
ImageBase< itkGetStaticConstMacro(ImageDimension) > ImageBaseType
Definition: itkTransformToSpatialJacobianSource.h:112
itk::TransformToSpatialJacobianSource::GenerateOutputInformation
virtual void GenerateOutputInformation(void)
itk::TransformToSpatialJacobianSource::IndexType
OutputImageType::IndexType IndexType
Definition: itkTransformToSpatialJacobianSource.h:105
itk::TransformToSpatialJacobianSource::ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkTransformToSpatialJacobianSource.h:76
itk::TransformToSpatialJacobianSource::SpatialJacobianType
TransformType::SpatialJacobianType SpatialJacobianType
Definition: itkTransformToSpatialJacobianSource.h:98
itk::TransformToSpatialJacobianSource::operator=
void operator=(const Self &)
itk::TransformToSpatialJacobianSource::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const
itk::TransformToSpatialJacobianSource::m_OutputDirection
DirectionType m_OutputDirection
Definition: itkTransformToSpatialJacobianSource.h:212
itk::TransformToSpatialJacobianSource::GetOutputSize
virtual const SizeType & GetOutputSize()
itk
Definition: itkAdvancedImageToImageMetric.h:40
itk::TransformToSpatialJacobianSource::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: itkTransformToSpatialJacobianSource.h:81
itk::TransformToSpatialJacobianSource::LinearGenerateData
void LinearGenerateData(void)
itk::TransformToSpatialJacobianSource::NonlinearThreadedGenerateData
void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
itk::TransformToSpatialJacobianSource::GetMTime
unsigned long GetMTime(void) const
itk::TransformToSpatialJacobianSource::BeforeThreadedGenerateData
virtual void BeforeThreadedGenerateData(void)
itk::AdvancedTransform::ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAdvancedTransform.h:96
itk::TransformToSpatialJacobianSource::OutputImagePointer
OutputImageType::Pointer OutputImagePointer
Definition: itkTransformToSpatialJacobianSource.h:79
itk::TransformToSpatialJacobianSource::PixelType
OutputImageType::PixelType PixelType
Definition: itkTransformToSpatialJacobianSource.h:101
itk::TransformToSpatialJacobianSource::PointType
OutputImageType::PointType PointType
Definition: itkTransformToSpatialJacobianSource.h:106


Generated on OURCE_DATE_EPOCH for elastix by doxygen 1.8.18 elastix logo