diff --git a/src/SMDS/SMDS_Mesh.cxx b/src/SMDS/SMDS_Mesh.cxx index 6b0ca83d5..07385dc90 100644 --- a/src/SMDS/SMDS_Mesh.cxx +++ b/src/SMDS/SMDS_Mesh.cxx @@ -731,11 +731,7 @@ SMDS_MeshFace * SMDS_Mesh::createQuadrangle(SMDS_MeshNode * node1, void SMDS_Mesh::RemoveNode(const SMDS_MeshNode * node) { - SMDS_Iterator * it= - node->GetInverseElementIterator(); - while(it->more()) RemoveElement(it->next(),true); - myNodeIDFactory->ReleaseID(node->GetID()); - myNodes.erase(const_cast(node)); + RemoveElement(node, true); } /////////////////////////////////////////////////////////////////////////////// @@ -744,10 +740,7 @@ void SMDS_Mesh::RemoveNode(const SMDS_MeshNode * node) void SMDS_Mesh::RemoveEdge(const SMDS_MeshEdge * edge) { - /** @todo to be fix */ - myEdges.erase(const_cast(edge)); - //removeElementDependencies(edge); - delete edge; + RemoveElement(edge,true); } /////////////////////////////////////////////////////////////////////////////// @@ -756,10 +749,7 @@ void SMDS_Mesh::RemoveEdge(const SMDS_MeshEdge * edge) void SMDS_Mesh::RemoveFace(const SMDS_MeshFace * face) { - /** @todo to be fix */ - myFaces.erase(const_cast(face)); - //removeElementDependencies(face); - delete face; + RemoveElement(face, true); } /////////////////////////////////////////////////////////////////////////////// @@ -768,64 +758,7 @@ void SMDS_Mesh::RemoveFace(const SMDS_MeshFace * face) void SMDS_Mesh::RemoveVolume(const SMDS_MeshVolume * volume) { - /** @todo to be fix */ - myVolumes.erase(const_cast(volume)); - //removeElementDependencies(volume); - delete volume; -} - -/////////////////////////////////////////////////////////////////////////////// -/// Remove no longer used sub element of an element. Unbind the element ID -/////////////////////////////////////////////////////////////////////////////// - -void SMDS_Mesh::removeElementDependencies(SMDS_MeshElement * element) -{ - /** @todo to be fix */ - myElementIDFactory->ReleaseID(element->GetID()); - SMDS_Iterator * it=element->nodesIterator(); - while(it->more()) - { - SMDS_MeshNode * node=static_cast( - const_cast(it->next())); - node->RemoveInverseElement(element); - if(node->emptyInverseElements()) RemoveNode(node); - } -} - -//======================================================================= -//function : RemoveElement -//purpose : -//======================================================================= - -void SMDS_Mesh::RemoveElement(const SMDS_MeshElement * elem, - const bool removenodes) -{ - /** @todo to be fix */ - switch(elem->GetType()) - { - case SMDSAbs_Node: - RemoveNode((const SMDS_MeshNode*)elem); - return; - case SMDSAbs_Edge: - RemoveEdge((const SMDS_MeshEdge*)elem); - break; - case SMDSAbs_Face: - RemoveFace((const SMDS_MeshFace*)elem); - break; - case SMDSAbs_Volume: - RemoveVolume((const SMDS_MeshVolume*)elem); - break; - default : - MESSAGE("remove function : unknown type"); - return; - } -/* - SMDS_Iterator * it=elem->nodesIterator(); - while(it->more()) - { - const SMDS_MeshNode * node=it->next(); - - }*/ + RemoveElement(volume, true); } //======================================================================= @@ -1440,3 +1373,201 @@ SMDS_Iterator * SMDS_Mesh::volumesIterator() const return new MyIterator(myVolumes); } +/////////////////////////////////////////////////////////////////////////////// +/// Do intersection of sets (more than 2) +/////////////////////////////////////////////////////////////////////////////// +set * intersectionOfSets( + set vs[], int numberOfSets) +{ + set* rsetA=new set(vs[0]); + set* rsetB; + + for(int i=0; i(); + set_intersection( + rsetA->begin(), rsetA->end(), + vs[i+1].begin(), vs[i+1].end(), + inserter(*rsetB, rsetB->begin())); + delete rsetA; + rsetA=rsetB; + } + return rsetA; +} + +/////////////////////////////////////////////////////////////////////////////// +/// Return the list of finit elements owning the given element +/////////////////////////////////////////////////////////////////////////////// +set * getFinitElements(const SMDS_MeshElement * element) +{ + int numberOfSets=element->NbNodes(); + set initSet[numberOfSets]; + + SMDS_Iterator * itNodes=element->nodesIterator(); + + int i=0; + while(itNodes->more()) + { + const SMDS_MeshNode * n=static_cast(itNodes->next()); + SMDS_Iterator * itFe = n->GetInverseElementIterator(); + + //initSet[i]=set(); + while(itFe->more()) initSet[i].insert(itFe->next()); + + i++; + delete itFe; + } + delete itNodes; + + return intersectionOfSets(initSet, numberOfSets); +} + +/////////////////////////////////////////////////////////////////////////////// +/// Return the list of nodes used only by the given elements +/////////////////////////////////////////////////////////////////////////////// +set * getExclusiveNodes( + set& elements) +{ + set * toReturn=new set(); + set::iterator itElements=elements.begin(); + + while(itElements!=elements.end()) + { + SMDS_Iterator * itNodes= + (*itElements)->nodesIterator(); + itElements++; + + while(itNodes->more()) + { + const SMDS_MeshNode * n=static_cast(itNodes->next()); + SMDS_Iterator * itFe = n->GetInverseElementIterator(); + set s; + while(itFe->more()) s.insert(itFe->next()); + delete itFe; + if(s==elements) toReturn->insert(n); + } + delete itNodes; + } + return toReturn; +} + +/////////////////////////////////////////////////////////////////////////////// +///Find the children of an element that are made of given nodes +///@param setOfChildren The set in which matching children will be inserted +///@param element The element were to search matching children +///@param nodes The nodes that the children must have to be selected +/////////////////////////////////////////////////////////////////////////////// +void SMDS_Mesh::addChildrenWithNodes(set& setOfChildren, + const SMDS_MeshElement * element, set& nodes) +{ + + switch(element->GetType()) + { + case SMDSAbs_Node: + MESSAGE("Internal Error: This should not append"); + break; + case SMDSAbs_Edge: + { + SMDS_Iterator * itn=element->nodesIterator(); + while(itn->more()) + { + const SMDS_MeshElement * e=itn->next(); + if(nodes.find(e)!=nodes.end()) setOfChildren.insert(element); + } + delete itn; + } break; + case SMDSAbs_Face: + { + SMDS_Iterator * itn=element->nodesIterator(); + while(itn->more()) + { + const SMDS_MeshElement * e=itn->next(); + if(nodes.find(e)!=nodes.end()) setOfChildren.insert(element); + } + delete itn; + if(hasConstructionEdges()) + { + SMDS_Iterator* ite=element->edgesIterator(); + while(ite->more()) + addChildrenWithNodes(setOfChildren, ite->next(), nodes); + delete ite; + } + } break; + case SMDSAbs_Volume: + { + if(hasConstructionFaces()) + { + SMDS_Iterator * ite=element->facesIterator(); + while(ite->more()) + addChildrenWithNodes(setOfChildren, ite->next(), nodes); + delete ite; + } + else if(hasConstructionEdges()) + { + SMDS_Iterator * ite=element->edgesIterator(); + while(ite->more()) + addChildrenWithNodes(setOfChildren, ite->next(), nodes); + delete ite; + } + } + } +} + +/////////////////////////////////////////////////////////////////////////////// +///@param elem The element to delete +///@param removenodes if true remaining nodes will be removed +/////////////////////////////////////////////////////////////////////////////// +void SMDS_Mesh::RemoveElement(const SMDS_MeshElement * elem, + const bool removenodes) +{ + set * s1=getFinitElements(elem); + + set * s2=getExclusiveNodes(*s1); + set s3; + set::iterator it=s1->begin(); + while(it!=s1->end()) + { + addChildrenWithNodes(s3, *it ,*s2); + s3.insert(*it); + it++; + } + if(elem->GetType()!=SMDSAbs_Node) s3.insert(elem); + it=s3.begin(); + while(it!=s3.end()) + { + switch((*it)->GetType()) + { + case SMDSAbs_Node: + MESSAGE("Internal Error: This should not happen"); + break; + case SMDSAbs_Edge: + myEdges.erase(static_cast( + const_cast(*it))); + break; + case SMDSAbs_Face: + myFaces.erase(static_cast( + const_cast(*it))); + break; + case SMDSAbs_Volume: + myVolumes.erase(static_cast( + const_cast(*it))); + break; + } + delete (*it); + it++; + } + if(removenodes) + { + it=s2->begin(); + while(it!=s2->end()) + { + myNodes.erase(static_cast( + const_cast(*it))); + delete *it; + it++; + } + } + + delete s2; + delete s1; +} diff --git a/src/SMDS/SMDS_Mesh.hxx b/src/SMDS/SMDS_Mesh.hxx index 5af0795fd..f869cd705 100644 --- a/src/SMDS/SMDS_Mesh.hxx +++ b/src/SMDS/SMDS_Mesh.hxx @@ -190,7 +190,6 @@ class SMDS_Mesh:public SMDS_MeshObject SMDS_MeshNode * node2, SMDS_MeshNode * node3); SMDS_MeshFace * createQuadrangle(SMDS_MeshNode * node1, SMDS_MeshNode * node2, SMDS_MeshNode * node3, SMDS_MeshNode * node4); - void removeElementDependencies(SMDS_MeshElement * element); const SMDS_MeshEdge* FindEdge(const SMDS_MeshNode * n1, const SMDS_MeshNode * n2) const; SMDS_MeshEdge* FindEdgeOrCreate(const SMDS_MeshNode * n1, @@ -209,6 +208,8 @@ class SMDS_Mesh:public SMDS_MeshObject const SMDS_MeshNode *n2, const SMDS_MeshNode *n3, const SMDS_MeshNode *n4); + void addChildrenWithNodes(set& setOfChildren, + const SMDS_MeshElement * element, set& nodes); // Fields PRIVATE typedef set SetOfNodes;