ITK  4.12.0
Insight Segmentation and Registration Toolkit
itkWarpImageFilter.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 itkWarpImageFilter_h
19 #define itkWarpImageFilter_h
20 #include "itkImageBase.h"
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
83 template<
84  typename TInputImage,
85  typename TOutputImage,
86  typename TDisplacementField
87  >
88 class ITK_TEMPLATE_EXPORT WarpImageFilter:
89  public ImageToImageFilter< TInputImage, TOutputImage >
90 {
91 public:
97 
99  itkNewMacro(Self);
100 
103 
105  typedef typename TOutputImage::RegionType OutputImageRegionType;
106 
108  typedef typename Superclass::InputImageType InputImageType;
109  typedef typename Superclass::InputImagePointer InputImagePointer;
110  typedef typename Superclass::OutputImageType OutputImageType;
111  typedef typename Superclass::OutputImagePointer OutputImagePointer;
112  typedef typename Superclass::InputImageConstPointer InputImageConstPointer;
113  typedef typename OutputImageType::IndexType IndexType;
115  typedef typename OutputImageType::SizeType SizeType;
116  typedef typename OutputImageType::PixelType PixelType;
117  typedef typename OutputImageType::InternalPixelType PixelComponentType;
118  typedef typename OutputImageType::SpacingType SpacingType;
119 
121  itkStaticConstMacro(ImageDimension, unsigned int,
122  TOutputImage::ImageDimension);
123  itkStaticConstMacro(InputImageDimension, unsigned int,
124  TInputImage::ImageDimension);
125  itkStaticConstMacro(DisplacementFieldDimension, unsigned int,
126  TDisplacementField::ImageDimension);
127 
130 
132  typedef TDisplacementField DisplacementFieldType;
134  typedef typename DisplacementFieldType::PixelType DisplacementType;
135 
136 #ifdef ITKV3_COMPATIBILITY
137  typedef TDisplacementField DeformationFieldType;
138  typedef typename DeformationFieldType::Pointer DeformationFieldPointer;
139  typedef typename DeformationFieldType::PixelType DeformationType;
140 #endif
141 
143  typedef double CoordRepType;
148 
151 
153  typedef typename TOutputImage::DirectionType DirectionType;
154 
156  void SetDisplacementField(const DisplacementFieldType *field);
157 
159  DisplacementFieldType * GetDisplacementField();
160 
161 #ifdef ITKV3_COMPATIBILITY
162  void SetDeformationField(const DisplacementFieldType *field)
163  {
164  this->SetDisplacementField(field);
165  }
166  DeformationFieldType * GetDeformationField(void)
167  {
168  return static_cast<DeformationFieldType *> (GetDisplacementField());
169  }
170 #endif
171 
173  itkSetObjectMacro(Interpolator, InterpolatorType);
174  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
176 
178  itkSetMacro(OutputSpacing, SpacingType);
179  virtual void SetOutputSpacing(const double *values);
181 
183  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
184 
186  itkSetMacro(OutputOrigin, PointType);
187  virtual void SetOutputOrigin(const double *values);
189 
191  itkGetConstReferenceMacro(OutputOrigin, PointType);
192 
194  itkSetMacro(OutputDirection, DirectionType);
195  itkGetConstReferenceMacro(OutputDirection, DirectionType);
197 
199  void SetOutputParametersFromImage(const ImageBaseType *image);
200 
203  itkSetMacro(OutputStartIndex, IndexType);
204 
206  itkGetConstReferenceMacro(OutputStartIndex, IndexType);
207 
209  itkSetMacro(OutputSize, SizeType);
210 
212  itkGetConstReferenceMacro(OutputSize, SizeType);
213 
215  itkSetMacro(EdgePaddingValue, PixelType);
216 
218  itkGetConstMacro(EdgePaddingValue, PixelType);
219 
225  virtual void GenerateOutputInformation() ITK_OVERRIDE;
226 
233  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
234 
237  virtual void BeforeThreadedGenerateData() ITK_OVERRIDE;
238 
241  virtual void AfterThreadedGenerateData() ITK_OVERRIDE;
242 
243 #ifdef ITK_USE_CONCEPT_CHECKING
244  // Begin concept checking
245  itkConceptMacro( SameDimensionCheck1,
246  ( Concept::SameDimension< ImageDimension, InputImageDimension > ) );
247  itkConceptMacro( SameDimensionCheck2,
248  ( Concept::SameDimension< ImageDimension, DisplacementFieldDimension > ) );
249  itkConceptMacro( InputHasNumericTraitsCheck,
250  ( Concept::HasNumericTraits< typename TInputImage::InternalPixelType > ) );
251  itkConceptMacro( DisplacementFieldHasNumericTraitsCheck,
252  ( Concept::HasNumericTraits< typename TDisplacementField::PixelType::ValueType > ) );
253  // End concept checking
254 #endif
255 
256 protected:
257  WarpImageFilter();
258  // ~WarpImageFilter() {} default implementation is ok
259 
260  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
261 
265  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
266  ThreadIdType threadId) ITK_OVERRIDE;
267 
274  virtual void VerifyInputInformation() ITK_OVERRIDE;
275 
284  void EvaluateDisplacementAtPhysicalPoint(const PointType & p, DisplacementType & output);
285 
289  void EvaluateDisplacementAtPhysicalPoint(const PointType & p, const DisplacementFieldType * fieldPtr,
290  DisplacementType & output);
291 
292  bool m_DefFieldSameInformation;
293  // variables for deffield interpoator
294  IndexType m_StartIndex, m_EndIndex;
295 
296 private:
297  ITK_DISALLOW_COPY_AND_ASSIGN(WarpImageFilter);
298 
299  PixelType m_EdgePaddingValue;
300  SpacingType m_OutputSpacing;
301  PointType m_OutputOrigin;
302  DirectionType m_OutputDirection;
303 
304  InterpolatorPointer m_Interpolator;
305  SizeType m_OutputSize; // Size of the output image
306  IndexType m_OutputStartIndex; // output image start index
307 
308 };
309 } // end namespace itk
310 
311 #ifndef ITK_MANUAL_INSTANTIATION
312 #include "itkWarpImageFilter.hxx"
313 #endif
314 
315 #endif
virtual void PrintSelf(std::ostream &os, Indent indent) const override
Superclass::InputImageType InputImageType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
TDisplacementField DisplacementFieldType
OutputImageType::InternalPixelType PixelComponentType
OutputImageType::SizeType SizeType
signed long IndexValueType
Definition: itkIntTypes.h:150
InterpolateImageFunction< InputImageType, CoordRepType > InterpolatorType
SmartPointer< const Self > ConstPointer
TOutputImage::RegionType OutputImageRegionType
OutputImageType::IndexType IndexType
Base class for all process objects that output image data.
SmartPointer< Self > Pointer
ImageBase< itkGetStaticConstMacro(ImageDimension) > ImageBaseType
Superclass::InputImagePointer InputImagePointer
DisplacementFieldType::PixelType DisplacementType
SmartPointer< Self > Pointer
OutputImageType::PixelType PixelType
TOutputImage::DirectionType DirectionType
OutputImageType::IndexValueType IndexValueType
OutputImageType::SpacingType SpacingType
Warps an image using an input displacement field.
DisplacementFieldType::Pointer DisplacementFieldPointer
Linearly interpolate an image at specified positions.
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Point< CoordRepType, itkGetStaticConstMacro(ImageDimension) > PointType
Base class for all image interpolaters.
Base class for templated image classes.
Definition: itkImageBase.h:114
InterpolatorType::Pointer InterpolatorPointer
Superclass::OutputImagePointer OutputImagePointer
LinearInterpolateImageFunction< InputImageType, CoordRepType > DefaultInterpolatorType
Base class for filters that take an image as input and produce an image as output.
Superclass::OutputImageType OutputImageType
#define itkConceptMacro(name, concept)
Superclass::InputImageConstPointer InputImageConstPointer