]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
0021459: EDF 1495 SMESH: Manipulation of discrete elements with attributes
authoreap <eap@opencascade.com>
Thu, 19 Jul 2012 13:14:28 +0000 (13:14 +0000)
committereap <eap@opencascade.com>
Thu, 19 Jul 2012 13:14:28 +0000 (13:14 +0000)
+  void AllocateDiameters( vtkIdType maxVtkID );
+  void SetBallDiameter( vtkIdType vtkID, double diameter );
+  double GetBallDiameter( vtkIdType vtkID ) const;

src/SMDS/SMDS_UnstructuredGrid.cxx
src/SMDS/SMDS_UnstructuredGrid.hxx

index 099b8cbe85235b332d2cdfb89df257bde51bfde0..0d9b1e6c3c458775feb84294a542778c70e26bf4 100644 (file)
@@ -27,7 +27,9 @@
 #include "utilities.h"
 
 #include <vtkCellArray.h>
+#include <vtkCellData.h>
 #include <vtkCellLinks.h>
+#include <vtkDoubleArray.h>
 #include <vtkIdTypeArray.h>
 #include <vtkUnsignedCharArray.h>
 
@@ -218,10 +220,12 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
       ++i;
     int endBloc = i;
     if ( endBloc > startBloc )
-      copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell,
-               alreadyCopied, startBloc, endBloc);
+      copyBloc(newTypes,
+               idCellsOldToNew, idNodesOldToNew,
+               newConnectivity, newLocations,
+               pointsCell, alreadyCopied,
+               startBloc, endBloc);
   }
-
   newConnectivity->Squeeze();
 
   if (1/*newNodeSize*/)
@@ -231,6 +235,19 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
       MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
     }
 
+  if (vtkDoubleArray* diameters =
+      vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
+  {
+    for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
+    {
+      if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
+        continue;
+      int newCellId = idCellsOldToNew[ oldCellID ];
+      if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
+        diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
+    }
+  }
+
   if (this->FaceLocations)
     {
       vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
@@ -275,7 +292,9 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
       newFaces->Delete();
     }
   else
+  {
     this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
+  }
 
   newPoints->Delete();
   newTypes->Delete();
@@ -299,10 +318,15 @@ void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& id
     }
 }
 
-void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew,
-                                     std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity,
-                                     vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
-                                     int start, int end)
+void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
+                                     std::vector<int>&     idCellsOldToNew,
+                                     std::vector<int>&     idNodesOldToNew,
+                                     vtkCellArray*         newConnectivity,
+                                     vtkIdTypeArray*       newLocations,
+                                     vtkIdType*            pointsCell,
+                                     int&                  alreadyCopied,
+                                     int                   start,
+                                     int                   end)
 {
   MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
   for (int j = start; j < end; j++)
@@ -1072,3 +1096,50 @@ SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
   // TODO update sub-shape list of elements and nodes
   return 0;
 }
+
+//================================================================================
+/*!
+ * \brief Allocates data array for ball diameters
+ *  \param MaxVtkID - max ID of a ball element
+ */
+//================================================================================
+
+void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
+{
+  SetBallDiameter( MaxVtkID, 0 );
+}
+
+//================================================================================
+/*!
+ * \brief Sets diameter of a ball element
+ *  \param vtkID - vtk id of the ball element
+ *  \param diameter - diameter of the ball element
+ */
+//================================================================================
+
+void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
+{
+  vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
+  if ( !array )
+  {
+    array = vtkDoubleArray::New();
+    array->SetNumberOfComponents(1);
+    vtkDataSet::CellData->SetScalars( array );
+  }
+  array->InsertValue( vtkID, diameter );
+}
+
+//================================================================================
+/*!
+ * \brief Returns diameter of a ball element
+ *  \param vtkID - vtk id of the ball element
+ */
+//================================================================================
+
+double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
+{
+  if ( vtkDataSet::CellData )
+    return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );
+  return 0;
+}
+
index 6cf909263b5f0816926342962623dfc3d0ee8246..4c1cc7b80934ba768daedd1f73e4287421c11831 100644 (file)
@@ -66,9 +66,10 @@ class SMDS_EXPORT SMDS_UnstructuredGrid: public vtkUnstructuredGrid
 {
 public:
   void setSMDS_mesh(SMDS_Mesh *mesh);
-  void compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize, std::vector<int>& idCellsOldToNew,
-                   int newCellSize);
-
+  void compactGrid(std::vector<int>& idNodesOldToNew,
+                   int               newNodeSize,
+                   std::vector<int>& idCellsOldToNew,
+                   int               newCellSize);
   virtual unsigned long GetMTime();
   virtual void Update();
   virtual void UpdateInformation();
@@ -89,7 +90,8 @@ public:
   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
   int getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes);
   void BuildLinks();
-  SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2, std::set<int>& originalNodes,
+  SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2,
+                                       std::set<int>&                      originalNodes,
                                        std::map<int, std::map<int, int> >& nodeDomains,
                                        std::map<int, std::map<long,int> >& nodeQuadDomains);
   vtkCellLinks* GetLinks()
@@ -100,6 +102,10 @@ public:
   {
     return _downArray[vtkType];
   }
+  void AllocateDiameters( vtkIdType maxVtkID );
+  void SetBallDiameter( vtkIdType vtkID, double diameter );
+  double GetBallDiameter( vtkIdType vtkID ) const;
+
   static SMDS_UnstructuredGrid* New();
   SMDS_Mesh *_mesh;
 protected: