]> SALOME platform Git repositories - modules/gui.git/commitdiff
Salome HOME
Junctions imp.
authorenk <enk@opencascade.com>
Wed, 24 Oct 2007 06:19:29 +0000 (06:19 +0000)
committerenk <enk@opencascade.com>
Wed, 24 Oct 2007 06:19:29 +0000 (06:19 +0000)
selection

src/SVTK/SALOME_Actor.cxx
src/SVTK/SALOME_Actor.h
src/SVTK/SVTK_Actor.cxx

index eb9d4ebcfd3cae1aa5e2631cdafc29fc43dbce54..ad6d31369bd33618568b411a12356c314451fd3f 100644 (file)
@@ -58,6 +58,7 @@
 
 #include <vtkInteractorStyle.h>
 #include <vtkRenderWindowInteractor.h>
+#include <vtkUnstructuredGrid.h>
 
 #include <TColStd_MapOfInteger.hxx>
 #include <TColStd_IndexedMapOfInteger.hxx>
@@ -831,3 +832,16 @@ SALOME_Actor
 {
   myHighlightActor->SetProperty(theProperty);
 }
+
+void
+SALOME_Actor
+::MapCells(const TColStd_IndexedMapOfInteger& theMapIndex,
+                vtkUnstructuredGrid* theUG)
+{
+  int aNbOfParts = theMapIndex.Extent();
+  for(int ind = 1; ind <= aNbOfParts; ind++){
+    int aPartId = theMapIndex( ind );
+    if(vtkCell* aCell = GetElemCell(aPartId))
+      theUG->InsertNextCell(aCell->GetCellType(),aCell->GetPointIds());
+  }
+}
\ No newline at end of file
index a44acc9fd35a369ded98aec4110e9889dc4572a3..a5677703b181c83c7feb1ffc5ddf7200f38c16aa 100644 (file)
@@ -49,6 +49,7 @@ class vtkCellPicker;
 class vtkOutlineSource;
 class vtkInteractorStyle;
 class vtkRenderWindowInteractor;
+class vtkUnstructuredGrid;
 
 class SVTK_Actor;
 class SVTK_RectPicker;
@@ -202,6 +203,9 @@ class SVTK_EXPORT SALOME_Actor : public VTKViewer_Actor
   void
   SetHighlightProperty(vtkProperty* theProperty);
 
+  virtual void MapCells(const TColStd_IndexedMapOfInteger& theMapIndex,
+                        vtkUnstructuredGrid* theUG);
+
  protected:
   //----------------------------------------------------------------------------
   vtkRenderWindowInteractor* myInteractor;
index 1d32a1c2b1f099d19fb54268ed880cc76671e4ad..dd34ab95eb2cf7756a04849134286438433e9dcb 100644 (file)
@@ -116,12 +116,7 @@ SVTK_Actor
   vtkDataSet *aSourceDataSet = theMapActor->GetInput();
   CopyPoints(GetSource(),aSourceDataSet);
 
-  int aNbOfParts = theMapIndex.Extent();
-  for(int ind = 1; ind <= aNbOfParts; ind++){
-    int aPartId = theMapIndex( ind );
-    if(vtkCell* aCell = theMapActor->GetElemCell(aPartId))
-      myUnstructuredGrid->InsertNextCell(aCell->GetCellType(),aCell->GetPointIds());
-  }
+  theMapActor->MapCells(theMapIndex,myUnstructuredGrid.GetPointer());
 
   UnShrink();
   if(theMapActor->IsShrunk()){