File: itkGPUImage.h

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (309 lines) | stat: -rw-r--r-- 8,354 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/

#ifndef itkGPUImage_h
#define itkGPUImage_h

#include "itkImage.h"
#include "itkGPUImageDataManager.h"
#include "itkVersion.h"
#include "itkObjectFactoryBase.h"

namespace itk
{
/** \class GPUImage
 *  \brief Templated n-dimensional image class for the GPU.
 *
 * Derived from itk Image class to use with GPU image filters.
 * This class manages both CPU and GPU memory implicitly, and
 * can be used with non-GPU itk filters as well. Memory transfer
 * between CPU and GPU is done automatically and implicitly.
 *
 * \ingroup ITKGPUCommon
 */
template <typename TPixel, unsigned int VImageDimension = 2>
class ITK_TEMPLATE_EXPORT GPUImage : public Image<TPixel, VImageDimension>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(GPUImage);

  using Self = GPUImage;
  using Superclass = Image<TPixel, VImageDimension>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;
  using ConstWeakPointer = WeakPointer<const Self>;

  itkNewMacro(Self);

  itkOverrideGetNameOfClassMacro(GPUImage);

  static constexpr unsigned int ImageDimension = VImageDimension;

  using typename Superclass::PixelType;
  using typename Superclass::ValueType;
  using typename Superclass::InternalPixelType;
  using typename Superclass::IOPixelType;
  using typename Superclass::DirectionType;
  using typename Superclass::SpacingType;
  using typename Superclass::PixelContainer;
  using typename Superclass::SizeType;
  using typename Superclass::IndexType;
  using typename Superclass::OffsetType;
  using typename Superclass::RegionType;
  using PixelContainerPointer = typename PixelContainer::Pointer;
  using PixelContainerConstPointer = typename PixelContainer::ConstPointer;
  using typename Superclass::AccessorType;

  using AccessorFunctorType = DefaultPixelAccessorFunctor<Self>;

  using NeighborhoodAccessorFunctorType = NeighborhoodAccessorFunctor<Self>;
  // NeighborhoodAccessorFunctorType;

  //
  // Allocate CPU and GPU memory space
  //
  void
  Allocate(bool initialize = false) override;

  void
  Initialize() override;

  void
  FillBuffer(const TPixel & value);

  void
  SetPixel(const IndexType & index, const TPixel & value);

  const TPixel &
  GetPixel(const IndexType & index) const;

  TPixel &
  GetPixel(const IndexType & index);

  const TPixel & operator[](const IndexType & index) const;

  TPixel & operator[](const IndexType & index);

  /** Explicit synchronize CPU/GPU buffers */
  void
  UpdateBuffers();

  //
  // Get CPU buffer pointer
  //
  TPixel *
  GetBufferPointer() override;

  const TPixel *
  GetBufferPointer() const override;

  /** Return the Pixel Accessor object */
  AccessorType
  GetPixelAccessor()
  {
    m_DataManager->SetGPUBufferDirty();
    return Superclass::GetPixelAccessor();
  }

  /** Return the Pixel Accesor object */
  const AccessorType
  GetPixelAccessor() const
  {
    m_DataManager->UpdateCPUBuffer();
    return Superclass::GetPixelAccessor();
  }

  /** Return the NeighborhoodAccessor functor */
  NeighborhoodAccessorFunctorType
  GetNeighborhoodAccessor()
  {
    m_DataManager->SetGPUBufferDirty();
    // return Superclass::GetNeighborhoodAccessor();
    return NeighborhoodAccessorFunctorType();
  }

  /** Return the NeighborhoodAccessor functor */
  const NeighborhoodAccessorFunctorType
  GetNeighborhoodAccessor() const
  {
    m_DataManager->UpdateCPUBuffer();
    // return Superclass::GetNeighborhoodAccessor();
    return NeighborhoodAccessorFunctorType();
  }

  void
  SetPixelContainer(PixelContainer * container);

  /** Return a pointer to the container. */
  PixelContainer *
  GetPixelContainer()
  {
    m_DataManager->SetGPUBufferDirty();
    return Superclass::GetPixelContainer();
  }

  const PixelContainer *
  GetPixelContainer() const
  {
    m_DataManager->UpdateCPUBuffer();
    return Superclass::GetPixelContainer();
  }

  void
  SetCurrentCommandQueue(int queueid)
  {
    m_DataManager->SetCurrentCommandQueue(queueid);
  }

  int
  GetCurrentCommandQueueID()
  {
    return m_DataManager->GetCurrentCommandQueueID();
  }

  itkGetModifiableObjectMacro(DataManager, GPUImageDataManager<GPUImage>);
  GPUDataManager *
  GetGPUDataManager();

  /* Override DataHasBeenGenerated() in DataObject class.
   * We need this because CPU time stamp is always bigger
   * than GPU's. That is because Modified() is called at
   * the end of each filter in the pipeline so although we
   * increment GPU's time stamp in GPUGenerateData() the
   * CPU's time stamp will be increased after that.
   */
  void
  DataHasBeenGenerated() override
  {
    Superclass::DataHasBeenGenerated();
    if (m_DataManager->IsCPUBufferDirty())
    {
      m_DataManager->Modified();
    }
  }

  /** Graft the data and information from one GPUImage to another. */
  virtual void
  Graft(const Self * data);

protected:
  void
  Graft(const DataObject * data) override;
  GPUImage();
  ~GPUImage() override;
  using Superclass::Graft;

private:
  typename GPUImageDataManager<GPUImage>::Pointer m_DataManager{};
};

class ITK_TEMPLATE_EXPORT GPUImageFactory : public itk::ObjectFactoryBase
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(GPUImageFactory);

  using Self = GPUImageFactory;
  using Superclass = itk::ObjectFactoryBase;
  using Pointer = itk::SmartPointer<Self>;
  using ConstPointer = itk::SmartPointer<const Self>;

  /** Class methods used to interface with the registered factories. */
  const char *
  GetITKSourceVersion() const override
  {
    return ITK_SOURCE_VERSION;
  }
  const char *
  GetDescription() const override
  {
    return "A Factory for GPUImage";
  }

  /** Method for class instantiation. */
  itkFactorylessNewMacro(Self);

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(GPUImageFactory);

  /** Register one factory of this type  */
  static void
  RegisterOneFactory()
  {
    auto factory = GPUImageFactory::New();

    itk::ObjectFactoryBase::RegisterFactory(factory);
  }

private:
#define OverrideImageTypeMacro(pt, dm)                         \
  this->RegisterOverride(typeid(itk::Image<pt, dm>).name(),    \
                         typeid(itk::GPUImage<pt, dm>).name(), \
                         "GPU Image Override",                 \
                         true,                                 \
                         itk::CreateObjectFunction<GPUImage<pt, dm>>::New())

  GPUImageFactory()
  {
    if (IsGPUAvailable())
    {
      // 1/2/3D
      OverrideImageTypeMacro(unsigned char, 1);
      OverrideImageTypeMacro(signed char, 1);
      OverrideImageTypeMacro(int, 1);
      OverrideImageTypeMacro(unsigned int, 1);
      OverrideImageTypeMacro(float, 1);
      OverrideImageTypeMacro(double, 1);

      OverrideImageTypeMacro(unsigned char, 2);
      OverrideImageTypeMacro(signed char, 2);
      OverrideImageTypeMacro(int, 2);
      OverrideImageTypeMacro(unsigned int, 2);
      OverrideImageTypeMacro(float, 2);
      OverrideImageTypeMacro(double, 2);

      OverrideImageTypeMacro(unsigned char, 3);
      OverrideImageTypeMacro(signed char, 3);
      OverrideImageTypeMacro(int, 3);
      OverrideImageTypeMacro(unsigned int, 3);
      OverrideImageTypeMacro(float, 3);
      OverrideImageTypeMacro(double, 3);
    }
  }
};

template <typename T>
class ITK_TEMPLATE_EXPORT GPUTraits
{
public:
  using Type = T;
};

template <typename TPixelType, unsigned int VDimension>
class ITK_TEMPLATE_EXPORT GPUTraits<Image<TPixelType, VDimension>>
{
public:
  using Type = GPUImage<TPixelType, VDimension>;
};

} // end namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkGPUImage.hxx"
#endif

#endif