ITK  6.0.0
Insight Toolkit
itkPointSetToPointSetMetricWithIndexv4.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 itkPointSetToPointSetMetricWithIndexv4_h
19#define itkPointSetToPointSetMetricWithIndexv4_h
20
22
23#include "itkFixedArray.h"
24#include "itkPointsLocator.h"
25#include "itkPointSet.h"
26
27namespace itk
28{
69template <typename TFixedPointSet,
70 typename TMovingPointSet = TFixedPointSet,
71 class TInternalComputationValueType = double>
72class ITK_TEMPLATE_EXPORT PointSetToPointSetMetricWithIndexv4
73 : public ObjectToObjectMetric<TFixedPointSet::PointDimension,
74 TMovingPointSet::PointDimension,
75 Image<TInternalComputationValueType, TFixedPointSet::PointDimension>,
76 TInternalComputationValueType>
77{
78public:
79 ITK_DISALLOW_COPY_AND_MOVE(PointSetToPointSetMetricWithIndexv4);
80
83 using Superclass = ObjectToObjectMetric<TFixedPointSet::PointDimension,
84 TMovingPointSet::PointDimension,
86 TInternalComputationValueType>;
89
91 itkOverrideGetNameOfClassMacro(PointSetToPointSetMetricWithIndexv4);
92
94 using typename Superclass::MeasureType;
95
97 using typename Superclass::ParametersType;
98 using typename Superclass::ParametersValueType;
99 using typename Superclass::NumberOfParametersType;
100
102 using typename Superclass::DerivativeType;
103
105 using typename Superclass::FixedTransformType;
106 using typename Superclass::FixedTransformPointer;
107 using typename Superclass::FixedInputPointType;
108 using typename Superclass::FixedOutputPointType;
109 using typename Superclass::FixedTransformParametersType;
110
111 using typename Superclass::MovingTransformType;
112 using typename Superclass::MovingTransformPointer;
113 using typename Superclass::MovingInputPointType;
114 using typename Superclass::MovingOutputPointType;
115 using typename Superclass::MovingTransformParametersType;
116
117 using typename Superclass::JacobianType;
118 using typename Superclass::FixedTransformJacobianType;
119 using typename Superclass::MovingTransformJacobianType;
120
121 using DisplacementFieldTransformType = typename Superclass::MovingDisplacementFieldTransformType;
122
123 using ObjectType = typename Superclass::ObjectType;
124
126 using typename Superclass::DimensionType;
127
129 using FixedPointSetType = TFixedPointSet;
131 using FixedPixelType = typename TFixedPointSet::PixelType;
132 using FixedPointsContainer = typename TFixedPointSet::PointsContainer;
133
134 static constexpr DimensionType FixedPointDimension = Superclass::FixedDimension;
135
137 using MovingPointSetType = TMovingPointSet;
139 using MovingPixelType = typename TMovingPointSet::PixelType;
140 using MovingPointsContainer = typename TMovingPointSet::PointsContainer;
141
142 static constexpr DimensionType MovingPointDimension = Superclass::MovingDimension;
143
150 static constexpr DimensionType PointDimension = Superclass::FixedDimension;
151
154 using CoordinateType = typename PointType::CoordinateType;
155#ifndef ITK_FUTURE_LEGACY_REMOVE
156 using CoordRepType ITK_FUTURE_DEPRECATED(
157 "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
158#endif
160 using PointsConstIterator = typename PointsContainer::ConstIterator;
161 using PointIdentifier = typename PointsContainer::ElementIdentifier;
162
166
169
172
174 using VirtualImageType = typename Superclass::VirtualImageType;
175 using typename Superclass::VirtualImagePointer;
176 using typename Superclass::VirtualPixelType;
177 using typename Superclass::VirtualRegionType;
178 using typename Superclass::VirtualSizeType;
179 using typename Superclass::VirtualSpacingType;
180 using VirtualOriginType = typename Superclass::VirtualPointType;
181 using typename Superclass::VirtualPointType;
182 using typename Superclass::VirtualDirectionType;
183 using VirtualRadiusType = typename Superclass::VirtualSizeType;
184 using typename Superclass::VirtualIndexType;
185 using typename Superclass::VirtualPointSetType;
186 using typename Superclass::VirtualPointSetPointer;
187
189 void
190 SetFixedObject(const ObjectType * object) override
191 {
192 auto * pointSet = dynamic_cast<FixedPointSetType *>(const_cast<ObjectType *>(object));
193 if (pointSet != nullptr)
194 {
195 this->SetFixedPointSet(pointSet);
196 }
197 else
198 {
199 itkExceptionMacro("Incorrect object type. Should be a point set.");
200 }
201 }
205 void
206 SetMovingObject(const ObjectType * object) override
207 {
208 auto * pointSet = dynamic_cast<MovingPointSetType *>(const_cast<ObjectType *>(object));
209 if (pointSet != nullptr)
210 {
211 this->SetMovingPointSet(pointSet);
212 }
213 else
214 {
215 itkExceptionMacro("Incorrect object type. Should be a point set.");
216 }
217 }
221 itkSetConstObjectMacro(FixedPointSet, FixedPointSetType);
222 itkGetConstObjectMacro(FixedPointSet, FixedPointSetType);
226 itkGetModifiableObjectMacro(FixedTransformedPointSet, FixedTransformedPointSetType);
227
229 itkSetConstObjectMacro(MovingPointSet, MovingPointSetType);
230 itkGetConstObjectMacro(MovingPointSet, MovingPointSetType);
234 itkGetModifiableObjectMacro(MovingTransformedPointSet, MovingTransformedPointSetType);
235
241
253 GetValue() const override;
254
265 void
266 GetDerivative(DerivativeType &) const override;
267
278 void
280
285 const VirtualPointSetType *
287
292 void
293 Initialize() override;
294
295 bool
297 {
298 /* An arbitrary point in the virtual domain will not always
299 * correspond to a point within either point set. */
300 return false;
301 }
302
314 itkSetMacro(StoreDerivativeAsSparseFieldForLocalSupportTransforms, bool);
315 itkGetConstMacro(StoreDerivativeAsSparseFieldForLocalSupportTransforms, bool);
316 itkBooleanMacro(StoreDerivativeAsSparseFieldForLocalSupportTransforms);
322 itkSetMacro(CalculateValueAndDerivativeInTangentSpace, bool);
323 itkGetConstMacro(CalculateValueAndDerivativeInTangentSpace, bool);
324 itkBooleanMacro(CalculateValueAndDerivativeInTangentSpace);
327protected:
330 void
331 PrintSelf(std::ostream & os, Indent indent) const override;
332
333 typename FixedPointSetType::ConstPointer m_FixedPointSet{};
334 mutable typename FixedTransformedPointSetType::Pointer m_FixedTransformedPointSet{};
335
336 mutable typename PointsLocatorType::Pointer m_FixedTransformedPointsLocator{};
337
338 typename MovingPointSetType::ConstPointer m_MovingPointSet{};
339 mutable typename MovingTransformedPointSetType::Pointer m_MovingTransformedPointSet{};
340
341 mutable typename PointsLocatorType::Pointer m_MovingTransformedPointsLocator{};
342
344 mutable VirtualPointSetPointer m_VirtualTransformedPointSet{};
345
350 bool m_UsePointSetData{};
351
359 bool m_CalculateValueAndDerivativeInTangentSpace{};
360
363 virtual void
365
370 virtual void
372
378 virtual SizeValueType
380
383 void
384 CalculateValueAndDerivative(MeasureType & calculatedValue, DerivativeType & derivative, bool calculateValue) const;
385
392 void
394
401 void
403
408 void
410
415 void
417
418 using typename Superclass::MetricCategoryType;
419
422 GetMetricCategory() const override
423 {
424 return MetricCategoryType::POINT_SET_METRIC;
425 }
426
427 virtual bool
429 {
430 return true;
431 }
432
433 virtual bool
435 {
436 return true;
437 }
438
444 virtual MeasureType
445 GetLocalNeighborhoodValueWithIndex(const PointIdentifier &, const PointType &, const PixelType & pixel) const = 0;
446
453 virtual LocalDerivativeType
455
461 virtual void
463 const PointType &,
464 MeasureType &,
466 const PixelType & pixel) const = 0;
467
468private:
469 mutable bool m_MovingTransformPointLocatorsNeedInitialization{};
470 mutable bool m_FixedTransformPointLocatorsNeedInitialization{};
471
472 // Flag to keep track of whether a warning has already been issued
473 // regarding the number of valid points.
474 mutable bool m_HaveWarnedAboutNumberOfValidPoints{};
475
476 // Flag to store derivatives at fixed point locations with the rest being zero gradient
477 // (default = true).
478 bool m_StoreDerivativeAsSparseFieldForLocalSupportTransforms{};
479
480 mutable ModifiedTimeType m_MovingTransformedPointSetTime{};
481 mutable ModifiedTimeType m_FixedTransformedPointSetTime{};
482
483 // Create ranges over the point set for multithreaded computation of value and derivatives
484 using PointIdentifierPair = std::pair<PointIdentifier, PointIdentifier>;
485 using PointIdentifierRanges = std::vector<PointIdentifierPair>;
488};
489} // end namespace itk
490
491#ifndef ITK_MANUAL_INSTANTIATION
492# include "itkPointSetToPointSetMetricWithIndexv4.hxx"
493#endif
494
495#endif
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:54
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Computes similarity between regions of two objects.
const PointIdentifierRanges CreateRanges() const
typename PointsContainer::ElementIdentifier PointIdentifier
virtual void InitializeForIteration() const
virtual SizeValueType CalculateNumberOfValidFixedPoints() const
void GetDerivative(DerivativeType &) const override
~PointSetToPointSetMetricWithIndexv4() override=default
const VirtualPointSetType * GetVirtualTransformedPointSet() const
typename TMovingPointSet::PointsContainer MovingPointsContainer
virtual void GetLocalNeighborhoodValueAndDerivativeWithIndex(const PointIdentifier &, const PointType &, MeasureType &, LocalDerivativeType &, const PixelType &pixel) const =0
void CalculateValueAndDerivative(MeasureType &calculatedValue, DerivativeType &derivative, bool calculateValue) const
virtual LocalDerivativeType GetLocalNeighborhoodDerivativeWithIndex(const PointIdentifier &, const PointType &, const PixelType &pixel) const
MeasureType GetValue() const override
SizeValueType GetNumberOfComponents() const
void StorePointDerivative(const VirtualPointType &, const DerivativeType &, DerivativeType &) const
typename Superclass::MovingDisplacementFieldTransformType DisplacementFieldTransformType
std::pair< PointIdentifier, PointIdentifier > PointIdentifierPair
void GetValueAndDerivative(MeasureType &, DerivativeType &) const override
virtual MeasureType GetLocalNeighborhoodValueWithIndex(const PointIdentifier &, const PointType &, const PixelType &pixel) const =0
void PrintSelf(std::ostream &os, Indent indent) const override
typename PointsContainer::ConstIterator PointsConstIterator
typename TFixedPointSet::PointsContainer FixedPointsContainer
typename PointsLocatorType::NeighborsIdentifierType NeighborsIdentifierType
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
Accelerate geometric searches for points.
typename TreeType::InstanceIdentifierVectorType NeighborsIdentifierType
SmartPointer< const Self > ConstPointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:105
unsigned long SizeValueType
Definition: itkIntTypes.h:86