Salome HOME
Implementation of PAL13462 (Complete MergeNodes).
authorrnv <rnv@opencascade.com>
Fri, 2 Mar 2007 08:03:40 +0000 (08:03 +0000)
committerrnv <rnv@opencascade.com>
Fri, 2 Mar 2007 08:03:40 +0000 (08:03 +0000)
idl/SMESH_MeshEditor.idl
src/SMESHGUI/SMESHGUI_MergeNodesDlg.cxx
src/SMESHGUI/SMESHGUI_MergeNodesDlg.h
src/SMESH_I/SMESH_MeshEditor_i.cxx
src/SMESH_I/SMESH_MeshEditor_i.hxx
src/SMESH_SWIG/smesh.py

index 81c2bf9903489941f12277d0c9f6988470726895..492e95440bfcca176d74a8f1b5cf266c51eeac3d 100644 (file)
@@ -329,6 +329,10 @@ module SMESH
     void FindCoincidentNodes (in  double              Tolerance,
                               out array_of_long_array GroupsOfNodes);
 
+    void FindCoincidentNodesOnPart (in  SMESH_IDSource      SubMeshOrGroup,
+                                   in  double              Tolerance,
+                                   out array_of_long_array GroupsOfNodes);
+
     void MergeNodes (in array_of_long_array GroupsOfNodes);
 
     void MergeEqualElements();
index 203d28328bb615a67666721a49da81fae3f617cc..aea515e4039e2490469f86c0d4f91dd04db582ea 100644 (file)
@@ -36,6 +36,8 @@
 
 #include "SMESH_Actor.h"
 #include "SMESH_TypeFilter.hxx"
+#include "SMESH_LogicalFilter.hxx"
+#include "SMESHGUI_MeshUtils.h"
 #include "SMDS_Mesh.hxx"
 
 #include "GEOMBase.h"
@@ -54,8 +56,6 @@
 
 #include "utilities.h"
 
-#include CORBA_SERVER_HEADER(SMESH_MeshEditor)
-
 // OCCT Includes
 #include <TColStd_MapOfInteger.hxx>
 
@@ -74,6 +74,9 @@
 #include <qpixmap.h>
 #include <qheader.h>
 
+//IDL Headers
+#include CORBA_SERVER_HEADER(SMESH_Group)
+
 using namespace std;
 
 //=================================================================================
@@ -160,7 +163,7 @@ SMESHGUI_MergeNodesDlg::SMESHGUI_MergeNodesDlg( SMESHGUI* theModule, const char*
 
   // Controls for mesh defining
   GroupMesh = new QGroupBox(this, "GroupMesh");
-  GroupMesh->setTitle(tr("SMESH_MESH"));
+  GroupMesh->setTitle(tr("SMESH_SELECT_WHOLE_MESH"));
   GroupMesh->setColumnLayout(0, Qt::Vertical);
   GroupMesh->layout()->setSpacing(0);
   GroupMesh->layout()->setMargin(0);
@@ -256,12 +259,24 @@ SMESHGUI_MergeNodesDlg::SMESHGUI_MergeNodesDlg( SMESHGUI* theModule, const char*
   myEditCurrentArgument = (QWidget*)LineEditMesh; 
 
   myActor = 0;
+  mySubMeshOrGroup = SMESH::SMESH_subMesh::_nil();
 
   mySelector = (SMESH::GetViewWindow( mySMESHGUI ))->GetSelector();
 
   mySMESHGUI->SetActiveDialogBox((QDialog*)this);
-
-  myMeshOrSubMeshFilter = new SMESH_TypeFilter (MESHorSUBMESH);
+  
+  // Costruction of the logical filter
+  SMESH_TypeFilter* aMeshOrSubMeshFilter = new SMESH_TypeFilter (MESHorSUBMESH);
+  SMESH_TypeFilter* aSmeshGroupFilter    = new SMESH_TypeFilter (GROUP);
+  
+  QPtrList<SUIT_SelectionFilter> aListOfFilters;
+  if (aMeshOrSubMeshFilter) aListOfFilters.append(aMeshOrSubMeshFilter);
+  if (aSmeshGroupFilter)    aListOfFilters.append(aSmeshGroupFilter);
+
+  myMeshOrSubMeshOrGroupFilter =
+    new SMESH_LogicalFilter (aListOfFilters, SMESH_LogicalFilter::LO_OR);
+  
+  //myMeshOrSubMeshFilter = new SMESH_TypeFilter (MESHorSUBMESH);
 
   /* signals and slots connections */
   connect(buttonOk, SIGNAL(clicked()),     this, SLOT(ClickOnOk()));
@@ -464,7 +479,10 @@ void SMESHGUI_MergeNodesDlg::onDetect()
     ListEdit->clear();
 
     SMESH::array_of_long_array_var aNodeGroups;
-    aMeshEditor->FindCoincidentNodes(SpinBoxTolerance->GetValue(), aNodeGroups);
+    if(!mySubMeshOrGroup->_is_nil())
+      aMeshEditor->FindCoincidentNodesOnPart(mySubMeshOrGroup, SpinBoxTolerance->GetValue(), aNodeGroups);
+    else
+      aMeshEditor->FindCoincidentNodes(SpinBoxTolerance->GetValue(), aNodeGroups);
 
     for (int i = 0; i < aNodeGroups->length(); i++) {
       SMESH::long_array& aGroup = aNodeGroups[i];
@@ -626,7 +644,7 @@ void SMESHGUI_MergeNodesDlg::SetEditCurrentArgument()
     SMESH::SetPointRepresentation(false);
     if ( SVTK_ViewWindow* aViewWindow = SMESH::GetViewWindow( mySMESHGUI ))
       aViewWindow->SetSelectionMode(ActorSelection);
-    mySelectionMgr->installFilter(myMeshOrSubMeshFilter);
+    mySelectionMgr->installFilter(myMeshOrSubMeshOrGroupFilter);
   }
 
   myEditCurrentArgument->setFocus();
@@ -643,23 +661,37 @@ void SMESHGUI_MergeNodesDlg::SelectionIntoArgument()
   if (myEditCurrentArgument == (QWidget*)LineEditMesh) {
     QString aString = "";
     LineEditMesh->setText(aString);
-
+    
     ListCoincident->clear();
     ListEdit->clear();
-
+    myActor = 0;
+    
     int nbSel = SMESH::GetNameOfSelectedIObjects(mySelectionMgr, aString);
     if (nbSel != 1)
       return;
 
     SALOME_ListIO aList;
     mySelectionMgr->selectedObjects(aList,SVTK_Viewer::Type());
-
+    
     Handle(SALOME_InteractiveObject) IO = aList.First();
-    myMesh = SMESH::IObjectToInterface<SMESH::SMESH_Mesh>(IO);
-    myActor = SMESH::FindActorByEntry(aList.First()->getEntry());
-    if (myMesh->_is_nil() || !myActor)
+    myMesh = SMESH::GetMeshByIO(IO);
+    
+    if (myMesh->_is_nil())
       return;
-
+    
+    myActor = SMESH::FindActorByEntry(IO->getEntry());
+    if (!myActor)
+      myActor = SMESH::FindActorByObject(myMesh);
+    if(!myActor)
+      return;
+    
+    mySubMeshOrGroup = SMESH::SMESH_IDSource::_nil();
+    
+    if ((!SMESH::IObjectToInterface<SMESH::SMESH_subMesh>(IO)->_is_nil() || //SUBMESH OR GROUP
+         !SMESH::IObjectToInterface<SMESH::SMESH_GroupBase>(IO)->_is_nil()) &&
+        !SMESH::IObjectToInterface<SMESH::SMESH_IDSource>(IO)->_is_nil())
+      mySubMeshOrGroup = SMESH::IObjectToInterface<SMESH::SMESH_IDSource>(IO);
+     
     LineEditMesh->setText(aString);
   }
 }
index eda961e8c0cc58341743cf5756b24b13281a0b53..4177984e6d35a820ffd67badea6dc025fa802a35 100644 (file)
@@ -56,7 +56,7 @@ class SVTK_Selector;
 
 // IDL Headers
 #include <SALOMEconfig.h>
-#include CORBA_SERVER_HEADER(SMESH_Mesh)
+#include CORBA_SERVER_HEADER(SMESH_MeshEditor)
 
 
 //=================================================================================
@@ -89,9 +89,10 @@ private:
     QWidget*                      myEditCurrentArgument;
 
     SMESH::SMESH_Mesh_var         myMesh;
+    SMESH::SMESH_IDSource_var     mySubMeshOrGroup;
     SMESH_Actor*                  myActor;
     //Handle(SMESH_TypeFilter)      myMeshOrSubMeshFilter;
-    SUIT_SelectionFilter*         myMeshOrSubMeshFilter;
+    SUIT_SelectionFilter*         myMeshOrSubMeshOrGroupFilter;
 
     QButtonGroup*     GroupConstructors;
     QRadioButton*     RadioButton1;
index 80bcaa23ddcfff82b29ce8d46186a3dc6c247f3a..13aa0d7143e342e67b45a5112a8e48a0fac0edec 100644 (file)
@@ -1776,6 +1776,64 @@ void SMESH_MeshEditor_i::FindCoincidentNodes (CORBA::Double                  Tol
                 << Tolerance << " )";
 }
 
+//=======================================================================
+//function : FindCoincidentNodesOnPart
+//purpose  :
+//=======================================================================
+void SMESH_MeshEditor_i::FindCoincidentNodesOnPart(SMESH::SMESH_IDSource_ptr      theObject,
+                                                   CORBA::Double                  Tolerance,
+                                                   SMESH::array_of_long_array_out GroupsOfNodes)
+{
+  initData();
+  SMESH::long_array_var aElementsId = theObject->GetIDs();
+
+  SMESHDS_Mesh* aMesh = GetMeshDS();
+  set<const SMDS_MeshNode*> nodes;
+
+  if ( !CORBA::is_nil(SMESH::SMESH_GroupBase::_narrow(theObject)) &&
+      SMESH::SMESH_GroupBase::_narrow(theObject)->GetType() == SMESH::NODE) {
+    for(int i = 0; i < aElementsId->length(); i++) {
+      CORBA::Long ind = aElementsId[i];
+      const SMDS_MeshNode * elem = aMesh->FindNode(ind);
+      if(elem)
+        nodes.insert(elem);
+    }
+  }
+  else {
+    for(int i = 0; i < aElementsId->length(); i++) {
+      CORBA::Long ind = aElementsId[i];
+      const SMDS_MeshElement * elem = aMesh->FindElement(ind);
+      if(elem) {
+        SMDS_ElemIteratorPtr nIt = elem->nodesIterator();
+        while ( nIt->more() )
+          nodes.insert( nodes.end(),static_cast<const SMDS_MeshNode*>(nIt->next()));
+      }
+    }
+  }
+    
+  
+  ::SMESH_MeshEditor::TListOfListOfNodes aListOfListOfNodes;
+  ::SMESH_MeshEditor anEditor( myMesh );
+  if(!nodes.empty())
+    anEditor.FindCoincidentNodes( nodes, Tolerance, aListOfListOfNodes );
+  
+  GroupsOfNodes = new SMESH::array_of_long_array;
+  GroupsOfNodes->length( aListOfListOfNodes.size() );
+  ::SMESH_MeshEditor::TListOfListOfNodes::iterator llIt = aListOfListOfNodes.begin();
+  for ( CORBA::Long i = 0; llIt != aListOfListOfNodes.end(); llIt++, i++ ) {
+    list< const SMDS_MeshNode* >& aListOfNodes = *llIt;
+    list< const SMDS_MeshNode* >::iterator lIt = aListOfNodes.begin();;
+    SMESH::long_array& aGroup = GroupsOfNodes[ i ];
+    aGroup.length( aListOfNodes.size() );
+    for ( int j = 0; lIt != aListOfNodes.end(); lIt++, j++ )
+      aGroup[ j ] = (*lIt)->GetID();
+  }
+  // Update Python script
+  TPythonDump() << "coincident_nodes_on_part = " << this << ".FindCoincidentNodesOnPart( "
+                <<theObject<<", "
+                << Tolerance << " )";
+}
+
 //=======================================================================
 //function : MergeNodes
 //purpose  :
index 24a91212fa9a2873e49fbec1b16bce57faa7f472..331177bbdf32dbd3aa51920a124ad9128bc08112 100644 (file)
@@ -204,6 +204,9 @@ class SMESH_MeshEditor_i: public POA_SMESH::SMESH_MeshEditor
 
   void FindCoincidentNodes (CORBA::Double                  Tolerance,
                             SMESH::array_of_long_array_out GroupsOfNodes);
+  void FindCoincidentNodesOnPart(SMESH::SMESH_IDSource_ptr      theObject,
+                                 CORBA::Double                  Tolerance,
+                                 SMESH::array_of_long_array_out GroupsOfNodes);
   void MergeNodes (const SMESH::array_of_long_array& GroupsOfNodes);
   void MergeEqualElements();
   CORBA::Long MoveClosestNodeToPoint(CORBA::Double x,
index d25d317e63991774d537977e553e82fee3b13628..d4df27ac668050b7ca76d89756331b58483edb36 100644 (file)
@@ -2397,6 +2397,13 @@ class Mesh:
     def FindCoincidentNodes (self, Tolerance):
         return self.editor.FindCoincidentNodes(Tolerance)
 
+    ## Find group of nodes close to each other within Tolerance.
+    #  @param Tolerance tolerance value
+    #  @param SubMeshOrGroup SubMesh or Group
+    #  @param list of group of nodes
+    def FindCoincidentNodesOnPart (self, SubMeshOrGroup, Tolerance):
+        return self.editor.FindCoincidentNodesOnPart(SubMeshOrGroup, Tolerance)
+
     ## Merge nodes
     #  @param list of group of nodes
     def MergeNodes (self, GroupsOfNodes):