go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkAdvancedImageToImageMetric.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
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 * http://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 itkAdvancedImageToImageMetric_h
19#define itkAdvancedImageToImageMetric_h
20
21#include "itkImageToImageMetric.h"
22
23#include "itkImageSamplerBase.h"
24#include "itkGradientImageFilter.h"
25#include "itkBSplineInterpolateImageFunction.h"
29#include "itkFixedArray.h"
31
32#include "itkImageMaskSpatialObject.h"
33
34// Needed for checking for B-spline for faster implementation
37
38#include <cassert>
39#include <memory> // For unique_ptr.
40#include <typeinfo>
41
42namespace itk
43{
44
81template <class TFixedImage, class TMovingImage>
82class ITK_TEMPLATE_EXPORT AdvancedImageToImageMetric : public ImageToImageMetric<TFixedImage, TMovingImage>
84public:
85 ITK_DISALLOW_COPY_AND_MOVE(AdvancedImageToImageMetric);
89 using Superclass = ImageToImageMetric<TFixedImage, TMovingImage>;
90 using Pointer = SmartPointer<Self>;
91 using ConstPointer = SmartPointer<const Self>;
94 itkTypeMacro(AdvancedImageToImageMetric, ImageToImageMetric);
97 itkStaticConstMacro(MovingImageDimension, unsigned int, TMovingImage::ImageDimension);
98 itkStaticConstMacro(FixedImageDimension, unsigned int, TFixedImage::ImageDimension);
101 using typename Superclass::CoordinateRepresentationType;
102 using typename Superclass::MovingImageType;
103 using typename Superclass::MovingImagePixelType;
104 using MovingImagePointer = typename MovingImageType::Pointer;
105 using typename Superclass::MovingImageConstPointer;
106 using typename Superclass::FixedImageType;
107 using FixedImagePointer = typename FixedImageType::Pointer;
108 using typename Superclass::FixedImageConstPointer;
109 using typename Superclass::FixedImageRegionType;
110 using typename Superclass::TransformType;
111 using typename Superclass::TransformPointer;
112 using typename Superclass::InputPointType;
113 using typename Superclass::OutputPointType;
114 using typename Superclass::TransformParametersType;
115 using typename Superclass::TransformJacobianType;
116 using typename Superclass::InterpolatorType;
117 using typename Superclass::InterpolatorPointer;
118 using typename Superclass::RealType;
119 using typename Superclass::GradientPixelType;
120 using typename Superclass::GradientImageType;
121 using typename Superclass::GradientImagePointer;
123 // Overrule the mask type from its base class, ITK ImageToImageMetric.
124 using FixedImageMaskType = ImageMaskSpatialObject<Self::FixedImageDimension>;
125 using FixedImageMaskPointer = SmartPointer<FixedImageMaskType>;
126 using FixedImageMaskConstPointer = SmartPointer<const FixedImageMaskType>;
127 using MovingImageMaskType = ImageMaskSpatialObject<Self::MovingImageDimension>;
128 using MovingImageMaskPointer = SmartPointer<MovingImageMaskType>;
129 using MovingImageMaskConstPointer = SmartPointer<const MovingImageMaskType>;
131 using typename Superclass::MeasureType;
132 using typename Superclass::DerivativeType;
133 using DerivativeValueType = typename DerivativeType::ValueType;
134 using typename Superclass::ParametersType;
137 using FixedImagePixelType = typename FixedImageType::PixelType;
138 using MovingImageRegionType = typename MovingImageType::RegionType;
139 using MovingImageDerivativeScalesType = FixedArray<double, Self::MovingImageDimension>;
150 using FixedImageLimiterOutputType = typename FixedImageLimiterType::OutputType;
153 using MovingImageLimiterOutputType = typename MovingImageLimiterType::OutputType;
156 using ScalarType = typename TransformType::ScalarType;
158 using NumberOfParametersType = typename AdvancedTransformType::NumberOfParametersType;
170 using ThreadInfoType = MultiThreaderBase::WorkUnitInfo;
174 virtual void
175 SetFixedImageMask(const FixedImageMaskType * const arg)
177 assert(arg == nullptr || typeid(*arg) == typeid(FixedImageMaskType));
178 Superclass::SetFixedImageMask(arg);
181 virtual void
182 SetMovingImageMask(const MovingImageMaskType * const arg)
184 assert(arg == nullptr || typeid(*arg) == typeid(MovingImageMaskType));
185 Superclass::SetMovingImageMask(arg);
189 GetFixedImageMask() const override
191 const auto * const mask = Superclass::GetFixedImageMask();
192 assert(mask == nullptr || typeid(*mask) == typeid(FixedImageMaskType));
193 return static_cast<const FixedImageMaskType *>(mask);
194 }
195
196 const MovingImageMaskType *
197 GetMovingImageMask() const override
198 {
199 const auto * const mask = Superclass::GetMovingImageMask();
200 assert(mask == nullptr || typeid(*mask) == typeid(MovingImageMaskType));
201 return static_cast<const MovingImageMaskType *>(mask);
202 }
203
205 virtual void
207 {
208 this->Superclass::SetTransform(arg);
209 if (m_AdvancedTransform != arg)
210 {
211 m_AdvancedTransform = arg;
212 this->Modified();
213 }
214 }
215
216
218 AdvancedTransformType *
219 GetTransform() override
220 {
221 return m_AdvancedTransform.GetPointer();
222 }
223
225 const AdvancedTransformType *
226 GetTransform() const override
227 {
228 return m_AdvancedTransform.GetPointer();
231
233 itkSetObjectMacro(ImageSampler, ImageSamplerType);
237 return m_ImageSampler.GetPointer();
239
240
243 itkGetConstMacro(UseImageSampler, bool);
244
248 itkSetMacro(RequiredRatioOfValidSamples, double);
249 itkGetConstMacro(RequiredRatioOfValidSamples, double);
250
253 itkSetObjectMacro(MovingImageLimiter, MovingImageLimiterType);
254 itkGetConstObjectMacro(MovingImageLimiter, MovingImageLimiterType);
255 itkSetObjectMacro(FixedImageLimiter, FixedImageLimiterType);
256 itkGetConstObjectMacro(FixedImageLimiter, FixedImageLimiterType);
257
264 itkSetMacro(MovingLimitRangeRatio, double);
265 itkGetConstMacro(MovingLimitRangeRatio, double);
266 itkSetMacro(FixedLimitRangeRatio, double);
267 itkGetConstMacro(FixedLimitRangeRatio, double);
268
271 itkGetConstMacro(UseFixedImageLimiter, bool);
272 itkGetConstMacro(UseMovingImageLimiter, bool);
273
281 itkSetMacro(UseMovingImageDerivativeScales, bool);
282 itkGetConstMacro(UseMovingImageDerivativeScales, bool);
284 itkSetMacro(ScaleGradientWithRespectToMovingImageOrientation, bool);
285 itkGetConstMacro(ScaleGradientWithRespectToMovingImageOrientation, bool);
286
287 itkSetMacro(MovingImageDerivativeScales, MovingImageDerivativeScalesType);
288 itkGetConstReferenceMacro(MovingImageDerivativeScales, MovingImageDerivativeScalesType);
289
298 void
299 Initialize() override;
300
302 itkSetMacro(UseMetricSingleThreaded, bool);
303 itkGetConstReferenceMacro(UseMetricSingleThreaded, bool);
304 itkBooleanMacro(UseMetricSingleThreaded);
305
307 // \todo: maybe these can be united, check base class.
308 itkSetMacro(UseMultiThread, bool);
309 itkGetConstReferenceMacro(UseMultiThread, bool);
310 itkBooleanMacro(UseMultiThread);
311
317 virtual void
318 BeforeThreadedGetValueAndDerivative(const TransformParametersType & parameters) const;
319
320protected:
323
325 ~AdvancedImageToImageMetric() override = default;
326
328 void
329 PrintSelf(std::ostream & os, Indent indent) const override;
330
334 using FixedImageIndexType = typename FixedImageType::IndexType;
335 using FixedImageIndexValueType = typename FixedImageIndexType::IndexValueType;
336 using MovingImageIndexType = typename MovingImageType::IndexType;
337 using FixedImagePointType = typename TransformType::InputPointType;
338 using MovingImagePointType = typename TransformType::OutputPointType;
339 using MovingImageContinuousIndexType = typename InterpolatorType::ContinuousIndexType;
340
343 BSplineInterpolateImageFunction<MovingImageType, CoordinateRepresentationType, double>;
344 using BSplineInterpolatorPointer = typename BSplineInterpolatorType::Pointer;
346 BSplineInterpolateImageFunction<MovingImageType, CoordinateRepresentationType, float>;
347 using BSplineInterpolatorFloatPointer = typename BSplineInterpolatorFloatType::Pointer;
353 using MovingImageDerivativeType = typename BSplineInterpolatorType::CovariantVectorType;
354
357
363 mutable ImageSamplerPointer m_ImageSampler{ nullptr };
366 typename AdvancedTransformType::Pointer m_AdvancedTransform{ nullptr };
369 mutable bool m_TransformIsBSpline{ false };
372 FixedImagePixelType m_FixedImageTrueMin{ 0 };
373 FixedImagePixelType m_FixedImageTrueMax{ 1 };
374 MovingImagePixelType m_MovingImageTrueMin{ 0 };
375 MovingImagePixelType m_MovingImageTrueMax{ 1 };
376 FixedImageLimiterOutputType m_FixedImageMinLimit{ 0 };
377 FixedImageLimiterOutputType m_FixedImageMaxLimit{ 1 };
378 MovingImageLimiterOutputType m_MovingImageMinLimit{ 0 };
379 MovingImageLimiterOutputType m_MovingImageMaxLimit{ 1 };
380
384 virtual void ThreadedGetValue(ThreadIdType) const {}
385
387 virtual void
388 AfterThreadedGetValue(MeasureType &) const
389 {}
390
392 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
394
396 void
398
400 virtual void ThreadedGetValueAndDerivative(ThreadIdType) const {}
401
403 virtual void
404 AfterThreadedGetValueAndDerivative(MeasureType &, DerivativeType &) const
405 {}
406
408 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
412 void
414
416 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
418
420 bool m_UseMetricSingleThreaded{ true };
421 bool m_UseMultiThread{ false };
422
427 {
428 // To give the threads access to all members.
430 // Used for accumulating derivatives
433 };
434 mutable MultiThreaderParameterType m_ThreaderMetricParameters{};
435
445 // test per thread struct with padding and alignment
447 {
449 MeasureType st_Value;
450 DerivativeType st_Derivative;
451 };
452 itkPadStruct(ITK_CACHE_LINE_ALIGNMENT,
454 PaddedGetValueAndDerivativePerThreadStruct);
455 itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT,
456 PaddedGetValueAndDerivativePerThreadStruct,
457 AlignedGetValueAndDerivativePerThreadStruct);
458 mutable std::unique_ptr<AlignedGetValueAndDerivativePerThreadStruct[]> m_GetValueAndDerivativePerThreadVariables{
459 nullptr
460 };
461 mutable ThreadIdType m_GetValueAndDerivativePerThreadVariablesSize{ 0 };
462
464 virtual void
466
472 virtual void
474
477 itkSetMacro(UseImageSampler, bool);
478
481 void
482 CheckNumberOfSamples(unsigned long wanted, unsigned long found) const;
483
488 void
490
499 virtual bool
501 RealType & movingImageValue,
502 MovingImageDerivativeType * gradient) const
503 {
504 return EvaluateMovingImageValueAndDerivativeWithOptionalThreadId(mappedPoint, movingImageValue, gradient);
505 }
506
507 /* A faster version of `EvaluateMovingImageValueAndDerivative`: Non-virtual, using multithreading, and doing less
508 * dynamic memory allocation/decallocation operations, internally. */
509 bool
511 RealType & movingImageValue,
512 MovingImageDerivativeType * gradient,
513 const ThreadIdType threadId) const
514 {
515 return EvaluateMovingImageValueAndDerivativeWithOptionalThreadId(mappedPoint, movingImageValue, gradient, threadId);
516 }
517
522 virtual void
523 EvaluateTransformJacobianInnerProduct(const TransformJacobianType & jacobian,
524 const MovingImageDerivativeType & movingImageDerivative,
525 DerivativeType & imageJacobian) const;
526
533 void
535
537 void
539
542 TransformPoint(const FixedImagePointType & fixedImagePoint) const;
543
550 bool
552 TransformJacobianType & jacobian,
553 NonZeroJacobianIndicesType & nzji) const;
554
556 virtual bool
558
561 void
563
566 itkSetMacro(UseFixedImageLimiter, bool);
567 itkSetMacro(UseMovingImageLimiter, bool);
568
569 double m_FixedLimitRangeRatio{ 0.01 };
570 double m_MovingLimitRangeRatio{ 0.01 };
571
572 // Prevent accidentally calling SetFixedImageMask or SetMovingImageMask through the ITK ImageToImageMetric interface.
573 void
574 SetFixedImageMask(typename Superclass::FixedImageMaskType *) final
575 {
576 itkExceptionMacro("Intentionally left unimplemented!");
577 }
578
579 void
580 SetFixedImageMask(const typename Superclass::FixedImageMaskType *) final
581 {
582 itkExceptionMacro("Intentionally left unimplemented!");
583 }
584
585 void
586 SetMovingImageMask(typename Superclass::MovingImageMaskType *) final
587 {
588 itkExceptionMacro("Intentionally left unimplemented!");
589 }
590
591 void
592 SetMovingImageMask(const typename Superclass::MovingImageMaskType *) final
593 {
594 itkExceptionMacro("Intentionally left unimplemented!");
595 }
596
597 // Protected using-declaration, to avoid `-Woverloaded-virtual` warnings from GCC (GCC 11.4) or clang (macos-12).
598 using Superclass::SetTransform;
599
600private:
601 template <typename... TOptionalThreadId>
602 bool
604 RealType & movingImageValue,
605 MovingImageDerivativeType * gradient,
606 const TOptionalThreadId... optionalThreadId) const;
607
609 FixedImageLimiterPointer m_FixedImageLimiter{ nullptr };
610 MovingImageLimiterPointer m_MovingImageLimiter{ nullptr };
611 LinearInterpolatorPointer m_LinearInterpolator{ nullptr };
612 BSplineInterpolatorPointer m_BSplineInterpolator{ nullptr };
613 BSplineInterpolatorFloatPointer m_BSplineInterpolatorFloat{ nullptr };
614 ReducedBSplineInterpolatorPointer m_ReducedBSplineInterpolator{ nullptr };
615
617 bool m_UseImageSampler{ false };
618 bool m_UseFixedImageLimiter{ false };
619 bool m_UseMovingImageLimiter{ false };
620 double m_RequiredRatioOfValidSamples{ 0.25 };
621 bool m_UseMovingImageDerivativeScales{ false };
622 bool m_ScaleGradientWithRespectToMovingImageOrientation{ false };
623
624 MovingImageDerivativeScalesType m_MovingImageDerivativeScales{ MovingImageDerivativeScalesType::Filled(1.0) };
625
626 // Private using-declarations, to avoid `-Woverloaded-virtual` warnings from GCC (GCC 11.4) or clang (macos-12).
627 using Superclass::TransformPoint;
628
630 {};
631
632 // Prevent accidentally accessing m_FixedImageMask or m_MovingImageMask directly from the ITK's ImageToImageMetric.
633 static constexpr DummyMask m_FixedImageMask{};
634 static constexpr DummyMask m_MovingImageMask{};
635};
636
637} // end namespace itk
638
639#ifndef ITK_MANUAL_INSTANTIATION
640# include "itkAdvancedImageToImageMetric.hxx"
641#endif
642
643#endif // end #ifndef itkAdvancedImageToImageMetric_h
Deformable transform using a B-spline representation.
This class combines two transforms: an 'initial transform' with a 'current transform'.
An extension of the ITK ImageToImageMetric. It is the intended base class for all elastix metrics.
virtual void AfterThreadedGetValue(MeasureType &) const
virtual bool IsInsideMovingMask(const MovingImagePointType &point) const
virtual void ThreadedGetValueAndDerivative(ThreadIdType) const
const MovingImageMaskType * GetMovingImageMask() const override
typename TransformType::OutputPointType MovingImagePointType
void PrintSelf(std::ostream &os, Indent indent) const override
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
virtual bool EvaluateMovingImageValueAndDerivative(const MovingImagePointType &mappedPoint, RealType &movingImageValue, MovingImageDerivativeType *gradient) const
typename FixedImageType::PixelType FixedImagePixelType
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, float > BSplineInterpolatorFloatType
typename DerivativeType::ValueType DerivativeValueType
typename FixedImageLimiterType::Pointer FixedImageLimiterPointer
void CheckNumberOfSamples(unsigned long wanted, unsigned long found) const
virtual void InitializeImageSampler()
void SetFixedImageMask(const typename Superclass::FixedImageMaskType *) final
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename MovingImageType::IndexType MovingImageIndexType
virtual void SetTransform(AdvancedTransformType *arg)
void LaunchGetValueAndDerivativeThreaderCallback() const
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION GetValueAndDerivativeThreaderCallback(void *arg)
bool FastEvaluateMovingImageValueAndDerivative(const MovingImagePointType &mappedPoint, RealType &movingImageValue, MovingImageDerivativeType *gradient, const ThreadIdType threadId) const
void SetFixedImageMask(typename Superclass::FixedImageMaskType *) final
typename BSplineOrder3TransformType::Pointer BSplineOrder3TransformPointer
void LaunchGetValueThreaderCallback() const
virtual void InitializeThreadingParameters() const
~AdvancedImageToImageMetric() override=default
typename ReducedBSplineInterpolatorType::Pointer ReducedBSplineInterpolatorPointer
virtual void BeforeThreadedGetValueAndDerivative(const TransformParametersType &parameters) const
typename BSplineOrder2TransformType::Pointer BSplineOrder2TransformPointer
typename FixedImageType::IndexType FixedImageIndexType
virtual void ThreadedGetValue(ThreadIdType) const
itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedGetValueAndDerivativePerThreadStruct, AlignedGetValueAndDerivativePerThreadStruct)
typename BSplineInterpolatorType::Pointer BSplineInterpolatorPointer
bool EvaluateMovingImageValueAndDerivativeWithOptionalThreadId(const MovingImagePointType &mappedPoint, RealType &movingImageValue, MovingImageDerivativeType *gradient, const TOptionalThreadId... optionalThreadId) const
typename BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
virtual void AfterThreadedGetValueAndDerivative(MeasureType &, DerivativeType &) const
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION AccumulateDerivativesThreaderCallback(void *arg)
virtual void EvaluateTransformJacobianInnerProduct(const TransformJacobianType &jacobian, const MovingImageDerivativeType &movingImageDerivative, DerivativeType &imageJacobian) const
const FixedImageMaskType * GetFixedImageMask() const override
MovingImagePointType TransformPoint(const FixedImagePointType &fixedImagePoint) const
typename MovingImageLimiterType::OutputType MovingImageLimiterOutputType
typename TransformType::InputPointType FixedImagePointType
const AdvancedTransformType * GetTransform() const override
void SetMovingImageMask(const typename Superclass::MovingImageMaskType *) final
typename FixedImageLimiterType::OutputType FixedImageLimiterOutputType
void SetMovingImageMask(typename Superclass::MovingImageMaskType *) final
typename BSplineInterpolatorFloatType::Pointer BSplineInterpolatorFloatPointer
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, double > BSplineInterpolatorType
AdvancedTransformType * GetTransform() override
typename TransformType::ScalarType ScalarType
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION GetValueThreaderCallback(void *arg)
SmartPointer< const FixedImageMaskType > FixedImageMaskConstPointer
typename ImageSamplerType::Pointer ImageSamplerPointer
bool EvaluateTransformJacobian(const FixedImagePointType &fixedImagePoint, TransformJacobianType &jacobian, NonZeroJacobianIndicesType &nzji) const
typename InterpolatorType::ContinuousIndexType MovingImageContinuousIndexType
typename AdvancedTransformType::NumberOfParametersType NumberOfParametersType
itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, GetValueAndDerivativePerThreadStruct, PaddedGetValueAndDerivativePerThreadStruct)
ImageMaskSpatialObject< Self::MovingImageDimension > MovingImageMaskType
typename LinearInterpolatorType::Pointer LinearInterpolatorPointer
Linearly interpolate an image at specified positions.
Transform maps points, vectors and covariant vectors from an input space to an output space.
SmartPointer< Self > Pointer
std::vector< unsigned long > NonZeroJacobianIndicesType
This class is a base class for any image sampler.
Base class for all ITK limiter function objects.
ImageMaskSpatialObject< Self::FixedImageDimension > FixedImageMaskType
Evaluates the B-Spline interpolation of an image. Spline order may be from 0 to 5.
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.


Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo