ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkSpecialCoordinatesImage.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 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_ASSIGN(SpecialCoordinatesImage);
99 
106 
108  itkNewMacro(Self);
109 
111  itkTypeMacro(SpecialCoordinatesImage, ImageBase);
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 
145 
147  using OffsetType = typename Superclass::OffsetType;
148 
150  using SizeType = typename Superclass::SizeType;
151 
154 
158 
163  using SpacingType = typename Superclass::SpacingType;
164 
168 
172 
175  void Allocate(bool initialize=false) override;
176 
179  void Initialize() override;
180 
183  void FillBuffer(const TPixel & value);
184 
190  void SetPixel(const IndexType & index, const TPixel & value)
191  {
192  OffsetValueType offset = this->FastComputeOffset(index);
193  ( *m_Buffer )[offset] = value;
194  }
195 
200  const TPixel & GetPixel(const IndexType & index) const
201  {
202  OffsetValueType offset = this->FastComputeOffset(index);
203  return ( ( *m_Buffer )[offset] );
204  }
205 
210  TPixel & GetPixel(const IndexType & index)
211  {
212  OffsetValueType offset = this->FastComputeOffset(index);
213  return ( ( *m_Buffer )[offset] );
214  }
215 
220  TPixel & operator[](const IndexType & index) { return this->GetPixel(index); }
221 
226  const TPixel & operator[](const IndexType & index) const { return this->GetPixel(index); }
227 
230  TPixel * GetBufferPointer() { return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr; }
231  const TPixel * GetBufferPointer() const { return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr; }
233 
235  PixelContainer * GetPixelContainer() { return m_Buffer.GetPointer(); }
236 
237  const PixelContainer * GetPixelContainer() const { return m_Buffer.GetPointer(); }
238 
241  void SetPixelContainer(PixelContainer *container);
242 
245 
247  const AccessorType GetPixelAccessor() const { return AccessorType(); }
248 
254  void SetSpacing(const SpacingType &) override {}
255  void SetSpacing(const double[VImageDimension]) override {}
256  void SetSpacing(const float[VImageDimension]) override {}
257  void SetOrigin(const PointType) override {}
258  void SetOrigin(const double[VImageDimension]) override {}
259  void SetOrigin(const float[VImageDimension]) override {}
261 
262  /* It is ILLEGAL in C++ to make a templated member function virtual! */
263  /* Therefore, we must just let templates take care of everything. */
264  /*
265  template<typename TCoordRep>
266  virtual bool TransformPhysicalPointToContinuousIndex(
267  const Point<TCoordRep, VImageDimension>& point,
268  ContinuousIndex<TCoordRep, VImageDimension>& index ) const = 0;
269 
270  template<typename TCoordRep>
271  virtual bool TransformPhysicalPointToIndex(
272  const Point<TCoordRep, VImageDimension>&,
273  IndexType & index ) const = 0;
274 
275  template<typename TCoordRep>
276  virtual void TransformContinuousIndexToPhysicalPoint(
277  const ContinuousIndex<TCoordRep, VImageDimension>& index,
278  Point<TCoordRep, VImageDimension>& point ) const = 0;
279 
280  template<typename TCoordRep>
281  virtual void TransformIndexToPhysicalPoint(
282  const IndexType & index,
283  Point<TCoordRep, VImageDimension>& point ) const = 0;
284  */
285 
286 protected:
288  void PrintSelf(std::ostream & os, Indent indent) const override;
289 
290  ~SpecialCoordinatesImage() override = default;
291 
292 private:
295 };
296 } // end namespace itk
297 
298 #ifndef ITK_MANUAL_INSTANTIATION
299 #include "itkSpecialCoordinatesImage.hxx"
300 #endif
301 
302 #endif
void SetSpacing(const float[VImageDimension]) override
void SetOrigin(const PointType) override
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
An image region represents a structured region of data.
Templated n-dimensional nonrectilinear-coordinate image base class.
Implements a weak reference to an object.
const AccessorType GetPixelAccessor() const
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image...
Definition: itkOffset.h:67
const TPixel * GetBufferPointer() const
void SetSpacing(const double[VImageDimension]) override
Provides a common API for pixel accessors for Image and VectorImage.
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
Base class for templated image classes.
Definition: itkImageBase.h:105
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Control indentation during Print() invocation.
Definition: itkIndent.h:49
const PixelContainer * GetPixelContainer() const
void SetSpacing(const SpacingType &) override
Give access to partial aspects a type.
Base class for most ITK classes.
Definition: itkObject.h:60
typename PixelContainer::ConstPointer PixelContainerConstPointer
void SetOrigin(const float[VImageDimension]) override
void SetOrigin(const double[VImageDimension]) override
Base class for all data objects in ITK.
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
Defines an itk::Image front-end to a standard C-array.
typename PixelContainer::Pointer PixelContainerPointer