ITK  6.0.0
Insight Toolkit
itkGPUImage.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
19#ifndef itkGPUImage_h
20#define itkGPUImage_h
21
22#include "itkImage.h"
24#include "itkVersion.h"
26
27namespace itk
28{
39template <typename TPixel, unsigned int VImageDimension = 2>
40class ITK_TEMPLATE_EXPORT GPUImage : public Image<TPixel, VImageDimension>
41{
42public:
43 ITK_DISALLOW_COPY_AND_MOVE(GPUImage);
44
45 using Self = GPUImage;
50
51 itkNewMacro(Self);
52
53 itkOverrideGetNameOfClassMacro(GPUImage);
54
55 static constexpr unsigned int ImageDimension = VImageDimension;
56
57 using typename Superclass::PixelType;
58 using typename Superclass::ValueType;
59 using typename Superclass::InternalPixelType;
60 using typename Superclass::IOPixelType;
61 using typename Superclass::DirectionType;
62 using typename Superclass::SpacingType;
63 using typename Superclass::PixelContainer;
64 using typename Superclass::SizeType;
65 using typename Superclass::IndexType;
66 using typename Superclass::OffsetType;
67 using typename Superclass::RegionType;
70 using typename Superclass::AccessorType;
71
73
75 // NeighborhoodAccessorFunctorType;
76
77 //
78 // Allocate CPU and GPU memory space
79 //
80 void
81 Allocate(bool initialize = false) override;
82
83 void
84 Initialize() override;
85
86 void
87 FillBuffer(const TPixel & value);
88
89 void
90 SetPixel(const IndexType & index, const TPixel & value);
91
92 const TPixel &
93 GetPixel(const IndexType & index) const;
94
95 TPixel &
96 GetPixel(const IndexType & index);
97
98 const TPixel &
99 operator[](const IndexType & index) const;
100
101 TPixel &
102 operator[](const IndexType & index);
103
105 void
107
108 //
109 // Get CPU buffer pointer
110 //
111 TPixel *
113
114 const TPixel *
115 GetBufferPointer() const override;
116
120 {
121 m_DataManager->SetGPUBufferDirty();
122 return Superclass::GetPixelAccessor();
123 }
127 const AccessorType
129 {
130 m_DataManager->UpdateCPUBuffer();
131 return Superclass::GetPixelAccessor();
132 }
136 NeighborhoodAccessorFunctorType
138 {
139 m_DataManager->SetGPUBufferDirty();
140 // return Superclass::GetNeighborhoodAccessor();
142 }
146 const NeighborhoodAccessorFunctorType
148 {
149 m_DataManager->UpdateCPUBuffer();
150 // return Superclass::GetNeighborhoodAccessor();
152 }
155 void
157
161 {
162 m_DataManager->SetGPUBufferDirty();
163 return Superclass::GetPixelContainer();
164 }
167 const PixelContainer *
169 {
170 m_DataManager->UpdateCPUBuffer();
171 return Superclass::GetPixelContainer();
172 }
173
174 void
176 {
177 m_DataManager->SetCurrentCommandQueue(queueid);
178 }
179
180 int
182 {
183 return m_DataManager->GetCurrentCommandQueueID();
184 }
185
186 itkGetModifiableObjectMacro(DataManager, GPUImageDataManager<GPUImage>);
189
190 /* Override DataHasBeenGenerated() in DataObject class.
191 * We need this because CPU time stamp is always bigger
192 * than GPU's. That is because Modified() is called at
193 * the end of each filter in the pipeline so although we
194 * increment GPU's time stamp in GPUGenerateData() the
195 * CPU's time stamp will be increased after that.
196 */
197 void
199 {
200 Superclass::DataHasBeenGenerated();
201 if (m_DataManager->IsCPUBufferDirty())
202 {
203 m_DataManager->Modified();
204 }
205 }
206
208 virtual void
209 Graft(const Self * data);
210
211protected:
212 void
213 Graft(const DataObject * data) override;
215 ~GPUImage() override;
216 using Superclass::Graft;
217
218private:
220};
221
222class ITK_TEMPLATE_EXPORT GPUImageFactory : public itk::ObjectFactoryBase
223{
224public:
225 ITK_DISALLOW_COPY_AND_MOVE(GPUImageFactory);
226
231
233 const char *
234 GetITKSourceVersion() const override
235 {
236 return ITK_SOURCE_VERSION;
237 }
238 const char *
239 GetDescription() const override
240 {
241 return "A Factory for GPUImage";
242 }
246 itkFactorylessNewMacro(Self);
247
249 itkOverrideGetNameOfClassMacro(GPUImageFactory);
250
252 static void
254 {
255 auto factory = GPUImageFactory::New();
256
258 }
259
260private:
261#define OverrideImageTypeMacro(pt, dm) \
262 this->RegisterOverride(typeid(itk::Image<pt, dm>).name(), \
263 typeid(itk::GPUImage<pt, dm>).name(), \
264 "GPU Image Override", \
265 true, \
266 itk::CreateObjectFunction<GPUImage<pt, dm>>::New())
267
269 {
270 if (IsGPUAvailable())
271 {
272 // 1/2/3D
273 OverrideImageTypeMacro(unsigned char, 1);
274 OverrideImageTypeMacro(signed char, 1);
276 OverrideImageTypeMacro(unsigned int, 1);
277 OverrideImageTypeMacro(float, 1);
278 OverrideImageTypeMacro(double, 1);
279
280 OverrideImageTypeMacro(unsigned char, 2);
281 OverrideImageTypeMacro(signed char, 2);
283 OverrideImageTypeMacro(unsigned int, 2);
284 OverrideImageTypeMacro(float, 2);
285 OverrideImageTypeMacro(double, 2);
286
287 OverrideImageTypeMacro(unsigned char, 3);
288 OverrideImageTypeMacro(signed char, 3);
290 OverrideImageTypeMacro(unsigned int, 3);
291 OverrideImageTypeMacro(float, 3);
292 OverrideImageTypeMacro(double, 3);
293 }
294 }
295};
296
297template <typename T>
298class ITK_TEMPLATE_EXPORT GPUTraits
299{
300public:
301 using Type = T;
302};
303
304template <typename TPixelType, unsigned int VDimension>
305class ITK_TEMPLATE_EXPORT GPUTraits<Image<TPixelType, VDimension>>
306{
307public:
309};
310
311} // end namespace itk
312
313#ifndef ITK_MANUAL_INSTANTIATION
314# include "itkGPUImage.hxx"
315#endif
316
317#endif
Base class for all data objects in ITK.
Provides a common API for pixel accessors for Image and VectorImage.
Give access to partial aspects a type.
GPU memory manager implemented using OpenCL. Required by GPUImage class.
const char * GetDescription() const override
Definition: itkGPUImage.h:239
static Pointer New()
const char * GetITKSourceVersion() const override
Definition: itkGPUImage.h:234
static void RegisterOneFactory()
Definition: itkGPUImage.h:253
Templated n-dimensional image class for the GPU.
Definition: itkGPUImage.h:41
void UpdateBuffers()
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition: itkGPUImage.h:147
void Graft(const DataObject *data) override
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition: itkGPUImage.h:69
const TPixel & GetPixel(const IndexType &index) const
int GetCurrentCommandQueueID()
Definition: itkGPUImage.h:181
TPixel * GetBufferPointer() override
const TPixel * GetBufferPointer() const override
const AccessorType GetPixelAccessor() const
Definition: itkGPUImage.h:128
const PixelContainer * GetPixelContainer() const
Definition: itkGPUImage.h:168
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition: itkGPUImage.h:137
PixelContainer * GetPixelContainer()
Definition: itkGPUImage.h:160
void SetPixel(const IndexType &index, const TPixel &value)
void SetCurrentCommandQueue(int queueid)
Definition: itkGPUImage.h:175
virtual void Graft(const Self *data)
void Allocate(bool initialize=false) override
GPUDataManager * GetGPUDataManager()
TPixel & operator[](const IndexType &index)
~GPUImage() override
void DataHasBeenGenerated() override
Definition: itkGPUImage.h:198
void Initialize() override
AccessorType GetPixelAccessor()
Definition: itkGPUImage.h:119
const TPixel & operator[](const IndexType &index) const
void SetPixelContainer(PixelContainer *container)
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkGPUImage.h:68
TPixel & GetPixel(const IndexType &index)
void FillBuffer(const TPixel &value)
Templated n-dimensional image class.
Definition: itkImage.h:89
Defines an itk::Image front-end to a standard C-array.
Light weight base class for most itk classes.
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
Create instances of classes using an object factory.
static bool RegisterFactory(ObjectFactoryBase *, InsertionPositionEnum where=InsertionPositionEnum::INSERT_AT_BACK, vcl_size_t position=0)
Base class for most ITK classes.
Definition: itkObject.h:62
Implements a weak reference to an object.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
#define OverrideImageTypeMacro(pt, dm)
Definition: itkGPUImage.h:261
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:39
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
bool IsGPUAvailable()