]> Creatis software - clitk.git/blob - registration/clitkGenericMetric.h
Moved from repository clitk to clitk.private/tests_dav
[clitk.git] / registration / clitkGenericMetric.h
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to:
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://www.centreleonberard.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17 ===========================================================================**/
18
19 #ifndef __clitkGenericMetric_h
20 #define __clitkGenericMetric_h
21
22 //clitk include
23 #include "clitkNormalizedCorrelationImageToImageMetric.h"
24 #include "clitkCorrelationRatioImageToImageMetric.h"
25 #include "itkMeanSquaresImageToImageMetricFor3DBLUTFFD.h"
26 #include "itkMattesMutualInformationImageToImageMetricFor3DBLUTFFD.h"
27 #include "clitkNormalizedCorrelationImageToImageMetricFor3DBLUTFFD.h"
28
29 //itk include
30 #include "itkSpatialObject.h"
31 #include "itkNormalizeImageFilter.h"
32 #include "itkImageToImageMetric.h"
33 #include "itkMeanSquaresImageToImageMetric.h"
34 #include "itkCorrelationCoefficientHistogramImageToImageMetric.h"
35 #include "itkGradientDifferenceImageToImageMetric.h"
36 #include "itkMutualInformationImageToImageMetric.h"
37 #include "itkMutualInformationHistogramImageToImageMetric.h"
38 #include "itkMattesMutualInformationImageToImageMetric.h"
39 #include "itkNormalizedMutualInformationHistogramImageToImageMetric.h"
40
41
42
43 /*
44
45 Requires at least the following section is the .ggo file. Adapt the defaults to the application
46
47
48 section "Metric (optimized, threaded versions are available for *, compile ITK with REVIEW and OPT_REGISTRATION enabled. Further optimized versions ** for  BLUT FFD optimizing a !3D! vector field  )"
49
50 option "metric"         -       "Type: 0=SSD*, 1=Normalized CC*, 2=Histogram CC, 3=Gradient-Difference, 4=Viola-Wells MI, 5=Histogram MI, 6=Mattes' MI*, 7=Normalized MI, 8=CR, 9=SSD for BLUT FFD**, 10=CC for BLUT FFD**, 11=Mattes' MI for BLUT FFD**"       int     no      default="0"
51 option "samples"        -       "Specify fraction [0, 1] of samples of the reference image used for the metric (* only). Use high fraction for detailed images (eg. 0.2, 0.5), for smooth images 0.01 might be enough." float   no      default="1"
52 option "intThreshold"   -       "Fixed image samples intensity threshold (* only)"                                                      float   no
53 option "subtractMean"   -       "1: Subtract mean for NCC calculation (narrows optimal)"                                                flag    on
54 option "bins"           -       "2,5-8: Number of histogram bins"                                                                       int     no      default="50"
55 option "random"         -       "4,6: Samples should be taken randomly, otherwise uniformly"                                            flag    off
56 option "stdDev"         -       "4: specify the standard deviation in mm of the gaussian kernels for both PDF estimations"              float   no      default="0.4"
57 option "explicitPDFDerivatives" -       "6: Calculate PDF derivatives explicitly (rigid=true; FFD=false)"                               flag    on
58
59
60 The use will look something like
61
62 typedef clitk::GenericMetric<args_info_type, FixedImageType, MovingImageType> GenericMetricType;
63 typename GenericMetricType::Pointer genericMetric=GenericMetricType::New();
64 genericMetric->SetArgsInfo(m_ArgsInfo);
65 genericMetric->SetFixedImageRegion(fixedImage->GetLargestPossibleRegion());
66 genericMetric->SetFixedImage(fixedImage);
67 genericMetric->SetFixedImageMask(fixedImageMask);
68 typedef itk::ImageToImageMetric<FixedImageType, MovingImageType> MetricType;
69 typename MetricType::Pointer metric=genericMetric->GetMetricPointer();
70
71 */
72
73
74 namespace clitk
75 {
76
77 template <  class args_info_type, class FixedImageType,  class MovingImageType >
78 class GenericMetric : public itk::LightObject
79 {
80 public:
81   //==============================================
82   typedef GenericMetric     Self;
83   typedef itk::LightObject     Superclass;
84   typedef itk::SmartPointer<Self>            Pointer;
85   typedef itk::SmartPointer<const Self>      ConstPointer;
86
87   /** Method for creation through the object factory. */
88   itkNewMacro(Self);
89
90   /** Determine the image dimension. */
91   itkStaticConstMacro(FixedImageDimension, unsigned int,
92                       FixedImageType::ImageDimension );
93   itkStaticConstMacro(MovingImageDimension, unsigned int,
94                       MovingImageType::ImageDimension );
95
96   // Typedef
97   typedef itk::ImageToImageMetric<FixedImageType, MovingImageType> MetricType;
98   typedef typename MetricType::Pointer MetricPointer;
99   typedef typename FixedImageType::RegionType FixedImageRegionType;
100   typedef itk::SpatialObject<itkGetStaticConstMacro(FixedImageDimension)> FixedImageMaskType;
101   typedef typename FixedImageMaskType::Pointer MaskPointer;
102   typedef typename FixedImageType::PixelType FixedImagePixelType;
103   typedef typename FixedImageType::IndexType FixedImageIndexType;
104   typedef typename FixedImageType::PointType FixedImagePointType;
105
106   //==============================================
107   //Set members
108   void SetArgsInfo(args_info_type args_info) {
109     m_ArgsInfo= args_info;
110     m_Verbose=m_ArgsInfo.verbose_flag;
111   }
112   void SetFixedImageRegion(const FixedImageRegionType f) {
113     m_FixedImageRegion=f;
114     m_FixedImageRegionGiven=true;
115   }
116   void SetFixedImage(typename FixedImageType::Pointer f) {
117     m_FixedImage=f;
118   }
119   void SetFixedImageMask( const FixedImageMaskType* f) {
120     m_FixedImageMask=f;
121   }
122
123   //==============================================
124   //Get members
125   MetricPointer GetMetricPointer(void);
126   bool GetMaximize(void) {
127     return m_Maximize;
128   }
129
130
131   //==============================================
132 protected:
133   GenericMetric();
134   ~GenericMetric() {};
135
136 private:
137   args_info_type m_ArgsInfo;
138   MetricPointer m_Metric;
139   bool m_Maximize;
140   bool m_Verbose;
141   bool m_FixedImageRegionGiven;
142   FixedImageRegionType m_FixedImageRegion;
143   typename FixedImageType::Pointer m_FixedImage;
144   typename FixedImageMaskType::ConstPointer m_FixedImageMask;
145
146 #if defined(ITK_USE_OPTIMIZED_REGISTRATION_METHODS) || ITK_VERSION_MAJOR >= 4
147   FixedImagePixelType m_FixedImageSamplesIntensityThreshold;
148   bool m_UseFixedImageSamplesIntensityThreshold;
149 #endif
150
151 };
152
153 } // end namespace clitk
154 #ifndef ITK_MANUAL_INSTANTIATION
155 #include "clitkGenericMetric.txx"
156 #endif
157
158 #endif // #define __clitkGenericMetric_h