]> SALOME platform Git repositories - modules/smesh.git/blob - src/SMDS/SMDS_UnstructuredGrid.hxx
Salome HOME
0021422: EDF 1963 SMESH: Viscous layer algorithm fails in some cases
[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 "SMESH_SMDS.hxx"
28
29 #include <vtkUnstructuredGrid.h>
30 #include <vtkCellLinks.h>
31 #include "chrono.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 class SMDS_Downward;
50 class SMDS_Mesh;
51 class SMDS_MeshVolume;
52
53 class SMDS_EXPORT SMDS_CellLinks: public vtkCellLinks
54 {
55 public:
56   vtkCellLinks::Link* ResizeL(vtkIdType sz);
57   vtkIdType GetLinksSize();
58   static SMDS_CellLinks* New();
59 protected:
60   SMDS_CellLinks();
61   ~SMDS_CellLinks();
62 };
63
64 class SMDS_EXPORT SMDS_UnstructuredGrid: public vtkUnstructuredGrid
65 {
66 public:
67   void setSMDS_mesh(SMDS_Mesh *mesh);
68   void compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize, std::vector<int>& idCellsOldToNew,
69                    int newCellSize);
70
71   virtual unsigned long GetMTime();
72   virtual void Update();
73   virtual void UpdateInformation();
74   virtual vtkPoints *GetPoints();
75
76   //#ifdef VTK_HAVE_POLYHEDRON
77   int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
78   //#endif
79
80   int CellIdToDownId(int vtkCellId);
81   void setCellIdToDownId(int vtkCellId, int downId);
82   void CleanDownwardConnectivity();
83   void BuildDownwardConnectivity(bool withEdges);
84   int GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId);
85   int GetParentVolumes(int* volVtkIds, int vtkId);
86   int GetParentVolumes(int* volVtkIds, int downId, unsigned char downType);
87   void GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType);
88   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
89   int getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes);
90   void BuildLinks();
91   SMDS_MeshVolume* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2, std::set<int>& originalNodes,
92                                          std::map<int, std::map<int, int> >& nodeDomains,
93                                          std::map<int, std::map<long,int> >& nodeQuadDomains);
94   vtkCellLinks* GetLinks()
95   {
96     return Links;
97   }
98   SMDS_Downward* getDownArray(unsigned char vtkType)
99   {
100     return _downArray[vtkType];
101   }
102   static SMDS_UnstructuredGrid* New();
103   SMDS_Mesh *_mesh;
104 protected:
105   SMDS_UnstructuredGrid();
106   ~SMDS_UnstructuredGrid();
107   void copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied, int start, int end);
108   void copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, std::vector<int>& idNodesOldToNew,
109                 vtkCellArray* newConnectivity, vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
110                 int start, int end);
111
112   std::vector<int> _cellIdToDownId; //!< convert vtk Id to downward[vtkType] id, initialized with -1
113   std::vector<unsigned char> _downTypes;
114   std::vector<SMDS_Downward*> _downArray;
115 };
116
117 #endif  /* _SMDS_UNSTRUCTUREDGRID_HXX */
118