]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
To restore picking functionality for Gauss Points
authorapo <apo@opencascade.com>
Wed, 17 Jan 2007 08:15:22 +0000 (08:15 +0000)
committerapo <apo@opencascade.com>
Wed, 17 Jan 2007 08:15:22 +0000 (08:15 +0000)
src/CONVERTOR/VISUConvertor.cxx
src/CONVERTOR/VISU_Convertor_impl.cxx
src/CONVERTOR/VISU_Structures_impl.cxx
src/CONVERTOR/VISU_Structures_impl.hxx
src/OBJECT/VISU_GaussPtsAct.cxx
src/PIPELINE/VISUPipeLine.cxx

index 3a07760d01b741100c544e89073e08fb642ac6ef..c9517fdabe08b28a2bf26ff5f7aa952f07d8f075 100644 (file)
@@ -46,6 +46,7 @@ static int MYDEBUG = 1;
 static int MYDEBUG = 0;
 #endif
 
+//#define _DEBUG_ID_MAPPING_
 #define _DEXCEPT_
 
 typedef vtkUnstructuredGrid TOutput;
@@ -91,26 +92,32 @@ void parseFile(const char* theFileName)
            if(anEntity != VISU::NODE_ENTITY){
              VISU::PGaussPtsIDMapper aGaussMesh = 
                aCon->GetTimeStampOnGaussPts(aMeshName,anEntity,aFieldName,aTimeStamp);
-             /*
-             VISU::TVTKOutput* aDataSet = aGaussMesh->GetVTKOutput();
+#ifdef _DEBUG_ID_MAPPING_
+             vtkDataSet* aDataSet = aGaussMesh->GetOutput();
+             aDataSet->Update();
              int aNbCells = aDataSet->GetNumberOfCells();
+             cout<<"aNbCells = "<<aNbCells<<endl;
              for(int anCellId = 0; anCellId < aNbCells; anCellId++){
                VISU::TGaussPointID anObjID = aGaussMesh->GetObjID(anCellId);
-               cout<<anObjID.first<<"; "<<anObjID.second<<endl;
+               cout<<anObjID.first<<"; "<<anObjID.second<<"; "<<aGaussMesh->GetNodeVTKID(anObjID.first)<<endl;
+               vtkFloatingPointType* aCoord = aGaussMesh->GetNodeCoord(anCellId);
+               cout<<aCoord[0]<<"; "<<aCoord[1]<<"; "<<aCoord[2]<<endl;
              }
-             */
+#endif
            }else{
+             continue;
              VISU::PIDMapper anIDMapper = 
                aCon->GetTimeStampOnMesh(aMeshName,anEntity,aFieldName,aTimeStamp);
-             /*
-             VISU::TVTKOutput* aDataSet = anIDMapper->GetVTKOutput();
+#ifdef _DEBUG_ID_MAPPING_
+             vtkDataSet* aDataSet = anIDMapper->GetOutput();
+             aDataSet->Update();
              int aNbCells = aDataSet->GetNumberOfCells();
              for(int anCellId = 0; anCellId < aNbCells; anCellId++){
                int anObjID = anIDMapper->GetElemObjID(anCellId);
                int aVTKID  = anIDMapper->GetElemVTKID(anObjID);
                cout<<anObjID<<"; "<<aVTKID<<endl;
              }
-             */
+#endif
            }
            //goto OK;
          }
@@ -132,18 +139,16 @@ void parseFile(const char* theFileName)
       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
        const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
        VISU::PIDMapper anIDMapper = aCon->GetMeshOnEntity(aMeshName,anEntity);
-       {
-         /*
-         VISU::TVTKOutput* aDataSet = anIDMapper->GetVTKOutput();
-         int aNbCells, anCellId, anObjID, aVTKID;
-         aNbCells = aDataSet->GetNumberOfCells();
-         for(anCellId = 0; anCellId < aNbCells; anCellId++){
-           anObjID = anIDMapper->GetElemObjID(anCellId);
-           aVTKID  = anIDMapper->GetElemVTKID(anObjID);
-           cout<<anObjID<<"; "<<aVTKID<<endl;
-         }
-         */
+#ifdef _DEBUG_ID_MAPPING_
+       vtkDataSet* aDataSet = anIDMapper->GetOutput();
+       int aNbCells, anCellId, anObjID, aVTKID;
+       aNbCells = aDataSet->GetNumberOfCells();
+       for(anCellId = 0; anCellId < aNbCells; anCellId++){
+         anObjID = anIDMapper->GetElemObjID(anCellId);
+         aVTKID  = anIDMapper->GetElemVTKID(anObjID);
+         cout<<anObjID<<"; "<<aVTKID<<endl;
        }
+#endif
       }
 
       //Import families
index 323216f9b491fa457c46ee80314ed025b532f231..5054f0a4ed158058cb4bec1cbf54c3c4bc6aeac1 100644 (file)
@@ -905,11 +905,10 @@ VISU_Convertor_impl
       GetMeshOnProfile(aMesh, aVTKMeshOnEntity, aProfile);
 
       VISU::PGaussMeshImpl aGaussMesh = aValForTime->myGaussMesh;
-      VISU::TPolyDataHolder& aGaussPtsSource = aGaussMesh->mySource;
-      if(!aGaussPtsSource.myIsVTKDone){
+      if(!aGaussMesh->myIsVTKDone){
        BuildGaussMesh(aMesh, aVTKMeshOnEntity, aGaussMesh);
        aGaussMesh->myParent = aProfile.get();
-       aGaussPtsSource.myIsVTKDone = true;
+       aGaussMesh->myIsVTKDone = true;
       }
 
       aGaussPtsIDFilter->myIDMapper = aGaussMesh;
index 2a65b133e8bba79bff9771da0ce2a662b5871f1b..100bcf1493eb9cfe65a11b988fc3943798d6282b 100644 (file)
@@ -682,19 +682,12 @@ namespace VISU
 
     return aSubMeshImpl.GetObjID(anInputID,aStartId);
   }
-  
+
   vtkPolyData* 
   TGaussMeshImpl
   ::GetPolyDataOutput()
   {
-    return mySource.GetPolyDataOutput();
-  }
-
-  vtkDataSet* 
-  TGaussMeshImpl
-  ::GetOutput()
-  {
-    return GetPolyDataOutput();
+    return TAppendPolyDataHolder::GetPolyDataOutput();
   }
 
   unsigned long int
@@ -702,7 +695,6 @@ namespace VISU
   ::GetMemorySize()
   {
     size_t aSize = TAppendPolyDataHolder::GetMemorySize();
-    aSize += mySource.GetMemorySize();
     TGeom2GaussSubMesh::const_iterator anIter = myGeom2GaussSubMesh.begin();
     TGeom2GaussSubMesh::const_iterator anIterEnd = myGeom2GaussSubMesh.end();
     for(; anIter != anIterEnd; anIter++){
index dc2e7b9bd6ab5d6f619ab2bf8456329018a0a035..7783b66c9fd6f01fef9bd48a251723067e600cc6 100644 (file)
@@ -481,11 +481,6 @@ namespace VISU
     vtkPolyData* 
     GetPolyDataOutput();
 
-    //! Reimplement the TIDMapper::GetOutput
-    virtual
-    vtkDataSet* 
-    GetOutput();
-
     //! Gets memory size used by the instance (bytes).
     virtual
     unsigned long int
@@ -496,7 +491,6 @@ namespace VISU
     TNamedIDMapper*
     GetParent();
 
-    TPolyDataHolder mySource; //!< Keeps VTK representation of the Gauss Points
     TNamedIDMapper* myParent; //!< Refer to parent mesh
     TGaussSubMeshArr myGaussSubMeshArr; //!< Keeps sequence of TGaussSubMesh as they were added into TAppendFilterHolder
     TGeom2GaussSubMesh myGeom2GaussSubMesh; //!< Keeps TGaussSubMesh according to their geometrical type
index 4ad22030f115a098eeac9f2cd1d43dc2cf1895c8..dd2e4ab9b789058ba89eb756d48849aaa5d377ec 100644 (file)
 #include <vtkTextMapper.h>
 #include <vtkTextProperty.h>
 
-#include <vtkCellData.h>
 #include <vtkPointData.h>
-
 #include <vtkDataArray.h>
-#include <vtkFloatArray.h>
 
 #include <vtkSphereSource.h>
 #include <vtkPolyDataMapper.h>
@@ -714,8 +711,8 @@ VISU_GaussPtsAct
        if(anIsChanged){
          vtkFloatingPointType* aNodeCoord = GetNodeCoord(anObjId);
          vtkDataSet* aDataSet = GetInput();
-         vtkCellData* aCellData = aDataSet->GetCellData();
-         if(vtkDataArray *aScalarArray = aCellData->GetScalars()){
+         vtkPointData* aPointData = aDataSet->GetPointData();
+         if(vtkDataArray *aScalarArray = aPointData->GetScalars()){
            vtkFloatingPointType aPyramidHeight = myPickingSettings->GetPyramidHeight();
            aPyramidHeight = aPyramidHeight*myGaussPointsPL->GetMaxPointSize();
            //vtkFloatingPointType aColor[3] = myPreHighlightActor->GetProperty()->GetColor();
@@ -927,7 +924,7 @@ VISU_GaussPtsAct
     vtkFloatingPointType aWorldCoord[4] = {aNodeCoord[0], aNodeCoord[1], aNodeCoord[2], 1.};
     //
     vtkDataSet* aDataSet = GetInput();
-    vtkCellData* aDataSetAttributes = aDataSet->GetCellData();
+    vtkPointData* aDataSetAttributes = aDataSet->GetPointData();
     //
     if(vtkDataArray* aScalarArray = aDataSetAttributes->GetScalars()){
       vtkFloatingPointType aVal = aScalarArray->GetTuple1(aVtkId);
@@ -964,18 +961,19 @@ VISU_GaussPtsAct
     }
 
     if(vtkDataArray* aFieldArray = aDataSetAttributes->GetArray("VISU_FIELD")){
-      if(vtkFloatArray *aFloatArray = dynamic_cast<vtkFloatArray*>(aFieldArray)){
-       int aNbComp = aFloatArray->GetNumberOfComponents();
-       aStr<<"\nData: {";
-       int anId = 0;
-       while(anId < aNbComp){
-         vtkFloatingPointType aComp = aFloatArray->GetComponent(aVtkId,anId++);
-         aStr<<aComp;
-         if(anId < aNbComp)
-           aStr<<"; ";
-       }
-       aStr<<"}";
+      int aNbComp = aFieldArray->GetNumberOfComponents();
+      std::vector<vtkFloatingPointType> aTuple(aNbComp);
+      aFieldArray->GetTuple(aVtkId, &aTuple[0]);
+      
+      aStr<<"\nData: {";
+      int anId = 0;
+      while(anId < aNbComp){
+       vtkFloatingPointType aComp = aTuple[anId++];
+       aStr<<aComp;
+       if(anId < aNbComp)
+         aStr<<"; ";
       }
+      aStr<<"}";
     }
     //
     // myTextActor
index 6ad4721037e1629091f5b2e971af749349cbb7e4..e192a47e6bbc84ea0e21128feed6d6ea340a9969 100644 (file)
@@ -38,7 +38,7 @@
 
 #include "VISU_Convertor.hxx"
 
-typedef VISU_DeformedShapePL TPresent;
+typedef VISU_GaussPointsPL TPresent;
 
 #include <vtkUnstructuredGrid.h>
 #include <vtkDataSetMapper.h>
@@ -58,6 +58,7 @@ using namespace std;
 
 static int isOnlyMesh = false;
 
+//#define _DEBUG_ID_MAPPING_
 
 //----------------------------------------------------------------------------
 template<class TPipeLine>
@@ -100,6 +101,20 @@ CreateColoredPL<VISU_GaussPointsPL>(VISU_Convertor* theConvertor,
 
   aPresent->Update();
 
+#ifdef _DEBUG_ID_MAPPING_
+  vtkDataSet* aDataSet = aPresent->GetOutput();
+  aDataSet->Update();
+  int aNbCells = aDataSet->GetNumberOfCells();
+  cout<<"aNbCells = "<<aNbCells<<endl;
+  for(int anCellId = 0; anCellId < aNbCells; anCellId++){
+    vtkIdType anObjID = aPresent->GetNodeObjID(anCellId);
+    vtkIdType aVtkID = aPresent->GetNodeVTKID(anObjID);
+    cout<<anObjID<<"; "<<aVtkID<<"; - ";
+    vtkFloatingPointType* aCoord = aPresent->GetNodeCoord(anObjID);
+    cout<<aCoord[0]<<"; "<<aCoord[1]<<"; "<<aCoord[2]<<endl;
+  }
+#endif
+
   return aPresent;
 }
 
@@ -177,10 +192,10 @@ main(int argc, char** argv)
        VISU::TFieldMap::const_iterator aFieldMapIter = aFieldMap.begin();
        for(; aFieldMapIter != aFieldMap.end(); aFieldMapIter++){
          const VISU::PField aField = aFieldMapIter->second;
-         ///*
+         /*
          if(aField->myNbComp == 1) 
            continue;
-         //*/
+         */
          const string& aFieldName = aFieldMapIter->first;
          const VISU::TValField& aValField = aField->myValField;
          VISU::TValField::const_iterator aValFieldIter = aValField.begin();
@@ -195,7 +210,7 @@ main(int argc, char** argv)
                                                 aFieldName,
                                                 aTimeStamp);
          }else{
-           //continue;
+           continue;
            aPresent = CreateColoredPL<TPresent>(aConvertor,
                                                 aMeshName,
                                                 anEntity,