ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkTransformToDeformationFieldSource.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 itkTransformToDeformationFieldSource_h
19 #define itkTransformToDeformationFieldSource_h
20 
21 #include "itkTransform.h"
22 #include "itkImageSource.h"
23 
24 #ifndef ITKV3_COMPATIBILITY
25 #error "This file is only valid when ITKV3_COMPATIBILITY is turned on. Users are encouraged to convert to itkTransformToDisplacementFieldSource.h in ITKv4"
26 #endif
27 
28 namespace itk
29 {
62 template< typename TOutputImage,
63  typename TTransformPrecisionType = double >
65  public ImageSource< TOutputImage >
66 {
67 public:
73 
74  typedef TOutputImage OutputImageType;
75  typedef typename OutputImageType::Pointer OutputImagePointer;
76  typedef typename OutputImageType::ConstPointer OutputImageConstPointer;
77  typedef typename OutputImageType::RegionType OutputImageRegionType;
78 
80  itkNewMacro(Self);
81 
84 
86  itkStaticConstMacro(ImageDimension, unsigned int,
87  TOutputImage::ImageDimension);
88 
90  typedef Transform< TTransformPrecisionType,
91  itkGetStaticConstMacro(ImageDimension),
92  itkGetStaticConstMacro(ImageDimension) > TransformType;
95 
97  typedef typename OutputImageType::PixelType PixelType;
98  typedef typename PixelType::ValueType PixelValueType;
99  typedef typename OutputImageType::RegionType RegionType;
100  typedef typename RegionType::SizeType SizeType;
101  typedef typename OutputImageType::IndexType IndexType;
102  typedef typename OutputImageType::PointType PointType;
103  typedef typename OutputImageType::SpacingType SpacingType;
104  typedef typename OutputImageType::PointType OriginType;
105  typedef typename OutputImageType::DirectionType DirectionType;
106 
109 
117  itkSetConstObjectMacro(Transform, TransformType);
118 
120  itkGetConstObjectMacro(Transform, TransformType);
121 
123  virtual void SetOutputSize(const SizeType & size);
124 
126  virtual const SizeType & GetOutputSize();
127 
130  virtual void SetOutputIndex(const IndexType & index);
131 
133  virtual const IndexType & GetOutputIndex();
134 
136  itkSetMacro(OutputRegion, OutputImageRegionType);
137 
139  itkGetConstReferenceMacro(OutputRegion, OutputImageRegionType);
140 
142  itkSetMacro(OutputSpacing, SpacingType);
143  virtual void SetOutputSpacing(const double *values);
145 
147  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
148 
150  itkSetMacro(OutputOrigin, OriginType);
151  virtual void SetOutputOrigin(const double *values);
153 
155  itkGetConstReferenceMacro(OutputOrigin, OriginType);
156 
158  itkSetMacro(OutputDirection, DirectionType);
159  itkGetConstReferenceMacro(OutputDirection, DirectionType);
161 
163  void SetOutputParametersFromImage(const ImageBaseType *image);
164 
166  virtual void GenerateOutputInformation();
167 
169  virtual void BeforeThreadedGenerateData();
170 
172  unsigned long GetMTime() const;
173 
174 #ifdef ITK_USE_CONCEPT_CHECKING
175  // Begin concept checking
176  itkStaticConstMacro(PixelDimension, unsigned int,
177  PixelType::Dimension);
178  itkConceptMacro( SameDimensionCheck,
180  // End concept checking
181 #endif
182 
183 protected:
186 
187  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
188 
193  const OutputImageRegionType & outputRegionForThread,
194  ThreadIdType threadId);
195 
200  const OutputImageRegionType & outputRegionForThread,
201  ThreadIdType threadId);
202 
207  const OutputImageRegionType & outputRegionForThread,
208  ThreadIdType threadId);
209 
210 private:
211 
212  TransformToDeformationFieldSource(const Self &) ITK_DELETE_FUNCTION;
213  void operator=(const Self &) ITK_DELETE_FUNCTION;
214 
216  RegionType m_OutputRegion; // region of the output image
217  TransformPointerType m_Transform; // Coordinate transform to use
218  SpacingType m_OutputSpacing; // output image spacing
219  OriginType m_OutputOrigin; // output image origin
220  DirectionType m_OutputDirection; // output image direction cosines
221 }; // end class
222  // TransformToDeformationFieldSource
223 } // end namespace itk
224 
225 #ifndef ITK_MANUAL_INSTANTIATION
226 #include "itkTransformToDeformationFieldSource.hxx"
227 #endif
228 
229 #endif // end #ifndef itkTransformToDeformationFieldSource_h
virtual void SetOutputSize(const SizeType &size)
virtual void SetOutputOrigin(OriginType _arg)
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
void PrintSelf(std::ostream &os, Indent indent) const override
virtual const SizeType & GetOutputSize()
virtual void SetOutputSpacing(SpacingType _arg)
Base class for all process objects that output image data.
void SetOutputParametersFromImage(const ImageBaseType *image)
ImageBase< itkGetStaticConstMacro(ImageDimension) > ImageBaseType
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
void LinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
Base class for templated image classes.
Definition: itkImageBase.h:112
Generate a deformation field from a coordinate transform.
virtual const IndexType & GetOutputIndex()
OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Transform< TTransformPrecisionType, itkGetStaticConstMacro(ImageDimension), itkGetStaticConstMacro(ImageDimension) > TransformType
virtual void SetOutputIndex(const IndexType &index)
#define itkConceptMacro(name, concept)