Salome HOME
Merge from V6_3_BR 06/06/2011
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.hxx
1 // Copyright (C) 2010-2011  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.
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 <vtkUnstructuredGrid.h>
28 #include <vtkCellLinks.h>
29 #include "chrono.hxx"
30
31 #include <vector>
32 #include <set>
33 #include <map>
34
35 //#define VTK_HAVE_POLYHEDRON
36 //#ifdef VTK_HAVE_POLYHEDRON
37 #define VTK_MAXTYPE VTK_POLYHEDRON
38 //#else
39 //  #define VTK_MAXTYPE VTK_QUADRATIC_PYRAMID
40 //#endif
41
42 #define NBMAXNEIGHBORS 100
43
44 // allow very huge polyhedrons in tests
45 #define NBMAXNODESINCELL 5000
46
47 class SMDS_Downward;
48 class SMDS_Mesh;
49 class SMDS_MeshVolume;
50
51 class SMDS_CellLinks: public vtkCellLinks
52 {
53 public:
54   vtkCellLinks::Link* ResizeL(vtkIdType sz);
55   vtkIdType GetLinksSize();
56   static SMDS_CellLinks* New();
57 protected:
58   SMDS_CellLinks();
59   ~SMDS_CellLinks();
60 };
61
62 class SMDS_UnstructuredGrid: public vtkUnstructuredGrid
63 {
64 public:
65   void setSMDS_mesh(SMDS_Mesh *mesh);
66   void compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize, std::vector<int>& idCellsOldToNew,
67                    int newCellSize);
68
69   virtual unsigned long GetMTime();
70   virtual void Update();
71   virtual void UpdateInformation();
72   virtual vtkPoints *GetPoints();
73
74   //#ifdef VTK_HAVE_POLYHEDRON
75   int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
76   //#endif
77
78   int CellIdToDownId(int vtkCellId);
79   void setCellIdToDownId(int vtkCellId, int downId);
80   void CleanDownwardConnectivity();
81   void BuildDownwardConnectivity(bool withEdges);
82   int GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId);
83   int GetParentVolumes(int* volVtkIds, int vtkId);
84   int GetParentVolumes(int* volVtkIds, int downId, unsigned char downType);
85   void GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType);
86   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
87   int getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes);
88   void BuildLinks();
89   SMDS_MeshVolume* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2, std::set<int>& originalNodes,
90                                          std::map<int, std::map<int, int> >& nodeDomains,
91                                          std::map<int, std::map<long,int> >& nodeQuadDomains);
92   vtkCellLinks* GetLinks()
93   {
94     return Links;
95   }
96   SMDS_Downward* getDownArray(unsigned char vtkType)
97   {
98     return _downArray[vtkType];
99   }
100   static SMDS_UnstructuredGrid* New();
101   SMDS_Mesh *_mesh;
102 protected:
103   SMDS_UnstructuredGrid();
104   ~SMDS_UnstructuredGrid();
105   void copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied, int start, int end);
106   void copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, std::vector<int>& idNodesOldToNew,
107                 vtkCellArray* newConnectivity, vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
108                 int start, int end);
109
110   std::vector<int> _cellIdToDownId; //!< convert vtk Id to downward[vtkType] id, initialized with -1
111   std::vector<unsigned char> _downTypes;
112   std::vector<SMDS_Downward*> _downArray;
113 };
114
115 #endif  /* _SMDS_UNSTRUCTUREDGRID_HXX */
116