ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkWarpVectorImageFilter.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 itkWarpVectorImageFilter_h
19 #define itkWarpVectorImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 #include "itkPoint.h"
24 #include "itkFixedArray.h"
25 
26 namespace itk
27 {
87 template<
88  typename TInputImage,
89  typename TOutputImage,
90  typename TDisplacementField
91  >
92 class ITK_TEMPLATE_EXPORT WarpVectorImageFilter:
93  public ImageToImageFilter< TInputImage, TOutputImage >
94 {
95 public:
101 
103  itkNewMacro(Self);
104 
107 
109  typedef typename TOutputImage::RegionType OutputImageRegionType;
110 
112  typedef typename Superclass::InputImageType InputImageType;
113  typedef typename Superclass::InputImagePointer InputImagePointer;
114  typedef typename Superclass::OutputImageType OutputImageType;
115  typedef typename Superclass::OutputImagePointer OutputImagePointer;
116  typedef typename Superclass::InputImageConstPointer InputImageConstPointer;
117 
120  typedef typename OutputImageType::PixelType PixelType;
121  typedef typename OutputImageType::SpacingType SpacingType;
122  typedef typename OutputImageType::PixelType::ValueType ValueType;
123 
125  itkStaticConstMacro(ImageDimension, unsigned int,
126  TOutputImage::ImageDimension);
127 
129  itkStaticConstMacro(PixelDimension, unsigned int,
131 
133  typedef TDisplacementField DisplacementFieldType;
134  typedef typename DisplacementFieldType::Pointer DisplacementFieldPointer;
135  typedef typename DisplacementFieldType::PixelType DisplacementType;
136 
137 #ifdef ITKV3_COMPATIBILITY
138  typedef TDisplacementField DeformationFieldType;
139  typedef typename DeformationFieldType::Pointer DeformationFieldPointer;
140  typedef typename DeformationFieldType::PixelType DeformationType;
141 #endif
142 
144  typedef double CoordRepType;
149 
152 
155 
157  void SetDisplacementField(const DisplacementFieldType *field);
158 
160  void SetDisplacementField(DisplacementFieldType *field);
161 
163  DisplacementFieldType * GetDisplacementField();
164 
165 #ifdef ITKV3_COMPATIBILITY
166  void SetDeformationField(const DeformationFieldType *field)
167  {
168  this->SetDisplacementField(field);
169  }
170 
171  void SetDeformationField(DeformationFieldType *field)
172  {
173  this->SetDisplacementField(field);
174  }
175 
176  DeformationFieldType * GetDeformationField(void)
177  {
178  return static_cast<DeformationFieldType *> (GetDisplacementField());
179  }
180 #endif
181 
183  itkSetObjectMacro(Interpolator, InterpolatorType);
184  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
186 
188  itkSetMacro(OutputSpacing, SpacingType);
189  virtual void SetOutputSpacing(const double *values);
191 
193  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
194 
196  itkSetMacro(OutputOrigin, PointType);
197  virtual void SetOutputOrigin(const double *values);
199 
201  itkGetConstReferenceMacro(OutputOrigin, PointType);
202 
204  itkSetMacro(OutputDirection, DirectionType);
205  itkGetConstReferenceMacro(OutputDirection, DirectionType);
207 
209  itkSetMacro(EdgePaddingValue, PixelType);
210 
212  itkGetConstMacro(EdgePaddingValue, PixelType);
213 
219  virtual void GenerateOutputInformation() ITK_OVERRIDE;
220 
227  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
228 
231  virtual void BeforeThreadedGenerateData() ITK_OVERRIDE;
232 
233 #ifdef ITK_USE_CONCEPT_CHECKING
234  // Begin concept checking
235  itkConceptMacro( InputHasNumericTraitsCheck,
236  ( Concept::HasNumericTraits< typename TInputImage::PixelType::ValueType > ) );
237  itkConceptMacro( OutputHasNumericTraitsCheck,
238  ( Concept::HasNumericTraits< ValueType > ) );
239  itkConceptMacro( DisplacementFieldHasNumericTraitsCheck,
240  ( Concept::HasNumericTraits< typename TDisplacementField::PixelType::ValueType > ) );
241  // End concept checking
242 #endif
243 
244 protected:
245  WarpVectorImageFilter();
246  ~WarpVectorImageFilter() ITK_OVERRIDE {}
247  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
248 
252  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
253  ThreadIdType threadId) ITK_OVERRIDE;
254 
255 private:
256  ITK_DISALLOW_COPY_AND_ASSIGN(WarpVectorImageFilter);
257 
262 
264 };
265 } // end namespace itk
266 
267 #ifndef ITK_MANUAL_INSTANTIATION
268 #include "itkWarpVectorImageFilter.hxx"
269 #endif
270 
271 #endif
Superclass::InputImageConstPointer InputImageConstPointer
InterpolatorType::Pointer InterpolatorPointer
Superclass::InputImageType InputImageType
OutputImageType::IndexType IndexType
Base class for all process objects that output image data.
Superclass::OutputImageType OutputImageType
TOutputImage::DirectionType DirectionType
Point< CoordRepType, itkGetStaticConstMacro(ImageDimension) > PointType
TOutputImage::RegionType OutputImageRegionType
Warps an image using an input displacement field.
VectorInterpolateImageFunction< InputImageType, CoordRepType > InterpolatorType
DisplacementFieldType::Pointer DisplacementFieldPointer
OutputImageType::PixelType::ValueType ValueType
Superclass::OutputImagePointer OutputImagePointer
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Superclass::InputImagePointer InputImagePointer
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Base class for all vector image interpolaters.
OutputImageType::PixelType PixelType
VectorLinearInterpolateImageFunction< InputImageType, CoordRepType > DefaultInterpolatorType
OutputImageType::SizeType SizeType
#define itkConceptMacro(name, concept)
SmartPointer< const Self > ConstPointer
DisplacementFieldType::PixelType DisplacementType
OutputImageType::SpacingType SpacingType
Linearly interpolate a vector image at specified positions.