ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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;
67 using OutputImageRegionType = typename OutputImageType::RegionType;
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;
85 using RegionType = typename OutputImageType::RegionType;
86 using SizeType = typename RegionType::SizeType;
87 using IndexType = typename OutputImageType::IndexType;
88 using PointType = typename OutputImageType::PointType;
89 using SpacingType = typename OutputImageType::SpacingType;
90 using OriginType = typename OutputImageType::PointType;
91 using DirectionType = typename OutputImageType::DirectionType;
92
95
101 virtual void
103 const TransformInputType *
104 GetInput() const;
105 itkSetGetDecoratedObjectInputMacro(Transform, TransformType);
107
110 itkSetMacro(OutputStartIndex, IndexType);
111 itkGetConstReferenceMacro(OutputStartIndex, IndexType);
113
115 itkSetMacro(Size, SizeType);
116 itkGetConstReferenceMacro(Size, SizeType);
118
120 itkSetMacro(OutputSpacing, SpacingType);
121 virtual void
124
126 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
127
129 itkSetMacro(OutputOrigin, OriginType);
130 virtual void
133
135 itkGetConstReferenceMacro(OutputOrigin, OriginType);
136
138 itkSetMacro(OutputDirection, DirectionType);
139 itkGetConstReferenceMacro(OutputDirection, DirectionType);
141
149
152
155 itkSetMacro(UseReferenceImage, bool);
156 itkBooleanMacro(UseReferenceImage);
157 itkGetConstMacro(UseReferenceImage, bool);
159
160 static constexpr unsigned int PixelDimension = PixelType::Dimension;
162
163protected:
166
168 void
170
172 void
173 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
174
175
179 void
181
185 void
186 LinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread);
187
188 void
189 PrintSelf(std::ostream & os, Indent indent) const override;
190
191private:
193 SizeType m_Size{}; // size of the output region
194 IndexType m_OutputStartIndex{}; // start index of the output region
195 SpacingType m_OutputSpacing{}; // output image spacing
196 OriginType m_OutputOrigin{}; // output image origin
197 DirectionType m_OutputDirection{}; // output image direction cosines
198 bool m_UseReferenceImage{ false };
199};
200} // end namespace itk
201
202#ifndef ITK_MANUAL_INSTANTIATION
203# include "itkTransformToDisplacementFieldFilter.hxx"
204#endif
205
206#endif
Decorates any subclass of itkObject with a DataObject API.
Base class for templated image classes.
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implements transparent reference counting.
virtual void SetInput(const TransformInputType *input)
typename OutputImageType::DirectionType DirectionType
const TransformInputType * GetInput() const
virtual void SetOutputOrigin(const SpacePrecisionType *origin)
Transform< TParametersValueType, ImageDimension, ImageDimension > TransformType
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.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
double SpacePrecisionType
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70