ITK  5.4.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 
27 namespace itk
28 {
94 template <typename TPixel, unsigned int VImageDimension = 2>
95 class ITK_TEMPLATE_EXPORT SpecialCoordinatesImage : public ImageBase<VImageDimension>
96 {
97 public:
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  OffsetValueType offset = this->FastComputeOffset(index);
197  (*m_Buffer)[offset] = value;
198  }
199 
204  const TPixel &
205  GetPixel(const IndexType & index) const
206  {
207  OffsetValueType offset = this->FastComputeOffset(index);
208  return ((*m_Buffer)[offset]);
209  }
210 
215  TPixel &
216  GetPixel(const IndexType & index)
217  {
218  OffsetValueType offset = this->FastComputeOffset(index);
219  return ((*m_Buffer)[offset]);
220  }
221 
226  TPixel & operator[](const IndexType & index) { return this->GetPixel(index); }
227 
232  const TPixel & operator[](const IndexType & index) const { return this->GetPixel(index); }
233 
236  TPixel *
238  {
239  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
240  }
241  const TPixel *
243  {
244  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
245  }
249  PixelContainer *
251  {
252  return m_Buffer.GetPointer();
253  }
254 
255  const PixelContainer *
257  {
258  return m_Buffer.GetPointer();
259  }
260 
263  void
264  SetPixelContainer(PixelContainer * container);
265 
267  AccessorType
269  {
270  return AccessorType();
271  }
272 
274  const AccessorType
276  {
277  return AccessorType();
278  }
279 
285  void
286  SetSpacing(const SpacingType &) override
287  {}
288  void
289  SetSpacing(const double[VImageDimension]) override
290  {}
291  void
292  SetSpacing(const float[VImageDimension]) override
293  {}
294  void
295  SetOrigin(const PointType) override
296  {}
297  void
298  SetOrigin(const double[VImageDimension]) override
299  {}
300  void
301  SetOrigin(const float[VImageDimension]) override
302  {}
305  /* It is ILLEGAL in C++ to make a templated member function virtual! */
306  /* Therefore, we must just let templates take care of everything. */
307  /*
308  template<typename TCoordRep>
309  virtual bool TransformPhysicalPointToContinuousIndex(
310  const Point<TCoordRep, VImageDimension>& point,
311  ContinuousIndex<TCoordRep, VImageDimension>& index ) const = 0;
312 
313  template<typename TCoordRep>
314  virtual bool TransformPhysicalPointToIndex(
315  const Point<TCoordRep, VImageDimension>&,
316  IndexType & index ) const = 0;
317 
318  template<typename TCoordRep>
319  virtual void TransformContinuousIndexToPhysicalPoint(
320  const ContinuousIndex<TCoordRep, VImageDimension>& index,
321  Point<TCoordRep, VImageDimension>& point ) const = 0;
322 
323  template<typename TCoordRep>
324  virtual void TransformIndexToPhysicalPoint(
325  const IndexType & index,
326  Point<TCoordRep, VImageDimension>& point ) const = 0;
327  */
328 
329 protected:
330  SpecialCoordinatesImage() = default;
331  void
332  PrintSelf(std::ostream & os, Indent indent) const override;
333 
334  ~SpecialCoordinatesImage() override = default;
335 
336 private:
339 };
340 } // end namespace itk
341 
342 #ifndef ITK_MANUAL_INSTANTIATION
343 # include "itkSpecialCoordinatesImage.hxx"
344 #endif
345 
346 #endif
itk::SpecialCoordinatesImage::operator[]
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
Definition: itkSpecialCoordinatesImage.h:226
itk::SpecialCoordinatesImage< TPixel, 3 >::PixelContainerPointer
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkSpecialCoordinatesImage.h:170
itk::ImageBase::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:147
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::SpecialCoordinatesImage::GetPixelAccessor
AccessorType GetPixelAccessor()
Definition: itkSpecialCoordinatesImage.h:268
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::Index< VImageDimension >
itkContinuousIndex.h
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:114
itk::SpecialCoordinatesImage::GetPixelContainer
PixelContainer * GetPixelContainer()
Definition: itkSpecialCoordinatesImage.h:250
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::Vector< SpacingValueType, VImageDimension >
itk::SpecialCoordinatesImage::SetSpacing
void SetSpacing(const double[VImageDimension]) override
Definition: itkSpecialCoordinatesImage.h:289
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SpecialCoordinatesImage< TPixel, 3 >::IOPixelType
PixelType IOPixelType
Definition: itkSpecialCoordinatesImage.h:126
itk::SpecialCoordinatesImage::GetPixelAccessor
const AccessorType GetPixelAccessor() const
Definition: itkSpecialCoordinatesImage.h:275
itk::SmartPointer< Self >
itk::DefaultPixelAccessor
Give access to partial aspects a type.
Definition: itkDefaultPixelAccessor.h:54
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::SpecialCoordinatesImage
Templated n-dimensional nonrectilinear-coordinate image base class.
Definition: itkSpecialCoordinatesImage.h:95
itk::SpecialCoordinatesImage< TPixel, 3 >::ValueType
TPixel ValueType
Definition: itkSpecialCoordinatesImage.h:118
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itkDefaultPixelAccessor.h
itk::DefaultPixelAccessorFunctor
Provides a common API for pixel accessors for Image and VectorImage.
Definition: itkDefaultPixelAccessorFunctor.h:47
itk::SpecialCoordinatesImage::SetSpacing
void SetSpacing(const SpacingType &) override
Definition: itkSpecialCoordinatesImage.h:286
itk::SpecialCoordinatesImage::operator[]
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Definition: itkSpecialCoordinatesImage.h:232
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::SpecialCoordinatesImage::GetBufferPointer
const TPixel * GetBufferPointer() const
Definition: itkSpecialCoordinatesImage.h:242
itk::SpecialCoordinatesImage< TPixel, 3 >::PixelContainerConstPointer
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition: itkSpecialCoordinatesImage.h:171
itkImportImageContainer.h
itk::SpecialCoordinatesImage::GetPixel
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
Definition: itkSpecialCoordinatesImage.h:216
itk::SpecialCoordinatesImage::GetPixelContainer
const PixelContainer * GetPixelContainer() const
Definition: itkSpecialCoordinatesImage.h:256
itk::ImportImageContainer
Defines an itk::Image front-end to a standard C-array.
Definition: itkImportImageContainer.h:45
itk::SpecialCoordinatesImage::GetBufferPointer
TPixel * GetBufferPointer()
Definition: itkSpecialCoordinatesImage.h:237
itk::WeakPointer
Implements a weak reference to an object.
Definition: itkWeakPointer.h:44
itk::SpecialCoordinatesImage::SetOrigin
void SetOrigin(const float[VImageDimension]) override
Definition: itkSpecialCoordinatesImage.h:301
itk::SpecialCoordinatesImage::SetPixel
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
Definition: itkSpecialCoordinatesImage.h:194
itk::SpecialCoordinatesImage::SetSpacing
void SetSpacing(const float[VImageDimension]) override
Definition: itkSpecialCoordinatesImage.h:292
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::SpecialCoordinatesImage::SetOrigin
void SetOrigin(const PointType) override
Definition: itkSpecialCoordinatesImage.h:295
itkDefaultPixelAccessorFunctor.h
itk::SpecialCoordinatesImage< TPixel, 3 >::PixelType
TPixel PixelType
Definition: itkSpecialCoordinatesImage.h:115
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::SpecialCoordinatesImage::SetOrigin
void SetOrigin(const double[VImageDimension]) override
Definition: itkSpecialCoordinatesImage.h:298
itk::Point< PointValueType, VImageDimension >
itk::SpecialCoordinatesImage::GetPixel
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Definition: itkSpecialCoordinatesImage.h:205
New
static Pointer New()
itk::SpecialCoordinatesImage< TPixel, 3 >::InternalPixelType
TPixel InternalPixelType
Definition: itkSpecialCoordinatesImage.h:124
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itkImageBase.h