go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkParzenWindowHistogramImageToImageMetric.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 __itkParzenWindowHistogramImageToImageMetric_H__
19 #define __itkParzenWindowHistogramImageToImageMetric_H__
20 
22 #include "itkKernelFunctionBase2.h"
23 
24 
25 namespace itk
26 {
73 template< class TFixedImage, class TMovingImage >
75  public AdvancedImageToImageMetric< TFixedImage, TMovingImage >
76 {
77 public:
78 
83  typedef SmartPointer< const Self > ConstPointer;
84 
87 
104  typedef typename Superclass::RealType RealType;
131 
133  itkStaticConstMacro( FixedImageDimension, unsigned int,
134  FixedImageType::ImageDimension );
135 
137  itkStaticConstMacro( MovingImageDimension, unsigned int,
138  MovingImageType::ImageDimension );
139 
146  void Initialize( void ) throw ( ExceptionObject );
147 
153  const ParametersType & parameters,
154  DerivativeType & Derivative ) const;
155 
161  void GetValueAndDerivative( const ParametersType & parameters,
162  MeasureType & value, DerivativeType & derivative ) const;
163 
170  itkSetClampMacro( NumberOfFixedHistogramBins, unsigned long,
171  4, NumericTraits< unsigned long >::max() );
172  itkGetMacro( NumberOfFixedHistogramBins, unsigned long );
173 
180  itkSetClampMacro( NumberOfMovingHistogramBins, unsigned long,
181  4, NumericTraits< unsigned long >::max() );
182  itkGetMacro( NumberOfMovingHistogramBins, unsigned long );
183 
185  itkSetClampMacro( FixedKernelBSplineOrder, unsigned int, 0, 3 );
186  itkGetConstMacro( FixedKernelBSplineOrder, unsigned int );
187 
189  itkSetClampMacro( MovingKernelBSplineOrder, unsigned int, 0, 3 );
190  itkGetConstMacro( MovingKernelBSplineOrder, unsigned int );
191 
195  itkSetMacro( UseExplicitPDFDerivatives, bool );
196  itkGetConstReferenceMacro( UseExplicitPDFDerivatives, bool );
197  itkBooleanMacro( UseExplicitPDFDerivatives );
198 
202  itkSetMacro( UseDerivative, bool );
203  itkGetConstMacro( UseDerivative, bool );
204 
208  itkSetMacro( UseFiniteDifferenceDerivative, bool );
209  itkGetConstMacro( UseFiniteDifferenceDerivative, bool );
210 
215  itkSetMacro( FiniteDifferencePerturbation, double );
216  itkGetConstMacro( FiniteDifferencePerturbation, double );
217 
218 protected:
219 
222 
225 
227  void PrintSelf( std::ostream & os, Indent indent ) const;
228 
234  typedef typename FixedImageType::OffsetValueType OffsetValueType;
243 
245  typedef double PDFValueType;
246  typedef float PDFDerivativeValueType;
249  typedef typename JointPDFType::Pointer JointPDFPointer;
251  typedef typename JointPDFDerivativesType::Pointer JointPDFDerivativesPointer;
253  typedef typename IncrementalMarginalPDFType::Pointer IncrementalMarginalPDFPointer;
254  typedef JointPDFType::IndexType JointPDFIndexType;
255  typedef JointPDFType::RegionType JointPDFRegionType;
256  typedef JointPDFType::SizeType JointPDFSizeType;
257  typedef JointPDFDerivativesType::IndexType JointPDFDerivativesIndexType;
258  typedef JointPDFDerivativesType::RegionType JointPDFDerivativesRegionType;
259  typedef JointPDFDerivativesType::SizeType JointPDFDerivativesSizeType;
260  typedef IncrementalMarginalPDFType::IndexType IncrementalMarginalPDFIndexType;
261  typedef IncrementalMarginalPDFType::RegionType IncrementalMarginalPDFRegionType;
262  typedef IncrementalMarginalPDFType::SizeType IncrementalMarginalPDFSizeType;
264 
268 
272  mutable double m_Alpha;
275 
287  mutable JointPDFRegionType m_JointPDFWindow; // no need for mutable anymore?
294 
299 
301  mutable std::vector< JointPDFPointer > m_ThreaderJointPDFs;
302 
306  struct ParzenWindowHistogramMultiThreaderParameterType // can't we use the one from AdvancedImageToImageMetric ?
307  {
309  };
310  ParzenWindowHistogramMultiThreaderParameterType m_ParzenWindowHistogramThreaderParameters;
311 
313  {
316  };
317  itkPadStruct( ITK_CACHE_LINE_ALIGNMENT, ParzenWindowHistogramGetValueAndDerivativePerThreadStruct,
318  PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct );
319  itkAlignedTypedef( ITK_CACHE_LINE_ALIGNMENT, PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct,
320  AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct );
321  mutable AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct * m_ParzenWindowHistogramGetValueAndDerivativePerThreadVariables;
323 
325  virtual void InitializeThreadingParameters( void ) const;
326 
328  inline void ThreadedComputePDFs( ThreadIdType threadId );
329 
331  inline void AfterThreadedComputePDFs( void ) const;
332 
334  static ITK_THREAD_RETURN_TYPE ComputePDFsThreaderCallback( void * arg );
335 
338 
346  double parzenWindowTerm, OffsetValueType parzenWindowIndex,
347  const KernelFunctionType * kernel,
348  ParzenValueContainerType & parzenValues ) const;
349 
354  const RealType & fixedImageValue,
355  const RealType & movingImageValue,
356  const DerivativeType * imageJacobian,
357  const NonZeroJacobianIndicesType * nzji,
358  JointPDFType * jointPDF ) const;
359 
371  RealType fixedImageValue, RealType movingImageValue, RealType movingMaskValue,
372  const DerivativeType & movingImageValuesRight,
373  const DerivativeType & movingImageValuesLeft,
374  const DerivativeType & movingMaskValuesRight,
375  const DerivativeType & movingMaskValuesLeft,
376  const NonZeroJacobianIndicesType & nzji ) const;
377 
384  const JointPDFIndexType & pdfIndex, double factor,
385  const DerivativeType & imageJacobian,
386  const NonZeroJacobianIndicesType & nzji ) const;
387 
389  virtual void NormalizeJointPDF(
390  JointPDFType * pdf, const double & factor ) const;
391 
394  JointPDFDerivativesType * pdf, const double & factor ) const;
395 
400  virtual void ComputeMarginalPDF(
401  const JointPDFType * jointPDF,
402  MarginalPDFType & marginalPDF,
403  const unsigned int & direction ) const;
404 
409  const JointPDFDerivativesType * incrementalPDF,
410  IncrementalMarginalPDFType * fixedIncrementalMarginalPDF,
411  IncrementalMarginalPDFType * movingIncrementalMarginalPDF ) const;
412 
422  virtual void ComputePDFsAndPDFDerivatives( const ParametersType & parameters ) const;
423 
447  virtual void ComputePDFsAndIncrementalPDFs( const ParametersType & parameters ) const;
448 
457  virtual void ComputePDFsSingleThreaded( const ParametersType & parameters ) const;
458 
459  virtual void ComputePDFs( const ParametersType & parameters ) const;
460 
462  virtual void InitializeHistograms( void );
463 
464  virtual void InitializeKernels( void );
465 
471  const ParametersType & itkNotUsed( parameters ),
472  MeasureType & itkNotUsed( value ),
473  DerivativeType & itkNotUsed( derivative ) ) const {}
474 
480  const ParametersType & itkNotUsed( parameters ),
481  MeasureType & itkNotUsed( value ),
482  DerivativeType & itkNotUsed( derivative ) ) const {}
483 
484 private:
485 
487  ParzenWindowHistogramImageToImageMetric( const Self & ); // purposely not implemented
489  void operator=( const Self & ); // purposely not implemented
490 
500 
501 };
502 
503 } // end namespace itk
504 
505 #ifndef ITK_MANUAL_INSTANTIATION
506 #include "itkParzenWindowHistogramImageToImageMetric.hxx"
507 #endif
508 
509 #endif // end #ifndef __itkParzenWindowHistogramImageToImageMetric_H__
itk::ParzenWindowHistogramImageToImageMetric::TransformPointer
Superclass::TransformPointer TransformPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:97
itk::ParzenWindowHistogramImageToImageMetric::FixedImageLimiterOutputType
Superclass::FixedImageLimiterOutputType FixedImageLimiterOutputType
Definition: itkParzenWindowHistogramImageToImageMetric.h:126
itk::AdvancedImageToImageMetric::TransformParametersType
Superclass::TransformParametersType TransformParametersType
Definition: itkAdvancedImageToImageMetric.h:113
itk::ParzenWindowHistogramImageToImageMetric::FixedImagePixelType
Superclass::FixedImagePixelType FixedImagePixelType
Definition: itkParzenWindowHistogramImageToImageMetric.h:118
itk::ParzenWindowHistogramImageToImageMetric::InterpolatorType
Superclass::InterpolatorType InterpolatorType
Definition: itkParzenWindowHistogramImageToImageMetric.h:102
itk::ParzenWindowHistogramImageToImageMetric::NormalizeJointPDFDerivatives
virtual void NormalizeJointPDFDerivatives(JointPDFDerivativesType *pdf, const double &factor) const
itk::ParzenWindowHistogramImageToImageMetric::Pointer
SmartPointer< Self > Pointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:82
itk::ParzenWindowHistogramImageToImageMetric::OutputPointType
Superclass::OutputPointType OutputPointType
Definition: itkParzenWindowHistogramImageToImageMetric.h:99
itk::ParzenWindowHistogramImageToImageMetric::JointPDFPointer
JointPDFType::Pointer JointPDFPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:249
itk::ParzenWindowHistogramImageToImageMetric::UpdateJointPDFAndIncrementalPDFs
virtual void UpdateJointPDFAndIncrementalPDFs(RealType fixedImageValue, RealType movingImageValue, RealType movingMaskValue, const DerivativeType &movingImageValuesRight, const DerivativeType &movingImageValuesLeft, const DerivativeType &movingMaskValuesRight, const DerivativeType &movingMaskValuesLeft, const NonZeroJacobianIndicesType &nzji) const
itk::ParzenWindowHistogramImageToImageMetric::DerivativeType
Superclass::DerivativeType DerivativeType
Definition: itkParzenWindowHistogramImageToImageMetric.h:115
itkKernelFunctionBase2.h
itk::ParzenWindowHistogramImageToImageMetric::InputPointType
Superclass::InputPointType InputPointType
Definition: itkParzenWindowHistogramImageToImageMetric.h:98
itk::AdvancedImageToImageMetric::FixedImageRegionType
Superclass::FixedImageRegionType FixedImageRegionType
Definition: itkAdvancedImageToImageMetric.h:108
itk::ParzenWindowHistogramImageToImageMetric::JointPDFType
Image< PDFValueType, 2 > JointPDFType
Definition: itkParzenWindowHistogramImageToImageMetric.h:248
itk::AdvancedImageToImageMetric::ImageSamplerPointer
ImageSamplerType::Pointer ImageSamplerPointer
Definition: itkAdvancedImageToImageMetric.h:139
itk::ParzenWindowHistogramImageToImageMetric::GetValueAndAnalyticDerivative
virtual void GetValueAndAnalyticDerivative(const ParametersType &, MeasureType &, DerivativeType &) const
Definition: itkParzenWindowHistogramImageToImageMetric.h:470
itk::ParzenWindowHistogramImageToImageMetric::m_NumberOfFixedHistogramBins
unsigned long m_NumberOfFixedHistogramBins
Definition: itkParzenWindowHistogramImageToImageMetric.h:492
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramImageToImageMetric
ParzenWindowHistogramImageToImageMetric()
itk::ParzenWindowHistogramImageToImageMetric::m_MovingImageMarginalPDF
MarginalPDFType m_MovingImageMarginalPDF
Definition: itkParzenWindowHistogramImageToImageMetric.h:278
itk::ParzenWindowHistogramImageToImageMetric::TransformParametersType
Superclass::TransformParametersType TransformParametersType
Definition: itkParzenWindowHistogramImageToImageMetric.h:100
itk::ParzenWindowHistogramImageToImageMetric::MovingImageLimiterOutputType
Superclass::MovingImageLimiterOutputType MovingImageLimiterOutputType
Definition: itkParzenWindowHistogramImageToImageMetric.h:127
itk::ParzenWindowHistogramImageToImageMetric::ComputePDFsAndPDFDerivatives
virtual void ComputePDFsAndPDFDerivatives(const ParametersType &parameters) const
itk::ParzenWindowHistogramImageToImageMetric::GetValueAndDerivative
void GetValueAndDerivative(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const
itk::AdvancedImageToImageMetric::GradientImageType
Superclass::GradientImageType GradientImageType
Definition: itkAdvancedImageToImageMetric.h:119
itk::ParzenWindowHistogramImageToImageMetric::m_MovingParzenTermToIndexOffset
double m_MovingParzenTermToIndexOffset
Definition: itkParzenWindowHistogramImageToImageMetric.h:293
SmartPointer< Self >
itk::ParzenWindowHistogramImageToImageMetric::JointPDFDerivativesIndexType
JointPDFDerivativesType::IndexType JointPDFDerivativesIndexType
Definition: itkParzenWindowHistogramImageToImageMetric.h:257
itk::AdvancedImageToImageMetric::DerivativeValueType
DerivativeType::ValueType DerivativeValueType
Definition: itkAdvancedImageToImageMetric.h:129
itk::ParzenWindowHistogramImageToImageMetric::m_FiniteDifferencePerturbation
double m_FiniteDifferencePerturbation
Definition: itkParzenWindowHistogramImageToImageMetric.h:499
itk::ParzenWindowHistogramImageToImageMetric::ParametersType
Superclass::ParametersType ParametersType
Definition: itkParzenWindowHistogramImageToImageMetric.h:117
itk::AdvancedImageToImageMetric
An extension of the ITK ImageToImageMetric. It is the intended base class for all elastix metrics.
Definition: itkAdvancedImageToImageMetric.h:81
itk::ParzenWindowHistogramImageToImageMetric::GetValueAndFiniteDifferenceDerivative
virtual void GetValueAndFiniteDifferenceDerivative(const ParametersType &, MeasureType &, DerivativeType &) const
Definition: itkParzenWindowHistogramImageToImageMetric.h:479
itk::ParzenWindowHistogramImageToImageMetric::KernelFunctionType
KernelFunctionBase2< PDFValueType > KernelFunctionType
Definition: itkParzenWindowHistogramImageToImageMetric.h:266
itk::ParzenWindowHistogramImageToImageMetric::FixedImageIndexType
Superclass::FixedImageIndexType FixedImageIndexType
Definition: itkParzenWindowHistogramImageToImageMetric.h:232
itk::ParzenWindowHistogramImageToImageMetric::MovingImageType
Superclass::MovingImageType MovingImageType
Definition: itkParzenWindowHistogramImageToImageMetric.h:90
itk::AdvancedImageToImageMetric::DerivativeType
Superclass::DerivativeType DerivativeType
Definition: itkAdvancedImageToImageMetric.h:128
itk::ParzenWindowHistogramImageToImageMetric::m_MovingIncrementalMarginalPDFLeft
IncrementalMarginalPDFPointer m_MovingIncrementalMarginalPDFLeft
Definition: itkParzenWindowHistogramImageToImageMetric.h:286
itk::ParzenWindowHistogramImageToImageMetric::ComputePDFs
virtual void ComputePDFs(const ParametersType &parameters) const
itk::ParzenWindowHistogramImageToImageMetric::m_UseDerivative
bool m_UseDerivative
Definition: itkParzenWindowHistogramImageToImageMetric.h:496
itk::ParzenWindowHistogramImageToImageMetric::IncrementalMarginalPDFRegionType
IncrementalMarginalPDFType::RegionType IncrementalMarginalPDFRegionType
Definition: itkParzenWindowHistogramImageToImageMetric.h:261
itk::ParzenWindowHistogramImageToImageMetric::InterpolatorPointer
Superclass::InterpolatorPointer InterpolatorPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:103
itk::ParzenWindowHistogramImageToImageMetric::JointPDFSizeType
JointPDFType::SizeType JointPDFSizeType
Definition: itkParzenWindowHistogramImageToImageMetric.h:256
itk::ParzenWindowHistogramImageToImageMetric::m_UseFiniteDifferenceDerivative
bool m_UseFiniteDifferenceDerivative
Definition: itkParzenWindowHistogramImageToImageMetric.h:498
Array< PDFValueType >
itk::ParzenWindowHistogramImageToImageMetric::m_JointPDF
JointPDFPointer m_JointPDF
Definition: itkParzenWindowHistogramImageToImageMetric.h:279
itk::ParzenWindowHistogramImageToImageMetric::m_DerivativeMovingKernel
KernelFunctionPointer m_DerivativeMovingKernel
Definition: itkParzenWindowHistogramImageToImageMetric.h:298
itk::ParzenWindowHistogramImageToImageMetric::GradientImageType
Superclass::GradientImageType GradientImageType
Definition: itkParzenWindowHistogramImageToImageMetric.h:106
itk::AdvancedImageToImageMetric::CoordinateRepresentationType
Superclass::CoordinateRepresentationType CoordinateRepresentationType
Definition: itkAdvancedImageToImageMetric.h:100
itk::ParzenWindowHistogramImageToImageMetric::PDFValueType
double PDFValueType
Definition: itkParzenWindowHistogramImageToImageMetric.h:245
itk::ParzenWindowHistogramImageToImageMetric::MovingImageIndexType
Superclass::MovingImageIndexType MovingImageIndexType
Definition: itkParzenWindowHistogramImageToImageMetric.h:235
itk::ParzenWindowHistogramImageToImageMetric::m_PerturbedAlphaLeft
DerivativeType m_PerturbedAlphaLeft
Definition: itkParzenWindowHistogramImageToImageMetric.h:274
itk::ParzenWindowHistogramImageToImageMetric::ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:83
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramGetValueAndDerivativePerThreadStruct::st_JointPDF
JointPDFPointer st_JointPDF
Definition: itkParzenWindowHistogramImageToImageMetric.h:315
itk::ParzenWindowHistogramImageToImageMetric::InitializeKernels
virtual void InitializeKernels(void)
itk::ParzenWindowHistogramImageToImageMetric::MovingImagePointType
Superclass::MovingImagePointType MovingImagePointType
Definition: itkParzenWindowHistogramImageToImageMetric.h:237
itk::AdvancedImageToImageMetric::MovingImageDerivativeScalesType
FixedArray< double, Self::MovingImageDimension > MovingImageDerivativeScalesType
Definition: itkAdvancedImageToImageMetric.h:135
itk::AdvancedImageToImageMetric::MovingImageLimiterOutputType
MovingImageLimiterType::OutputType MovingImageLimiterOutputType
Definition: itkAdvancedImageToImageMetric.h:149
itk::AdvancedImageToImageMetric::MovingImageDerivativeType
BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
Definition: itkAdvancedImageToImageMetric.h:323
Image
itk::ParzenWindowHistogramImageToImageMetric::ParzenValueContainerType
Array< PDFValueType > ParzenValueContainerType
Definition: itkParzenWindowHistogramImageToImageMetric.h:263
itk::ParzenWindowHistogramImageToImageMetric::m_FixedImageMarginalPDF
MarginalPDFType m_FixedImageMarginalPDF
Definition: itkParzenWindowHistogramImageToImageMetric.h:277
itk::ParzenWindowHistogramImageToImageMetric::m_FixedIncrementalMarginalPDFRight
IncrementalMarginalPDFPointer m_FixedIncrementalMarginalPDFRight
Definition: itkParzenWindowHistogramImageToImageMetric.h:283
itk::ParzenWindowHistogramImageToImageMetric::IncrementalMarginalPDFIndexType
IncrementalMarginalPDFType::IndexType IncrementalMarginalPDFIndexType
Definition: itkParzenWindowHistogramImageToImageMetric.h:260
itk::ParzenWindowHistogramImageToImageMetric::JointPDFDerivativesType
Image< PDFDerivativeValueType, 3 > JointPDFDerivativesType
Definition: itkParzenWindowHistogramImageToImageMetric.h:250
itk::ParzenWindowHistogramImageToImageMetric::IncrementalMarginalPDFPointer
IncrementalMarginalPDFType::Pointer IncrementalMarginalPDFPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:253
itk::ParzenWindowHistogramImageToImageMetric::GradientImagePointer
Superclass::GradientImagePointer GradientImagePointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:107
itk::ParzenWindowHistogramImageToImageMetric::NormalizeJointPDF
virtual void NormalizeJointPDF(JointPDFType *pdf, const double &factor) const
itk::ParzenWindowHistogramImageToImageMetric::UpdateJointPDFDerivatives
void UpdateJointPDFDerivatives(const JointPDFIndexType &pdfIndex, double factor, const DerivativeType &imageJacobian, const NonZeroJacobianIndicesType &nzji) const
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramGetValueAndDerivativePerThreadStruct
Definition: itkParzenWindowHistogramImageToImageMetric.h:313
itk::ParzenWindowHistogramImageToImageMetric::m_FixedParzenTermToIndexOffset
double m_FixedParzenTermToIndexOffset
Definition: itkParzenWindowHistogramImageToImageMetric.h:292
itk::AdvancedImageToImageMetric::OutputPointType
Superclass::OutputPointType OutputPointType
Definition: itkAdvancedImageToImageMetric.h:112
itk::AdvancedImageToImageMetric::ParametersType
Superclass::ParametersType ParametersType
Definition: itkAdvancedImageToImageMetric.h:130
itk::AdvancedImageToImageMetric::MovingImageLimiterType
LimiterFunctionBase< RealType, MovingImageDimension > MovingImageLimiterType
Definition: itkAdvancedImageToImageMetric.h:147
itk::ParzenWindowHistogramImageToImageMetric::m_NumberOfMovingHistogramBins
unsigned long m_NumberOfMovingHistogramBins
Definition: itkParzenWindowHistogramImageToImageMetric.h:493
itk::ParzenWindowHistogramImageToImageMetric::itkStaticConstMacro
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
itk::AdvancedImageToImageMetric::GradientImageFilterType
Superclass::GradientImageFilterType GradientImageFilterType
Definition: itkAdvancedImageToImageMetric.h:121
itk::ParzenWindowHistogramImageToImageMetric::ThreaderType
Superclass::ThreaderType ThreaderType
Definition: itkParzenWindowHistogramImageToImageMetric.h:129
itk::AdvancedImageToImageMetric::InterpolatorPointer
Superclass::InterpolatorPointer InterpolatorPointer
Definition: itkAdvancedImageToImageMetric.h:116
itk::AdvancedImageToImageMetric::GradientImagePointer
Superclass::GradientImagePointer GradientImagePointer
Definition: itkAdvancedImageToImageMetric.h:120
itk::ParzenWindowHistogramImageToImageMetric::m_MovingIncrementalMarginalPDFRight
IncrementalMarginalPDFPointer m_MovingIncrementalMarginalPDFRight
Definition: itkParzenWindowHistogramImageToImageMetric.h:284
itk::ParzenWindowHistogramImageToImageMetric::m_IncrementalJointPDFLeft
JointPDFDerivativesPointer m_IncrementalJointPDFLeft
Definition: itkParzenWindowHistogramImageToImageMetric.h:282
itk::AdvancedImageToImageMetric::ImageSampleContainerType
ImageSamplerType::OutputVectorContainerType ImageSampleContainerType
Definition: itkAdvancedImageToImageMetric.h:140
itk::AdvancedImageToImageMetric::InterpolatorType
Superclass::InterpolatorType InterpolatorType
Definition: itkAdvancedImageToImageMetric.h:115
itk::ParzenWindowHistogramImageToImageMetric::m_IncrementalJointPDFRight
JointPDFDerivativesPointer m_IncrementalJointPDFRight
Definition: itkParzenWindowHistogramImageToImageMetric.h:281
itk::ParzenWindowHistogramImageToImageMetric::itkPadStruct
itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, ParzenWindowHistogramGetValueAndDerivativePerThreadStruct, PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct)
itk::AdvancedImageToImageMetric::NonZeroJacobianIndicesType
AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
Definition: itkAdvancedImageToImageMetric.h:330
itk::AdvancedImageToImageMetric::TransformJacobianType
Superclass::TransformJacobianType TransformJacobianType
Definition: itkAdvancedImageToImageMetric.h:114
itk::ParzenWindowHistogramImageToImageMetric::IncrementalMarginalPDFSizeType
IncrementalMarginalPDFType::SizeType IncrementalMarginalPDFSizeType
Definition: itkParzenWindowHistogramImageToImageMetric.h:262
itk::AdvancedImageToImageMetric::FixedImageMaskType
Superclass::FixedImageMaskType FixedImageMaskType
Definition: itkAdvancedImageToImageMetric.h:123
itk::ParzenWindowHistogramImageToImageMetric::ThreadInfoType
Superclass::ThreadInfoType ThreadInfoType
Definition: itkParzenWindowHistogramImageToImageMetric.h:130
itk::AdvancedImageToImageMetric::FixedImageIndexType
FixedImageType::IndexType FixedImageIndexType
Definition: itkAdvancedImageToImageMetric.h:303
itk::AdvancedImageToImageMetric::MovingImageRegionType
MovingImageType::RegionType MovingImageRegionType
Definition: itkAdvancedImageToImageMetric.h:134
itk::ParzenWindowHistogramImageToImageMetric::MovingImageContinuousIndexType
Superclass::MovingImageContinuousIndexType MovingImageContinuousIndexType
Definition: itkParzenWindowHistogramImageToImageMetric.h:238
itk::ParzenWindowHistogramImageToImageMetric::m_FixedIncrementalMarginalPDFLeft
IncrementalMarginalPDFPointer m_FixedIncrementalMarginalPDFLeft
Definition: itkParzenWindowHistogramImageToImageMetric.h:285
itk::ParzenWindowHistogramImageToImageMetric::ThreadedComputePDFs
void ThreadedComputePDFs(ThreadIdType threadId)
itk::ParzenWindowHistogramImageToImageMetric::CentralDifferenceGradientFilterType
Superclass::CentralDifferenceGradientFilterType CentralDifferenceGradientFilterType
Definition: itkParzenWindowHistogramImageToImageMetric.h:241
itk::ParzenWindowHistogramImageToImageMetric::operator=
void operator=(const Self &)
itk::AdvancedImageToImageMetric::MovingImageMaskType
Superclass::MovingImageMaskType MovingImageMaskType
Definition: itkAdvancedImageToImageMetric.h:125
itk::ParzenWindowHistogramImageToImageMetric::JointPDFDerivativesRegionType
JointPDFDerivativesType::RegionType JointPDFDerivativesRegionType
Definition: itkParzenWindowHistogramImageToImageMetric.h:258
ThreadIdType
itk::ParzenWindowHistogramImageToImageMetric::m_ParzenWindowHistogramThreaderParameters
ParzenWindowHistogramMultiThreaderParameterType m_ParzenWindowHistogramThreaderParameters
Definition: itkParzenWindowHistogramImageToImageMetric.h:310
itk::AdvancedImageToImageMetric::GradientImageFilterPointer
Superclass::GradientImageFilterPointer GradientImageFilterPointer
Definition: itkAdvancedImageToImageMetric.h:122
itk::AdvancedImageToImageMetric::FixedImageType
Superclass::FixedImageType FixedImageType
Definition: itkAdvancedImageToImageMetric.h:105
itk::AdvancedImageToImageMetric::ThreaderType
itk::MultiThreader ThreaderType
Definition: itkAdvancedImageToImageMetric.h:171
itk::ParzenWindowHistogramImageToImageMetric::AfterThreadedComputePDFs
void AfterThreadedComputePDFs(void) const
itk::ParzenWindowHistogramImageToImageMetric::Self
ParzenWindowHistogramImageToImageMetric Self
Definition: itkParzenWindowHistogramImageToImageMetric.h:80
itk::ParzenWindowHistogramImageToImageMetric::m_MovingImageBinSize
double m_MovingImageBinSize
Definition: itkParzenWindowHistogramImageToImageMetric.h:291
itk::ParzenWindowHistogramImageToImageMetric::IncrementalMarginalPDFType
Image< PDFValueType, 2 > IncrementalMarginalPDFType
Definition: itkParzenWindowHistogramImageToImageMetric.h:252
itk::AdvancedImageToImageMetric::MovingImagePixelType
Superclass::MovingImagePixelType MovingImagePixelType
Definition: itkAdvancedImageToImageMetric.h:102
itk::ParzenWindowHistogramImageToImageMetric::ImageSamplerType
Superclass::ImageSamplerType ImageSamplerType
Definition: itkParzenWindowHistogramImageToImageMetric.h:120
itk::ParzenWindowHistogramImageToImageMetric::MovingImageDerivativeScalesType
Superclass::MovingImageDerivativeScalesType MovingImageDerivativeScalesType
Definition: itkParzenWindowHistogramImageToImageMetric.h:128
itk::AdvancedImageToImageMetric::FixedImagePixelType
FixedImageType::PixelType FixedImagePixelType
Definition: itkAdvancedImageToImageMetric.h:133
itk::ParzenWindowHistogramImageToImageMetric::JointPDFIndexType
JointPDFType::IndexType JointPDFIndexType
Definition: itkParzenWindowHistogramImageToImageMetric.h:254
itk::ParzenWindowHistogramImageToImageMetric::m_ParzenWindowHistogramGetValueAndDerivativePerThreadVariables
AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct * m_ParzenWindowHistogramGetValueAndDerivativePerThreadVariables
Definition: itkParzenWindowHistogramImageToImageMetric.h:321
itk::ParzenWindowHistogramImageToImageMetric::JointPDFDerivativesSizeType
JointPDFDerivativesType::SizeType JointPDFDerivativesSizeType
Definition: itkParzenWindowHistogramImageToImageMetric.h:259
itk::ParzenWindowHistogramImageToImageMetric::Superclass
AdvancedImageToImageMetric< TFixedImage, TMovingImage > Superclass
Definition: itkParzenWindowHistogramImageToImageMetric.h:81
itk::ParzenWindowHistogramImageToImageMetric::MovingImageConstPointer
Superclass::MovingImageConstPointer MovingImageConstPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:92
itk::ParzenWindowHistogramImageToImageMetric::RealType
Superclass::RealType RealType
Definition: itkParzenWindowHistogramImageToImageMetric.h:104
itk::ParzenWindowHistogramImageToImageMetric::FixedImageMaskType
Superclass::FixedImageMaskType FixedImageMaskType
Definition: itkParzenWindowHistogramImageToImageMetric.h:110
itk::AdvancedImageToImageMetric::FixedImagePointType
TransformType::InputPointType FixedImagePointType
Definition: itkAdvancedImageToImageMetric.h:306
itk::AdvancedImageToImageMetric::FixedImageLimiterOutputType
FixedImageLimiterType::OutputType FixedImageLimiterOutputType
Definition: itkAdvancedImageToImageMetric.h:146
itk::ParzenWindowHistogramImageToImageMetric::m_MovingImageNormalizedMin
double m_MovingImageNormalizedMin
Definition: itkParzenWindowHistogramImageToImageMetric.h:288
itk::ParzenWindowHistogramImageToImageMetric::MovingImageMaskPointer
Superclass::MovingImageMaskPointer MovingImageMaskPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:113
itk::ParzenWindowHistogramImageToImageMetric::GradientPixelType
Superclass::GradientPixelType GradientPixelType
Definition: itkParzenWindowHistogramImageToImageMetric.h:105
itk::ParzenWindowHistogramImageToImageMetric::FixedImageType
Superclass::FixedImageType FixedImageType
Definition: itkParzenWindowHistogramImageToImageMetric.h:93
itk::ParzenWindowHistogramImageToImageMetric::MovingImageMaskType
Superclass::MovingImageMaskType MovingImageMaskType
Definition: itkParzenWindowHistogramImageToImageMetric.h:112
itk::ParzenWindowHistogramImageToImageMetric::FixedImageMaskPointer
Superclass::FixedImageMaskPointer FixedImageMaskPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:111
itk::ParzenWindowHistogramImageToImageMetric::ComputeMarginalPDF
virtual void ComputeMarginalPDF(const JointPDFType *jointPDF, MarginalPDFType &marginalPDF, const unsigned int &direction) const
itk::ParzenWindowHistogramImageToImageMetric::BSplineInterpolatorType
Superclass::BSplineInterpolatorType BSplineInterpolatorType
Definition: itkParzenWindowHistogramImageToImageMetric.h:239
itk::ParzenWindowHistogramImageToImageMetric::m_FixedKernel
KernelFunctionPointer m_FixedKernel
Definition: itkParzenWindowHistogramImageToImageMetric.h:296
itk::ParzenWindowHistogramImageToImageMetric::OffsetValueType
FixedImageType::OffsetValueType OffsetValueType
Definition: itkParzenWindowHistogramImageToImageMetric.h:234
itk::ParzenWindowHistogramImageToImageMetric::ComputePDFsAndIncrementalPDFs
virtual void ComputePDFsAndIncrementalPDFs(const ParametersType &parameters) const
itk::ParzenWindowHistogramImageToImageMetric::m_ThreaderJointPDFs
std::vector< JointPDFPointer > m_ThreaderJointPDFs
Definition: itkParzenWindowHistogramImageToImageMetric.h:301
itk::ParzenWindowHistogramImageToImageMetric::m_FixedImageNormalizedMin
double m_FixedImageNormalizedMin
Definition: itkParzenWindowHistogramImageToImageMetric.h:289
itk::ParzenWindowHistogramImageToImageMetric::TransformType
Superclass::TransformType TransformType
Definition: itkParzenWindowHistogramImageToImageMetric.h:96
itk::ParzenWindowHistogramImageToImageMetric::m_Alpha
double m_Alpha
Definition: itkParzenWindowHistogramImageToImageMetric.h:272
itk::ParzenWindowHistogramImageToImageMetric::GradientImageFilterPointer
Superclass::GradientImageFilterPointer GradientImageFilterPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:109
itk::AdvancedImageToImageMetric::ImageSampleContainerPointer
ImageSamplerType::OutputVectorContainerPointer ImageSampleContainerPointer
Definition: itkAdvancedImageToImageMetric.h:141
itk::ParzenWindowHistogramImageToImageMetric::ComputePDFsThreaderCallback
static ITK_THREAD_RETURN_TYPE ComputePDFsThreaderCallback(void *arg)
itk::AdvancedImageToImageMetric::FixedImageLimiterType
LimiterFunctionBase< RealType, FixedImageDimension > FixedImageLimiterType
Definition: itkAdvancedImageToImageMetric.h:144
itk::AdvancedImageToImageMetric::RealType
Superclass::RealType RealType
Definition: itkAdvancedImageToImageMetric.h:117
itk::ParzenWindowHistogramImageToImageMetric::MeasureType
Superclass::MeasureType MeasureType
Definition: itkParzenWindowHistogramImageToImageMetric.h:114
itk::ParzenWindowHistogramImageToImageMetric::FixedImageConstPointer
Superclass::FixedImageConstPointer FixedImageConstPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:94
itk::ParzenWindowHistogramImageToImageMetric::CoordinateRepresentationType
Superclass::CoordinateRepresentationType CoordinateRepresentationType
Definition: itkParzenWindowHistogramImageToImageMetric.h:86
itk::ParzenWindowHistogramImageToImageMetric::MovingImageLimiterType
Superclass::MovingImageLimiterType MovingImageLimiterType
Definition: itkParzenWindowHistogramImageToImageMetric.h:125
itk::ParzenWindowHistogramImageToImageMetric::m_ParzenWindowHistogramGetValueAndDerivativePerThreadVariablesSize
ThreadIdType m_ParzenWindowHistogramGetValueAndDerivativePerThreadVariablesSize
Definition: itkParzenWindowHistogramImageToImageMetric.h:322
itk::ParzenWindowHistogramImageToImageMetric::EvaluateParzenValues
void EvaluateParzenValues(double parzenWindowTerm, OffsetValueType parzenWindowIndex, const KernelFunctionType *kernel, ParzenValueContainerType &parzenValues) const
itk::KernelFunctionBase2
Kernel used for density estimation and nonparameteric regression.
Definition: itkKernelFunctionBase2.h:41
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramMultiThreaderParameterType
Definition: itkParzenWindowHistogramImageToImageMetric.h:307
itk::ParzenWindowHistogramImageToImageMetric::JointPDFRegionType
JointPDFType::RegionType JointPDFRegionType
Definition: itkParzenWindowHistogramImageToImageMetric.h:255
itk::ParzenWindowHistogramImageToImageMetric::ImageSampleContainerPointer
Superclass::ImageSampleContainerPointer ImageSampleContainerPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:123
itk::ParzenWindowHistogramImageToImageMetric::m_MovingKernel
KernelFunctionPointer m_MovingKernel
Definition: itkParzenWindowHistogramImageToImageMetric.h:297
itk::ParzenWindowHistogramImageToImageMetric::TransformJacobianType
Superclass::TransformJacobianType TransformJacobianType
Definition: itkParzenWindowHistogramImageToImageMetric.h:101
itk::ParzenWindowHistogramImageToImageMetric::MarginalPDFType
Array< PDFValueType > MarginalPDFType
Definition: itkParzenWindowHistogramImageToImageMetric.h:247
itk::ParzenWindowHistogramImageToImageMetric::m_FixedImageBinSize
double m_FixedImageBinSize
Definition: itkParzenWindowHistogramImageToImageMetric.h:290
itk::AdvancedImageToImageMetric::MovingImageContinuousIndexType
InterpolatorType::ContinuousIndexType MovingImageContinuousIndexType
Definition: itkAdvancedImageToImageMetric.h:308
itk::AdvancedImageToImageMetric::MeasureType
Superclass::MeasureType MeasureType
Definition: itkAdvancedImageToImageMetric.h:127
itk::ParzenWindowHistogramImageToImageMetric::itkAlignedTypedef
itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct, AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct)
itk::ParzenWindowHistogramImageToImageMetric::GetDerivative
void GetDerivative(const ParametersType &parameters, DerivativeType &Derivative) const
itk::ParzenWindowHistogramImageToImageMetric::ImageSampleContainerType
Superclass::ImageSampleContainerType ImageSampleContainerType
Definition: itkParzenWindowHistogramImageToImageMetric.h:122
itk
Definition: itkAdvancedImageToImageMetric.h:40
itk::ParzenWindowHistogramImageToImageMetric::Initialize
void Initialize(void)
itk::ParzenWindowHistogramImageToImageMetric::ImageSamplerPointer
Superclass::ImageSamplerPointer ImageSamplerPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:121
itk::ParzenWindowHistogramImageToImageMetric::DerivativeValueType
Superclass::DerivativeValueType DerivativeValueType
Definition: itkParzenWindowHistogramImageToImageMetric.h:116
itk::ParzenWindowHistogramImageToImageMetric::FixedImageLimiterType
Superclass::FixedImageLimiterType FixedImageLimiterType
Definition: itkParzenWindowHistogramImageToImageMetric.h:124
itk::AdvancedImageToImageMetric::FixedImageMaskPointer
Superclass::FixedImageMaskPointer FixedImageMaskPointer
Definition: itkAdvancedImageToImageMetric.h:124
itk::ParzenWindowHistogramImageToImageMetric::m_UseExplicitPDFDerivatives
bool m_UseExplicitPDFDerivatives
Definition: itkParzenWindowHistogramImageToImageMetric.h:497
itk::ParzenWindowHistogramImageToImageMetric::m_PerturbedAlphaRight
DerivativeType m_PerturbedAlphaRight
Definition: itkParzenWindowHistogramImageToImageMetric.h:273
itk::ParzenWindowHistogramImageToImageMetric::m_MovingKernelBSplineOrder
unsigned int m_MovingKernelBSplineOrder
Definition: itkParzenWindowHistogramImageToImageMetric.h:495
itk::ParzenWindowHistogramImageToImageMetric::InitializeThreadingParameters
virtual void InitializeThreadingParameters(void) const
itk::AdvancedImageToImageMetric::MovingImagePointType
TransformType::OutputPointType MovingImagePointType
Definition: itkAdvancedImageToImageMetric.h:307
itk::AdvancedImageToImageMetric::MovingImageMaskPointer
Superclass::MovingImageMaskPointer MovingImageMaskPointer
Definition: itkAdvancedImageToImageMetric.h:126
itk::ParzenWindowHistogramImageToImageMetric::FixedImagePointType
Superclass::FixedImagePointType FixedImagePointType
Definition: itkParzenWindowHistogramImageToImageMetric.h:236
itk::ParzenWindowHistogramImageToImageMetric::m_JointPDFDerivatives
JointPDFDerivativesPointer m_JointPDFDerivatives
Definition: itkParzenWindowHistogramImageToImageMetric.h:280
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramGetValueAndDerivativePerThreadStruct::st_NumberOfPixelsCounted
SizeValueType st_NumberOfPixelsCounted
Definition: itkParzenWindowHistogramImageToImageMetric.h:314
itk::ParzenWindowHistogramImageToImageMetric::ComputeIncrementalMarginalPDFs
virtual void ComputeIncrementalMarginalPDFs(const JointPDFDerivativesType *incrementalPDF, IncrementalMarginalPDFType *fixedIncrementalMarginalPDF, IncrementalMarginalPDFType *movingIncrementalMarginalPDF) const
itk::ParzenWindowHistogramImageToImageMetric::LaunchComputePDFsThreaderCallback
void LaunchComputePDFsThreaderCallback(void) const
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramImageToImageMetric
ParzenWindowHistogramImageToImageMetric(const Self &)
itk::ParzenWindowHistogramImageToImageMetric::KernelFunctionPointer
KernelFunctionType::Pointer KernelFunctionPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:267
itk::AdvancedImageToImageMetric::MovingImageType
Superclass::MovingImageType MovingImageType
Definition: itkAdvancedImageToImageMetric.h:101
itk::ParzenWindowHistogramImageToImageMetric
A base class for image metrics based on a joint histogram computed using Parzen Windowing.
Definition: itkParzenWindowHistogramImageToImageMetric.h:76
itk::ParzenWindowHistogramImageToImageMetric::MovingImagePixelType
Superclass::MovingImagePixelType MovingImagePixelType
Definition: itkParzenWindowHistogramImageToImageMetric.h:91
itk::AdvancedImageToImageMetric::BSplineInterpolatorType
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, double > BSplineInterpolatorType
Definition: itkAdvancedImageToImageMetric.h:312
itk::ParzenWindowHistogramImageToImageMetric::FixedImageRegionType
Superclass::FixedImageRegionType FixedImageRegionType
Definition: itkParzenWindowHistogramImageToImageMetric.h:95
itk::ParzenWindowHistogramImageToImageMetric::ComputePDFsSingleThreaded
virtual void ComputePDFsSingleThreaded(const ParametersType &parameters) const
itk::AdvancedImageToImageMetric::GradientPixelType
Superclass::GradientPixelType GradientPixelType
Definition: itkAdvancedImageToImageMetric.h:118
itk::AdvancedImageToImageMetric::CentralDifferenceGradientFilterType
GradientImageFilter< MovingImageType, RealType, RealType > CentralDifferenceGradientFilterType
Definition: itkAdvancedImageToImageMetric.h:325
itk::ParzenWindowHistogramImageToImageMetric::m_JointPDFWindow
JointPDFRegionType m_JointPDFWindow
Definition: itkParzenWindowHistogramImageToImageMetric.h:287
itk::ParzenWindowHistogramImageToImageMetric::MovingImageRegionType
Superclass::MovingImageRegionType MovingImageRegionType
Definition: itkParzenWindowHistogramImageToImageMetric.h:119
itk::ParzenWindowHistogramImageToImageMetric::GradientImageFilterType
Superclass::GradientImageFilterType GradientImageFilterType
Definition: itkParzenWindowHistogramImageToImageMetric.h:108
itk::AdvancedImageToImageMetric::TransformPointer
Superclass::TransformPointer TransformPointer
Definition: itkAdvancedImageToImageMetric.h:110
itk::AdvancedImageToImageMetric::TransformType
Superclass::TransformType TransformType
Definition: itkAdvancedImageToImageMetric.h:109
itk::AdvancedImageToImageMetric::MovingImageIndexType
MovingImageType::IndexType MovingImageIndexType
Definition: itkAdvancedImageToImageMetric.h:305
itk::ParzenWindowHistogramImageToImageMetric::~ParzenWindowHistogramImageToImageMetric
virtual ~ParzenWindowHistogramImageToImageMetric()
itk::ParzenWindowHistogramImageToImageMetric::JointPDFDerivativesPointer
JointPDFDerivativesType::Pointer JointPDFDerivativesPointer
Definition: itkParzenWindowHistogramImageToImageMetric.h:251
itk::ParzenWindowHistogramImageToImageMetric::PDFDerivativeValueType
float PDFDerivativeValueType
Definition: itkParzenWindowHistogramImageToImageMetric.h:246
itk::ParzenWindowHistogramImageToImageMetric::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const
itk::AdvancedImageToImageMetric::FixedImageIndexValueType
FixedImageIndexType::IndexValueType FixedImageIndexValueType
Definition: itkAdvancedImageToImageMetric.h:304
itk::ParzenWindowHistogramImageToImageMetric::NonZeroJacobianIndicesType
Superclass::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
Definition: itkParzenWindowHistogramImageToImageMetric.h:242
itkAdvancedImageToImageMetric.h
itk::ParzenWindowHistogramImageToImageMetric::UpdateJointPDFAndDerivatives
virtual void UpdateJointPDFAndDerivatives(const RealType &fixedImageValue, const RealType &movingImageValue, const DerivativeType *imageJacobian, const NonZeroJacobianIndicesType *nzji, JointPDFType *jointPDF) const
itk::AdvancedImageToImageMetric::FixedImageConstPointer
Superclass::FixedImageConstPointer FixedImageConstPointer
Definition: itkAdvancedImageToImageMetric.h:107
itk::AdvancedImageToImageMetric::MovingImageConstPointer
Superclass::MovingImageConstPointer MovingImageConstPointer
Definition: itkAdvancedImageToImageMetric.h:104
itk::ParzenWindowHistogramImageToImageMetric::MovingImageDerivativeType
Superclass::MovingImageDerivativeType MovingImageDerivativeType
Definition: itkParzenWindowHistogramImageToImageMetric.h:240
itk::AdvancedImageToImageMetric::InputPointType
Superclass::InputPointType InputPointType
Definition: itkAdvancedImageToImageMetric.h:111
itk::AdvancedImageToImageMetric::ImageSamplerType
ImageSamplerBase< FixedImageType > ImageSamplerType
Definition: itkAdvancedImageToImageMetric.h:138
itk::ParzenWindowHistogramImageToImageMetric::InitializeHistograms
virtual void InitializeHistograms(void)
itk::ParzenWindowHistogramImageToImageMetric::ParzenWindowHistogramMultiThreaderParameterType::m_Metric
Self * m_Metric
Definition: itkParzenWindowHistogramImageToImageMetric.h:308
itk::ParzenWindowHistogramImageToImageMetric::FixedImageIndexValueType
Superclass::FixedImageIndexValueType FixedImageIndexValueType
Definition: itkParzenWindowHistogramImageToImageMetric.h:233
itk::ParzenWindowHistogramImageToImageMetric::itkStaticConstMacro
itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension)
itk::ParzenWindowHistogramImageToImageMetric::m_FixedKernelBSplineOrder
unsigned int m_FixedKernelBSplineOrder
Definition: itkParzenWindowHistogramImageToImageMetric.h:494
itk::AdvancedImageToImageMetric::ThreadInfoType
ThreaderType::ThreadInfoStruct ThreadInfoType
Definition: itkAdvancedImageToImageMetric.h:172


Generated on OURCE_DATE_EPOCH for elastix by doxygen 1.8.18 elastix logo