ITK  4.9.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 SpecialCoordinatesImage:public ImageBase< VImageDimension >
96 {
97 public:
104 
106  itkNewMacro(Self);
107 
110 
113  typedef TPixel PixelType;
114 
116  typedef TPixel ValueType;
117 
122  typedef TPixel InternalPixelType;
123 
125 
129 
134 
139  itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
140 
143 
146 
148  typedef typename Superclass::SizeType SizeType;
149 
152 
156 
162 
166 
170 
173  virtual void Allocate(bool initialize=false) ITK_OVERRIDE;
174 
177  virtual void Initialize() ITK_OVERRIDE;
178 
181  void FillBuffer(const TPixel & value);
182 
188  void SetPixel(const IndexType & index, const TPixel & value)
189  {
190  OffsetValueType offset = this->FastComputeOffset(index);
191  ( *m_Buffer )[offset] = value;
192  }
193 
198  const TPixel & GetPixel(const IndexType & index) const
199  {
200  OffsetValueType offset = this->FastComputeOffset(index);
201  return ( ( *m_Buffer )[offset] );
202  }
203 
208  TPixel & GetPixel(const IndexType & index)
209  {
210  OffsetValueType offset = this->FastComputeOffset(index);
211  return ( ( *m_Buffer )[offset] );
212  }
213 
218  TPixel & operator[](const IndexType & index) { return this->GetPixel(index); }
219 
224  const TPixel & operator[](const IndexType & index) const { return this->GetPixel(index); }
225 
228  TPixel * GetBufferPointer() { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
229  const TPixel * GetBufferPointer() const { return m_Buffer ? m_Buffer->GetBufferPointer() : ITK_NULLPTR; }
231 
234 
235  const PixelContainer * GetPixelContainer() const { return m_Buffer.GetPointer(); }
236 
239  void SetPixelContainer(PixelContainer *container);
240 
243 
245  const AccessorType GetPixelAccessor(void) const { return AccessorType(); }
246 
252  virtual void SetSpacing(const SpacingType &) ITK_OVERRIDE {}
253  virtual void SetSpacing(const double[VImageDimension]) ITK_OVERRIDE {}
254  virtual void SetSpacing(const float[VImageDimension]) ITK_OVERRIDE {}
255  virtual void SetOrigin(const PointType) ITK_OVERRIDE {}
256  virtual void SetOrigin(const double[VImageDimension]) ITK_OVERRIDE {}
257  virtual void SetOrigin(const float[VImageDimension]) ITK_OVERRIDE {}
259 
260  /* It is ILLEGAL in C++ to make a templated member function virtual! */
261  /* Therefore, we must just let templates take care of everything. */
262  /*
263  template<typename TCoordRep>
264  virtual bool TransformPhysicalPointToContinuousIndex(
265  const Point<TCoordRep, VImageDimension>& point,
266  ContinuousIndex<TCoordRep, VImageDimension>& index ) const = 0;
267 
268  template<typename TCoordRep>
269  virtual bool TransformPhysicalPointToIndex(
270  const Point<TCoordRep, VImageDimension>&,
271  IndexType & index ) const = 0;
272 
273  template<typename TCoordRep>
274  virtual void TransformContinuousIndexToPhysicalPoint(
275  const ContinuousIndex<TCoordRep, VImageDimension>& index,
276  Point<TCoordRep, VImageDimension>& point ) const = 0;
277 
278  template<typename TCoordRep>
279  virtual void TransformIndexToPhysicalPoint(
280  const IndexType & index,
281  Point<TCoordRep, VImageDimension>& point ) const = 0;
282  */
283 
284 protected:
286  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
287 
289 
290 private:
291  SpecialCoordinatesImage(const Self &) ITK_DELETE_FUNCTION;
292  void operator=(const Self &) ITK_DELETE_FUNCTION;
293 
296 };
297 } // end namespace itk
298 
299 #ifndef ITK_MANUAL_INSTANTIATION
300 #include "itkSpecialCoordinatesImage.hxx"
301 #endif
302 
303 #endif
virtual void Initialize() override
virtual void Allocate(bool initialize=false) override
Index< VImageDimension > IndexType
Definition: itkImageBase.h:137
DefaultPixelAccessor< PixelType > AccessorType
void PrintSelf(std::ostream &os, Indent indent) const override
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
Templated n-dimensional nonrectilinear-coordinate image base class.
static const unsigned int ImageDimension
ObjectType * GetPointer() const
DefaultPixelAccessorFunctor< Self > AccessorFunctorType
Implements a weak reference to an object.
Size< VImageDimension > SizeType
Definition: itkImageBase.h:146
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:162
ImageBase< VImageDimension > Superclass
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
virtual void SetSpacing(const SpacingType &) override
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:142
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:157
const TPixel * GetBufferPointer() const
void FillBuffer(const TPixel &value)
OffsetValueType FastComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:731
Provides a common API for pixel accessors for Image and VectorImage.
PixelContainer::ConstPointer PixelContainerConstPointer
virtual void SetSpacing(const float[VImageDimension]) override
OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:143
PixelContainer::Pointer PixelContainerPointer
virtual void SetOrigin(const double[VImageDimension]) override
Base class for templated image classes.
Definition: itkImageBase.h:112
virtual void SetOrigin(const float[VImageDimension]) override
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
SmartPointer< const Self > ConstPointer
void SetPixelContainer(PixelContainer *container)
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:150
Give access to partial aspects a type.
WeakPointer< const Self > ConstWeakPointer
virtual void SetSpacing(const double[VImageDimension]) override
virtual void SetOrigin(const PointType) override
Base class for all data objects in ITK.
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
const AccessorType GetPixelAccessor(void) const
Defines an itk::Image front-end to a standard C-array.
ImportImageContainer< SizeValueType, PixelType > PixelContainer