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 >
57 public SingleValuedPointSetToPointSetMetric< TFixedPointSet, TMovingPointSet >
58{
59public:
60
64 TFixedPointSet, TMovingPointSet > Superclass;
66 typedef SmartPointer< const Self > ConstPointer;
67
69 itkNewMacro( Self );
70
72 itkTypeMacro( StatisticalShapePointPenalty,
74
81
89
92
95
96 typedef typename OutputPointType::CoordRepType CoordRepType;
98 typedef vnl_matrix< CoordRepType > VnlMatrixType;
99 //typedef itk::Array<VnlVectorType *> ProposalDerivativeType;
100 typedef typename std::vector< VnlVectorType * > ProposalDerivativeType;
101 //typedef typename vnl_vector<VnlVectorType *> ProposalDerivativeType; //Cannot be linked
102 typedef vnl_svd_economy< CoordRepType > PCACovarianceType;
103
105 void Initialize( void ) override;
106
108 MeasureType GetValue( const TransformParametersType & parameters ) const override;
109
111 void GetDerivative( const TransformParametersType & parameters,
112 DerivativeType & Derivative ) const override;
113
116 MeasureType & Value, DerivativeType & Derivative ) const override;
117
119 itkSetClampMacro( ShrinkageIntensity, MeasureType, 0.0, 1.0 );
120 itkGetMacro( ShrinkageIntensity, MeasureType );
121
122 itkSetMacro( ShrinkageIntensityNeedsUpdate, bool );
123 itkBooleanMacro( ShrinkageIntensityNeedsUpdate );
124
126 itkSetClampMacro( BaseVariance, MeasureType,
127 -1.0, NumericTraits< MeasureType >::max() );
128 itkGetMacro( BaseVariance, MeasureType );
129
130 itkSetMacro( BaseVarianceNeedsUpdate, bool );
131 itkBooleanMacro( BaseVarianceNeedsUpdate );
132
133 itkSetClampMacro( CentroidXVariance, MeasureType,
134 -1.0, NumericTraits< MeasureType >::max() );
135 itkGetMacro( CentroidXVariance, MeasureType );
136
137 itkSetClampMacro( CentroidYVariance, MeasureType,
138 -1.0, NumericTraits< MeasureType >::max() );
139 itkGetMacro( CentroidYVariance, MeasureType );
140
141 itkSetClampMacro( CentroidZVariance, MeasureType,
142 -1.0, NumericTraits< MeasureType >::max() );
143 itkGetMacro( CentroidZVariance, MeasureType );
144
145 itkSetClampMacro( SizeVariance, MeasureType,
146 -1.0, NumericTraits< MeasureType >::max() );
147 itkGetMacro( SizeVariance, MeasureType );
148
149 itkSetMacro( VariancesNeedsUpdate, bool );
150 itkBooleanMacro( VariancesNeedsUpdate );
151
152 itkSetClampMacro( CutOffValue, MeasureType,
153 0.0, NumericTraits< MeasureType >::max() );
154 itkGetMacro( CutOffValue, MeasureType );
155
156 itkSetClampMacro( CutOffSharpness, MeasureType,
157 NumericTraits< MeasureType >::NonpositiveMin(), NumericTraits< MeasureType >::max() );
158 itkGetMacro( CutOffSharpness, MeasureType );
159
160 itkSetMacro( ShapeModelCalculation, int );
161 itkGetConstReferenceMacro( ShapeModelCalculation, int );
162
163 itkSetMacro( NormalizedShapeModel, bool );
164 itkGetConstReferenceMacro( NormalizedShapeModel, bool );
165 itkBooleanMacro( NormalizedShapeModel );
166
167 itkSetConstObjectMacro( EigenVectors, vnl_matrix< double > );
168 itkSetConstObjectMacro( EigenValues, vnl_vector< double > );
169 itkSetConstObjectMacro( MeanVector, vnl_vector< double > );
170
171 itkSetConstObjectMacro( CovarianceMatrix, vnl_matrix< double > );
172
173protected:
174
177
179 void PrintSelf( std::ostream & os, Indent indent ) const override;
180
181private:
182
183 StatisticalShapePointPenalty( const Self & ); // purposely not implemented
184 void operator=( const Self & ); // purposely not implemented
185
186 void FillProposalVector( const OutputPointType & fixedPoint,
187 const unsigned int vertexindex ) const;
188
189 void FillProposalDerivative( const OutputPointType & fixedPoint,
190 const unsigned int vertexindex ) const;
191
193 const unsigned int shapeLength ) const;
194
196 const unsigned int shapeLength ) const;
197
198 void UpdateL2( const unsigned int shapeLength ) const;
199
200 void NormalizeProposalVector( const unsigned int shapeLength ) const;
201
202 void UpdateL2AndNormalizeProposalDerivative( const unsigned int shapeLength ) const;
203
204 void CalculateValue( MeasureType & value, VnlVectorType & differenceVector,
205 VnlVectorType & centerrotated, VnlVectorType & eigrot ) const;
206
207 void CalculateDerivative( DerivativeType & derivative, const MeasureType & value,
208 const VnlVectorType & differenceVector, const VnlVectorType & centerrotated,
209 const VnlVectorType & eigrot, const unsigned int shapeLength ) const;
210
211 void CalculateCutOffValue( MeasureType & value ) const;
212
214 typename DerivativeType::element_type & derivativeElement,
215 const MeasureType &value ) const;
216
221
223
231 double m_SizeStd;
232
236
238
240 unsigned int m_ProposalLength;
245 double m_BaseStd;
248
251
252};
253
254} // end namespace itk
255
256#ifndef ITK_MANUAL_INSTANTIATION
257#include "itkStatisticalShapePointPenalty.hxx"
258#endif
259
260#endif
FixedPointSetType::PointDataContainer::ConstIterator PointDataIterator
TransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
AdvancedTransform< CoordinateRepresentationType, itkGetStaticConstMacro(FixedPointSetDimension), itkGetStaticConstMacro(MovingPointSetDimension) > TransformType
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
SingleValuedPointSetToPointSetMetric< TFixedPointSet, TMovingPointSet > Superclass
Superclass::PointDataIterator PointDataIterator
void NormalizeProposalVector(const unsigned int shapeLength) const
Superclass::DerivativeValueType DerivativeValueType
Superclass::MovingPointSetConstPointer MovingPointSetConstPointer
void UpdateCentroidAndAlignProposalDerivative(const unsigned int shapeLength) const
vnl_svd_economy< CoordRepType > PCACovarianceType
void UpdateL2AndNormalizeProposalDerivative(const unsigned int shapeLength) const
void CalculateCutOffDerivative(typename DerivativeType::element_type &derivativeElement, const MeasureType &value) const
Superclass::FixedPointSetType FixedPointSetType
Superclass::TransformParametersType TransformParametersType
StatisticalShapePointPenalty(const Self &)
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const override
std::vector< VnlVectorType * > ProposalDerivativeType
void Initialize(void) override
Superclass::TransformJacobianType TransformJacobianType
void CalculateDerivative(DerivativeType &derivative, const MeasureType &value, const VnlVectorType &differenceVector, const VnlVectorType &centerrotated, 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
void CalculateCutOffValue(MeasureType &value) const
Superclass::MovingPointSetType MovingPointSetType
void GetDerivative(const TransformParametersType &parameters, DerivativeType &Derivative) const override
void FillProposalDerivative(const OutputPointType &fixedPoint, const unsigned int vertexindex) const
Superclass::FixedPointSetConstPointer FixedPointSetConstPointer
void UpdateL2(const unsigned int shapeLength) const
Superclass::NonZeroJacobianIndicesType NonZeroJacobianIndicesType


Generated on 1667476801 for elastix by doxygen 1.9.4 elastix logo