ITK  6.0.0
Insight Toolkit
itkWarpVectorImageFilter.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
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 * https://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
23#include "itkPoint.h"
24#include "itkFixedArray.h"
25
26namespace itk
27{
88template <typename TInputImage, typename TOutputImage, typename TDisplacementField>
89class ITK_TEMPLATE_EXPORT WarpVectorImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
90{
91public:
92 ITK_DISALLOW_COPY_AND_MOVE(WarpVectorImageFilter);
93
99
101 itkNewMacro(Self);
102
104 itkOverrideGetNameOfClassMacro(WarpVectorImageFilter);
105
108
110 using typename Superclass::InputImageType;
111 using typename Superclass::InputImagePointer;
112 using typename Superclass::OutputImageType;
113 using typename Superclass::OutputImagePointer;
114 using typename Superclass::InputImageConstPointer;
115
118 using PixelType = typename OutputImageType::PixelType;
119 using SpacingType = typename OutputImageType::SpacingType;
120 using ValueType = typename OutputImageType::PixelType::ValueType;
121
123 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
124
126 static constexpr unsigned int PixelDimension = PixelType::Dimension;
127
129 using DisplacementFieldType = TDisplacementField;
131 using DisplacementType = typename DisplacementFieldType::PixelType;
132
134 using CoordRepType = double;
138
141
144
146 void
148
150 void
152
156
158 itkSetObjectMacro(Interpolator, InterpolatorType);
159 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
163 itkSetMacro(OutputSpacing, SpacingType);
164 virtual void
165 SetOutputSpacing(const double * spacing);
169 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
170
172 itkSetMacro(OutputOrigin, PointType);
173 virtual void
174 SetOutputOrigin(const double * origin);
178 itkGetConstReferenceMacro(OutputOrigin, PointType);
179
181 itkSetMacro(OutputDirection, DirectionType);
182 itkGetConstReferenceMacro(OutputDirection, DirectionType);
186 itkSetMacro(EdgePaddingValue, PixelType);
187
189 itkGetConstMacro(EdgePaddingValue, PixelType);
190
196 void
198
205 void
207
210 void
212
213#ifdef ITK_USE_CONCEPT_CHECKING
214 // Begin concept checking
216 itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<ValueType>));
217 itkConceptMacro(DisplacementFieldHasNumericTraitsCheck,
219 // End concept checking
220#endif
221
222protected:
224 ~WarpVectorImageFilter() override = default;
225 void
226 PrintSelf(std::ostream & os, Indent indent) const override;
227
231 void
232 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
233
234
235private:
236 PixelType m_EdgePaddingValue{};
237 SpacingType m_OutputSpacing{};
238 PointType m_OutputOrigin{};
239 DirectionType m_OutputDirection{};
240
241 InterpolatorPointer m_Interpolator{};
242};
243} // end namespace itk
244
245#ifndef ITK_MANUAL_INSTANTIATION
246# include "itkWarpVectorImageFilter.hxx"
247#endif
248
249#endif
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all vector image interpolators.
Linearly interpolate a vector image at specified positions.
Warps an image using an input displacement field.
typename OutputImageType::SizeType SizeType
void GenerateInputRequestedRegion() override
void PrintSelf(std::ostream &os, Indent indent) const override
DisplacementFieldType * GetDisplacementField()
typename DisplacementFieldType::PixelType DisplacementType
void BeforeThreadedGenerateData() override
typename InterpolatorType::Pointer InterpolatorPointer
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
typename OutputImageType::IndexType IndexType
void SetDisplacementField(const DisplacementFieldType *field)
typename OutputImageType::PixelType::ValueType ValueType
void SetDisplacementField(DisplacementFieldType *field)
typename TOutputImage::DirectionType DirectionType
typename DisplacementFieldType::Pointer DisplacementFieldPointer
typename OutputImageType::SpacingType SpacingType
virtual void SetOutputOrigin(const double *origin)
void GenerateOutputInformation() override
~WarpVectorImageFilter() override=default
typename OutputImageType::PixelType PixelType
virtual void SetOutputSpacing(const double *spacing)
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....