ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkDemonsRegistrationFunction.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 itkDemonsRegistrationFunction_h
19#define itkDemonsRegistrationFunction_h
20
22#include "itkPoint.h"
25#include <mutex>
26
27namespace itk
28{
53template <typename TFixedImage, typename TMovingImage, typename TDisplacementField>
54class ITK_TEMPLATE_EXPORT DemonsRegistrationFunction
55 : public PDEDeformableRegistrationFunction<TFixedImage, TMovingImage, TDisplacementField>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(DemonsRegistrationFunction);
59
65
67 itkNewMacro(Self);
68
70 itkOverrideGetNameOfClassMacro(DemonsRegistrationFunction);
71
73 using typename Superclass::MovingImageType;
74 using typename Superclass::MovingImagePointer;
75
77 using typename Superclass::FixedImageType;
78 using typename Superclass::FixedImagePointer;
79 using IndexType = typename FixedImageType::IndexType;
80 using SizeType = typename FixedImageType::SizeType;
81 using SpacingType = typename FixedImageType::SpacingType;
82
86
88 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
89
91 using typename Superclass::PixelType;
92 using typename Superclass::RadiusType;
93 using typename Superclass::NeighborhoodType;
94 using typename Superclass::FloatOffsetType;
95 using typename Superclass::TimeStepType;
96
98 using CoordinateType = double;
99#ifndef ITK_FUTURE_LEGACY_REMOVE
100 using CoordRepType ITK_FUTURE_DEPRECATED(
101 "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
102#endif
107
110
114
118
120 void
125
127 InterpolatorType *
132
134 TimeStepType
135 ComputeGlobalTimeStep(void * itkNotUsed(GlobalData)) const override
136 {
137 return m_TimeStep;
138 }
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 return global;
151 }
152
157 void
158 ReleaseGlobalDataPointer(void * gd) const override;
159
161 void
163
169 ComputeUpdate(const NeighborhoodType & it, void * gd, const FloatOffsetType & offset = FloatOffsetType(0.0)) override;
170
174 virtual double
175 GetMetric() const
176 {
177 return m_Metric;
178 }
179
181 virtual double
183 {
184 return m_RMSChange;
185 }
186
190 virtual void
192 {
194 }
195 virtual bool
200
201
206 virtual void
208
209 virtual double
211
212protected:
214 ~DemonsRegistrationFunction() override = default;
215 void
216 PrintSelf(std::ostream & os, Indent indent) const override;
217
220
229
230private:
232 // SpacingType m_FixedImageSpacing;
233 // PointType m_FixedImageOrigin;
235 double m_Normalizer{};
236
239
243
246
249
252
255
259 mutable double m_Metric{};
260 mutable double m_SumOfSquaredDifference{};
262 mutable double m_RMSChange{};
263 mutable double m_SumOfSquaredChange{};
264
266 mutable std::mutex m_MetricCalculationMutex{};
267};
268} // end namespace itk
269
270#ifndef ITK_MANUAL_INSTANTIATION
271# include "itkDemonsRegistrationFunction.hxx"
272#endif
273
274#endif
Calculate the derivative by central differencing.
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
A templated class holding a n-Dimensional covariant vector.
~DemonsRegistrationFunction() override=default
void ReleaseGlobalDataPointer(void *gd) const override
LinearInterpolateImageFunction< MovingImageType, CoordinateType > DefaultInterpolatorType
void InitializeIteration() override
CentralDifferenceImageFunction< MovingImageType, CoordinateType > MovingImageGradientCalculatorType
PixelType ComputeUpdate(const NeighborhoodType &it, void *gd, const FloatOffsetType &offset=FloatOffsetType(0.0)) override
void PrintSelf(std::ostream &os, Indent indent) const override
void SetMovingImageInterpolator(InterpolatorType *ptr)
PDEDeformableRegistrationFunction< FixedImageType, MovingImageType, DisplacementFieldType > Superclass
TimeStepType ComputeGlobalTimeStep(void *GlobalData) const override
virtual double GetIntensityDifferenceThreshold() const
virtual void SetIntensityDifferenceThreshold(double)
ConstNeighborhoodIterator< TImageType, DefaultBoundaryConditionType > NeighborhoodType
static constexpr unsigned int ImageDimension
Vector< float, Self::ImageDimension > FloatOffsetType
typename ImageType::PixelType PixelType
typename ConstNeighborhoodIterator< TImageType >::RadiusType RadiusType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Base class for all image interpolators.
Point< CoordinateType, Self::ImageDimension > PointType
Linearly interpolate an image at specified positions.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86