ITK  4.3.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  class TInputImage,
85  class TOutputImage,
86  class TDisplacementField
87  >
88 class ITK_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::SpacingType SpacingType;
118 
120  itkStaticConstMacro(ImageDimension, unsigned int,
121  TOutputImage::ImageDimension);
122  itkStaticConstMacro(InputImageDimension, unsigned int,
123  TInputImage::ImageDimension);
124  itkStaticConstMacro(DisplacementFieldDimension, unsigned int,
125  TDisplacementField::ImageDimension);
126 
129 
131  typedef TDisplacementField DisplacementFieldType;
132  typedef typename DisplacementFieldType::Pointer DisplacementFieldPointer;
133  typedef typename DisplacementFieldType::PixelType DisplacementType;
134 
135 #ifdef ITKV3_COMPATIBILITY
136  typedef TDisplacementField DeformationFieldType;
137  typedef typename DeformationFieldType::Pointer DeformationFieldPointer;
138  typedef typename DeformationFieldType::PixelType DeformationType;
139 #endif
140 
142  typedef double CoordRepType;
147 
150 
152  typedef typename TOutputImage::DirectionType DirectionType;
153 
155  void SetDisplacementField(const DisplacementFieldType *field);
156 
158  DisplacementFieldType * GetDisplacementField(void);
159 
160 #ifdef ITKV3_COMPATIBILITY
161  void SetDeformationField(const DisplacementFieldType *field)
162  {
163  this->SetDisplacementField(field);
164  }
165  DeformationFieldType * GetDeformationField(void)
166  {
167  return static_cast<DeformationFieldType *> (GetDisplacementField());
168  }
169 #endif
170 
172  itkSetObjectMacro(Interpolator, InterpolatorType);
173 
175  itkGetObjectMacro(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();
226 
233  virtual void GenerateInputRequestedRegion();
234 
237  virtual void BeforeThreadedGenerateData();
238 
241  virtual void AfterThreadedGenerateData();
242 
243 #ifdef ITK_USE_CONCEPT_CHECKING
244 
245  itkConceptMacro( SameDimensionCheck1,
246  ( Concept::SameDimension< ImageDimension, InputImageDimension > ) );
247  itkConceptMacro( SameDimensionCheck2,
248  ( Concept::SameDimension< ImageDimension, DisplacementFieldDimension > ) );
249  itkConceptMacro( InputHasNumericTraitsCheck,
250  ( Concept::HasNumericTraits< typename TInputImage::PixelType > ) );
251  itkConceptMacro( DisplacementFieldHasNumericTraitsCheck,
252  ( Concept::HasNumericTraits< typename TDisplacementField::PixelType::ValueType > ) );
253 
255 #endif
256 
257 protected:
258  WarpImageFilter();
259  // ~WarpImageFilter() {} default implementation is ok
260 
261  void PrintSelf(std::ostream & os, Indent indent) const;
262 
266  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
267  ThreadIdType threadId);
268 
274  virtual void VerifyInputInformation() {}
275 
276 private:
277  WarpImageFilter(const Self &); //purposely not implemented
278  void operator=(const Self &); //purposely not implemented
279 
283  DisplacementType EvaluateDisplacementAtPhysicalPoint(const PointType & p);
284 
289 
291  SizeType m_OutputSize; // Size of the output image
292  IndexType m_OutputStartIndex; // output image start index
294  // variables for deffield interpoator
295  IndexType m_StartIndex, m_EndIndex;
296 };
297 } // end namespace itk
298 
299 #ifndef ITK_MANUAL_INSTANTIATION
300 #include "itkWarpImageFilter.hxx"
301 #endif
302 
303 #endif
304