Salome HOME
IPAL22870: Incorrect numbering of volume elements
authorouv <ouv@opencascade.com>
Wed, 14 Mar 2012 08:59:47 +0000 (08:59 +0000)
committerouv <ouv@opencascade.com>
Wed, 14 Mar 2012 08:59:47 +0000 (08:59 +0000)
src/OBJECT/SMESH_CellLabelActor.cxx
src/OBJECT/SMESH_CellLabelActor.h
src/OBJECT/SMESH_DeviceActor.cxx
src/OBJECT/SMESH_NodeLabelActor.cxx
src/OBJECT/SMESH_NodeLabelActor.h

index 399bc869f78ebb824c0305ce75e71a9e4a7efc67..ad297451d5bb8ca0740b0caeefab5566011bc175 100644 (file)
@@ -38,7 +38,7 @@
 #include <vtkPointData.h>
 #include <vtkProperty2D.h>
 #include <vtkRenderer.h>
-#include <vtkPolyData.h>
+#include <vtkUnstructuredGrid.h>
 #include <vtkCellData.h>
 
 vtkStandardNewMacro(SMESH_CellLabelActor);
@@ -49,7 +49,7 @@ vtkStandardNewMacro(SMESH_CellLabelActor);
 SMESH_CellLabelActor::SMESH_CellLabelActor() {
     //Definition of cells numbering pipeline
   //---------------------------------------
-  myCellsNumDataSet = vtkPolyData::New();
+  myCellsNumDataSet = vtkUnstructuredGrid::New();
 
   myCellCenters = VTKViewer_CellCenters::New();
   myCellCenters->SetInput(myCellsNumDataSet);
@@ -125,19 +125,19 @@ SMESH_CellLabelActor::~SMESH_CellLabelActor() {
 
 void SMESH_CellLabelActor::SetCellsLabeled(bool theIsCellsLabeled) {
   myTransformFilter->Update();
-  vtkPolyData* aGrid = vtkPolyData::SafeDownCast(myTransformFilter->GetOutput());
+  vtkUnstructuredGrid* aGrid = vtkUnstructuredGrid::SafeDownCast(myTransformFilter->GetOutput());
   if(!aGrid)
     return;
 
   myIsCellsLabeled = theIsCellsLabeled && aGrid->GetNumberOfPoints();
   if(myIsCellsLabeled){
     myCellsNumDataSet->ShallowCopy(aGrid);
-    vtkDataSet *aDataSet = myCellsNumDataSet;
+    vtkUnstructuredGrid *aDataSet = myCellsNumDataSet;
     int aNbElem = aDataSet->GetNumberOfCells();
     vtkIntArray *anArray = vtkIntArray::New();
     anArray->SetNumberOfValues(aNbElem);
     for(int anId = 0; anId < aNbElem; anId++){
-      int aSMDSId = GetElemObjId(anId);
+      int aSMDSId = myVisualObj->GetElemObjId(anId);
       anArray->SetValue(anId,aSMDSId);
     }
     aDataSet->GetCellData()->SetScalars(anArray);
index dda47fea53a391ce448b66c4de7b9c8242dfc11e..3b2ba3b3b6d941c0ada608b74c22911b6b9c8951 100644 (file)
@@ -32,7 +32,7 @@ class vtkSelectVisiblePoints;
 class vtkLabeledDataMapper;
 class vtkActor2D;
 class vtkMaskPoints;
-class vtkPolyData;
+class vtkUnstructuredGrid;
 
 class VTKViewer_CellCenters;
 
@@ -65,7 +65,7 @@ protected:
   ~SMESH_CellLabelActor();
 
   bool myIsCellsLabeled;
-  vtkPolyData* myCellsNumDataSet;
+  vtkUnstructuredGrid* myCellsNumDataSet;
   vtkActor2D *myCellsLabels;
   vtkMaskPoints* myClsMaskPoints;
   VTKViewer_CellCenters* myCellCenters;
index 43bf340cbf9c164acc522a673a90d35b7a7ad050..113ae791f144dea86803576c0c2692f35ac77c72 100644 (file)
@@ -240,17 +240,17 @@ SMESH_DeviceActor
     myPassFilter[ anId + 1]->SetInput( myPassFilter[ anId ]->GetOutput() );
     
     anId++; // 1
-    myGeomFilter->SetInput( myPassFilter[ anId ]->GetOutput() );
+    myTransformFilter->SetInput( myPassFilter[ anId ]->GetOutput() );
 
     anId++; // 2
-    myPassFilter[ anId ]->SetInput( myGeomFilter->GetOutput() ); 
+    myPassFilter[ anId ]->SetInput( myTransformFilter->GetOutput() );
     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
 
     anId++; // 3
-    myTransformFilter->SetInput( myPassFilter[ anId ]->GetPolyDataOutput() );
+    myGeomFilter->SetInput( myPassFilter[ anId ]->GetOutput() );
 
     anId++; // 4
-    myPassFilter[ anId ]->SetInput( myTransformFilter->GetOutput() );
+    myPassFilter[ anId ]->SetInput( myGeomFilter->GetOutput() ); 
     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
 
     anId++; // 5
index 51c913585c4f634ce9504e93505e5091d4f38bd6..0f8f3e0bfded86c32ce53b33c09dff446971f29e 100644 (file)
@@ -37,7 +37,7 @@
 #include <vtkPointData.h>
 #include <vtkProperty2D.h>
 #include <vtkRenderer.h>
-#include <vtkPolyData.h>
+#include <vtkUnstructuredGrid.h>
 
 vtkStandardNewMacro(SMESH_NodeLabelActor);
 
@@ -47,7 +47,7 @@ vtkStandardNewMacro(SMESH_NodeLabelActor);
 SMESH_NodeLabelActor::SMESH_NodeLabelActor() {
   //Definition of points numbering pipeline
   //---------------------------------------
-  myPointsNumDataSet = vtkPolyData::New();
+  myPointsNumDataSet = vtkUnstructuredGrid::New();
 
   myPtsMaskPoints = vtkMaskPoints::New();
   myPtsMaskPoints->SetInput(myPointsNumDataSet);
@@ -114,7 +114,7 @@ SMESH_NodeLabelActor::~SMESH_NodeLabelActor() {
 
 void SMESH_NodeLabelActor::SetPointsLabeled(bool theIsPointsLabeled) {
   myTransformFilter->Update();
-  vtkDataSet* aGrid = vtkPolyData::SafeDownCast(myTransformFilter->GetOutput());
+  vtkDataSet* aGrid = vtkUnstructuredGrid::SafeDownCast(myTransformFilter->GetOutput());
 
   if(!aGrid)
     return;
@@ -124,7 +124,7 @@ void SMESH_NodeLabelActor::SetPointsLabeled(bool theIsPointsLabeled) {
   if ( myIsPointsLabeled )
   {
     myPointsNumDataSet->ShallowCopy(aGrid);
-    vtkDataSet *aDataSet = myPointsNumDataSet;
+    vtkUnstructuredGrid *aDataSet = myPointsNumDataSet;
     
     int aNbElem = aDataSet->GetNumberOfPoints();
     
@@ -133,7 +133,7 @@ void SMESH_NodeLabelActor::SetPointsLabeled(bool theIsPointsLabeled) {
     
     for ( vtkIdType anId = 0; anId < aNbElem; anId++ )
     {
-      int aSMDSId = GetNodeObjId( anId );
+      int aSMDSId = myVisualObj->GetNodeObjId( anId );
       anArray->SetValue( anId, aSMDSId );
     }
     
index 0a9ffaa70bc75a495a48f36787bd8fa5a1f9d0b6..406bad68296d30df4ab5d7bbb0a1ab49d45835d3 100644 (file)
@@ -32,7 +32,7 @@ class vtkSelectVisiblePoints;
 class vtkLabeledDataMapper;
 class vtkActor2D;
 class vtkMaskPoints;
-class vtkPolyData;
+class vtkUnstructuredGrid;
 
 
 class SMESHOBJECT_EXPORT SMESH_NodeLabelActor : public SMESH_DeviceActor {
@@ -63,7 +63,7 @@ protected:
   ~SMESH_NodeLabelActor();
   
   bool myIsPointsLabeled;
-  vtkPolyData* myPointsNumDataSet;
+  vtkUnstructuredGrid* myPointsNumDataSet;
   vtkActor2D *myPointLabels;
   vtkMaskPoints* myPtsMaskPoints;
   vtkLabeledDataMapper* myPtsLabeledDataMapper;