Developer Documentation
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
MultiGridOctreeData.h
1 /*
2 Copyright (c) 2006, Michael Kazhdan and Matthew Bolitho
3 All rights reserved.
4 
5 Redistribution and use in source and binary forms, with or without modification,
6 are permitted provided that the following conditions are met:
7 
8 Redistributions of source code must retain the above copyright notice, this list of
9 conditions and the following disclaimer. Redistributions in binary form must reproduce
10 the above copyright notice, this list of conditions and the following disclaimer
11 in the documentation and/or other materials provided with the distribution.
12 
13 Neither the name of the Johns Hopkins University nor the names of its contributors
14 may be used to endorse or promote products derived from this software without specific
15 prior written permission.
16 
17 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY
18 EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO THE IMPLIED WARRANTIES
19 OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT
20 SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
21 INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
22 TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR
23 BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
24 CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
25 ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH
26 DAMAGE.
27 */
28 
29 #ifndef MULTI_GRID_OCTREE_DATA_INCLUDED
30 #define MULTI_GRID_OCTREE_DATA_INCLUDED
31 
32 #define GRADIENT_DOMAIN_SOLUTION 1 // Given the constraint vector-field V(p), there are two ways to solve for the coefficients, x, of the indicator function
33  // with respect to the B-spline basis {B_i(p)}
34  // 1] Find x minimizing:
35  // || V(p) - \sum_i \nabla x_i B_i(p) ||^2
36  // which is solved by the system A_1x = b_1 where:
37  // A_1[i,j] = < \nabla B_i(p) , \nabla B_j(p) >
38  // b_1[i] = < \nabla B_i(p) , V(p) >
39  // 2] Formulate this as a Poisson equation:
40  // \sum_i x_i \Delta B_i(p) = \nabla \cdot V(p)
41  // which is solved by the system A_2x = b_2 where:
42  // A_2[i,j] = - < \Delta B_i(p) , B_j(p) >
43  // b_2[i] = - < B_i(p) , \nabla \cdot V(p) >
44  // Although the two system matrices should be the same (assuming that the B_i satisfy dirichlet/neumann boundary conditions)
45  // the constraint vectors can differ when V does not satisfy the Neumann boundary conditions:
46  // A_1[i,j] = \int_R < \nabla B_i(p) , \nabla B_j(p) >
47  // = \int_R [ \nabla \cdot ( B_i(p) \nabla B_j(p) ) - B_i(p) \Delta B_j(p) ]
48  // = \int_dR < N(p) , B_i(p) \nabla B_j(p) > + A_2[i,j]
49  // and the first integral is zero if either f_i is zero on the boundary dR or the derivative of B_i across the boundary is zero.
50  // However, for the constraints we have:
51  // b_1(i) = \int_R < \nabla B_i(p) , V(p) >
52  // = \int_R [ \nabla \cdot ( B_i(p) V(p) ) - B_i(p) \nabla \cdot V(p) ]
53  // = \int_dR < N(p) , B_i(p) V(p) > + b_2[i]
54  // In particular, this implies that if the B_i satisfy the Neumann boundary conditions (rather than Dirichlet),
55  // and V is not zero across the boundary, then the two constraints are different.
56  // Forcing the < V(p) , N(p) > = 0 on the boundary, by killing off the component of the vector-field in the normal direction
57  // (FORCE_NEUMANN_FIELD), makes the two systems equal, and the value of this flag should be immaterial.
58  // Note that under interpretation 1, we have:
59  // \sum_i b_1(i) = < \nabla \sum_ i B_i(p) , V(p) > = 0
60  // because the B_i's sum to one. However, in general, we could have
61  // \sum_i b_2(i) \neq 0.
62  // This could cause trouble because the constant functions are in the kernel of the matrix A, so CG will misbehave if the constraint
63  // has a non-zero DC term. (Again, forcing < V(p) , N(p) > = 0 along the boundary resolves this problem.)
64 
65 #define FORCE_NEUMANN_FIELD 1 // This flag forces the normal component across the boundary of the integration domain to be zero.
66  // This should be enabled if GRADIENT_DOMAIN_SOLUTION is not, so that CG doesn't run into trouble.
67 
68 #define ROBERTO_TOLDO_FIX 1
69 
70 #if !FORCE_NEUMANN_FIELD
71 #pragma message( "[WARNING] Not zeroing out normal component on boundary" )
72 #endif // !FORCE_NEUMANN_FIELD
73 
74 #include "Hash.h"
75 #include "BSplineData.h"
76 
77 char* outputFile=NULL;
78 int echoStdout=0;
79 void DumpOutput( const char* format , ... )
80 {
81  if( outputFile )
82  {
83  FILE* fp = fopen( outputFile , "a" );
84  va_list args;
85  va_start( args , format );
86  vfprintf( fp , format , args );
87  fclose( fp );
88  va_end( args );
89  }
90  if( echoStdout )
91  {
92  va_list args;
93  va_start( args , format );
94  vprintf( format , args );
95  va_end( args );
96  }
97 }
98 void DumpOutput2( char* str , const char* format , ... )
99 {
100  if( outputFile )
101  {
102  FILE* fp = fopen( outputFile , "a" );
103  va_list args;
104  va_start( args , format );
105  vfprintf( fp , format , args );
106  fclose( fp );
107  va_end( args );
108  }
109  if( echoStdout )
110  {
111  va_list args;
112  va_start( args , format );
113  vprintf( format , args );
114  va_end( args );
115  }
116  va_list args;
117  va_start( args , format );
118  vsprintf( str , format , args );
119  va_end( args );
120  if( str[strlen(str)-1]=='\n' ) str[strlen(str)-1] = 0;
121 }
122 
123 
124 
125 typedef float Real;
126 typedef float MatrixReal;
128 
129 
130 class RootInfo
131 {
132 public:
133  const TreeOctNode* node;
134  int edgeIndex;
135  long long key;
136 };
137 
139 {
140 public:
141  static long long EdgeIndex( const TreeOctNode* node , int eIndex , int maxDepth , int index[DIMENSION] );
142  static long long EdgeIndex( const TreeOctNode* node , int eIndex , int maxDepth );
143  static long long FaceIndex( const TreeOctNode* node , int fIndex , int maxDepth,int index[DIMENSION] );
144  static long long FaceIndex( const TreeOctNode* node , int fIndex , int maxDepth );
145  static long long CornerIndex( int depth , const int offSet[DIMENSION] , int cIndex , int maxDepth , int index[DIMENSION] );
146  static long long CornerIndex( const TreeOctNode* node , int cIndex , int maxDepth , int index[DIMENSION] );
147  static long long CornerIndex( const TreeOctNode* node , int cIndex , int maxDepth );
148  static long long CenterIndex( int depth , const int offSet[DIMENSION] , int maxDepth , int index[DIMENSION] );
149  static long long CenterIndex( const TreeOctNode* node , int maxDepth , int index[DIMENSION] );
150  static long long CenterIndex( const TreeOctNode* node , int maxDepth );
151  static long long CornerIndexKey( const int index[DIMENSION] );
152 };
154 {
155 public:
156  Pointer( TreeOctNode* ) treeNodes;
157  int *nodeCount;
158  int maxDepth;
159  SortedTreeNodes( void );
160  ~SortedTreeNodes( void );
161  void set( TreeOctNode& root );
163  {
164  int idx[Cube::CORNERS];
165  CornerIndices( void ) { memset( idx , -1 , sizeof( int ) * Cube::CORNERS ); }
166  int& operator[] ( int i ) { return idx[i]; }
167  const int& operator[] ( int i ) const { return idx[i]; }
168  };
170  {
171  CornerTableData( void ) { cCount=0; }
172  ~CornerTableData( void ) { clear(); }
173  void clear( void ) { cTable.clear() ; cCount = 0; }
174  CornerIndices& operator[] ( const TreeOctNode* node );
175  const CornerIndices& operator[] ( const TreeOctNode* node ) const;
176  CornerIndices& cornerIndices( const TreeOctNode* node );
177  const CornerIndices& cornerIndices( const TreeOctNode* node ) const;
178  int cCount;
179  std::vector< CornerIndices > cTable;
180  std::vector< int > offsets;
181  };
182  void setCornerTable( CornerTableData& cData , const TreeOctNode* rootNode , int depth , int threads ) const;
183  void setCornerTable( CornerTableData& cData , const TreeOctNode* rootNode , int threads ) const { setCornerTable( cData , rootNode , maxDepth-1 , threads ); }
184  void setCornerTable( CornerTableData& cData , int threads ) const { setCornerTable( cData , NULL , maxDepth-1 , threads ); }
185  int getMaxCornerCount( int depth , int maxDepth , int threads ) const ;
186  struct EdgeIndices
187  {
188  int idx[Cube::EDGES];
189  EdgeIndices( void ) { memset( idx , -1 , sizeof( int ) * Cube::EDGES ); }
190  int& operator[] ( int i ) { return idx[i]; }
191  const int& operator[] ( int i ) const { return idx[i]; }
192  };
194  {
195  EdgeTableData( void ) { eCount=0; }
196  ~EdgeTableData( void ) { clear(); }
197  void clear( void ) { eTable.clear() , eCount=0; }
198  EdgeIndices& operator[] ( const TreeOctNode* node );
199  const EdgeIndices& operator[] ( const TreeOctNode* node ) const;
200  EdgeIndices& edgeIndices( const TreeOctNode* node );
201  const EdgeIndices& edgeIndices( const TreeOctNode* node ) const;
202  int eCount;
203  std::vector< EdgeIndices > eTable;
204  std::vector< int > offsets;
205  };
206  void setEdgeTable( EdgeTableData& eData , const TreeOctNode* rootNode , int depth , int threads );
207  void setEdgeTable( EdgeTableData& eData , const TreeOctNode* rootNode , int threads ) { setEdgeTable( eData , rootNode , maxDepth-1 , threads ); }
208  void setEdgeTable( EdgeTableData& eData , int threads ) { setEdgeTable( eData , NULL , maxDepth-1 , threads ); }
209  int getMaxEdgeCount( const TreeOctNode* rootNode , int depth , int threads ) const ;
210 };
211 
213 {
214 public:
215  int nodeIndex;
216  union
217  {
218  int mcIndex;
219  struct
220  {
221  Real centerWeightContribution;
222  int normalIndex;
223  };
224  };
225  Real constraint , solution;
226  int pointIndex;
227 
228  TreeNodeData(void);
229  ~TreeNodeData(void);
230 };
231 
232 
233 template< int Degree >
234 class Octree
235 {
236  SortedTreeNodes _sNodes;
237  int _minDepth;
238  bool _constrainValues;
239  int _boundaryType;
240  Real _scale;
241  Point3D< Real > _center;
242  std::vector< int > _pointCount;
243  struct PointData
244  {
245  Point3D< Real > position;
246  Real coarserValue;
247  Real weight;
248  PointData( Point3D< Real > p=Point3D< Real >() , Real w=0 ):position(p),coarserValue(Real(0)),weight(w){}
249  };
250  std::vector< PointData > _points;
251 
252  bool _inBounds( Point3D< Real > ) const;
253 
254  Real radius;
255  int width;
256  Real GetLaplacian( const int index[DIMENSION] ) const;
257  // Note that this is a slight misnomer. We're only taking the diveregence/Laplacian in the weak sense, so there is a change of sign.
258  Real GetLaplacian( const TreeOctNode* node1 , const TreeOctNode* node2 ) const;
259  Real GetDivergence( const TreeOctNode* node1 , const TreeOctNode* node2 , const Point3D<Real>& normal1 ) const;
260  Real GetDivergenceMinusLaplacian( const TreeOctNode* node1 , const TreeOctNode* node2 , Real value1 , const Point3D<Real>& normal1 ) const;
261 
263  {
264  public:
265  int adjacencyCount;
266  void Function(const TreeOctNode* node1,const TreeOctNode* node2);
267  };
269  public:
270  int *adjacencies,adjacencyCount;
271  void Function(const TreeOctNode* node1,const TreeOctNode* node2);
272  };
273 
275  public:
276  int depth;
277  void Function(TreeOctNode* node1,const TreeOctNode* node2);
278  };
280  {
281  public:
282  int fIndex , maxDepth;
283  std::vector< std::pair< RootInfo , RootInfo > >* edges;
284  hash_map< long long , std::pair< RootInfo , int > >* vertexCount;
285  void Function( const TreeOctNode* node1 , const TreeOctNode* node2 );
286  };
287 
288  int _SolveFixedDepthMatrix( int depth , const SortedTreeNodes& sNodes , Real* subConstraints , bool showResidual , int minIters , double accuracy , bool noSolve = false , int fixedIters=-1 );
289  int _SolveFixedDepthMatrix( int depth , const SortedTreeNodes& sNodes , Real* subConstraints , int startingDepth , bool showResidual , int minIters , double accuracy , bool noSolve = false , int fixedIters=-1 );
290 
291  void SetMatrixRowBounds( const TreeOctNode* node , int rDepth , const int rOff[3] , int& xStart , int& xEnd , int& yStart , int& yEnd , int& zStart , int& zEnd ) const;
292  int GetMatrixRowSize( const TreeOctNode::Neighbors5& neighbors5 ) const;
293  int GetMatrixRowSize( const TreeOctNode::Neighbors5& neighbors5 , int xStart , int xEnd , int yStart , int yEnd , int zStart , int zEnd ) const;
294  int SetMatrixRow( const TreeOctNode::Neighbors5& neighbors5 , Pointer( MatrixEntry< MatrixReal > ) row , int offset , const double stencil[5][5][5] ) const;
295  int SetMatrixRow( const TreeOctNode::Neighbors5& neighbors5 , Pointer( MatrixEntry< MatrixReal > ) row , int offset , const double stencil[5][5][5] , int xStart , int xEnd , int yStart , int yEnd , int zStart , int zEnd ) const;
296  void SetDivergenceStencil( int depth , Point3D< double > stencil[5][5][5] , bool scatter ) const;
297  void SetLaplacianStencil( int depth , double stencil[5][5][5] ) const;
298  template< class C , int N > struct Stencil{ C values[N][N][N]; };
299  void SetLaplacianStencils( int depth , Stencil< double , 5 > stencil[2][2][2] ) const;
300  void SetDivergenceStencils( int depth , Stencil< Point3D< double > , 5 > stencil[2][2][2] , bool scatter ) const;
301  void SetEvaluationStencils( int depth , Stencil< Real , 3 > stencil1[8] , Stencil< Real , 3 > stencil2[8][8] ) const;
302 
303  static void UpdateCoarserSupportBounds( const TreeOctNode* node , int& startX , int& endX , int& startY , int& endY , int& startZ , int& endZ );
304  void UpdateConstraintsFromCoarser( const TreeOctNode::NeighborKey5& neighborKey5 , TreeOctNode* node , Real* metSolution , const Stencil< double , 5 >& stencil ) const;
305  void SetCoarserPointValues( int depth , const SortedTreeNodes& sNodes , Real* metSolution );
306  Real WeightedCoarserFunctionValue( const TreeOctNode::NeighborKey3& neighborKey3 , const TreeOctNode* node , Real* metSolution ) const;
307  void UpSampleCoarserSolution( int depth , const SortedTreeNodes& sNodes , PoissonVector< Real >& solution ) const;
308  void DownSampleFinerConstraints( int depth , SortedTreeNodes& sNodes ) const;
309  template< class C > void DownSample( int depth , const SortedTreeNodes& sNodes , C* constraints ) const;
310  template< class C > void UpSample( int depth , const SortedTreeNodes& sNodes , C* coefficients ) const;
311  int GetFixedDepthLaplacian( SparseSymmetricMatrix< Real >& matrix , int depth , const SortedTreeNodes& sNodes , Real* subConstraints );
312  int GetRestrictedFixedDepthLaplacian( SparseSymmetricMatrix< Real >& matrix , int depth , const int* entries , int entryCount , const TreeOctNode* rNode, Real radius , const SortedTreeNodes& sNodes , Real* subConstraints );
313 
314  void SetIsoCorners( Real isoValue , TreeOctNode* leaf , SortedTreeNodes::CornerTableData& cData , Pointer( char ) valuesSet , Pointer( Real ) values , TreeOctNode::ConstNeighborKey3& nKey , const Real* metSolution , const Stencil< Real , 3 > stencil1[8] , const Stencil< Real , 3 > stencil2[8][8] );
315  static int IsBoundaryFace( const TreeOctNode* node , int faceIndex , int subdivideDepth );
316  static int IsBoundaryEdge( const TreeOctNode* node , int edgeIndex , int subdivideDepth );
317  static int IsBoundaryEdge( const TreeOctNode* node , int dir , int x , int y , int subidivideDepth );
318 
319  // For computing the iso-surface there is a lot of re-computation of information across shared geometry.
320  // For function values we don't care so much.
321  // For edges we need to be careful so that the mesh remains water-tight
323  {
324  // Edge to iso-vertex map
325  hash_map< long long , int > boundaryRoots;
326  // Vertex to ( value , normal ) map
327  hash_map< long long , std::pair< Real , Point3D< Real > > > *boundaryValues;
328  Pointer( int ) interiorRoots;
329  Pointer( Real ) cornerValues;
330  Pointer( Point3D< Real > ) cornerNormals;
331  Pointer( char ) cornerValuesSet;
332  Pointer( char ) cornerNormalsSet;
333  Pointer( char ) edgesSet;
334  };
335 
336  int SetBoundaryMCRootPositions( int sDepth , Real isoValue , RootData& rootData , CoredMeshData* mesh , int nonLinearFit );
337  int SetMCRootPositions( TreeOctNode* node , int sDepth , Real isoValue , TreeOctNode::ConstNeighborKey5& neighborKey5 , RootData& rootData ,
338  std::vector< Point3D< Real > >* interiorPositions , CoredMeshData* mesh , const Real* metSolution , int nonLinearFit );
339  int GetMCIsoTriangles( TreeOctNode* node , CoredMeshData* mesh , RootData& rootData ,
340  std::vector< Point3D< Real > >* interiorPositions , int offSet , int sDepth , bool polygonMesh , std::vector< Point3D< Real > >* barycenters );
341  static int AddTriangles( CoredMeshData* mesh , std::vector<CoredPointIndex>& edges , std::vector< Point3D< Real > >* interiorPositions , int offSet , bool polygonMesh , std::vector< Point3D< Real > >* barycenters );
342 
343 
344  void GetMCIsoEdges( TreeOctNode* node , int sDepth , std::vector< std::pair< RootInfo , RootInfo > >& edges );
345  static int GetEdgeLoops( std::vector< std::pair< RootInfo , RootInfo > >& edges , std::vector< std::vector< std::pair< RootInfo , RootInfo > > >& loops);
346  static int InteriorFaceRootCount( const TreeOctNode* node , const int &faceIndex , int maxDepth );
347  static int EdgeRootCount( const TreeOctNode* node , int edgeIndex , int maxDepth );
348  static void GetRootSpan( const RootInfo& ri , Point3D< Real >& start , Point3D< Real >& end );
349  int GetRoot( const RootInfo& ri , Real isoValue , TreeOctNode::ConstNeighborKey5& neighborKey5 , Point3D<Real> & position , RootData& rootData , int sDepth , const Real* metSolution , int nonLinearFit );
350  static int GetRootIndex( const TreeOctNode* node , int edgeIndex , int maxDepth , RootInfo& ri );
351  static int GetRootIndex( const TreeOctNode* node , int edgeIndex , int maxDepth , int sDepth , RootInfo& ri );
352  static int GetRootIndex( const RootInfo& ri , RootData& rootData , CoredPointIndex& index );
353  static int GetRootPair( const RootInfo& root , int maxDepth , RootInfo& pair );
354 
355  int UpdateWeightContribution( TreeOctNode* node , const Point3D<Real>& position , TreeOctNode::NeighborKey3& neighborKey , Real weight=Real(1.0) );
356  Real GetSampleWeight( TreeOctNode* node , const Point3D<Real>& position , TreeOctNode::NeighborKey3& neighborKey );
357  void GetSampleDepthAndWeight( TreeOctNode* node , const Point3D<Real>& position , TreeOctNode::NeighborKey3& neighborKey , Real samplesPerNode , Real& depth , Real& weight );
358  int SplatOrientedPoint( TreeOctNode* node , const Point3D<Real>& point , const Point3D<Real>& normal , TreeOctNode::NeighborKey3& neighborKey );
359  int SplatOrientedPoint( TreeOctNode* node , const Point3D<Real>& point , const Point3D<Real>& normal , TreeOctNode::NeighborKey5& neighborKey );
360  Real SplatOrientedPoint( const Point3D<Real>& point , const Point3D<Real>& normal , TreeOctNode::NeighborKey3& neighborKey , int kernelDepth , Real samplesPerNode , int minDepth , int maxDepth );
361  Real SplatOrientedPoint( const Point3D<Real>& point , const Point3D<Real>& normal , TreeOctNode::NeighborKey3& neighborKey3 , TreeOctNode::NeighborKey5& neighborKey5 , int kernelDepth , Real samplesPerNode , int minDepth , int maxDepth );
362 
363  int HasNormals(TreeOctNode* node,Real epsilon);
364  Real getCornerValue( const TreeOctNode::ConstNeighborKey3& neighborKey3 , const TreeOctNode* node , int corner , const Real* metSolution );
365  Point3D< Real > getCornerNormal( const TreeOctNode::ConstNeighborKey5& neighborKey5 , const TreeOctNode* node , int corner , const Real* metSolution );
366  Real getCornerValue( const TreeOctNode::ConstNeighborKey3& neighborKey3 , const TreeOctNode* node , int corner , const Real* metSolution , const Real stencil1[3][3][3] , const Real stencil2[3][3][3] );
367  Real getCenterValue( const TreeOctNode::ConstNeighborKey3& neighborKey3 , const TreeOctNode* node );
368  static bool _IsInset( const TreeOctNode* node );
369  static bool _IsInsetSupported( const TreeOctNode* node );
370 public:
371  int threads;
372  static double maxMemoryUsage;
373  static double MemoryUsage( void );
374  std::vector< Point3D<Real> >* normals;
375  Real postDerivativeSmooth;
376  TreeOctNode tree;
378  Octree( void );
379 
380  void setBSplineData( int maxDepth , int boundaryType=BSplineElements< Degree >::NONE );
381  void finalize( int subdivisionDepth );
382  int refineBoundary( int subdivisionDepth );
383  Pointer( Real ) GetSolutionGrid( int& res , Real isoValue=0.f , int depth=-1 );
384  int setTree( char* fileName , int maxDepth , int minDepth , int kernelDepth , Real samplesPerNode ,
385  Real scaleFactor , int useConfidence , Real constraintWeight , int adaptiveExponent , XForm4x4< Real > xForm=XForm4x4< Real >::Identity() );
386  int setTreeMemory( std::vector< Real >& _pts_stream, int maxDepth , int minDepth ,
387  int splatDepth , Real samplesPerNode , Real scaleFactor ,
388  int useConfidence , Real constraintWeight , int adaptiveExponent , XForm4x4< Real > xForm=XForm4x4< Real >::Identity() );
389  void SetLaplacianConstraints(void);
390  void ClipTree(void);
391  int LaplacianMatrixIteration( int subdivideDepth , bool showResidual , int minIters , double accuracy , int maxSolveDepth , int fixedIters );
392 
393  Real GetIsoValue( void );
394  void GetMCIsoTriangles( Real isoValue , int subdivideDepth , CoredMeshData* mesh , int fullDepthIso=0 , int nonLinearFit=1 , bool addBarycenter=false , bool polygonMesh=false );
395 };
396 
397 #ifndef DOXY_IGNORE_THIS
398 #include "MultiGridOctreeData.inl"
399 #endif
400 #endif // MULTI_GRID_OCTREE_DATA_INCLUDED
STL namespace.