ITK  6.0.0
Insight Toolkit
itkLevelSetMotionRegistrationFunction.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 itkLevelSetMotionRegistrationFunction_h
19#define itkLevelSetMotionRegistrationFunction_h
20
22#include "itkPoint.h"
23#include "itkCovariantVector.h"
26#include <mutex>
27
28namespace itk
29{
54template <typename TFixedImage, typename TMovingImage, typename TDisplacementField>
55class ITK_TEMPLATE_EXPORT LevelSetMotionRegistrationFunction
56 : public PDEDeformableRegistrationFunction<TFixedImage, TMovingImage, TDisplacementField>
57{
58public:
59 ITK_DISALLOW_COPY_AND_MOVE(LevelSetMotionRegistrationFunction);
60
64
67
69 itkNewMacro(Self);
70
72 itkOverrideGetNameOfClassMacro(LevelSetMotionRegistrationFunction);
73
75 using typename Superclass::MovingImageType;
76 using typename Superclass::MovingImagePointer;
77 using MovingSpacingType = typename MovingImageType::SpacingType;
78
80 using typename Superclass::FixedImageType;
81 using typename Superclass::FixedImagePointer;
84 using SpacingType = typename FixedImageType::SpacingType;
85
87 using typename Superclass::DisplacementFieldType;
88 using typename Superclass::DisplacementFieldTypePointer;
89
91 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
92
94 using typename Superclass::PixelType;
95 using typename Superclass::RadiusType;
96 using typename Superclass::NeighborhoodType;
97 using typename Superclass::FloatOffsetType;
98 using typename Superclass::TimeStepType;
99
101 using CoordinateType = double;
102#ifndef ITK_FUTURE_LEGACY_REMOVE
103 using CoordRepType ITK_FUTURE_DEPRECATED(
104 "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
105#endif
110
114
118
120 void
122 {
123 m_MovingImageInterpolator = ptr;
124 }
125
127 InterpolatorType *
129 {
130 return m_MovingImageInterpolator;
131 }
132
137 TimeStepType
138 ComputeGlobalTimeStep(void * GlobalData) const override;
139
142 void *
143 GetGlobalDataPointer() const override
144 {
145 auto * global = new GlobalDataStruct();
146
147 global->m_SumOfSquaredDifference = 0.0;
148 global->m_NumberOfPixelsProcessed = 0L;
149 global->m_SumOfSquaredChange = 0;
150 global->m_MaxL1Norm = NumericTraits<double>::NonpositiveMin();
151 return global;
152 }
153
158 void
159 ReleaseGlobalDataPointer(void * gd) const override;
160
162 void
164
170 ComputeUpdate(const NeighborhoodType & it, void * gd, const FloatOffsetType & offset = FloatOffsetType(0.0)) override;
171
175 virtual double
176 GetMetric() const
177 {
178 return m_Metric;
179 }
180
182 virtual double
184 {
185 return m_RMSChange;
186 }
187
194 virtual void
195 SetAlpha(double);
196
197 virtual double
198 GetAlpha() const;
199
204 virtual void
206
207 virtual double
209
212 virtual void
214
215 virtual double
217
220 virtual void
222
223 virtual double
225
229 void
231
232 bool
234
235protected:
238 void
239 PrintSelf(std::ostream & os, Indent indent) const override;
240
243
247 {
252 };
253
254private:
256 SpacingType m_FixedImageSpacing{};
257 PointType m_FixedImageOrigin{};
258
260 MovingImageSmoothingFilterPointer m_MovingImageSmoothingFilter{};
261
263 InterpolatorPointer m_MovingImageInterpolator{};
264 InterpolatorPointer m_SmoothMovingImageInterpolator{};
265
268 double m_Alpha{};
269
271 double m_GradientMagnitudeThreshold{};
272
274 double m_IntensityDifferenceThreshold{};
275
277 double m_GradientSmoothingStandardDeviations{};
278
282 mutable double m_Metric{};
283 mutable double m_SumOfSquaredDifference{};
284 mutable SizeValueType m_NumberOfPixelsProcessed{};
285 mutable double m_RMSChange{};
286 mutable double m_SumOfSquaredChange{};
287
289 mutable std::mutex m_MetricCalculationMutex{};
290
291 bool m_UseImageSpacing{ true };
292};
293} // end namespace itk
294
295#ifndef ITK_MANUAL_INSTANTIATION
296# include "itkLevelSetMotionRegistrationFunction.hxx"
297#endif
298
299#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
A templated class holding a n-Dimensional covariant vector.
typename ConstNeighborhoodIterator< TDisplacementField >::RadiusType RadiusType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for all image interpolators.
PixelType ComputeUpdate(const NeighborhoodType &it, void *gd, const FloatOffsetType &offset=FloatOffsetType(0.0)) override
void PrintSelf(std::ostream &os, Indent indent) const override
typename MovingImageSmoothingFilterType::Pointer MovingImageSmoothingFilterPointer
typename MovingImageType::SpacingType MovingSpacingType
void ReleaseGlobalDataPointer(void *gd) const override
virtual double GetGradientMagnitudeThreshold() const
virtual double GetGradientSmoothingStandardDeviations() const
virtual void SetIntensityDifferenceThreshold(double)
virtual void SetGradientSmoothingStandardDeviations(double)
virtual void SetGradientMagnitudeThreshold(double)
virtual double GetIntensityDifferenceThreshold() const
~LevelSetMotionRegistrationFunction() override=default
TimeStepType ComputeGlobalTimeStep(void *GlobalData) const override
Light weight base class for most itk classes.
Linearly interpolate an image at specified positions.
static constexpr T NonpositiveMin()
Computes the smoothing of an image by convolution with the Gaussian kernels implemented as IIR filter...
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:63
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86