ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkMattesMutualInformationImageToImageMetric.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 itkMattesMutualInformationImageToImageMetric_h
19#define itkMattesMutualInformationImageToImageMetric_h
20
22#include "itkPoint.h"
23#include "itkIndex.h"
25#include "itkArray2D.h"
26
27#include <memory> // For unique_ptr.
28#include <mutex>
29
30
31namespace itk
32{
103template <typename TFixedImage, typename TMovingImage>
105 : public ImageToImageMetric<TFixedImage, TMovingImage>
106{
107public:
108 ITK_DISALLOW_COPY_AND_MOVE(MattesMutualInformationImageToImageMetric);
109
115
117 itkNewMacro(Self);
118
120 itkOverrideGetNameOfClassMacro(MattesMutualInformationImageToImageMetric);
121
123 using typename Superclass::TransformType;
124 using typename Superclass::TransformPointer;
126 using typename Superclass::InterpolatorType;
127 using typename Superclass::MeasureType;
128 using typename Superclass::DerivativeType;
129 using typename Superclass::ParametersType;
130 using typename Superclass::FixedImageType;
131 using typename Superclass::MovingImageType;
137
141 using typename Superclass::WeightsValueType;
142 using typename Superclass::IndexValueType;
143
144 using OffsetValueType = typename FixedImageType::OffsetValueType;
145
147 static constexpr unsigned int MovingImageDimension = MovingImageType::ImageDimension;
148
156 void
157 Initialize() override;
158
161 GetValue(const ParametersType & parameters) const override;
162
164 void
165 GetDerivative(const ParametersType & parameters, DerivativeType & derivative) const override;
166
168 void
170 MeasureType & value,
171 DerivativeType & derivative) const override;
172
179 itkSetClampMacro(NumberOfHistogramBins, SizeValueType, 5, NumericTraits<SizeValueType>::max());
180 itkGetConstReferenceMacro(NumberOfHistogramBins, SizeValueType);
182
207 itkSetMacro(UseExplicitPDFDerivatives, bool);
208 itkGetConstReferenceMacro(UseExplicitPDFDerivatives, bool);
209 itkBooleanMacro(UseExplicitPDFDerivatives);
211
213 using PDFValueType = double; // NOTE: floating point precision is not as stable. Double precision proves faster and
214 // more robust in real-world testing.
215
219
224 const typename JointPDFType::Pointer
226 {
227 if (this->m_MMIMetricPerThreadVariables == nullptr)
228 {
229 return JointPDFType::Pointer(nullptr);
230 }
231 return this->m_MMIMetricPerThreadVariables[0].JointPDF;
232 }
233
234
241 const typename JointPDFDerivativesType::Pointer
243 {
244 if (this->m_MMIMetricPerThreadVariables == nullptr)
245 {
246 return JointPDFDerivativesType::Pointer(nullptr);
247 }
248 return this->m_MMIMetricPerThreadVariables[0].JointPDFDerivatives;
249 }
250
251
252protected:
255 void
256 PrintSelf(std::ostream & os, Indent indent) const override;
257
258private:
267
271
273 void
275
277 void
279
281 void
283 unsigned int sampleNumber,
284 int pdfMovingIndex,
285 const ImageDerivativesType & movingImageGradientValue,
286 PDFValueType cubicBSplineDerivativeValue) const;
287
288 void
289 GetValueThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override;
290 void
291 GetValueThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override;
292 // NOTE: The signature in base class requires that movingImageValue is of type double
293 bool
295 SizeValueType fixedImageSample,
296 const MovingImagePointType & mappedPoint,
297 double movingImageValue) const override;
298
299 void
300 GetValueAndDerivativeThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override;
301 void
302 GetValueAndDerivativeThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override;
303 // NOTE: The signature in base class requires that movingImageValue is of type double
304 bool
306 SizeValueType fixedImageSample,
307 const MovingImagePointType & mappedPoint,
308 double movingImageValue,
309 const ImageDerivativesType & movingImageGradientValue) const override;
310
321
325
327
329 using MarginalPDFType = std::vector<PDFValueType>;
331
350
351#if !defined(ITK_WRAPPING_PARSER)
352 itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, MMIMetricPerThreadStruct, PaddedMMIMetricPerThreadStruct);
353 itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedMMIMetricPerThreadStruct, AlignedMMIMetricPerThreadStruct);
354 // Due to a bug in older version of Visual Studio where std::vector resize
355 // uses a value instead of a const reference, this must be a pointer.
356 // See
357 // https://thetweaker.wordpress.com/2010/05/05/stdvector-of-aligned-elements/
358 // https://connect.microsoft.com/VisualStudio/feedback/details/692988
359 std::unique_ptr<AlignedMMIMetricPerThreadStruct[]> m_MMIMetricPerThreadVariables;
360#endif
361
363 mutable bool m_ImplicitDerivativesSecondPass{ false };
364};
365} // end namespace itk
366
367#ifndef ITK_MANUAL_INSTANTIATION
368# include "itkMattesMutualInformationImageToImageMetric.hxx"
369#endif
370
371#endif
Array2D class representing a 2D array.
Definition itkArray2D.h:43
Derivative of a BSpline kernel used for density estimation and nonparametric regression.
BSpline kernel used for density estimation and nonparametric regression.
typename FixedImageType::ConstPointer FixedImageConstPointer
typename TransformType::Pointer TransformPointer
typename BSplineTransformWeightsType::ValueType WeightsValueType
Array< ParametersValueType > DerivativeType
typename MovingImageType::ConstPointer MovingImageConstPointer
typename BSplineTransformType::ParameterIndexArrayType BSplineTransformIndexArrayType
CovariantVector< double, Self::MovingImageDimension > ImageDerivativesType
typename BSplineTransformIndexArrayType::ValueType IndexValueType
typename BSplineTransformType::WeightsType BSplineTransformWeightsType
typename TransformType::OutputPointType MovingImagePointType
InterpolateImageFunction< MovingImageType, CoordinateRepresentationType > InterpolatorType
typename Superclass::ParametersValueType CoordinateRepresentationType
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
Transform< CoordinateRepresentationType, Self::MovingImageDimension, Self::FixedImageDimension > TransformType
Superclass::ParametersType ParametersType
typename TransformType::JacobianType TransformJacobianType
Templated n-dimensional image class.
Definition itkImage.h:89
Size< VImageDimension > SizeType
ImageRegion< VImageDimension > RegionType
SmartPointer< Self > Pointer
Definition itkImage.h:96
Index< VImageDimension > IndexType
Control indentation during Print() invocation.
Definition itkIndent.h:50
void GetValueAndDerivative(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
void GetDerivative(const ParametersType &parameters, DerivativeType &derivative) const override
void ComputeFixedImageParzenWindowIndices(FixedImageSampleContainer &samples)
CovariantVector< double, Self::MovingImageDimension > ImageDerivativesType
void ComputePDFDerivatives(ThreadIdType threadId, unsigned int sampleNumber, int pdfMovingIndex, const ImageDerivativesType &movingImageGradientValue, PDFValueType cubicBSplineDerivativeValue) const
bool GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample, const MovingImagePointType &mappedPoint, double movingImageValue, const ImageDerivativesType &movingImageGradientValue) const override
const JointPDFDerivativesType::Pointer GetJointPDFDerivatives() const
void GetValueAndDerivativeThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override
std::unique_ptr< AlignedMMIMetricPerThreadStruct[]> m_MMIMetricPerThreadVariables
typename TransformType::OutputPointType MovingImagePointType
BSplineDerivativeKernelFunction< 3, PDFValueType > CubicBSplineDerivativeFunctionType
itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, MMIMetricPerThreadStruct, PaddedMMIMetricPerThreadStruct)
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
~MattesMutualInformationImageToImageMetric() override=default
void GetValueThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override
itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedMMIMetricPerThreadStruct, AlignedMMIMetricPerThreadStruct)
void GetValueThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override
MeasureType GetValue(const ParametersType &parameters) const override
void PrintSelf(std::ostream &os, Indent indent) const override
bool GetValueThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample, const MovingImagePointType &mappedPoint, double movingImageValue) const override
void GetValueAndDerivativeThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override
static constexpr T max(const T &)
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
unsigned long SizeValueType
Definition itkIntTypes.h:86