ITK  6.0.0
Insight Toolkit
itkTransformToDisplacementFieldFilter.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 itkTransformToDisplacementFieldFilter_h
19#define itkTransformToDisplacementFieldFilter_h
20
22#include "itkTransform.h"
23#include "itkImageSource.h"
24
25namespace itk
26{
54template <typename TOutputImage, typename TParametersValueType = double>
55class ITK_TEMPLATE_EXPORT TransformToDisplacementFieldFilter : public ImageSource<TOutputImage>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(TransformToDisplacementFieldFilter);
59
65
66 using OutputImageType = TOutputImage;
68
70 itkNewMacro(Self);
71
73 itkOverrideGetNameOfClassMacro(TransformToDisplacementFieldFilter);
74
76 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
77
81
83 using PixelType = typename OutputImageType::PixelType;
84 using PixelValueType = typename PixelType::ValueType;
89 using SpacingType = typename OutputImageType::SpacingType;
92
95
100 using Superclass::SetInput;
101 virtual void
103 const TransformInputType *
104 GetInput() const;
105 itkSetGetDecoratedObjectInputMacro(Transform, TransformType);
110 itkSetMacro(OutputStartIndex, IndexType);
111 itkGetConstReferenceMacro(OutputStartIndex, IndexType);
115 itkSetMacro(Size, SizeType);
116 itkGetConstReferenceMacro(Size, SizeType);
120 itkSetMacro(OutputSpacing, SpacingType);
121 virtual void
126 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
127
129 itkSetMacro(OutputOrigin, OriginType);
130 virtual void
135 itkGetConstReferenceMacro(OutputOrigin, OriginType);
136
138 itkSetMacro(OutputDirection, DirectionType);
139 itkGetConstReferenceMacro(OutputDirection, DirectionType);
149
152
155 itkSetMacro(UseReferenceImage, bool);
156 itkBooleanMacro(UseReferenceImage);
157 itkGetConstMacro(UseReferenceImage, bool);
160#ifdef ITK_USE_CONCEPT_CHECKING
161 // Begin concept checking
162 static constexpr unsigned int PixelDimension = PixelType::Dimension;
164 // End concept checking
165#endif
166
167protected:
170
172 void
174
176 void
177 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
178
179
183 void
185
189 void
190 LinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread);
191
192 void
193 PrintSelf(std::ostream & os, Indent indent) const override;
194
195private:
197 SizeType m_Size{}; // size of the output region
198 IndexType m_OutputStartIndex{}; // start index of the output region
199 SpacingType m_OutputSpacing{}; // output image spacing
200 OriginType m_OutputOrigin{}; // output image origin
201 DirectionType m_OutputDirection{}; // output image direction cosines
202 bool m_UseReferenceImage{ false };
203};
204} // end namespace itk
205
206#ifndef ITK_MANUAL_INSTANTIATION
207# include "itkTransformToDisplacementFieldFilter.hxx"
208#endif
209
210#endif
Decorates any subclass of itkObject with a DataObject API.
Base class for templated image classes.
Definition: itkImageBase.h:115
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
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...
Generate a displacement field from a coordinate transform.
virtual void SetInput(const TransformInputType *input)
typename OutputImageType::DirectionType DirectionType
const TransformInputType * GetInput() const
virtual void SetOutputOrigin(const SpacePrecisionType *origin)
itkGetInputMacro(ReferenceImage, ReferenceImageBaseType)
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void SetOutputSpacing(const SpacePrecisionType *spacing)
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
void LinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread)
itkSetInputMacro(ReferenceImage, ReferenceImageBaseType)
void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread)
~TransformToDisplacementFieldFilter() override=default
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:84
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
double SpacePrecisionType
Definition: itkFloatTypes.h:30
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:70