Salome HOME
PR: synchro V7_main tag mergefrom_V6_main_06Mar13
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
index 099b8cbe85235b332d2cdfb89df257bde51bfde0..d3ea0b6209b38f50eb2cf7eec3763b69caecc765 100644 (file)
@@ -17,7 +17,6 @@
 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 
-#define CHRONODEF
 #include "SMDS_UnstructuredGrid.hxx"
 #include "SMDS_Mesh.hxx"
 #include "SMDS_MeshInfo.hxx"
@@ -27,7 +26,9 @@
 #include "utilities.h"
 
 #include <vtkCellArray.h>
+#include <vtkCellData.h>
 #include <vtkCellLinks.h>
+#include <vtkDoubleArray.h>
 #include <vtkIdTypeArray.h>
 #include <vtkUnsignedCharArray.h>
 
@@ -42,14 +43,14 @@ SMDS_CellLinks* SMDS_CellLinks::New()
   return new SMDS_CellLinks();
 }
 
-vtkCellLinks::Link* SMDS_CellLinks::ResizeL(vtkIdType sz)
+void SMDS_CellLinks::ResizeForPoint(vtkIdType vtkID)
 {
-  return vtkCellLinks::Resize(sz);
-}
-
-vtkIdType SMDS_CellLinks::GetLinksSize()
-{
-  return this->Size;
+  if ( vtkID > this->MaxId )
+  {
+    this->MaxId = vtkID;
+    if ( vtkID >= this->Size ) 
+      vtkCellLinks::Resize( vtkID+SMDS_Mesh::chunkSize );
+  }
 }
 
 SMDS_CellLinks::SMDS_CellLinks() :
@@ -86,7 +87,8 @@ unsigned long SMDS_UnstructuredGrid::GetMTime()
   MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
   return mtime;
 }
-
+// OUV_PORTING_VTK6: seems to be useless
+/*
 void SMDS_UnstructuredGrid::Update()
 {
   MESSAGE("SMDS_UnstructuredGrid::Update");
@@ -98,7 +100,7 @@ void SMDS_UnstructuredGrid::UpdateInformation()
   MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
   return vtkUnstructuredGrid::UpdateInformation();
 }
-
+*/
 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
 {
   // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
@@ -152,7 +154,7 @@ void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
                                         std::vector<int>& idCellsOldToNew, int newCellSize)
 {
-  MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);CHRONO(1);
+  MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);//CHRONO(1);
   int alreadyCopied = 0;
 
   // --- if newNodeSize, create a new compacted vtkPoints
@@ -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++)
@@ -739,7 +763,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
  * @param vtkId the vtk id of the cell
  * @return number of neighbors
  */
-int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId)
+int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
 {
   int vtkType = this->GetCellType(vtkId);
   int cellDim = SMDS_Downward::getCellDimension(vtkType);
@@ -774,9 +798,27 @@ int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsi
           downIds[nb] = downId;
           downTypes[nb] = cellType;
           nb++;
+          if (nb >= NBMAXNEIGHBORS)
+            {
+              INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
+              return nb;
+            }
+        }
+      if (getSkin)
+        {
+          if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
+            {
+              neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
+              downIds[nb] = downId;
+              downTypes[nb] = cellType;
+              nb++;
+              if (nb >= NBMAXNEIGHBORS)
+                {
+                  INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
+                  return nb;
+                }
+            }
         }
-      if (nb >= NBMAXNEIGHBORS)
-        assert(0);
     }
   return nb;
 }
@@ -1072,3 +1114,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;
+}
+