]> SALOME platform Git repositories - modules/gui.git/commitdiff
Salome HOME
PR: SMDS memory improvement, VTK version from Paraview 3.9
authorprascle <prascle>
Mon, 29 Nov 2010 12:51:19 +0000 (12:51 +0000)
committerprascle <prascle>
Mon, 29 Nov 2010 12:51:19 +0000 (12:51 +0000)
src/SVTK/Makefile.am
src/SVTK/SVTK_Actor.cxx
src/Session/SALOME_Session_Server.cxx
src/VTKViewer/VTKViewer_ExtractUnstructuredGrid.cxx
src/VTKViewer/VTKViewer_GeometryFilter.cxx

index 232e9ef8354ae17caf9e5ca925764cabbdc1a120..bae773e7353268e35c421969e495c8b4312d577d 100755 (executable)
@@ -162,7 +162,7 @@ libSVTK_la_LDFLAGS =                        \
 libSVTK_la_LIBADD = ../Qtx/libqtx.la ../SUIT/libsuit.la ../ViewerTools/libViewerTools.la \
                    ../OBJECT/libSalomeObject.la ../Prs/libSalomePrs.la \
                    ../VTKViewer/libVTKViewer.la ../OpenGLUtils/libOpenGLUtils.la \
-                   -lSALOMELocalTrace
+                   -lSALOMELocalTrace -lOpUtil
 
 # Executable
 bin_PROGRAMS = SVTK
index 66d71c86673f86e4623903be73d9b2e31268bb9e..e513649b5f1ab95968f4656bd5d5fe443aae1db6 100644 (file)
 #include <vtkRenderer.h>
 
 #include <vtkCell.h>
+#include <vtkPolyhedron.h>
 #include <vtkPolyData.h>
 
+#include "Utils_SALOME_Exception.hxx"
+#include "utilities.h"
+
 using namespace std;
 
 
@@ -125,7 +129,18 @@ SVTK_Actor
   for(int ind = 1; ind <= aNbOfParts; ind++){
     int aPartId = theMapIndex( ind );
     if(vtkCell* aCell = theMapActor->GetElemCell(aPartId))
-      myUnstructuredGrid->InsertNextCell(aCell->GetCellType(),aCell->GetPointIds());
+      {
+      if (aCell->GetCellType() != VTK_POLYHEDRON)
+        myUnstructuredGrid->InsertNextCell(aCell->GetCellType(),aCell->GetPointIds());
+      else
+        {
+          vtkPolyhedron *polyhedron = dynamic_cast<vtkPolyhedron*>(aCell);
+          if (!polyhedron)
+            throw SALOME_Exception(LOCALIZED ("not a polyhedron"));
+          vtkIdType *pts = polyhedron->GetFaces();
+          myUnstructuredGrid->InsertNextCell(aCell->GetCellType(),pts[0], pts+1);
+        }
+      }
   }
 
   UnShrink();
index d3d61e2e12449a752ca7b141c3b8f6f3818588b5..507c399921759d796b8326bf0bcca2a24caad040 100755 (executable)
@@ -96,7 +96,7 @@ void MessageOutput( QtMsgType type, const char* msg )
   switch ( type )
   {
   case QtDebugMsg:
-    MESSAGE( "Debug: " << msg );
+    //MESSAGE( "Debug: " << msg );
     break;
   case QtWarningMsg:
     MESSAGE( "Warning: " << msg );
index d2ee28678a2afc8a2f20ce2f0bcc2fc667540c78..2ce0a7165c8ff53f605e6097abfb027d899538e2 100755 (executable)
@@ -37,6 +37,8 @@
 #include <vtkInformation.h>
 #include <vtkInformationVector.h>
 
+#include "utilities.h"
+
 using namespace std;
 
 #if defined __GNUC__
@@ -66,6 +68,7 @@ void VTKViewer_ExtractUnstructuredGrid::RegisterCell(vtkIdType theCellId){
 void VTKViewer_ExtractUnstructuredGrid::RegisterCellsWithType(vtkIdType theCellType){
 //  if(0 && MYDEBUG) MESSAGE("RegisterCellsWithType - theCellType = "<<theCellType);
   myCellTypes.insert(theCellType);
+  //MESSAGE("myCellTypes.insert " << theCellType);
   Modified();
 }
 
@@ -96,6 +99,8 @@ vtkIdType VTKViewer_ExtractUnstructuredGrid::GetOutputId(int theInId) const{
 inline void InsertCell(vtkUnstructuredGrid *theInput,
                        vtkCellArray *theConnectivity, 
                        vtkUnsignedCharArray* theCellTypesArray,
+                       vtkIdTypeArray*& theFaces,
+                       vtkIdTypeArray*& theFaceLocations,
                        vtkIdType theCellId, 
                        vtkIdList *theIdList,
                        bool theStoreMapping,
@@ -110,9 +115,40 @@ inline void InsertCell(vtkUnstructuredGrid *theInput,
   for(vtkIdType i = 0; i < aNbIds; i++){
     theIdList->SetId(i,aPntIds->GetId(i));
   }
-  theConnectivity->InsertNextCell(theIdList);
-
   vtkIdType aCellType = aCell->GetCellType();
+  if (aCellType != VTK_POLYHEDRON)
+    {
+      theConnectivity->InsertNextCell(theIdList);
+      if (theFaceLocations)
+        theFaceLocations->InsertNextValue(-1);
+    }
+  else
+    {
+      //MESSAGE("InsertCell type VTK_POLYHEDRON " << theStoreMapping);
+      if (!theFaces)
+        {
+          theFaces = vtkIdTypeArray::New();
+          theFaces->Allocate(theCellTypesArray->GetSize());
+          theFaceLocations = vtkIdTypeArray::New();
+          theFaceLocations->Allocate(theCellTypesArray->GetSize());
+          // FaceLocations must be padded until the current position
+          for (vtkIdType i = 0; i <= theCellTypesArray->GetMaxId(); i++)
+            {
+              theFaceLocations->InsertNextValue(-1);
+            }
+        }
+      // insert face location
+      theFaceLocations->InsertNextValue(theFaces->GetMaxId() + 1);
+
+      // insert cell connectivity and faces stream
+      vtkIdType nfaces;
+      vtkIdType* face;
+      vtkIdType realnpts;
+      theInput->GetFaceStream(theCellId, nfaces, face);
+      vtkUnstructuredGrid::DecomposeAPolyhedronCell(
+          nfaces, face, realnpts, theConnectivity, theFaces);
+    }
+
   theCellTypesArray->InsertNextValue(aCellType);
   if(theStoreMapping){
     theOut2InId.push_back(theCellId);
@@ -122,7 +158,7 @@ inline void InsertCell(vtkUnstructuredGrid *theInput,
 
 inline void InsertPointCell(vtkCellArray *theConnectivity, 
                             vtkUnsignedCharArray* theCellTypesArray,
-                            vtkIdType theCellId, 
+                            vtkIdType theCellId,
                             vtkIdList *theIdList,
                             bool theStoreMapping,
                             vtkIdType theOutId, 
@@ -190,19 +226,23 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
       aCellTypesArray->SetNumberOfComponents(1);
       aCellTypesArray->Allocate(aNbElems*aCellTypesArray->GetNumberOfComponents());
+
+      vtkIdTypeArray *newFaces = 0;
+      vtkIdTypeArray *newFaceLocations = 0;
+
       if(!myCellIds.empty() && myCellTypes.empty()){
         if(myStoreMapping) myOut2InId.reserve(myCellIds.size());
         if(myChangeMode == eAdding){
           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
             if(myCellIds.find(aCellId) != myCellIds.end()){
-              InsertCell(anInput,aConnectivity,aCellTypesArray,aCellId,anIdList,
+              InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
                          myStoreMapping,anOutId,myOut2InId,myIn2OutId);
             }
           }
         }else{
           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
             if(myCellIds.find(aCellId) == myCellIds.end()){
-              InsertCell(anInput,aConnectivity,aCellTypesArray,aCellId,anIdList,
+              InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
                          myStoreMapping,anOutId,myOut2InId,myIn2OutId);
             }
           }
@@ -212,7 +252,7 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
             vtkIdType aType = anInput->GetCellType(aCellId);
             if(myCellTypes.find(aType) != myCellTypes.end()){
-              InsertCell(anInput,aConnectivity,aCellTypesArray,aCellId,anIdList,
+              InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
                          myStoreMapping,anOutId,myOut2InId,myIn2OutId);
             }
           }
@@ -220,7 +260,7 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
             vtkIdType aType = anInput->GetCellType(aCellId);
             if(myCellTypes.find(aType) == myCellTypes.end()){
-              InsertCell(anInput,aConnectivity,aCellTypesArray,aCellId,anIdList,
+              InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
                          myStoreMapping,anOutId,myOut2InId,myIn2OutId);
             }
           }
@@ -231,7 +271,7 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
             vtkIdType aType = anInput->GetCellType(aCellId);
             if(myCellTypes.find(aType) != myCellTypes.end()){
               if(myCellIds.find(aCellId) != myCellIds.end()){
-                InsertCell(anInput,aConnectivity,aCellTypesArray,aCellId,anIdList,
+                InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
                            myStoreMapping,anOutId,myOut2InId,myIn2OutId);
               }
             }
@@ -241,7 +281,7 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
             vtkIdType aType = anInput->GetCellType(aCellId);
             if(myCellTypes.find(aType) == myCellTypes.end()){
               if(myCellIds.find(aCellId) == myCellIds.end()){
-                InsertCell(anInput,aConnectivity,aCellTypesArray,aCellId,anIdList,
+                InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
                            myStoreMapping,anOutId,myOut2InId,myIn2OutId);
               }
             }
@@ -256,7 +296,7 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
         for(vtkIdType i = 0, *pts, npts; aConnectivity->GetNextCell(npts,pts); i++){
           aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
         }
-        anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
+        anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity,newFaceLocations,newFaces);
         anOutput->SetPoints(anInput->GetPoints());
         aCellLocationsArray->Delete();
       }
@@ -348,7 +388,7 @@ void VTKViewer_ExtractUnstructuredGrid::Execute()
       for(vtkIdType i = 0, *pts, npts; aConnectivity->GetNextCell(npts,pts); i++){
         aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
       }
-      anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
+      anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity,0, 0);
       anOutput->SetPoints(anInput->GetPoints());
       aCellLocationsArray->Delete();
     }
index ae44393d4e170814106bb27e205d6a97c5ec8824..0052072e9562d837c7efb95797f65d7735bd25a5 100755 (executable)
@@ -55,6 +55,8 @@
 #include <map>
 #include <set>
 
+#include "utilities.h"
+
 #if defined __GNUC__
   #if __GNUC__ == 2
     #define __GNUC_2__
@@ -480,6 +482,53 @@ VTKViewer_GeometryFilter
             }
           break;
         }
+
+        case VTK_POLYHEDRON:
+          {
+            //MESSAGE("VTK_POLYHEDRON geometry filter");
+            vtkIdType nFaces = 0;
+            vtkIdType* ptIds = 0;
+            int idp = 0;
+            input->GetFaceStream(cellId, nFaces, ptIds);
+            for (faceId = 0; faceId < nFaces; faceId++)
+              {
+                faceIds->Reset();
+                numFacePts = ptIds[idp];
+                //MESSAGE("numFacePts="<< numFacePts);
+                int pt0 = ++idp;
+                for (int i = 0; i < numFacePts; i++)
+                  {
+                    //MESSAGE("ptIds[" << idp + i << "]=" << ptIds[idp + i]);
+                    faceIds->InsertNextId(ptIds[idp + i]);
+                  }
+                idp += numFacePts;
+                switch (numFacePts)
+                  {
+                  case 3:
+                    aCellType = VTK_TRIANGLE;
+                    break;
+                  case 4:
+                    aCellType = VTK_QUAD;
+                    break;
+                  default:
+                    aCellType = VTK_POLYGON;
+                    break;
+                  }
+                // TODO understand and fix display of several polyhedrons
+                input->GetCellNeighbors(cellId, faceIds, cellIds);
+                if (cellIds->GetNumberOfIds() <= 0 || myShowInside
+                    || (!allVisible && !cellVis[cellIds->GetId(0)]))
+                  {
+                    for (i = 0; i < numFacePts; i++)
+                      aNewPts[i] = ptIds[pt0 + i];
+                    newCellId = output->InsertNextCell(aCellType, numFacePts, aNewPts);
+                    if (myStoreMapping)
+                      myVTK2ObjIds.push_back(cellId);
+                    outputCD->CopyData(cd, cellId, newCellId);
+                  }
+              }
+            break;
+          }
         //Quadratic cells
         case VTK_QUADRATIC_EDGE:
         case VTK_QUADRATIC_TRIANGLE: