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