ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
69
70template <
71 typename TFixedImage,
72 typename TMovingImage = TFixedImage,
74 TFixedImage::ImageDimension>,
75 class TDisplacements =
78class ITK_TEMPLATE_EXPORT BlockMatchingImageFilter : public MeshToMeshFilter<TFeatures, TDisplacements>
79{
80public:
81 ITK_DISALLOW_COPY_AND_MOVE(BlockMatchingImageFilter);
82
83 static constexpr unsigned int ImageDimension = TFixedImage::ImageDimension;
84
89
91 using FixedImageType = TFixedImage;
92 using FixedImageConstPointer = typename FixedImageType::ConstPointer;
93 using FixedImagePixelType = typename FixedImageType::PixelType;
94
96 using MovingImageType = TMovingImage;
97 using MovingImageConstPointer = typename MovingImageType::ConstPointer;
98
100 using FeaturePointsType = TFeatures;
101 using FeaturePointsPointer = typename FeaturePointsType::Pointer;
102 using FeaturePointsConstPointer = typename FeaturePointsType::ConstPointer;
103 using FeaturePointsPhysicalCoordinates = typename FeaturePointsType::PointType;
104
106 using DisplacementsType = TDisplacements;
107 using DisplacementsPointer = typename DisplacementsType::Pointer;
108 using DisplacementsConstPointer = typename DisplacementsType::ConstPointer;
109 using DisplacementsVector = typename DisplacementsType::PixelType;
110
112 using SimilaritiesType = TSimilarities;
113 using SimilaritiesPointer = typename SimilaritiesType::Pointer;
114 using SimilaritiesConstPointer = typename SimilaritiesType::ConstPointer;
115 using SimilaritiesValue = typename SimilaritiesType::PixelType;
116
122
124 itkNewMacro(Self);
125
127 itkOverrideGetNameOfClassMacro(BlockMatchingImageFilter);
128
130 itkSetMacro(BlockRadius, ImageSizeType);
131 itkGetConstMacro(BlockRadius, ImageSizeType);
133
135 itkSetMacro(SearchRadius, ImageSizeType);
136 itkGetConstMacro(SearchRadius, ImageSizeType);
138
143
148
153
154 inline DisplacementsType *
156 {
157 return dynamic_cast<DisplacementsType *>(this->ProcessObject::GetOutput(0));
158 }
159
160 inline SimilaritiesType *
162 {
163 return dynamic_cast<SimilaritiesType *>(this->ProcessObject::GetOutput(1));
164 }
165
167 using Superclass::MakeOutput;
170
171protected:
176 void
178
180 void
182
184 virtual void
186
187 virtual void
189
191 virtual void
193
195 void
196 GenerateData() override;
197
199 ~BlockMatchingImageFilter() override = default;
200
201 void
202 PrintSelf(std::ostream & os, Indent indent) const override;
203
210 ThreaderCallback(void * arg);
211
215 {
217 };
218
219private:
220 // algorithm parameters
223
224 // temporary dynamic arrays for storing threads outputs
226 std::unique_ptr<DisplacementsVector[]> m_DisplacementsVectorsArray;
227 std::unique_ptr<SimilaritiesValue[]> m_SimilaritiesValuesArray;
228};
229} // end namespace itk
230
231#ifndef ITK_MANUAL_INSTANTIATION
232# include "itkBlockMatchingImageFilter.hxx"
233#endif
234
235#endif
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
static constexpr unsigned int ImageDimension
MeshToMeshFilter< FeaturePointsType, DisplacementsType > Superclass
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)
ImageRegion< ImageDimension > ImageRegionType
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.
SmartPointer< Self > Pointer
An image region represents a structured region of data.
Control indentation during Print() invocation.
Definition itkIndent.h:50
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition itkPointSet.h:82
DataObject * GetOutput(const DataObjectIdentifierType &key)
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
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
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70