go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkCUDAResampleImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright UMC Utrecht and contributors
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 #ifndef __itkCUDAResamplerImageFilter_h
19 #define __itkCUDAResamplerImageFilter_h
20 
21 #include "itkImage.h"
22 #include "itkResampleImageFilter.h"
25 #include "itkBSplineDeformableTransform.h"
26 #include "cudaResampleImageFilter.cuh"
27 
28 namespace itk
29 {
30 
43 template< typename TInputImage, typename TOutputImage, typename TInterpolatorPrecisionType = float >
44 class ITK_EXPORT itkCUDAResampleImageFilter :
45  public ResampleImageFilter< TInputImage, TOutputImage, TInterpolatorPrecisionType >
46 {
47 public:
48 
51  typedef ResampleImageFilter<
52  TInputImage, TOutputImage, TInterpolatorPrecisionType > Superclass;
54  typedef SmartPointer< const Self > ConstPointer;
55 
57  itkNewMacro( Self );
58 
60  itkTypeMacro( itkCUDAResampleImageFilter, ResampleImageFilter );
61 
63  typedef typename Superclass::InputImageType InputImageType;
64  typedef typename Superclass::OutputImageType OutputImageType;
65  typedef typename Superclass::InputImagePointer InputImagePointer;
66  typedef typename Superclass::InputImageConstPointer InputImageConstPointer;
67  typedef typename Superclass::OutputImagePointer OutputImagePointer;
68  typedef typename Superclass::InputImageRegionType InputImageRegionType;
69 
70  typedef typename Superclass::TransformType TransformType;
71  typedef typename Superclass::TransformPointerType TransformPointerType;
72  typedef typename Superclass::InterpolatorType InterpolatorType;
73  typedef typename Superclass::InterpolatorPointerType InterpolatorPointerType;
74 
75  typedef typename Superclass::SizeType SizeType;
76  typedef typename Superclass::IndexType IndexType;
77  typedef typename Superclass::PointType PointType;
78  typedef typename Superclass::PixelType PixelType;
79  typedef typename Superclass::InputPixelType InputPixelType;
80  typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
81  typedef typename Superclass::SpacingType SpacingType;
82  typedef typename Superclass::OriginPointType OriginPointType;
83  typedef typename Superclass::DirectionType DirectionType;
84  typedef typename Superclass::ImageBaseType ImageBaseType;
85 
88  TInterpolatorPrecisionType, 3 > InternalComboTransformType;
90  TInterpolatorPrecisionType, 3, 3 > InternalAdvancedBSplineTransformType;
93  typedef BSplineDeformableTransform<
94  TInterpolatorPrecisionType, 3, 3 > InternalBSplineTransformType;
95  typedef cuda::CUDAResampleImageFilter<
96  typename InternalBSplineTransformType::ParametersValueType,
97  typename TInputImage::PixelType, float > CudaResampleImageFilterType;
98 
100  itkSetMacro( UseCuda, bool );
101  itkGetConstMacro( UseCuda, bool );
102  itkBooleanMacro( UseCuda );
103 
105  itkSetMacro( UseGPUToCastData, bool );
106  itkGetConstMacro( UseGPUToCastData, bool );
107  itkBooleanMacro( UseGPUToCastData );
108 
110  itkSetMacro( UseFastCUDAKernel, bool );
111  itkGetConstMacro( UseFastCUDAKernel, bool );
112  itkBooleanMacro( UseFastCUDAKernel );
113 
115  virtual void GenerateData( void );
116 
119  {
120 public:
121 
122  std::vector< std::string > m_Warnings;
123 
124  void ResetWarningReport( void )
125  {
126  this->m_Warnings.resize( 0 );
127  }
128 
129 
130  std::string GetWarningReportAsString( void ) const
131  {
132  std::string warnings = "\n---------------------------------\n";
133  for( std::size_t i = 0; i < this->m_Warnings.size(); i++ )
134  {
135  warnings += "itkCUDAResampleImageFilter: " + this->m_Warnings[ i ];
136  warnings += "\n---------------------------------\n";
137  }
138  return warnings;
139  }
140 
141 
142  };
143 
144  //itkGetConstReferenceMacro( WarningReport, WarningReportType );
145  virtual const WarningReportType & GetWarningReport( void ) const
146  {
147  return this->m_WarningReport;
148  }
149 
150 
151 protected:
152 
155 
156  virtual void CheckForValidConfiguration( ValidTransformPointer & bSplineTransform );
157 
158 private:
159 
161  bool m_UseCuda;
164 
165  CudaResampleImageFilterType m_CudaResampleImageFilter;
166  WarningReportType m_WarningReport;
167 
172  bool CheckForValidTransform( ValidTransformPointer & bSplineTransform ) const;
173 
177  bool CheckForValidInterpolator( void ) const;
178 
182  bool CheckForValidDirectionCosines( ValidTransformPointer bSplineTransform ); //const;
183 
184  // NOTE: const can be added again in ITK4. It's due to GetInput() being not const-correct.
185 
187  void CopyParameters( ValidTransformPointer bSplineTransform );
188 
189 };
190 
191 // end class itkCUDAResampleImageFilter
192 
193 } // end namespace itk
194 
195 #ifndef ITK_MANUAL_INSTANTIATION
196 #include "itkCUDAResampleImageFilter.hxx"
197 #endif
198 
199 #endif // end #ifndef __itkCUDAResamplerImageFilter_h
This class combines two transforms: an 'initial transform' with a 'current transform'.
Superclass::OutputImageType OutputImageType
Superclass::TransformPointerType TransformPointerType
Deformable transform using a B-spline representation.
Superclass::InterpolatorPointerType InterpolatorPointerType
Resample an image on the GPU via a coordinate transform.
Superclass::InputImageRegionType InputImageRegionType
BSplineDeformableTransform< TInterpolatorPrecisionType, 3, 3 > InternalBSplineTransformType
Superclass::InputImagePointer InputImagePointer
virtual const WarningReportType & GetWarningReport(void) const
CudaResampleImageFilterType m_CudaResampleImageFilter
InternalAdvancedBSplineTransformType::ConstPointer ValidTransformConstPointer
AdvancedCombinationTransform< TInterpolatorPrecisionType, 3 > InternalComboTransformType
InternalAdvancedBSplineTransformType::Pointer ValidTransformPointer
Superclass::OutputImagePointer OutputImagePointer
Superclass::InputImageType InputImageType
Superclass::InputImageConstPointer InputImageConstPointer
Superclass::OriginPointType OriginPointType
Superclass::OutputImageRegionType OutputImageRegionType
SmartPointer< const Self > ConstPointer
Superclass::InputPixelType InputPixelType
AdvancedBSplineDeformableTransform< TInterpolatorPrecisionType, 3, 3 > InternalAdvancedBSplineTransformType
Superclass::InterpolatorType InterpolatorType
cuda::CUDAResampleImageFilter< typename InternalBSplineTransformType::ParametersValueType, typename TInputImage::PixelType, float > CudaResampleImageFilterType
ResampleImageFilter< TInputImage, TOutputImage, TInterpolatorPrecisionType > Superclass


Generated on 04-09-2015 for elastix by doxygen 1.8.9.1 elastix logo