ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkJointHistogramMutualInformationImageToImageMetricv4.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
19#ifndef itkJointHistogramMutualInformationImageToImageMetricv4_h
20#define itkJointHistogramMutualInformationImageToImageMetricv4_h
21
23#include "itkImage.h"
25
28
29namespace itk
30{
37template <typename TFixedImage,
38 typename TMovingImage,
39 typename TVirtualImage = TFixedImage,
40 typename TInternalComputationValueType = double,
41 typename TMetricTraits =
44 : public ImageToImageMetricv4<TFixedImage, TMovingImage, TVirtualImage, TInternalComputationValueType, TMetricTraits>
45{
46public:
48
51 using Superclass =
55
57 itkNewMacro(Self);
58
61
64
66
68 using InternalComputationValueType = TInternalComputationValueType;
69
71 using typename Superclass::ParametersType;
74
76 using typename Superclass::MeasureType;
77 using typename Superclass::DerivativeType;
84
87
89 using typename Superclass::VirtualIndexType;
90 using typename Superclass::VirtualPointType;
92
93 /* Image dimension accessors */
94 static constexpr typename TVirtualImage::ImageDimensionType VirtualImageDimension = TVirtualImage::ImageDimension;
95 static constexpr typename TMovingImage::ImageDimensionType MovingImageDimension = TMovingImage::ImageDimension;
96
98 using PDFValueType = TInternalComputationValueType;
99
108
110 itkGetModifiableObjectMacro(JointPDF, JointPDFType);
111
112 // Declare the type for the derivative calculation
115 using JPDFGradientImagePointer = typename JPDFGradientImageType::Pointer;
116
119 using MarginalGradientImagePointer = typename MarginalGradientImageType::Pointer;
120
126
132
133
135 itkSetClampMacro(NumberOfHistogramBins, SizeValueType, 5, NumericTraits<SizeValueType>::max());
136 itkGetConstReferenceMacro(NumberOfHistogramBins, SizeValueType);
138
140 itkSetMacro(VarianceForJointPDFSmoothing, TInternalComputationValueType);
141 itkGetMacro(VarianceForJointPDFSmoothing, TInternalComputationValueType);
143
145 void
146 Initialize() override;
147
149 GetValue() const override;
150
151protected:
154
158 void
159 InitializeForIteration() const override;
160
164
167 inline void
169 const MovingImagePixelType movingImageValue,
170 JointPDFPointType & jointPDFpoint) const;
171
174 Self>;
175 friend class JointHistogramMutualInformationComputeJointPDFThreaderBase<ThreadedIndexedContainerPartitioner, Self>;
177 ThreadedImageRegionPartitioner<Self::VirtualImageDimension>,
178 Self>;
179 friend class JointHistogramMutualInformationComputeJointPDFThreader<ThreadedIndexedContainerPartitioner, Self>;
180
182 JointHistogramMutualInformationComputeJointPDFThreader<ThreadedImageRegionPartitioner<Self::VirtualImageDimension>,
183 Self>;
185 JointHistogramMutualInformationComputeJointPDFThreader<ThreadedIndexedContainerPartitioner, Self>;
186
187 typename JointHistogramMutualInformationDenseComputeJointPDFThreaderType::Pointer
189 typename JointHistogramMutualInformationSparseComputeJointPDFThreaderType::Pointer
191
193 ThreadedImageRegionPartitioner<Superclass::VirtualImageDimension>,
194 Superclass,
195 Self>;
196 friend class JointHistogramMutualInformationGetValueAndDerivativeThreader<ThreadedIndexedContainerPartitioner,
197 Superclass,
198 Self>;
199
202 ThreadedImageRegionPartitioner<Superclass::VirtualImageDimension>,
203 Superclass,
204 Self>;
206 JointHistogramMutualInformationGetValueAndDerivativeThreader<ThreadedIndexedContainerPartitioner, Superclass, Self>;
207
209 void
210 PrintSelf(std::ostream & os, Indent indent) const override;
211
214
215private:
218
221
224
226 TInternalComputationValueType m_VarianceForJointPDFSmoothing{};
227
230 TInternalComputationValueType m_FixedImageTrueMin{};
231 TInternalComputationValueType m_FixedImageTrueMax{};
232 TInternalComputationValueType m_MovingImageTrueMin{};
233 TInternalComputationValueType m_MovingImageTrueMax{};
234 TInternalComputationValueType m_FixedImageBinSize{};
235 TInternalComputationValueType m_MovingImageBinSize{};
236
237 TInternalComputationValueType m_JointPDFSum{};
239
240 TInternalComputationValueType m_Log2{};
242};
243
244} // end namespace itk
245
246#ifndef ITK_MANUAL_INSTANTIATION
247# include "itkJointHistogramMutualInformationImageToImageMetricv4.hxx"
248#endif
249
250#endif
OptimizerParameters< TInternalComputationValueType > ParametersType
TInternalComputationValueType ParametersValueType
A simple structure holding type information for ImageToImageMetricv4 classes.
Computes the gradient of an image by convolution with the first derivative of a Gaussian.
Image< CovariantVector< typename NumericTraits< typename JointPDFType::PixelType >::RealType, JointPDFType::ImageDimension >, JointPDFType::ImageDimension > OutputImageType
Templated n-dimensional image class.
Definition itkImage.h:89
Size< VImageDimension > SizeType
Point< PointValueType, VImageDimension > PointType
ImageRegion< VImageDimension > RegionType
typename IndexType::IndexValueType IndexValueType
Vector< SpacingValueType, VImageDimension > SpacingType
SmartPointer< Self > Pointer
Definition itkImage.h:96
Index< VImageDimension > IndexType
Provide a threaded computation of the joint PDF for JointHistogramMutualInformationImageToImageMetric...
Processes points for JointHistogramMutualInformationImageToImageMetricv4 GetValueAndDerivative().
JointHistogramMutualInformationGetValueAndDerivativeThreader< ThreadedImageRegionPartitioner< Superclass::VirtualImageDimension >, Superclass, Self > JointHistogramMutualInformationDenseGetValueAndDerivativeThreaderType
JointHistogramMutualInformationGetValueAndDerivativeThreader< ThreadedIndexedContainerPartitioner, Superclass, Self > JointHistogramMutualInformationSparseGetValueAndDerivativeThreaderType
JointHistogramMutualInformationSparseComputeJointPDFThreaderType::Pointer m_JointHistogramMutualInformationSparseComputeJointPDFThreader
void ComputeJointPDFPoint(const FixedImagePixelType fixedImageValue, const MovingImagePixelType movingImageValue, JointPDFPointType &jointPDFpoint) const
void PrintSelf(std::ostream &os, Indent indent) const override
ImageToImageMetricv4< TFixedImage, TMovingImage, TVirtualImage, TInternalComputationValueType, TMetricTraits > Superclass
itk::GradientRecursiveGaussianImageFilter< MarginalPDFType > MarginalGradientFilterType
JointHistogramMutualInformationComputeJointPDFThreader< ThreadedIndexedContainerPartitioner, Self > JointHistogramMutualInformationSparseComputeJointPDFThreaderType
JointHistogramMutualInformationComputeJointPDFThreader< ThreadedImageRegionPartitioner< Self::VirtualImageDimension >, Self > JointHistogramMutualInformationDenseComputeJointPDFThreaderType
JointHistogramMutualInformationDenseComputeJointPDFThreaderType::Pointer m_JointHistogramMutualInformationDenseComputeJointPDFThreader
Linearly interpolate an image at specified positions.
static constexpr T max(const T &)
Array< TInternalComputationValueType > DerivativeType
Implements transparent reference counting.
Class for partitioning of an ImageRegion.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86