go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkSumSquaredTissueVolumeDifferenceImageToImageMetric.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 itkSumSquaredTissueVolumeDifferenceImageToImageMetric_h
19#define itkSumSquaredTissueVolumeDifferenceImageToImageMetric_h
20
22
23namespace itk
24{
25
60template <class TFixedImage, class TMovingImage>
62 : public AdvancedImageToImageMetric<TFixedImage, TMovingImage>
63{
64public:
66
70 using Pointer = SmartPointer<Self>;
71 using ConstPointer = SmartPointer<const Self>;
72
74 itkNewMacro(Self);
75
78
80 using typename Superclass::CoordinateRepresentationType;
81 using typename Superclass::MovingImageType;
82 using typename Superclass::MovingImagePixelType;
83 using typename Superclass::MovingImageConstPointer;
84 using typename Superclass::FixedImageType;
85 using typename Superclass::FixedImageConstPointer;
86 using typename Superclass::FixedImageRegionType;
87 using typename Superclass::InputPointType;
88 using typename Superclass::OutputPointType;
89 using typename Superclass::TransformParametersType;
90 using typename Superclass::TransformJacobianType;
91 using typename Superclass::InterpolatorType;
92 using typename Superclass::InterpolatorPointer;
93 using typename Superclass::RealType;
94 using typename Superclass::GradientPixelType;
95 using typename Superclass::GradientImageType;
96 using typename Superclass::GradientImagePointer;
97 using typename Superclass::FixedImageMaskType;
101 using typename Superclass::MeasureType;
102 using typename Superclass::DerivativeType;
103 using typename Superclass::DerivativeValueType;
104 using typename Superclass::ParametersType;
105 using typename Superclass::FixedImagePixelType;
107 using typename Superclass::ImageSamplerType;
108 using typename Superclass::ImageSamplerPointer;
116
119 using SpatialJacobianType = typename TransformType::SpatialJacobianType;
120 using JacobianOfSpatialJacobianType = typename TransformType::JacobianOfSpatialJacobianType;
121 using SpatialHessianType = typename TransformType::SpatialHessianType;
122 using JacobianOfSpatialHessianType = typename TransformType::JacobianOfSpatialHessianType;
123 using InternalMatrixType = typename TransformType::InternalMatrixType;
124
126 itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension);
127
129 itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension);
130
132 virtual MeasureType
133 GetValueSingleThreaded(const TransformParametersType & parameters) const;
134
135 MeasureType
136 GetValue(const TransformParametersType & parameters) const override;
137
139 void
140 GetDerivative(const TransformParametersType & parameters, DerivativeType & derivative) const override;
141
143 void
144 GetValueAndDerivative(const TransformParametersType & parameters,
145 MeasureType & Value,
146 DerivativeType & Derivative) const override;
147
149 void
150 GetValueAndDerivativeSingleThreaded(const TransformParametersType & parameters,
151 MeasureType & measure,
152 DerivativeType & derivative) const;
153
155 itkSetMacro(AirValue, RealType);
156 itkGetConstMacro(AirValue, RealType);
157
159 itkSetMacro(TissueValue, RealType);
160 itkGetConstMacro(TissueValue, RealType);
161
162protected:
165
166 void
167 PrintSelf(std::ostream & os, Indent indent) const override;
168
172 using typename Superclass::FixedImageIndexType;
175 using typename Superclass::FixedImagePointType;
181
185 void
186 EvaluateTransformJacobianInnerProduct(const TransformJacobianType & jacobian,
187 const MovingImageDerivativeType & movingImageDerivative,
188 DerivativeType & imageJacobian) const override;
189
192 void
193 UpdateValueAndDerivativeTerms(const RealType fixedImageValue,
194 const RealType movingImageValue,
195 const DerivativeType & imageJacobian,
196 const NonZeroJacobianIndicesType & nzji,
197 const RealType spatialJacobianDeterminant,
198 const DerivativeType & jacobianOfSpatialJacobianDeterminant,
199 MeasureType & measure,
200 DerivativeType & deriv) const;
201
207 void
209 const JacobianOfSpatialJacobianType & jacobianOfSpatialJacobian,
210 const SpatialJacobianType & inverseSpatialJacobian,
211 DerivativeType & jacobianOfSpatialJacobianDeterminant) const;
212
214 void
215 ThreadedGetValue(ThreadIdType threadID) const override;
216
218 void
219 AfterThreadedGetValue(MeasureType & value) const override;
220
222 void
223 ThreadedGetValueAndDerivative(ThreadIdType threadId) const override;
224
226 void
227 AfterThreadedGetValueAndDerivative(MeasureType & measure, DerivativeType & derivative) const override;
228
229private:
231 RealType m_AirValue{ -1000.0 };
232
234 RealType m_TissueValue{ 55.0 };
235
236}; // end class SumSquaredTissueVolumeDifferenceImageToImageMetric
237
238} // end namespace itk
239
240#ifndef ITK_MANUAL_INSTANTIATION
241# include "itkSumSquaredTissueVolumeDifferenceImageToImageMetric.hxx"
242#endif
243
244#endif // end #ifndef itkSumSquaredTissueVolumeDifferenceImageToImageMetric_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
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename MovingImageType::IndexType MovingImageIndexType
typename ImageSamplerType::OutputVectorContainerType ImageSampleContainerType
ImageMaskSpatialObject< Self::FixedImageDimension > FixedImageMaskType
typename FixedImageType::IndexType FixedImageIndexType
SmartPointer< MovingImageMaskType > MovingImageMaskPointer
typename BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
typename MovingImageLimiterType::OutputType MovingImageLimiterOutputType
typename TransformType::InputPointType FixedImagePointType
typename FixedImageLimiterType::OutputType FixedImageLimiterOutputType
SmartPointer< FixedImageMaskType > FixedImageMaskPointer
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, double > BSplineInterpolatorType
typename ImageSamplerType::Pointer ImageSamplerPointer
typename InterpolatorType::ContinuousIndexType MovingImageContinuousIndexType
ImageMaskSpatialObject< Self::MovingImageDimension > MovingImageMaskType
Transform maps points, vectors and covariant vectors from an input space to an output space.
This class is a base class for any image sampler.
Base class for all ITK limiter function objects.
Compute sum of square tissue volume difference between two images.
MeasureType GetValue(const TransformParametersType &parameters) const override
virtual MeasureType GetValueSingleThreaded(const TransformParametersType &parameters) const
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension)
void EvaluateJacobianOfSpatialJacobianDeterminantInnerProduct(const JacobianOfSpatialJacobianType &jacobianOfSpatialJacobian, const SpatialJacobianType &inverseSpatialJacobian, DerivativeType &jacobianOfSpatialJacobianDeterminant) const
ITK_DISALLOW_COPY_AND_MOVE(SumSquaredTissueVolumeDifferenceImageToImageMetric)
void UpdateValueAndDerivativeTerms(const RealType fixedImageValue, const RealType movingImageValue, const DerivativeType &imageJacobian, const NonZeroJacobianIndicesType &nzji, const RealType spatialJacobianDeterminant, const DerivativeType &jacobianOfSpatialJacobianDeterminant, MeasureType &measure, DerivativeType &deriv) const
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
void AfterThreadedGetValue(MeasureType &value) const override
void ThreadedGetValue(ThreadIdType threadID) const override
void EvaluateTransformJacobianInnerProduct(const TransformJacobianType &jacobian, const MovingImageDerivativeType &movingImageDerivative, DerivativeType &imageJacobian) const override
void ThreadedGetValueAndDerivative(ThreadIdType threadId) const override
void PrintSelf(std::ostream &os, Indent indent) const override
void GetDerivative(const TransformParametersType &parameters, DerivativeType &derivative) const override
typename BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
void GetValueAndDerivativeSingleThreaded(const TransformParametersType &parameters, MeasureType &measure, DerivativeType &derivative) const
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const override
void AfterThreadedGetValueAndDerivative(MeasureType &measure, DerivativeType &derivative) const override


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