Salome HOME
Merge from V6_3_BR 06/06/2011
[modules/smesh.git] / src / SMDS / SMDS_Mesh.hxx
1 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  SMESH SMDS : implementaion of Salome mesh data structure
24 //  File   : SMDS_Mesh.hxx
25 //  Module : SMESH
26 //
27 #ifndef _SMDS_Mesh_HeaderFile
28 #define _SMDS_Mesh_HeaderFile
29
30 #include "SMESH_SMDS.hxx"
31
32 #include "SMDS_MeshNode.hxx"
33 #include "SMDS_MeshCell.hxx"
34 #include "SMDS_Mesh0DElement.hxx"
35 #include "SMDS_MeshEdge.hxx"
36 #include "SMDS_MeshFace.hxx"
37 #include "SMDS_MeshVolume.hxx"
38 #include "SMDS_MeshNodeIDFactory.hxx"
39 #include "SMDS_MeshElementIDFactory.hxx"
40 #include "SMDS_MeshInfo.hxx"
41 #include "SMDS_ElemIterator.hxx"
42 #include "SMDS_VolumeOfNodes.hxx"
43 #include "SMDS_VtkEdge.hxx"
44 #include "SMDS_VtkFace.hxx"
45 #include "SMDS_VtkVolume.hxx"
46 #include "ObjectPool.hxx"
47 #include "SMDS_UnstructuredGrid.hxx"
48
49 #include <boost/shared_ptr.hpp>
50 #include <set>
51 #include <list>
52 #include <vector>
53 #include <vtkSystemIncludes.h>
54
55 #include "Utils_SALOME_Exception.hxx"
56 #define MYASSERT(val) if (!(val)) throw SALOME_Exception(LOCALIZED("assertion not verified"));
57
58 class SMDS_EXPORT SMDS_Mesh:public SMDS_MeshObject{
59 public:
60   friend class SMDS_MeshIDFactory;
61   friend class SMDS_MeshNodeIDFactory;
62   friend class SMDS_MeshElementIDFactory;
63   friend class SMDS_MeshVolumeVtkNodes;
64   friend class SMDS_MeshNode;
65
66   SMDS_Mesh();
67   
68   //! to retreive this SMDS_Mesh instance from its elements (index stored in SMDS_Elements)
69   static std::vector<SMDS_Mesh*> _meshList;
70
71   //! actual nodes coordinates, cells definition and reverse connectivity are stored in a vtkUnstructuredGrid
72   inline SMDS_UnstructuredGrid* getGrid() {return myGrid; }
73   inline int getMeshId() {return myMeshId; }
74
75   SMDS_NodeIteratorPtr nodesIterator          (bool idInceasingOrder=false) const;
76   SMDS_0DElementIteratorPtr elements0dIterator(bool idInceasingOrder=false) const;
77   SMDS_EdgeIteratorPtr edgesIterator          (bool idInceasingOrder=false) const;
78   SMDS_FaceIteratorPtr facesIterator          (bool idInceasingOrder=false) const;
79   SMDS_VolumeIteratorPtr volumesIterator      (bool idInceasingOrder=false) const;
80   SMDS_ElemIteratorPtr elementsIterator(SMDSAbs_ElementType type=SMDSAbs_All) const;
81
82   SMDSAbs_ElementType GetElementType( const int id, const bool iselem ) const;
83
84   SMDS_Mesh *AddSubMesh();
85
86   virtual SMDS_MeshNode* AddNodeWithID(double x, double y, double z, int ID);
87   virtual SMDS_MeshNode* AddNode(double x, double y, double z);
88
89   virtual SMDS_Mesh0DElement* Add0DElementWithID(int n, int ID);
90   virtual SMDS_Mesh0DElement* Add0DElementWithID(const SMDS_MeshNode * n, int ID);
91   virtual SMDS_Mesh0DElement* Add0DElement      (const SMDS_MeshNode * n);
92
93   virtual SMDS_MeshEdge* AddEdgeWithID(int n1, int n2, int ID);
94   virtual SMDS_MeshEdge* AddEdgeWithID(const SMDS_MeshNode * n1,
95                                        const SMDS_MeshNode * n2,
96                                        int ID);
97   virtual SMDS_MeshEdge* AddEdge(const SMDS_MeshNode * n1,
98                                  const SMDS_MeshNode * n2);
99
100   // 2d order edge with 3 nodes: n12 - node between n1 and n2
101   virtual SMDS_MeshEdge* AddEdgeWithID(int n1, int n2, int n12, int ID);
102   virtual SMDS_MeshEdge* AddEdgeWithID(const SMDS_MeshNode * n1,
103                                        const SMDS_MeshNode * n2,
104                                        const SMDS_MeshNode * n12,
105                                        int ID);
106   virtual SMDS_MeshEdge* AddEdge(const SMDS_MeshNode * n1,
107                                  const SMDS_MeshNode * n2,
108                                  const SMDS_MeshNode * n12);
109
110   virtual SMDS_MeshFace* AddFaceWithID(int n1, int n2, int n3, int ID);
111   virtual SMDS_MeshFace* AddFaceWithID(const SMDS_MeshNode * n1,
112                                        const SMDS_MeshNode * n2,
113                                        const SMDS_MeshNode * n3,
114                                        int ID);
115   virtual SMDS_MeshFace* AddFace(const SMDS_MeshNode * n1,
116                                  const SMDS_MeshNode * n2,
117                                  const SMDS_MeshNode * n3);
118
119   virtual SMDS_MeshFace* AddFaceWithID(int n1, int n2, int n3, int n4, int ID);
120   virtual SMDS_MeshFace* AddFaceWithID(const SMDS_MeshNode * n1,
121                                        const SMDS_MeshNode * n2,
122                                        const SMDS_MeshNode * n3,
123                                        const SMDS_MeshNode * n4,
124                                        int ID);
125   virtual SMDS_MeshFace* AddFace(const SMDS_MeshNode * n1,
126                                  const SMDS_MeshNode * n2,
127                                  const SMDS_MeshNode * n3,
128                                  const SMDS_MeshNode * n4);
129
130   virtual SMDS_MeshFace* AddFaceWithID(const SMDS_MeshEdge * e1,
131                                        const SMDS_MeshEdge * e2,
132                                        const SMDS_MeshEdge * e3, int ID);
133   virtual SMDS_MeshFace* AddFace(const SMDS_MeshEdge * e1,
134                                  const SMDS_MeshEdge * e2,
135                                  const SMDS_MeshEdge * e3);
136
137   virtual SMDS_MeshFace* AddFaceWithID(const SMDS_MeshEdge * e1,
138                                        const SMDS_MeshEdge * e2,
139                                        const SMDS_MeshEdge * e3,
140                                        const SMDS_MeshEdge * e4, int ID);
141   virtual SMDS_MeshFace* AddFace(const SMDS_MeshEdge * e1,
142                                  const SMDS_MeshEdge * e2,
143                                  const SMDS_MeshEdge * e3,
144                                  const SMDS_MeshEdge * e4);
145
146   // 2d order triangle of 6 nodes
147   virtual SMDS_MeshFace* AddFaceWithID(int n1, int n2, int n3,
148                                        int n12,int n23,int n31, int ID);
149   virtual SMDS_MeshFace* AddFaceWithID(const SMDS_MeshNode * n1,
150                                        const SMDS_MeshNode * n2,
151                                        const SMDS_MeshNode * n3,
152                                        const SMDS_MeshNode * n12,
153                                        const SMDS_MeshNode * n23,
154                                        const SMDS_MeshNode * n31,
155                                        int ID);
156   virtual SMDS_MeshFace* AddFace(const SMDS_MeshNode * n1,
157                                  const SMDS_MeshNode * n2,
158                                  const SMDS_MeshNode * n3,
159                                  const SMDS_MeshNode * n12,
160                                  const SMDS_MeshNode * n23,
161                                  const SMDS_MeshNode * n31);
162
163   // 2d order quadrangle
164   virtual SMDS_MeshFace* AddFaceWithID(int n1, int n2, int n3, int n4,
165                                        int n12,int n23,int n34,int n41, int ID);
166   virtual SMDS_MeshFace* AddFaceWithID(const SMDS_MeshNode * n1,
167                                        const SMDS_MeshNode * n2,
168                                        const SMDS_MeshNode * n3,
169                                        const SMDS_MeshNode * n4,
170                                        const SMDS_MeshNode * n12,
171                                        const SMDS_MeshNode * n23,
172                                        const SMDS_MeshNode * n34,
173                                        const SMDS_MeshNode * n41,
174                                        int ID);
175   virtual SMDS_MeshFace* AddFace(const SMDS_MeshNode * n1,
176                                  const SMDS_MeshNode * n2,
177                                  const SMDS_MeshNode * n3,
178                                  const SMDS_MeshNode * n4,
179                                  const SMDS_MeshNode * n12,
180                                  const SMDS_MeshNode * n23,
181                                  const SMDS_MeshNode * n34,
182                                  const SMDS_MeshNode * n41);
183
184   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4, int ID);
185   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
186                                            const SMDS_MeshNode * n2,
187                                            const SMDS_MeshNode * n3,
188                                            const SMDS_MeshNode * n4,
189                                            int ID);
190   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
191                                      const SMDS_MeshNode * n2,
192                                      const SMDS_MeshNode * n3,
193                                      const SMDS_MeshNode * n4);
194
195   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4,
196                                            int n5, int ID);
197   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
198                                            const SMDS_MeshNode * n2,
199                                            const SMDS_MeshNode * n3,
200                                            const SMDS_MeshNode * n4,
201                                            const SMDS_MeshNode * n5,
202                                            int ID);
203   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
204                                      const SMDS_MeshNode * n2,
205                                      const SMDS_MeshNode * n3,
206                                      const SMDS_MeshNode * n4,
207                                      const SMDS_MeshNode * n5);
208
209   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4,
210                                            int n5, int n6, int ID);
211   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
212                                            const SMDS_MeshNode * n2,
213                                            const SMDS_MeshNode * n3,
214                                            const SMDS_MeshNode * n4,
215                                            const SMDS_MeshNode * n5,
216                                            const SMDS_MeshNode * n6,
217                                            int ID);
218   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
219                                      const SMDS_MeshNode * n2,
220                                      const SMDS_MeshNode * n3,
221                                      const SMDS_MeshNode * n4,
222                                      const SMDS_MeshNode * n5,
223                                      const SMDS_MeshNode * n6);
224
225   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4,
226                                            int n5, int n6, int n7, int n8, int ID);
227   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
228                                            const SMDS_MeshNode * n2,
229                                            const SMDS_MeshNode * n3,
230                                            const SMDS_MeshNode * n4,
231                                            const SMDS_MeshNode * n5,
232                                            const SMDS_MeshNode * n6,
233                                            const SMDS_MeshNode * n7,
234                                            const SMDS_MeshNode * n8,
235                                            int ID);
236   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
237                                      const SMDS_MeshNode * n2,
238                                      const SMDS_MeshNode * n3,
239                                      const SMDS_MeshNode * n4,
240                                      const SMDS_MeshNode * n5,
241                                      const SMDS_MeshNode * n6,
242                                      const SMDS_MeshNode * n7,
243                                      const SMDS_MeshNode * n8);
244
245   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshFace * f1,
246                                            const SMDS_MeshFace * f2,
247                                            const SMDS_MeshFace * f3,
248                                            const SMDS_MeshFace * f4, int ID);
249   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshFace * f1,
250                                      const SMDS_MeshFace * f2,
251                                      const SMDS_MeshFace * f3,
252                                      const SMDS_MeshFace * f4);
253
254   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshFace * f1,
255                                            const SMDS_MeshFace * f2,
256                                            const SMDS_MeshFace * f3,
257                                            const SMDS_MeshFace * f4,
258                                            const SMDS_MeshFace * f5, int ID);
259   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshFace * f1,
260                                      const SMDS_MeshFace * f2,
261                                      const SMDS_MeshFace * f3,
262                                      const SMDS_MeshFace * f4,
263                                      const SMDS_MeshFace * f5);
264
265   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshFace * f1,
266                                            const SMDS_MeshFace * f2,
267                                            const SMDS_MeshFace * f3,
268                                            const SMDS_MeshFace * f4,
269                                            const SMDS_MeshFace * f5,
270                                            const SMDS_MeshFace * f6, int ID);
271   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshFace * f1,
272                                      const SMDS_MeshFace * f2,
273                                      const SMDS_MeshFace * f3,
274                                      const SMDS_MeshFace * f4,
275                                      const SMDS_MeshFace * f5,
276                                      const SMDS_MeshFace * f6);
277
278   // 2d order tetrahedron of 10 nodes
279   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4,
280                                            int n12,int n23,int n31,
281                                            int n14,int n24,int n34, int ID);
282   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
283                                            const SMDS_MeshNode * n2,
284                                            const SMDS_MeshNode * n3,
285                                            const SMDS_MeshNode * n4,
286                                            const SMDS_MeshNode * n12,
287                                            const SMDS_MeshNode * n23,
288                                            const SMDS_MeshNode * n31,
289                                            const SMDS_MeshNode * n14,
290                                            const SMDS_MeshNode * n24,
291                                            const SMDS_MeshNode * n34,
292                                            int ID);
293   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
294                                      const SMDS_MeshNode * n2,
295                                      const SMDS_MeshNode * n3,
296                                      const SMDS_MeshNode * n4,
297                                      const SMDS_MeshNode * n12,
298                                      const SMDS_MeshNode * n23,
299                                      const SMDS_MeshNode * n31,
300                                      const SMDS_MeshNode * n14,
301                                      const SMDS_MeshNode * n24,
302                                      const SMDS_MeshNode * n34);
303
304   // 2d order pyramid of 13 nodes
305   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4, int n5,
306                                            int n12,int n23,int n34,int n41,
307                                            int n15,int n25,int n35,int n45,
308                                            int ID);
309   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
310                                            const SMDS_MeshNode * n2,
311                                            const SMDS_MeshNode * n3,
312                                            const SMDS_MeshNode * n4,
313                                            const SMDS_MeshNode * n5,
314                                            const SMDS_MeshNode * n12,
315                                            const SMDS_MeshNode * n23,
316                                            const SMDS_MeshNode * n34,
317                                            const SMDS_MeshNode * n41,
318                                            const SMDS_MeshNode * n15,
319                                            const SMDS_MeshNode * n25,
320                                            const SMDS_MeshNode * n35,
321                                            const SMDS_MeshNode * n45,
322                                            int ID);
323   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
324                                      const SMDS_MeshNode * n2,
325                                      const SMDS_MeshNode * n3,
326                                      const SMDS_MeshNode * n4,
327                                      const SMDS_MeshNode * n5,
328                                      const SMDS_MeshNode * n12,
329                                      const SMDS_MeshNode * n23,
330                                      const SMDS_MeshNode * n34,
331                                      const SMDS_MeshNode * n41,
332                                      const SMDS_MeshNode * n15,
333                                      const SMDS_MeshNode * n25,
334                                      const SMDS_MeshNode * n35,
335                                      const SMDS_MeshNode * n45);
336
337   // 2d order Pentahedron with 15 nodes
338   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3,
339                                            int n4, int n5, int n6,
340                                            int n12,int n23,int n31,
341                                            int n45,int n56,int n64,
342                                            int n14,int n25,int n36,
343                                            int ID);
344   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
345                                            const SMDS_MeshNode * n2,
346                                            const SMDS_MeshNode * n3,
347                                            const SMDS_MeshNode * n4,
348                                            const SMDS_MeshNode * n5,
349                                            const SMDS_MeshNode * n6,
350                                            const SMDS_MeshNode * n12,
351                                            const SMDS_MeshNode * n23,
352                                            const SMDS_MeshNode * n31,
353                                            const SMDS_MeshNode * n45,
354                                            const SMDS_MeshNode * n56,
355                                            const SMDS_MeshNode * n64,
356                                            const SMDS_MeshNode * n14,
357                                            const SMDS_MeshNode * n25,
358                                            const SMDS_MeshNode * n36,
359                                            int ID);
360   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
361                                      const SMDS_MeshNode * n2,
362                                      const SMDS_MeshNode * n3,
363                                      const SMDS_MeshNode * n4,
364                                      const SMDS_MeshNode * n5,
365                                      const SMDS_MeshNode * n6,
366                                      const SMDS_MeshNode * n12,
367                                      const SMDS_MeshNode * n23,
368                                      const SMDS_MeshNode * n31,
369                                      const SMDS_MeshNode * n45,
370                                      const SMDS_MeshNode * n56,
371                                      const SMDS_MeshNode * n64,
372                                      const SMDS_MeshNode * n14,
373                                      const SMDS_MeshNode * n25,
374                                      const SMDS_MeshNode * n36);
375
376   // 2d oreder Hexahedrons with 20 nodes
377   virtual SMDS_MeshVolume* AddVolumeWithID(int n1, int n2, int n3, int n4,
378                                            int n5, int n6, int n7, int n8,
379                                            int n12,int n23,int n34,int n41,
380                                            int n56,int n67,int n78,int n85,
381                                            int n15,int n26,int n37,int n48,
382                                            int ID);
383   virtual SMDS_MeshVolume* AddVolumeWithID(const SMDS_MeshNode * n1,
384                                            const SMDS_MeshNode * n2,
385                                            const SMDS_MeshNode * n3,
386                                            const SMDS_MeshNode * n4,
387                                            const SMDS_MeshNode * n5,
388                                            const SMDS_MeshNode * n6,
389                                            const SMDS_MeshNode * n7,
390                                            const SMDS_MeshNode * n8,
391                                            const SMDS_MeshNode * n12,
392                                            const SMDS_MeshNode * n23,
393                                            const SMDS_MeshNode * n34,
394                                            const SMDS_MeshNode * n41,
395                                            const SMDS_MeshNode * n56,
396                                            const SMDS_MeshNode * n67,
397                                            const SMDS_MeshNode * n78,
398                                            const SMDS_MeshNode * n85,
399                                            const SMDS_MeshNode * n15,
400                                            const SMDS_MeshNode * n26,
401                                            const SMDS_MeshNode * n37,
402                                            const SMDS_MeshNode * n48,
403                                            int ID);
404   virtual SMDS_MeshVolume* AddVolume(const SMDS_MeshNode * n1,
405                                      const SMDS_MeshNode * n2,
406                                      const SMDS_MeshNode * n3,
407                                      const SMDS_MeshNode * n4,
408                                      const SMDS_MeshNode * n5,
409                                      const SMDS_MeshNode * n6,
410                                      const SMDS_MeshNode * n7,
411                                      const SMDS_MeshNode * n8,
412                                      const SMDS_MeshNode * n12,
413                                      const SMDS_MeshNode * n23,
414                                      const SMDS_MeshNode * n34,
415                                      const SMDS_MeshNode * n41,
416                                      const SMDS_MeshNode * n56,
417                                      const SMDS_MeshNode * n67,
418                                      const SMDS_MeshNode * n78,
419                                      const SMDS_MeshNode * n85,
420                                      const SMDS_MeshNode * n15,
421                                      const SMDS_MeshNode * n26,
422                                      const SMDS_MeshNode * n37,
423                                      const SMDS_MeshNode * n48);
424
425   virtual SMDS_MeshFace* AddPolygonalFaceWithID (std::vector<int> nodes_ids,
426                                                  const int        ID);
427
428   virtual SMDS_MeshFace* AddPolygonalFaceWithID (std::vector<const SMDS_MeshNode*> nodes,
429                                                  const int                         ID);
430
431   virtual SMDS_MeshFace* AddPolygonalFace (std::vector<const SMDS_MeshNode*> nodes);
432
433   virtual SMDS_MeshVolume* AddPolyhedralVolumeWithID
434                            (std::vector<int> nodes_ids,
435                             std::vector<int> quantities,
436                             const int        ID);
437
438   virtual SMDS_MeshVolume* AddPolyhedralVolumeWithID
439                            (std::vector<const SMDS_MeshNode*> nodes,
440                             std::vector<int>                  quantities,
441                             const int                         ID);
442
443   virtual SMDS_MeshVolume* AddPolyhedralVolume
444                            (std::vector<const SMDS_MeshNode*> nodes,
445                             std::vector<int>                  quantities);
446
447   virtual SMDS_MeshVolume* AddVolumeFromVtkIds(const std::vector<vtkIdType>& vtkNodeIds);
448
449   virtual SMDS_MeshVolume* AddVolumeFromVtkIdsWithID(const std::vector<vtkIdType>& vtkNodeIds,
450                                                      const int ID);
451
452   virtual void RemoveElement(const SMDS_MeshElement *        elem,
453                              std::list<const SMDS_MeshElement *>& removedElems,
454                              std::list<const SMDS_MeshElement *>& removedNodes,
455                              const bool                      removenodes = false);
456   virtual void RemoveElement(const SMDS_MeshElement * elem, bool removenodes = false);
457   virtual void RemoveNode(const SMDS_MeshNode * node);
458   virtual void Remove0DElement(const SMDS_Mesh0DElement * elem0d);
459   virtual void RemoveEdge(const SMDS_MeshEdge * edge);
460   virtual void RemoveFace(const SMDS_MeshFace * face);
461   virtual void RemoveVolume(const SMDS_MeshVolume * volume);
462
463   /*! Remove only the given element and only if it is free.
464    *  Method does not work for meshes with descendants.
465    *  Implemented for fast cleaning of meshes.
466    */
467   virtual void RemoveFreeElement(const SMDS_MeshElement * elem);
468
469   virtual void Clear();
470
471   virtual bool RemoveFromParent();
472   virtual bool RemoveSubMesh(const SMDS_Mesh * aMesh);
473
474   bool ChangeElementNodes(const SMDS_MeshElement * elem,
475                           const SMDS_MeshNode    * nodes[],
476                           const int                nbnodes);
477   bool ChangePolyhedronNodes(const SMDS_MeshElement *                 elem,
478                              const std::vector<const SMDS_MeshNode*>& nodes,
479                              const std::vector<int> &                 quantities);
480
481   virtual void Renumber (const bool isNodes, const int startID = 1, const int deltaID = 1);
482   // Renumber all nodes or elements.
483   virtual void compactMesh();
484
485   const SMDS_MeshNode *FindNode(int idnode) const;
486   const SMDS_MeshNode *FindNodeVtk(int idnode) const;
487   const SMDS_Mesh0DElement* Find0DElement(int idnode) const;
488   const SMDS_MeshEdge *FindEdge(int idnode1, int idnode2) const;
489   const SMDS_MeshEdge *FindEdge(int idnode1, int idnode2, int idnode3) const;
490   const SMDS_MeshFace *FindFace(int idnode1, int idnode2, int idnode3) const;
491   const SMDS_MeshFace *FindFace(int idnode1, int idnode2, int idnode3, int idnode4) const;
492   const SMDS_MeshFace *FindFace(int idnode1, int idnode2, int idnode3,
493                                 int idnode4, int idnode5, int idnode6) const;
494   const SMDS_MeshFace *FindFace(int idnode1, int idnode2, int idnode3, int idnode4,
495                                 int idnode5, int idnode6, int idnode7, int idnode8) const;
496   const SMDS_MeshElement *FindElement(int IDelem) const;
497   static const SMDS_Mesh0DElement* Find0DElement(const SMDS_MeshNode * n);
498   static const SMDS_MeshEdge* FindEdge(const SMDS_MeshNode * n1,
499                                        const SMDS_MeshNode * n2);
500   static const SMDS_MeshEdge* FindEdge(const SMDS_MeshNode * n1,
501                                        const SMDS_MeshNode * n2,
502                                        const SMDS_MeshNode * n3);
503   static const SMDS_MeshFace* FindFace(const SMDS_MeshNode *n1,
504                                        const SMDS_MeshNode *n2,
505                                        const SMDS_MeshNode *n3);
506   static const SMDS_MeshFace* FindFace(const SMDS_MeshNode *n1,
507                                        const SMDS_MeshNode *n2,
508                                        const SMDS_MeshNode *n3,
509                                        const SMDS_MeshNode *n4);
510   static const SMDS_MeshFace* FindFace(const SMDS_MeshNode *n1,
511                                        const SMDS_MeshNode *n2,
512                                        const SMDS_MeshNode *n3,
513                                        const SMDS_MeshNode *n4,
514                                        const SMDS_MeshNode *n5,
515                                        const SMDS_MeshNode *n6);
516   static const SMDS_MeshFace* FindFace(const SMDS_MeshNode *n1,
517                                        const SMDS_MeshNode *n2,
518                                        const SMDS_MeshNode *n3,
519                                        const SMDS_MeshNode *n4,
520                                        const SMDS_MeshNode *n5,
521                                        const SMDS_MeshNode *n6,
522                                        const SMDS_MeshNode *n7,
523                                        const SMDS_MeshNode *n8);
524
525   const SMDS_MeshFace *FindFace(const std::vector<int>& nodes_ids) const;
526   static const SMDS_MeshFace* FindFace(const std::vector<const SMDS_MeshNode *>& nodes);
527   static const SMDS_MeshElement* FindElement(const std::vector<const SMDS_MeshNode *>& nodes,
528                                              const SMDSAbs_ElementType                 type=SMDSAbs_All,
529                                              const bool                                noMedium=true);
530
531   /*!
532    * \brief Raise an exception if free memory (ram+swap) too low
533     * \param doNotRaise - if true, suppres exception, just return free memory size
534     * \retval int - amount of available memory in MB or negative number in failure case
535    */
536   static int CheckMemory(const bool doNotRaise=false) throw (std::bad_alloc);
537
538   int MaxNodeID() const;
539   int MinNodeID() const;
540   int MaxElementID() const;
541   int MinElementID() const;
542
543   const SMDS_MeshInfo& GetMeshInfo() const { return myInfo; }
544
545   int NbNodes() const;
546   int Nb0DElements() const;
547   int NbEdges() const;
548   int NbFaces() const;
549   int NbVolumes() const;
550   int NbSubMesh() const;
551   void DumpNodes() const;
552   void Dump0DElements() const;
553   void DumpEdges() const;
554   void DumpFaces() const;
555   void DumpVolumes() const;
556   void DebugStats() const;
557   SMDS_Mesh *boundaryFaces();
558   SMDS_Mesh *boundaryEdges();
559   virtual ~SMDS_Mesh();
560   bool hasConstructionEdges();
561   bool hasConstructionFaces();
562   bool hasInverseElements();
563   void setConstructionEdges(bool);
564   void setConstructionFaces(bool);
565   void setInverseElements(bool);
566
567   /*!
568    * Checks if the element is present in mesh.
569    * Useful to determine dead pointers.
570    * Use this function for debug purpose only! Do not check in the code
571    * using it even in _DEBUG_ mode
572    */
573   bool Contains (const SMDS_MeshElement* elem) const;
574
575   typedef std::vector<SMDS_MeshNode *> SetOfNodes;
576   typedef std::vector<SMDS_MeshCell *> SetOfCells;
577
578   void updateNodeMinMax();
579   void updateBoundingBox();
580   double getMaxDim();
581   int fromVtkToSmds(int vtkid);
582
583   void incrementNodesCapacity(int nbNodes);
584   void incrementCellsCapacity(int nbCells);
585   void adjustStructure();
586   void dumpGrid(string ficdump="dumpGrid");
587   static int chunkSize;
588
589   //! low level modification: add, change or remove node or element
590   inline void setMyModified() { this->myModified = true; }
591
592   void Modified();
593   unsigned long GetMTime();
594   bool isCompacted();
595
596 protected:
597   SMDS_Mesh(SMDS_Mesh * parent);
598
599   SMDS_MeshFace * createTriangle(const SMDS_MeshNode * node1,
600                                  const SMDS_MeshNode * node2,
601                                  const SMDS_MeshNode * node3,
602                                  int ID);
603   SMDS_MeshFace * createQuadrangle(const SMDS_MeshNode * node1,
604                                    const SMDS_MeshNode * node2,
605                                    const SMDS_MeshNode * node3,
606                                    const SMDS_MeshNode * node4,
607                                    int ID);
608 //  SMDS_Mesh0DElement* Find0DElementOrCreate(const SMDS_MeshNode * n);
609   SMDS_MeshEdge* FindEdgeOrCreate(const SMDS_MeshNode * n1,
610                                   const SMDS_MeshNode * n2);
611   SMDS_MeshFace* FindFaceOrCreate(const SMDS_MeshNode *n1,
612                                   const SMDS_MeshNode *n2,
613                                   const SMDS_MeshNode *n3);
614   SMDS_MeshFace* FindFaceOrCreate(const SMDS_MeshNode *n1,
615                                   const SMDS_MeshNode *n2,
616                                   const SMDS_MeshNode *n3,
617                                   const SMDS_MeshNode *n4);
618
619   bool registerElement(int ID, SMDS_MeshElement * element);
620
621   void addChildrenWithNodes(std::set<const SMDS_MeshElement*>& setOfChildren,
622                             const SMDS_MeshElement * element,
623                             std::set<const SMDS_MeshElement*>& nodes);
624
625   inline void adjustmyCellsCapacity(int ID)
626   {
627     assert(ID >= 0);
628     myElementIDFactory->adjustMaxId(ID);
629     if (ID >= myCells.size())
630       myCells.resize(ID+SMDS_Mesh::chunkSize,0);
631   }
632
633   inline void adjustBoundingBox(double x, double y, double z)
634   {
635     if (x > xmax) xmax = x;
636     else if (x < xmin) xmin = x;
637     if (y > ymax) ymax = y;
638     else if (y < ymin) ymin = y;
639     if (z > zmax) zmax = z;
640     else if (z < zmin) zmin = z;
641   }
642
643   // Fields PRIVATE
644
645   //! index of this SMDS_mesh in the static vector<SMDS_Mesh*> _meshList
646   int myMeshId;
647
648   //! actual nodes coordinates, cells definition and reverse connectivity are stored in a vtkUnstructuredGrid
649   SMDS_UnstructuredGrid*      myGrid;
650
651   //! Small objects like SMDS_MeshNode are allocated by chunks to limit memory costs of new
652   ObjectPool<SMDS_MeshNode>* myNodePool;
653
654   //! Small objects like SMDS_VtkVolume are allocated by chunks to limit memory costs of new
655   ObjectPool<SMDS_VtkVolume>* myVolumePool;
656   ObjectPool<SMDS_VtkFace>* myFacePool;
657   ObjectPool<SMDS_VtkEdge>* myEdgePool;
658
659   //! SMDS_MeshNodes refer to vtk nodes (vtk id = index in myNodes),store reference to this mesh, and subshape
660   SetOfNodes             myNodes;
661
662   //! SMDS_MeshCells refer to vtk cells (vtk id != index in myCells),store reference to this mesh, and subshape
663   SetOfCells             myCells;
664
665   //! for cells only: index = ID for SMDS users, value = ID in vtkUnstructuredGrid
666   //std::vector<int>       myCellIdSmdsToVtk;
667
668   //! for cells only: index = ID in vtkUnstructuredGrid, value = ID for SMDS users
669   std::vector<int>       myCellIdVtkToSmds;
670
671   SMDS_Mesh *            myParent;
672   std::list<SMDS_Mesh *> myChildren;
673   SMDS_MeshNodeIDFactory *myNodeIDFactory;
674   SMDS_MeshElementIDFactory *myElementIDFactory;
675   SMDS_MeshInfo          myInfo;
676
677   //! use a counter to keep track of modifications
678   unsigned long myModifTime, myCompactTime;
679
680   int myNodeMin;
681   int myNodeMax;
682
683   bool myHasConstructionEdges;
684   bool myHasConstructionFaces;
685   bool myHasInverseElements;
686
687   //! any add, remove or change of node or cell
688   bool myModified;
689
690   double xmin;
691   double xmax;
692   double ymin;
693   double ymax;
694   double zmin;
695   double zmax;
696 };
697
698
699 #endif