go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkAdvancedNormalizedCorrelationImageToImageMetric.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 itkAdvancedNormalizedCorrelationImageToImageMetric_h
19#define itkAdvancedNormalizedCorrelationImageToImageMetric_h
20
22#include <vector>
23
24namespace itk
25{
90template <class TFixedImage, class TMovingImage>
92 : public AdvancedImageToImageMetric<TFixedImage, TMovingImage>
93{
94public:
96
100 using Pointer = SmartPointer<Self>;
101 using ConstPointer = SmartPointer<const Self>;
102
104 itkNewMacro(Self);
105
108
110 using typename Superclass::CoordinateRepresentationType;
111 using typename Superclass::MovingImageType;
112 using typename Superclass::MovingImagePixelType;
113 using typename Superclass::MovingImageConstPointer;
114 using typename Superclass::FixedImageType;
115 using typename Superclass::FixedImageConstPointer;
116 using typename Superclass::FixedImageRegionType;
117 using typename Superclass::TransformType;
118 using typename Superclass::TransformPointer;
119 using typename Superclass::InputPointType;
120 using typename Superclass::OutputPointType;
121 using typename Superclass::TransformParametersType;
122 using typename Superclass::TransformJacobianType;
124 using typename Superclass::InterpolatorType;
125 using typename Superclass::InterpolatorPointer;
126 using typename Superclass::RealType;
127 using typename Superclass::GradientPixelType;
128 using typename Superclass::GradientImageType;
129 using typename Superclass::GradientImagePointer;
130 using typename Superclass::GradientImageFilterType;
131 using typename Superclass::GradientImageFilterPointer;
132 using typename Superclass::FixedImageMaskType;
133 using typename Superclass::FixedImageMaskPointer;
134 using typename Superclass::MovingImageMaskType;
135 using typename Superclass::MovingImageMaskPointer;
136 using typename Superclass::MeasureType;
137 using typename Superclass::DerivativeType;
138 using typename Superclass::DerivativeValueType;
139 using typename Superclass::ParametersType;
140 using typename Superclass::FixedImagePixelType;
142 using typename Superclass::ImageSamplerType;
143 using typename Superclass::ImageSamplerPointer;
151 using typename Superclass::ThreaderType;
152 using typename Superclass::ThreadInfoType;
153
155 itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension);
156
158 itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension);
159
161 MeasureType
162 GetValue(const TransformParametersType & parameters) const override;
163
165 void
166 GetDerivative(const TransformParametersType & parameters, DerivativeType & derivative) const override;
167
169 void
170 GetValueAndDerivativeSingleThreaded(const TransformParametersType & parameters,
171 MeasureType & value,
172 DerivativeType & derivative) const;
173
174 void
175 GetValueAndDerivative(const TransformParametersType & parameters,
176 MeasureType & value,
177 DerivativeType & derivative) const override;
178
184 itkSetMacro(SubtractMean, bool);
185 itkGetConstReferenceMacro(SubtractMean, bool);
186 itkBooleanMacro(SubtractMean);
187
188protected:
191
192 void
193 PrintSelf(std::ostream & os, Indent indent) const override;
194
198 using typename Superclass::FixedImageIndexType;
201 using typename Superclass::FixedImagePointType;
208
212 void
213 UpdateDerivativeTerms(const RealType & fixedImageValue,
214 const RealType & movingImageValue,
215 const DerivativeType & imageJacobian,
216 const NonZeroJacobianIndicesType & nzji,
217 DerivativeType & derivativeF,
218 DerivativeType & derivativeM,
219 DerivativeType & differential) const;
220
225 void
227
229 void
230 ThreadedGetValueAndDerivative(ThreadIdType threadID) override;
231
233 void
234 AfterThreadedGetValueAndDerivative(MeasureType & value, DerivativeType & derivative) const override;
235
237 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
239
240private:
241 mutable bool m_SubtractMean;
242
243 using AccumulateType = typename NumericTraits<MeasureType>::AccumulateType;
244
249 {
251
257 };
258
260 {
267 DerivativeType st_DerivativeF;
268 DerivativeType st_DerivativeM;
269 DerivativeType st_Differential;
270 };
271 itkPadStruct(ITK_CACHE_LINE_ALIGNMENT,
273 PaddedCorrelationGetValueAndDerivativePerThreadStruct);
274 itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT,
275 PaddedCorrelationGetValueAndDerivativePerThreadStruct,
276 AlignedCorrelationGetValueAndDerivativePerThreadStruct);
277 mutable std::vector<AlignedCorrelationGetValueAndDerivativePerThreadStruct>
279};
280
281} // end namespace itk
282
283#ifndef ITK_MANUAL_INSTANTIATION
284# include "itkAdvancedNormalizedCorrelationImageToImageMetric.hxx"
285#endif
286
287#endif // end #ifndef itkAdvancedNormalizedCorrelationImageToImageMetric_h
An extension of the ITK ImageToImageMetric. It is the intended base class for all elastix metrics.
typename TransformType::OutputPointType MovingImagePointType
typename ImageSamplerType::OutputVectorContainerPointer ImageSampleContainerPointer
typename MovingImageType::RegionType MovingImageRegionType
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
typename FixedImageType::PixelType FixedImagePixelType
typename DerivativeType::ValueType DerivativeValueType
GradientImageFilter< MovingImageType, RealType, RealType > CentralDifferenceGradientFilterType
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename MovingImageType::IndexType MovingImageIndexType
typename ImageSamplerType::OutputVectorContainerType ImageSampleContainerType
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, double > BSplineInterpolatorType
typename ThreaderType::WorkUnitInfo ThreadInfoType
typename FixedImageType::IndexType FixedImageIndexType
typename BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
typename MovingImageLimiterType::OutputType MovingImageLimiterOutputType
typename TransformType::InputPointType FixedImagePointType
typename FixedImageLimiterType::OutputType FixedImageLimiterOutputType
typename ImageSamplerType::Pointer ImageSamplerPointer
typename InterpolatorType::ContinuousIndexType MovingImageContinuousIndexType
typename AdvancedTransformType::NumberOfParametersType NumberOfParametersType
Computes normalized correlation between two images, based on AdvancedImageToImageMetric....
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION AccumulateDerivativesThreaderCallback(void *arg)
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedCorrelationGetValueAndDerivativePerThreadStruct, AlignedCorrelationGetValueAndDerivativePerThreadStruct)
itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension)
MeasureType GetValue(const TransformParametersType &parameters) const override
void ThreadedGetValueAndDerivative(ThreadIdType threadID) override
void GetDerivative(const TransformParametersType &parameters, DerivativeType &derivative) const override
void UpdateDerivativeTerms(const RealType &fixedImageValue, const RealType &movingImageValue, const DerivativeType &imageJacobian, const NonZeroJacobianIndicesType &nzji, DerivativeType &derivativeF, DerivativeType &derivativeM, DerivativeType &differential) const
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, CorrelationGetValueAndDerivativePerThreadStruct, PaddedCorrelationGetValueAndDerivativePerThreadStruct)
void GetValueAndDerivativeSingleThreaded(const TransformParametersType &parameters, MeasureType &value, DerivativeType &derivative) const
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
void PrintSelf(std::ostream &os, Indent indent) const override
ITK_DISALLOW_COPY_AND_MOVE(AdvancedNormalizedCorrelationImageToImageMetric)
void AfterThreadedGetValueAndDerivative(MeasureType &value, DerivativeType &derivative) const override
std::vector< AlignedCorrelationGetValueAndDerivativePerThreadStruct > m_CorrelationGetValueAndDerivativePerThreadVariables
This class is a base class for any image sampler.
Base class for all ITK limiter function objects.


Generated on 2023-01-13 for elastix by doxygen 1.9.6 elastix logo