ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkGPUBinaryThresholdImageFilter.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 itkGPUBinaryThresholdImageFilter_h
19 #define itkGPUBinaryThresholdImageFilter_h
20 
21 #include "itkOpenCLUtil.h"
22 #include "itkGPUFunctorBase.h"
23 #include "itkGPUKernelManager.h"
26 
27 namespace itk
28 {
29 
30 namespace Functor
31 {
32 template< typename TInput, typename TOutput >
34 {
35 public:
37  {
42  }
43 
45  }
46 
47  void SetLowerThreshold(const TInput & thresh)
48  {
49  m_LowerThreshold = thresh;
50  }
51  void SetUpperThreshold(const TInput & thresh)
52  {
53  m_UpperThreshold = thresh;
54  }
55  void SetInsideValue(const TOutput & value)
56  {
57  m_InsideValue = value;
58  }
59  void SetOutsideValue(const TOutput & value)
60  {
61  m_OutsideValue = value;
62  }
63 
66  int SetGPUKernelArguments(GPUKernelManager::Pointer KernelManager, int KernelHandle)
67  {
68  KernelManager->SetKernelArg(KernelHandle, 0, sizeof(TInput), &(m_LowerThreshold) );
69  KernelManager->SetKernelArg(KernelHandle, 1, sizeof(TInput), &(m_UpperThreshold) );
70  KernelManager->SetKernelArg(KernelHandle, 2, sizeof(TOutput), &(m_InsideValue) );
71  KernelManager->SetKernelArg(KernelHandle, 3, sizeof(TOutput), &(m_OutsideValue) );
72  return 4;
73  }
75 
76 private:
79  TOutput m_InsideValue;
80  TOutput m_OutsideValue;
81 };
82 } // end of namespace Functor
83 
85 itkGPUKernelClassMacro(GPUBinaryThresholdImageFilterKernel);
86 
94 template< typename TInputImage, typename TOutputImage >
96  public
97  GPUUnaryFunctorImageFilter< TInputImage, TOutputImage,
98  Functor::GPUBinaryThreshold<
99  typename TInputImage::PixelType,
100  typename TOutputImage::PixelType >,
101  BinaryThresholdImageFilter<TInputImage, TOutputImage> >
102 {
103 public:
104 
107  typedef GPUUnaryFunctorImageFilter< TInputImage, TOutputImage,
109  typename TInputImage::PixelType,
110  typename TOutputImage::PixelType >,
115 
117  itkNewMacro(Self);
118 
121 
123  typedef typename TInputImage::PixelType InputPixelType;
124  typedef typename TOutputImage::PixelType OutputPixelType;
125 
128 
130  itkGetOpenCLSourceFromKernelMacro(GPUBinaryThresholdImageFilterKernel);
131 
132 protected:
135  }
136 
139  //virtual void BeforeThreadedGenerateData();
140 
143  virtual void GPUGenerateData() ITK_OVERRIDE;
144 
145 private:
146  GPUBinaryThresholdImageFilter(const Self &) ITK_DELETE_FUNCTION;
147  void operator=(const Self &) ITK_DELETE_FUNCTION;
148 
149 };
150 
158 {
159 public:
161  typedef ObjectFactoryBase Superclass;
164 
166  virtual const char* GetITKSourceVersion() const ITK_OVERRIDE
167  {
168  return ITK_SOURCE_VERSION;
169  }
170  const char* GetDescription() const ITK_OVERRIDE
171  {
172  return "A Factory for GPUBinaryThresholdImageFilter";
173  }
175 
177  itkFactorylessNewMacro(Self);
178 
181 
183  static void RegisterOneFactory(void)
184  {
186 
188  }
189 
190 private:
191  GPUBinaryThresholdImageFilterFactory(const Self&) ITK_DELETE_FUNCTION;
192  void operator=(const Self&) ITK_DELETE_FUNCTION;
193 
194 #define OverrideThresholdFilterTypeMacro(ipt,opt,dm) \
195  { \
196  typedef itk::Image<ipt,dm> InputImageType; \
197  typedef itk::Image<opt,dm> OutputImageType; \
198  this->RegisterOverride( \
199  typeid(itk::BinaryThresholdImageFilter<InputImageType,OutputImageType>).name(), \
200  typeid(itk::GPUBinaryThresholdImageFilter<InputImageType,OutputImageType>).name(), \
201  "GPU Binary Threshold Image Filter Override", \
202  true, \
203  itk::CreateObjectFunction<GPUBinaryThresholdImageFilter<InputImageType,OutputImageType> >::New() ); \
204  }
205 
207  {
208  if( IsGPUAvailable() )
209  {
210  OverrideThresholdFilterTypeMacro(unsigned char, unsigned char, 1);
211  OverrideThresholdFilterTypeMacro(char, char, 1);
212  OverrideThresholdFilterTypeMacro(float,float,1);
214  OverrideThresholdFilterTypeMacro(unsigned int,unsigned int,1);
215  OverrideThresholdFilterTypeMacro(double,double,1);
216 
217  OverrideThresholdFilterTypeMacro(unsigned char, unsigned char, 2);
218  OverrideThresholdFilterTypeMacro(char, char, 2);
219  OverrideThresholdFilterTypeMacro(float,float,2);
221  OverrideThresholdFilterTypeMacro(unsigned int,unsigned int,2);
222  OverrideThresholdFilterTypeMacro(double,double,2);
223 
224  OverrideThresholdFilterTypeMacro(unsigned char, unsigned char, 3);
225  OverrideThresholdFilterTypeMacro(unsigned short, unsigned short, 3);
226  OverrideThresholdFilterTypeMacro(char, char, 3);
227  OverrideThresholdFilterTypeMacro(float,float,3);
229  OverrideThresholdFilterTypeMacro(unsigned int,unsigned int,3);
230  OverrideThresholdFilterTypeMacro(double,double,3);
231  }
232  }
233 
234 };
235 
236 } // end of namespace itk
237 
238 #ifndef ITK_MANUAL_INSTANTIATION
239 #include "itkGPUBinaryThresholdImageFilter.hxx"
240 #endif
241 
242 #endif
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:40
SimpleDataObjectDecorator< InputPixelType > InputPixelObjectType
int SetGPUKernelArguments(GPUKernelManager::Pointer KernelManager, int KernelHandle)
Base class for all process objects that output image data.
Create instances of classes using an object factory.
itkGPUKernelClassMacro(GPUImageOpsKernel)
GPU version of binary threshold image filter.
virtual const char * GetITKSourceVersion() const override
Binarize an input image by thresholding.
GPUUnaryFunctorImageFilter< TInputImage, TOutputImage, Functor::GPUBinaryThreshold< typename TInputImage::PixelType, typename TOutputImage::PixelType >, BinaryThresholdImageFilter< TInputImage, TOutputImage > > GPUSuperclass
Base functor class for GPU functor image filters.
Decorates any &quot;simple&quot; data type (data types without smart pointers) with a DataObject API...
static T max(const T &)
#define OverrideThresholdFilterTypeMacro(ipt, opt, dm)
BinaryThresholdImageFilter< TInputImage, TOutputImage > CPUSuperclass
static bool RegisterFactory(ObjectFactoryBase *, InsertionPositionType where=INSERT_AT_BACK, vcl_size_t position=0)
Implements pixel-wise generic operation on one image using the GPU.
bool IsGPUAvailable()
itkGetOpenCLSourceFromKernelMacro(GPUBinaryThresholdImageFilterKernel)
virtual void GPUGenerateData() override