ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkTransformToDisplacementFieldFilter.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 __itkTransformToDisplacementFieldFilter_h
19 #define __itkTransformToDisplacementFieldFilter_h
20 
21 #include "itkDataObjectDecorator.h"
22 #include "itkTransform.h"
23 #include "itkImageSource.h"
24 
25 namespace itk
26 {
54 template< typename TOutputImage,
55  typename TScalar = double >
57  public ImageSource< TOutputImage >
58 {
59 public:
65 
66  typedef TOutputImage OutputImageType;
67  typedef typename OutputImageType::RegionType OutputImageRegionType;
68 
70  itkNewMacro(Self);
71 
74 
76  itkStaticConstMacro(ImageDimension, unsigned int, TOutputImage::ImageDimension);
77 
81 
83  typedef typename OutputImageType::PixelType PixelType;
84  typedef typename PixelType::ValueType PixelValueType;
85  typedef typename OutputImageType::RegionType RegionType;
86  typedef typename RegionType::SizeType SizeType;
87  typedef typename OutputImageType::IndexType IndexType;
88  typedef typename OutputImageType::PointType PointType;
89  typedef typename OutputImageType::SpacingType SpacingType;
90  typedef typename OutputImageType::PointType OriginType;
91  typedef typename OutputImageType::DirectionType DirectionType;
92 
95 
100  using Superclass::SetInput;
101  virtual void SetInput( const TransformInputType * transform );
102  const TransformInputType * GetInput() const;
105 
108  itkSetMacro(OutputStartIndex, IndexType);
109  itkGetConstReferenceMacro(OutputStartIndex, IndexType);
111 
113  itkSetMacro(Size, SizeType);
114  itkGetConstReferenceMacro(Size, SizeType);
116 
118  itkSetMacro(OutputSpacing, SpacingType);
119  virtual void SetOutputSpacing(const SpacePrecisionType *values);
121 
123  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
124 
126  itkSetMacro(OutputOrigin, OriginType);
127  virtual void SetOutputOrigin(const SpacePrecisionType *values);
129 
131  itkGetConstReferenceMacro(OutputOrigin, OriginType);
132 
134  itkSetMacro(OutputDirection, DirectionType);
135  itkGetConstReferenceMacro(OutputDirection, DirectionType);
137 
144  itkSetInputMacro(ReferenceImage, ReferenceImageBaseType);
145 
147  itkGetInputMacro(ReferenceImage, ReferenceImageBaseType);
148 
151  itkSetMacro(UseReferenceImage, bool);
152  itkBooleanMacro(UseReferenceImage);
153  itkGetConstMacro(UseReferenceImage, bool);
155 
156 #ifdef ITK_USE_CONCEPT_CHECKING
157  // Begin concept checking
158  itkStaticConstMacro(PixelDimension, unsigned int,
159  PixelType::Dimension);
160  itkConceptMacro( SameDimensionCheck,
162  // End concept checking
163 #endif
164 
165 protected:
168 
170  virtual void GenerateOutputInformation();
171 
175  virtual void ThreadedGenerateData( const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId );
176 
180  void NonlinearThreadedGenerateData( const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId );
181 
185  void LinearThreadedGenerateData( const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId );
186 
187  virtual void PrintSelf(std::ostream & os, Indent indent) const;
188 
189 private:
190  TransformToDisplacementFieldFilter( const Self & ); //purposely not implemented
191  void operator=( const Self & ); //purposely not implemented
192 
194  SizeType m_Size; // size of the output region
195  IndexType m_OutputStartIndex; // start index of the output region
196  SpacingType m_OutputSpacing; // output image spacing
197  OriginType m_OutputOrigin; // output image origin
198  DirectionType m_OutputDirection; // output image direction cosines
200 
201 };
202 } // end namespace itk
203 
204 #ifndef ITK_MANUAL_INSTANTIATION
205 #include "itkTransformToDisplacementFieldFilter.hxx"
206 #endif
207 
208 #endif
itkSetGetDecoratedObjectInputMacro(Transform, TransformType)
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Generate a displacement field from a coordinate transform.
void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
Base class for all process objects that output image data.
virtual void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
virtual void SetOutputOrigin(OriginType _arg)
Decorates any subclass of itkObject with a DataObject API.
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Transform< TScalar, ImageDimension, ImageDimension > TransformType
virtual void SetOutputSpacing(SpacingType _arg)
virtual void SetInput(const TransformInputType *transform)
const TransformInputType * GetInput() const
Base class for templated image classes.
Definition: itkImageBase.h:112
void LinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
Control indentation during Print() invocation.
Definition: itkIndent.h:49
double SpacePrecisionType
Definition: itkFloatTypes.h:30
#define itkConceptMacro(name, concept)
virtual void PrintSelf(std::ostream &os, Indent indent) const
unsigned int ThreadIdType
Definition: itkIntTypes.h:159