ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkMutualInformationImageToImageMetric.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 itkMutualInformationImageToImageMetric_h
19#define itkMutualInformationImageToImageMetric_h
20
22#include "itkPoint.h"
23
24#include "itkIndex.h"
26
27namespace itk
28{
90template <typename TFixedImage, typename TMovingImage>
91class ITK_TEMPLATE_EXPORT MutualInformationImageToImageMetric : public ImageToImageMetric<TFixedImage, TMovingImage>
92{
93public:
94 ITK_DISALLOW_COPY_AND_MOVE(MutualInformationImageToImageMetric);
95
101
103 itkNewMacro(Self);
104
106 itkOverrideGetNameOfClassMacro(MutualInformationImageToImageMetric);
107
109 using typename Superclass::TransformType;
110 using typename Superclass::TransformPointer;
112 using typename Superclass::InterpolatorType;
113 using typename Superclass::MeasureType;
114 using typename Superclass::DerivativeType;
115 using typename Superclass::ParametersType;
116 using typename Superclass::FixedImageType;
117 using typename Superclass::MovingImageType;
120
122 using FixedImageIndexType = typename FixedImageType::IndexType;
123 using FixedImageIndexValueType = typename FixedImageIndexType::IndexValueType;
124 using MovingImageIndexType = typename MovingImageType::IndexType;
125 using FixedImagePointType = typename TransformType::InputPointType;
126 using MovingImagePointType = typename TransformType::OutputPointType;
127
129
131 static constexpr unsigned int MovingImageDimension = MovingImageType::ImageDimension;
132
134 void
135 GetDerivative(const ParametersType & parameters, DerivativeType & derivative) const override;
136
139 GetValue(const ParametersType & parameters) const override;
140
142 void
144 MeasureType & value,
145 DerivativeType & derivative) const override;
146
151 void
152 SetNumberOfSpatialSamples(unsigned int num);
153
155 itkGetConstReferenceMacro(NumberOfSpatialSamples, unsigned int);
156
162 itkSetClampMacro(MovingImageStandardDeviation,
163 double,
166 itkGetConstReferenceMacro(MovingImageStandardDeviation, double);
168
174 itkSetClampMacro(FixedImageStandardDeviation,
175 double,
178 itkGetConstMacro(FixedImageStandardDeviation, double);
180
183 itkSetObjectMacro(KernelFunction, KernelFunctionType);
184 itkGetModifiableObjectMacro(KernelFunction, KernelFunctionType);
186
187
188protected:
191 void
192 PrintSelf(std::ostream & os, Indent indent) const override;
193
194private:
201 {
202 public:
204 ~SpatialSample() = default;
205
207 double FixedImageValue{ 0.0 };
208 double MovingImageValue{ 0.0 };
209 };
210
212 using SpatialSampleContainer = std::vector<SpatialSample>;
213
217
221
226
228
240 virtual void
242
243 /*
244 * Calculate derivatives of the image intensity at the specified point with respect to the transform parameters.
245 *
246 * \todo This should really be done by the mapper.
247 *
248 * \todo This is a temporary solution until this feature is implemented
249 * in the mapper. This solution only works for any transform
250 * that support ComputeJacobianWithRespectToParameters()
251 */
252 void
254
257
259};
260} // end namespace itk
261
262#ifndef ITK_MANUAL_INSTANTIATION
263# include "itkMutualInformationImageToImageMetric.hxx"
264#endif
265
266#endif
Calculate the derivative by central differencing.
typename FixedImageType::ConstPointer FixedImageConstPointer
typename TransformType::Pointer TransformPointer
Array< ParametersValueType > DerivativeType
typename MovingImageType::ConstPointer MovingImageConstPointer
InterpolateImageFunction< MovingImageType, CoordinateRepresentationType > InterpolatorType
typename Superclass::ParametersValueType CoordinateRepresentationType
Transform< CoordinateRepresentationType, Self::MovingImageDimension, Self::FixedImageDimension > TransformType
Superclass::ParametersType ParametersType
typename TransformType::JacobianType TransformJacobianType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Kernel used for density estimation and nonparametric regression.
typename Superclass::MovingImageConstPointer MovingImageCosntPointer
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename TransformType::OutputPointType MovingImagePointType
virtual void SampleFixedImageDomain(SpatialSampleContainer &samples) const
~MutualInformationImageToImageMetric() override=default
MeasureType GetValue(const ParametersType &parameters) const override
CentralDifferenceImageFunction< MovingImageType, CoordinateRepresentationType > DerivativeFunctionType
void GetValueAndDerivative(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
void SetNumberOfSpatialSamples(unsigned int num)
void CalculateDerivatives(const FixedImagePointType &, DerivativeType &, TransformJacobianType &) const
ImageToImageMetric< TFixedImage, TMovingImage > Superclass
void PrintSelf(std::ostream &os, Indent indent) const override
typename TransformType::JacobianType TransformJacobianType
void GetDerivative(const ParametersType &parameters, DerivativeType &derivative) const override
static constexpr T NonpositiveMin()
static constexpr T max(const T &)
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....