ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkGPUImage.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 
19 #ifndef itkGPUImage_h
20 #define itkGPUImage_h
21 
22 #include "itkImage.h"
23 #include "itkGPUImageDataManager.h"
24 #include "itkVersion.h"
25 #include "itkObjectFactoryBase.h"
26 
27 namespace itk
28 {
39 template <typename TPixel, unsigned int VImageDimension = 2>
40 class ITK_TEMPLATE_EXPORT GPUImage : public Image<TPixel,VImageDimension>
41 {
42 public:
43  ITK_DISALLOW_COPY_AND_ASSIGN(GPUImage);
44 
45  using Self = GPUImage;
50 
51  itkNewMacro(Self);
52 
53  itkTypeMacro(GPUImage, Image);
54 
55  static constexpr unsigned int ImageDimension = VImageDimension;
56 
57  using PixelType = typename Superclass::PixelType;
58  using ValueType = typename Superclass::ValueType;
59  using InternalPixelType = typename Superclass::InternalPixelType;
60  using IOPixelType = typename Superclass::IOPixelType;
62  using SpacingType = typename Superclass::SpacingType;
63  using PixelContainer = typename Superclass::PixelContainer;
64  using SizeType = typename Superclass::SizeType;
65  using IndexType = typename Superclass::IndexType;
66  using OffsetType = typename Superclass::OffsetType;
70  using AccessorType = typename Superclass::AccessorType;
71 
73 
75  // NeighborhoodAccessorFunctorType;
76 
77  //
78  // Allocate CPU and GPU memory space
79  //
80  void Allocate(bool initialize=false) override;
81 
82  void Initialize() override;
83 
84  void FillBuffer(const TPixel & value);
85 
86  void SetPixel(const IndexType & index, const TPixel & value);
87 
88  const TPixel & GetPixel(const IndexType & index) const;
89 
90  TPixel & GetPixel(const IndexType & index);
91 
92  const TPixel & operator[](const IndexType & index) const;
93 
94  TPixel & operator[](const IndexType & index);
95 
97  void UpdateBuffers();
98 
99  //
100  // Get CPU buffer pointer
101  //
102  TPixel* GetBufferPointer() override;
103 
104  const TPixel * GetBufferPointer() const override;
105 
108  {
109  m_DataManager->SetGPUBufferDirty();
110  return Superclass::GetPixelAccessor();
111  }
113 
116  {
117  m_DataManager->UpdateCPUBuffer();
118  return Superclass::GetPixelAccessor();
119  }
121 
124  {
125  m_DataManager->SetGPUBufferDirty();
126  //return Superclass::GetNeighborhoodAccessor();
128  }
130 
133  {
134  m_DataManager->UpdateCPUBuffer();
135  //return Superclass::GetNeighborhoodAccessor();
137  }
139 
140  void SetPixelContainer(PixelContainer *container);
141 
144  {
145  m_DataManager->SetGPUBufferDirty(); return Superclass::GetPixelContainer();
146  }
147 
149  {
150  m_DataManager->UpdateCPUBuffer();
151  return Superclass::GetPixelContainer();
152  }
153 
154  void SetCurrentCommandQueue( int queueid )
155  {
156  m_DataManager->SetCurrentCommandQueue( queueid );
157  }
158 
160  return m_DataManager->GetCurrentCommandQueueID();
161  }
162 
163  itkGetModifiableObjectMacro(DataManager, GPUImageDataManager< GPUImage >);
164  GPUDataManager * GetGPUDataManager();
165 
166  /* Override DataHasBeenGenerated() in DataObject class.
167  * We need this because CPU time stamp is always bigger
168  * than GPU's. That is because Modified() is called at
169  * the end of each filter in the pipeline so although we
170  * increment GPU's time stamp in GPUGenerateData() the
171  * CPU's time stamp will be increased after that.
172  */
173  void DataHasBeenGenerated() override
174  {
175  Superclass::DataHasBeenGenerated();
176  if( m_DataManager->IsCPUBufferDirty() )
177  {
178  m_DataManager->Modified();
179  }
180  }
181 
183  virtual void Graft(const Self *data);
184 
185 protected:
186  void Graft(const DataObject *data) override;
187  GPUImage();
188  ~GPUImage() override;
189  using Superclass::Graft;
190 
191 private:
193 };
194 
195 class ITK_TEMPLATE_EXPORT GPUImageFactory : public itk::ObjectFactoryBase
196 {
197 public:
198  ITK_DISALLOW_COPY_AND_ASSIGN(GPUImageFactory);
199 
204 
206  const char* GetITKSourceVersion() const override {
207  return ITK_SOURCE_VERSION;
208  }
209  const char* GetDescription() const override {
210  return "A Factory for GPUImage";
211  }
213 
215  itkFactorylessNewMacro(Self);
216 
219 
221  static void RegisterOneFactory()
222  {
224 
226  }
227 
228 private:
229 #define OverrideImageTypeMacro(pt,dm) this->RegisterOverride( \
230  typeid(itk::Image<pt,dm>).name(), \
231  typeid(itk::GPUImage<pt,dm>).name(), \
232  "GPU Image Override", \
233  true, \
234  itk::CreateObjectFunction<GPUImage<pt,dm> >::New() )
235 
237  {
238  if( IsGPUAvailable() )
239  {
240  // 1/2/3D
241  OverrideImageTypeMacro(unsigned char, 1);
242  OverrideImageTypeMacro(signed char, 1);
243  OverrideImageTypeMacro(int, 1);
244  OverrideImageTypeMacro(unsigned int, 1);
245  OverrideImageTypeMacro(float, 1);
246  OverrideImageTypeMacro(double, 1);
247 
248  OverrideImageTypeMacro(unsigned char, 2);
249  OverrideImageTypeMacro(signed char, 2);
250  OverrideImageTypeMacro(int, 2);
251  OverrideImageTypeMacro(unsigned int, 2);
252  OverrideImageTypeMacro(float, 2);
253  OverrideImageTypeMacro(double, 2);
254 
255  OverrideImageTypeMacro(unsigned char, 3);
256  OverrideImageTypeMacro(signed char, 3);
257  OverrideImageTypeMacro(int, 3);
258  OverrideImageTypeMacro(unsigned int, 3);
259  OverrideImageTypeMacro(float, 3);
260  OverrideImageTypeMacro(double, 3);
261  }
262  }
263 
264 };
265 
266 template <typename T>
267 class ITK_TEMPLATE_EXPORT GPUTraits
268 {
269 public:
270  using Type = T;
271 };
272 
273 template <typename TPixelType, unsigned int NDimension>
274 class ITK_TEMPLATE_EXPORT GPUTraits< Image< TPixelType, NDimension > >
275 {
276 public:
278 };
279 
280 } // end namespace itk
281 
282 #ifndef ITK_MANUAL_INSTANTIATION
283 #include "itkGPUImage.hxx"
284 #endif
285 
286 #endif
static void RegisterOneFactory()
Definition: itkGPUImage.h:221
void DataHasBeenGenerated() override
Definition: itkGPUImage.h:173
TPixel PixelType
Definition: itkImage.h:95
Light weight base class for most itk classes.
PixelType IOPixelType
Definition: itkImage.h:106
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:40
static Pointer New()
const char * GetDescription() const override
Definition: itkGPUImage.h:209
Create instances of classes using an object factory.
Implements a weak reference to an object.
void SetCurrentCommandQueue(int queueid)
Definition: itkGPUImage.h:154
GPU memory manager implemented using OpenCL. Required by GPUImage class.
TPixel ValueType
Definition: itkImage.h:98
const AccessorType GetPixelAccessor() const
Definition: itkGPUImage.h:115
PixelContainer * GetPixelContainer()
Definition: itkGPUImage.h:143
Templated n-dimensional image class for the GPU.
Definition: itkGPUImage.h:40
AccessorType GetPixelAccessor()
Definition: itkGPUImage.h:107
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
const PixelContainer * GetPixelContainer() const
Definition: itkGPUImage.h:148
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition: itkGPUImage.h:132
typename Superclass::RegionType RegionType
Definition: itkImage.h:139
typename Superclass::SizeType SizeType
Definition: itkImage.h:128
typename Superclass::SpacingType SpacingType
Definition: itkImage.h:143
GPUImageDataManager< GPUImage >::Pointer m_DataManager
Definition: itkGPUImage.h:192
Provides a common API for pixel accessors for Image and VectorImage.
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition: itkGPUImage.h:123
typename Superclass::IndexType IndexType
Definition: itkImage.h:121
#define OverrideImageTypeMacro(pt, dm)
Definition: itkGPUImage.h:229
int GetCurrentCommandQueueID()
Definition: itkGPUImage.h:159
static bool RegisterFactory(ObjectFactoryBase *, InsertionPositionType where=INSERT_AT_BACK, vcl_size_t position=0)
typename Superclass::DirectionType DirectionType
Definition: itkImage.h:135
bool IsGPUAvailable()
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkImage.h:151
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition: itkImage.h:152
const char * GetITKSourceVersion() const override
Definition: itkGPUImage.h:206
Give access to partial aspects a type.
typename Superclass::OffsetType OffsetType
Definition: itkImage.h:125
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
Defines an itk::Image front-end to a standard C-array.
TPixel InternalPixelType
Definition: itkImage.h:104