]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
To introduce the common functionality which can be used for mapping ids purposes
authorapo <apo@opencascade.com>
Wed, 28 Mar 2007 07:11:11 +0000 (07:11 +0000)
committerapo <apo@opencascade.com>
Wed, 28 Mar 2007 07:11:11 +0000 (07:11 +0000)
(In context of fix for Bug NPAL15278 - "EDF 347 : ScalarMaponDeformedShape")

src/CONVERTOR/VISU_ConvertorUtils.cxx
src/CONVERTOR/VISU_ConvertorUtils.hxx

index a926628fa340fe2d3652b4cd2a41a62bf8643f9e..d1334e0213c903e89fe484c05b3ac343725aadf1 100644 (file)
@@ -32,6 +32,9 @@
 #include <vtkCellData.h>
 #include <vtkDataSet.h>
 
+#include <vtkIntArray.h>
+#include <algorithm>
+
 #ifdef _DEBUG_
 static int MYDEBUG = 0;
 #else
@@ -41,6 +44,7 @@ static int MYDEBUG = 0;
 namespace VISU
 {
 
+  //---------------------------------------------------------------
   void 
   WriteToFile(vtkUnstructuredGrid* theDataSet, const std::string& theFileName)
   {
@@ -53,26 +57,99 @@ namespace VISU
   }
 
 
+  //---------------------------------------------------------------
   bool 
   IsDataOnPoints(vtkDataSet* theDataSet)
   {
     theDataSet->Update();
-    if(theDataSet->GetPointData()->GetArray("VISU_POINTS_MAPPER"))
-      return theDataSet->GetPointData()->GetNumberOfArrays() > 1;
-    return theDataSet->GetPointData()->GetNumberOfArrays() > 0;
+    vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
+    if(aDataSetAttributes->GetArray("VISU_POINTS_MAPPER"))
+      return aDataSetAttributes->GetNumberOfArrays() > 1;
+    return aDataSetAttributes->GetNumberOfArrays() > 0;
   }
 
 
+  //---------------------------------------------------------------
   bool 
   IsDataOnCells(vtkDataSet* theDataSet)
   {
     theDataSet->Update();
-    if(theDataSet->GetCellData()->GetArray("VISU_CELLS_MAPPER"))
-      return theDataSet->GetCellData()->GetNumberOfArrays() > 1;
-    return theDataSet->GetCellData()->GetNumberOfArrays() > 0;
+    vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
+    if(aDataSetAttributes->GetArray("VISU_CELLS_MAPPER"))
+      return aDataSetAttributes->GetNumberOfArrays() > 1;
+    return aDataSetAttributes->GetNumberOfArrays() > 0;
+  }
+
+
+  //---------------------------------------------------------------
+  vtkIdType
+  GetElemVTKID(vtkDataSet *theDataSet, vtkIdType theID)
+  {
+    theDataSet->Update();
+    vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
+    if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_CELLS_MAPPER")){
+      if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
+       int aMaxId = anIntArray->GetMaxId();
+       int* aPointer = anIntArray->GetPointer(0);
+       int* anEndPointer = anIntArray->GetPointer(aMaxId + 1);
+       int* aPtr = std::find(aPointer, anEndPointer, theID);
+       return aPtr - aPointer;
+      }
+    }
+    return -1;
+  }
+
+
+  //---------------------------------------------------------------
+  vtkIdType
+  GetElemObjID(vtkDataSet *theDataSet, vtkIdType theID)
+  {
+    theDataSet->Update();
+    vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
+    if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_CELLS_MAPPER")){
+      if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
+       return anIntArray->GetValue(theID);
+      }
+    }
+    return -1;
+  }
+
+
+  //---------------------------------------------------------------
+  vtkIdType
+  GetNodeVTKID(vtkDataSet *theDataSet, vtkIdType theID)
+  {
+    theDataSet->Update();
+    vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
+    if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_POINTS_MAPPER")){
+      if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
+       int aMaxId = anIntArray->GetMaxId();
+       int* aPointer = anIntArray->GetPointer(0);
+       int* anEndPointer = anIntArray->GetPointer(aMaxId + 1);
+       int* aPtr = std::find(aPointer, anEndPointer, theID);
+       return aPtr - aPointer;
+      }
+    }
+    return -1;
+  }
+
+
+  //---------------------------------------------------------------
+  vtkIdType
+  GetNodeObjID(vtkDataSet *theDataSet, vtkIdType theID)
+  {
+    theDataSet->Update();
+    vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
+    if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_POINTS_MAPPER")){
+      if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
+       return anIntArray->GetValue(theID);
+      }
+    }
+    return -1;
   }
 
 
+  //---------------------------------------------------------------
   TTimerLog
   ::TTimerLog(int theIsDebug,
              const std::string& theName):
@@ -85,6 +162,7 @@ namespace VISU
     BEGMSG(myIsDebug > 1,"{\n");
   }
 
+  //---------------------------------------------------------------
   TTimerLog
   ::~TTimerLog()
   {
index 2a86f810b01bf2277effaaf664016a779f0f299f..6c952733c00c0835fe77381e70ddf4942fe3b2ac 100644 (file)
@@ -30,6 +30,7 @@
 #include <string>
 
 #include <vtkCellType.h>
+#include <vtkSystemIncludes.h>
 
 #include "MED_Utilities.hxx"
 
@@ -54,6 +55,18 @@ namespace VISU
   bool 
   IsDataOnPoints(vtkDataSet* theDataSet);
 
+  vtkIdType
+  GetElemVTKID(vtkDataSet *theDataSet, vtkIdType theID);
+
+  vtkIdType
+  GetElemObjID(vtkDataSet *theDataSet, vtkIdType theID);
+
+  vtkIdType
+  GetNodeVTKID(vtkDataSet *theDataSet, vtkIdType theID);
+
+  vtkIdType
+  GetNodeObjID(vtkDataSet *theDataSet, vtkIdType theID);
+
   class TTimerLog
   {
     int myIsDebug;