go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkStatisticalShapePointPenalty.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 itkStatisticalShapePointPenalty_h
19#define itkStatisticalShapePointPenalty_h
20
22
23#include "itkPoint.h"
24#include "itkPointSet.h"
25#include "itkImage.h"
26#include "itkArray.h"
27#include <itkVariableSizeMatrix.h>
28
29#include <vnl/vnl_matrix.h>
30#include <vnl/vnl_math.h>
31#include <vnl/vnl_vector.h>
32#include <vnl/algo/vnl_real_eigensystem.h>
33#include <vnl/algo/vnl_symmetric_eigensystem.h>
34//#include <vnl/algo/vnl_svd.h>
35#include <vnl/algo/vnl_svd_economy.h>
36
37#include <string>
38
39namespace itk
40{
55template <class TFixedPointSet, class TMovingPointSet>
56class ITK_TEMPLATE_EXPORT StatisticalShapePointPenalty
57 : public SingleValuedPointSetToPointSetMetric<TFixedPointSet, TMovingPointSet>
58{
59public:
61
65 using Pointer = SmartPointer<Self>;
66 using ConstPointer = SmartPointer<const Self>;
67
69 itkNewMacro(Self);
70
73
75 using typename Superclass::TransformType;
76 using typename Superclass::TransformPointer;
80
81 using typename Superclass::MeasureType;
82 using typename Superclass::DerivativeType;
84 using typename Superclass::FixedPointSetType;
85 using typename Superclass::MovingPointSetType;
88
89 using typename Superclass::PointIterator;
90 using typename Superclass::PointDataIterator;
91
92 using typename Superclass::InputPointType;
93 using typename Superclass::OutputPointType;
94
95 using CoordRepType = typename OutputPointType::CoordRepType;
96 using VnlVectorType = vnl_vector<CoordRepType>;
97 using VnlMatrixType = vnl_matrix<CoordRepType>;
98 // typedef itk::Array<VnlVectorType *> ProposalDerivativeType;
99 using ProposalDerivativeType = typename std::vector<VnlVectorType *>;
100 // typedef typename vnl_vector<VnlVectorType *> ProposalDerivativeType; //Cannot be linked
101 using PCACovarianceType = vnl_svd_economy<CoordRepType>;
102
104 void
105 Initialize() override;
106
108 MeasureType
109 GetValue(const TransformParametersType & parameters) const override;
110
112 void
113 GetDerivative(const TransformParametersType & parameters, DerivativeType & Derivative) const override;
114
116 void
118 MeasureType & Value,
119 DerivativeType & Derivative) const override;
120
122 itkSetClampMacro(ShrinkageIntensity, MeasureType, 0.0, 1.0);
123 itkGetConstMacro(ShrinkageIntensity, MeasureType);
124
125 itkSetMacro(ShrinkageIntensityNeedsUpdate, bool);
126 itkBooleanMacro(ShrinkageIntensityNeedsUpdate);
127
129 itkSetClampMacro(BaseVariance, MeasureType, -1.0, NumericTraits<MeasureType>::max());
130 itkGetConstMacro(BaseVariance, MeasureType);
131
132 itkSetMacro(BaseVarianceNeedsUpdate, bool);
133 itkBooleanMacro(BaseVarianceNeedsUpdate);
134
135 itkSetClampMacro(CentroidXVariance, MeasureType, -1.0, NumericTraits<MeasureType>::max());
136 itkGetConstMacro(CentroidXVariance, MeasureType);
137
138 itkSetClampMacro(CentroidYVariance, MeasureType, -1.0, NumericTraits<MeasureType>::max());
139 itkGetConstMacro(CentroidYVariance, MeasureType);
140
141 itkSetClampMacro(CentroidZVariance, MeasureType, -1.0, NumericTraits<MeasureType>::max());
142 itkGetConstMacro(CentroidZVariance, MeasureType);
143
144 itkSetClampMacro(SizeVariance, MeasureType, -1.0, NumericTraits<MeasureType>::max());
145 itkGetConstMacro(SizeVariance, MeasureType);
146
147 itkSetMacro(VariancesNeedsUpdate, bool);
148 itkBooleanMacro(VariancesNeedsUpdate);
149
150 itkSetClampMacro(CutOffValue, MeasureType, 0.0, NumericTraits<MeasureType>::max());
151 itkGetConstMacro(CutOffValue, MeasureType);
152
153 itkSetClampMacro(CutOffSharpness,
154 MeasureType,
155 NumericTraits<MeasureType>::NonpositiveMin(),
156 NumericTraits<MeasureType>::max());
157 itkGetConstMacro(CutOffSharpness, MeasureType);
158
159 itkSetMacro(ShapeModelCalculation, int);
160 itkGetConstReferenceMacro(ShapeModelCalculation, int);
161
162 itkSetMacro(NormalizedShapeModel, bool);
163 itkGetConstReferenceMacro(NormalizedShapeModel, bool);
164 itkBooleanMacro(NormalizedShapeModel);
165
166 itkSetConstObjectMacro(EigenVectors, vnl_matrix<double>);
167 itkSetConstObjectMacro(EigenValues, vnl_vector<double>);
168 itkSetConstObjectMacro(MeanVector, vnl_vector<double>);
169
170 itkSetConstObjectMacro(CovarianceMatrix, vnl_matrix<double>);
171
172protected:
175
177 void
178 PrintSelf(std::ostream & os, Indent indent) const override;
179
180private:
181 void
182 FillProposalVector(const OutputPointType & fixedPoint, const unsigned int vertexindex) const;
183
184 void
185 FillProposalDerivative(const OutputPointType & fixedPoint, const unsigned int vertexindex) const;
186
187 void
188 UpdateCentroidAndAlignProposalVector(const unsigned int shapeLength) const;
189
190 void
191 UpdateCentroidAndAlignProposalDerivative(const unsigned int shapeLength) const;
192
193 void
194 UpdateL2(const unsigned int shapeLength) const;
195
196 void
197 NormalizeProposalVector(const unsigned int shapeLength) const;
198
199 void
200 UpdateL2AndNormalizeProposalDerivative(const unsigned int shapeLength) const;
201
202 void
203 CalculateValue(MeasureType & value,
204 VnlVectorType & differenceVector,
205 VnlVectorType & centerrotated,
206 VnlVectorType & eigrot) const;
207
208 void
209 CalculateDerivative(DerivativeType & derivative,
210 const MeasureType & value,
211 const VnlVectorType & differenceVector,
212 const VnlVectorType & eigrot,
213 const unsigned int shapeLength) const;
214
215 void
216 CalculateCutOffValue(MeasureType & value) const;
217
218 void
219 CalculateCutOffDerivative(typename DerivativeType::element_type & derivativeElement, const MeasureType & value) const;
220
221 const VnlVectorType * m_MeanVector{};
222 const VnlMatrixType * m_CovarianceMatrix{};
223 const VnlMatrixType * m_EigenVectors{};
224 const VnlVectorType * m_EigenValues{};
225
226 VnlMatrixType * m_InverseCovarianceMatrix{};
227
228 double m_CentroidXVariance{};
229 double m_CentroidXStd{};
230 double m_CentroidYVariance{};
231 double m_CentroidYStd{};
232 double m_CentroidZVariance{};
233 double m_CentroidZStd{};
234 double m_SizeVariance{};
235 double m_SizeStd{};
236
237 bool m_ShrinkageIntensityNeedsUpdate{};
238 bool m_BaseVarianceNeedsUpdate{};
239 bool m_VariancesNeedsUpdate{};
240
241 VnlVectorType * m_EigenValuesRegularized{};
242
243 mutable ProposalDerivativeType * m_ProposalDerivative{};
244 unsigned int m_ProposalLength{};
245 bool m_NormalizedShapeModel{};
246 int m_ShapeModelCalculation{};
247 double m_ShrinkageIntensity{};
248 double m_BaseVariance{};
249 double m_BaseStd{};
250 mutable VnlVectorType m_ProposalVector{};
251 mutable VnlVectorType m_MeanValues{};
252
253 double m_CutOffValue{};
254 double m_CutOffSharpness{};
255};
256
257} // end namespace itk
258
259#ifndef ITK_MANUAL_INSTANTIATION
260# include "itkStatisticalShapePointPenalty.hxx"
261#endif
262
263#endif
Transform maps points, vectors and covariant vectors from an input space to an output space.
typename FixedPointSetType::PointDataContainer::ConstIterator PointDataIterator
typename TransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
typename MovingPointSetType::ConstPointer MovingPointSetConstPointer
typename FixedPointSetType::ConstPointer FixedPointSetConstPointer
typename FixedPointSetType::PointsContainer::ConstIterator PointIterator
Computes the Mahalanobis distance between the transformed shape and a mean shape. A model mean and co...
void UpdateCentroidAndAlignProposalVector(const unsigned int shapeLength) const
void NormalizeProposalVector(const unsigned int shapeLength) const
void UpdateCentroidAndAlignProposalDerivative(const unsigned int shapeLength) const
void UpdateL2AndNormalizeProposalDerivative(const unsigned int shapeLength) const
vnl_svd_economy< CoordRepType > PCACovarianceType
void CalculateCutOffDerivative(typename DerivativeType::element_type &derivativeElement, const MeasureType &value) const
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const override
void CalculateDerivative(DerivativeType &derivative, const MeasureType &value, const VnlVectorType &differenceVector, const VnlVectorType &eigrot, const unsigned int shapeLength) const
void FillProposalVector(const OutputPointType &fixedPoint, const unsigned int vertexindex) const
void CalculateValue(MeasureType &value, VnlVectorType &differenceVector, VnlVectorType &centerrotated, VnlVectorType &eigrot) const
void PrintSelf(std::ostream &os, Indent indent) const override
MeasureType GetValue(const TransformParametersType &parameters) const override
typename TransformType::ParametersType TransformParametersType
void CalculateCutOffValue(MeasureType &value) const
typename OutputPointType::CoordRepType CoordRepType
typename TransformType::OutputPointType OutputPointType
ITK_DISALLOW_COPY_AND_MOVE(StatisticalShapePointPenalty)
void GetDerivative(const TransformParametersType &parameters, DerivativeType &Derivative) const override
void FillProposalDerivative(const OutputPointType &fixedPoint, const unsigned int vertexindex) const
void UpdateL2(const unsigned int shapeLength) const
typename std::vector< VnlVectorType * > ProposalDerivativeType


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