Salome HOME
Update copyrights 2014.
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.hxx
1 // Copyright (C) 2010-2014  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 "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_MeshCell;
52 class SMDS_MeshVolume;
53
54 class SMDS_EXPORT SMDS_CellLinks: public vtkCellLinks
55 {
56 public:
57   void ResizeForPoint(vtkIdType vtkID);
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,
69                    int               newNodeSize,
70                    std::vector<int>& idCellsOldToNew,
71                    int               newCellSize);
72   virtual unsigned long GetMTime();
73   // OUV_PORTING_VTK6: seems to be useless
74   //virtual void Update();
75   //virtual void UpdateInformation();
76   virtual vtkPoints *GetPoints();
77
78   //#ifdef VTK_HAVE_POLYHEDRON
79   int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
80   //#endif
81
82   int CellIdToDownId(int vtkCellId);
83   void setCellIdToDownId(int vtkCellId, int downId);
84   void CleanDownwardConnectivity();
85   void BuildDownwardConnectivity(bool withEdges);
86   int GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin=false);
87   int GetParentVolumes(int* volVtkIds, int vtkId);
88   int GetParentVolumes(int* volVtkIds, int downId, unsigned char downType);
89   void GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType);
90   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
91   int getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes);
92   void BuildLinks();
93   SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2,
94                                        std::set<int>&                      originalNodes,
95                                        std::map<int, std::map<int, int> >& nodeDomains,
96                                        std::map<int, std::map<long,int> >& nodeQuadDomains);
97   vtkCellLinks* GetLinks()
98   {
99     return Links;
100   }
101   SMDS_Downward* getDownArray(unsigned char vtkType)
102   {
103     return _downArray[vtkType];
104   }
105   void AllocateDiameters( vtkIdType maxVtkID );
106   void SetBallDiameter( vtkIdType vtkID, double diameter );
107   double GetBallDiameter( vtkIdType vtkID ) const;
108
109   static SMDS_UnstructuredGrid* New();
110   SMDS_Mesh *_mesh;
111
112 protected:
113   SMDS_UnstructuredGrid();
114   ~SMDS_UnstructuredGrid();
115   void copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied, int start, int end);
116   void copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, std::vector<int>& idNodesOldToNew,
117                 vtkCellArray* newConnectivity, vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
118                 int start, int end);
119
120   std::vector<int> _cellIdToDownId; //!< convert vtk Id to downward[vtkType] id, initialized with -1
121   std::vector<unsigned char> _downTypes;
122   std::vector<SMDS_Downward*> _downArray;
123 };
124
125 #endif  /* _SMDS_UNSTRUCTUREDGRID_HXX */
126