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 & operator[](const IndexType & index) const;
99
100 TPixel & operator[](const IndexType & index);
101
103 void
105
106 //
107 // Get CPU buffer pointer
108 //
109 TPixel *
111
112 const TPixel *
113 GetBufferPointer() const override;
114
118 {
119 m_DataManager->SetGPUBufferDirty();
120 return Superclass::GetPixelAccessor();
121 }
125 const AccessorType
127 {
128 m_DataManager->UpdateCPUBuffer();
129 return Superclass::GetPixelAccessor();
130 }
134 NeighborhoodAccessorFunctorType
136 {
137 m_DataManager->SetGPUBufferDirty();
138 // return Superclass::GetNeighborhoodAccessor();
140 }
144 const NeighborhoodAccessorFunctorType
146 {
147 m_DataManager->UpdateCPUBuffer();
148 // return Superclass::GetNeighborhoodAccessor();
150 }
153 void
155
159 {
160 m_DataManager->SetGPUBufferDirty();
161 return Superclass::GetPixelContainer();
162 }
165 const PixelContainer *
167 {
168 m_DataManager->UpdateCPUBuffer();
169 return Superclass::GetPixelContainer();
170 }
171
172 void
174 {
175 m_DataManager->SetCurrentCommandQueue(queueid);
176 }
177
178 int
180 {
181 return m_DataManager->GetCurrentCommandQueueID();
182 }
183
184 itkGetModifiableObjectMacro(DataManager, GPUImageDataManager<GPUImage>);
187
188 /* Override DataHasBeenGenerated() in DataObject class.
189 * We need this because CPU time stamp is always bigger
190 * than GPU's. That is because Modified() is called at
191 * the end of each filter in the pipeline so although we
192 * increment GPU's time stamp in GPUGenerateData() the
193 * CPU's time stamp will be increased after that.
194 */
195 void
197 {
198 Superclass::DataHasBeenGenerated();
199 if (m_DataManager->IsCPUBufferDirty())
200 {
201 m_DataManager->Modified();
202 }
203 }
204
206 virtual void
207 Graft(const Self * data);
208
209protected:
210 void
211 Graft(const DataObject * data) override;
213 ~GPUImage() override;
214 using Superclass::Graft;
215
216private:
218};
219
220class ITK_TEMPLATE_EXPORT GPUImageFactory : public itk::ObjectFactoryBase
221{
222public:
223 ITK_DISALLOW_COPY_AND_MOVE(GPUImageFactory);
224
229
231 const char *
232 GetITKSourceVersion() const override
233 {
234 return ITK_SOURCE_VERSION;
235 }
236 const char *
237 GetDescription() const override
238 {
239 return "A Factory for GPUImage";
240 }
244 itkFactorylessNewMacro(Self);
245
247 itkOverrideGetNameOfClassMacro(GPUImageFactory);
248
250 static void
252 {
253 auto factory = GPUImageFactory::New();
254
256 }
257
258private:
259#define OverrideImageTypeMacro(pt, dm) \
260 this->RegisterOverride(typeid(itk::Image<pt, dm>).name(), \
261 typeid(itk::GPUImage<pt, dm>).name(), \
262 "GPU Image Override", \
263 true, \
264 itk::CreateObjectFunction<GPUImage<pt, dm>>::New())
265
267 {
268 if (IsGPUAvailable())
269 {
270 // 1/2/3D
271 OverrideImageTypeMacro(unsigned char, 1);
272 OverrideImageTypeMacro(signed char, 1);
274 OverrideImageTypeMacro(unsigned int, 1);
275 OverrideImageTypeMacro(float, 1);
276 OverrideImageTypeMacro(double, 1);
277
278 OverrideImageTypeMacro(unsigned char, 2);
279 OverrideImageTypeMacro(signed char, 2);
281 OverrideImageTypeMacro(unsigned int, 2);
282 OverrideImageTypeMacro(float, 2);
283 OverrideImageTypeMacro(double, 2);
284
285 OverrideImageTypeMacro(unsigned char, 3);
286 OverrideImageTypeMacro(signed char, 3);
288 OverrideImageTypeMacro(unsigned int, 3);
289 OverrideImageTypeMacro(float, 3);
290 OverrideImageTypeMacro(double, 3);
291 }
292 }
293};
294
295template <typename T>
296class ITK_TEMPLATE_EXPORT GPUTraits
297{
298public:
299 using Type = T;
300};
301
302template <typename TPixelType, unsigned int VDimension>
303class ITK_TEMPLATE_EXPORT GPUTraits<Image<TPixelType, VDimension>>
304{
305public:
307};
308
309} // end namespace itk
310
311#ifndef ITK_MANUAL_INSTANTIATION
312# include "itkGPUImage.hxx"
313#endif
314
315#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:237
static Pointer New()
const char * GetITKSourceVersion() const override
Definition: itkGPUImage.h:232
static void RegisterOneFactory()
Definition: itkGPUImage.h:251
Templated n-dimensional image class for the GPU.
Definition: itkGPUImage.h:41
void UpdateBuffers()
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition: itkGPUImage.h:145
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:179
TPixel * GetBufferPointer() override
const TPixel * GetBufferPointer() const override
const AccessorType GetPixelAccessor() const
Definition: itkGPUImage.h:126
const PixelContainer * GetPixelContainer() const
Definition: itkGPUImage.h:166
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition: itkGPUImage.h:135
PixelContainer * GetPixelContainer()
Definition: itkGPUImage.h:158
void SetPixel(const IndexType &index, const TPixel &value)
void SetCurrentCommandQueue(int queueid)
Definition: itkGPUImage.h:173
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:196
void Initialize() override
AccessorType GetPixelAccessor()
Definition: itkGPUImage.h:117
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:259
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:39
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
bool IsGPUAvailable()