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  virtual ~PimpleImageBase() = default;
42 
43  virtual std::unique_ptr<PimpleImageBase>
44  ProxyCopy() = 0;
45 
46  virtual PixelIDValueEnum
47  GetPixelID() const = 0;
48  virtual unsigned int
49  GetDimension() const = 0;
50  virtual uint64_t
51  GetNumberOfPixels() const = 0;
52  virtual unsigned int
53  GetNumberOfComponentsPerPixel() const = 0;
54 
55  virtual std::unique_ptr<PimpleImageBase>
56  ShallowCopy() const = 0;
57  virtual std::unique_ptr<PimpleImageBase>
58  DeepCopy() const = 0;
59  virtual itk::DataObject *
60  GetDataBase() = 0;
61  virtual const itk::DataObject *
62  GetDataBase() const = 0;
63 
64  virtual unsigned int
65  GetWidth() const
66  {
67  return this->GetSize(0);
68  }
69  virtual unsigned int
70  GetHeight() const
71  {
72  return this->GetSize(1);
73  }
74  virtual unsigned int
75  GetDepth() const
76  {
77  return this->GetSize(2);
78  }
79 
80  virtual std::vector<unsigned int>
81  GetSize() const = 0;
82  virtual unsigned int
83  GetSize(unsigned int dimension) const = 0;
84 
85  virtual bool
86  IsCongruentImageGeometry(const PimpleImageBase * otherImage,
87  double coordinateTolerance,
88  double directionTolerance) const = 0;
89 
90  virtual bool
91  IsSameImageGeometryAs(const PimpleImageBase * otherImage,
92  double coordinateTolerance,
93  double directionTolerance) const = 0;
94 
95  virtual std::vector<double>
96  GetOrigin() const = 0;
97  virtual void
98  SetOrigin(const std::vector<double> & orgn) = 0;
99  virtual std::vector<double>
100  GetSpacing() const = 0;
101  virtual void
102  SetSpacing(const std::vector<double> & spc) = 0;
103 
104  virtual std::vector<double>
105  GetDirection() const = 0;
106  virtual void
107  SetDirection(const std::vector<double> & direction) = 0;
108 
109  virtual std::vector<int64_t>
110  TransformPhysicalPointToIndex(const std::vector<double> & pt) const = 0;
111  virtual std::vector<double>
112  TransformIndexToPhysicalPoint(const std::vector<int64_t> & idx) const = 0;
113  virtual std::vector<double>
114  TransformPhysicalPointToContinuousIndex(const std::vector<double> & pt) const = 0;
115  virtual std::vector<double>
116  TransformContinuousIndexToPhysicalPoint(const std::vector<double> & idx) const = 0;
117 
118  virtual std::vector<double>
119  EvaluateAtContinuousIndex(const std::vector<double> & index, InterpolatorEnum interp) const = 0;
120 
121  virtual std::string
122  ToString() const = 0;
123 
124  virtual int
125  GetReferenceCountOfImage() const = 0;
126 
127  virtual int8_t
128  GetPixelAsInt8(const std::vector<uint32_t> & idx) const = 0;
129  virtual uint8_t
130  GetPixelAsUInt8(const std::vector<uint32_t> & idx) const = 0;
131  virtual int16_t
132  GetPixelAsInt16(const std::vector<uint32_t> & idx) const = 0;
133  virtual uint16_t
134  GetPixelAsUInt16(const std::vector<uint32_t> & idx) const = 0;
135  virtual int32_t
136  GetPixelAsInt32(const std::vector<uint32_t> & idx) const = 0;
137  virtual uint32_t
138  GetPixelAsUInt32(const std::vector<uint32_t> & idx) const = 0;
139  virtual int64_t
140  GetPixelAsInt64(const std::vector<uint32_t> & idx) const = 0;
141  virtual uint64_t
142  GetPixelAsUInt64(const std::vector<uint32_t> & idx) const = 0;
143  virtual float
144  GetPixelAsFloat(const std::vector<uint32_t> & idx) const = 0;
145  virtual double
146  GetPixelAsDouble(const std::vector<uint32_t> & idx) const = 0;
147 
148  virtual std::vector<int8_t>
149  GetPixelAsVectorInt8(const std::vector<uint32_t> & idx) const = 0;
150  virtual std::vector<uint8_t>
151  GetPixelAsVectorUInt8(const std::vector<uint32_t> & idx) const = 0;
152  virtual std::vector<int16_t>
153  GetPixelAsVectorInt16(const std::vector<uint32_t> & idx) const = 0;
154  virtual std::vector<uint16_t>
155  GetPixelAsVectorUInt16(const std::vector<uint32_t> & idx) const = 0;
156  virtual std::vector<int32_t>
157  GetPixelAsVectorInt32(const std::vector<uint32_t> & idx) const = 0;
158  virtual std::vector<uint32_t>
159  GetPixelAsVectorUInt32(const std::vector<uint32_t> & idx) const = 0;
160  virtual std::vector<int64_t>
161  GetPixelAsVectorInt64(const std::vector<uint32_t> & idx) const = 0;
162  virtual std::vector<uint64_t>
163  GetPixelAsVectorUInt64(const std::vector<uint32_t> & idx) const = 0;
164  virtual std::vector<float>
165  GetPixelAsVectorFloat32(const std::vector<uint32_t> & idx) const = 0;
166  virtual std::vector<double>
167  GetPixelAsVectorFloat64(const std::vector<uint32_t> & idx) const = 0;
168 
169  virtual std::complex<float>
170  GetPixelAsComplexFloat32(const std::vector<uint32_t> & idx) const = 0;
171  virtual std::complex<double>
172  GetPixelAsComplexFloat64(const std::vector<uint32_t> & idx) const = 0;
173 
174  virtual void
175  SetPixelAsInt8(const std::vector<uint32_t> & idx, int8_t v) = 0;
176  virtual void
177  SetPixelAsUInt8(const std::vector<uint32_t> & idx, uint8_t v) = 0;
178  virtual void
179  SetPixelAsInt16(const std::vector<uint32_t> & idx, int16_t v) = 0;
180  virtual void
181  SetPixelAsUInt16(const std::vector<uint32_t> & idx, uint16_t v) = 0;
182  virtual void
183  SetPixelAsInt32(const std::vector<uint32_t> & idx, int32_t v) = 0;
184  virtual void
185  SetPixelAsUInt32(const std::vector<uint32_t> & idx, uint32_t v) = 0;
186  virtual void
187  SetPixelAsInt64(const std::vector<uint32_t> & idx, int64_t v) = 0;
188  virtual void
189  SetPixelAsUInt64(const std::vector<uint32_t> & idx, uint64_t v) = 0;
190  virtual void
191  SetPixelAsFloat(const std::vector<uint32_t> & idx, float v) = 0;
192  virtual void
193  SetPixelAsDouble(const std::vector<uint32_t> & idx, double v) = 0;
194 
195  virtual void
196  SetPixelAsVectorInt8(const std::vector<uint32_t> & idx, const std::vector<int8_t> & v) = 0;
197  virtual void
198  SetPixelAsVectorUInt8(const std::vector<uint32_t> & idx, const std::vector<uint8_t> & v) = 0;
199  virtual void
200  SetPixelAsVectorInt16(const std::vector<uint32_t> & idx, const std::vector<int16_t> & v) = 0;
201  virtual void
202  SetPixelAsVectorUInt16(const std::vector<uint32_t> & idx, const std::vector<uint16_t> & v) = 0;
203  virtual void
204  SetPixelAsVectorInt32(const std::vector<uint32_t> & idx, const std::vector<int32_t> & v) = 0;
205  virtual void
206  SetPixelAsVectorUInt32(const std::vector<uint32_t> & idx, const std::vector<uint32_t> & v) = 0;
207  virtual void
208  SetPixelAsVectorInt64(const std::vector<uint32_t> & idx, const std::vector<int64_t> & v) = 0;
209  virtual void
210  SetPixelAsVectorUInt64(const std::vector<uint32_t> & idx, const std::vector<uint64_t> & v) = 0;
211  virtual void
212  SetPixelAsVectorFloat32(const std::vector<uint32_t> & idx, const std::vector<float> & v) = 0;
213  virtual void
214  SetPixelAsVectorFloat64(const std::vector<uint32_t> & idx, const std::vector<double> & v) = 0;
215 
216  virtual void
217  SetPixelAsComplexFloat32(const std::vector<uint32_t> & idx, const std::complex<float> v) = 0;
218  virtual void
219  SetPixelAsComplexFloat64(const std::vector<uint32_t> & idx, const std::complex<double> v) = 0;
220 
221 
222  virtual int8_t *
223  GetBufferAsInt8() = 0;
224  virtual uint8_t *
225  GetBufferAsUInt8() = 0;
226  virtual int16_t *
227  GetBufferAsInt16() = 0;
228  virtual uint16_t *
229  GetBufferAsUInt16() = 0;
230  virtual int32_t *
231  GetBufferAsInt32() = 0;
232  virtual uint32_t *
233  GetBufferAsUInt32() = 0;
234  virtual int64_t *
235  GetBufferAsInt64() = 0;
236  virtual uint64_t *
237  GetBufferAsUInt64() = 0;
238  virtual float *
239  GetBufferAsFloat() = 0;
240  virtual double *
241  GetBufferAsDouble() = 0;
242  virtual void *
243  GetBufferAsVoid() = 0;
244 
245  virtual const int8_t *
246  GetBufferAsInt8() const = 0;
247  virtual const uint8_t *
248  GetBufferAsUInt8() const = 0;
249  virtual const int16_t *
250  GetBufferAsInt16() const = 0;
251  virtual const uint16_t *
252  GetBufferAsUInt16() const = 0;
253  virtual const int32_t *
254  GetBufferAsInt32() const = 0;
255  virtual const uint32_t *
256  GetBufferAsUInt32() const = 0;
257  virtual const int64_t *
258  GetBufferAsInt64() const = 0;
259  virtual const uint64_t *
260  GetBufferAsUInt64() const = 0;
261  virtual const float *
262  GetBufferAsFloat() const = 0;
263  virtual const double *
264  GetBufferAsDouble() const = 0;
265  virtual const void *
266  GetBufferAsVoid() const = 0;
267 };
268 
269 } // namespace itk::simple
270 
271 
272 #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:100
SITKCommon_HIDDEN
#define SITKCommon_HIDDEN
Definition: sitkCommon.h:44
sitkPixelIDTokens.h
itk::simple::PimpleImageBase::GetWidth
virtual unsigned int GetWidth() const
Definition: sitkPimpleImageBase.h:65
itk::simple::PimpleImageBase::GetHeight
virtual unsigned int GetHeight() const
Definition: sitkPimpleImageBase.h:70
itk::simple
Definition: sitkAdditionalProcedures.h:28
itk::simple::PimpleImageBase::GetDepth
virtual unsigned int GetDepth() const
Definition: sitkPimpleImageBase.h:75
itk::DataObject