ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkANTSNeighborhoodCorrelationImageToImageMetricv4.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 itkANTSNeighborhoodCorrelationImageToImageMetricv4_h
19 #define itkANTSNeighborhoodCorrelationImageToImageMetricv4_h
20 
23 
24 namespace itk {
25 
93 template<typename TFixedImage, typename TMovingImage, typename TVirtualImage = TFixedImage,
94  typename TInternalComputationValueType = double,
95  typename TMetricTraits = DefaultImageToImageMetricTraitsv4<TFixedImage,TMovingImage,TVirtualImage,TInternalComputationValueType>
96  >
98  public ImageToImageMetricv4<TFixedImage, TMovingImage, TVirtualImage, TInternalComputationValueType, TMetricTraits>
99 {
100 public:
101 
104  typedef ImageToImageMetricv4<TFixedImage, TMovingImage, TVirtualImage,
105  TInternalComputationValueType,TMetricTraits> Superclass;
108 
110  itkNewMacro(Self);
111 
113  itkTypeMacro(Self, Superclass);
114 
116  typedef typename Superclass::MeasureType MeasureType;
117  typedef typename Superclass::DerivativeType DerivativeType;
118  typedef typename Superclass::DerivativeValueType DerivativeValueType;
119  typedef typename Superclass::VirtualPointType VirtualPointType;
120  typedef typename Superclass::FixedImagePointType FixedImagePointType;
121  typedef typename Superclass::FixedImagePixelType FixedImagePixelType;
122  typedef typename Superclass::FixedTransformType FixedTransformType;
123  typedef typename Superclass::FixedImageGradientType FixedImageGradientType;
124  typedef typename FixedTransformType::JacobianType FixedImageJacobianType;
125 
126  typedef typename Superclass::MovingImagePointType MovingImagePointType;
127  typedef typename Superclass::MovingImagePixelType MovingImagePixelType;
128  typedef typename Superclass::MovingImageGradientType MovingImageGradientType;
129  typedef typename Superclass::MovingTransformType MovingTransformType;
130  typedef typename MovingTransformType::JacobianType MovingImageJacobianType;
131  typedef typename Superclass::JacobianType JacobianType;
132 
133  typedef typename Superclass::VirtualImageGradientType VirtualImageGradientType;
134 
135  typedef typename Superclass::FixedImageType FixedImageType;
136  typedef typename Superclass::MovingImageType MovingImageType;
138  typedef typename Superclass::FixedOutputPointType FixedOutputPointType;
139  typedef typename Superclass::MovingOutputPointType MovingOutputPointType;
140 
141  typedef typename Superclass::FixedTransformType::JacobianType
143  typedef typename Superclass::MovingTransformType::JacobianType
145 
146  typedef typename Superclass::NumberOfParametersType NumberOfParametersType;
148 
149  typedef typename VirtualImageType::RegionType ImageRegionType;
152 
153  /* Image dimension accessors */
154  itkStaticConstMacro(FixedImageDimension, ImageDimensionType,
155  FixedImageType::ImageDimension);
156 
157  itkStaticConstMacro(MovingImageDimension, ImageDimensionType,
158  MovingImageType::ImageDimension);
159 
160  itkStaticConstMacro(VirtualImageDimension, ImageDimensionType,
161  VirtualImageType::ImageDimension);
162 
163  // Set the radius of the neighborhood window centered at each pixel.
164  // See the note above about using a radius less than 2.
165  itkSetMacro(Radius, RadiusType);
166 
167  // Get the Radius of the neighborhood window centered at each pixel
168  itkGetMacro(Radius, RadiusType);
169  itkGetConstMacro(Radius, RadiusType);
170 
171  void Initialize(void) ITK_OVERRIDE;
172 
173 protected:
175  virtual ~ANTSNeighborhoodCorrelationImageToImageMetricv4() ITK_OVERRIDE;
176 
180 
182  typedef ANTSNeighborhoodCorrelationImageToImageMetricv4GetValueAndDerivativeThreader< ThreadedIndexedContainerPartitioner, Superclass, Self >
184 
185  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
186 
187 private:
188  ITK_DISALLOW_COPY_AND_ASSIGN(ANTSNeighborhoodCorrelationImageToImageMetricv4);
189 
190  // Radius of the neighborhood window centered at each pixel
191  RadiusType m_Radius;
192 };
193 
194 } // end namespace itk
195 
196 #ifndef ITK_MANUAL_INSTANTIATION
197 #include "itkANTSNeighborhoodCorrelationImageToImageMetricv4.hxx"
198 #endif
199 
200 #endif
Light weight base class for most itk classes.
Computes normalized cross correlation using a small neighborhood for each voxel between two images...
Class for partitioning of an ImageRegion.
Superclass::DimensionType ImageDimensionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
ImageToImageMetricv4< TFixedImage, TMovingImage, TVirtualImage, TInternalComputationValueType, TMetricTraits > Superclass
Superclass::VirtualImageType VirtualImageType