ITK  5.2.0
Insight Toolkit
itkGPUMeanImageFilter.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  * 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 itkGPUMeanImageFilter_h
19 #define itkGPUMeanImageFilter_h
20 
21 #include "itkMeanImageFilter.h"
22 #include "itkGPUBoxImageFilter.h"
23 #include "itkVersion.h"
24 #include "itkObjectFactoryBase.h"
25 #include "itkOpenCLUtil.h"
26 
27 namespace itk
28 {
40 itkGPUKernelClassMacro(GPUMeanImageFilterKernel);
41 
42 template <typename TInputImage, typename TOutputImage>
43 class ITK_TEMPLATE_EXPORT GPUMeanImageFilter
44  : // public GPUImageToImageFilter<
45  // TInputImage, TOutputImage,
46  // MeanImageFilter< TInputImage,
47  // TOutputImage > >
48  public GPUBoxImageFilter<TInputImage, TOutputImage, MeanImageFilter<TInputImage, TOutputImage>>
49 {
50 public:
51  ITK_DISALLOW_COPY_AND_MOVE(GPUMeanImageFilter);
52 
58 
59  itkNewMacro(Self);
60 
63 
65  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
66  using OutputImagePixelType = typename Superclass::OutputImagePixelType;
67 
69  using InputImageType = TInputImage;
70  using InputImagePointer = typename InputImageType::Pointer;
71  using InputImageConstPointer = typename InputImageType::ConstPointer;
73  using InputImagePixelType = typename InputImageType::PixelType;
74 
76  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
77  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
78 
80  itkGetOpenCLSourceFromKernelMacro(GPUMeanImageFilterKernel);
81 
82 protected:
84  ~GPUMeanImageFilter() override;
85 
86  void
87  PrintSelf(std::ostream & os, Indent indent) const override;
88 
89  void
90  GPUGenerateData() override;
91 
92 private:
94 };
95 
103 {
104 public:
105  ITK_DISALLOW_COPY_AND_MOVE(GPUMeanImageFilterFactory);
106 
111 
113  const char *
114  GetITKSourceVersion() const override
115  {
116  return ITK_SOURCE_VERSION;
117  }
118  const char *
119  GetDescription() const override
120  {
121  return "A Factory for GPUMeanImageFilter";
122  }
124 
126  itkFactorylessNewMacro(Self);
127 
130 
132  static void
134  {
136 
138  }
139 
140 private:
141 #define OverrideMeanFilterTypeMacro(ipt, opt, dm) \
142  { \
143  using InputImageType = Image<ipt, dm>; \
144  using OutputImageType = Image<opt, dm>; \
145  this->RegisterOverride(typeid(MeanImageFilter<InputImageType, OutputImageType>).name(), \
146  typeid(GPUMeanImageFilter<InputImageType, OutputImageType>).name(), \
147  "GPU Mean Image Filter Override", \
148  true, \
149  CreateObjectFunction<GPUMeanImageFilter<InputImageType, OutputImageType>>::New()); \
150  } \
151  ITK_MACROEND_NOOP_STATEMENT
152 
154  {
155  if (IsGPUAvailable())
156  {
157  OverrideMeanFilterTypeMacro(unsigned char, unsigned char, 1);
158  OverrideMeanFilterTypeMacro(char, char, 1);
159  OverrideMeanFilterTypeMacro(float, float, 1);
160  OverrideMeanFilterTypeMacro(int, int, 1);
161  OverrideMeanFilterTypeMacro(unsigned int, unsigned int, 1);
162  OverrideMeanFilterTypeMacro(double, double, 1);
163 
164  OverrideMeanFilterTypeMacro(unsigned char, unsigned char, 2);
165  OverrideMeanFilterTypeMacro(char, char, 2);
166  OverrideMeanFilterTypeMacro(float, float, 2);
167  OverrideMeanFilterTypeMacro(int, int, 2);
168  OverrideMeanFilterTypeMacro(unsigned int, unsigned int, 2);
169  OverrideMeanFilterTypeMacro(double, double, 2);
170 
171  OverrideMeanFilterTypeMacro(unsigned char, unsigned char, 3);
172  OverrideMeanFilterTypeMacro(char, char, 3);
173  OverrideMeanFilterTypeMacro(float, float, 3);
174  OverrideMeanFilterTypeMacro(int, int, 3);
175  OverrideMeanFilterTypeMacro(unsigned int, unsigned int, 3);
176  OverrideMeanFilterTypeMacro(double, double, 3);
177  }
178  }
179 };
180 
181 } // end namespace itk
182 
183 #ifndef ITK_MANUAL_INSTANTIATION
184 # include "itkGPUMeanImageFilter.hxx"
185 #endif
186 
187 #endif
itk::ObjectFactoryBase
Create instances of classes using an object factory.
Definition: itkObjectFactoryBase.h:62
OverrideMeanFilterTypeMacro
#define OverrideMeanFilterTypeMacro(ipt, opt, dm)
Definition: itkGPUMeanImageFilter.h:141
itk::GPUMeanImageFilterFactory::GetITKSourceVersion
const char * GetITKSourceVersion() const override
Definition: itkGPUMeanImageFilter.h:114
itk::GPUMeanImageFilterFactory::GetDescription
const char * GetDescription() const override
Definition: itkGPUMeanImageFilter.h:119
itk::GPUMeanImageFilterFactory
Object Factory implementation for GPUMeanImageFilter.
Definition: itkGPUMeanImageFilter.h:102
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itkOpenCLUtil.h
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
ITK_SOURCE_VERSION
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:39
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itkObjectFactoryBase.h
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::GPUMeanImageFilterFactory::Self
GPUMeanImageFilterFactory Self
Definition: itkGPUMeanImageFilter.h:107
itk::ObjectFactoryBase::ObjectFactoryBase
ObjectFactoryBase()
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::GPUMeanImageFilterFactory::New
static Pointer New()
itk::GPUMeanImageFilterFactory::RegisterOneFactory
static void RegisterOneFactory()
Definition: itkGPUMeanImageFilter.h:133
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkVersion.h
itk::itkGPUKernelClassMacro
itkGPUKernelClassMacro(GPUImageOpsKernel)
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::GPUMeanImageFilterFactory::GPUMeanImageFilterFactory
GPUMeanImageFilterFactory()
Definition: itkGPUMeanImageFilter.h:153
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::GPUMeanImageFilter::m_MeanFilterGPUKernelHandle
int m_MeanFilterGPUKernelHandle
Definition: itkGPUMeanImageFilter.h:93
itk::GPUMeanImageFilter
GPU-enabled implementation of the MeanImageFilter.
Definition: itkGPUMeanImageFilter.h:43
itk::ObjectFactoryBase::RegisterFactory
static bool RegisterFactory(ObjectFactoryBase *, InsertionPositionEnum where=InsertionPositionEnum::INSERT_AT_BACK, vcl_size_t position=0)
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::GPUBoxImageFilter
A base class for all the GPU filters working on a box neighborhood.
Definition: itkGPUBoxImageFilter.h:42
itkMeanImageFilter.h
itkGPUBoxImageFilter.h
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::IsGPUAvailable
bool IsGPUAvailable()