SimpleITK  
sitkPimpleImageBase.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 * 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 sitkPimpleImageBase_h
19 #define sitkPimpleImageBase_h
20 
21 #include <vector>
22 #include "sitkPixelIDTokens.h"
23 #include "sitkTemplateFunctions.h"
24 
25 namespace itk::simple
26  {
27 
39  {
40  public:
41 
42  virtual ~PimpleImageBase( ) = default;
43 
44  virtual std::unique_ptr<PimpleImageBase> ProxyCopy() = 0;
45 
46  virtual PixelIDValueEnum GetPixelID() const = 0;
47  virtual unsigned int GetDimension( ) const = 0;
48  virtual uint64_t GetNumberOfPixels( ) const = 0;
49  virtual unsigned int GetNumberOfComponentsPerPixel( ) const = 0;
50 
51  virtual std::unique_ptr<PimpleImageBase> ShallowCopy() const = 0;
52  virtual std::unique_ptr<PimpleImageBase> DeepCopy() const = 0;
53  virtual itk::DataObject* GetDataBase( ) = 0;
54  virtual const itk::DataObject* GetDataBase( ) const = 0;
55 
56  virtual unsigned int GetWidth( ) const { return this->GetSize( 0 ); }
57  virtual unsigned int GetHeight( ) const { return this->GetSize( 1 ); }
58  virtual unsigned int GetDepth( ) const { return this->GetSize( 2 ); }
59 
60  virtual std::vector< unsigned int > GetSize( ) const = 0;
61  virtual unsigned int GetSize( unsigned int dimension ) const = 0;
62 
63 
64  virtual std::vector<double> GetOrigin( ) const = 0;
65  virtual void SetOrigin( const std::vector<double> &orgn ) = 0;
66  virtual std::vector<double> GetSpacing( ) const = 0;
67  virtual void SetSpacing( const std::vector<double> &spc ) = 0;
68 
69  virtual std::vector< double > GetDirection( ) const = 0;
70  virtual void SetDirection( const std::vector< double > &direction ) = 0;
71 
72  virtual std::vector<int64_t> TransformPhysicalPointToIndex( const std::vector<double> &pt) const = 0;
73  virtual std::vector<double> TransformIndexToPhysicalPoint( const std::vector<int64_t> &idx) const = 0;
74  virtual std::vector<double> TransformPhysicalPointToContinuousIndex( const std::vector<double> &pt) const = 0;
75  virtual std::vector<double> TransformContinuousIndexToPhysicalPoint( const std::vector<double> &idx) const = 0;
76 
77  virtual std::vector<double> EvaluateAtContinuousIndex( const std::vector<double> &index, InterpolatorEnum interp) const = 0;
78 
79  virtual std::string ToString() const = 0;
80 
81  virtual int GetReferenceCountOfImage() const = 0;
82 
83  virtual int8_t GetPixelAsInt8( const std::vector<uint32_t> &idx) const = 0;
84  virtual uint8_t GetPixelAsUInt8( const std::vector<uint32_t> &idx) const = 0;
85  virtual int16_t GetPixelAsInt16( const std::vector<uint32_t> &idx ) const = 0;
86  virtual uint16_t GetPixelAsUInt16( const std::vector<uint32_t> &idx ) const = 0;
87  virtual int32_t GetPixelAsInt32( const std::vector<uint32_t> &idx ) const = 0;
88  virtual uint32_t GetPixelAsUInt32( const std::vector<uint32_t> &idx ) const = 0;
89  virtual int64_t GetPixelAsInt64( const std::vector<uint32_t> &idx ) const = 0;
90  virtual uint64_t GetPixelAsUInt64( const std::vector<uint32_t> &idx ) const = 0;
91  virtual float GetPixelAsFloat( const std::vector<uint32_t> &idx ) const = 0;
92  virtual double GetPixelAsDouble( const std::vector<uint32_t> &idx ) const = 0;
93 
94  virtual std::vector<int8_t> GetPixelAsVectorInt8( const std::vector<uint32_t> &idx) const = 0;
95  virtual std::vector<uint8_t> GetPixelAsVectorUInt8( const std::vector<uint32_t> &idx) const = 0;
96  virtual std::vector<int16_t> GetPixelAsVectorInt16( const std::vector<uint32_t> &idx ) const = 0;
97  virtual std::vector<uint16_t> GetPixelAsVectorUInt16( const std::vector<uint32_t> &idx ) const = 0;
98  virtual std::vector<int32_t> GetPixelAsVectorInt32( const std::vector<uint32_t> &idx ) const = 0;
99  virtual std::vector<uint32_t> GetPixelAsVectorUInt32( const std::vector<uint32_t> &idx ) const = 0;
100  virtual std::vector<int64_t> GetPixelAsVectorInt64( const std::vector<uint32_t> &idx ) const = 0;
101  virtual std::vector<uint64_t> GetPixelAsVectorUInt64( const std::vector<uint32_t> &idx ) const = 0;
102  virtual std::vector<float> GetPixelAsVectorFloat32( const std::vector<uint32_t> &idx ) const = 0;
103  virtual std::vector<double> GetPixelAsVectorFloat64( const std::vector<uint32_t> &idx ) const = 0;
104 
105  virtual std::complex<float> GetPixelAsComplexFloat32( const std::vector<uint32_t> &idx ) const = 0;
106  virtual std::complex<double> GetPixelAsComplexFloat64( const std::vector<uint32_t> &idx ) const = 0;
107 
108  virtual void SetPixelAsInt8( const std::vector<uint32_t> &idx, int8_t v ) = 0;
109  virtual void SetPixelAsUInt8( const std::vector<uint32_t> &idx, uint8_t v ) = 0;
110  virtual void SetPixelAsInt16( const std::vector<uint32_t> &idx, int16_t v ) = 0;
111  virtual void SetPixelAsUInt16( const std::vector<uint32_t> &idx, uint16_t v ) = 0;
112  virtual void SetPixelAsInt32( const std::vector<uint32_t> &idx, int32_t v ) = 0;
113  virtual void SetPixelAsUInt32( const std::vector<uint32_t> &idx, uint32_t v ) = 0;
114  virtual void SetPixelAsInt64( const std::vector<uint32_t> &idx, int64_t v ) = 0;
115  virtual void SetPixelAsUInt64( const std::vector<uint32_t> &idx, uint64_t v ) = 0;
116  virtual void SetPixelAsFloat( const std::vector<uint32_t> &idx, float v ) = 0;
117  virtual void SetPixelAsDouble( const std::vector<uint32_t> &idx, double v ) = 0;
118 
119  virtual void SetPixelAsVectorInt8( const std::vector<uint32_t> &idx, const std::vector<int8_t> &v ) = 0;
120  virtual void SetPixelAsVectorUInt8( const std::vector<uint32_t> &idx, const std::vector<uint8_t> &v ) = 0;
121  virtual void SetPixelAsVectorInt16( const std::vector<uint32_t> &idx, const std::vector<int16_t> &v ) = 0;
122  virtual void SetPixelAsVectorUInt16( const std::vector<uint32_t> &idx, const std::vector<uint16_t> &v ) = 0;
123  virtual void SetPixelAsVectorInt32( const std::vector<uint32_t> &idx, const std::vector<int32_t> &v ) = 0;
124  virtual void SetPixelAsVectorUInt32( const std::vector<uint32_t> &idx, const std::vector<uint32_t> &v ) = 0;
125  virtual void SetPixelAsVectorInt64( const std::vector<uint32_t> &idx, const std::vector<int64_t> &v ) = 0;
126  virtual void SetPixelAsVectorUInt64( const std::vector<uint32_t> &idx, const std::vector<uint64_t> &v ) = 0;
127  virtual void SetPixelAsVectorFloat32( const std::vector<uint32_t> &idx, const std::vector<float> &v ) = 0;
128  virtual void SetPixelAsVectorFloat64( const std::vector<uint32_t> &idx, const std::vector<double> &v ) = 0;
129 
130  virtual void SetPixelAsComplexFloat32( const std::vector<uint32_t> &idx, const std::complex<float> v ) = 0;
131  virtual void SetPixelAsComplexFloat64( const std::vector<uint32_t> &idx, const std::complex<double> v ) = 0;
132 
133 
134  virtual int8_t *GetBufferAsInt8( ) = 0;
135  virtual uint8_t *GetBufferAsUInt8( ) = 0;
136  virtual int16_t *GetBufferAsInt16( ) = 0;
137  virtual uint16_t *GetBufferAsUInt16( ) = 0;
138  virtual int32_t *GetBufferAsInt32( ) = 0;
139  virtual uint32_t *GetBufferAsUInt32( ) = 0;
140  virtual int64_t *GetBufferAsInt64( ) = 0;
141  virtual uint64_t *GetBufferAsUInt64( ) = 0;
142  virtual float *GetBufferAsFloat( ) = 0;
143  virtual double *GetBufferAsDouble( ) = 0;
144  virtual void *GetBufferAsVoid( ) = 0;
145 
146  virtual const int8_t *GetBufferAsInt8( ) const = 0;
147  virtual const uint8_t *GetBufferAsUInt8( ) const = 0;
148  virtual const int16_t *GetBufferAsInt16( ) const = 0;
149  virtual const uint16_t *GetBufferAsUInt16( ) const = 0;
150  virtual const int32_t *GetBufferAsInt32( ) const = 0;
151  virtual const uint32_t *GetBufferAsUInt32( ) const = 0;
152  virtual const int64_t *GetBufferAsInt64( ) const = 0;
153  virtual const uint64_t *GetBufferAsUInt64( ) const = 0;
154  virtual const float *GetBufferAsFloat( ) const = 0;
155  virtual const double *GetBufferAsDouble( ) const = 0;
156  virtual const void *GetBufferAsVoid( ) const = 0;
157  };
158 
159  }
160 
161 
162 #endif // sitkPimpleImageBase_h
itk::simple::InterpolatorEnum
InterpolatorEnum
Definition: sitkInterpolator.h:28
sitkTemplateFunctions.h
itk::simple::PimpleImageBase
Private implementation idiom image base class.
Definition: sitkPimpleImageBase.h:38
itk::simple::PixelIDValueEnum
PixelIDValueEnum
Enumerated values of pixelIDs.
Definition: sitkPixelIDValues.h:91
SITKCommon_HIDDEN
#define SITKCommon_HIDDEN
Definition: sitkCommon.h:44
sitkPixelIDTokens.h
itk::simple::PimpleImageBase::GetWidth
virtual unsigned int GetWidth() const
Definition: sitkPimpleImageBase.h:56
itk::simple::PimpleImageBase::GetHeight
virtual unsigned int GetHeight() const
Definition: sitkPimpleImageBase.h:57
itk::simple
Definition: sitkAdditionalProcedures.h:28
itk::simple::PimpleImageBase::GetDepth
virtual unsigned int GetDepth() const
Definition: sitkPimpleImageBase.h:58
itk::DataObject