ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkInvertDisplacementFieldImageFilter.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 itkInvertDisplacementFieldImageFilter_h
19#define itkInvertDisplacementFieldImageFilter_h
20
24#include <mutex>
25
26namespace itk
27{
28
103
104template <typename TInputImage, typename TOutputImage = TInputImage>
105class ITK_TEMPLATE_EXPORT InvertDisplacementFieldImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
106{
107public:
108 ITK_DISALLOW_COPY_AND_MOVE(InvertDisplacementFieldImageFilter);
109
114
116 itkNewMacro(Self);
117
119 itkOverrideGetNameOfClassMacro(InvertDisplacementFieldImageFilter);
120
122 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
123
124 using InputFieldType = TInputImage;
125 using OutputFieldType = TOutputImage;
126
129
131 using PixelType = typename OutputFieldType::PixelType;
132 using VectorType = typename OutputFieldType::PixelType;
133 using RegionType = typename OutputFieldType::RegionType;
135 using IndexType = typename OutputFieldType::IndexType;
136
137 using PointType = typename OutputFieldType::PointType;
138 using SpacingType = typename OutputFieldType::SpacingType;
139 using OriginType = typename OutputFieldType::PointType;
140 using SizeType = typename OutputFieldType::SizeType;
141 using DirectionType = typename OutputFieldType::DirectionType;
142
144 using RealType = typename VectorType::ComponentType;
148
150 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
151
153 void
155 {
156 itkDebugMacro("setting deformation field to " << field);
157 if (field != this->GetInput(0))
158 {
159 this->SetInput(0, field);
160 this->Modified();
161 if (!this->m_Interpolator.IsNull())
162 {
163 this->m_Interpolator->SetInputImage(field);
164 }
165 }
166 }
167
171 const InputFieldType *
173 {
174 return this->GetInput(0);
175 }
176
179 itkSetInputMacro(InverseFieldInitialEstimate, InverseDisplacementFieldType);
180 itkGetInputMacro(InverseFieldInitialEstimate, InverseDisplacementFieldType);
182
183 /* Set the interpolator. */
184 virtual void
186
187 /* Set/Get the number of iterations */
188 itkSetMacro(MaximumNumberOfIterations, unsigned int);
189 itkGetConstMacro(MaximumNumberOfIterations, unsigned int);
190
191 /* Set/Get the mean stopping criterion */
192 itkSetMacro(MeanErrorToleranceThreshold, RealType);
193 itkGetConstMacro(MeanErrorToleranceThreshold, RealType);
194
195 /* Set/Get the max stopping criterion */
196 itkSetMacro(MaxErrorToleranceThreshold, RealType);
197 itkGetConstMacro(MaxErrorToleranceThreshold, RealType);
198
199 /* Get the max norm */
200 itkGetConstMacro(MaxErrorNorm, RealType);
201
202 /* Get the mean norm */
203 itkGetConstMacro(MeanErrorNorm, RealType);
204
205 /* Should we force the boundary to have zero displacement? */
206 itkSetMacro(EnforceBoundaryCondition, bool);
207 itkGetMacro(EnforceBoundaryCondition, bool);
208 itkBooleanMacro(EnforceBoundaryCondition);
209
210protected:
213
216
218 void
219 PrintSelf(std::ostream & os, Indent indent) const override;
220
222 void
223 GenerateData() override;
224
226 void
228
229private:
232
233 unsigned int m_MaximumNumberOfIterations{ 20 };
234
237
238 // internal ivars necessary for multithreading basic operations
239
240 typename DisplacementFieldType::Pointer m_ComposedField{};
242
249 std::mutex m_Mutex{};
250};
251
252} // end namespace itk
253
254#ifndef ITK_MANUAL_INSTANTIATION
255# include "itkInvertDisplacementFieldImageFilter.hxx"
256#endif
257
258#endif
virtual void SetInput(const InputImageType *input)
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
void PrintSelf(std::ostream &os, Indent indent) const override
itkSetInputMacro(InverseFieldInitialEstimate, InverseDisplacementFieldType)
virtual void SetInterpolator(InterpolatorType *interpolator)
VectorLinearInterpolateImageFunction< InputFieldType, RealType > DefaultInterpolatorType
itkGetInputMacro(InverseFieldInitialEstimate, InverseDisplacementFieldType)
~InvertDisplacementFieldImageFilter() override=default
VectorInterpolateImageFunction< InputFieldType, RealType > InterpolatorType
void DynamicThreadedGenerateData(const RegionType &) override
ImageToImageFilter< TInputImage, TOutputImage > Superclass
virtual void Modified() const
Implements transparent reference counting.
Base class for all vector image interpolators.
Linearly interpolate a vector image at specified positions.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....