// ========================================================
-// class : SortableElement
-// purpose : allow sorting elements basing on their nodes
+// class : ComparableElement
+// purpose : allow comparing elements basing on their nodes
// ========================================================
-class SortableElement : public set <const SMDS_MeshElement*>
+
+class ComparableElement : public boost::container::flat_set< int >
{
+ typedef boost::container::flat_set< int > int_set;
+
+ const SMDS_MeshElement* myElem;
+ int mySumID;
+ mutable int myGroupID;
+
public:
- SortableElement( const SMDS_MeshElement* theElem )
+ ComparableElement( const SMDS_MeshElement* theElem ):
+ myElem ( theElem ), mySumID( 0 ), myGroupID( -1 )
{
- myElem = theElem;
- SMDS_ElemIteratorPtr nodeIt = theElem->nodesIterator();
- while ( nodeIt->more() )
- this->insert( nodeIt->next() );
+ this->reserve( theElem->NbNodes() );
+ for ( SMDS_ElemIteratorPtr nodeIt = theElem->nodesIterator(); nodeIt->more(); )
+ {
+ int id = nodeIt->next()->GetID();
+ mySumID += id;
+ this->insert( id );
+ }
}
- const SMDS_MeshElement* Get() const
- { return myElem; }
+ const SMDS_MeshElement* GetElem() const { return myElem; }
+
+ int& GroupID() const { return myGroupID; }
+ //int& GroupID() const { return const_cast< int& >( myGroupID ); }
+
+ ComparableElement( const ComparableElement& theSource ) // move copy
+ {
+ ComparableElement& src = const_cast< ComparableElement& >( theSource );
+ (int_set&) (*this ) = boost::move( src );
+ myElem = src.myElem;
+ mySumID = src.mySumID;
+ myGroupID = src.myGroupID;
+ }
+
+ static int HashCode(const ComparableElement& se, int limit )
+ {
+ return ::HashCode( se.mySumID, limit );
+ }
+ static Standard_Boolean IsEqual(const ComparableElement& se1, const ComparableElement& se2 )
+ {
+ return ( se1 == se2 );
+ }
-private:
- mutable const SMDS_MeshElement* myElem;
};
//=======================================================================
// Search among theElements or in the whole mesh if theElements is empty
//=======================================================================
-void SMESH_MeshEditor::FindEqualElements(TIDSortedElemSet & theElements,
- TListOfListOfElementsID & theGroupsOfElementsID)
+void SMESH_MeshEditor::FindEqualElements( TIDSortedElemSet & theElements,
+ TListOfListOfElementsID & theGroupsOfElementsID )
{
ClearLastCreated();
- typedef map< SortableElement, int > TMapOfNodeSet;
- typedef list<int> TGroupOfElems;
-
SMDS_ElemIteratorPtr elemIt;
if ( theElements.empty() ) elemIt = GetMeshDS()->elementsIterator();
else elemIt = SMESHUtils::elemSetIterator( theElements );
- vector< TGroupOfElems > arrayOfGroups;
- TGroupOfElems groupOfElems;
- TMapOfNodeSet mapOfNodeSet;
+ typedef NCollection_Map< ComparableElement, ComparableElement > TMapOfElements;
+ typedef std::list<int> TGroupOfElems;
+ TMapOfElements mapOfElements;
+ std::vector< TGroupOfElems > arrayOfGroups;
+ TGroupOfElems groupOfElems;
- for ( int iGroup = 0; elemIt->more(); )
+ while ( elemIt->more() )
{
const SMDS_MeshElement* curElem = elemIt->next();
- SortableElement SE(curElem);
+ ComparableElement compElem = curElem;
// check uniqueness
- pair< TMapOfNodeSet::iterator, bool> pp = mapOfNodeSet.insert(make_pair(SE, iGroup));
- if ( !pp.second ) { // one more coincident elem
- TMapOfNodeSet::iterator& itSE = pp.first;
- int iG = itSE->second;
+ const ComparableElement& elemInSet = mapOfElements.Added( compElem );
+ if ( elemInSet.GetElem() != curElem ) // coincident elem
+ {
+ int& iG = elemInSet.GroupID();
+ if ( iG < 0 )
+ {
+ iG = arrayOfGroups.size();
+ arrayOfGroups.push_back( groupOfElems );
+ arrayOfGroups[ iG ].push_back( elemInSet.GetElem()->GetID() );
+ }
arrayOfGroups[ iG ].push_back( curElem->GetID() );
}
- else {
- arrayOfGroups.push_back( groupOfElems );
- arrayOfGroups.back().push_back( curElem->GetID() );
- iGroup++;
- }
}
groupOfElems.clear();
- vector< TGroupOfElems >::iterator groupIt = arrayOfGroups.begin();
+ std::vector< TGroupOfElems >::iterator groupIt = arrayOfGroups.begin();
for ( ; groupIt != arrayOfGroups.end(); ++groupIt )
{
if ( groupIt->size() > 1 ) {
- //groupOfElems.sort(); -- theElements is sorted already
- theGroupsOfElementsID.push_back( groupOfElems );
- theGroupsOfElementsID.back().splice( theGroupsOfElementsID.back().end(), *groupIt );
+ //groupOfElems.sort(); -- theElements are sorted already
+ theGroupsOfElementsID.emplace_back( *groupIt );
}
}
}
TIDSortedElemSet aMeshElements; /* empty input ==
to merge equal elements in the whole mesh */
TListOfListOfElementsID aGroupsOfElementsID;
- FindEqualElements(aMeshElements, aGroupsOfElementsID);
- MergeElements(aGroupsOfElementsID);
+ FindEqualElements( aMeshElements, aGroupsOfElementsID );
+ MergeElements( aGroupsOfElementsID );
}
//=======================================================================
myOctreeNode->NodesAround( thePnt.Coord(), dist2Nodes, myHalfLeafSize );
if ( !dist2Nodes.empty() )
return dist2Nodes.begin()->second;
- std::list<const SMDS_MeshNode*> nodes;
+
+ std::vector<const SMDS_MeshNode*> nodes;
//myOctreeNode->NodesAround( &tgtNode, &nodes, myHalfLeafSize );
double minSqDist = DBL_MAX;
if ( nodes.empty() ) // get all nodes of OctreeNode's closest to thePnt
{
// sort leafs by their distance from thePnt
- typedef std::map< double, SMESH_OctreeNode* > TDistTreeMap;
+ typedef std::multimap< double, SMESH_OctreeNode* > TDistTreeMap;
TDistTreeMap treeMap;
std::list< SMESH_OctreeNode* > treeList;
std::list< SMESH_OctreeNode* >::iterator trIt;
{
const Bnd_B3d& box = *tree->getBox();
double sqDist = thePnt.SquareDistance( 0.5 * ( box.CornerMin() + box.CornerMax() ));
- std::pair<TDistTreeMap::iterator,bool> it_in =
- treeMap.insert( std::make_pair( sqDist, tree ));
- if ( !it_in.second ) // not unique distance to box center
- treeMap.insert( it_in.first, std::make_pair( sqDist + 1e-13*treeMap.size(), tree ));
+ treeMap.insert( std::make_pair( sqDist, tree ));
}
}
// find distance after which there is no sense to check tree's
if ( sqDist_tree->first > sqLimit )
break;
SMESH_OctreeNode* tree = sqDist_tree->second;
- tree->NodesAround( tree->GetNodeIterator()->next(), &nodes );
+ tree->AllNodesAround( tree->GetNodeIterator()->next(), &nodes );
}
}
// find closest among nodes
minSqDist = DBL_MAX;
const SMDS_MeshNode* closestNode = 0;
- std::list<const SMDS_MeshNode*>::iterator nIt = nodes.begin();
- for ( ; nIt != nodes.end(); ++nIt ) {
- double sqDist = thePnt.SquareDistance( SMESH_TNodeXYZ( *nIt ) );
+ for ( size_t i = 0; i < nodes.size(); ++i )
+ {
+ double sqDist = thePnt.SquareDistance( SMESH_NodeXYZ( nodes[ i ]));
if ( minSqDist > sqDist ) {
- closestNode = *nIt;
+ closestNode = nodes[ i ];
minSqDist = sqDist;
}
}
//---------------------------------------------------------------------
/*!
- * \brief Finds nodes located within a tolerance near a point
+ * \brief Finds nodes located within a tolerance near a point
*/
int FindNearPoint(const gp_Pnt& point,
const double tolerance,
{
public:
- typedef boost::container::flat_set< const SMDS_MeshElement* > TElemSeq;
+ typedef boost::container::flat_set< const SMDS_MeshElement*, TIDCompare > TElemSeq;
ElementBndBoxTree(const SMDS_Mesh& mesh,
SMDSAbs_ElementType elemType,
#include "SMESH_OctreeNode.hxx"
#include "SMDS_SetIterator.hxx"
+#include "SMESH_MeshAlgos.hxx"
#include "SMESH_TypeDefs.hxx"
+
#include <gp_Pnt.hxx>
using namespace std;
SMESH_OctreeNode::SMESH_OctreeNode (const TIDSortedNodeSet & theNodes, const int maxLevel,
const int maxNbNodes , const double minBoxSize )
:SMESH_Octree( new Limit( maxLevel,minBoxSize,maxNbNodes)),
- myNodes(theNodes)
+ myNodes( theNodes.begin(), theNodes.end() )
{
compute();
}
Bnd_B3d* SMESH_OctreeNode::buildRootBox()
{
Bnd_B3d* box = new Bnd_B3d;
- TIDSortedNodeSet::iterator it = myNodes.begin();
- for (; it != myNodes.end(); it++) {
- const SMDS_MeshNode* n1 = *it;
- gp_XYZ p1( n1->X(), n1->Y(), n1->Z() );
- box->Add(p1);
- }
+ for ( size_t i = 0; i < myNodes.size(); ++i )
+ box->Add( SMESH_NodeXYZ( myNodes[ i ]));
+
if ((int) myNodes.size() <= getMaxNbNodes() )
myIsLeaf = true;
*/
//====================================================================================
-const bool SMESH_OctreeNode::isInside (const gp_XYZ& p, const double precision)
+const bool SMESH_OctreeNode::isInside ( const gp_XYZ& p, const double precision )
{
- if (precision <= 0.)
- return !(getBox()->IsOut(p));
+ if ( precision <= 0.)
+ return !( getBox()->IsOut(p) );
+
Bnd_B3d BoxWithPrecision = *getBox();
- BoxWithPrecision.Enlarge(precision);
+ BoxWithPrecision.Enlarge( precision );
return ! BoxWithPrecision.IsOut(p);
}
* Shares the father's data with each of his child
*/
//================================================
+
void SMESH_OctreeNode::buildChildrenData()
{
gp_XYZ min = getBox()->CornerMin();
gp_XYZ max = getBox()->CornerMax();
gp_XYZ mid = (min + max)/2.;
- TIDSortedNodeSet::iterator it = myNodes.begin();
- while (it != myNodes.end())
+ for ( int i = 0; i < 8; i++ )
{
- const SMDS_MeshNode* n1 = *it;
- int ChildBoxNum = getChildIndex( n1->X(), n1->Y(), n1->Z(), mid );
- SMESH_OctreeNode* myChild = dynamic_cast<SMESH_OctreeNode*> (myChildren[ChildBoxNum]);
- myChild->myNodes.insert(myChild->myNodes.end(),n1);
- myNodes.erase( it );
- it = myNodes.begin();
+ SMESH_OctreeNode* myChild = static_cast<SMESH_OctreeNode*>( myChildren[ i ]);
+ myChild->myNodes.reserve( myNodes.size() / 8 );
}
- for (int i = 0; i < 8; i++)
+
+ for ( size_t i = 0; i < myNodes.size(); ++i )
{
- SMESH_OctreeNode* myChild = dynamic_cast<SMESH_OctreeNode*> (myChildren[i]);
+ SMESH_NodeXYZ n = myNodes[ i ];
+ int ChildBoxNum = getChildIndex( n.X(), n.Y(), n.Z(), mid );
+ SMESH_OctreeNode* myChild = static_cast<SMESH_OctreeNode*>( myChildren[ ChildBoxNum ]);
+ myChild->myNodes.push_back( myNodes[ i ]);
+ }
+ SMESHUtils::FreeVector( myNodes );
+
+ for ( int i = 0; i < 8; i++ )
+ {
+ SMESH_OctreeNode* myChild = static_cast<SMESH_OctreeNode*>( myChildren[ i ]);
if ((int) myChild->myNodes.size() <= getMaxNbNodes() )
+ {
myChild->myIsLeaf = true;
+ if ( myChild->myNodes.empty() )
+ SMESHUtils::FreeVector( myChild->myNodes );
+ }
}
}
* \param Result - list of Nodes potentials to be near Node
*/
//====================================================================
-void SMESH_OctreeNode::NodesAround (const SMDS_MeshNode * Node,
- list<const SMDS_MeshNode*>* Result,
- const double precision)
+
+void SMESH_OctreeNode::AllNodesAround (const SMDS_MeshNode * Node,
+ std::vector<const SMDS_MeshNode*>* Result,
+ const double precision)
{
- SMESH_TNodeXYZ p(Node);
- if (isInside(p, precision))
+ SMESH_NodeXYZ p = Node;
+ if ( isInside( p, precision ))
{
- if (isLeaf())
+ if ( isLeaf() )
{
- Result->insert(Result->end(), myNodes.begin(), myNodes.end());
+ Result->insert( Result->end(), myNodes.begin(), myNodes.end() );
}
else
{
- for (int i = 0; i < 8; i++)
+ for ( int i = 0; i < 8; i++ )
{
- SMESH_OctreeNode* myChild = dynamic_cast<SMESH_OctreeNode*> (myChildren[i]);
- myChild->NodesAround(Node, Result, precision);
+ SMESH_OctreeNode* myChild = static_cast<SMESH_OctreeNode*> (myChildren[i]);
+ myChild->AllNodesAround( Node, Result, precision );
}
}
}
*/
//================================================================================
-bool SMESH_OctreeNode::NodesAround(const gp_XYZ &node,
+bool SMESH_OctreeNode::NodesAround(const gp_XYZ & node,
map<double, const SMDS_MeshNode*>& dist2Nodes,
double precision)
{
else if ( precision == 0. )
precision = maxSize() / 2;
- if (isInside(node, precision))
+ if ( isInside( node, precision ))
{
- if (!isLeaf())
+ if ( !isLeaf() )
{
// first check a child containing node
gp_XYZ mid = (getBox()->CornerMin() + getBox()->CornerMax()) / 2.;
int nodeChild = getChildIndex( node.X(), node.Y(), node.Z(), mid );
if ( ((SMESH_OctreeNode*) myChildren[nodeChild])->NodesAround(node, dist2Nodes, precision))
return true;
-
+
for (int i = 0; i < 8; i++)
if ( i != nodeChild )
if (((SMESH_OctreeNode*) myChildren[i])->NodesAround(node, dist2Nodes, precision))
else if ( NbNodes() > 0 )
{
double minDist = precision * precision;
- TIDSortedNodeSet::iterator nIt = myNodes.begin();
- for ( ; nIt != myNodes.end(); ++nIt )
+ for ( size_t i = 0; i < myNodes.size(); ++i )
{
- SMESH_TNodeXYZ p2( *nIt );
- double dist2 = ( node - p2 ).SquareModulus();
+ SMESH_NodeXYZ p2 = myNodes[ i ];
+ double dist2 = ( node - p2 ).SquareModulus();
if ( dist2 < minDist )
- dist2Nodes.insert( make_pair( minDist = dist2, p2._node ));
+ dist2Nodes.insert( std::make_pair( minDist = dist2, myNodes[ i ] ));
}
-// if ( dist2Nodes.size() > 1 ) // leave only closest node in dist2Nodes
-// dist2Nodes.erase( ++dist2Nodes.begin(), dist2Nodes.end());
+ // if ( dist2Nodes.size() > 1 ) // leave only closest node in dist2Nodes
+ // dist2Nodes.erase( ++dist2Nodes.begin(), dist2Nodes.end());
// true if an exact overlapping found
return ( sqrt( minDist ) <= precision * 1e-12 );
if ( isLeaf() && NbNodes() )
{
double minDist2 = precision * precision;
- TIDSortedNodeSet::iterator nIt = myNodes.begin();
- for ( ; nIt != myNodes.end(); ++nIt )
+ for ( size_t i = 0; i < myNodes.size(); ++i )
{
- SMESH_TNodeXYZ p2( *nIt );
+ SMESH_NodeXYZ p2 = myNodes[ i ];
double dist2 = ( point - p2 ).SquareModulus();
if ( dist2 <= minDist2 )
- nodes.push_back( p2._node );
+ nodes.push_back( myNodes[ i ] );
}
}
else if ( myChildren )
{
for (int i = 0; i < 8; i++)
{
- SMESH_OctreeNode* myChild = dynamic_cast<SMESH_OctreeNode*> (myChildren[i]);
- myChild->NodesAround( point, nodes, precision);
+ SMESH_OctreeNode* myChild = static_cast<SMESH_OctreeNode*>( myChildren[ i ]);
+ myChild->NodesAround( point, nodes, precision );
}
}
}
* \param maxNbNodes - maximum Nodes in a Leaf of the SMESH_OctreeNode constructed, default value is 5
*/
//=============================
+
void SMESH_OctreeNode::FindCoincidentNodes (TIDSortedNodeSet& theSetOfNodes,
- list< list< const SMDS_MeshNode*> >* theGroupsOfNodes,
- const double theTolerance,
- const int maxLevel,
- const int maxNbNodes)
+ TListOfNodeLists* theGroupsOfNodes,
+ const double theTolerance,
+ const int maxLevel,
+ const int maxNbNodes)
{
- // VSR 14/10/2011: limit max number of the levels in order to avoid endless recursing
+ // VSR 14/10/2011: limit max number of the levels in order to avoid endless recursion
const int MAX_LEVEL = 10;
- SMESH_OctreeNode theOctreeNode(theSetOfNodes, maxLevel < 0 ? MAX_LEVEL : maxLevel, maxNbNodes, theTolerance);
+ SMESH_OctreeNode theOctreeNode(theSetOfNodes,
+ maxLevel < 0 ? MAX_LEVEL : maxLevel,
+ maxNbNodes,
+ theTolerance);
theOctreeNode.FindCoincidentNodes (&theSetOfNodes, theTolerance, theGroupsOfNodes);
}
* \param theGroupsOfNodes - list of nodes closed to each other returned
*/
//=============================
-void SMESH_OctreeNode::FindCoincidentNodes ( TIDSortedNodeSet* theSetOfNodes,
- const double theTolerance,
- list< list< const SMDS_MeshNode*> >* theGroupsOfNodes)
+
+void SMESH_OctreeNode::FindCoincidentNodes ( TIDSortedNodeSet* theSetOfNodes,
+ const double theTolerance,
+ TListOfNodeLists* theGroupsOfNodes )
{
- TIDSortedNodeSet::iterator it1 = theSetOfNodes->begin();
- list<const SMDS_MeshNode*>::iterator it2;
+ // un-mark all nodes; we mark nodes added to theGroupsOfNodes
+ SMESH_MeshAlgos::MarkElems( SMESHUtils::elemSetIterator( *theSetOfNodes ), false );
- list<const SMDS_MeshNode*> ListOfCoincidentNodes;
+ vector<const SMDS_MeshNode*> coincidentNodes;
TIDCompare idLess;
- while (it1 != theSetOfNodes->end())
+ TIDSortedNodeSet::iterator it1 = theSetOfNodes->begin();
+ for ( ; it1 != theSetOfNodes->end(); ++it1 )
{
const SMDS_MeshNode * n1 = *it1;
+ if ( n1->isMarked() )
+ continue;
+ n1->setIsMarked( true );
- // Searching for Nodes around n1 and put them in ListofCoincidentNodes.
+ // Searching for Nodes around n1 and put them in coincidentNodes.
// Found nodes are also erased from theSetOfNodes
- FindCoincidentNodes(n1, theSetOfNodes, &ListOfCoincidentNodes, theTolerance);
+ coincidentNodes.clear();
+ findCoincidentNodes( n1, theSetOfNodes, &coincidentNodes, theTolerance );
- if ( !ListOfCoincidentNodes.empty() )
+ if ( !coincidentNodes.empty() )
{
// We build a list {n1 + his neighbors} and add this list in theGroupsOfNodes
- if ( idLess( n1, ListOfCoincidentNodes.front() )) ListOfCoincidentNodes.push_front( n1 );
- else ListOfCoincidentNodes.push_back ( n1 );
- ListOfCoincidentNodes.sort( idLess );
- theGroupsOfNodes->push_back( list<const SMDS_MeshNode*>() );
- theGroupsOfNodes->back().splice( theGroupsOfNodes->back().end(), ListOfCoincidentNodes );
- }
+ std::sort( coincidentNodes.begin(), coincidentNodes.end(), idLess );
+ list<const SMDS_MeshNode*> newGroup;
+ newGroup.push_back( n1 );
+ newGroup.insert( newGroup.end(), coincidentNodes.begin(), coincidentNodes.end() );
- theSetOfNodes->erase(it1);
- it1 = theSetOfNodes->begin();
+ theGroupsOfNodes->emplace_back( newGroup );
+ }
}
}
* \param precision - Precision used
*/
//======================================================================================
-void SMESH_OctreeNode::FindCoincidentNodes (const SMDS_MeshNode * Node,
- TIDSortedNodeSet* SetOfNodes,
- list<const SMDS_MeshNode*>* Result,
- const double precision)
+
+void SMESH_OctreeNode::findCoincidentNodes (const SMDS_MeshNode * Node,
+ TIDSortedNodeSet* SetOfNodes,
+ std::vector<const SMDS_MeshNode*>* Result,
+ const double precision)
{
- gp_Pnt p1 (Node->X(), Node->Y(), Node->Z());
- bool isInsideBool = isInside( p1.XYZ(), precision );
+ SMESH_NodeXYZ p1 = Node;
- if (isInsideBool)
+ if ( isInside( p1, precision ))
{
// I'm only looking in the leaves, since all the nodes are stored there.
- if (isLeaf())
+ if ( isLeaf() )
{
- TIDSortedNodeSet::iterator it = myNodes.begin();
const double tol2 = precision * precision;
- bool squareBool;
- while (it != myNodes.end())
+ for ( size_t i = 0; i < myNodes.size(); ++i )
{
- const SMDS_MeshNode* n2 = *it;
- squareBool = false;
- // We're only looking at nodes with a superior Id.
- // JFA: Why?
- //if (Node->GetID() < n2->GetID())
- if (Node->GetID() != n2->GetID()) // JFA: for bug 0020185
- {
- gp_Pnt p2 (n2->X(), n2->Y(), n2->Z());
- // Distance optimized computation
- squareBool = (p1.SquareDistance( p2 ) <= tol2);
+ if ( myNodes[ i ]->isMarked() ) // coincident node already found
+ continue;
- // If n2 inside the SquareDistance, we add it in Result and remove it from SetOfNodes and myNodes
- if (squareBool)
+ //if ( Node != myNodes[ i ]) // JFA: for bug 0020185
+ {
+ // If n2 inside the SquareDistance, we add it in Result
+ bool coincide = ( p1.SquareDistance( myNodes[ i ]) <= tol2 );
+ if ( coincide )
{
- Result->insert(Result->begin(), n2);
- SetOfNodes->erase( n2 );
- myNodes.erase( *it++ ); // it++ goes forward and returns it's previous position
+ Result->push_back ( myNodes[ i ]);
+ myNodes[ i ]->setIsMarked( true );
}
}
- if ( !squareBool )
- it++;
}
- if ( !Result->empty() )
- myNodes.erase(Node); // JFA: for bug 0020185
}
else
{
// If I'm not a leaf, I'm going to see my children !
- for (int i = 0; i < 8; i++)
+ for ( int i = 0; i < 8; i++ )
{
- SMESH_OctreeNode* myChild = dynamic_cast<SMESH_OctreeNode*> (myChildren[i]);
- myChild->FindCoincidentNodes(Node, SetOfNodes, Result, precision);
+ SMESH_OctreeNode* myChild = static_cast<SMESH_OctreeNode*> (myChildren[i]);
+ myChild->findCoincidentNodes( Node, SetOfNodes, Result, precision );
}
}
}
{
if ( isLeaf() )
{
- TIDSortedNodeSet::iterator pNode = myNodes.find( node );
- bool nodeInMe = ( pNode != myNodes.end() );
+ std::vector< const SMDS_MeshNode* >::iterator pNode =
+ std::find( myNodes.begin(), myNodes.end(), node );
+ bool nodeInMe = ( pNode != myNodes.end() );
bool pointInMe = isInside( toPnt.Coord(), 1e-10 );
if ( pointInMe != nodeInMe )
{
if ( pointInMe )
- myNodes.insert( node );
+ myNodes.push_back( node );
else
- myNodes.erase( node );
+ myNodes.erase( pNode );
}
}
else if ( myChildren )
* \brief Return iterator over children
*/
//================================================================================
+
SMESH_OctreeNodeIteratorPtr SMESH_OctreeNode::GetChildrenIterator()
{
return SMESH_OctreeNodeIteratorPtr
* \brief Return nodes iterator
*/
//================================================================================
+
SMDS_NodeIteratorPtr SMESH_OctreeNode::GetNodeIterator()
{
- return SMDS_NodeIteratorPtr
- ( new SMDS_SetIterator< SMDS_pNode, TIDSortedNodeSet::const_iterator >
- ( myNodes.begin(), myNodes.size() ? myNodes.end() : myNodes.begin()));
+ return boost::make_shared< SMDS_NodeVectorIterator >( myNodes.begin(), myNodes.end());
}
typedef SMDS_Iterator<SMESH_OctreeNode*> SMESH_OctreeNodeIterator;
typedef boost::shared_ptr<SMESH_OctreeNodeIterator> SMESH_OctreeNodeIteratorPtr;
typedef std::set< const SMDS_MeshNode*, TIDCompare > TIDSortedNodeSet;
+typedef std::list< std::list< const SMDS_MeshNode*> > TListOfNodeLists;
class SMESHUtils_EXPORT SMESH_OctreeNode : public SMESH_Octree
{
virtual const bool isInside(const gp_XYZ& p, const double precision = 0.);
// Return in Result a list of Nodes potentials to be near Node
- void NodesAround(const SMDS_MeshNode * node,
- std::list<const SMDS_MeshNode*>* result,
- const double precision = 0.);
+ void AllNodesAround(const SMDS_MeshNode * node,
+ std::vector<const SMDS_MeshNode*>* result,
+ const double precision = 0.);
// Return in dist2Nodes nodes mapped to their square distance from Node
bool NodesAround(const gp_XYZ& point,
// Return in theGroupsOfNodes a list of group of nodes close to each other within theTolerance
// Search for all the nodes in nodes
- void FindCoincidentNodes ( TIDSortedNodeSet* nodes,
- const double theTolerance,
- std::list< std::list< const SMDS_MeshNode*> >* theGroupsOfNodes);
+ void FindCoincidentNodes ( TIDSortedNodeSet* nodes,
+ const double theTolerance,
+ TListOfNodeLists * theGroupsOfNodes);
// Static method that return in theGroupsOfNodes a list of group of nodes close to each other within
// theTolerance search for all the nodes in nodes
- static void FindCoincidentNodes ( TIDSortedNodeSet& nodes,
- std::list< std::list< const SMDS_MeshNode*> >* theGroupsOfNodes,
- const double theTolerance = 0.00001,
- const int maxLevel = -1,
- const int maxNbNodes = 5);
+ static void FindCoincidentNodes ( TIDSortedNodeSet& nodes,
+ TListOfNodeLists* theGroupsOfNodes,
+ const double theTolerance = 0.00001,
+ const int maxLevel = -1,
+ const int maxNbNodes = 5);
/*!
* \brief Update data according to node movement
*/
- void UpdateByMoveNode( const SMDS_MeshNode* node, const gp_Pnt& toPnt );
+ void UpdateByMoveNode( const SMDS_MeshNode* node, const gp_Pnt& toPnt );
/*!
* \brief Return iterator over children
*/
// Construct an empty SMESH_OctreeNode used by SMESH_Octree::buildChildren()
virtual SMESH_Octree* newChild() const;
- // Return in result a list of nodes closed to Node and remove it from SetOfNodes
- void FindCoincidentNodes( const SMDS_MeshNode * Node,
- TIDSortedNodeSet* SetOfNodes,
- std::list<const SMDS_MeshNode*>* Result,
- const double precision);
+ // Return in result a list of nodes closed to Node
+ void findCoincidentNodes( const SMDS_MeshNode * Node,
+ TIDSortedNodeSet* SetOfNodes,
+ std::vector<const SMDS_MeshNode*>* Result,
+ const double precision);
// The set of nodes inside the box of the Octree (Empty if Octree is not a leaf)
- TIDSortedNodeSet myNodes;
+ std::vector< const SMDS_MeshNode* > myNodes;
};
CORBA::Boolean theCommonGroups)
throw ( SALOME::SALOME_Exception )
{
- typedef list<SMESH::SMESH_Group_var> TListOfNewGroups;
- typedef map< pair<string, SMESH::ElementType>, TListOfNewGroups > TGroupsMap;
-
- TPythonDump* pPythonDump = new TPythonDump;
- TPythonDump& aPythonDump = *pPythonDump; // prevent dump of called methods
+ std::unique_ptr< TPythonDump > pPythonDump( new TPythonDump );
+ TPythonDump& pythonDump = *pPythonDump; // prevent dump of called methods
// create mesh
- SMESH::SMESH_Mesh_var aNewMesh = CreateEmptyMesh();
-
- if ( aNewMesh->_is_nil() )
- return aNewMesh._retn();
-
- SMESH_Mesh_i* aNewImpl = SMESH::DownCast<SMESH_Mesh_i*>( aNewMesh );
- if ( !aNewImpl )
- return aNewMesh._retn();
+ SMESH::SMESH_Mesh_var newMesh = CreateEmptyMesh();
+ SMESH_Mesh_i* newImpl = SMESH::DownCast<SMESH_Mesh_i*>( newMesh );
+ if ( !newImpl ) return newMesh._retn();
- ::SMESH_Mesh& aLocMesh = aNewImpl->GetImpl();
- SMESHDS_Mesh* aNewMeshDS = aLocMesh.GetMeshDS();
+ ::SMESH_Mesh& locMesh = newImpl->GetImpl();
+ SMESHDS_Mesh* newMeshDS = locMesh.GetMeshDS();
- TGroupsMap aGroupsMap;
- TListOfNewGroups aListOfNewGroups;
- ::SMESH_MeshEditor aNewEditor(&aLocMesh);
- SMESH::ListOfGroups_var aListOfGroups;
+ typedef std::list<SMESH::SMESH_Group_var> TListOfNewGroups;
+ typedef std::pair<string, SMESH::ElementType > TNameAndType;
+ typedef std::map< TNameAndType, TListOfNewGroups > TGroupsMap;
+ TGroupsMap groupsMap;
+ TListOfNewGroups listOfNewGroups;
- ::SMESH_MeshEditor::ElemFeatures elemType;
- std::vector<const SMDS_MeshNode*> aNodesArray;
+ ::SMESH_MeshEditor newEditor( &locMesh );
+ ::SMESH_MeshEditor::ElemFeatures elemType;
// loop on sub-meshes
- for ( CORBA::ULong i = 0; i < theMeshesArray.length(); i++)
+ for ( CORBA::ULong i = 0; i < theMeshesArray.length(); i++ )
{
if ( CORBA::is_nil( theMeshesArray[i] )) continue;
- SMESH::SMESH_Mesh_var anInitMesh = theMeshesArray[i]->GetMesh();
- if ( anInitMesh->_is_nil() ) continue;
- SMESH_Mesh_i* anInitImpl = SMESH::DownCast<SMESH_Mesh_i*>( anInitMesh );
- if ( !anInitImpl ) continue;
- anInitImpl->Load();
-
- //::SMESH_Mesh& aInitLocMesh = anInitImpl->GetImpl();
- //SMESHDS_Mesh* anInitMeshDS = aInitLocMesh.GetMeshDS();
+ SMESH::SMESH_Mesh_var initMesh = theMeshesArray[i]->GetMesh();
+ SMESH_Mesh_i* initImpl = SMESH::DownCast<SMESH_Mesh_i*>( initMesh );
+ if ( !initImpl ) continue;
+ initImpl->Load();
+
+ // assure that IDs increments by one during iteration
+ ::SMESH_Mesh& initLocMesh = initImpl->GetImpl();
+ SMESHDS_Mesh* initMeshDS = initLocMesh.GetMeshDS();
+ if ( initMeshDS->MaxNodeID() != initMeshDS->NbNodes() ||
+ initMeshDS->MaxElementID() != initMeshDS->NbElements() )
+ initMeshDS->CompactMesh();
// remember nb of elements before filling in
- SMESH::long_array_var prevState = aNewMesh->GetNbElementsByType();
+ SMESH::long_array_var prevState = newMesh->GetNbElementsByType();
+
+ // copy nodes
- typedef std::map<const SMDS_MeshElement*, const SMDS_MeshElement*, TIDCompare > TEEMap;
- TEEMap elemsMap, nodesMap;
+ std::vector< const SMDS_MeshElement* > newNodes( initMeshDS->NbNodes() + 1, 0 );
+ SMDS_ElemIteratorPtr elemIt = initImpl->GetElements( theMeshesArray[i], SMESH::NODE );
+ while ( elemIt->more() )
+ {
+ SMESH_NodeXYZ node = elemIt->next();
+ newNodes[ node->GetID() ] = newMeshDS->AddNode( node.X(), node.Y(), node.Z() );
+ }
- // loop on elements of a sub-mesh
- SMDS_ElemIteratorPtr itElems = anInitImpl->GetElements( theMeshesArray[i], SMESH::ALL );
- const SMDS_MeshElement* anElem;
- const SMDS_MeshElement* aNewElem;
- const SMDS_MeshNode* aNode;
- const SMDS_MeshNode* aNewNode;
- int anElemNbNodes;
+ // copy elements
- while ( itElems->more() )
+ std::vector< const SMDS_MeshElement* > newElems( initMeshDS->NbElements() + 1, 0 );
+ elemIt = initImpl->GetElements( theMeshesArray[i], SMESH::ALL );
+ while ( elemIt->more() )
{
- anElem = itElems->next();
- anElemNbNodes = anElem->NbNodes();
- aNodesArray.resize( anElemNbNodes );
+ const SMDS_MeshElement* elem = elemIt->next();
+ elemType.myNodes.resize( elem->NbNodes() );
- // loop on nodes of an element
- SMDS_ElemIteratorPtr itNodes = anElem->nodesIterator();
+ SMDS_NodeIteratorPtr itNodes = elem->nodeIterator();
for ( int k = 0; itNodes->more(); k++)
{
- aNode = static_cast<const SMDS_MeshNode*>( itNodes->next() );
- TEEMap::iterator n2nnIt = nodesMap.find( aNode );
- if ( n2nnIt == nodesMap.end() )
- {
- aNewNode = aNewMeshDS->AddNode(aNode->X(), aNode->Y(), aNode->Z());
- nodesMap.insert( make_pair( aNode, aNewNode ));
- }
- else
- {
- aNewNode = static_cast<const SMDS_MeshNode*>( n2nnIt->second );
- }
- aNodesArray[k] = aNewNode;
+ const SMDS_MeshNode* node = itNodes->next();
+ elemType.myNodes[ k ] = static_cast< const SMDS_MeshNode*> ( newNodes[ node->GetID() ]);
}
// creates a corresponding element on existent nodes in new mesh
- if ( anElem->GetType() == SMDSAbs_Node )
- aNewElem = 0;
- else
- aNewElem =
- aNewEditor.AddElement( aNodesArray, elemType.Init( anElem, /*basicOnly=*/false ));
-
- if ( aNewElem )
- elemsMap.insert( make_pair( anElem, aNewElem ));
-
- } //elems loop
-
- aNewEditor.ClearLastCreated(); // forget the history
+ newElems[ elem->GetID() ] =
+ newEditor.AddElement( elemType.myNodes, elemType.Init( elem, /*basicOnly=*/false ));
+ }
+ newEditor.ClearLastCreated(); // forget the history
// create groups of just added elements
- SMESH::SMESH_Group_var aNewGroup;
- SMESH::ElementType aGroupType;
+ SMESH::SMESH_Group_var newGroup;
+ SMESH::ElementType groupType;
if ( theCommonGroups )
{
- SMESH::long_array_var curState = aNewMesh->GetNbElementsByType();
+ // type names
+ const char* typeNames[] = { "All","Nodes","Edges","Faces","Volumes","0DElems","Balls" };
+ { // check of typeNames: compilation failure mains that NB_ELEMENT_TYPES changed:
+ const int nbNames = sizeof(typeNames) / sizeof(const char*);
+ int _assert[( nbNames == SMESH::NB_ELEMENT_TYPES ) ? 2 : -1 ]; _assert[0]=_assert[1]=0;
+ }
+
+ SMESH::long_array_var curState = newMesh->GetNbElementsByType();
- for( aGroupType = SMESH::NODE;
- aGroupType < SMESH::NB_ELEMENT_TYPES;
- aGroupType = (SMESH::ElementType)( aGroupType + 1 ))
+ for( groupType = SMESH::NODE;
+ groupType < SMESH::NB_ELEMENT_TYPES;
+ groupType = (SMESH::ElementType)( groupType + 1 ))
{
- if ( curState[ aGroupType ] <= prevState[ aGroupType ])
- continue;
+ if ( curState[ groupType ] <= prevState[ groupType ])
+ continue; // no elements of groupType added from the i-th mesh
// make a group name
- const char* typeNames[] = { "All","Nodes","Edges","Faces","Volumes","0DElems","Balls" };
- { // check of typeNames: compilation failure mains that NB_ELEMENT_TYPES changed:
- const int nbNames = sizeof(typeNames) / sizeof(const char*);
- int _assert[( nbNames == SMESH::NB_ELEMENT_TYPES ) ? 2 : -1 ]; _assert[0]=_assert[1]=0;
- }
- string groupName = "Gr";
- SALOMEDS::SObject_wrap aMeshSObj = ObjectToSObject( myCurrentStudy, theMeshesArray[i] );
- if ( aMeshSObj ) {
- CORBA::String_var name = aMeshSObj->GetName();
+ std::string groupName = "Gr";
+ SALOMEDS::SObject_wrap meshSO = ObjectToSObject( myCurrentStudy, theMeshesArray[i] );
+ if ( meshSO ) {
+ CORBA::String_var name = meshSO->GetName();
groupName += name;
}
groupName += "_";
- groupName += typeNames[ aGroupType ];
+ groupName += typeNames[ groupType ];
// make and fill a group
- TEEMap & e2neMap = ( aGroupType == SMESH::NODE ) ? nodesMap : elemsMap;
- aNewGroup = aNewImpl->CreateGroup( aGroupType, groupName.c_str() );
- if ( SMESH_Group_i* grp_i = SMESH::DownCast<SMESH_Group_i*>( aNewGroup ))
+ newGroup = newImpl->CreateGroup( groupType, groupName.c_str() );
+ std::vector< const SMDS_MeshElement* > & elemVec =
+ ( groupType == SMESH::NODE ) ? newNodes : newElems;
+ if ( SMESH_Group_i* grp_i = SMESH::DownCast<SMESH_Group_i*>( newGroup ))
{
if ( SMESHDS_Group* grpDS = dynamic_cast<SMESHDS_Group*>( grp_i->GetGroupDS() ))
{
- TEEMap::iterator e2neIt = e2neMap.begin();
- for ( ; e2neIt != e2neMap.end(); ++e2neIt )
+ for ( size_t j = 0; j < elemVec.size(); ++j )
{
- aNewElem = e2neIt->second;
- if ( aNewElem->GetType() == grpDS->GetType() )
- {
- grpDS->Add( aNewElem );
-
- if ( prevState[ aGroupType ]++ >= curState[ aGroupType ] )
- break;
- }
+ if ( elemVec[j] && elemVec[j]->GetType() == grpDS->GetType() )
+ grpDS->Add( elemVec[j] );
}
}
}
- aListOfNewGroups.clear();
- aListOfNewGroups.push_back(aNewGroup);
- aGroupsMap.insert(make_pair( make_pair(groupName, aGroupType), aListOfNewGroups ));
+ listOfNewGroups.clear();
+ listOfNewGroups.push_back( newGroup );
+ groupsMap.insert( std::make_pair( TNameAndType( groupName, groupType ),
+ listOfNewGroups ));
}
}
- if ( SMESH_Mesh_i* anSrcImpl = SMESH::DownCast<SMESH_Mesh_i*>( theMeshesArray[i] ))
+ if ( SMESH_Mesh_i* initImpl = SMESH::DownCast<SMESH_Mesh_i*>( theMeshesArray[i] ))
{
- // copy orphan nodes
- if ( anSrcImpl->NbNodes() > (int)nodesMap.size() )
- {
- SMDS_ElemIteratorPtr itNodes = anInitImpl->GetElements( theMeshesArray[i], SMESH::NODE );
- while ( itNodes->more() )
- {
- const SMDS_MeshNode* aNode = static_cast< const SMDS_MeshNode* >( itNodes->next() );
- if ( aNode->NbInverseElements() == 0 )
- {
- aNewNode = aNewMeshDS->AddNode(aNode->X(), aNode->Y(), aNode->Z());
- nodesMap.insert( make_pair( aNode, aNewNode ));
- }
- }
- }
-
// copy groups
- SMESH::SMESH_GroupBase_ptr aGroup;
- CORBA::String_var aGroupName;
- SMESH::long_array_var anNewIDs = new SMESH::long_array();
+ SMESH::SMESH_GroupBase_ptr group;
+ CORBA::String_var groupName;
+ SMESH::long_array_var newIDs = new SMESH::long_array();
// loop on groups of a source mesh
- aListOfGroups = anSrcImpl->GetGroups();
- for ( CORBA::ULong iG = 0; iG < aListOfGroups->length(); iG++ )
+ SMESH::ListOfGroups_var listOfGroups = initImpl->GetGroups();
+ for ( CORBA::ULong iG = 0; iG < listOfGroups->length(); iG++ )
{
- aGroup = aListOfGroups[iG];
- aGroupType = aGroup->GetType();
- aGroupName = aGroup->GetName();
- string aName = aGroupName.in();
+ group = listOfGroups[iG];
+ groupType = group->GetType();
+ groupName = group->GetName();
+ std::string name = groupName.in();
// convert a list of IDs
- anNewIDs->length( aGroup->Size() );
- TEEMap & e2neMap = ( aGroupType == SMESH::NODE ) ? nodesMap : elemsMap;
- SMDS_ElemIteratorPtr itGrElems = anSrcImpl->GetElements( aGroup, SMESH::ALL );
- int iElem = 0;
+ newIDs->length( group->Size() );
+ std::vector< const SMDS_MeshElement* > & elemVec =
+ ( groupType == SMESH::NODE ) ? newNodes : newElems;
+ SMDS_ElemIteratorPtr itGrElems = initImpl->GetElements( group, SMESH::ALL );
+ int nbElems = 0;
while ( itGrElems->more() )
{
- anElem = itGrElems->next();
- TEEMap::iterator e2neIt = e2neMap.find( anElem );
- if ( e2neIt != e2neMap.end() )
- anNewIDs[ iElem++ ] = e2neIt->second->GetID();
+ const SMDS_MeshElement* elem = itGrElems->next();
+ const SMDS_MeshElement* newElem = elemVec[ elem->GetID() ];
+ if ( newElem )
+ newIDs[ nbElems++ ] = newElem->GetID();
}
- anNewIDs->length( iElem );
+ newIDs->length( nbElems );
- // check a current group name and type don't have identical ones in final mesh
- aListOfNewGroups.clear();
- TGroupsMap::iterator anIter = aGroupsMap.find( make_pair( aName, aGroupType ));
- if ( anIter == aGroupsMap.end() ) {
+ // check that a current group name and type don't have identical ones in final mesh
+ listOfNewGroups.clear();
+ TNameAndType nameAndType( name, groupType );
+ TGroupsMap::iterator anIter = groupsMap.find( nameAndType );
+ if ( anIter == groupsMap.end() )
+ {
// add a new group in the mesh
- aNewGroup = aNewImpl->CreateGroup( aGroupType, aGroupName.in() );
- // add elements into new group
- aNewGroup->Add( anNewIDs );
+ newGroup = newImpl->CreateGroup( groupType, groupName.in() );
+ newGroup->Add( newIDs );
- aListOfNewGroups.push_back(aNewGroup);
- aGroupsMap.insert(make_pair( make_pair(aName, aGroupType), aListOfNewGroups ));
+ listOfNewGroups.push_back( newGroup );
+ groupsMap.insert( std::make_pair( nameAndType, listOfNewGroups ));
}
-
- else if ( theUniteIdenticalGroups ) {
+ else if ( theUniteIdenticalGroups )
+ {
// unite identical groups
TListOfNewGroups& aNewGroups = anIter->second;
- aNewGroups.front()->Add( anNewIDs );
+ aNewGroups.front()->Add( newIDs );
}
-
- else {
+ else
+ {
// rename identical groups
- aNewGroup = aNewImpl->CreateGroup(aGroupType, aGroupName.in());
- aNewGroup->Add( anNewIDs );
+ newGroup = newImpl->CreateGroup( groupType, groupName );
+ newGroup->Add( newIDs );
- TListOfNewGroups& aNewGroups = anIter->second;
- string aNewGroupName;
- if (aNewGroups.size() == 1) {
- aNewGroupName = aName + "_1";
- aNewGroups.front()->SetName(aNewGroupName.c_str());
+ TListOfNewGroups& newGroups = anIter->second;
+ std::string newGroupName;
+ if ( newGroups.size() == 1 )
+ {
+ newGroupName = name + "_1";
+ newGroups.front()->SetName( newGroupName.c_str() );
}
- char aGroupNum[128];
- sprintf(aGroupNum, "%u", (unsigned int)aNewGroups.size()+1);
- aNewGroupName = aName + "_" + string(aGroupNum);
- aNewGroup->SetName(aNewGroupName.c_str());
- aNewGroups.push_back(aNewGroup);
+ newGroupName = name + "_" + SMESH_Comment( newGroups.size() + 1 );
+ newGroup->SetName( newGroupName.c_str() );
+ newGroups.push_back( newGroup );
}
- } //groups loop
+ } // loop on groups
} // if an IDSource is a mesh
} //meshes loop
- if (theMergeNodesAndElements) // merge nodes
+ if ( theMergeNodesAndElements ) // merge nodes
{
- TIDSortedNodeSet aMeshNodes; // no input nodes
- SMESH_MeshEditor::TListOfListOfNodes aGroupsOfNodes;
- aNewEditor.FindCoincidentNodes( aMeshNodes, theMergeTolerance, aGroupsOfNodes,
- /*SeparateCornersAndMedium=*/ false );
- aNewEditor.MergeNodes( aGroupsOfNodes );
+ TIDSortedNodeSet meshNodes; // no input nodes == treat all
+ SMESH_MeshEditor::TListOfListOfNodes groupsOfNodes;
+ newEditor.FindCoincidentNodes( meshNodes, theMergeTolerance, groupsOfNodes,
+ /*SeparateCornersAndMedium=*/ false );
+ newEditor.MergeNodes( groupsOfNodes );
// merge elements
- aNewEditor.MergeEqualElements();
+ newEditor.MergeEqualElements();
}
// Update Python script
- aPythonDump << aNewMesh << " = " << this << "."
- << ( theCommonGroups ? "ConcatenateWithGroups" : "Concatenate" )
- << "([";
- for ( CORBA::ULong i = 0; i < theMeshesArray.length(); i++) {
- if (i > 0) aPythonDump << ", ";
- aPythonDump << theMeshesArray[i];
- }
- aPythonDump << "], ";
- aPythonDump << theUniteIdenticalGroups << ", "
- << theMergeNodesAndElements << ", "
- << TVar( theMergeTolerance ) << ")";
+ pythonDump << newMesh << " = " << this
+ << "." << ( theCommonGroups ? "ConcatenateWithGroups" : "Concatenate" ) << "("
+ << theMeshesArray << ", "
+ << theUniteIdenticalGroups << ", "
+ << theMergeNodesAndElements << ", "
+ << TVar( theMergeTolerance ) << ")";
- delete pPythonDump; // enable python dump from GetGroups()
+ pPythonDump.reset(); // enable python dump from GetGroups()
// 0020577: EDF 1164 SMESH: Bad dump of concatenate with create common groups
- if ( !aNewMesh->_is_nil() )
+ if ( !newMesh->_is_nil() )
{
- SMESH::ListOfGroups_var groups = aNewMesh->GetGroups();
+ SMESH::ListOfGroups_var groups = newMesh->GetGroups();
}
// IPAL21468 Change icon of compound because it need not be computed.
- SALOMEDS::SObject_wrap aMeshSObj = ObjectToSObject( myCurrentStudy, aNewMesh );
- SetPixMap( aMeshSObj, "ICON_SMESH_TREE_MESH" );
+ SALOMEDS::SObject_wrap meshSO = ObjectToSObject( myCurrentStudy, newMesh );
+ SetPixMap( meshSO, "ICON_SMESH_TREE_MESH" );
- if (aNewMeshDS)
- aNewMeshDS->Modified();
+ newMeshDS->Modified();
- return aNewMesh._retn();
+ return newMesh._retn();
}
//================================================================================
{
const SMDS_MeshElement* res = _node;
_node = 0;
- while (( _elemIter->more() || _nodeIter->more() ) && !_node )
+ while ( !_node && ( _elemIter->more() || _nodeIter->more() ))
{
if ( _nodeIter->more() )
{