ITK  5.4.0
Insight Toolkit
itkGPUImageToImageFilter.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#ifndef itkGPUImageToImageFilter_h
19#define itkGPUImageToImageFilter_h
20
22#include "itkGPUKernelManager.h"
23
24namespace itk
25{
38template <typename TInputImage,
39 typename TOutputImage,
40 typename TParentImageFilter = ImageToImageFilter<TInputImage, TOutputImage>>
41class ITK_TEMPLATE_EXPORT GPUImageToImageFilter : public TParentImageFilter
42{
43public:
44 ITK_DISALLOW_COPY_AND_MOVE(GPUImageToImageFilter);
45
48 using Superclass = TParentImageFilter;
51
52 itkNewMacro(Self);
53
55 itkOverrideGetNameOfClassMacro(GPUImageToImageFilter);
56
59 using typename Superclass::OutputImageRegionType;
60 using typename Superclass::OutputImagePixelType;
61
63 using InputImageType = TInputImage;
67 using InputImagePixelType = typename InputImageType::PixelType;
68
70 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
71 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
72
73 // macro to set if GPU is used
74 itkSetMacro(GPUEnabled, bool);
75 itkGetConstMacro(GPUEnabled, bool);
76 itkBooleanMacro(GPUEnabled);
77
78 void
79 GenerateData() override;
80 virtual void
82 virtual void
84
85protected:
86 void
87 GraftOutput(DataObject * output) override;
88 void
89 GraftOutput(const DataObjectIdentifierType & key, DataObject * output) override;
92
93 void
94 PrintSelf(std::ostream & os, Indent indent) const override;
95
96 virtual void
98 {}
99
100 // GPU kernel manager
101 typename GPUKernelManager::Pointer m_GPUKernelManager{};
102
103 // GPU kernel handle - kernel should be defined in specific filter (not in the
104 // base class)
105 // int m_KernelHandle;
106
107private:
108 bool m_GPUEnabled{ true };
109};
110
111} // end namespace itk
112
113#ifndef ITK_MANUAL_INSTANTIATION
114# include "itkGPUImageToImageFilter.hxx"
115#endif
116
117#endif
Base class for all data objects in ITK.
class to abstract the behaviour of the GPU filters.
void GenerateData() override
void GraftOutput(DataObject *output) override
virtual void GraftOutput(const DataObjectIdentifierType &key, typename itk::GPUTraits< TOutputImage >::Type *output)
void GraftOutput(const DataObjectIdentifierType &key, DataObject *output) override
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void GraftOutput(typename itk::GPUTraits< TOutputImage >::Type *output)
Base class for all process objects that output image data.
Superclass::DataObjectIdentifierType DataObjectIdentifierType
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::Pointer InputImagePointer
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::RegionType InputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
DataObject::DataObjectIdentifierType DataObjectIdentifierType
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....