ITK  6.0.0
Insight Toolkit
itkBlockMatchingImageFilter.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 itkBlockMatchingImageFilter_h
19#define itkBlockMatchingImageFilter_h
20
21#include "itkMeshToMeshFilter.h"
22#include "itkImage.h"
23#include "itkPointSet.h"
24#include "itkVector.h"
26
27#include <memory> // For unique_ptr.
28
29namespace itk
30{
72template <
73 typename TFixedImage,
74 typename TMovingImage = TFixedImage,
75 typename TFeatures = PointSet<Matrix<SpacePrecisionType, TFixedImage::ImageDimension, TFixedImage::ImageDimension>,
76 TFixedImage::ImageDimension>,
77 class TDisplacements =
78 PointSet<Vector<typename TFeatures::PointType::ValueType, TFeatures::PointDimension>, TFeatures::PointDimension>,
79 class TSimilarities = PointSet<SpacePrecisionType, TDisplacements::PointDimension>>
80class ITK_TEMPLATE_EXPORT BlockMatchingImageFilter : public MeshToMeshFilter<TFeatures, TDisplacements>
81{
82public:
83 ITK_DISALLOW_COPY_AND_MOVE(BlockMatchingImageFilter);
84
85 static constexpr unsigned int ImageDimension = TFixedImage::ImageDimension;
86
91
93 using FixedImageType = TFixedImage;
95 using FixedImagePixelType = typename FixedImageType::PixelType;
96
98 using MovingImageType = TMovingImage;
100
102 using FeaturePointsType = TFeatures;
106
108 using DisplacementsType = TDisplacements;
111 using DisplacementsVector = typename DisplacementsType::PixelType;
112
114 using SimilaritiesType = TSimilarities;
117 using SimilaritiesValue = typename SimilaritiesType::PixelType;
118
124
126 itkNewMacro(Self);
127
129 itkOverrideGetNameOfClassMacro(BlockMatchingImageFilter);
130
132 itkSetMacro(BlockRadius, ImageSizeType);
133 itkGetConstMacro(BlockRadius, ImageSizeType);
137 itkSetMacro(SearchRadius, ImageSizeType);
138 itkGetConstMacro(SearchRadius, ImageSizeType);
156 inline DisplacementsType *
158 {
159 return dynamic_cast<DisplacementsType *>(this->ProcessObject::GetOutput(0));
160 }
161
162 inline SimilaritiesType *
164 {
165 return dynamic_cast<SimilaritiesType *>(this->ProcessObject::GetOutput(1));
166 }
167
169 using Superclass::MakeOutput;
172
173protected:
178 void
180
182 void
184
186 virtual void
188
189 virtual void
191
193 virtual void
195
197 void
198 GenerateData() override;
199
201 ~BlockMatchingImageFilter() override = default;
202
203 void
204 PrintSelf(std::ostream & os, Indent indent) const override;
205
212 ThreaderCallback(void * arg);
213
217 {
219 };
220
221private:
222 // algorithm parameters
223 ImageSizeType m_BlockRadius{};
224 ImageSizeType m_SearchRadius{};
225
226 // temporary dynamic arrays for storing threads outputs
227 SizeValueType m_PointsCount{};
228 std::unique_ptr<DisplacementsVector[]> m_DisplacementsVectorsArray;
229 std::unique_ptr<SimilaritiesValue[]> m_SimilaritiesValuesArray;
230};
231} // end namespace itk
232
233#ifndef ITK_MANUAL_INSTANTIATION
234# include "itkBlockMatchingImageFilter.hxx"
235#endif
236
237#endif
Computes displacements of given points from a fixed image in a floating image.
void PrintSelf(std::ostream &os, Indent indent) const override
itkSetInputMacro(MovingImage, MovingImageType)
itkGetInputMacro(FeaturePoints, FeaturePointsType)
typename FixedImageType::ConstPointer FixedImageConstPointer
typename DisplacementsType::ConstPointer DisplacementsConstPointer
itkGetInputMacro(FixedImage, FixedImageType)
typename FixedImageType::PixelType FixedImagePixelType
virtual void BeforeThreadedGenerateData()
typename SimilaritiesType::PixelType SimilaritiesValue
typename MovingImageType::ConstPointer MovingImageConstPointer
typename SimilaritiesType::Pointer SimilaritiesPointer
void EnlargeOutputRequestedRegion(DataObject *output) override
std::unique_ptr< SimilaritiesValue[]> m_SimilaritiesValuesArray
typename SimilaritiesType::ConstPointer SimilaritiesConstPointer
virtual void ThreadedGenerateData(ThreadIdType threadId)
typename DisplacementsType::PixelType DisplacementsVector
void GenerateOutputInformation() override
itkSetInputMacro(FixedImage, FixedImageType)
typename FeaturePointsType::ConstPointer FeaturePointsConstPointer
typename FeaturePointsType::PointType FeaturePointsPhysicalCoordinates
typename DisplacementsType::Pointer DisplacementsPointer
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreaderCallback(void *arg)
DataObject::Pointer MakeOutput(ProcessObject::DataObjectPointerArraySizeType idx) override
Make a DataObject of the correct type to used as the specified output.
itkSetInputMacro(FeaturePoints, FeaturePointsType)
typename FeaturePointsType::Pointer FeaturePointsPointer
itkGetInputMacro(MovingImage, MovingImageType)
virtual void AfterThreadedGenerateData()
std::unique_ptr< DisplacementsVector[]> m_DisplacementsVectorsArray
~BlockMatchingImageFilter() override=default
Base class for all data objects in ITK.
An image region represents a structured region of data.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
MeshToMeshFilter is the base class for all process objects that output mesh data, and require mesh da...
DataObject * GetOutput(const DataObjectIdentifierType &key)
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:69