ITK  5.0.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 ITK_TEMPLATE_EXPORT SparseImage:public Image< TNode *, VImageDimension >
68 {
69 public:
70  ITK_DISALLOW_COPY_AND_ASSIGN(SparseImage);
71 
73  using Self = SparseImage;
78 
80  itkNewMacro(Self);
81 
83  itkTypeMacro(SparseImage, Image);
84 
86  static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
87 
89  using NodeType = TNode;
90 
92  using IndexType = typename Superclass::IndexType;
93 
97 
98  using IOPixelType = typename Superclass::IOPixelType;
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  void Initialize() override;
136 
137 protected:
138  SparseImage();
139  ~SparseImage() override = default;
140 
141  void PrintSelf(std::ostream & os, Indent indent) const override;
142 
143 private:
146 
148 
149 };
150 } // end namespace itk
151 
152 #ifndef ITK_MANUAL_INSTANTIATION
153 #include "itkSparseImage.hxx"
154 #endif
155 
156 #endif
NodeStoreType::Pointer m_NodeStore
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
NodeListType::Pointer m_NodeList
Implements a weak reference to an object.
A storage type for sparse image data.
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
NodeType * AddNode(const IndexType &index)
typename Superclass::IOPixelType IOPixelType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Base class for most ITK classes.
Definition: itkObject.h:60
A very simple linked list that is used to manage nodes in a layer of a sparse field level-set solver...
A specialized memory management object for allocating and destroying contiguous blocks of objects...
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
NodeListType * GetNodeList()