ITK  5.4.0
Insight Toolkit
itkOctree.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 itkOctree_h
19 #define itkOctree_h
20 
21 #include "itkOctreeNode.h"
22 #include "itkImage.h"
23 #include "itkCommonEnums.h"
24 
25 namespace itk
26 {
27 enum
28 {
32 };
33 
40 class ITKCommon_EXPORT OctreeBase : public Object
41 {
42 public:
43 
45  using Self = OctreeBase;
47 
49 
50 #if !defined(ITK_LEGACY_REMOVE)
51  // We need to expose the enum values at the class level
52  // for backwards compatibility
53  static constexpr OctreeEnum UNKNOWN_PLANE = OctreeEnum::UNKNOWN_PLANE;
54  static constexpr OctreeEnum SAGITAL_PLANE = OctreeEnum::SAGITAL_PLANE;
55  static constexpr OctreeEnum CORONAL_PLANE = OctreeEnum::CORONAL_PLANE;
56  static constexpr OctreeEnum TRANSVERSE_PLANE = OctreeEnum::TRANSVERSE_PLANE;
57 #endif
58 
63  virtual OctreeNode *
64  GetTree() = 0;
65 
70  virtual unsigned int
71  GetDepth() = 0;
72 
78  virtual unsigned int
79  GetWidth() = 0;
80 
82  virtual void
83  SetDepth(unsigned int depth) = 0;
84 
86  virtual void
87  SetWidth(unsigned int width) = 0;
88 
94  virtual void
95  BuildFromBuffer(const void * buffer,
96  const unsigned int xsize,
97  const unsigned int ysize,
98  const unsigned int zsize) = 0;
99 
110  virtual const OctreeNodeBranch *
111  GetColorTable() const = 0;
112 
114  virtual int
115  GetColorTableSize() const = 0;
116 };
117 
126 template <typename TPixel, unsigned int ColorTableSize, typename MappingFunctionType>
127 class ITK_TEMPLATE_EXPORT Octree : public OctreeBase
128 {
129 public:
130  ITK_DISALLOW_COPY_AND_MOVE(Octree);
131 
133  using Self = Octree;
138 
140  itkNewMacro(Self);
141 
143  itkOverrideGetNameOfClassMacro(Octree);
144 
146  GetImage();
147 
148  void
149  BuildFromBuffer(const void * frombuffer,
150  const unsigned int xsize,
151  const unsigned int ysize,
152  const unsigned int zsize) override;
153 
154  void BuildFromImage(Image<TPixel, 3> * fromImage);
155 
156  Octree();
157  ~Octree() override;
158  void
159  SetColor(unsigned int color)
160  {
161  m_Tree.SetColor(color);
162  }
163  void
165  {
166  m_Tree.SetBranch(branch);
167  }
168  void
169  SetTrueDims(const unsigned int Dim0, const unsigned int Dim1, const unsigned int Dim2);
170 
171  int
172  GetValue(const unsigned int Dim0, const unsigned int Dim1, const unsigned int Dim2);
173 
174  void
175  SetWidth(unsigned int width) override;
176 
177  void
178  SetDepth(unsigned int depth) override;
179 
180  unsigned int
181  GetWidth() override;
182 
183  unsigned int
184  GetDepth() override;
185 
186  /***
187  * Exposes enum values for backwards compatibility
188  * */
189 #if !defined(ITK_LEGACY_REMOVE)
190  static constexpr OctreeEnum UNKNOWN_PLANE = OctreeEnum::UNKNOWN_PLANE;
191  static constexpr OctreeEnum SAGITAL_PLANE = OctreeEnum::SAGITAL_PLANE;
192  static constexpr OctreeEnum CORONAL_PLANE = OctreeEnum::CORONAL_PLANE;
193  static constexpr OctreeEnum TRANSVERSE_PLANE = OctreeEnum::TRANSVERSE_PLANE;
194 #endif
195 
196  OctreeNode *
197  GetTree() override;
198 
199  const OctreeNodeBranch *
200  GetColorTable() const override;
201 
202  int
203  GetColorTableSize() const override;
204 
205 private:
207  maskToOctree(const TPixel * Mask,
208  unsigned int width,
209  unsigned int x,
210  unsigned int y,
211  unsigned int z,
212  unsigned int xsize,
213  unsigned int ysize,
214  unsigned int zsize);
215 
216  OctreeEnum m_Plane{ OctreeEnum::UNKNOWN_PLANE }; // The orientation of the plane for this octree
217 
218  // The width of the Octree. This is always a power of 2,
219  // and large enough to contain MAX(DIMS[1,2,3])
220  unsigned int m_Width{ 0 };
221 
222  unsigned int m_Depth{ 0 }; // The depth of the Octree
223  unsigned int m_TrueDims[3]{}; // The true dimensions of the image
224  OctreeNodeBranch m_ColorTable[ColorTableSize];
225  OctreeNode m_Tree{};
226  MappingFunctionType m_MappingFunction{};
227 };
228 } // namespace itk
229 
230 #ifndef ITK_MANUAL_INSTANTIATION
231 # include "itkOctree.hxx"
232 #endif
233 
234 #endif // itkOctree_h
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkOctreeNode.h
itk::OctreeEnums::Octree
Octree
Definition: itkCommonEnums.h:220
itk::B2_MASKFILE_GRAY
Definition: itkOctree.h:31
itk::OctreeBase
class ITK_FORWARD_EXPORT OctreeBase
Definition: itkOctreeNode.h:28
itk::Octree
Represent a 3D Image with an Octree data structure.
Definition: itkOctree.h:127
itkImage.h
itk::Octree::SetColor
void SetColor(unsigned int color)
Definition: itkOctree.h:159
itk::SmartPointer< Self >
itk::OctreeBase
Provides non-templated access to templated instances of Octree.
Definition: itkOctree.h:40
itk::B2_MASKFILE_BLACK
Definition: itkOctree.h:29
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::OctreeNodeBranch
Definition: itkOctreeNode.h:146
itk::B2_MASKFILE_WHITE
Definition: itkOctree.h:30
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itkCommonEnums.h
itk::OctreeNode
A data structure representing a node in an Octree.
Definition: itkOctreeNode.h:46
itk::Octree::ImageTypePointer
typename ImageType::Pointer ImageTypePointer
Definition: itkOctree.h:137
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::Octree::SetTree
void SetTree(OctreeNodeBranch *branch)
Definition: itkOctree.h:164