mkazhdan / PoissonRecon

Poisson Surface Reconstruction
MIT License
1.55k stars 424 forks source link

Some Questions About the Code of Screened Poisson Surface Reconstruction #43

Open Jacklwln opened 6 years ago

Jacklwln commented 6 years ago

Hello Misha, There are some questions about the code of "Screened Poisson Surface Reconstruction " . They are: (1) Why do you normalize the input data as mapping the data into the 0-1 range when you construct Octree? (2) I have found that the functions SetSliceIsoCorners(...) , SetSliceIsoVertices(...) and so on are executed on the slices that located between 0-0.5 range at all depth, but not 0-1. I mean corners and iso-vertices are gotten by going through a part of slices, but not all slices, and just process the nodes whose returned value of function _IsValidNode< 0 >( ...) is true, which makes me confused. Could you please tell me why set node flags? And why do not go through all slices or node to get corners and iso-vertices? (3) The reconstructed surfaces are different when using Dirichlet and Neumann boundary conditions, why does Neumann boundary condition could result in an open surface, but the Dirichlet boundary condition leads to a closed surface? And I have other two questions that are not related to the code, what are the differences for constructing surface between triangular mesh and quadrilateral mesh? And quadrilateral mesh could be a wider application?

Cheers, Jack

mkazhdan commented 6 years ago
  1. It easier to define an octree if it partitions a fixed/canonical (e.g. [0,1]^3) cube. 

  2. I'm not sure about this question. However, there are two reasons not to go through all nodes: -- First: If you are using odd-degree elements, then you need to pad out the nodes for indexing reasons (coefficients are stored at cell corners instead of centers). This means that you effectively have an octree over the domain [0,2]^3 for indexing but only use a subset of that for representing implicit functions. -- Second: For reconstructions at depth d, input point data is stored in nodes of depth d, to even if those are not used for representing the implicit function (e.g. if the sampling density is too low).

  3. Dirichlet constraints force the implicit function to have a value of 0 on the boundary of the [0,1]^3 cube. Since we extract the iso-surface a t a value of 0.5, this implies that the implicit surface cannot pass through the boundary of the cube, and hence has to be watertight. 

On December 14, 2017 8:56:00 PM EST, Jack notifications@github.com wrote:

Hello Misha, There are some questions about the code of "Screened Poisson Surface Reconstruction " . They are: (1) Why do you normalize the input data as mapping the data into the 0~1 range when you construct Octree? (2) I have found that the functions SetSliceIsoCorners(...) , SetSliceIsoVertices(...) and so on are executed on the slices that located between 0~0.5 range at all depth, but not 0~1. I mean corners and iso-vertices are gotten by going through a part of slices, but not all slices, and just process the nodes whose returned value of function _IsValidNode< 0 >( ...) is true, which makes me confused. Could you please tell me why set node flags? And why do not go through all slices or node to get corners and iso-vertices? (3) The reconstructed surfaces are different when using Dirichlet and Neumann boundary conditions, why does Neumann boundary condition could result in an open surface, but the Dirichlet boundary condition leads to a closed surface? And I have other two questions that are not related to the code, what are the differences for constructing surface between triangular mesh and quadrilateral mesh? And quadrilateral mesh could be a wider application?

Cheers, Jack

-- You are receiving this because you are subscribed to this thread. Reply to this email directly or view it on GitHub: https://github.com/mkazhdan/PoissonRecon/issues/43

Jacklwln commented 6 years ago

Thank you very much! However, the functions as follow: template< class Real > template< int FEMDegree > bool Octree< Real >::IsValidNode( const TreeOctNode* node , bool dirichlet ) { if( !node || node->depth()<1 ) return false; int d , off[3]; node->depthAndOffset( d , off ); int dim = BSplineData< FEMDegree >::Dimension( d-1 ); if (FEMDegree & 1 && dirichlet) return !(off[0] <= 0 || off[1] <= 0 || off[2] <= 0 || off[0] >= dim - 1 || off[1] >= dim - 1 || off[2] >= dim - 1); else return !( off[0]< 0 || off[1]< 0 || off[2]< 0 || off[0]> dim-1 || off[1]> dim-1 || off[2]> dim-1 ); }

template< class Real > void Octree< Real >::_SetValidityFlags( void ) { for( int i=0 ; i<_sNodes.end( _sNodes.levels()-1 ) ; i++ ) { _sNodes.treeNodes[i]->nodeData.flags = 0; if( IsValidNode< 0 >( _sNodes.treeNodes[i] , _dirichlet ) ) _sNodes.treeNodes[i]->nodeData.flags |= (1<<0); if( IsValidNode< 1 >( _sNodes.treeNodes[i] , _dirichlet ) ) _sNodes.treeNodes[i]->nodeData.flags |= (1<<1); } }

and

template< int FEMDegree > bool _IsValidNode(const TreeOctNode* node) const { return node && (node->nodeData.flags & (1 << (FEMDegree & 1))); } at last, if the offset of a node is bigger than "dim-1", the function __IsValidNode(...) as above could set the node is invalid, even if the node exists. The functions above make node invalid, which confuses me, I do not know why a node would be set invalid when the offset of the node is bigger than "dim-1", could you please tell why? And I believe this is one reason not to go through all nodes.

mkazhdan commented 6 years ago

In that case the node is only used for indexing a basis function, but the interior of the node is outside the unit cube.

On December 15, 2017 9:50:17 AM EST, Jack notifications@github.com wrote:

Thank you very much! However, the functions as follow: template< class Real > template< int FEMDegree > bool Octree< Real >::IsValidNode( const TreeOctNode* node , bool dirichlet ) { if( !node || node->depth()<1 ) return false; int d , off[3]; node->depthAndOffset( d , off ); int dim = BSplineData< FEMDegree >::Dimension( d-1 ); if (FEMDegree & 1 && dirichlet) return !(off[0] <= 0 || off[1] <= 0 || off[2] <= 0 || off[0] >= dim - 1 || off[1] >= dim - 1 || off[2] >= dim - 1); else return !( off[0]< 0 || off[1]< 0 || off[2]< 0 || off[0]> dim-1 || off[1]> dim-1 || off[2]> dim-1 ); }

template< class Real > void Octree< Real >::_SetValidityFlags( void ) { for( int i=0 ; i<_sNodes.end( _sNodes.levels()-1 ) ; i++ ) { _sNodes.treeNodes[i]->nodeData.flags = 0; if( IsValidNode< 0 >( _sNodes.treeNodes[i] , _dirichlet ) ) _sNodes.treeNodes[i]->nodeData.flags |= (1<<0); if( IsValidNode< 1 >( _sNodes.treeNodes[i] , _dirichlet ) ) _sNodes.treeNodes[i]->nodeData.flags |= (1<<1); } }

and

template< int FEMDegree > bool _IsValidNode(const TreeOctNode* node) const { return node && (node->nodeData.flags & (1 << (FEMDegree & 1))); } at last, if the offset of a node is bigger than "dim-1", the function __IsValidNode(...) as above could set the node is invalid, even if the node exists. The functions above make node invalid, which confuses me, I do not know why a node would be set invalid when the offset of the node is bigger than "dim-1", could you please tell why? And I believe this is one reason not to go through all nodes.

-- You are receiving this because you commented. Reply to this email directly or view it on GitHub: https://github.com/mkazhdan/PoissonRecon/issues/43#issuecomment-352024091

Jacklwln commented 6 years ago

Thank you for your reply.