ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
25namespace itk
26{
27enum
28{
32};
33
40class ITKCommon_EXPORT OctreeBase : public Object
41{
42public:
46
48
49#if !defined(ITK_LEGACY_REMOVE)
50 // We need to expose the enum values at the class level
51 // for backwards compatibility
52 static constexpr OctreeEnum UNKNOWN_PLANE = OctreeEnum::UNKNOWN_PLANE;
53 static constexpr OctreeEnum SAGITTAL_PLANE = OctreeEnum::SAGITTAL_PLANE;
54 [[deprecated("Use SAGITTAL_PLANE instead")]] static constexpr OctreeEnum SAGITAL_PLANE = OctreeEnum::SAGITTAL_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
126template <typename TPixel, unsigned int ColorTableSize, typename MappingFunctionType>
127class ITK_TEMPLATE_EXPORT Octree : public OctreeBase
128{
129public:
130 ITK_DISALLOW_COPY_AND_MOVE(Octree);
131
133 using Self = Octree;
138
140 itkNewMacro(Self);
141
143 itkOverrideGetNameOfClassMacro(Octree);
144
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
156
158 ~Octree() override;
159 void
160 SetColor(unsigned int color)
161 {
162 m_Tree.SetColor(color);
163 }
164 void
166 {
167 m_Tree.SetBranch(branch);
168 }
169 void
170 SetTrueDims(const unsigned int Dim0, const unsigned int Dim1, const unsigned int Dim2);
171
172 int
173 GetValue(const unsigned int Dim0, const unsigned int Dim1, const unsigned int Dim2);
174
175 void
176 SetWidth(unsigned int width) override;
177
178 void
179 SetDepth(unsigned int depth) override;
180
181 unsigned int
182 GetWidth() override;
183
184 unsigned int
185 GetDepth() override;
186
187#if !defined(ITK_LEGACY_REMOVE)
188 /***
189 * Exposes enum values for backwards compatibility
190 * */
192 static constexpr OctreeEnum UNKNOWN_PLANE = OctreeEnum::UNKNOWN_PLANE;
193 static constexpr OctreeEnum SAGITTAL_PLANE = OctreeEnum::SAGITTAL_PLANE;
194 [[deprecated("Use SAGITTAL_PLANE instead")]] static constexpr OctreeEnum SAGITAL_PLANE = OctreeEnum::SAGITTAL_PLANE;
195 static constexpr OctreeEnum CORONAL_PLANE = OctreeEnum::CORONAL_PLANE;
196 static constexpr OctreeEnum TRANSVERSE_PLANE = OctreeEnum::TRANSVERSE_PLANE;
198#endif
199
200 OctreeNode *
201 GetTree() override;
202
203 const OctreeNodeBranch *
204 GetColorTable() const override;
205
206 int
207 GetColorTableSize() const override;
208
209private:
211 maskToOctree(const TPixel * Mask,
212 unsigned int width,
213 unsigned int x,
214 unsigned int y,
215 unsigned int z,
216 unsigned int xsize,
217 unsigned int ysize,
218 unsigned int zsize);
219
220 OctreeEnum m_Plane{ OctreeEnum::UNKNOWN_PLANE }; // The orientation of the plane for this octree
221
222 // The width of the Octree. This is always a power of 2,
223 // and large enough to contain MAX(DIMS[1,2,3])
224 unsigned int m_Width{ 0 };
225
226 unsigned int m_Depth{ 0 }; // The depth of the Octree
227 unsigned int m_TrueDims[3]{}; // The true dimensions of the image
230 MappingFunctionType m_MappingFunction{};
231};
232} // namespace itk
233
234#ifndef ITK_MANUAL_INSTANTIATION
235# include "itkOctree.hxx"
236#endif
237
238#endif // itkOctree_h
Templated n-dimensional image class.
Definition itkImage.h:89
SmartPointer< Self > Pointer
Definition itkImage.h:96
Provides non-templated access to templated instances of Octree.
Definition itkOctree.h:41
virtual void SetWidth(unsigned int width)=0
SmartPointer< Self > Pointer
Definition itkOctree.h:45
virtual int GetColorTableSize() const =0
virtual OctreeNode * GetTree()=0
virtual unsigned int GetWidth()=0
virtual const OctreeNodeBranch * GetColorTable() const =0
virtual void SetDepth(unsigned int depth)=0
OctreeEnums::Octree OctreeEnum
Definition itkOctree.h:47
OctreeBase Self
Definition itkOctree.h:44
virtual unsigned int GetDepth()=0
virtual void BuildFromBuffer(const void *buffer, const unsigned int xsize, const unsigned int ysize, const unsigned int zsize)=0
A data structure representing a node in an Octree.
void BuildFromImage(Image< TPixel, 3 > *fromImage)
int GetValue(const unsigned int Dim0, const unsigned int Dim1, const unsigned int Dim2)
unsigned int m_Width
Definition itkOctree.h:224
const OctreeNodeBranch * GetColorTable() const override
Image< TPixel, 3 > ImageType
Definition itkOctree.h:136
~Octree() override
OctreeNode m_Tree
Definition itkOctree.h:229
void SetWidth(unsigned int width) override
typename ImageType::Pointer ImageTypePointer
Definition itkOctree.h:137
SmartPointer< Self > Pointer
Definition itkOctree.h:135
OctreeNodeBranch m_ColorTable[ColorTableSize]
Definition itkOctree.h:228
OctreeEnum m_Plane
Definition itkOctree.h:220
OctreeNodeBranch * maskToOctree(const TPixel *Mask, unsigned int width, unsigned int x, unsigned int y, unsigned int z, unsigned int xsize, unsigned int ysize, unsigned int zsize)
OctreeBase Superclass
Definition itkOctree.h:134
int GetColorTableSize() const override
OctreeNode * GetTree() override
unsigned int GetDepth() override
unsigned int m_Depth
Definition itkOctree.h:226
ImageTypePointer GetImage()
void SetDepth(unsigned int depth) override
MappingFunctionType m_MappingFunction
Definition itkOctree.h:230
unsigned int GetWidth() override
void SetTree(OctreeNodeBranch *branch)
Definition itkOctree.h:165
Octree Self
Definition itkOctree.h:133
void SetTrueDims(const unsigned int Dim0, const unsigned int Dim1, const unsigned int Dim2)
void BuildFromBuffer(const void *frombuffer, const unsigned int xsize, const unsigned int ysize, const unsigned int zsize) override
void SetColor(unsigned int color)
Definition itkOctree.h:160
unsigned int m_TrueDims[3]
Definition itkOctree.h:227
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
@ B2_MASKFILE_GRAY
Definition itkOctree.h:31
@ B2_MASKFILE_BLACK
Definition itkOctree.h:29
@ B2_MASKFILE_WHITE
Definition itkOctree.h:30