Salome HOME
Keep compatibility with ParaView-5.0.1
[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
49 #ifndef VTK_HAS_MTIME_TYPE
50 #define VTK_HAS_MTIME_TYPE
51 typedef unsigned long int vtkMTimeType;
52 #endif
53
54
55 class SMDS_Downward;
56 class SMDS_Mesh;
57 class SMDS_MeshCell;
58 class SMDS_MeshVolume;
59
60 class SMDS_EXPORT SMDS_CellLinks: public vtkCellLinks
61 {
62 public:
63   void ResizeForPoint(vtkIdType vtkID);
64   void BuildLinks(vtkDataSet *data, vtkCellArray *Connectivity, vtkUnsignedCharArray* types);
65   static SMDS_CellLinks* New();
66 protected:
67   SMDS_CellLinks();
68   ~SMDS_CellLinks();
69 };
70
71 class SMDS_EXPORT SMDS_UnstructuredGrid: public vtkUnstructuredGrid
72 {
73 public:
74   void setSMDS_mesh(SMDS_Mesh *mesh);
75   void compactGrid(std::vector<int>& idNodesOldToNew,
76                    int               newNodeSize,
77                    std::vector<int>& idCellsOldToNew,
78                    int               newCellSize);
79   virtual vtkMTimeType GetMTime();
80   virtual vtkPoints *GetPoints();
81
82   int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
83
84   int CellIdToDownId(int vtkCellId);
85   void setCellIdToDownId(int vtkCellId, int downId);
86   void CleanDownwardConnectivity();
87   void BuildDownwardConnectivity(bool withEdges);
88   int GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin=false);
89   int GetParentVolumes(int* volVtkIds, int vtkId);
90   int GetParentVolumes(int* volVtkIds, int downId, unsigned char downType);
91   void GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType);
92   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
93   int getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes);
94   SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2,
95                                        std::set<int>&                      originalNodes,
96                                        std::map<int, std::map<int, int> >& nodeDomains,
97                                        std::map<int, std::map<long,int> >& nodeQuadDomains);
98   void BuildLinks();
99   void DeleteLinks();
100   SMDS_CellLinks* GetLinks();
101   bool HasLinks() const { return this->Links; }
102
103   SMDS_Downward* getDownArray(unsigned char vtkType)
104   {
105     return _downArray[vtkType];
106   }
107   void AllocateDiameters( vtkIdType maxVtkID );
108   void SetBallDiameter( vtkIdType vtkID, double diameter );
109   double GetBallDiameter( vtkIdType vtkID ) const;
110
111   static SMDS_UnstructuredGrid* New();
112   SMDS_Mesh *_mesh;
113
114 protected:
115   SMDS_UnstructuredGrid();
116   ~SMDS_UnstructuredGrid();
117   void copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied, int start, int end);
118   void copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, std::vector<int>& idNodesOldToNew,
119                 vtkCellArray* newConnectivity, vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
120                 int start, int end);
121
122   std::vector<int> _cellIdToDownId; //!< convert vtk Id to downward[vtkType] id, initialized with -1
123   std::vector<unsigned char> _downTypes;
124   std::vector<SMDS_Downward*> _downArray;
125 };
126
127 #endif  /* _SMDS_UNSTRUCTUREDGRID_HXX */
128