From d5df7d32732e309b147b231b28caa0fa4672027e Mon Sep 17 00:00:00 2001 From: eap Date: Wed, 3 Mar 2010 07:48:07 +0000 Subject: [PATCH] 0020725: EDF 1242 SMESH : Crash avec Convert lin--> quad avec BLSURF/GHS3D on 64bits 0020721: EDF 1233 SMESH : Crash/bad behavior of 'Convert linear Quadratic with Medium Nodes on Geometry' feature with BLSurf/Ghs3D * Fix GetNodeUV() for the case of surface both U and V periodic. * Protect QFace::GetBoundaryLink() from infinite recursion. * Protect QFace::GetLinkChain() from stack overflow. --- src/SMESH/SMESH_MeshEditor.cxx | 600 +++++++++++++++++++++++++++++---- src/SMESH/SMESH_MeshEditor.hxx | 12 +- 2 files changed, 554 insertions(+), 58 deletions(-) diff --git a/src/SMESH/SMESH_MeshEditor.cxx b/src/SMESH/SMESH_MeshEditor.cxx index 1032898a3..04c644c34 100644 --- a/src/SMESH/SMESH_MeshEditor.cxx +++ b/src/SMESH/SMESH_MeshEditor.cxx @@ -38,11 +38,12 @@ #include "SMESHDS_Group.hxx" #include "SMESHDS_Mesh.hxx" -#include "SMESH_subMesh.hxx" +#include "SMESH_Algo.hxx" #include "SMESH_ControlsDef.hxx" +#include "SMESH_Group.hxx" #include "SMESH_MesherHelper.hxx" #include "SMESH_OctreeNode.hxx" -#include "SMESH_Group.hxx" +#include "SMESH_subMesh.hxx" #include "utilities.h" @@ -51,11 +52,17 @@ #include #include #include +#include #include +#include #include +#include #include #include +#include #include +#include +#include #include #include #include @@ -75,6 +82,7 @@ #include #include #include + #include #include @@ -1102,6 +1110,7 @@ bool SMESH_MeshEditor::QuadToTri (TIDSortedElemSet & theElems, //function : BestSplit //purpose : Find better diagonal for cutting. //======================================================================= + int SMESH_MeshEditor::BestSplit (const SMDS_MeshElement* theQuad, SMESH::Controls::NumericalFunctorPtr theCrit) { @@ -1143,6 +1152,164 @@ int SMESH_MeshEditor::BestSplit (const SMDS_MeshElement* theQuad, return -1; } +namespace +{ + // Methods of splitting volumes into tetra + + const int theHexTo5[5*4] = + { + 0, 1, 5, 2, + 0, 4, 5, 7, + 0, 3, 7, 2, + 5, 6, 7, 2, + 0, 2, 5, 7 + }; + const int theHexTo6[6*4] = + { + 0, 1, 5, 2, + 0, 4, 5, 7, + 0, 3, 7, 2, + 5, 6, 7, 2, + 0, 2, 5, 7 + }; + const int thePyraTo2[2*4] = + { + 0, 1, 2, 4, + 0, 2, 3, 4 + }; + + const int thePentaTo8[8*4] = + { + 0, 1, 2, 6, + 3, 5, 4, 6, + 0, 3, 4, 6, + 0, 4, 1, 6, + 1, 4, 5, 6, + 1, 5, 2, 6, + 2, 5, 3, 6, + 2, 3, 0, 6 + }; + + struct TSplitMethod + { + int _nbTetra; + const int* _connectivity; + bool _addNode; // additional node is to be created + TSplitMethod( int nbTet=0, const int* conn=0, bool addNode=false) + : _nbTetra(nbTet), _connectivity(conn), _addNode(addNode) {} + }; + + /*! + * \brief return TSplitMethod for the given element + */ + TSplitMethod getSplitMethod( const SMDS_MeshElement* vol, const int theMethodFlags) + { + TSplitMethod method; + if ( vol->GetType() == SMDSAbs_Volume && !vol->IsPoly()) + switch ( vol->NbNodes() ) + { + case 8: + case 20: + if ( theMethodFlags & SMESH_MeshEditor::HEXA_TO_5 ) + method = TSplitMethod( 5, theHexTo5 ); + else + method = TSplitMethod( 6, theHexTo6 ); + break; + case 5: + case 13: + method = TSplitMethod( 2, thePyraTo2 ); + break; + case 6: + case 15: + method = TSplitMethod( 8, thePentaTo8, /*addNode=*/true ); + break; + default:; + } + return method; + } +} + +//======================================================================= +//function : SplitVolumesIntoTetra +//purpose : Split volumic elements into tetrahedra. +//======================================================================= + +void SMESH_MeshEditor::SplitVolumesIntoTetra (const TIDSortedElemSet & theElems, + const int theMethodFlags) +{ + // sdt-like iterator on coordinates of nodes of mesh element + typedef SMDS_StdIterator< TNodeXYZ, SMDS_ElemIteratorPtr > NXyzIterator; + NXyzIterator xyzEnd; + + SMESH_MesherHelper helper( *GetMesh()); + + TIDSortedElemSet::const_iterator elem = theElems.begin(); + for ( ; elem != theElems.end(); ++elem ) + { + SMDSAbs_EntityType geomType = (*elem)->GetEntityType(); + if ( geomType <= SMDSEntity_Quad_Tetra ) + continue; // tetra or face or edge + + if ( (*elem)->IsQuadratic() ) + { + // add quadratic links to the helper + SMDS_VolumeTool vol( *elem ); + for ( int iF = 0; iF < vol.NbFaces(); ++iF ) + { + const SMDS_MeshNode** fNodes = vol.GetFaceNodes( iF ); + for ( int iN = 0; iN < vol.NbFaceNodes( iF ); iN += 2) + helper.AddTLinkNode( fNodes[iF], fNodes[iF+2], fNodes[iF+1] ); + } + helper.SetIsQuadratic( true ); + } + else + { + helper.SetIsQuadratic( false ); + } + + vector tetras; // splits of a volume + + if ( geomType == SMDSEntity_Polyhedra ) + { + // Each face of a polyhedron is split into triangles and + // each of triangles and a cell barycenter form a tetrahedron. + + SMDS_VolumeTool vol( *elem ); + + // make a node at barycenter + gp_XYZ gc = std::accumulate( NXyzIterator((*elem)->nodesIterator()), xyzEnd,gp_XYZ(0,0,0)); + gc /= vol.NbNodes(); + SMDS_MeshNode* gcNode = GetMeshDS()->AddNode( gc.X(), gc.Y(), gc.Z() ); + + for ( int iF = 0; iF < vol.NbFaces(); ++iF ) + { + const SMDS_MeshNode** fNodes = vol.GetFaceNodes( iF ); + int nbFNodes = vol.NbFaceNodes( iF ); + int nbTria = nbFNodes - 2; + bool extFace = vol.IsFaceExternal( iF ); + SMDS_MeshElement* tet; + for ( int i = 0; i < nbTria; ++i ) + { + if ( extFace ) + tet = helper.AddVolume( fNodes[0], fNodes[i+1], fNodes[i+2], gcNode ); + else + tet = helper.AddVolume( fNodes[0], fNodes[i+2], fNodes[i+1], gcNode ); + tetras.push_back( tet ); + } + } + + } + else + { + + TSplitMethod splitMethod = getSplitMethod( *elem, theMethodFlags ); + if ( splitMethod._nbTetra < 1 ) continue; + + vector volNodes( (*elem)->begin_nodes(), (*elem)->end_nodes()); + } + } +} + //======================================================================= //function : AddToSameGroups //purpose : add elemToAdd to the groups the elemInGroups belongs to @@ -5584,6 +5751,7 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint() ElementBndBoxTree(const SMDS_Mesh& mesh, SMDSAbs_ElementType elemType); void getElementsNearPoint( const gp_Pnt& point, TIDSortedElemSet& foundElems); + void getElementsNearLine ( const gp_Ax1& line, TIDSortedElemSet& foundElems); ~ElementBndBoxTree(); protected: @@ -5709,6 +5877,31 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint() } } + //================================================================================ + /*! + * \brief Return elements which can be intersected by the line + */ + //================================================================================ + + void ElementBndBoxTree::getElementsNearLine( const gp_Ax1& line, + TIDSortedElemSet& foundElems) + { + if ( level() && getBox().IsOut( line )) + return; + + if ( isLeaf() ) + { + for ( int i = 0; i < _elements.size(); ++i ) + if ( !_elements[i]->IsOut( line )) + foundElems.insert( _elements[i]->_element ); + } + else + { + for (int i = 0; i < 8; i++) + ((ElementBndBoxTree*) myChildren[i])->getElementsNearLine( line, foundElems ); + } + } + //================================================================================ /*! * \brief Construct the element box @@ -5729,60 +5922,78 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint() //======================================================================= /*! - * \brief Implementation of search for the elements by point + * \brief Implementation of search for the elements by point and + * of classification of point in 2D mesh */ //======================================================================= struct SMESH_ElementSearcherImpl: public SMESH_ElementSearcher { - SMESHDS_Mesh* _mesh; - ElementBndBoxTree* _ebbTree; - SMESH_NodeSearcherImpl* _nodeSearcher; - SMDSAbs_ElementType _elementType; - - SMESH_ElementSearcherImpl( SMESHDS_Mesh& mesh ): _mesh(&mesh),_ebbTree(0),_nodeSearcher(0) {} + SMESHDS_Mesh* _mesh; + ElementBndBoxTree* _ebbTree; + SMESH_NodeSearcherImpl* _nodeSearcher; + SMDSAbs_ElementType _elementType; + double _tolerance; + set _internalFaces; + + SMESH_ElementSearcherImpl( SMESHDS_Mesh& mesh ) + : _mesh(&mesh),_ebbTree(0),_nodeSearcher(0), _tolerance(-1) {} ~SMESH_ElementSearcherImpl() { if ( _ebbTree ) delete _ebbTree; _ebbTree = 0; if ( _nodeSearcher ) delete _nodeSearcher; _nodeSearcher = 0; } + virtual int FindElementsByPoint(const gp_Pnt& point, + SMDSAbs_ElementType type, + vector< const SMDS_MeshElement* >& foundElements); + virtual TopAbs_State GetPointState(const gp_Pnt& point); - /*! - * \brief Find elements of given type where the given point is IN or ON. - * Returns nb of found elements and elements them-selves. - * - * 'ALL' type means elements of any type excluding nodes and 0D elements - */ - int FindElementsByPoint(const gp_Pnt& point, - SMDSAbs_ElementType type, - vector< const SMDS_MeshElement* >& foundElements) + struct TInters //!< data of intersection of the line and the mesh face { - foundElements.clear(); + const SMDS_MeshElement* _face; + gp_Vec _faceNorm; + bool _coincides; //!< the line lays in face plane + TInters(const SMDS_MeshElement* face, const gp_Vec& faceNorm, bool coinc=false) + : _face(face), _faceNorm( faceNorm ), _coincides( coinc ) {} + }; + double getTolerance(); + bool getIntersParamOnLine(const gp_Lin& line, const SMDS_MeshElement* face, + const double tolerance, double & param); + void findOuterBoundary(); + bool isOuterBoundary(const SMDS_MeshElement* face) const { return !_internalFaces.count(face);} +}; +//======================================================================= +/*! + * \brief define tolerance for search + */ +//======================================================================= + +double SMESH_ElementSearcherImpl::getTolerance() +{ + if ( _tolerance < 0 ) + { const SMDS_MeshInfo& meshInfo = _mesh->GetMeshInfo(); - // ----------------- - // define tolerance - // ----------------- - double tolerance = 0; + _tolerance = 0; if ( _nodeSearcher && meshInfo.NbNodes() > 1 ) { double boxSize = _nodeSearcher->getTree()->maxSize(); - tolerance = 1e-8 * boxSize/* / meshInfo.NbNodes()*/; + _tolerance = 1e-8 * boxSize/* / meshInfo.NbNodes()*/; } else if ( _ebbTree && meshInfo.NbElements() > 0 ) { double boxSize = _ebbTree->maxSize(); - tolerance = 1e-8 * boxSize/* / meshInfo.NbElements()*/; + _tolerance = 1e-8 * boxSize/* / meshInfo.NbElements()*/; } - if ( tolerance == 0 ) + if ( _tolerance == 0 ) { // define tolerance by size of a most complex element int complexType = SMDSAbs_Volume; while ( complexType > SMDSAbs_All && meshInfo.NbElements( SMDSAbs_ElementType( complexType )) < 1 ) --complexType; - if ( complexType == SMDSAbs_All ) return foundElements.size(); // empty mesh + if ( complexType == SMDSAbs_All ) return 0; // empty mesh double elemSize; if ( complexType == int( SMDSAbs_Node )) @@ -5804,50 +6015,325 @@ struct SMESH_ElementSearcherImpl: public SMESH_ElementSearcher elemSize = max( dist, elemSize ); } } - tolerance = 1e-6 * elemSize; + _tolerance = 1e-6 * elemSize; + } + } + return _tolerance; +} + +//================================================================================ +/*! + * \brief Find intersection of the line and an edge of face and return parameter on line + */ +//================================================================================ + +bool SMESH_ElementSearcherImpl::getIntersParamOnLine(const gp_Lin& line, + const SMDS_MeshElement* face, + const double tol, + double & param) +{ + int nbInts = 0; + param = 0; + + GeomAPI_ExtremaCurveCurve anExtCC; + Handle(Geom_Curve) lineCurve = new Geom_Line( line ); + + int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes(); + for ( int i = 0; i < nbNodes && nbInts < 2; ++i ) + { + GC_MakeSegment edge( SMESH_MeshEditor::TNodeXYZ( face->GetNode( i )), + SMESH_MeshEditor::TNodeXYZ( face->GetNode( (i+1)%nbNodes) )); + anExtCC.Init( lineCurve, edge); + if ( anExtCC.NbExtrema() > 0 && anExtCC.LowerDistance() <= tol) + { + Quantity_Parameter pl, pe; + anExtCC.LowerDistanceParameters( pl, pe ); + param += pl; + if ( ++nbInts == 2 ) + break; } + } + if ( nbInts > 0 ) param /= nbInts; + return nbInts > 0; +} +//================================================================================ +/*! + * \brief Find all faces belonging to the outer boundary of mesh + */ +//================================================================================ - // ================================================================================= - if ( type == SMDSAbs_Node || type == SMDSAbs_0DElement ) +void SMESH_ElementSearcherImpl::findOuterBoundary() +{ + +} + +//======================================================================= +/*! + * \brief Find elements of given type where the given point is IN or ON. + * Returns nb of found elements and elements them-selves. + * + * 'ALL' type means elements of any type excluding nodes and 0D elements + */ +//======================================================================= + +int SMESH_ElementSearcherImpl:: +FindElementsByPoint(const gp_Pnt& point, + SMDSAbs_ElementType type, + vector< const SMDS_MeshElement* >& foundElements) +{ + foundElements.clear(); + + double tolerance = getTolerance(); + + // ================================================================================= + if ( type == SMDSAbs_Node || type == SMDSAbs_0DElement ) + { + if ( !_nodeSearcher ) + _nodeSearcher = new SMESH_NodeSearcherImpl( _mesh ); + + const SMDS_MeshNode* closeNode = _nodeSearcher->FindClosestTo( point ); + if ( !closeNode ) return foundElements.size(); + + if ( point.Distance( SMESH_MeshEditor::TNodeXYZ( closeNode )) > tolerance ) + return foundElements.size(); // to far from any node + + if ( type == SMDSAbs_Node ) + { + foundElements.push_back( closeNode ); + } + else { - if ( !_nodeSearcher ) - _nodeSearcher = new SMESH_NodeSearcherImpl( _mesh ); + SMDS_ElemIteratorPtr elemIt = closeNode->GetInverseElementIterator( SMDSAbs_0DElement ); + while ( elemIt->more() ) + foundElements.push_back( elemIt->next() ); + } + } + // ================================================================================= + else // elements more complex than 0D + { + if ( !_ebbTree || _elementType != type ) + { + if ( _ebbTree ) delete _ebbTree; + _ebbTree = new ElementBndBoxTree( *_mesh, _elementType = type ); + } + TIDSortedElemSet suspectElems; + _ebbTree->getElementsNearPoint( point, suspectElems ); + TIDSortedElemSet::iterator elem = suspectElems.begin(); + for ( ; elem != suspectElems.end(); ++elem ) + if ( !SMESH_MeshEditor::isOut( *elem, point, tolerance )) + foundElements.push_back( *elem ); + } + return foundElements.size(); +} - const SMDS_MeshNode* closeNode = _nodeSearcher->FindClosestTo( point ); - if ( !closeNode ) return foundElements.size(); +//================================================================================ +/*! + * \brief Classify the given point in the closed 2D mesh + */ +//================================================================================ - if ( point.Distance( SMESH_MeshEditor::TNodeXYZ( closeNode )) > tolerance ) - return foundElements.size(); // to far from any node +TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point) +{ + double tolerance = getTolerance(); + if ( !_ebbTree || _elementType != SMDSAbs_Face ) + { + if ( _ebbTree ) delete _ebbTree; + _ebbTree = new ElementBndBoxTree( *_mesh, _elementType = SMDSAbs_Face ); + } + // algo: analyse transition of a line starting at the point through mesh boundary; + // try several lines, if none of attemps gives a clear answer, we give up as the + // task can be too complex including internal boundaries, concave surfaces etc. + + const int nbAxes = 3; + gp_Dir axisDir[ nbAxes ] = { gp::DX(), gp::DY(), gp::DZ() }; + map< double, TInters > paramOnLine2TInters[ nbAxes ]; + list< TInters > tangentInters[ nbAxes ]; // of faces whose plane includes the line + multimap< int, int > nbInt2Axis; // to find the simplest case + for ( int axis = 0; axis < nbAxes; ++axis ) + { + gp_Ax1 lineAxis( point, axisDir[axis]); + gp_Lin line ( lineAxis ); + + TIDSortedElemSet suspectFaces; // faces possibly intersecting the line + _ebbTree->getElementsNearLine( lineAxis, suspectFaces ); - if ( type == SMDSAbs_Node ) + // Intersect faces with the line + + map< double, TInters > & u2inters = paramOnLine2TInters[ axis ]; + TIDSortedElemSet::iterator face = suspectFaces.begin(); + for ( ; face != suspectFaces.end(); ++face ) + { + // get face plane + gp_XYZ fNorm; + if ( !SMESH_Algo::FaceNormal( *face, fNorm, /*normalized=*/false)) continue; + gp_Pln facePlane( SMESH_MeshEditor::TNodeXYZ( (*face)->GetNode(0)), fNorm ); + + // intersection + IntAna_IntConicQuad intersection( line, IntAna_Quadric( facePlane )); + if ( !intersection.IsDone() ) + continue; + if ( intersection.IsInQuadric() ) { - foundElements.push_back( closeNode ); + tangentInters[ axis ].push_back( TInters( *face, fNorm, true )); } - else + else if ( ! intersection.IsParallel() && intersection.NbPoints() > 0 ) { - SMDS_ElemIteratorPtr elemIt = closeNode->GetInverseElementIterator( SMDSAbs_0DElement ); - while ( elemIt->more() ) - foundElements.push_back( elemIt->next() ); + gp_Pnt intersectionPoint = intersection.Point(1); + if ( !SMESH_MeshEditor::isOut( *face, intersectionPoint, tolerance )) + u2inters.insert(make_pair( intersection.ParamOnConic(1), TInters( *face, fNorm ))); } } - // ================================================================================= - else // elements more complex than 0D + // Analyse intersections roughly + + int nbInter = u2inters.size(); + if ( nbInter == 0 ) + return TopAbs_OUT; + + double f = u2inters.begin()->first, l = u2inters.rbegin()->first; + if ( nbInter == 1 ) + return fabs( f ) < tolerance ? TopAbs_ON : TopAbs_UNKNOWN; + + if ( fabs( f ) < tolerance || fabs( l ) < tolerance ) + return TopAbs_ON; + + if ( (f<0) == (l<0) ) + return TopAbs_OUT; + + int nbIntBeforePoint = std::distance( u2inters.begin(), u2inters.lower_bound(0)); + int nbIntAfterPoint = nbInter - nbIntBeforePoint; + if ( nbIntBeforePoint == 1 || nbIntAfterPoint == 1 ) + return TopAbs_IN; + + nbInt2Axis.insert( make_pair( min( nbIntBeforePoint, nbIntAfterPoint ), axis )); + + } // three attempts - loop on CS axes + + // Analyse intersections thoroughly + // We make two loops, on the first one we correctly exclude touching intersections, + // on the second, we additionally just throw away intersections with small angles + + for ( int angleCheck = 0; angleCheck < 2; ++angleCheck ) + { + multimap< int, int >::const_iterator nb_axis = nbInt2Axis.begin(); + for ( ; nb_axis != nbInt2Axis.end(); ++nb_axis ) { - if ( !_ebbTree || _elementType != type ) + int axis = nb_axis->second; + map< double, TInters > & u2inters = paramOnLine2TInters[ axis ]; + + gp_Ax1 lineAxis( point, axisDir[axis]); + gp_Lin line ( lineAxis ); + + // add tangent intersections to u2inters + double param; + list< TInters >::const_iterator tgtInt = tangentInters[ axis ].begin(); + for ( ; tgtInt != tangentInters[ axis ].end(); ++tgtInt ) + if ( getIntersParamOnLine( line, tgtInt->_face, tolerance, param )) + u2inters.insert(make_pair( param, *tgtInt )); + tangentInters[ axis ].clear(); + + // Count intersections before and after the point excluding touching ones. + + int nbIntBeforePoint = 0, nbIntAfterPoint = 0; + double f = numeric_limits::max(), l = -numeric_limits::max(); + map< double, TInters >::iterator u_int2 = u2inters.begin(), u_int1 = u_int2++; + bool ok = ! u_int1->second._coincides; + while ( ok && u_int1 != u2inters.end() ) + { + // skip intersections at the same point (if line pass through edge or node) + int nbSamePnt = 0; + double u = u_int1->first; + while ( u_int2 != u2inters.end() && fabs( u_int2->first - u ) < tolerance ) + { + ++nbSamePnt; + ++u_int2; + } + + // skip tangent intersections + int nbTgt = 0; + const SMDS_MeshElement* prevFace = u_int1->second._face; + while ( ok && u_int2->second._coincides ) + { + if ( SMESH_Algo::GetCommonNodes(prevFace , u_int2->second._face).empty() ) + ok = false; + else + { + nbTgt++; + u_int2++; + ok = ( u_int2 != u2inters.end() ); + } + } + if ( !ok ) break; + + // skip intersections at the same point after tangent intersections + if ( nbTgt > 0 ) + { + double u = u_int2->first; + ++u_int2; + while ( u_int2 != u2inters.end() && fabs( u_int2->first - u ) < tolerance ) + { + ++nbSamePnt; + ++u_int2; + } + } + + bool touchingInt = false; + if ( nbSamePnt + nbTgt > 0 ) + { + double minDot = numeric_limits::max(), maxDot = -numeric_limits::max(); + map< double, TInters >::iterator u_int = u_int1; + for ( ; u_int != u_int2; ++u_int ) + { + if ( u_int->second._coincides ) continue; + double dot = u_int->second._faceNorm * line.Direction(); + if ( dot > maxDot ) maxDot = dot; + if ( dot < minDot ) minDot = dot; + } + touchingInt = ( minDot*maxDot < 0 ); + } + // throw away intersection with lower angles + if ( !touchingInt && angleCheck ) + { + const double angTol = 2 * Standard_PI180, normAng = Standard_PI / 2; + double angle = u_int1->second._faceNorm.Angle( line.Direction() ); + touchingInt = ( fabs( angle - normAng ) < angTol ); + } + if ( !touchingInt ) + { + if ( u < 0 ) + ++nbIntBeforePoint; + else + ++nbIntAfterPoint; + + if ( u < f ) f = u; + if ( u > l ) l = u; + } + + u_int1 = u_int2++; // to next intersection + + } // loop on intersections with one line + + if ( ok ) { - if ( _ebbTree ) delete _ebbTree; - _ebbTree = new ElementBndBoxTree( *_mesh, _elementType = type ); + if ( nbIntBeforePoint == 0 || nbIntAfterPoint == 0) + return TopAbs_OUT; + + if ( nbIntBeforePoint + nbIntAfterPoint == 1 ) + return fabs( f ) < tolerance ? TopAbs_ON : TopAbs_UNKNOWN; + + if ( nbIntBeforePoint == 1 || nbIntAfterPoint == 1 ) + return fabs( f ) < tolerance ? TopAbs_ON : TopAbs_UNKNOWN; + + if ( fabs( f ) < tolerance || fabs( l ) < tolerance ) + return TopAbs_ON; + + if ( (f<0) == (l<0) ) + return TopAbs_OUT; } - TIDSortedElemSet suspectElems; - _ebbTree->getElementsNearPoint( point, suspectElems ); - TIDSortedElemSet::iterator elem = suspectElems.begin(); - for ( ; elem != suspectElems.end(); ++elem ) - if ( !SMESH_MeshEditor::isOut( *elem, point, tolerance )) - foundElements.push_back( *elem ); - } - return foundElements.size(); - } -}; // struct SMESH_ElementSearcherImpl + } // loop on intersections of the tree lines - thorough analysis + } // two attempts - with and w/o angleCheck + + return TopAbs_UNKNOWN; +} //======================================================================= /*! diff --git a/src/SMESH/SMESH_MeshEditor.hxx b/src/SMESH/SMESH_MeshEditor.hxx index db3c262e1..cd6041943 100644 --- a/src/SMESH/SMESH_MeshEditor.hxx +++ b/src/SMESH/SMESH_MeshEditor.hxx @@ -78,6 +78,7 @@ struct SMESH_NodeSearcher /*! * \brief Find elements of given type where the given point is IN or ON. * Returns nb of found elements and elements them-selves. + * Another task is to find out if the given point is out of closed 2D mesh. * * 'ALL' type means elements of any type excluding nodes and 0D elements */ @@ -88,6 +89,8 @@ struct SMESH_ElementSearcher virtual int FindElementsByPoint(const gp_Pnt& point, SMDSAbs_ElementType type, std::vector< const SMDS_MeshElement* >& foundElems)=0; + + virtual TopAbs_State GetPointState(const gp_Pnt& point) = 0; }; //======================================================================= @@ -124,7 +127,7 @@ public: struct TNodeXYZ : public gp_XYZ { const SMDS_MeshNode* _node; - TNodeXYZ( const SMDS_MeshElement* e):_node(0) { + TNodeXYZ( const SMDS_MeshElement* e):gp_XYZ(0,0,0),_node(0) { if (e) { ASSERT( e->GetType() == SMDSAbs_Node ); _node = static_cast(e); @@ -221,6 +224,13 @@ public: SMESH::Controls::NumericalFunctorPtr theCriterion); + enum SplitVolumToTetraFlags { HEXA_TO_5 = 1, HEXA_TO_6 = 2 };//!