ITK  6.0.0
Insight Toolkit
itkGPUBinaryThresholdImageFilter.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 itkGPUBinaryThresholdImageFilter_h
19#define itkGPUBinaryThresholdImageFilter_h
20
21#include "itkOpenCLUtil.h"
22#include "itkGPUFunctorBase.h"
23#include "itkGPUKernelManager.h"
26
27namespace itk
28{
29
30namespace Functor
31{
32template <typename TInput, typename TOutput>
33class ITK_TEMPLATE_EXPORT GPUBinaryThreshold : public GPUFunctorBase
34{
35public:
37 {
38 m_LowerThreshold = NumericTraits<TInput>::NonpositiveMin();
39 m_UpperThreshold = NumericTraits<TInput>::max();
40 m_OutsideValue = TOutput{};
41 m_InsideValue = NumericTraits<TOutput>::max();
42 }
43
44 ~GPUBinaryThreshold() override = default;
45
46 void
47 SetLowerThreshold(const TInput & thresh)
48 {
49 m_LowerThreshold = thresh;
50 }
51 void
52 SetUpperThreshold(const TInput & thresh)
53 {
54 m_UpperThreshold = thresh;
55 }
56 void
57 SetInsideValue(const TOutput & value)
58 {
59 m_InsideValue = value;
60 }
61 void
62 SetOutsideValue(const TOutput & value)
63 {
64 m_OutsideValue = value;
65 }
66
69 int
70 SetGPUKernelArguments(GPUKernelManager::Pointer KernelManager, int KernelHandle) override
71 {
72 KernelManager->SetKernelArg(KernelHandle, 0, sizeof(TInput), &(m_LowerThreshold));
73 KernelManager->SetKernelArg(KernelHandle, 1, sizeof(TInput), &(m_UpperThreshold));
74 KernelManager->SetKernelArg(KernelHandle, 2, sizeof(TOutput), &(m_InsideValue));
75 KernelManager->SetKernelArg(KernelHandle, 3, sizeof(TOutput), &(m_OutsideValue));
76 return 4;
77 }
80private:
81 TInput m_LowerThreshold{};
82 TInput m_UpperThreshold{};
83 TOutput m_InsideValue{};
84 TOutput m_OutsideValue{};
85};
86} // end of namespace Functor
87
89itkGPUKernelClassMacro(GPUBinaryThresholdImageFilterKernel);
90
98template <typename TInputImage, typename TOutputImage>
99class ITK_TEMPLATE_EXPORT GPUBinaryThresholdImageFilter
101 TInputImage,
102 TOutputImage,
103 Functor::GPUBinaryThreshold<typename TInputImage::PixelType, typename TOutputImage::PixelType>,
104 BinaryThresholdImageFilter<TInputImage, TOutputImage>>
105{
106public:
107 ITK_DISALLOW_COPY_AND_MOVE(GPUBinaryThresholdImageFilter);
108
112 TInputImage,
113 TOutputImage,
119
121 itkNewMacro(Self);
122
124 itkOverrideGetNameOfClassMacro(GPUBinaryThresholdImageFilter);
125
127 using InputPixelType = typename TInputImage::PixelType;
128 using OutputPixelType = typename TOutputImage::PixelType;
129
132
134 itkGetOpenCLSourceFromKernelMacro(GPUBinaryThresholdImageFilterKernel);
135
136protected:
138 ~GPUBinaryThresholdImageFilter() override = default;
139
142 // virtual void BeforeThreadedGenerateData();
143
146 void
147 GPUGenerateData() override;
148};
149
157{
158public:
159 ITK_DISALLOW_COPY_AND_MOVE(GPUBinaryThresholdImageFilterFactory);
160
165
167 const char *
168 GetITKSourceVersion() const override
169 {
170 return ITK_SOURCE_VERSION;
171 }
172 const char *
173 GetDescription() const override
174 {
175 return "A Factory for GPUBinaryThresholdImageFilter";
176 }
180 itkFactorylessNewMacro(Self);
181
183 itkOverrideGetNameOfClassMacro(GPUBinaryThresholdImageFilterFactory);
184
186 static void
188 {
190
192 }
193
194private:
195#define OverrideThresholdFilterTypeMacro(ipt, opt, dm) \
196 { \
197 using InputImageType = itk::Image<ipt, dm>; \
198 using OutputImageType = itk::Image<opt, dm>; \
199 this->RegisterOverride( \
200 typeid(itk::BinaryThresholdImageFilter<InputImageType, OutputImageType>).name(), \
201 typeid(itk::GPUBinaryThresholdImageFilter<InputImageType, OutputImageType>).name(), \
202 "GPU Binary Threshold Image Filter Override", \
203 true, \
204 itk::CreateObjectFunction<GPUBinaryThresholdImageFilter<InputImageType, OutputImageType>>::New()); \
205 } \
206 ITK_MACROEND_NOOP_STATEMENT
207
209 {
210 if (IsGPUAvailable())
211 {
212 OverrideThresholdFilterTypeMacro(unsigned char, unsigned char, 1);
214 OverrideThresholdFilterTypeMacro(float, float, 1);
216 OverrideThresholdFilterTypeMacro(unsigned int, unsigned int, 1);
217 OverrideThresholdFilterTypeMacro(double, double, 1);
218
219 OverrideThresholdFilterTypeMacro(unsigned char, unsigned char, 2);
221 OverrideThresholdFilterTypeMacro(float, float, 2);
223 OverrideThresholdFilterTypeMacro(unsigned int, unsigned int, 2);
224 OverrideThresholdFilterTypeMacro(double, double, 2);
225
226 OverrideThresholdFilterTypeMacro(unsigned char, unsigned char, 3);
227 OverrideThresholdFilterTypeMacro(unsigned short, unsigned short, 3);
229 OverrideThresholdFilterTypeMacro(float, float, 3);
231 OverrideThresholdFilterTypeMacro(unsigned int, unsigned int, 3);
232 OverrideThresholdFilterTypeMacro(double, double, 3);
233 }
234 }
235};
236
237} // end of namespace itk
238
239#ifndef ITK_MANUAL_INSTANTIATION
240# include "itkGPUBinaryThresholdImageFilter.hxx"
241#endif
242
243#endif
Binarize an input image by thresholding.
typename TInputImage::PixelType InputPixelType
typename TOutputImage::PixelType OutputPixelType
~GPUBinaryThreshold() override=default
int SetGPUKernelArguments(GPUKernelManager::Pointer KernelManager, int KernelHandle) override
Base functor class for GPU functor image filters.
GPU version of binary threshold image filter.
~GPUBinaryThresholdImageFilter() override=default
itkGetOpenCLSourceFromKernelMacro(GPUBinaryThresholdImageFilterKernel)
Implements pixel-wise generic operation on one image using the GPU.
Base class for all process objects that output image data.
Light weight base class for most itk classes.
static constexpr T NonpositiveMin()
static constexpr T max(const T &)
Create instances of classes using an object factory.
static bool RegisterFactory(ObjectFactoryBase *, InsertionPositionEnum where=InsertionPositionEnum::INSERT_AT_BACK, vcl_size_t position=0)
Decorates any "simple" data type (data types without smart pointers) with a DataObject API.
#define OverrideThresholdFilterTypeMacro(ipt, opt, dm)
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:39
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
itkGPUKernelClassMacro(GPUImageOpsKernel)
bool IsGPUAvailable()