Salome HOME
b5a69d8b074a6b7e5037d0d650dca067179b3388
[tools/solverlab.git] / CDMATH / mesh / inc / Mesh.hxx
1 /*
2  * mesh.hxx
3  *
4  *  Created on: 22 janv. 2012
5  *      Authors: CDMAT
6  */
7
8 #ifndef MESH_HXX_
9 #define MESH_HXX_
10
11 #include "MEDCouplingUMesh.hxx"
12 #include "MEDCouplingIMesh.hxx"
13 #include "MEDCouplingFieldDouble.hxx"
14
15 /**
16  * Mesh class is defined by
17  * - case 1: file name of mesh med file (general unstructured)
18  * - case 2: 1D cartesian, xmin and xmax and number of cells
19  * - case 3: 2D cartesian, xmin, xmax, ymin and ymax and numbers of cells in x direction and y direction
20  * - case 4: 3D cartesian, xmin, xmax, ymin, ymax, zmin and zmax and numbers of cells in x direction, y direction and z direction
21  * - case 5: 2D regular triangular mesh
22  * - case 6: 3D regular hexahedral mesh
23  * - case 7: 1D unstructured
24  */
25
26 namespace MEDCoupling
27 {
28 class MEDFileUMesh;
29 class MEDCouplingMesh;
30 class MEDCouplingIMesh;
31 class MEDCouplingUMesh;
32 class DataArrayIdType;
33 }
34 namespace ParaMEDMEM
35 {
36 class DataArrayIdType;
37 }
38 #include <MCAuto.hxx>
39 #include "NormalizedGeometricTypes"
40
41 class Node;
42 class Cell;
43 class Face;
44
45 typedef enum
46   {
47     CELLS = 0,
48     NODES = 1,
49     FACES = 2,
50   } EntityType;
51
52 #include <vector>
53 #include <string>
54 #include <map>
55
56 class Mesh
57 {
58
59 public: //----------------------------------------------------------------
60         /**
61          * \brief default constructor
62          */
63         Mesh ( void ) ;
64
65         /**
66          * \brief constructor with data to load a general unstructured mesh
67          * @param filename name of mesh file
68          * @param meshLevel : relative mesh dimension : 0->cells, 1->Faces etc
69          */
70         Mesh ( const std::string filename, int meshLevel=0 ) ;
71
72         /**
73          * \brief constructor with data for a regular 1D grid 
74          * @param xmin : minimum x
75          * @param xmax : maximum x
76          * @param nx : Number of cells in x direction
77          */
78         Mesh( double xmin, double xmax, int nx, std::string meshName="MESH1D_Regular_Grid" ) ;
79
80         /**
81          * \brief constructor with data for an unstructured 1D mesh
82          * @param points : abscissas of the mesh nodes
83          */
84         Mesh( std::vector<double> points, std::string meshName="MESH1D_unstructured" ) ;
85
86         /**
87          * \brief constructor with data for a regular 2D grid 
88          * @param xmin : minimum x
89          * @param xmax : maximum x
90          * @param ymin : minimum y
91          * @param ymax : maximum y
92          * @param nx : Number of cells in x direction
93          * @param ny : Number of cells in y direction
94      * @param split_to_triangles_policy : each rectangle will be split into 2 triangles with orientation of the cut depending if value is 0 or 1
95          */
96         Mesh( double xmin, double xmax, int nx, double ymin, double ymax, int ny, int split_to_triangles_policy=-1, std::string meshName="MESH2D_Regular_Rectangle_Grid") ;
97
98         /**
99          * \brief constructor with data for a regular 3D grid 
100          * @param xmin : minimum x
101          * @param xmax : maximum x
102          * @param ymin : minimum y
103          * @param ymax : maximum y
104          * @param zmin : minimum z
105          * @param zmax : maximum z
106          * @param nx : Number of cells in x direction
107          * @param ny : Number of cells in y direction
108          * @param nz : Number of cells in z direction
109      * @param split_to_tetrahedra_policy : each cuboid will be split into 5 tetrahedra if value is INTERP_KERNEL::PLANAR_FACE_5 or 6 tetrahedra if the value is INTERP_KERNEL::PLANAR_FACE_6
110          */
111         Mesh( double xmin, double xmax, int nx, double ymin, double ymax, int ny, double zmin, double zmax, int nz, int split_to_tetrahedra_policy=-1, std::string meshName="MESH3D_Regular_Cuboid_Grid") ;
112
113         Mesh( const MEDCoupling::MEDCouplingIMesh* mesh ) ;
114         Mesh( const MEDCoupling::MEDCouplingUMesh* mesh ) ;
115
116         /**
117          * \brief constructor with data
118          * @param filename : file name of mesh med file
119          * @param meshLevel : relative mesh dimension : 0->cells, 1->Faces etc
120          */
121         void readMeshMed( const std::string filename, int meshLevel=0 ) ;
122
123         /**
124          * \brief constructor by copy
125          * @param mesh : The Mesh object to be copied
126          */
127         Mesh ( const Mesh & mesh ) ;
128
129         /**
130          * \brief destructor
131          */
132         ~Mesh( void ) ;
133
134         /**
135          * \brief return mesh name
136          * @return _name
137          */
138         std::string getName( void ) const ;
139
140         /**
141          * \brief return Space dimension
142          * @return _spaceDim
143          */
144         int getSpaceDimension( void ) const ;
145
146         /**
147          * \brief Mesh dimension
148          * @return _meshDim
149          */
150         int getMeshDimension( void ) const ;
151
152         /**
153          * \brief return The nodes in this mesh
154          * @return _nodes
155          */
156         Node* getNodes ( void ) const ;
157
158         /**
159          * \brief return The cells in this mesh
160          * @return _vertices
161          */
162         Cell* getCells ( void ) const ;
163
164         /**
165          * \brief return the faces in this mesh
166          * @return _vertices
167          */
168         Face* getFaces ( void ) const ;
169
170         /**
171          * \brief return the number of nodes in this mesh
172          * @return _numberOfNodes
173          */
174         int getNumberOfNodes ( void )  const ;
175
176         /**
177          * \brief return the number of faces in this mesh
178          * @return _numberOfFaces
179          */
180         int getNumberOfFaces ( void )  const ;
181
182         /**
183          * \brief return the number of cells in this mesh
184          * @return _numberOfCells
185          */
186         int getNumberOfCells ( void )  const ;
187
188         /**
189          * \brief return the number of edges in this mesh
190          * @return _numberOfEdges
191          */
192         int getNumberOfEdges ( void )  const ;
193
194         /**
195          * \brief return the cell i in this mesh
196          * @return _cells[i]
197          */
198         Cell& getCell ( int i )  const ;
199
200         /**
201          * return The face i in this mesh
202          * @return _faces[i]
203          */
204         Face& getFace ( int i )  const ;
205
206         /**
207          * \brief return The node i in this mesh
208          * @return _nodes[i]
209          */
210         Node& getNode ( int i )  const ;
211
212         /**
213          * \brief return number of cell in x direction (structured mesh)
214          * @return _nX
215          */
216         int getNx( void )  const ;
217
218         /**
219          * \brief return number of cell in y direction (structured mesh)
220          * @return _nY
221          */
222         int getNy( void )  const ;
223
224         /**
225          * \brief return number of cell in z direction (structured mesh)
226          * @return _nZ
227          */
228         int getNz( void )  const ;
229
230         double getXMin( void )  const ;
231
232         double getXMax( void )  const ;
233
234         double getYMin( void )  const ;
235
236         double getYMax( void )  const ;
237
238         double getZMin( void )  const ;
239
240         double getZMax( void )  const ;
241
242         std::vector<double> getDXYZ() const ;// for structured meshes
243
244         std::vector<mcIdType> getCellGridStructure() const;// for structured meshes
245
246         /**
247          * \brief surcharge operator =
248          * @param mesh : The Mesh object to be copied
249          */
250         const Mesh& operator= ( const Mesh& mesh ) ;
251
252         /**
253          * \brief return the mesh MEDCoupling
254          * @return _mesh
255          */
256         MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> getMEDCouplingMesh ( void )  const ;
257
258         /**
259          * \brief computes the skin surrounding the mesh
260          */
261         Mesh getBoundaryMesh ( void )  const ;
262
263         /**
264          * \brief return the list of face group names
265          * return _faceGroupNames
266          */
267         std::vector<std::string> getNameOfFaceGroups( void )  const ;
268
269         /**
270          * \brief return the list of node group names
271          * return _nodeGroupNames
272          */
273         std::vector<std::string> getNameOfNodeGroups( void )  const ;
274
275         /**
276          * \brief return the list of face groups
277          * @return _faceGroups
278          */
279         std::vector<MEDCoupling::MEDCouplingUMesh *> getFaceGroups( void )  const ;
280
281         /**
282          * \brief return the list of node groups
283          * @return _nodeGroups
284          */
285         std::vector<MEDCoupling::DataArrayIdType *> getNodeGroups( void )  const ;
286
287     /*
288      * Functions to extract boundary nodes and faces Ids
289      */
290      /**
291       *  \brief return the list of boundary faces Ids
292       *  @return _boundaryFaceIds
293       */
294     std::vector< int > getBoundaryFaceIds() const;
295     /**
296      * \brief list of boundary nodes Ids
297      * @return _boundaryNodeIds
298      */
299     std::vector< int > getBoundaryNodeIds() const;
300     /*
301      * Functions to extract group nodes and faces ids
302      */
303      /** 
304       * @return list of face group Ids
305       */
306     std::vector< int > getGroupFaceIds(std::string groupName) const;
307     /**
308      * @return list of node group Ids
309      * */
310     std::vector< int > getGroupNodeIds(std::string groupName) const;
311  
312         /**
313          * \brief write mesh in the VTK format
314          */
315         void writeVTK ( const std::string fileName ) const ;
316
317         /**
318          * \brief write mesh in the MED format
319          */
320         void writeMED ( const std::string fileName ) const ;
321
322         void setGroupAtPlan(double value, int direction, double eps, std::string groupName) ;
323
324         void setGroupAtFaceByCoords(double x, double y, double z, double eps, std::string groupName) ;
325
326         void setFaceGroupByIds(std::vector< int > faceIds, std::string groupName) ;
327
328         void setNodeGroupByIds(std::vector< int > nodeIds, std::string groupName) ;
329
330         /*
331      * Functions to manage periodic boundary condition in square/cubic geometries 
332      */
333     void setPeriodicFaces(bool check_groups= false, bool use_central_inversion=false) ;
334     int getIndexFacePeriodic(int indexFace, bool check_groups= false, bool use_central_inversion=false);
335     void setBoundaryNodesFromFaces();
336     std::map<int,int> getIndexFacePeriodic( void ) const;
337     bool isIndexFacePeriodicSet() const ;
338     
339         bool isBorderNode(int nodeid) const ;
340         bool isBorderFace(int faceid) const ;
341         
342         bool isTriangular() const ;
343         bool isTetrahedral() const ;
344         bool isQuadrangular() const ;
345         bool isHexahedral() const ;
346     bool isStructured() const ;
347     std::vector< std::string > getElementTypesNames() const ;
348         /**
349          * \brief Compute the minimum value over all cells of the ratio cell perimeter/cell vaolume
350          */
351     double minRatioVolSurf() const;
352     
353         /**
354          * \brief Compute the maximum number of neighbours around an element (cells around a cell or nodes around a node)
355          */
356     int getMaxNbNeighbours(EntityType type) const;
357     
358 private: //----------------------------------------------------------------
359
360         MEDCoupling::MEDCouplingUMesh*  setMesh( void ) ;
361
362         void setGroups( const MEDCoupling::MEDFileUMesh* medmesh, MEDCoupling::MEDCouplingUMesh*  mu) ;
363
364     std::string _name;
365     
366         /**
367          * \brief Space dimension
368          */
369         int _spaceDim ;
370
371         /**
372          * \brief Mesh dimension
373          */
374         int _meshDim ;
375     
376     /*
377      * Structured mesh parameters
378      */
379
380     bool _isStructured;
381     
382         double _xMin;
383
384         double _xMax;
385
386         double _yMin;
387
388         double _yMax;
389
390         double _zMin;
391
392         double _zMax;
393
394         std::vector<mcIdType> _nxyz;
395
396         std::vector<double> _dxyz;
397         /*
398          * The nodes in this mesh.
399          */
400         Node *_nodes;
401
402         /*
403          * The number of nodes in this mesh.
404          */
405         int _numberOfNodes;
406
407         /*
408          * The faces in this mesh.
409          */
410         Face *_faces;
411
412         /*
413          * The numbers of faces in this mesh.
414          */
415         int _numberOfFaces;
416
417         /*
418          * The cells in this mesh.
419          */
420         Cell *_cells;
421
422         /*
423          * The number of cells in this mesh.
424          */
425         int _numberOfCells;
426
427         /*
428          * The number of edges in this mesh.
429          */
430         int _numberOfEdges;//Useful to deduce the number of non zero coefficients in the finite element matrix 
431
432         /*
433          * The names of face groups.
434          */
435         std::vector<std::string> _faceGroupNames;
436
437         /*
438          * The names of node groups.
439          */
440         std::vector<std::string> _nodeGroupNames;
441
442         /*
443          * The list of face groups.
444          */
445         std::vector<MEDCoupling::MEDCouplingUMesh *> _faceGroups;
446         /*
447          * The list of node groups.
448          */
449         std::vector<MEDCoupling::DataArrayIdType *> _nodeGroups;
450         /*
451          * The mesh MEDCoupling
452          */
453         MEDCoupling::MCAuto<MEDCoupling::MEDCouplingMesh> _mesh;
454
455         std::vector< INTERP_KERNEL::NormalizedCellType > _eltsTypes;//List of cell types contained in the mesh
456         std::vector< std::string > _eltsTypesNames;//List of cell types contained in the mesh
457     std::vector< INTERP_KERNEL::NormalizedCellType > getElementTypes() const;    
458     
459     /*
460      * Tools to manage periodic boundary conditions in square/cube geometries
461      */
462      bool _indexFacePeriodicSet;
463      std::map<int,int> _indexFacePeriodicMap;
464     
465     /* List of boundary faces*/
466     std::vector< int > _boundaryFaceIds;
467     /* List of boundary nodes*/
468     std::vector< int > _boundaryNodeIds;
469     /* Boundary mesh */
470     const MEDCoupling::MEDCouplingUMesh * _boundaryMesh;
471 };
472
473 #endif /* MESH_HXX_ */