ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkSparseImage.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 __itkSparseImage_h
19 #define __itkSparseImage_h
20 
21 #include "itkImage.h"
22 #include "itkSparseFieldLayer.h"
23 #include "itkObjectStore.h"
24 
25 namespace itk
26 {
66 template< typename TNode, unsigned int VImageDimension = 2 >
67 class SparseImage:public Image< TNode *, VImageDimension >
68 {
69 public:
71  typedef SparseImage Self;
76 
78  itkNewMacro(Self);
79 
81  itkTypeMacro(SparseImage, Image);
82 
84  itkStaticConstMacro(ImageDimension, unsigned int,
86 
88  typedef TNode NodeType;
89 
91  typedef typename Superclass::IndexType IndexType;
92 
97 
99 
103 
107  { return NeighborhoodAccessorFunctorType(); }
108 
112  { return NeighborhoodAccessorFunctorType(); }
113 
116  NodeType * AddNode(const IndexType & index)
117  {
118  m_NodeList->PushFront( m_NodeStore->Borrow() );
119  NodeType *node = m_NodeList->Front();
120  node->m_Index = index;
121  this->SetPixel(index, node);
122  return node;
123  }
125 
129  {
130  return m_NodeList;
131  }
132 
135  virtual void Initialize();
136 
137 protected:
138  SparseImage();
140 
141  void PrintSelf(std::ostream & os, Indent indent) const;
142 
143 private:
146 
148 
149  SparseImage(const Self &); //purposely not implemented
150  void operator=(const Self &); //purposely not implemented
151 };
152 } // end namespace itk
153 
154 #ifndef ITK_MANUAL_INSTANTIATION
155 #include "itkSparseImage.hxx"
156 #endif
157 
158 #endif
NodeStoreType::Pointer m_NodeStore
SparseImage Self
SparseFieldLayer< NodeType > NodeListType
WeakPointer< const Self > ConstWeakPointer
void PrintSelf(std::ostream &os, Indent indent) const
Image< TNode *, VImageDimension > Superclass
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
NodeListType::Pointer m_NodeList
SmartPointer< Self > Pointer
Implements a weak reference to an object.
A storage type for sparse image data.
Superclass::IndexType IndexType
virtual void Initialize()
ObjectStore< NodeType > NodeStoreType
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
SmartPointer< const Self > ConstPointer
Superclass::IndexType IndexType
Definition: itkImage.h:122
Superclass::IOPixelType IOPixelType
NodeType * AddNode(const IndexType &index)
NeighborhoodAccessorFunctor< Self > NeighborhoodAccessorFunctorType
static const unsigned int ImageDimension
Definition: itkImage.h:119
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void operator=(const Self &)
A very simple linked list that is used to manage nodes in a layer of a sparse field level-set solver...
static const unsigned int ImageDimension
A specialized memory management object for allocating and destroying contiguous blocks of objects...
void SetPixel(const IndexType &index, const TNode *&value)
Set a pixel value.
Definition: itkImage.h:187
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
NodeListType * GetNodeList()