Salome HOME
Merge branch 'V8_4_BR'
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.hxx
1 // Copyright (C) 2010-2016  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // File:    SMDS_UnstructuredGrid.hxx
21 // Author:  prascle
22 // Created: September 16, 2009, 10:28 PM
23
24 #ifndef _SMDS_UNSTRUCTUREDGRID_HXX
25 #define _SMDS_UNSTRUCTUREDGRID_HXX
26
27 #include "SMESH_SMDS.hxx"
28
29 #include <vtkUnstructuredGrid.h>
30 #include <vtkCellLinks.h>
31
32 #include <vector>
33 #include <set>
34 #include <map>
35
36 //#define VTK_HAVE_POLYHEDRON
37 //#ifdef VTK_HAVE_POLYHEDRON
38 #define VTK_MAXTYPE VTK_POLYHEDRON
39 //#else
40 //  #define VTK_MAXTYPE VTK_QUADRATIC_PYRAMID
41 //#endif
42
43 #define NBMAXNEIGHBORS 100
44
45 // allow very huge polyhedrons in tests
46 #define NBMAXNODESINCELL 5000
47
48 // Keep compatibility with paraview 5.0.1 on Linux
49 #ifndef WIN32
50   #ifndef VTK_HAS_MTIME_TYPE
51   #define VTK_HAS_MTIME_TYPE
52   typedef unsigned long int vtkMTimeType;
53   #endif
54 #endif
55
56 class SMDS_Downward;
57 class SMDS_Mesh;
58 class SMDS_MeshCell;
59 class SMDS_MeshVolume;
60
61 class SMDS_EXPORT SMDS_CellLinks: public vtkCellLinks
62 {
63 public:
64   void ResizeForPoint(vtkIdType vtkID);
65   void BuildLinks(vtkDataSet *data, vtkCellArray *Connectivity, vtkUnsignedCharArray* types);
66   static SMDS_CellLinks* New();
67 protected:
68   SMDS_CellLinks();
69   ~SMDS_CellLinks();
70 };
71
72 class SMDS_EXPORT SMDS_UnstructuredGrid: public vtkUnstructuredGrid
73 {
74 public:
75   void setSMDS_mesh(SMDS_Mesh *mesh);
76   void compactGrid(std::vector<int>& idNodesOldToNew,
77                    int               newNodeSize,
78                    std::vector<int>& idCellsOldToNew,
79                    int               newCellSize);
80   virtual vtkMTimeType GetMTime();
81   virtual vtkPoints *GetPoints();
82
83   int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
84
85   int CellIdToDownId(int vtkCellId);
86   void setCellIdToDownId(int vtkCellId, int downId);
87   void CleanDownwardConnectivity();
88   void BuildDownwardConnectivity(bool withEdges);
89   int GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin=false);
90   int GetParentVolumes(int* volVtkIds, int vtkId);
91   int GetParentVolumes(int* volVtkIds, int downId, unsigned char downType);
92   void GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType);
93   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
94   int getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes);
95   SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2,
96                                        std::set<int>&                      originalNodes,
97                                        std::map<int, std::map<int, int> >& nodeDomains,
98                                        std::map<int, std::map<long,int> >& nodeQuadDomains);
99   void BuildLinks();
100   void DeleteLinks();
101   SMDS_CellLinks* GetLinks();
102   bool HasLinks() const { return this->Links; }
103
104   SMDS_Downward* getDownArray(unsigned char vtkType)
105   {
106     return _downArray[vtkType];
107   }
108   void AllocateDiameters( vtkIdType maxVtkID );
109   void SetBallDiameter( vtkIdType vtkID, double diameter );
110   double GetBallDiameter( vtkIdType vtkID ) const;
111
112   static SMDS_UnstructuredGrid* New();
113   SMDS_Mesh *_mesh;
114
115 protected:
116   SMDS_UnstructuredGrid();
117   ~SMDS_UnstructuredGrid();
118   void copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied, int start, int end);
119   void copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, std::vector<int>& idNodesOldToNew,
120                 vtkCellArray* newConnectivity, vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
121                 int start, int end);
122
123   std::vector<int> _cellIdToDownId; //!< convert vtk Id to downward[vtkType] id, initialized with -1
124   std::vector<unsigned char> _downTypes;
125   std::vector<SMDS_Downward*> _downArray;
126 };
127
128 #endif  /* _SMDS_UNSTRUCTUREDGRID_HXX */
129