ITK  6.0.0
Insight Toolkit
itkSpecialCoordinatesImage.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 itkSpecialCoordinatesImage_h
19#define itkSpecialCoordinatesImage_h
20
21#include "itkImageBase.h"
25#include "itkContinuousIndex.h"
26
27namespace itk
28{
94template <typename TPixel, unsigned int VImageDimension = 2>
95class ITK_TEMPLATE_EXPORT SpecialCoordinatesImage : public ImageBase<VImageDimension>
96{
97public:
98 ITK_DISALLOW_COPY_AND_MOVE(SpecialCoordinatesImage);
99
106
108 itkNewMacro(Self);
109
111 itkOverrideGetNameOfClassMacro(SpecialCoordinatesImage);
112
115 using PixelType = TPixel;
116
118 using ValueType = TPixel;
119
124 using InternalPixelType = TPixel;
125
127
131
136
141 static constexpr unsigned int ImageDimension = VImageDimension;
142
144 using typename Superclass::IndexType;
145
147 using typename Superclass::OffsetType;
148
150 using typename Superclass::SizeType;
151
154
157 using typename Superclass::RegionType;
158
163 using typename Superclass::SpacingType;
164
167 using typename Superclass::PointType;
168
172
175 void
176 Allocate(bool initialize = false) override;
177
180 void
181 Initialize() override;
182
185 void
186 FillBuffer(const TPixel & value);
187
193 void
194 SetPixel(const IndexType & index, const TPixel & value)
195 {
196 const OffsetValueType offset = this->FastComputeOffset(index);
197 (*m_Buffer)[offset] = value;
198 }
199
204 const TPixel &
205 GetPixel(const IndexType & index) const
206 {
207 const OffsetValueType offset = this->FastComputeOffset(index);
208 return ((*m_Buffer)[offset]);
209 }
210
215 TPixel &
216 GetPixel(const IndexType & index)
217 {
218 const OffsetValueType offset = this->FastComputeOffset(index);
219 return ((*m_Buffer)[offset]);
220 }
221
226 TPixel &
227 operator[](const IndexType & index)
228 {
229 return this->GetPixel(index);
230 }
231
236 const TPixel &
237 operator[](const IndexType & index) const
238 {
239 return this->GetPixel(index);
240 }
241
244 TPixel *
246 {
247 return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
248 }
249 const TPixel *
251 {
252 return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
253 }
257 PixelContainer *
259 {
260 return m_Buffer.GetPointer();
261 }
262
263 const PixelContainer *
265 {
266 return m_Buffer.GetPointer();
267 }
268
271 void
273
277 {
278 return AccessorType();
279 }
280
282 const AccessorType
284 {
285 return AccessorType();
286 }
287
293 void
294 SetSpacing(const SpacingType &) override
295 {}
296 void
297 SetSpacing(const double[VImageDimension]) override
298 {}
299 void
300 SetSpacing(const float[VImageDimension]) override
301 {}
302 void
303 SetOrigin(const PointType) override
304 {}
305 void
306 SetOrigin(const double[VImageDimension]) override
307 {}
308 void
309 SetOrigin(const float[VImageDimension]) override
310 {}
313 /* It is ILLEGAL in C++ to make a templated member function virtual! */
314 /* Therefore, we must just let templates take care of everything. */
315 /*
316 template<typename TCoordinate>
317 virtual bool TransformPhysicalPointToContinuousIndex(
318 const Point<TCoordinate, VImageDimension>& point,
319 ContinuousIndex<TCoordinate, VImageDimension>& index ) const = 0;
320
321 template<typename TCoordinate>
322 virtual bool TransformPhysicalPointToIndex(
323 const Point<TCoordinate, VImageDimension>&,
324 IndexType & index ) const = 0;
325
326 template<typename TCoordinate>
327 virtual void TransformContinuousIndexToPhysicalPoint(
328 const ContinuousIndex<TCoordinate, VImageDimension>& index,
329 Point<TCoordinate, VImageDimension>& point ) const = 0;
330
331 template<typename TCoordinate>
332 virtual void TransformIndexToPhysicalPoint(
333 const IndexType & index,
334 Point<TCoordinate, VImageDimension>& point ) const = 0;
335 */
336
337protected:
339 void
340 PrintSelf(std::ostream & os, Indent indent) const override;
341
342 ~SpecialCoordinatesImage() override = default;
343
344private:
347};
348} // end namespace itk
349
350#ifndef ITK_MANUAL_INSTANTIATION
351# include "itkSpecialCoordinatesImage.hxx"
352#endif
353
354#endif
Base class for all data objects in ITK.
Provides a common API for pixel accessors for Image and VectorImage.
Give access to partial aspects a type.
Base class for templated image classes.
Definition: itkImageBase.h:115
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:147
Defines an itk::Image front-end to a standard C-array.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for most ITK classes.
Definition: itkObject.h:62
Templated n-dimensional nonrectilinear-coordinate image base class.
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
void SetOrigin(const double[VImageDimension]) override
void SetSpacing(const double[VImageDimension]) override
~SpecialCoordinatesImage() override=default
typename PixelContainer::Pointer PixelContainerPointer
void SetSpacing(const SpacingType &) override
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
void SetOrigin(const float[VImageDimension]) override
typename PixelContainer::ConstPointer PixelContainerConstPointer
void PrintSelf(std::ostream &os, Indent indent) const override
void SetSpacing(const float[VImageDimension]) override
void SetPixelContainer(PixelContainer *container)
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
void SetOrigin(const PointType) override
void FillBuffer(const TPixel &value)
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
void Allocate(bool initialize=false) override
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
const PixelContainer * GetPixelContainer() const
const AccessorType GetPixelAccessor() const
Implements a weak reference to an object.
SmartPointer< const Self > ConstPointer
static Pointer New()
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....