]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
To introduce ids mapping for mesh nodes
authorapo <apo@opencascade.com>
Wed, 11 Apr 2007 08:58:07 +0000 (08:58 +0000)
committerapo <apo@opencascade.com>
Wed, 11 Apr 2007 08:58:07 +0000 (08:58 +0000)
(In context of fix for Bug NPAL14169 - EDF VISU 256 : No continuum on a scalarmap)

src/CONVERTOR/VISU_AppendFilter.cxx
src/CONVERTOR/VISU_AppendFilter.hxx
src/CONVERTOR/VISU_Convertor_impl.cxx
src/CONVERTOR/VISU_Convertor_impl.hxx
src/PIPELINE/VISU_PrsMergerPL.cxx

index f5e9993995050ab0a2baaebba787a1acd1e182ea..b99e7d869c9222f239cf89e19e640e864d79616d 100644 (file)
@@ -57,21 +57,21 @@ VISU_AppendFilter
 
 void
 VISU_AppendFilter
-::SetPoints(vtkPoints* thePoints)
+::SetSharedPointsDataSet(vtkPointSet* thePointsDataSet)
 {
-  if(GetPoints() == thePoints)
+  if(GetSharedPointsDataSet() == thePointsDataSet)
     return;
 
-  myPoints = thePoints;
+  mySharedPointsDataSet = thePointsDataSet;
 
   Modified();
 }
 
-vtkPoints*
+vtkPointSet*
 VISU_AppendFilter
-::GetPoints()
+::GetSharedPointsDataSet()
 {
-  return myPoints.GetPointer();
+  return mySharedPointsDataSet.GetPointer();
 }
 
 void
@@ -219,8 +219,9 @@ void
 VISU_AppendFilter
 ::Execute()
 {
-  if(myPoints.GetPointer()){
-    if(myPoints->GetNumberOfPoints() < 1)
+  if(GetSharedPointsDataSet()){
+    vtkPoints* aPoints = GetSharedPointsDataSet()->GetPoints();
+    if(aPoints->GetNumberOfPoints() < 1)
       return;
   
     vtkUnstructuredGrid *anOutput = this->GetOutput(); 
@@ -241,7 +242,8 @@ VISU_AppendFilter
       
       // Append each input dataset together
       // 1.points
-      anOutput->SetPoints(myPoints.GetPointer());
+      anOutput->SetPoints(GetSharedPointsDataSet()->GetPoints());
+      anOutput->GetPointData()->PassData(GetSharedPointsDataSet()->GetPointData());
       
       // 2.cells
       vtkIdList *anIdList = vtkIdList::New(); 
@@ -280,7 +282,8 @@ VISU_AppendFilter
       
       // Append each input dataset together
       // 1.points
-      anOutput->SetPoints(myPoints.GetPointer());
+      anOutput->SetPoints(GetSharedPointsDataSet()->GetPoints());
+      anOutput->GetPointData()->PassData(GetSharedPointsDataSet()->GetPointData());
       
       // 2.cells
       vtkIdList *anIdList = vtkIdList::New(); 
index 91c7674c07e8a9268150255aac0ea9edd5f79fe2..fd9a827dfbc21d897d6f2b1337ee028fe166411d 100644 (file)
@@ -22,7 +22,7 @@
 #include <vtkAppendFilter.h>
 #include <vtkSmartPointer.h>
 
-class vtkPoints;
+class vtkPointSet;
 
 /*! \brief This class used same as vtkAppendFilter. See documentation on VTK for more information.
  */
@@ -39,10 +39,10 @@ public:
   vtkTypeRevisionMacro(VISU_AppendFilter, vtkAppendFilter);
 
   void
-  SetPoints(vtkPoints* thePoints);
+  SetSharedPointsDataSet(vtkPointSet* thePointsDataSet);
 
-  vtkPoints*
-  GetPoints();
+  vtkPointSet*
+  GetSharedPointsDataSet();
 
   void
   SetMergingInputs(bool theIsMergingInputs);
@@ -65,7 +65,7 @@ protected:
   virtual void Execute();
 
   bool myIsMergingInputs;
-  vtkSmartPointer<vtkPoints> myPoints;
+  vtkSmartPointer<vtkPointSet> mySharedPointsDataSet;
 };
 
 #endif
index aa52822c4a54bfd369db3cf7f3d019a30087aa40..1a91c358a37ad2bad568cdf67a19a1bbfe232f05 100644 (file)
@@ -252,10 +252,14 @@ namespace VISU
 
   //---------------------------------------------------------------
   TMeshImpl::TMeshImpl():
-    myPoints(vtkPoints::New()),
+    myPointsSource(vtkUnstructuredGrid::New()),
     myNbPoints(0) 
   {
-    myPoints->Delete();
+    vtkPoints* aPoints = vtkPoints::New();
+    myPointsSource->SetPoints(aPoints);
+    aPoints->Delete();
+
+    myPointsSource->Delete();
   }
 
 
@@ -1151,13 +1155,12 @@ namespace
   
 
   //---------------------------------------------------------------
-  vtkPoints*
-  GetPoints(const PMeshImpl& theMesh) 
+  vtkUnstructuredGrid*
+  GetPointsSource(const PMeshImpl& theMesh) 
   {
-    TVTKPoints& aPoints = theMesh->myPoints;
-    const TNamedPointCoords& aCoords = theMesh->myNamedPointCoords;
-
+    TVTKSource& aPointSource = theMesh->myPointsSource;
     if(!theMesh->myIsVTKDone){
+      const TNamedPointCoords& aCoords = theMesh->myNamedPointCoords;
       TCoordHelperPtr aCoordHelperPtr;
       {
        int aMeshDimension = theMesh->myDim;
@@ -1197,10 +1200,11 @@ namespace
        }
       }
 
+      vtkPoints* aPoints = aPointSource->GetPoints();
       vtkIdType aNbPoints = aCoords.GetNbPoints();
       aPoints->SetNumberOfPoints(aNbPoints);
       
-      INITMSG(MYDEBUG,"GetPoints - aNbPoints = "<<aNbPoints<<
+      INITMSG(MYDEBUG,"GetPointsSource - aNbPoints = "<<aNbPoints<<
              "; aDim = "<<theMesh->myDim<<
              endl);
 
@@ -1211,13 +1215,25 @@ namespace
                          aCoordHelperPtr->GetCoord(aCoordSlice,eY),
                          aCoordHelperPtr->GetCoord(aCoordSlice,eZ));
       }
-      
+
+      vtkIdType aNbTuples = aNbPoints;
+      vtkIntArray *aDataArray = vtkIntArray::New();
+      aDataArray->SetName("VISU_POINTS_MAPPER");
+      aDataArray->SetNumberOfComponents(1);
+      aDataArray->SetNumberOfTuples(aNbTuples);
+      for(vtkIdType aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
+       vtkIdType anObjID = aCoords.GetObjID(aTupleId);
+       aDataArray->SetValue(aTupleId, anObjID);
+       }
+      aPointSource->GetPointData()->AddArray(aDataArray);
+      aDataArray->Delete();
+
       theMesh->myIsVTKDone = true;
 
       if(MYVTKDEBUG) aPoints->DebugOn();
     }
 
-    return aPoints.GetPointer();
+    return aPointSource.GetPointer();
   }
   
   
@@ -1568,7 +1584,7 @@ namespace
     if(theSubProfile->myIsVTKDone)
       return;
     
-    aSource->SetPoints(GetPoints(theMesh));
+    aSource->ShallowCopy(GetPointsSource(theMesh));
     INITMSGA(MYDEBUG,0,"GetNumberOfPoints - "<<aSource->GetNumberOfPoints()<<endl);
     GetCells(aSource,theSubProfile,theProfile,theMeshOnEntity);
     BEGMSG(MYDEBUG,"GetNumberOfCells - "<<aSource->GetNumberOfCells()<<endl);
@@ -1592,7 +1608,7 @@ namespace
       return true;
    
     const TVTKAppendFilter& anAppendFilter = theProfile->GetFilter();
-    anAppendFilter->SetPoints(GetPoints(theMesh));
+    anAppendFilter->SetSharedPointsDataSet(GetPointsSource(theMesh));
 
     if(theProfile->myIsAll){
       TVTKOutput* aDataSet = theMeshOnEntity->GetVTKOutput();
@@ -1971,7 +1987,7 @@ VISU_Convertor_impl
       if(MYVTKDEBUG) anAppendFilter->DebugOn();
 
       LoadMeshOnEntity(aMesh,aMeshOnEntity);
-      anAppendFilter->SetPoints(GetPoints(aMesh));
+      anAppendFilter->SetSharedPointsDataSet(GetPointsSource(aMesh));
       
       const TGeom2SubMesh& aGeom2SubMesh = aMeshOnEntity->myGeom2SubMesh;
       TGeom2SubMesh::const_iterator anIter = aGeom2SubMesh.begin();
@@ -1989,7 +2005,7 @@ VISU_Convertor_impl
        //ENK: 23.11.2006
 
        const TVTKSource& aSource = aSubMesh->GetSource();
-       aSource->SetPoints(GetPoints(aMesh));
+       aSource->ShallowCopy(GetPointsSource(aMesh));
 
        aSubMesh->myStartID = aCellID;
        GetCellsOnSubMesh(aSource, aMeshOnEntity, aSubMesh, aVGeom);
@@ -2060,7 +2076,7 @@ VISU_Convertor_impl
       GetMeshOnEntity(theMeshName,theEntity);
 
       LoadFamilyOnEntity(aMesh,aMeshOnEntity,aFamily);
-      aSource->SetPoints(GetPoints(aMesh));
+      aSource->ShallowCopy(GetPointsSource(aMesh));
       GetCellsOnFamily(aSource,aMeshOnEntity,aFamily);
 
       aFamily->myNamedPointCoords = aMesh->myNamedPointCoords;
@@ -2122,7 +2138,7 @@ VISU_Convertor_impl
       const VISU::TFamilySet& aFamilySet = aGroup->myFamilySet;
 
       LoadMeshOnGroup(aMesh,aFamilySet);
-      anAppendFilter->SetPoints(GetPoints(aMesh));
+      anAppendFilter->SetSharedPointsDataSet(GetPointsSource(aMesh));
 
       TFamilySet::const_iterator anIter = aFamilySet.begin();
 
index f9f6b15722a43002db88357675a85667fdc16c41..fd1aaf052683466d28cba6a0ca6097e2237d9901 100644 (file)
@@ -279,7 +279,7 @@ namespace VISU
   {
     PNamedPointCoords myNamedPointCoords; //!< Keeps intermediate representation of the nodes
 
-    TVTKPoints myPoints; //!< Keeps VTK representation of the nodes
+    TVTKSource myPointsSource; //!< Keeps VTK representation of the nodes
     vtkIdType myNbPoints; //!< Keeps number of the nodes
 
     TMeshImpl();
index 529dc3a2d659c2ce5d55747fdf342ed5f9271112..f9d21555570a00b7a0c9cf40c41a8522d13ec9fe 100644 (file)
@@ -310,7 +310,7 @@ VISU_PrsMergerPL
     VISU::TVTKOutput* aScalarsOutput = aScalarsIDMapper->GetVTKOutput();
 
     // copy points to output from input scalar map
-    myAppendFilter->SetPoints(aScalarsOutput->GetPoints());
+    myAppendFilter->SetSharedPointsDataSet(aScalarsOutput);
       
     vtkIdType aNbGeoms = this->GetNbGeometry();
     for(vtkIdType aGeomId = 0; aGeomId < aNbGeoms; aGeomId++){