ITK  5.3.0
Insight Toolkit
itkANTSNeighborhoodCorrelationImageToImageMetricv4.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 {
26 
94 template <typename TFixedImage,
95  typename TMovingImage,
96  typename TVirtualImage = TFixedImage,
97  typename TInternalComputationValueType = double,
98  typename TMetricTraits =
99  DefaultImageToImageMetricTraitsv4<TFixedImage, TMovingImage, TVirtualImage, TInternalComputationValueType>>
101  : public ImageToImageMetricv4<TFixedImage, TMovingImage, TVirtualImage, TInternalComputationValueType, TMetricTraits>
102 {
103 public:
104  ITK_DISALLOW_COPY_AND_MOVE(ANTSNeighborhoodCorrelationImageToImageMetricv4);
105 
108  using Superclass =
112 
114  itkNewMacro(Self);
115 
118 
120  using typename Superclass::MeasureType;
121  using typename Superclass::DerivativeType;
122  using typename Superclass::DerivativeValueType;
123  using typename Superclass::VirtualPointType;
124  using typename Superclass::FixedImagePointType;
125  using typename Superclass::FixedImagePixelType;
126  using typename Superclass::FixedTransformType;
127  using typename Superclass::FixedImageGradientType;
129 
130  using typename Superclass::MovingImagePointType;
131  using typename Superclass::MovingImagePixelType;
132  using typename Superclass::MovingImageGradientType;
133  using typename Superclass::MovingTransformType;
135  using typename Superclass::JacobianType;
136 
137  using typename Superclass::VirtualImageGradientType;
138 
139  using typename Superclass::FixedImageType;
140  using typename Superclass::MovingImageType;
141  using VirtualImageType = typename Superclass::VirtualImageType;
142  using typename Superclass::FixedOutputPointType;
143  using typename Superclass::MovingOutputPointType;
144 
145  using FixedTransformJacobianType = typename Superclass::FixedTransformType::JacobianType;
146  using MovingTransformJacobianType = typename Superclass::MovingTransformType::JacobianType;
147 
148  using typename Superclass::NumberOfParametersType;
149  using typename Superclass::ImageDimensionType;
150 
154 
155  /* Image dimension accessors */
156  static constexpr ImageDimensionType FixedImageDimension = FixedImageType::ImageDimension;
157 
158  static constexpr ImageDimensionType MovingImageDimension = MovingImageType::ImageDimension;
159 
160  static constexpr ImageDimensionType VirtualImageDimension = VirtualImageType::ImageDimension;
161 
162  // Set the radius of the neighborhood window centered at each pixel.
163  // See the note above about using a radius less than 2.
164  itkSetMacro(Radius, RadiusType);
165 
166  // Get the Radius of the neighborhood window centered at each pixel
167  itkGetMacro(Radius, RadiusType);
168  itkGetConstMacro(Radius, RadiusType);
169 
170  void
171  Initialize() override;
172 
173 protected:
176 
178  ThreadedImageRegionPartitioner<VirtualImageDimension>,
179  Superclass,
184  Superclass,
186 
189  Superclass,
193  Superclass,
195 
196  void
197  PrintSelf(std::ostream & os, Indent indent) const override;
198 
199 private:
200  // Radius of the neighborhood window centered at each pixel
202 };
203 
204 } // end namespace itk
205 
206 #ifndef ITK_MANUAL_INSTANTIATION
207 # include "itkANTSNeighborhoodCorrelationImageToImageMetricv4.hxx"
208 #endif
209 
210 #endif
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::ImageRegionType
typename VirtualImageType::RegionType ImageRegionType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:151
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4
Computes normalized cross correlation using a small neighborhood for each voxel between two images,...
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:100
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::FixedTransformJacobianType
typename Superclass::FixedTransformType::JacobianType FixedTransformJacobianType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:145
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::IndexType
typename VirtualImageType::IndexType IndexType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:153
itk::ThreadedImageRegionPartitioner
Class for partitioning of an ImageRegion.
Definition: itkThreadedImageRegionPartitioner.h:45
itk::ThreadedIndexedContainerPartitioner
Partitions an indexed container.
Definition: itkThreadedIndexedContainerPartitioner.h:45
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::MovingImageJacobianType
typename MovingTransformType::JacobianType MovingImageJacobianType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:134
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::VirtualImageType
typename Superclass::VirtualImageType VirtualImageType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:141
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkImageToImageMetricv4.h
itkANTSNeighborhoodCorrelationImageToImageMetricv4GetValueAndDerivativeThreader.h
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::MovingTransformJacobianType
typename Superclass::MovingTransformType::JacobianType MovingTransformJacobianType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:146
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImageToImageMetricv4
Definition: itkImageToImageMetricv4.h:174
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::FixedImageJacobianType
typename FixedTransformType::JacobianType FixedImageJacobianType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:128
itk::Array2D
Array2D class representing a 2D array with size defined at construction time.
Definition: itkArray2D.h:45
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::RadiusType
typename VirtualImageType::SizeType RadiusType
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:152
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition: itkAddImageFilter.h:89
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4::m_Radius
RadiusType m_Radius
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4.h:201
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4GetValueAndDerivativeThreader
Threading implementation for ANTS CC metric ANTSNeighborhoodCorrelationImageToImageMetricv4 ....
Definition: itkANTSNeighborhoodCorrelationImageToImageMetricv4GetValueAndDerivativeThreader.h:58