Salome HOME
Fix for the "0021861: EDF 2226 : Documentation of numeric functor option in split...
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
index 2902856bc49b077481cc8fc40b878159cacc8967..1eb5c445d60e5e13ed1af797359a7f258ceb53d5 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2010-2011  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2010-2012  CEA/DEN, EDF R&D, OPEN CASCADE
 //
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
@@ -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++)
@@ -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;
 }
@@ -929,15 +971,17 @@ void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> loc
  * @param orderedNodes list of nodes to reorder (in out)
  * @return size of the list
  */
-int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes)
+int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
 {
   int vtkType = this->GetCellType(vtkVolId);
-  int cellDim = SMDS_Downward::getCellDimension(vtkType);
-  if (cellDim != 3)
-    return 0;
-  SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
-  int downVolId = this->_cellIdToDownId[vtkVolId];
-  downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
+  dim = SMDS_Downward::getCellDimension(vtkType);
+  if (dim == 3)
+    {
+      SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
+      int downVolId = this->_cellIdToDownId[vtkVolId];
+      downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
+    }
+  // else nothing to do;
   return orderedNodes.size();
 }
 
@@ -970,7 +1014,7 @@ void SMDS_UnstructuredGrid::BuildLinks()
  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
  * @return ok if success.
  */
-SMDS_MeshVolume* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
+SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
                                                   int domain1,
                                                   int domain2,
                                                   std::set<int>& originalNodes,
@@ -984,65 +1028,136 @@ SMDS_MeshVolume* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
   for (; it != originalNodes.end(); ++it)
     orderedOriginals.push_back(*it);
 
-  int nbNodes = this->getOrderedNodesOfFace(vtkVolId, orderedOriginals);
+  int dim = 0;
+  int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
   vector<vtkIdType> orderedNodes;
 
+  bool isQuadratic = false;
   switch (orderedOriginals.size())
   {
     case 3:
-    case 4:
-      for (int i = 0; i < nbNodes; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
-      for (int i = 0; i < nbNodes; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+      if (dim == 2)
+        isQuadratic = true;
       break;
     case 6:
     case 8:
-      {
-        long dom1 = domain1;
-        long dom2 = domain2;
-        long dom1_2; // for nodeQuadDomains
-        if (domain1 < domain2)
-          dom1_2 = dom1 + INT_MAX * dom2;
-        else
-          dom1_2 = dom2 + INT_MAX * dom1;
-        //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
-        int ima = orderedOriginals.size();
-        int mid = orderedOriginals.size() / 2;
-        //cerr << "ima=" << ima << " mid=" << mid << endl;
-        for (int i = 0; i < mid; i++)
-          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
-        for (int i = 0; i < mid; i++)
-          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
-        for (int i = mid; i < ima; i++)
-          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
-        for (int i = mid; i < ima; i++)
-          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
-        for (int i = 0; i < mid; i++)
-          {
-            int oldId = orderedOriginals[i];
-            int newId;
-            if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
-              newId = nodeQuadDomains[oldId][dom1_2];
-            else
-              {
-                double *coords = this->GetPoint(oldId);
-                SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
-                newId = newNode->getVtkId();
-                std::map<long, int> emptyMap;
-                nodeQuadDomains[oldId] = emptyMap;
-                nodeQuadDomains[oldId][dom1_2] = newId;
-              }
-            orderedNodes.push_back(newId);
-          }
-      }
+      isQuadratic = true;
       break;
     default:
-      ASSERT(0);
+      isQuadratic = false;
+      break;
   }
 
-  SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
+  if (isQuadratic)
+    {
+      long dom1 = domain1;
+      long dom2 = domain2;
+      long dom1_2; // for nodeQuadDomains
+      if (domain1 < domain2)
+        dom1_2 = dom1 + INT_MAX * dom2;
+      else
+        dom1_2 = dom2 + INT_MAX * dom1;
+      //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
+      int ima = orderedOriginals.size();
+      int mid = orderedOriginals.size() / 2;
+      //cerr << "ima=" << ima << " mid=" << mid << endl;
+      for (int i = 0; i < mid; i++)
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
+      for (int i = 0; i < mid; i++)
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+      for (int i = mid; i < ima; i++)
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
+      for (int i = mid; i < ima; i++)
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+      for (int i = 0; i < mid; i++)
+        {
+          int oldId = orderedOriginals[i];
+          int newId;
+          if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
+            newId = nodeQuadDomains[oldId][dom1_2];
+          else
+            {
+              double *coords = this->GetPoint(oldId);
+              SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
+              newId = newNode->getVtkId();
+              std::map<long, int> emptyMap;
+              nodeQuadDomains[oldId] = emptyMap;
+              nodeQuadDomains[oldId][dom1_2] = newId;
+            }
+          orderedNodes.push_back(newId);
+        }
+    }
+  else
+    {
+      for (int i = 0; i < nbNodes; i++)
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
+      if (dim == 3)
+        for (int i = 0; i < nbNodes; i++)
+          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+      else
+        for (int i = nbNodes-1; i >= 0; i--)
+          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+
+    }
+
+  if (dim == 3)
+    {
+      SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
+      return vol;
+    }
+  else if (dim == 2)
+    {
+      SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
+      return face;
+    }
 
   // TODO update sub-shape list of elements and nodes
-  return vol;
+  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;
+}
+