ACloudViewer  3.9.4
A Modern Library for 3D Data Processing
MultiGridOctreeData.h
Go to the documentation of this file.
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 // [COMMENTS]
29 // -- Throughout the code, should make a distinction between indices and offsets
30 // -- Make an instance of _Evaluate that samples the finite-elements correctly (specifically, to handle the boundaries)
31 // -- Make functions like depthAndOffset parity dependent (ideally all "depth"s should be relative to the B-Slpline resolution
32 // -- Make all points relative to the unit-cube, regardless of degree parity
33 // -- It's possible that for odd degrees, the iso-surfacing will fail because the leaves in the SortedTreeNodes do not form a partition of space
34 // -- [MAYBE] Treat normal field as a sum of delta functions, rather than a smoothed signal (again, so that high degrees aren't forced to generate smooth reconstructions)
35 // -- [MAYBE] Make the degree of the B-Spline with which the normals are splatted independent of the degree of the FEM system. (This way, higher degree systems aren't forced to generate smoother normal fields.)
36 
37 // [TODO]
38 // -- Currently, the implementation assumes that the boundary constraints are the same for vector fields and scalar fields
39 // -- Fix up the ordering in the divergence evaluation
40 
41 #ifndef MULTI_GRID_OCTREE_DATA_INCLUDED
42 #define MULTI_GRID_OCTREE_DATA_INCLUDED
43 
44 #define NEW_CODE 1
45 #define NEW_NEW_CODE 0 // Enabling this ensures that all the nodes contained in the support of the normal field are in the tree
46 
47 #define DATA_DEGREE 1 // The order of the B-Spline used to splat in data for color interpolation
48 #define WEIGHT_DEGREE 2 // The order of the B-Spline used to splat in the weights for density estimation
49 #define NORMAL_DEGREE 2 // The order of the B-Spline used to splat int the normals for constructing the Laplacian constraints
50 
51 //#define MAX_MEMORY_GB 15
52 #define MAX_MEMORY_GB 0
53 
54 #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
55  // with respect to the B-spline basis {B_i(p)}
56  // 1] Find x minimizing:
57  // || V(p) - \sum_i \nabla x_i B_i(p) ||^2
58  // which is solved by the system A_1x = b_1 where:
59  // A_1[i,j] = < \nabla B_i(p) , \nabla B_j(p) >
60  // b_1[i] = < \nabla B_i(p) , V(p) >
61  // 2] Formulate this as a Poisson equation:
62  // \sum_i x_i \Delta B_i(p) = \nabla \cdot V(p)
63  // which is solved by the system A_2x = b_2 where:
64  // A_2[i,j] = - < \Delta B_i(p) , B_j(p) >
65  // b_2[i] = - < B_i(p) , \nabla \cdot V(p) >
66  // Although the two system matrices should be the same (assuming that the B_i satisfy dirichlet/neumann boundary conditions)
67  // the constraint vectors can differ when V does not satisfy the Neumann boundary conditions:
68  // A_1[i,j] = \int_R < \nabla B_i(p) , \nabla B_j(p) >
69  // = \int_R [ \nabla \cdot ( B_i(p) \nabla B_j(p) ) - B_i(p) \Delta B_j(p) ]
70  // = \int_dR < N(p) , B_i(p) \nabla B_j(p) > + A_2[i,j]
71  // 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.
72  // However, for the constraints we have:
73  // b_1(i) = \int_R < \nabla B_i(p) , V(p) >
74  // = \int_R [ \nabla \cdot ( B_i(p) V(p) ) - B_i(p) \nabla \cdot V(p) ]
75  // = \int_dR < N(p) , B_i(p) V(p) > + b_2[i]
76  // In particular, this implies that if the B_i satisfy the Neumann boundary conditions (rather than Dirichlet),
77  // and V is not zero across the boundary, then the two constraints are different.
78  // Forcing the < V(p) , N(p) > = 0 on the boundary, by killing off the component of the vector-field in the normal direction
79  // (FORCE_NEUMANN_FIELD), makes the two systems equal, and the value of this flag should be immaterial.
80  // Note that under interpretation 1, we have:
81  // \sum_i b_1(i) = < \nabla \sum_ i B_i(p) , V(p) > = 0
82  // because the B_i's sum to one. However, in general, we could have
83  // \sum_i b_2(i) \neq 0.
84  // This could cause trouble because the constant functions are in the kernel of the matrix A, so CG will misbehave if the constraint
85  // has a non-zero DC term. (Again, forcing < V(p) , N(p) > = 0 along the boundary resolves this problem.)
86 
87 #define FORCE_NEUMANN_FIELD 1 // This flag forces the normal component across the boundary of the integration domain to be zero.
88  // This should be enabled if GRADIENT_DOMAIN_SOLUTION is not, so that CG doesn't run into trouble.
89 
90 #if !FORCE_NEUMANN_FIELD
91 #pragma message( "[WARNING] Not zeroing out normal component on boundary" )
92 #endif // !FORCE_NEUMANN_FIELD
93 
94 #include "Hash.h"
95 #include "BSplineData.h"
96 #include "PointStream.h"
97 
98 #ifndef _OPENMP
99 int omp_get_num_procs( void ){ return 1; }
100 int omp_get_thread_num( void ){ return 0; }
101 #endif // _OPENMP
102 
104 {
105 public:
106  static size_t NodeCount;
108  char flags;
109 
110  TreeNodeData( void );
111  ~TreeNodeData( void );
112 };
113 
115 {
117 public:
118  static const int VERTEX_COORDINATE_SHIFT = ( sizeof( long long ) * 8 ) / 3;
119  static long long EdgeIndex( const TreeOctNode* node , int eIndex , int maxDepth , int index[DIMENSION] );
120  static long long EdgeIndex( const TreeOctNode* node , int eIndex , int maxDepth );
121  static long long FaceIndex( const TreeOctNode* node , int fIndex , int maxDepth,int index[DIMENSION] );
122  static long long FaceIndex( const TreeOctNode* node , int fIndex , int maxDepth );
123  static long long CornerIndex( const TreeOctNode* node , int cIndex , int maxDepth , int index[DIMENSION] );
124  static long long CornerIndex( const TreeOctNode* node , int cIndex , int maxDepth );
125  static long long CenterIndex( const TreeOctNode* node , int maxDepth , int index[DIMENSION] );
126  static long long CenterIndex( const TreeOctNode* node , int maxDepth );
127  static long long CornerIndex( int depth , const int offSet[DIMENSION] , int cIndex , int maxDepth , int index[DIMENSION] );
128  static long long CenterIndex( int depth , const int offSet[DIMENSION] , int maxDepth , int index[DIMENSION] );
129  static long long CornerIndexKey( const int index[DIMENSION] );
130 };
131 
132 // This class stores the octree nodes, sorted by depth and then by z-slice.
133 // To support primal representations, the initializer takes a function that
134 // determines if a node should be included/indexed in the sorted list.
136 {
138 protected:
140  int _levels;
141 public:
143  int begin( int depth ) const{ return _sliceStart[depth][0]; }
144  int end( int depth ) const{ return _sliceStart[depth][(size_t)1<<depth]; }
145  int begin( int depth , int slice ) const{ return _sliceStart[depth][slice ] ; }
146  int end( int depth , int slice ) const{ if(depth<0||depth>=_levels||slice<0||slice>=(1<<depth)) printf( "uh oh\n" ) ; return _sliceStart[depth][slice+1]; }
147  int size( void ) const { return _sliceStart[_levels-1][(size_t)1<<(_levels-1)]; }
148  int size( int depth ) const { if(depth<0||depth>=_levels) printf( "uhoh\n" ); return _sliceStart[depth][(size_t)1<<depth] - _sliceStart[depth][0]; }
149  int size( int depth , int slice ) const { return _sliceStart[depth][slice+1] - _sliceStart[depth][slice]; }
150  int levels( void ) const { return _levels; }
151 
154  void set( TreeOctNode& root , std::vector< int >* map );
155  void set( TreeOctNode& root );
156 
157  template< int Indices >
158  struct _Indices
159  {
160  int idx[Indices];
161  _Indices( void ){ memset( idx , -1 , sizeof( int ) * Indices ); }
162  int& operator[] ( int i ) { return idx[i]; }
163  const int& operator[] ( int i ) const { return idx[i]; }
164  };
168 
170  {
176  ~SliceTableData( void ){ clear(); }
180  const SquareCornerIndices& cornerIndices( const TreeOctNode* node ) const;
181  const SquareCornerIndices& cornerIndices( int idx ) const;
184  const SquareEdgeIndices& edgeIndices( const TreeOctNode* node ) const;
185  const SquareEdgeIndices& edgeIndices( int idx ) const;
188  const SquareFaceIndices& faceIndices( const TreeOctNode* node ) const;
189  const SquareFaceIndices& faceIndices( int idx ) const;
190  protected:
191  Pointer( int ) _cMap;
192  Pointer( int ) _eMap;
193  Pointer( int ) _fMap;
194  friend class SortedTreeNodes;
195  };
197  {
202  ~XSliceTableData( void ){ clear(); }
203  void clear( void ) { DeletePointer( fTable ) ; DeletePointer( eTable ) ; fCount = eCount = 0; }
206  const SquareCornerIndices& edgeIndices( const TreeOctNode* node ) const;
207  const SquareCornerIndices& edgeIndices( int idx ) const;
210  const SquareEdgeIndices& faceIndices( const TreeOctNode* node ) const;
211  const SquareEdgeIndices& faceIndices( int idx ) const;
212  protected:
213  Pointer( int ) _eMap;
214  Pointer( int ) _fMap;
215  friend class SortedTreeNodes;
216  };
217  void setSliceTableData ( SliceTableData& sData , int depth , int offset , int threads ) const;
218  void setXSliceTableData( XSliceTableData& sData , int depth , int offset , int threads ) const;
219 };
220 
221 template< int Degree >
222 struct PointSupportKey : public OctNode< TreeNodeData >::NeighborKey< BSplineEvaluationData< Degree >::SupportEnd , -BSplineEvaluationData< Degree >::SupportStart >
223 {
226  static const int Size = LeftRadius + RightRadius + 1;
227 };
228 template< int Degree >
229 struct ConstPointSupportKey : public OctNode< TreeNodeData >::ConstNeighborKey< BSplineEvaluationData< Degree >::SupportEnd , -BSplineEvaluationData< Degree >::SupportStart >
230 {
233  static const int Size = LeftRadius + RightRadius + 1;
234 };
235 
236 template< class Real >
237 struct PointData
238 {
241  Real weight;
242  PointData( Point3D< Real > p=Point3D< Real >() , Real w=0 ) { position = p , weight = w , weightedCoarserDValue = Real(0); }
243 };
244 template< class Data , int Degree >
246 {
247  std::vector< int > indices;
248  std::vector< Data > data;
249  template< class TreeNodeData >
250  int index( const OctNode< TreeNodeData >* node ) const { return ( !node || node->nodeData.nodeIndex<0 || node->nodeData.nodeIndex>=(int)indices.size() ) ? -1 : indices[ node->nodeData.nodeIndex ]; }
251 #if NEW_NEW_CODE
252  int index( int nodeIndex ) const { return ( nodeIndex<0 || nodeIndex>=(int)indices.size() ) ? -1 : indices[ nodeIndex ]; }
253 #endif // NEW_NEW_CODE
254  void resize( size_t sz ){ indices.resize( sz , -1 ); }
255  void remapIndices( const std::vector< int >& map )
256  {
257  std::vector< int > temp = indices;
258  indices.resize( map.size() );
259  for( size_t i=0 ; i<map.size() ; i++ )
260  if( map[i]<(int)temp.size() ) indices[i] = temp[ map[i] ];
261  else indices[i] = -1;
262  }
263 };
264 template< class Data , int Degree >
266 {
267  Pointer( Data ) data;
268  DenseNodeData( void ) { data = NullPointer( Data ); }
269  DenseNodeData( size_t sz ){ if( sz ) data = NewPointer< Data >( sz ) ; else data = NullPointer( Data ); }
270  void resize( size_t sz ){ DeletePointer( data ) ; if( sz ) data = NewPointer< Data >( sz ) ; else data = NullPointer( Data ); }
271  Data& operator[] ( int idx ) { return data[idx]; }
272  const Data& operator[] ( int idx ) const { return data[idx]; }
273 };
274 
275 template< class C , int N > struct Stencil{ C values[N][N][N]; };
276 
277 template< int Degree1 , int Degree2 >
279 {
281  static const int OverlapSize = BSplineIntegrationData< Degree1 , Degree2 >::OverlapSize;
282  static const int OverlapStart = BSplineIntegrationData< Degree1 , Degree2 >::OverlapStart;
283  static const int OverlapEnd = BSplineIntegrationData< Degree1 , Degree2 >::OverlapEnd;
284 public:
285  static double GetLaplacian ( const typename FunctionIntegrator:: Integrator& integrator , const int off1[3] , const int off2[3] );
286  static double GetLaplacian ( const typename FunctionIntegrator::ChildIntegrator& integrator , const int off1[3] , const int off2[3] );
287  static double GetDivergence1( const typename FunctionIntegrator:: Integrator& integrator , const int off1[3] , const int off2[3] , Point3D< double > normal1 );
288  static double GetDivergence1( const typename FunctionIntegrator::ChildIntegrator& integrator , const int off1[3] , const int off2[3] , Point3D< double > normal1 );
289  static double GetDivergence2( const typename FunctionIntegrator:: Integrator& integrator , const int off1[3] , const int off2[3] , Point3D< double > normal2 );
290  static double GetDivergence2( const typename FunctionIntegrator::ChildIntegrator& integrator , const int off1[3] , const int off2[3] , Point3D< double > normal2 );
291  static Point3D< double > GetDivergence1 ( const typename FunctionIntegrator:: Integrator& integrator , const int off1[3] , const int off2[3] );
292  static Point3D< double > GetDivergence1 ( const typename FunctionIntegrator::ChildIntegrator& integrator , const int off1[3] , const int off2[3] );
293  static Point3D< double > GetDivergence2 ( const typename FunctionIntegrator:: Integrator& integrator , const int off1[3] , const int off2[3] );
294  static Point3D< double > GetDivergence2 ( const typename FunctionIntegrator::ChildIntegrator& integrator , const int off1[3] , const int off2[3] );
295  static void SetCentralDivergenceStencil ( const typename FunctionIntegrator:: Integrator& integrator , Stencil< Point3D< double > , OverlapSize >& stencil , bool scatter );
296  static void SetCentralDivergenceStencils( const typename FunctionIntegrator::ChildIntegrator& integrator , Stencil< Point3D< double > , OverlapSize > stencil[2][2][2] , bool scatter );
297  static void SetCentralLaplacianStencil ( const typename FunctionIntegrator:: Integrator& integrator , Stencil< double , OverlapSize >& stencil );
298  static void SetCentralLaplacianStencils ( const typename FunctionIntegrator::ChildIntegrator& integrator , Stencil< double , OverlapSize > stencil[2][2][2] );
299 };
300 
301 // Note that throughout this code, the "depth" parameter refers to the depth in the octree, not the corresponding depth
302 // of the B-Spline element
303 template< class Real >
304 class Octree
305 {
307 public:
308  template< int FEMDegree > static void FunctionIndex( const TreeOctNode* node , int idx[3] );
309 
310  typedef typename TreeOctNode:: NeighborKey< 1 , 1 > AdjacenctNodeKey;
311  typedef typename TreeOctNode::ConstNeighborKey< 1 , 1 > ConstAdjacenctNodeKey;
312 
313  template< class V >
315  {
316  V v;
317  Real w;
318  ProjectiveData( V vv=V(0) , Real ww=Real(0) ) : v(vv) , w(ww) { }
319  operator V (){ return w!=0 ? v/w : v*w; }
320  ProjectiveData& operator += ( const ProjectiveData& p ){ v += p.v , w += p.w ; return *this; }
321  ProjectiveData& operator -= ( const ProjectiveData& p ){ v -= p.v , w -= p.w ; return *this; }
322  ProjectiveData& operator *= ( Real s ){ v *= s , w *= s ; return *this; }
323  ProjectiveData& operator /= ( Real s ){ v /= s , w /= s ; return *this; }
324  ProjectiveData operator + ( const ProjectiveData& p ) const { return ProjectiveData( v+p.v , w+p.w ); }
325  ProjectiveData operator - ( const ProjectiveData& p ) const { return ProjectiveData( v-p.v , w-p.w ); }
326  ProjectiveData operator * ( Real s ) const { return ProjectiveData( v*s , w*s ); }
327  ProjectiveData operator / ( Real s ) const { return ProjectiveData( v/s , w/s ); }
328  };
329  template< int FEMDegree > static bool IsValidNode( const TreeOctNode* node , bool dirichlet );
330 protected:
331  template< int FEMDegree > bool _IsValidNode( const TreeOctNode* node ) const { return node && ( node->nodeData.flags & ( 1<<( FEMDegree&1 ) ) ) ; }
332 
342  Real _scale;
345 
346  bool _InBounds( Point3D< Real > ) const;
347  template< int FEMDegree > static int _Dimension( int depth ){ return BSplineData< FEMDegree >::Dimension( depth-1 ); }
348  static int _Resolution( int depth ){ return 1<<(depth-1); }
349  template< int FEMDegree > static bool _IsInteriorlySupported( int d , int x , int y , int z )
350  {
351  if( d-1>=0 )
352  {
353  int begin , end;
355  return ( x>=begin && x<end && y>=begin && y<end && z>=begin && z<end );
356  }
357  else return false;
358  }
359  template< int FEMDegree > static bool _IsInteriorlySupported( const TreeOctNode* node )
360  {
361  if( !node ) return false;
362  int d , off[3];
363  node->depthAndOffset( d , off );
364  return _IsInteriorlySupported< FEMDegree >( d , off[0] , off[1] , off[2] );
365  }
366  template< int FEMDegree1 , int FEMDegree2 > static bool _IsInteriorlyOverlapped( int d , int x , int y , int z )
367  {
368  if( d-1>=0 )
369  {
370  int begin , end;
372  return ( x>=begin && x<end && y>=begin && y<end && z>=begin && z<end );
373  }
374  else return false;
375  }
376  template< int FEMDegree1 , int FEMDegree2 > static bool _IsInteriorlyOverlapped( const TreeOctNode* node )
377  {
378  if( !node ) return false;
379  int d , off[3];
380  node->depthAndOffset( d , off );
381  return _IsInteriorlyOverlapped< FEMDegree1 , FEMDegree2 >( d , off[0] , off[1] , off[2] );
382  }
383  static void _DepthAndOffset( const TreeOctNode* node , int& d , int off[3] ){ node->depthAndOffset( d , off ) ; d -= 1; }
384  static int _Depth( const TreeOctNode* node ){ return node->depth()-1; }
385  static void _StartAndWidth( const TreeOctNode* node , Point3D< Real >& start , Real& width )
386  {
387  int d , off[3];
388  _DepthAndOffset( node , d , off );
389  if( d>=0 ) width = Real( 1.0 / (1<< d ) );
390  else width = Real( 1.0 * (1<<(-d)) );
391  for( int dd=0 ; dd<DIMENSION ; dd++ ) start[dd] = Real( off[dd] ) * width;
392  }
393  static void _CenterAndWidth( const TreeOctNode* node , Point3D< Real >& center , Real& width )
394  {
395  int d , off[3];
396  _DepthAndOffset( node , d , off );
397  width = Real( 1.0 / (1<<d) );
398  for( int dd=0 ; dd<DIMENSION ; dd++ ) center[dd] = Real( off[dd] + 0.5 ) * width;
399  }
400  template< int LeftRadius , int RightRadius >
401  static typename TreeOctNode::ConstNeighbors< LeftRadius + RightRadius + 1 >& _Neighbors( TreeOctNode::ConstNeighborKey< LeftRadius , RightRadius >& key , int depth ){ return key.neighbors[ depth + 1 ]; }
402  template< int LeftRadius , int RightRadius >
403  static typename TreeOctNode::Neighbors< LeftRadius + RightRadius + 1 >& _Neighbors( TreeOctNode::NeighborKey< LeftRadius , RightRadius >& key , int depth ){ return key.neighbors[ depth + 1 ]; }
404  template< int LeftRadius , int RightRadius >
405  static const typename TreeOctNode::template Neighbors< LeftRadius + RightRadius + 1 >& _Neighbors( const typename TreeOctNode::template NeighborKey< LeftRadius , RightRadius >& key , int depth ){ return key.neighbors[ depth + 1 ]; }
406  template< int LeftRadius , int RightRadius >
407  static const typename TreeOctNode::template ConstNeighbors< LeftRadius + RightRadius + 1 >& _Neighbors( const typename TreeOctNode::template ConstNeighborKey< LeftRadius , RightRadius >& key , int depth ){ return key.neighbors[ depth + 1 ]; }
408 
409  static void _SetFullDepth( TreeOctNode* node , int depth );
410  void _setFullDepth( int depth );
411 
413  // System construction code //
414  // MultiGridOctreeData.System.inl //
416  template< int FEMDegree >
417  void _setMultiColorIndices( int start , int end , std::vector< std::vector< int > >& indices ) const;
418  template< int FEMDegree >
419  int _SolveSystemGS( const BSplineData< FEMDegree >& bsData , SparseNodeData< PointData< Real > , 0 >& pointInfo , int depth , DenseNodeData< Real , FEMDegree >& solution , DenseNodeData< Real , FEMDegree >& constraints , DenseNodeData< Real , FEMDegree >& metSolutionConstraints , int iters , bool coarseToFine , bool showResidual=false , double* bNorm2=NULL , double* inRNorm2=NULL , double* outRNorm2=NULL , bool forceSilent=false );
420  template< int FEMDegree >
421  int _SolveSystemCG( const BSplineData< FEMDegree >& bsData , SparseNodeData< PointData< Real > , 0 >& pointInfo , int depth , DenseNodeData< Real , FEMDegree >& solution , DenseNodeData< Real , FEMDegree >& constraints , DenseNodeData< Real , FEMDegree >& metSolutionConstraints , int iters , bool coarseToFine , bool showResidual=false , double* bNorm2=NULL , double* inRNorm2=NULL , double* outRNorm2=NULL , double accuracy=0 );
422  template< int FEMDegree >
424  template< int FEMDegree >
425  int _GetMatrixRowSize( const typename TreeOctNode::Neighbors< BSplineIntegrationData< FEMDegree , FEMDegree >::OverlapSize >& neighbors ) const;
426 
427  template< int FEMDegree1 , int FEMDegree2 > static void _SetParentOverlapBounds( const TreeOctNode* node , int& startX , int& endX , int& startY , int& endY , int& startZ , int& endZ );
428  template< int FEMDegree >
430  // Updates the constraints @(depth-1) based on the solution coefficients @(depth)
431  template< int FEMDegree >
433  // Evaluate the points @(depth) using coefficients @(depth-1)
434  template< int FEMDegree >
435  void _SetPointValuesFromCoarser( SparseNodeData< PointData< Real > , 0 >& pointInfo , int highDepth , const BSplineData< FEMDegree >& bsData , const DenseNodeData< Real , FEMDegree >& upSampledCoefficients );
436  // Evalutes the solution @(depth) at the points @(depth-1) and updates the met constraints @(depth-1)
437  template< int FEMDegree >
438  void _SetPointConstraintsFromFiner( const SparseNodeData< PointData< Real > , 0 >& pointInfo , int highDepth , const BSplineData< FEMDegree >& bsData , const DenseNodeData< Real , FEMDegree >& finerCoefficients , DenseNodeData< Real , FEMDegree >& metConstraints ) const;
439  template< int FEMDegree >
440  Real _CoarserFunctionValue( Point3D< Real > p , const PointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , const BSplineData< FEMDegree >& bsData , const DenseNodeData< Real , FEMDegree >& upSampledCoefficients ) const;
441  template< int FEMDegree >
442  Real _FinerFunctionValue ( Point3D< Real > p , const PointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , const BSplineData< FEMDegree >& bsData , const DenseNodeData< Real , FEMDegree >& coefficients ) const;
443  template< int FEMDegree >
445  template< int FEMDegree >
447 
448  // Down samples constraints @(depth) to constraints @(depth-1)
449  template< class C , int FEMDegree > void _DownSample( int highDepth , DenseNodeData< C , FEMDegree >& constraints ) const;
450  // Up samples coefficients @(depth-1) to coefficients @(depth)
451  template< class C , int FEMDegree > void _UpSample( int highDepth , DenseNodeData< C , FEMDegree >& coefficients ) const;
452  template< class C , int FEMDegree > static void _UpSample( int highDepth , ConstPointer( C ) lowCoefficients , Pointer( C ) highCoefficients , bool dirichlet , int threads );
453 
455  // Code for splatting point-sample data //
456  // MultiGridOctreeData.WeightedSamples.inl //
458  template< int WeightDegree >
460  template< int WeightDegree >
462  template< int WeightDegree >
464  template< int WeightDegree >
465  void _GetSampleDepthAndWeight( const SparseNodeData< Real , WeightDegree >& densityWeights , const TreeOctNode* node , Point3D< Real > position , ConstPointSupportKey< WeightDegree >& weightKey , Real& depth , Real& weight ) const;
466  template< int WeightDegree >
467  void _GetSampleDepthAndWeight( const SparseNodeData< Real , WeightDegree >& densityWeights , TreeOctNode* node , Point3D< Real > position , PointSupportKey< WeightDegree >& weightKey , Real& depth , Real& weight );
468 public:
469  template< int WeightDegree >
470  void _GetSampleDepthAndWeight( const SparseNodeData< Real , WeightDegree >& densityWeights , Point3D< Real > position , PointSupportKey< WeightDegree >& weightKey , Real& depth , Real& weight );
471  template< int WeightDegree >
473 protected:
474  template< int DataDegree , class V > void _SplatPointData( TreeOctNode* node , Point3D< Real > point , V v , SparseNodeData< V , DataDegree >& data , PointSupportKey< DataDegree >& dataKey );
475  template< int WeightDegree , int DataDegree , class V > Real _SplatPointData( const SparseNodeData< Real , WeightDegree >& densityWeights , Point3D< Real > point , V v , SparseNodeData< V , DataDegree >& data , PointSupportKey< WeightDegree >& weightKey , PointSupportKey< DataDegree >& dataKey , int minDepth , int maxDepth , int dim=DIMENSION );
476  template< int WeightDegree , int DataDegree , class V > void _MultiSplatPointData( const SparseNodeData< Real , WeightDegree >* densityWeights , Point3D< Real > point , V v , SparseNodeData< V , DataDegree >& data , PointSupportKey< WeightDegree >& weightKey , PointSupportKey< DataDegree >& dataKey , int maxDepth , int dim=DIMENSION );
477  template< class V , int DataDegree > V _Evaluate( const DenseNodeData< V , DataDegree >& coefficients , Point3D< Real > p , const BSplineData< DataDegree >& bsData , const ConstPointSupportKey< DataDegree >& neighborKey ) const;
478  template< class V , int DataDegree > V _Evaluate( const SparseNodeData< V , DataDegree >& coefficients , Point3D< Real > p , const BSplineData< DataDegree >& bsData , const ConstPointSupportKey< DataDegree >& dataKey ) const;
479 public:
480  template< class V , int DataDegree > V Evaluate( const DenseNodeData< V , DataDegree >& coefficients , Point3D< Real > p , const BSplineData< DataDegree >& bsData ) const;
481  template< class V , int DataDegree > V Evaluate( const SparseNodeData< V , DataDegree >& coefficients , Point3D< Real > p , const BSplineData< DataDegree >& bsData ) const;
482  template< class V , int DataDegree > Pointer( V ) Evaluate( const DenseNodeData< V , DataDegree >& coefficients , int& res , Real isoValue=0.f , int depth=-1 , bool primal=false );
483 
484  template< int NormalDegree > int _HasNormals( TreeOctNode* node , const SparseNodeData< Point3D< Real > , NormalDegree >& normalInfo );
485  void _MakeComplete( void );
486  void _SetValidityFlags( void );
487  template< int NormalDegree > void _ClipTree( const SparseNodeData< Point3D< Real > , NormalDegree >& normalInfo );
488 
490  // Evaluation Methods //
491  // MultiGridOctreeData.Evaluation //
493  static const int CHILDREN = Cube::CORNERS;
494  template< int FEMDegree >
495  struct _Evaluator
496  {
507 
516  void set( int depth , bool dirichlet );
517  };
518  template< class V , int FEMDegree >
519  V _getCenterValue( const ConstPointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , const DenseNodeData< V , FEMDegree >& solution , const DenseNodeData< V , FEMDegree >& metSolution , const _Evaluator< FEMDegree >& evaluator , bool isInterior ) const;
520  template< class V , int FEMDegree >
521  V _getCornerValue( const ConstPointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , int corner , const DenseNodeData< V , FEMDegree >& solution , const DenseNodeData< V , FEMDegree >& metSolution , const _Evaluator< FEMDegree >& evaluator , bool isInterior ) const;
522  template< class V , int FEMDegree >
523  V _getEdgeValue ( const ConstPointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , int edge , const DenseNodeData< V , FEMDegree >& solution , const DenseNodeData< V , FEMDegree >& metSolution , const _Evaluator< FEMDegree >& evaluator , bool isInterior ) const;
524 
525  template< int FEMDegree >
526  std::pair< Real , Point3D< Real > > _getCornerValueAndGradient( const ConstPointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , int corner , const DenseNodeData< Real , FEMDegree >& solution , const DenseNodeData< Real , FEMDegree >& metSolution , const _Evaluator< FEMDegree >& evaluator , bool isInterior ) const;
527  template< int FEMDegree >
528  std::pair< Real , Point3D< Real > > _getEdgeValueAndGradient ( const ConstPointSupportKey< FEMDegree >& neighborKey , const TreeOctNode* node , int edge , const DenseNodeData< Real , FEMDegree >& solution , const DenseNodeData< Real , FEMDegree >& metSolution , const _Evaluator< FEMDegree >& evaluator , bool isInterior ) const;
529 
531  // Iso-Surfacing Methods //
532  // MultiGridOctreeData.IsoSurface.inl //
534  struct IsoEdge
535  {
536  long long edges[2];
537  IsoEdge( void ){ edges[0] = edges[1] = 0; }
538  IsoEdge( long long v1 , long long v2 ){ edges[0] = v1 , edges[1] = v2; }
539  long long& operator[]( int idx ){ return edges[idx]; }
540  const long long& operator[]( int idx ) const { return edges[idx]; }
541  };
542  struct FaceEdges
543  {
545  int count;
546  };
547  template< class Vertex >
548  struct SliceValues
549  {
552  Pointer( long long ) edgeKeys ; Pointer( char ) edgeSet;
555  hash_map< long long , std::vector< IsoEdge > > faceEdgeMap;
556  hash_map< long long , std::pair< int , Vertex > > edgeVertexMap;
557  hash_map< long long , long long > vertexPairMap;
558 
559  SliceValues( void );
560  ~SliceValues( void );
561  void reset( bool nonLinearFit );
562  protected:
564  };
565  template< class Vertex >
567  {
569  Pointer( long long ) edgeKeys ; Pointer( char ) edgeSet;
571  hash_map< long long , std::vector< IsoEdge > > faceEdgeMap;
572  hash_map< long long , std::pair< int , Vertex > > edgeVertexMap;
573  hash_map< long long , long long > vertexPairMap;
574 
575  XSliceValues( void );
576  ~XSliceValues( void );
577  void reset( void );
578  protected:
580  };
581  template< class Vertex >
582  struct SlabValues
583  {
586  SliceValues< Vertex >& sliceValues( int idx ){ return _sliceValues[idx&1]; }
587  const SliceValues< Vertex >& sliceValues( int idx ) const { return _sliceValues[idx&1]; }
588  XSliceValues< Vertex >& xSliceValues( int idx ){ return _xSliceValues[idx&1]; }
589  const XSliceValues< Vertex >& xSliceValues( int idx ) const { return _xSliceValues[idx&1]; }
590  };
591  template< class Vertex , int FEMDegree >
592  void SetSliceIsoCorners( const DenseNodeData< Real , FEMDegree >& solution , const DenseNodeData< Real , FEMDegree >& coarseSolution , Real isoValue , int depth , int slice , std::vector< SlabValues< Vertex > >& sValues , const _Evaluator< FEMDegree >& evaluator , int threads );
593  template< class Vertex , int FEMDegree >
594  void SetSliceIsoCorners( const DenseNodeData< Real , FEMDegree >& solution , const DenseNodeData< Real , FEMDegree >& coarseSolution , Real isoValue , int depth , int slice , int z , std::vector< SlabValues< Vertex > >& sValues , const _Evaluator< FEMDegree >& evaluator , int threads );
595  template< int WeightDegree , int ColorDegree , class Vertex >
596  void SetSliceIsoVertices( const BSplineData< ColorDegree >* colorBSData , const SparseNodeData< Real , WeightDegree >* densityWeights , const SparseNodeData< ProjectiveData< Point3D< Real > > , ColorDegree >* colorData , Real isoValue , int depth , int slice , int& vOffset , CoredMeshData< Vertex >& mesh , std::vector< SlabValues< Vertex > >& sValues , int threads );
597  template< int WeightDegree , int ColorDegree , class Vertex >
598  void SetSliceIsoVertices( const BSplineData< ColorDegree >* colorBSData , const SparseNodeData< Real , WeightDegree >* densityWeights , const SparseNodeData< ProjectiveData< Point3D< Real > > , ColorDegree >* colorData , Real isoValue , int depth , int slice , int z , int& vOffset , CoredMeshData< Vertex >& mesh , std::vector< SlabValues< Vertex > >& sValues , int threads );
599  template< int WeightDegree , int ColorDegree , class Vertex >
600  void SetXSliceIsoVertices( const BSplineData< ColorDegree >* colorBSData , const SparseNodeData< Real , WeightDegree >* densityWeights , const SparseNodeData< ProjectiveData< Point3D< Real > > , ColorDegree >* colorData , Real isoValue , int depth , int slab , int& vOffset , CoredMeshData< Vertex >& mesh , std::vector< SlabValues< Vertex > >& sValues , int threads );
601  template< class Vertex >
602  void CopyFinerSliceIsoEdgeKeys( int depth , int slice , std::vector< SlabValues< Vertex > >& sValues , int threads );
603  template< class Vertex >
604  void CopyFinerSliceIsoEdgeKeys( int depth , int slice , int z , std::vector< SlabValues< Vertex > >& sValues , int threads );
605  template< class Vertex >
606  void CopyFinerXSliceIsoEdgeKeys( int depth , int slab , std::vector< SlabValues< Vertex > >& sValues , int threads );
607  template< class Vertex >
608  void SetSliceIsoEdges( int depth , int slice , std::vector< SlabValues< Vertex > >& slabValues , int threads );
609  template< class Vertex >
610  void SetSliceIsoEdges( int depth , int slice , int z , std::vector< SlabValues< Vertex > >& slabValues , int threads );
611  template< class Vertex >
612  void SetXSliceIsoEdges( int depth , int slice , std::vector< SlabValues< Vertex > >& slabValues , int threads );
613 
614  template< class Vertex >
615  void SetIsoSurface( int depth , int offset , const SliceValues< Vertex >& bValues , const SliceValues< Vertex >& fValues , const XSliceValues< Vertex >& xValues , CoredMeshData< Vertex >& mesh , bool polygonMesh , bool addBarycenter , int& vOffset , int threads );
616 
617  template< class Vertex >
618  static int AddIsoPolygons( CoredMeshData< Vertex >& mesh , std::vector< std::pair< int , Vertex > >& polygon , bool polygonMesh , bool addBarycenter , int& vOffset );
619 
620  template< int WeightDegree , int ColorDegree , class Vertex >
621  bool GetIsoVertex( const BSplineData< ColorDegree >* colorBSData , const SparseNodeData< Real , WeightDegree >* densityWeights , const SparseNodeData< ProjectiveData< Point3D< Real > > , ColorDegree >* colorData , Real isoValue , ConstPointSupportKey< WeightDegree >& weightKey , ConstPointSupportKey< ColorDegree >& colorKey , const TreeOctNode* node , int edgeIndex , int z , const SliceValues< Vertex >& sValues , Vertex& vertex );
622  template< int WeightDegree , int ColorDegree , class Vertex >
623  bool GetIsoVertex( const BSplineData< ColorDegree >* colorBSData , const SparseNodeData< Real , WeightDegree >* densityWeights , const SparseNodeData< ProjectiveData< Point3D< Real > > , ColorDegree >* colorData , Real isoValue , ConstPointSupportKey< WeightDegree >& weightKey , ConstPointSupportKey< ColorDegree >& colorKey , const TreeOctNode* node , int cornerIndex , const SliceValues< Vertex >& bValues , const SliceValues< Vertex >& fValues , Vertex& vertex );
624 
625 public:
626  static double maxMemoryUsage;
627  int threads;
628 
629  static double MemoryUsage( void );
630  Octree( void );
631 
632  // After calling set tree, the indices of the octree node will be stored by depth, and within depth they will be sorted by slice
633  template< class PointReal , int NormalDegree , int WeightDegree , int DataDegree , class Data , class _Data >
634  int SetTree( OrientedPointStream< PointReal >* pointStream , int minDepth , int maxDepth , int fullDepth , int splatDepth , Real samplesPerNode ,
635  Real scaleFactor , bool useConfidence , bool useNormalWeight ,
636  Real constraintWeight , int adaptiveExponent ,
637  SparseNodeData< Real , WeightDegree >& densityWeights , SparseNodeData< PointData< Real > , 0 >& pointInfo , SparseNodeData< Point3D< Real > , NormalDegree >& normalInfo , SparseNodeData< Real , NormalDegree >& nodeWeights ,
638  SparseNodeData< ProjectiveData< _Data > , DataDegree >* dataValues ,
639  XForm4x4< Real >& xForm , bool dirichlet=false , bool makeComplete=false );
640 
641  template< int FEMDegree > void EnableMultigrid( std::vector< int >* map );
642 
643  template< int FEMDegree , int NormalDegree >
645  template< int FEMDegree >
646  DenseNodeData< Real , FEMDegree > SolveSystem( SparseNodeData< PointData< Real > , 0 >& pointInfo , DenseNodeData< Real , FEMDegree >& constraints , bool showResidual , int iters , int maxSolveDepth , int cgDepth=0 , double cgAccuracy=0 );
647 
648  template< int FEMDegree , int NormalDegree >
650  template< int FEMDegree , int WeightDegree , int ColorDegree , class Vertex >
651  void GetMCIsoSurface( const SparseNodeData< Real , WeightDegree >* densityWeights , const SparseNodeData< ProjectiveData< Point3D< Real > > , ColorDegree >* colorData , const DenseNodeData< Real , FEMDegree >& solution , Real isoValue , CoredMeshData< Vertex >& mesh , bool nonLinearFit=true , bool addBarycenter=false , bool polygonMesh=false );
652 
653  const TreeOctNode& tree( void ) const{ return _tree; }
654  size_t leaves( void ) const { return _tree.leaves(); }
655  size_t nodes( void ) const { return _tree.nodes(); }
656 };
657 
658 template< class Real >
659 void Reset( void )
660 {
663 }
664 
665 #include "MultiGridOctreeData.inl"
666 #include "MultiGridOctreeData.SortedTreeNodes.inl"
667 #include "MultiGridOctreeData.WeightedSamples.inl"
668 #include "MultiGridOctreeData.System.inl"
669 #include "MultiGridOctreeData.IsoSurface.inl"
670 #include "MultiGridOctreeData.Evaluation.inl"
671 #endif // MULTI_GRID_OCTREE_DATA_INCLUDED
#define DeletePointer(...)
Definition: Array.h:90
#define ConstPointer(...)
Definition: Array.h:85
#define NullPointer(...)
Definition: Array.h:86
#define Pointer(...)
Definition: Array.h:84
int width
int offset
math::float3 position
int omp_get_num_procs(void)
void Reset(void)
int omp_get_thread_num(void)
#define DIMENSION
Definition: Octree.h:38
boost::geometry::model::polygon< point_xy > polygon
Definition: TreeIso.cpp:37
#define NULL
static int Dimension(int depth)
Definition: BSplineData.h:356
static void InteriorSupportedSpan(int depth, int &begin, int &end)
Definition: BSplineData.h:110
static void InteriorOverlappedSpan(int depth, int &begin, int &end)
Definition: BSplineData.h:281
static const unsigned int EDGES
Definition: MarchingCubes.h:54
static const unsigned int FACES
Definition: MarchingCubes.h:54
static const unsigned int CORNERS
Definition: MarchingCubes.h:54
size_t leaves(void) const
size_t nodes(void) const
void depthAndOffset(int &depth, int offset[3]) const
NodeData nodeData
Definition: Octree.h:61
int depth(void) const
static bool _IsInteriorlySupported(int d, int x, int y, int z)
static void _UpSample(int highDepth, const C *lowCoefficients, C *highCoefficients, bool dirichlet, int threads)
void _DownSample(int highDepth, DenseNodeData< C, FEMDegree > &constraints) const
int _SolveSystemGS(const BSplineData< FEMDegree > &bsData, SparseNodeData< PointData< Real >, 0 > &pointInfo, int depth, DenseNodeData< Real, FEMDegree > &solution, DenseNodeData< Real, FEMDegree > &constraints, DenseNodeData< Real, FEMDegree > &metSolutionConstraints, int iters, bool coarseToFine, bool showResidual=false, double *bNorm2=NULL, double *inRNorm2=NULL, double *outRNorm2=NULL, bool forceSilent=false)
static int _Resolution(int depth)
void _setMultiColorIndices(int start, int end, std::vector< std::vector< int > > &indices) const
void _UpdateConstraintsFromFiner(const typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::ChildIntegrator &childIntegrator, const BSplineData< FEMDegree > &bsData, int highDepth, const DenseNodeData< Real, FEMDegree > &fineSolution, DenseNodeData< Real, FEMDegree > &coarseConstraints) const
void _SetPointConstraintsFromFiner(const SparseNodeData< PointData< Real >, 0 > &pointInfo, int highDepth, const BSplineData< FEMDegree > &bsData, const DenseNodeData< Real, FEMDegree > &finerCoefficients, DenseNodeData< Real, FEMDegree > &metConstraints) const
static bool IsValidNode(const TreeOctNode *node, bool dirichlet)
static void FunctionIndex(const TreeOctNode *node, int idx[3])
static bool _IsInteriorlyOverlapped(const TreeOctNode *node)
static void _DepthAndOffset(const TreeOctNode *node, int &d, int off[3])
static double maxMemoryUsage
void SetXSliceIsoEdges(int depth, int slice, std::vector< SlabValues< Vertex > > &slabValues, int threads)
void GetMCIsoSurface(const SparseNodeData< Real, WeightDegree > *densityWeights, const SparseNodeData< ProjectiveData< Point3D< Real > >, ColorDegree > *colorData, const DenseNodeData< Real, FEMDegree > &solution, Real isoValue, CoredMeshData< Vertex > &mesh, bool nonLinearFit=true, bool addBarycenter=false, bool polygonMesh=false)
void CopyFinerXSliceIsoEdgeKeys(int depth, int slab, std::vector< SlabValues< Vertex > > &sValues, int threads)
static void _StartAndWidth(const TreeOctNode *node, Point3D< Real > &start, Real &width)
static double MemoryUsage(void)
Real _GetSamplesPerNode(const SparseNodeData< Real, WeightDegree > &densityWeights, const TreeOctNode *node, Point3D< Real > position, ConstPointSupportKey< WeightDegree > &weightKey) const
DenseNodeData< Real, FEMDegree > SetLaplacianConstraints(const SparseNodeData< Point3D< Real >, NormalDegree > &normalInfo)
int SetTree(OrientedPointStream< PointReal > *pointStream, int minDepth, int maxDepth, int fullDepth, int splatDepth, Real samplesPerNode, Real scaleFactor, bool useConfidence, bool useNormalWeight, Real constraintWeight, int adaptiveExponent, SparseNodeData< Real, WeightDegree > &densityWeights, SparseNodeData< PointData< Real >, 0 > &pointInfo, SparseNodeData< Point3D< Real >, NormalDegree > &normalInfo, SparseNodeData< Real, NormalDegree > &nodeWeights, SparseNodeData< ProjectiveData< _Data >, DataDegree > *dataValues, XForm4x4< Real > &xForm, bool dirichlet=false, bool makeComplete=false)
void _setFullDepth(int depth)
int _GetMatrixAndUpdateConstraints(const SparseNodeData< PointData< Real >, 0 > &pointInfo, SparseMatrix< Real > &matrix, DenseNodeData< Real, FEMDegree > &constraints, typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::Integrator &integrator, typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::ChildIntegrator &childIntegrator, const BSplineData< FEMDegree > &bsData, int depth, const DenseNodeData< Real, FEMDegree > *metSolution, bool coarseToFine)
bool _constrainValues
static void _SetFullDepth(TreeOctNode *node, int depth)
bool GetIsoVertex(const BSplineData< ColorDegree > *colorBSData, const SparseNodeData< Real, WeightDegree > *densityWeights, const SparseNodeData< ProjectiveData< Point3D< Real > >, ColorDegree > *colorData, Real isoValue, ConstPointSupportKey< WeightDegree > &weightKey, ConstPointSupportKey< ColorDegree > &colorKey, const TreeOctNode *node, int cornerIndex, const SliceValues< Vertex > &bValues, const SliceValues< Vertex > &fValues, Vertex &vertex)
void SetSliceIsoEdges(int depth, int slice, int z, std::vector< SlabValues< Vertex > > &slabValues, int threads)
Real GetIsoValue(const DenseNodeData< Real, FEMDegree > &solution, const SparseNodeData< Real, NormalDegree > &nodeWeights)
void _AddWeightContribution(SparseNodeData< Real, WeightDegree > &densityWeights, TreeOctNode *node, Point3D< Real > position, PointSupportKey< WeightDegree > &weightKey, Real weight=Real(1.0))
Real _GetSamplesPerNode(const SparseNodeData< Real, WeightDegree > &densityWeights, TreeOctNode *node, Point3D< Real > position, PointSupportKey< WeightDegree > &weightKey)
void SetSliceIsoEdges(int depth, int slice, std::vector< SlabValues< Vertex > > &slabValues, int threads)
void EnableMultigrid(std::vector< int > *map)
static bool _IsInteriorlyOverlapped(int d, int x, int y, int z)
void _MultiSplatPointData(const SparseNodeData< Real, WeightDegree > *densityWeights, Point3D< Real > point, V v, SparseNodeData< V, DataDegree > &data, PointSupportKey< WeightDegree > &weightKey, PointSupportKey< DataDegree > &dataKey, int maxDepth, int dim=DIMENSION)
V * Evaluate(const DenseNodeData< V, DataDegree > &coefficients, int &res, Real isoValue=0.f, int depth=-1, bool primal=false)
static TreeOctNode::ConstNeighbors< LeftRadius+RightRadius+1 > & _Neighbors(TreeOctNode::ConstNeighborKey< LeftRadius, RightRadius > &key, int depth)
bool _IsValidNode(const TreeOctNode *node) const
void _UpSample(int highDepth, DenseNodeData< C, FEMDegree > &coefficients) const
V Evaluate(const DenseNodeData< V, DataDegree > &coefficients, Point3D< Real > p, const BSplineData< DataDegree > &bsData) const
void _ClipTree(const SparseNodeData< Point3D< Real >, NormalDegree > &normalInfo)
Real _FinerFunctionValue(Point3D< Real > p, const PointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, const BSplineData< FEMDegree > &bsData, const DenseNodeData< Real, FEMDegree > &coefficients) const
int _HasNormals(TreeOctNode *node, const SparseNodeData< Point3D< Real >, NormalDegree > &normalInfo)
void _GetSampleDepthAndWeight(const SparseNodeData< Real, WeightDegree > &densityWeights, Point3D< Real > position, PointSupportKey< WeightDegree > &weightKey, Real &depth, Real &weight)
TreeOctNode::ConstNeighborKey< 1, 1 > ConstAdjacenctNodeKey
int _SolveSystemCG(const BSplineData< FEMDegree > &bsData, SparseNodeData< PointData< Real >, 0 > &pointInfo, int depth, DenseNodeData< Real, FEMDegree > &solution, DenseNodeData< Real, FEMDegree > &constraints, DenseNodeData< Real, FEMDegree > &metSolutionConstraints, int iters, bool coarseToFine, bool showResidual=false, double *bNorm2=NULL, double *inRNorm2=NULL, double *outRNorm2=NULL, double accuracy=0)
void CopyFinerSliceIsoEdgeKeys(int depth, int slice, int z, std::vector< SlabValues< Vertex > > &sValues, int threads)
int _GetMatrixRowSize(const typename TreeOctNode::Neighbors< BSplineIntegrationData< FEMDegree, FEMDegree >::OverlapSize > &neighbors) const
SortedTreeNodes _sNodes
int _SetMatrixRow(const SparseNodeData< PointData< Real >, 0 > &pointInfo, const typename TreeOctNode::Neighbors< BSplineIntegrationData< FEMDegree, FEMDegree >::OverlapSize > &neighbors, MatrixEntry< Real > *row, int offset, const typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::Integrator &integrator, const Stencil< double, BSplineIntegrationData< FEMDegree, FEMDegree >::OverlapSize > &stencil, const BSplineData< FEMDegree > &bsData) const
void _SplatPointData(TreeOctNode *node, Point3D< Real > point, V v, SparseNodeData< V, DataDegree > &data, PointSupportKey< DataDegree > &dataKey)
void _GetSampleDepthAndWeight(const SparseNodeData< Real, WeightDegree > &densityWeights, TreeOctNode *node, Point3D< Real > position, PointSupportKey< WeightDegree > &weightKey, Real &depth, Real &weight)
void SetSliceIsoCorners(const DenseNodeData< Real, FEMDegree > &solution, const DenseNodeData< Real, FEMDegree > &coarseSolution, Real isoValue, int depth, int slice, int z, std::vector< SlabValues< Vertex > > &sValues, const _Evaluator< FEMDegree > &evaluator, int threads)
V Evaluate(const SparseNodeData< V, DataDegree > &coefficients, Point3D< Real > p, const BSplineData< DataDegree > &bsData) const
void SetXSliceIsoVertices(const BSplineData< ColorDegree > *colorBSData, const SparseNodeData< Real, WeightDegree > *densityWeights, const SparseNodeData< ProjectiveData< Point3D< Real > >, ColorDegree > *colorData, Real isoValue, int depth, int slab, int &vOffset, CoredMeshData< Vertex > &mesh, std::vector< SlabValues< Vertex > > &sValues, int threads)
static int AddIsoPolygons(CoredMeshData< Vertex > &mesh, std::vector< std::pair< int, Vertex > > &polygon, bool polygonMesh, bool addBarycenter, int &vOffset)
void CopyFinerSliceIsoEdgeKeys(int depth, int slice, std::vector< SlabValues< Vertex > > &sValues, int threads)
void SetSliceIsoVertices(const BSplineData< ColorDegree > *colorBSData, const SparseNodeData< Real, WeightDegree > *densityWeights, const SparseNodeData< ProjectiveData< Point3D< Real > >, ColorDegree > *colorData, Real isoValue, int depth, int slice, int &vOffset, CoredMeshData< Vertex > &mesh, std::vector< SlabValues< Vertex > > &sValues, int threads)
const TreeOctNode & tree(void) const
static const TreeOctNode::template Neighbors< LeftRadius+RightRadius+1 > & _Neighbors(const typename TreeOctNode::template NeighborKey< LeftRadius, RightRadius > &key, int depth)
void SetSliceIsoCorners(const DenseNodeData< Real, FEMDegree > &solution, const DenseNodeData< Real, FEMDegree > &coarseSolution, Real isoValue, int depth, int slice, std::vector< SlabValues< Vertex > > &sValues, const _Evaluator< FEMDegree > &evaluator, int threads)
Octree(void)
void SetSliceIsoVertices(const BSplineData< ColorDegree > *colorBSData, const SparseNodeData< Real, WeightDegree > *densityWeights, const SparseNodeData< ProjectiveData< Point3D< Real > >, ColorDegree > *colorData, Real isoValue, int depth, int slice, int z, int &vOffset, CoredMeshData< Vertex > &mesh, std::vector< SlabValues< Vertex > > &sValues, int threads)
V _Evaluate(const DenseNodeData< V, DataDegree > &coefficients, Point3D< Real > p, const BSplineData< DataDegree > &bsData, const ConstPointSupportKey< DataDegree > &neighborKey) const
int _GetSliceMatrixAndUpdateConstraints(const SparseNodeData< PointData< Real >, 0 > &pointInfo, SparseMatrix< Real > &matrix, DenseNodeData< Real, FEMDegree > &constraints, typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::Integrator &integrator, typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::ChildIntegrator &childIntegrator, const BSplineData< FEMDegree > &bsData, int depth, int slice, const DenseNodeData< Real, FEMDegree > &metSolution, bool coarseToFine)
static int _Dimension(int depth)
void _SetValidityFlags(void)
V _getEdgeValue(const ConstPointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, int edge, const DenseNodeData< V, FEMDegree > &solution, const DenseNodeData< V, FEMDegree > &metSolution, const _Evaluator< FEMDegree > &evaluator, bool isInterior) const
void SetIsoSurface(int depth, int offset, const SliceValues< Vertex > &bValues, const SliceValues< Vertex > &fValues, const XSliceValues< Vertex > &xValues, CoredMeshData< Vertex > &mesh, bool polygonMesh, bool addBarycenter, int &vOffset, int threads)
Real _CoarserFunctionValue(Point3D< Real > p, const PointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, const BSplineData< FEMDegree > &bsData, const DenseNodeData< Real, FEMDegree > &upSampledCoefficients) const
size_t leaves(void) const
TreeOctNode * _spaceRoot
static int _Depth(const TreeOctNode *node)
TreeOctNode::NeighborKey< 1, 1 > AdjacenctNodeKey
static void _SetParentOverlapBounds(const TreeOctNode *node, int &startX, int &endX, int &startY, int &endY, int &startZ, int &endZ)
static TreeOctNode::Neighbors< LeftRadius+RightRadius+1 > & _Neighbors(TreeOctNode::NeighborKey< LeftRadius, RightRadius > &key, int depth)
static const int CHILDREN
static bool _IsInteriorlySupported(const TreeOctNode *node)
Real _SplatPointData(const SparseNodeData< Real, WeightDegree > &densityWeights, Point3D< Real > point, V v, SparseNodeData< V, DataDegree > &data, PointSupportKey< WeightDegree > &weightKey, PointSupportKey< DataDegree > &dataKey, int minDepth, int maxDepth, int dim=DIMENSION)
V _getCornerValue(const ConstPointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, int corner, const DenseNodeData< V, FEMDegree > &solution, const DenseNodeData< V, FEMDegree > &metSolution, const _Evaluator< FEMDegree > &evaluator, bool isInterior) const
void _GetSampleDepthAndWeight(const SparseNodeData< Real, WeightDegree > &densityWeights, const TreeOctNode *node, Point3D< Real > position, ConstPointSupportKey< WeightDegree > &weightKey, Real &depth, Real &weight) const
bool GetIsoVertex(const BSplineData< ColorDegree > *colorBSData, const SparseNodeData< Real, WeightDegree > *densityWeights, const SparseNodeData< ProjectiveData< Point3D< Real > >, ColorDegree > *colorData, Real isoValue, ConstPointSupportKey< WeightDegree > &weightKey, ConstPointSupportKey< ColorDegree > &colorKey, const TreeOctNode *node, int edgeIndex, int z, const SliceValues< Vertex > &sValues, Vertex &vertex)
static const TreeOctNode::template ConstNeighbors< LeftRadius+RightRadius+1 > & _Neighbors(const typename TreeOctNode::template ConstNeighborKey< LeftRadius, RightRadius > &key, int depth)
Point3D< Real > _center
std::pair< Real, Point3D< Real > > _getEdgeValueAndGradient(const ConstPointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, int edge, const DenseNodeData< Real, FEMDegree > &solution, const DenseNodeData< Real, FEMDegree > &metSolution, const _Evaluator< FEMDegree > &evaluator, bool isInterior) const
std::pair< Real, Point3D< Real > > _getCornerValueAndGradient(const ConstPointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, int corner, const DenseNodeData< Real, FEMDegree > &solution, const DenseNodeData< Real, FEMDegree > &metSolution, const _Evaluator< FEMDegree > &evaluator, bool isInterior) const
size_t nodes(void) const
void _UpdateConstraintsFromCoarser(const SparseNodeData< PointData< Real >, 0 > &pointInfo, const typename TreeOctNode::Neighbors< BSplineIntegrationData< FEMDegree, FEMDegree >::OverlapSize > &neighbors, const typename TreeOctNode::Neighbors< BSplineIntegrationData< FEMDegree, FEMDegree >::OverlapSize > &pNeighbors, TreeOctNode *node, DenseNodeData< Real, FEMDegree > &constraints, const DenseNodeData< Real, FEMDegree > &metSolution, const typename BSplineIntegrationData< FEMDegree, FEMDegree >::FunctionIntegrator::ChildIntegrator &childIntegrator, const Stencil< double, BSplineIntegrationData< FEMDegree, FEMDegree >::OverlapSize > &stencil, const BSplineData< FEMDegree > &bsData) const
void _SetPointValuesFromCoarser(SparseNodeData< PointData< Real >, 0 > &pointInfo, int highDepth, const BSplineData< FEMDegree > &bsData, const DenseNodeData< Real, FEMDegree > &upSampledCoefficients)
V _Evaluate(const SparseNodeData< V, DataDegree > &coefficients, Point3D< Real > p, const BSplineData< DataDegree > &bsData, const ConstPointSupportKey< DataDegree > &dataKey) const
bool _InBounds(Point3D< Real >) const
static void _CenterAndWidth(const TreeOctNode *node, Point3D< Real > &center, Real &width)
V _getCenterValue(const ConstPointSupportKey< FEMDegree > &neighborKey, const TreeOctNode *node, const DenseNodeData< V, FEMDegree > &solution, const DenseNodeData< V, FEMDegree > &metSolution, const _Evaluator< FEMDegree > &evaluator, bool isInterior) const
void _GetSampleDepthAndWeight(const SparseNodeData< Real, WeightDegree > &densityWeights, Point3D< Real > position, ConstPointSupportKey< WeightDegree > &weightKey, Real &depth, Real &weight)
void _MakeComplete(void)
DenseNodeData< Real, FEMDegree > SolveSystem(SparseNodeData< PointData< Real >, 0 > &pointInfo, DenseNodeData< Real, FEMDegree > &constraints, bool showResidual, int iters, int maxSolveDepth, int cgDepth=0, double cgAccuracy=0)
TreeOctNode _tree
void setXSliceTableData(XSliceTableData &sData, int depth, int offset, int threads) const
int begin(int depth) const
_Indices< Square::CORNERS > SquareCornerIndices
int size(int depth) const
_Indices< Square::EDGES > SquareEdgeIndices
int end(int depth) const
int end(int depth, int slice) const
TreeOctNode ** treeNodes
int size(int depth, int slice) const
void set(TreeOctNode &root, std::vector< int > *map)
_Indices< Square::FACES > SquareFaceIndices
~SortedTreeNodes(void)
int size(void) const
void set(TreeOctNode &root)
int begin(int depth, int slice) const
int levels(void) const
void setSliceTableData(SliceTableData &sData, int depth, int offset, int threads) const
static double GetLaplacian(const typename FunctionIntegrator::Integrator &integrator, const int off1[3], const int off2[3])
static void SetCentralLaplacianStencil(const typename FunctionIntegrator::Integrator &integrator, Stencil< double, OverlapSize > &stencil)
static Point3D< double > GetDivergence1(const typename FunctionIntegrator::ChildIntegrator &integrator, const int off1[3], const int off2[3])
static Point3D< double > GetDivergence1(const typename FunctionIntegrator::Integrator &integrator, const int off1[3], const int off2[3])
static void SetCentralDivergenceStencil(const typename FunctionIntegrator::Integrator &integrator, Stencil< Point3D< double >, OverlapSize > &stencil, bool scatter)
static void SetCentralLaplacianStencils(const typename FunctionIntegrator::ChildIntegrator &integrator, Stencil< double, OverlapSize > stencil[2][2][2])
static double GetDivergence2(const typename FunctionIntegrator::ChildIntegrator &integrator, const int off1[3], const int off2[3], Point3D< double > normal2)
static double GetDivergence1(const typename FunctionIntegrator::ChildIntegrator &integrator, const int off1[3], const int off2[3], Point3D< double > normal1)
static Point3D< double > GetDivergence2(const typename FunctionIntegrator::ChildIntegrator &integrator, const int off1[3], const int off2[3])
static Point3D< double > GetDivergence2(const typename FunctionIntegrator::Integrator &integrator, const int off1[3], const int off2[3])
static double GetDivergence2(const typename FunctionIntegrator::Integrator &integrator, const int off1[3], const int off2[3], Point3D< double > normal2)
static double GetDivergence1(const typename FunctionIntegrator::Integrator &integrator, const int off1[3], const int off2[3], Point3D< double > normal1)
static void SetCentralDivergenceStencils(const typename FunctionIntegrator::ChildIntegrator &integrator, Stencil< Point3D< double >, OverlapSize > stencil[2][2][2], bool scatter)
static double GetLaplacian(const typename FunctionIntegrator::ChildIntegrator &integrator, const int off1[3], const int off2[3])
TreeNodeData(void)
static size_t NodeCount
~TreeNodeData(void)
static long long CenterIndex(const TreeOctNode *node, int maxDepth, int index[DIMENSION])
static long long FaceIndex(const TreeOctNode *node, int fIndex, int maxDepth, int index[DIMENSION])
static long long CenterIndex(int depth, const int offSet[DIMENSION], int maxDepth, int index[DIMENSION])
static long long EdgeIndex(const TreeOctNode *node, int eIndex, int maxDepth, int index[DIMENSION])
static long long CornerIndex(const TreeOctNode *node, int cIndex, int maxDepth, int index[DIMENSION])
static long long CenterIndex(const TreeOctNode *node, int maxDepth)
static const int VERTEX_COORDINATE_SHIFT
static long long CornerIndex(int depth, const int offSet[DIMENSION], int cIndex, int maxDepth, int index[DIMENSION])
static long long EdgeIndex(const TreeOctNode *node, int eIndex, int maxDepth)
static long long CornerIndex(const TreeOctNode *node, int cIndex, int maxDepth)
static long long FaceIndex(const TreeOctNode *node, int fIndex, int maxDepth)
static long long CornerIndexKey(const int index[DIMENSION])
static const int RightRadius
static const int LeftRadius
DenseNodeData(size_t sz)
Data & operator[](int idx)
void resize(size_t sz)
IsoEdge(long long v1, long long v2)
long long & operator[](int idx)
const long long & operator[](int idx) const
ProjectiveData & operator/=(Real s)
ProjectiveData & operator*=(Real s)
ProjectiveData operator/(Real s) const
ProjectiveData operator*(Real s) const
ProjectiveData operator-(const ProjectiveData &p) const
ProjectiveData & operator-=(const ProjectiveData &p)
ProjectiveData(V vv=V(0), Real ww=Real(0))
ProjectiveData operator+(const ProjectiveData &p) const
ProjectiveData & operator+=(const ProjectiveData &p)
const XSliceValues< Vertex > & xSliceValues(int idx) const
SliceValues< Vertex > _sliceValues[2]
XSliceValues< Vertex > & xSliceValues(int idx)
SliceValues< Vertex > & sliceValues(int idx)
XSliceValues< Vertex > _xSliceValues[2]
const SliceValues< Vertex > & sliceValues(int idx) const
void reset(bool nonLinearFit)
SortedTreeNodes::SliceTableData sliceData
Point3D< Real > * cornerGradients
std::unordered_map< long long, std::pair< int, Vertex > > edgeVertexMap
std::unordered_map< long long, long long > vertexPairMap
std::unordered_map< long long, std::vector< IsoEdge > > faceEdgeMap
std::unordered_map< long long, std::pair< int, Vertex > > edgeVertexMap
SortedTreeNodes::XSliceTableData xSliceData
std::unordered_map< long long, std::vector< IsoEdge > > faceEdgeMap
std::unordered_map< long long, long long > vertexPairMap
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dCellStencils[CHILDREN]
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dEdgeStencils[CHILDREN][Cube::EDGES]
void set(int depth, bool dirichlet)
BSplineEvaluationData< FEMDegree >::Evaluator evaluator
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > cornerStencil[Cube::CORNERS]
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dFaceStencil[Cube::FACES]
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > cornerStencils[CHILDREN][Cube::CORNERS]
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dCellStencil
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > faceStencils[CHILDREN][Cube::FACES]
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > edgeStencils[CHILDREN][Cube::EDGES]
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dCornerStencils[CHILDREN][Cube::CORNERS]
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > cellStencils[CHILDREN]
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dEdgeStencil[Cube::EDGES]
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > cellStencil
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > edgeStencil[Cube::EDGES]
BSplineEvaluationData< FEMDegree >::ChildEvaluator childEvaluator
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dFaceStencils[CHILDREN][Cube::FACES]
Stencil< double, BSplineEvaluationData< FEMDegree >::SupportSize > faceStencil[Cube::FACES]
Stencil< Point3D< double >, BSplineEvaluationData< FEMDegree >::SupportSize > dCornerStencil[Cube::CORNERS]
Real weightedCoarserDValue
Point3D< Real > position
PointData(Point3D< Real > p=Point3D< Real >(), Real w=0)
static const int Size
static const int LeftRadius
static const int RightRadius
SquareFaceIndices & faceIndices(int idx)
SquareEdgeIndices & edgeIndices(const TreeOctNode *node)
SquareCornerIndices & cornerIndices(const TreeOctNode *node)
const SquareCornerIndices & cornerIndices(const TreeOctNode *node) const
SquareCornerIndices & cornerIndices(int idx)
const SquareFaceIndices & faceIndices(int idx) const
SquareFaceIndices & faceIndices(const TreeOctNode *node)
SquareEdgeIndices & edgeIndices(int idx)
const SquareEdgeIndices & edgeIndices(int idx) const
const SquareCornerIndices & cornerIndices(int idx) const
const SquareEdgeIndices & edgeIndices(const TreeOctNode *node) const
const SquareFaceIndices & faceIndices(const TreeOctNode *node) const
SquareEdgeIndices & faceIndices(int idx)
SquareEdgeIndices & faceIndices(const TreeOctNode *node)
SquareCornerIndices & edgeIndices(int idx)
const SquareCornerIndices & edgeIndices(const TreeOctNode *node) const
const SquareEdgeIndices & faceIndices(int idx) const
SquareCornerIndices & edgeIndices(const TreeOctNode *node)
const SquareEdgeIndices & faceIndices(const TreeOctNode *node) const
const SquareCornerIndices & edgeIndices(int idx) const
std::vector< Data > data
void resize(size_t sz)
void remapIndices(const std::vector< int > &map)
std::vector< int > indices
int index(const OctNode< TreeNodeData > *node) const
C values[N][N][N]
Definition: lsd.c:149