ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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();
123 }
124
126 const AccessorType
128 {
129 m_DataManager->UpdateCPUBuffer();
131 }
132
134 NeighborhoodAccessorFunctorType
136 {
137 m_DataManager->SetGPUBufferDirty();
138 // return Superclass::GetNeighborhoodAccessor();
140 }
141
143 const NeighborhoodAccessorFunctorType
145 {
146 m_DataManager->UpdateCPUBuffer();
147 // return Superclass::GetNeighborhoodAccessor();
149 }
150
151 void
153
158 {
159 m_DataManager->SetGPUBufferDirty();
161 }
162
163 const PixelContainer *
165 {
166 m_DataManager->UpdateCPUBuffer();
168 }
169
170
171 void
173 {
174 m_DataManager->SetCurrentCommandQueue(queueid);
175 }
176
177 int
179 {
180 return m_DataManager->GetCurrentCommandQueueID();
181 }
182
183 itkGetModifiableObjectMacro(DataManager, GPUImageDataManager<GPUImage>);
186
187 /* Override DataHasBeenGenerated() in DataObject class.
188 * We need this because CPU time stamp is always bigger
189 * than GPU's. That is because Modified() is called at
190 * the end of each filter in the pipeline so although we
191 * increment GPU's time stamp in GPUGenerateData() the
192 * CPU's time stamp will be increased after that.
193 */
194 void
196 {
198 if (m_DataManager->IsCPUBufferDirty())
199 {
200 m_DataManager->Modified();
201 }
202 }
203
205 virtual void
206 Graft(const Self * data);
207
208protected:
209 void
210 Graft(const DataObject * data) override;
212 ~GPUImage() override;
213 using Superclass::Graft;
214
215private:
217};
218
219class ITK_TEMPLATE_EXPORT GPUImageFactory : public itk::ObjectFactoryBase
220{
221public:
222 ITK_DISALLOW_COPY_AND_MOVE(GPUImageFactory);
223
228
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 }
241
242
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.
virtual void DataHasBeenGenerated()
Provides a common API for pixel accessors for Image and VectorImage.
GPU memory manager implemented using OpenCL. Required by GPUImage class.
const char * GetDescription() const override
itk::SmartPointer< Self > Pointer
static Pointer New()
itk::SmartPointer< const Self > ConstPointer
GPUImageFactory Self
const char * GetITKSourceVersion() const override
static void RegisterOneFactory()
itk::ObjectFactoryBase Superclass
Templated n-dimensional image class for the GPU.
Definition itkGPUImage.h:41
void UpdateBuffers()
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
void Graft(const DataObject *data) override
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition itkGPUImage.h:69
DefaultPixelAccessorFunctor< Self > AccessorFunctorType
Definition itkGPUImage.h:72
GPUImage Self
Definition itkGPUImage.h:45
const TPixel & GetPixel(const IndexType &index) const
NeighborhoodAccessorFunctor< Self > NeighborhoodAccessorFunctorType
Definition itkGPUImage.h:74
int GetCurrentCommandQueueID()
TPixel * GetBufferPointer() override
const TPixel * GetBufferPointer() const override
const AccessorType GetPixelAccessor() const
Index< VImageDimension > IndexType
const PixelContainer * GetPixelContainer() const
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
PixelContainer * GetPixelContainer()
SmartPointer< const Self > ConstPointer
Definition itkGPUImage.h:48
void SetPixel(const IndexType &index, const TPixel &value)
void SetCurrentCommandQueue(int queueid)
virtual void Graft(const Self *data)
void Allocate(bool initialize=false) override
GPUDataManager * GetGPUDataManager()
TPixel & operator[](const IndexType &index)
~GPUImage() override
void DataHasBeenGenerated() override
void Initialize() override
AccessorType GetPixelAccessor()
const TPixel & operator[](const IndexType &index) const
SmartPointer< Self > Pointer
Definition itkGPUImage.h:47
WeakPointer< const Self > ConstWeakPointer
Definition itkGPUImage.h:49
GPUImageDataManager< GPUImage >::Pointer m_DataManager
Image< TPixel, VImageDimension > Superclass
Definition itkGPUImage.h:46
void SetPixelContainer(PixelContainer *container)
typename PixelContainer::Pointer PixelContainerPointer
Definition itkGPUImage.h:68
TPixel & GetPixel(const IndexType &index)
static constexpr unsigned int ImageDimension
Definition itkGPUImage.h:55
void FillBuffer(const TPixel &value)
GPUImage< TPixelType, VDimension > Type
Templated n-dimensional image class.
Definition itkImage.h:89
AccessorType GetPixelAccessor()
Definition itkImage.h:303
TPixel ValueType
Definition itkImage.h:111
virtual void Graft(const Self *image)
Size< VImageDimension > SizeType
TPixel InternalPixelType
Definition itkImage.h:117
ImageRegion< VImageDimension > RegionType
TPixel PixelType
Definition itkImage.h:108
PixelType IOPixelType
Definition itkImage.h:119
DefaultPixelAccessor< PixelType > AccessorType
Definition itkImage.h:123
Offset< VImageDimension > OffsetType
ImportImageContainer< SizeValueType, PixelType > PixelContainer
Definition itkImage.h:145
Vector< SpacingValueType, VImageDimension > SpacingType
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
PixelContainer * GetPixelContainer()
Definition itkImage.h:272
Index< VImageDimension > IndexType
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)
Implements transparent reference counting.
Implements a weak reference to an object.
#define OverrideImageTypeMacro(pt, dm)
#define ITK_SOURCE_VERSION
Definition itkVersion.h:39
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
bool IsGPUAvailable()