Salome HOME
23207: EDF 11232 SMESH: viscou layers builder fails at the very fist inflation
[modules/smesh.git] / src / SMESH / SMESH_Mesh.cxx
1 // Copyright (C) 2007-2016  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, or (at your option) any later version.
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 //  File   : SMESH_Mesh.cxx
24 //  Author : Paul RASCLE, EDF
25 //  Module : SMESH
26 //
27 #include "SMESH_Mesh.hxx"
28 #include "SMESH_MesherHelper.hxx"
29 #include "SMDS_MeshVolume.hxx"
30 #include "SMDS_SetIterator.hxx"
31 #include "SMESHDS_Document.hxx"
32 #include "SMESHDS_Group.hxx"
33 #include "SMESHDS_GroupOnGeom.hxx"
34 #include "SMESHDS_Script.hxx"
35 #include "SMESHDS_TSubMeshHolder.hxx"
36 #include "SMESH_Gen.hxx"
37 #include "SMESH_Group.hxx"
38 #include "SMESH_HypoFilter.hxx"
39 #include "SMESH_Hypothesis.hxx"
40 #include "SMESH_subMesh.hxx"
41
42 #include "utilities.h"
43
44 #include "DriverDAT_W_SMDS_Mesh.h"
45 #include "DriverGMF_Read.hxx"
46 #include "DriverGMF_Write.hxx"
47 #include "DriverMED_R_SMESHDS_Mesh.h"
48 #include "DriverMED_W_SMESHDS_Mesh.h"
49 #include "DriverSTL_R_SMDS_Mesh.h"
50 #include "DriverSTL_W_SMDS_Mesh.h"
51 #include "DriverUNV_R_SMDS_Mesh.h"
52 #include "DriverUNV_W_SMDS_Mesh.h"
53 #ifdef WITH_CGNS
54 #include "DriverCGNS_Read.hxx"
55 #include "DriverCGNS_Write.hxx"
56 #endif
57
58 #include <GEOMUtils.hxx>
59
60 #undef _Precision_HeaderFile
61 #include <BRepBndLib.hxx>
62 #include <BRepPrimAPI_MakeBox.hxx>
63 #include <Bnd_Box.hxx>
64 #include <TColStd_MapOfInteger.hxx>
65 #include <TopExp.hxx>
66 #include <TopExp_Explorer.hxx>
67 #include <TopTools_ListIteratorOfListOfShape.hxx>
68 #include <TopTools_ListOfShape.hxx>
69 #include <TopTools_MapOfShape.hxx>
70 #include <TopoDS_Iterator.hxx>
71
72 #include "SMESH_TryCatch.hxx" // include after OCCT headers!
73
74 #include "Utils_ExceptHandlers.hxx"
75
76 #ifndef WIN32
77 #include <boost/thread/thread.hpp>
78 #include <boost/bind.hpp>
79 #else 
80 #include <pthread.h>
81 #endif
82
83 using namespace std;
84
85 // maximum stored group name length in MED file
86 #define MAX_MED_GROUP_NAME_LENGTH 80
87
88 #ifdef _DEBUG_
89 static int MYDEBUG = 0;
90 #else
91 static int MYDEBUG = 0;
92 #endif
93
94 #define cSMESH_Hyp(h) static_cast<const SMESH_Hypothesis*>(h)
95
96 typedef SMESH_HypoFilter THypType;
97
98 class SMESH_Mesh::SubMeshHolder : public SMESHDS_TSubMeshHolder< SMESH_subMesh >
99 {
100 };
101
102 //=============================================================================
103 /*!
104  * 
105  */
106 //=============================================================================
107
108 SMESH_Mesh::SMESH_Mesh(int               theLocalId, 
109                        int               theStudyId, 
110                        SMESH_Gen*        theGen,
111                        bool              theIsEmbeddedMode,
112                        SMESHDS_Document* theDocument):
113   _groupId( 0 ), _nbSubShapes( 0 )
114 {
115   if(MYDEBUG) MESSAGE("SMESH_Mesh::SMESH_Mesh(int localId)");
116   _id            = theLocalId;
117   _studyId       = theStudyId;
118   _gen           = theGen;
119   _myDocument    = theDocument;
120   _myMeshDS      = theDocument->NewMesh(theIsEmbeddedMode,theLocalId);
121   _isShapeToMesh = false;
122   _isAutoColor   = false;
123   _isModified    = false;
124   _shapeDiagonal = 0.0;
125   _callUp        = NULL;
126   _myMeshDS->ShapeToMesh( PseudoShape() );
127   _subMeshHolder = new SubMeshHolder;
128 }
129
130 //================================================================================
131 /*!
132  * \brief Constructor of SMESH_Mesh being a base of some descendant class
133  */
134 //================================================================================
135
136 SMESH_Mesh::SMESH_Mesh():
137   _id(-1),
138   _studyId(-1),
139   _groupId( 0 ),
140   _nbSubShapes( 0 ),
141   _isShapeToMesh( false ),
142   _myDocument( 0 ),
143   _myMeshDS( 0 ),
144   _gen( 0 ),
145   _isAutoColor( false ),
146   _isModified( false ),
147   _shapeDiagonal( 0.0 ),
148   _callUp( 0 )
149 {
150   _subMeshHolder = new SubMeshHolder;
151 }
152
153 namespace
154 {
155 #ifndef WIN32
156   void deleteMeshDS(SMESHDS_Mesh* meshDS)
157   {
158     //cout << "deleteMeshDS( " << meshDS << endl;
159     delete meshDS;
160   }
161 #else
162   static void* deleteMeshDS(void* meshDS)
163   {
164     //cout << "deleteMeshDS( " << meshDS << endl;
165     SMESHDS_Mesh* m = (SMESHDS_Mesh*)meshDS;
166     if(m) {
167       delete m;
168     }
169     return 0;
170   }
171 #endif
172 }
173
174 //=============================================================================
175 /*!
176  *
177  */
178 //=============================================================================
179
180 SMESH_Mesh::~SMESH_Mesh()
181 {
182   if(MYDEBUG) MESSAGE("SMESH_Mesh::~SMESH_Mesh");
183
184   // avoid usual removal of elements while processing RemoveHypothesis( algo ) event
185   SMESHDS_SubMeshIteratorPtr smIt = _myMeshDS->SubMeshes();
186   while ( smIt->more() )
187     const_cast<SMESHDS_SubMesh*>( smIt->next() )->Clear();
188
189   // issue 0020340: EDF 1022 SMESH : Crash with FindNodeClosestTo in a second new study
190   //   Notify event listeners at least that something happens
191   if ( SMESH_subMesh * sm = GetSubMeshContaining(1))
192     sm->ComputeStateEngine( SMESH_subMesh::MESH_ENTITY_REMOVED );
193
194   // delete groups
195   map < int, SMESH_Group * >::iterator itg;
196   for (itg = _mapGroup.begin(); itg != _mapGroup.end(); itg++) {
197     SMESH_Group *aGroup = (*itg).second;
198     delete aGroup;
199   }
200   _mapGroup.clear();
201
202   // delete sub-meshes
203   delete _subMeshHolder;
204
205   if ( _callUp) delete _callUp;
206   _callUp = 0;
207
208   // remove self from studyContext
209   if ( _gen )
210   {
211     StudyContextStruct * studyContext = _gen->GetStudyContext( _studyId );
212     studyContext->mapMesh.erase( _id );
213   }
214   if ( _myDocument )
215     _myDocument->RemoveMesh( _id );
216   _myDocument = 0;
217
218   if ( _myMeshDS ) {
219     // delete _myMeshDS, in a thread in order not to block closing a study with large meshes
220 #ifndef WIN32
221     boost::thread aThread(boost::bind( & deleteMeshDS, _myMeshDS ));
222 #else
223     pthread_t thread;
224     int result=pthread_create(&thread, NULL, deleteMeshDS, (void*)_myMeshDS);
225 #endif
226   }
227 }
228
229 //================================================================================
230 /*!
231  * \brief Return true if a mesh with given id exists
232  */
233 //================================================================================
234
235 bool SMESH_Mesh::MeshExists( int meshId ) const
236 {
237   return _myDocument ? bool( _myDocument->GetMesh( meshId )) : false;
238 }
239
240 //================================================================================
241 /*!
242  * \brief Return a mesh by id
243  */
244 //================================================================================
245
246 SMESH_Mesh* SMESH_Mesh::FindMesh( int meshId ) const
247 {
248   if ( _id == meshId )
249     return (SMESH_Mesh*) this;
250
251   if ( StudyContextStruct *aStudyContext = _gen->GetStudyContext( _studyId ))
252   {
253     std::map < int, SMESH_Mesh * >::iterator i_m = aStudyContext->mapMesh.find( meshId );
254     if ( i_m != aStudyContext->mapMesh.end() )
255       return i_m->second;
256   }
257   return NULL;
258 }
259
260 //=============================================================================
261 /*!
262  * \brief Set geometry to be meshed
263  */
264 //=============================================================================
265
266 void SMESH_Mesh::ShapeToMesh(const TopoDS_Shape & aShape)
267 {
268   if(MYDEBUG) MESSAGE("SMESH_Mesh::ShapeToMesh");
269
270   if ( !aShape.IsNull() && _isShapeToMesh ) {
271     if ( aShape.ShapeType() != TopAbs_COMPOUND && // group contents is allowed to change
272          _myMeshDS->ShapeToMesh().ShapeType() != TopAbs_COMPOUND )
273       throw SALOME_Exception(LOCALIZED ("a shape to mesh has already been defined"));
274   }
275   // clear current data
276   if ( !_myMeshDS->ShapeToMesh().IsNull() )
277   {
278     // removal of a shape to mesh, delete objects referring to sub-shapes:
279     // - sub-meshes
280     _subMeshHolder->DeleteAll();
281     //  - groups on geometry
282     map <int, SMESH_Group *>::iterator i_gr = _mapGroup.begin();
283     while ( i_gr != _mapGroup.end() ) {
284       if ( dynamic_cast<SMESHDS_GroupOnGeom*>( i_gr->second->GetGroupDS() )) {
285         _myMeshDS->RemoveGroup( i_gr->second->GetGroupDS() );
286         delete i_gr->second;
287         _mapGroup.erase( i_gr++ );
288       }
289       else
290         i_gr++;
291     }
292     _mapAncestors.Clear();
293
294     // clear SMESHDS
295     TopoDS_Shape aNullShape;
296     _myMeshDS->ShapeToMesh( aNullShape );
297
298     _shapeDiagonal = 0.0;
299   }
300
301   // set a new geometry
302   if ( !aShape.IsNull() )
303   {
304     _myMeshDS->ShapeToMesh(aShape);
305     _isShapeToMesh = true;
306     _nbSubShapes = _myMeshDS->MaxShapeIndex();
307
308     // fill map of ancestors
309     fillAncestorsMap(aShape);
310   }
311   else
312   {
313     _isShapeToMesh = false;
314     _shapeDiagonal = 0.0;
315     _myMeshDS->ShapeToMesh( PseudoShape() );
316   }
317   _isModified = false;
318 }
319
320 //=======================================================================
321 /*!
322  * \brief Return geometry to be meshed. (It may be a PseudoShape()!)
323  */
324 //=======================================================================
325
326 TopoDS_Shape SMESH_Mesh::GetShapeToMesh() const
327 {
328   return _myMeshDS->ShapeToMesh();
329 }
330
331 //=======================================================================
332 /*!
333  * \brief Return a solid which is returned by GetShapeToMesh() if
334  *        a real geometry to be meshed was not set
335  */
336 //=======================================================================
337
338 const TopoDS_Solid& SMESH_Mesh::PseudoShape()
339 {
340   static TopoDS_Solid aSolid;
341   if ( aSolid.IsNull() )
342   {
343     aSolid = BRepPrimAPI_MakeBox(1,1,1);
344   }
345   return aSolid;
346 }
347
348 //=======================================================================
349 /*!
350  * \brief Return diagonal size of bounding box of a shape
351  */
352 //=======================================================================
353
354 double SMESH_Mesh::GetShapeDiagonalSize(const TopoDS_Shape & aShape)
355 {
356   if ( !aShape.IsNull() ) {
357     Bnd_Box Box;
358     // avoid too long waiting on large shapes. PreciseBoundingBox() was added
359     // to assure same result which else depends on presence of triangulation (IPAL52557).
360     const int maxNbFaces = 4000;
361     int nbFaces = 0;
362     for ( TopExp_Explorer f( aShape, TopAbs_FACE ); f.More() && nbFaces < maxNbFaces; f.Next() )
363       ++nbFaces;
364     bool isPrecise = false;
365     if ( nbFaces < maxNbFaces )
366       try {
367         GEOMUtils::PreciseBoundingBox( aShape, Box );
368         isPrecise = true;
369       }
370       catch (...) {
371         isPrecise = false;
372       }
373     if ( !isPrecise )
374     {
375       BRepBndLib::Add( aShape, Box );
376     }
377     if ( !Box.IsVoid() )
378       return sqrt( Box.SquareExtent() );
379   }
380   return 0;
381 }
382
383 //=======================================================================
384 /*!
385  * \brief Return diagonal size of bounding box of shape to mesh
386  */
387 //=======================================================================
388
389 double SMESH_Mesh::GetShapeDiagonalSize() const
390 {
391   if ( _shapeDiagonal == 0. && _isShapeToMesh )
392     const_cast<SMESH_Mesh*>(this)->_shapeDiagonal = GetShapeDiagonalSize( GetShapeToMesh() );
393
394   return _shapeDiagonal;
395 }
396
397 //================================================================================
398 /*!
399  * \brief Load mesh from study file
400  */
401 //================================================================================
402
403 void SMESH_Mesh::Load()
404 {
405   if (_callUp)
406     _callUp->Load();
407 }
408
409 //=======================================================================
410 /*!
411  * \brief Remove all nodes and elements
412  */
413 //=======================================================================
414
415 void SMESH_Mesh::Clear()
416 {
417   if ( HasShapeToMesh() ) // remove all nodes and elements
418   {
419     // clear mesh data
420     _myMeshDS->ClearMesh();
421
422     // update compute state of submeshes
423     if ( SMESH_subMesh *sm = GetSubMeshContaining( GetShapeToMesh() ) )
424     {
425       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
426       sm->ComputeSubMeshStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
427       sm->ComputeStateEngine( SMESH_subMesh::CLEAN ); // for event listeners (issue 0020918)
428       sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
429     }
430   }
431   else // remove only nodes/elements computed by algorithms
432   {
433     if ( SMESH_subMesh *sm = GetSubMeshContaining( GetShapeToMesh() ) )
434     {
435       sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
436       sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
437       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
438       sm->ComputeSubMeshStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
439     }
440   }
441   _isModified = false;
442 }
443
444 //=======================================================================
445 /*!
446  * \brief Remove all nodes and elements of indicated shape
447  */
448 //=======================================================================
449
450 void SMESH_Mesh::ClearSubMesh(const int theShapeId)
451 {
452   // clear sub-meshes; get ready to re-compute as a side-effect 
453   if ( SMESH_subMesh *sm = GetSubMeshContaining( theShapeId ) )
454   {
455     SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
456                                                              /*complexShapeFirst=*/false);
457     while ( smIt->more() )
458     {
459       sm = smIt->next();
460       TopAbs_ShapeEnum shapeType = sm->GetSubShape().ShapeType();      
461       if ( shapeType == TopAbs_VERTEX || shapeType < TopAbs_SOLID )
462         // all other shapes depends on vertices so they are already cleaned
463         sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
464       // to recompute even if failed
465       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
466     }
467   }
468 }
469
470 //=======================================================================
471 //function : UNVToMesh
472 //purpose  : 
473 //=======================================================================
474
475 int SMESH_Mesh::UNVToMesh(const char* theFileName)
476 {
477   if(MYDEBUG) MESSAGE("UNVToMesh - theFileName = "<<theFileName);
478   if(_isShapeToMesh)
479     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
480   _isShapeToMesh = false;
481   DriverUNV_R_SMDS_Mesh myReader;
482   myReader.SetMesh(_myMeshDS);
483   myReader.SetFile(theFileName);
484   myReader.SetMeshId(-1);
485   myReader.Perform();
486   if(MYDEBUG){
487     MESSAGE("UNVToMesh - _myMeshDS->NbNodes() = "<<_myMeshDS->NbNodes());
488     MESSAGE("UNVToMesh - _myMeshDS->NbEdges() = "<<_myMeshDS->NbEdges());
489     MESSAGE("UNVToMesh - _myMeshDS->NbFaces() = "<<_myMeshDS->NbFaces());
490     MESSAGE("UNVToMesh - _myMeshDS->NbVolumes() = "<<_myMeshDS->NbVolumes());
491   }
492   SMDS_MeshGroup* aGroup = (SMDS_MeshGroup*) myReader.GetGroup();
493   if (aGroup != 0) {
494     TGroupNamesMap aGroupNames = myReader.GetGroupNamesMap();
495     //const TGroupIdMap& aGroupId = myReader.GetGroupIdMap();
496     aGroup->InitSubGroupsIterator();
497     while (aGroup->MoreSubGroups()) {
498       SMDS_MeshGroup* aSubGroup = (SMDS_MeshGroup*) aGroup->NextSubGroup();
499       string aName = aGroupNames[aSubGroup];
500       int aId;
501
502       SMESH_Group* aSMESHGroup = AddGroup( aSubGroup->GetType(), aName.c_str(), aId );
503       if ( aSMESHGroup ) {
504         if(MYDEBUG) MESSAGE("UNVToMesh - group added: "<<aName);      
505         SMESHDS_Group* aGroupDS = dynamic_cast<SMESHDS_Group*>( aSMESHGroup->GetGroupDS() );
506         if ( aGroupDS ) {
507           aGroupDS->SetStoreName(aName.c_str());
508           aSubGroup->InitIterator();
509           const SMDS_MeshElement* aElement = 0;
510           while (aSubGroup->More()) {
511             aElement = aSubGroup->Next();
512             if (aElement) {
513               aGroupDS->SMDSGroup().Add(aElement);
514             }
515           }
516           if (aElement)
517             aGroupDS->SetType(aElement->GetType());
518         }
519       }
520     }
521   }
522   return 1;
523 }
524
525 //=======================================================================
526 //function : MEDToMesh
527 //purpose  : 
528 //=======================================================================
529
530 int SMESH_Mesh::MEDToMesh(const char* theFileName, const char* theMeshName)
531 {
532   if(MYDEBUG) MESSAGE("MEDToMesh - theFileName = "<<theFileName<<", mesh name = "<<theMeshName);
533   if(_isShapeToMesh)
534     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
535   _isShapeToMesh = false;
536   DriverMED_R_SMESHDS_Mesh myReader;
537   myReader.SetMesh(_myMeshDS);
538   myReader.SetMeshId(-1);
539   myReader.SetFile(theFileName);
540   myReader.SetMeshName(theMeshName);
541   Driver_Mesh::Status status = myReader.Perform();
542   if(MYDEBUG){
543     MESSAGE("MEDToMesh - _myMeshDS->NbNodes() = "<<_myMeshDS->NbNodes());
544     MESSAGE("MEDToMesh - _myMeshDS->NbEdges() = "<<_myMeshDS->NbEdges());
545     MESSAGE("MEDToMesh - _myMeshDS->NbFaces() = "<<_myMeshDS->NbFaces());
546     MESSAGE("MEDToMesh - _myMeshDS->NbVolumes() = "<<_myMeshDS->NbVolumes());
547   }
548 #ifdef _DEBUG_
549   SMESH_ComputeErrorPtr er = myReader.GetError();
550   if ( er && !er->IsOK() ) cout << er->myComment << endl;
551 #endif
552
553   // Reading groups (sub-meshes are out of scope of MED import functionality)
554   list<TNameAndType> aGroupNames = myReader.GetGroupNamesAndTypes();
555   if(MYDEBUG) MESSAGE("MEDToMesh - Nb groups = "<<aGroupNames.size()); 
556   int anId;
557   list<TNameAndType>::iterator name_type = aGroupNames.begin();
558   for ( ; name_type != aGroupNames.end(); name_type++ ) {
559     SMESH_Group* aGroup = AddGroup( name_type->second, name_type->first.c_str(), anId );
560     if ( aGroup ) {
561       if(MYDEBUG) MESSAGE("MEDToMesh - group added: "<<name_type->first.c_str());      
562       SMESHDS_Group* aGroupDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
563       if ( aGroupDS ) {
564         aGroupDS->SetStoreName( name_type->first.c_str() );
565         myReader.GetGroup( aGroupDS );
566       }
567     }
568   }
569   return (int) status;
570 }
571
572 //=======================================================================
573 //function : STLToMesh
574 //purpose  : 
575 //=======================================================================
576
577 int SMESH_Mesh::STLToMesh(const char* theFileName)
578 {
579   if(MYDEBUG) MESSAGE("STLToMesh - theFileName = "<<theFileName);
580   if(_isShapeToMesh)
581     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
582   _isShapeToMesh = false;
583   DriverSTL_R_SMDS_Mesh myReader;
584   myReader.SetMesh(_myMeshDS);
585   myReader.SetFile(theFileName);
586   myReader.SetMeshId(-1);
587   myReader.Perform();
588   if(MYDEBUG){
589     MESSAGE("STLToMesh - _myMeshDS->NbNodes() = "<<_myMeshDS->NbNodes());
590     MESSAGE("STLToMesh - _myMeshDS->NbEdges() = "<<_myMeshDS->NbEdges());
591     MESSAGE("STLToMesh - _myMeshDS->NbFaces() = "<<_myMeshDS->NbFaces());
592     MESSAGE("STLToMesh - _myMeshDS->NbVolumes() = "<<_myMeshDS->NbVolumes());
593   }
594   return 1;
595 }
596
597 //================================================================================
598 /*!
599  * \brief Reads the given mesh from the CGNS file
600  *  \param theFileName - name of the file
601  *  \retval int - Driver_Mesh::Status
602  */
603 //================================================================================
604
605 int SMESH_Mesh::CGNSToMesh(const char*  theFileName,
606                            const int    theMeshIndex,
607                            std::string& theMeshName)
608 {
609   int res = Driver_Mesh::DRS_FAIL;
610 #ifdef WITH_CGNS
611
612   DriverCGNS_Read myReader;
613   myReader.SetMesh(_myMeshDS);
614   myReader.SetFile(theFileName);
615   myReader.SetMeshId(theMeshIndex);
616   res = myReader.Perform();
617   theMeshName = myReader.GetMeshName();
618
619   // create groups
620   SynchronizeGroups();
621
622 #endif
623   return res;
624 }
625
626 //================================================================================
627 /*!
628  * \brief Fill its data by reading a GMF file
629  */
630 //================================================================================
631
632 SMESH_ComputeErrorPtr SMESH_Mesh::GMFToMesh(const char* theFileName,
633                                             bool        theMakeRequiredGroups)
634 {
635   DriverGMF_Read myReader;
636   myReader.SetMesh(_myMeshDS);
637   myReader.SetFile(theFileName);
638   myReader.SetMakeRequiredGroups( theMakeRequiredGroups );
639   myReader.Perform();
640   //theMeshName = myReader.GetMeshName();
641
642   // create groups
643   SynchronizeGroups();
644
645   return myReader.GetError();
646 }
647
648 //=============================================================================
649 /*!
650  * 
651  */
652 //=============================================================================
653
654 SMESH_Hypothesis::Hypothesis_Status
655 SMESH_Mesh::AddHypothesis(const TopoDS_Shape & aSubShape,
656                           int                  anHypId,
657                           std::string*         anError  ) throw(SALOME_Exception)
658 {
659   Unexpect aCatch(SalomeException);
660   if(MYDEBUG) MESSAGE("SMESH_Mesh::AddHypothesis");
661
662   if ( anError )
663     anError->clear();
664
665   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
666   if ( !subMesh || !subMesh->GetId())
667     return SMESH_Hypothesis::HYP_BAD_SUBSHAPE;
668
669   SMESH_Hypothesis *anHyp = GetHypothesis( anHypId );
670   if ( !anHyp )
671     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
672
673   bool isGlobalHyp = IsMainShape( aSubShape );
674
675   // NotConformAllowed can be only global
676   if ( !isGlobalHyp )
677   {
678     // NOTE: this is not a correct way to check a name of hypothesis,
679     // there should be an attribute of hypothesis saying that it can/can't
680     // be global/local
681     string hypName = anHyp->GetName();
682     if ( hypName == "NotConformAllowed" )
683     {
684       if(MYDEBUG) MESSAGE( "Hypotesis <NotConformAllowed> can be only global" );
685       return SMESH_Hypothesis::HYP_INCOMPATIBLE;
686     }
687   }
688
689   // shape
690
691   bool                     isAlgo = ( anHyp->GetType() != SMESHDS_Hypothesis::PARAM_ALGO );
692   SMESH_subMesh::algo_event event = isAlgo ? SMESH_subMesh::ADD_ALGO : SMESH_subMesh::ADD_HYP;
693
694   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
695
696   if ( anError && SMESH_Hypothesis::IsStatusFatal(ret) && subMesh->GetComputeError() )
697     *anError = subMesh->GetComputeError()->myComment;
698
699   // sub-shapes
700   if ( !SMESH_Hypothesis::IsStatusFatal(ret) &&
701        anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is added on father
702   {
703     event = isAlgo ? SMESH_subMesh::ADD_FATHER_ALGO : SMESH_subMesh::ADD_FATHER_HYP;
704
705     SMESH_Hypothesis::Hypothesis_Status ret2 =
706       subMesh->SubMeshesAlgoStateEngine(event, anHyp, /*exitOnFatal=*/true);
707     if (ret2 > ret)
708     {
709       ret = ret2;
710       if ( SMESH_Hypothesis::IsStatusFatal( ret ))
711       {
712         if ( anError && subMesh->GetComputeError() )
713           *anError = subMesh->GetComputeError()->myComment;
714         // remove anHyp
715         event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
716         subMesh->AlgoStateEngine(event, anHyp);
717       }
718     }
719
720     // check concurent hypotheses on ancestors
721     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !isGlobalHyp )
722     {
723       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
724       while ( smIt->more() ) {
725         SMESH_subMesh* sm = smIt->next();
726         if ( sm->IsApplicableHypotesis( anHyp )) {
727           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
728           if (ret2 > ret) {
729             ret = ret2;
730             break;
731           }
732         }
733       }
734     }
735   }
736   HasModificationsToDiscard(); // to reset _isModified flag if a mesh becomes empty
737
738   if(MYDEBUG) subMesh->DumpAlgoState(true);
739   if(MYDEBUG) SCRUTE(ret);
740   return ret;
741 }
742
743 //=============================================================================
744 /*!
745  * 
746  */
747 //=============================================================================
748
749 SMESH_Hypothesis::Hypothesis_Status
750   SMESH_Mesh::RemoveHypothesis(const TopoDS_Shape & aSubShape,
751                                int anHypId)throw(SALOME_Exception)
752 {
753   Unexpect aCatch(SalomeException);
754   if(MYDEBUG) MESSAGE("SMESH_Mesh::RemoveHypothesis");
755   
756   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
757   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
758     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
759   
760   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
761   if(MYDEBUG) {
762     SCRUTE(anHyp->GetType());
763   }
764   
765   // shape 
766   
767   bool                     isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
768   SMESH_subMesh::algo_event event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
769
770   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
771
772   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
773
774   // there may appear concurrent hyps that were covered by the removed hyp
775   if (ret < SMESH_Hypothesis::HYP_CONCURENT &&
776       subMesh->IsApplicableHypotesis( anHyp ) &&
777       subMesh->CheckConcurentHypothesis( anHyp->GetType() ) != SMESH_Hypothesis::HYP_OK)
778     ret = SMESH_Hypothesis::HYP_CONCURENT;
779
780   // sub-shapes
781   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
782       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is removed from father
783   {
784     event = isAlgo ? SMESH_subMesh::REMOVE_FATHER_ALGO : SMESH_subMesh::REMOVE_FATHER_HYP;
785
786     SMESH_Hypothesis::Hypothesis_Status ret2 =
787       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
788     if (ret2 > ret) // more severe
789       ret = ret2;
790
791     // check concurent hypotheses on ancestors
792     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !IsMainShape( aSubShape ) )
793     {
794       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
795       while ( smIt->more() ) {
796         SMESH_subMesh* sm = smIt->next();
797         if ( sm->IsApplicableHypotesis( anHyp )) {
798           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
799           if (ret2 > ret) {
800             ret = ret2;
801             break;
802           }
803         }
804       }
805     }
806   }
807
808   HasModificationsToDiscard(); // to reset _isModified flag if mesh become empty
809
810   if(MYDEBUG) subMesh->DumpAlgoState(true);
811   if(MYDEBUG) SCRUTE(ret);
812   return ret;
813 }
814
815 //=============================================================================
816 /*!
817  * 
818  */
819 //=============================================================================
820
821 const list<const SMESHDS_Hypothesis*>&
822 SMESH_Mesh::GetHypothesisList(const TopoDS_Shape & aSubShape) const
823   throw(SALOME_Exception)
824 {
825   return _myMeshDS->GetHypothesis(aSubShape);
826 }
827
828 //=======================================================================
829 /*!
830  * \brief Return the hypothesis assigned to the shape
831  *  \param aSubShape    - the shape to check
832  *  \param aFilter      - the hypothesis filter
833  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
834  *  \param assignedTo   - to return the shape the found hypo is assigned to
835  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
836  */
837 //=======================================================================
838
839 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const TopoDS_Shape &    aSubShape,
840                                                    const SMESH_HypoFilter& aFilter,
841                                                    const bool              andAncestors,
842                                                    TopoDS_Shape*           assignedTo) const
843 {
844   return GetHypothesis( const_cast< SMESH_Mesh* >(this)->GetSubMesh( aSubShape ),
845                         aFilter, andAncestors, assignedTo );
846 }
847
848 //=======================================================================
849 /*!
850  * \brief Return the hypothesis assigned to the shape of a sub-mesh
851  *  \param aSubMesh     - the sub-mesh to check
852  *  \param aFilter      - the hypothesis filter
853  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
854  *  \param assignedTo   - to return the shape the found hypo is assigned to
855  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
856  */
857 //=======================================================================
858
859 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const SMESH_subMesh *   aSubMesh,
860                                                    const SMESH_HypoFilter& aFilter,
861                                                    const bool              andAncestors,
862                                                    TopoDS_Shape*           assignedTo) const
863 {
864   if ( !aSubMesh ) return 0;
865
866   {
867     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
868     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
869     list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
870     for ( ; hyp != hypList.end(); hyp++ ) {
871       const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
872       if ( aFilter.IsOk( h, aSubShape)) {
873         if ( assignedTo ) *assignedTo = aSubShape;
874         return h;
875       }
876     }
877   }
878   if ( andAncestors )
879   {
880     // user sorted submeshes of ancestors, according to stored submesh priority
881     std::vector< SMESH_subMesh * > & ancestors =
882       const_cast< std::vector< SMESH_subMesh * > & > ( aSubMesh->GetAncestors() );
883     SortByMeshOrder( ancestors );
884
885     vector<SMESH_subMesh*>::const_iterator smIt = ancestors.begin(); 
886     for ( ; smIt != ancestors.end(); smIt++ )
887     {
888       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
889       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
890       list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
891       for ( ; hyp != hypList.end(); hyp++ ) {
892         const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
893         if (aFilter.IsOk( h, curSh )) {
894           if ( assignedTo ) *assignedTo = curSh;
895           return h;
896         }
897       }
898     }
899   }
900   return 0;
901 }
902
903 //================================================================================
904 /*!
905  * \brief Return hypotheses assigned to the shape
906   * \param aSubShape - the shape to check
907   * \param aFilter - the hypothesis filter
908   * \param aHypList - the list of the found hypotheses
909   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
910   * \retval int - number of unique hypos in aHypList
911  */
912 //================================================================================
913
914 int SMESH_Mesh::GetHypotheses(const TopoDS_Shape &                aSubShape,
915                               const SMESH_HypoFilter&             aFilter,
916                               list <const SMESHDS_Hypothesis * >& aHypList,
917                               const bool                          andAncestors,
918                               list< TopoDS_Shape > *              assignedTo/*=0*/) const
919 {
920   return GetHypotheses( const_cast< SMESH_Mesh* >(this)->GetSubMesh( aSubShape ),
921                         aFilter, aHypList, andAncestors, assignedTo );
922 }
923
924 //================================================================================
925 /*!
926  * \brief Return hypotheses assigned to the shape of a sub-mesh
927   * \param aSubShape - the sub-mesh to check
928   * \param aFilter - the hypothesis filter
929   * \param aHypList - the list of the found hypotheses
930   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
931   * \retval int - number of unique hypos in aHypList
932  */
933 //================================================================================
934
935 int SMESH_Mesh::GetHypotheses(const SMESH_subMesh *               aSubMesh,
936                               const SMESH_HypoFilter&             aFilter,
937                               list <const SMESHDS_Hypothesis * >& aHypList,
938                               const bool                          andAncestors,
939                               list< TopoDS_Shape > *              assignedTo/*=0*/) const
940 {
941   if ( !aSubMesh ) return 0;
942
943   set<string> hypTypes; // to exclude same type hypos from the result list
944   int nbHyps = 0;
945
946   // only one main hypothesis is allowed
947   bool mainHypFound = false;
948
949   // fill in hypTypes
950   list<const SMESHDS_Hypothesis*>::const_iterator hyp;
951   for ( hyp = aHypList.begin(); hyp != aHypList.end(); hyp++ ) {
952     if ( hypTypes.insert( (*hyp)->GetName() ).second )
953       nbHyps++;
954     if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
955       mainHypFound = true;
956   }
957
958   // get hypos from aSubShape
959   {
960     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
961     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
962     for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
963     {
964       const SMESH_Hypothesis* h = cSMESH_Hyp( *hyp );
965       if (( aFilter.IsOk( h, aSubShape )) &&
966           ( h->IsAuxiliary() || !mainHypFound ) &&
967           ( h->IsAuxiliary() || hypTypes.insert( h->GetName() ).second ))
968       {
969         aHypList.push_back( *hyp );
970         nbHyps++;
971         if ( !h->IsAuxiliary() )
972           mainHypFound = true;
973         if ( assignedTo ) assignedTo->push_back( aSubShape );
974       }
975     }
976   }
977
978   // get hypos from ancestors of aSubShape
979   if ( andAncestors )
980   {
981     // user sorted submeshes of ancestors, according to stored submesh priority
982     std::vector< SMESH_subMesh * > & ancestors =
983       const_cast< std::vector< SMESH_subMesh * > & > ( aSubMesh->GetAncestors() );
984     SortByMeshOrder( ancestors );
985
986     vector<SMESH_subMesh*>::const_iterator smIt = ancestors.begin();
987     for ( ; smIt != ancestors.end(); smIt++ )
988     {
989       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
990       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
991       for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
992       {
993         const SMESH_Hypothesis* h = cSMESH_Hyp( *hyp );
994         if (( aFilter.IsOk( h, curSh )) &&
995             ( h->IsAuxiliary() || !mainHypFound ) &&
996             ( h->IsAuxiliary() || hypTypes.insert( h->GetName() ).second ))
997         {
998           aHypList.push_back( *hyp );
999           nbHyps++;
1000           if ( !h->IsAuxiliary() )
1001             mainHypFound = true;
1002           if ( assignedTo ) assignedTo->push_back( curSh );
1003         }
1004       }
1005     }
1006   }
1007   return nbHyps;
1008 }
1009
1010 //================================================================================
1011 /*!
1012  * \brief Return a hypothesis by its ID
1013  */
1014 //================================================================================
1015
1016 SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const int anHypId) const
1017 {
1018   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
1019   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
1020     return NULL;
1021
1022   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
1023   return anHyp;
1024 }
1025
1026 //=============================================================================
1027 /*!
1028  * 
1029  */
1030 //=============================================================================
1031
1032 const list<SMESHDS_Command*> & SMESH_Mesh::GetLog() throw(SALOME_Exception)
1033 {
1034   Unexpect aCatch(SalomeException);
1035   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetLog");
1036   return _myMeshDS->GetScript()->GetCommands();
1037 }
1038
1039 //=============================================================================
1040 /*!
1041  * 
1042  */
1043 //=============================================================================
1044 void SMESH_Mesh::ClearLog() throw(SALOME_Exception)
1045 {
1046   Unexpect aCatch(SalomeException);
1047   if(MYDEBUG) MESSAGE("SMESH_Mesh::ClearLog");
1048   _myMeshDS->GetScript()->Clear();
1049 }
1050
1051 //=============================================================================
1052 /*!
1053  * Get or Create the SMESH_subMesh object implementation
1054  */
1055 //=============================================================================
1056
1057 SMESH_subMesh *SMESH_Mesh::GetSubMesh(const TopoDS_Shape & aSubShape)
1058   throw(SALOME_Exception)
1059 {
1060   int index = _myMeshDS->ShapeToIndex(aSubShape);
1061   if ( !index && aSubShape.IsNull() )
1062     return 0;
1063
1064   // for submeshes on GEOM Group
1065   if (( !index || index > _nbSubShapes ) && aSubShape.ShapeType() == TopAbs_COMPOUND ) {
1066     TopoDS_Iterator it( aSubShape );
1067     if ( it.More() )
1068     {
1069       index = _myMeshDS->AddCompoundSubmesh( aSubShape, it.Value().ShapeType() );
1070       // fill map of Ancestors
1071       while ( _nbSubShapes < index )
1072         fillAncestorsMap( _myMeshDS->IndexToShape( ++_nbSubShapes ));
1073     }
1074   }
1075   // if ( !index )
1076   //   return NULL; // neither sub-shape nor a group
1077
1078   SMESH_subMesh* aSubMesh = _subMeshHolder->Get( index );
1079   if ( !aSubMesh )
1080   {
1081     aSubMesh = new SMESH_subMesh(index, this, _myMeshDS, aSubShape);
1082     _subMeshHolder->Add( index, aSubMesh );
1083
1084     // include non-computable sub-meshes in SMESH_subMesh::_ancestors of sub-submeshes
1085     switch ( aSubShape.ShapeType() ) {
1086     case TopAbs_COMPOUND:
1087     case TopAbs_WIRE:
1088     case TopAbs_SHELL:
1089       for ( TopoDS_Iterator subIt( aSubShape ); subIt.More(); subIt.Next() )
1090       {
1091         SMESH_subMesh* sm = GetSubMesh( subIt.Value() );
1092         SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*inclideSelf=*/true);
1093         while ( smIt->more() )
1094           smIt->next()->ClearAncestors();
1095       }
1096     default:;
1097     }
1098   }
1099   return aSubMesh;
1100 }
1101
1102 //=============================================================================
1103 /*!
1104  * Get the SMESH_subMesh object implementation. Dont create it, return null
1105  * if it does not exist.
1106  */
1107 //=============================================================================
1108
1109 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const TopoDS_Shape & aSubShape) const
1110   throw(SALOME_Exception)
1111 {
1112   int index = _myMeshDS->ShapeToIndex(aSubShape);
1113   return GetSubMeshContaining( index );
1114 }
1115
1116 //=============================================================================
1117 /*!
1118  * Get the SMESH_subMesh object implementation. Dont create it, return null
1119  * if it does not exist.
1120  */
1121 //=============================================================================
1122
1123 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const int aShapeID) const
1124 throw(SALOME_Exception)
1125 {
1126   SMESH_subMesh *aSubMesh = _subMeshHolder->Get( aShapeID );
1127
1128   return aSubMesh;
1129 }
1130
1131 //================================================================================
1132 /*!
1133  * \brief Return sub-meshes of groups containing the given sub-shape
1134  */
1135 //================================================================================
1136
1137 list<SMESH_subMesh*>
1138 SMESH_Mesh::GetGroupSubMeshesContaining(const TopoDS_Shape & aSubShape) const
1139   throw(SALOME_Exception)
1140 {
1141   list<SMESH_subMesh*> found;
1142
1143   SMESH_subMesh * subMesh = GetSubMeshContaining(aSubShape);
1144   if ( !subMesh )
1145     return found;
1146
1147   // sub-meshes of groups have max IDs, so search from the map end
1148   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator( /*reverse=*/true ) );
1149   while ( smIt->more() ) {
1150     SMESH_subMesh*    sm = smIt->next();
1151     SMESHDS_SubMesh * ds = sm->GetSubMeshDS();
1152     if ( ds && ds->IsComplexSubmesh() ) {
1153       if ( SMESH_MesherHelper::IsSubShape( aSubShape, sm->GetSubShape() ))
1154       {
1155         found.push_back( sm );
1156         //break;
1157       }
1158     } else {
1159       break; // the rest sub-meshes are not those of groups
1160     }
1161   }
1162
1163   if ( found.empty() ) // maybe the main shape is a COMPOUND (issue 0021530)
1164   {
1165     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1) )
1166       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1167       {
1168         TopoDS_Iterator it( mainSM->GetSubShape() );
1169         if ( it.Value().ShapeType() == aSubShape.ShapeType() &&
1170              SMESH_MesherHelper::IsSubShape( aSubShape, mainSM->GetSubShape() ))
1171           found.push_back( mainSM );
1172       }
1173   }
1174   else // issue 0023068
1175   {
1176     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1) )
1177       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1178         found.push_back( mainSM );
1179   }
1180   return found;
1181 }
1182 //=======================================================================
1183 //function : IsUsedHypothesis
1184 //purpose  : Return True if anHyp is used to mesh aSubShape
1185 //=======================================================================
1186
1187 bool SMESH_Mesh::IsUsedHypothesis(SMESHDS_Hypothesis * anHyp,
1188                                   const SMESH_subMesh* aSubMesh)
1189 {
1190   SMESH_Hypothesis* hyp = static_cast<SMESH_Hypothesis*>(anHyp);
1191
1192   // check if anHyp can be used to mesh aSubMesh
1193   if ( !aSubMesh || !aSubMesh->IsApplicableHypotesis( hyp ))
1194     return false;
1195
1196   SMESH_Algo *algo = aSubMesh->GetAlgo();
1197
1198   // algorithm
1199   if (anHyp->GetType() > SMESHDS_Hypothesis::PARAM_ALGO)
1200     return ( anHyp == algo );
1201
1202   // algorithm parameter
1203   if (algo)
1204   {
1205     // look trough hypotheses used by algo
1206     const SMESH_HypoFilter* hypoKind;
1207     if (( hypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() ))) {
1208       list <const SMESHDS_Hypothesis * > usedHyps;
1209       if ( GetHypotheses( aSubMesh, *hypoKind, usedHyps, true ))
1210         return ( find( usedHyps.begin(), usedHyps.end(), anHyp ) != usedHyps.end() );
1211     }
1212   }
1213
1214   return false;
1215 }
1216
1217 //=============================================================================
1218 /*!
1219  *
1220  */
1221 //=============================================================================
1222
1223 // const list < SMESH_subMesh * >&
1224 // SMESH_Mesh::GetSubMeshUsingHypothesis(SMESHDS_Hypothesis * anHyp)
1225 //   throw(SALOME_Exception)
1226 // {
1227 //   _subMeshesUsingHypothesisList.clear();
1228 //   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1229 //   while ( smIt->more() )
1230 //   {
1231 //     SMESH_subMesh* aSubMesh = smIt->next();
1232 //     if ( IsUsedHypothesis ( anHyp, aSubMesh ))
1233 //       _subMeshesUsingHypothesisList.push_back( aSubMesh );
1234 //   }
1235 //   return _subMeshesUsingHypothesisList;
1236 // }
1237
1238 //=======================================================================
1239 //function : NotifySubMeshesHypothesisModification
1240 //purpose  : Say all submeshes using theChangedHyp that it has been modified
1241 //=======================================================================
1242
1243 void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* hyp)
1244 {
1245   Unexpect aCatch(SalomeException);
1246
1247   if ( !GetMeshDS()->IsUsedHypothesis( hyp ))
1248     return;
1249
1250   if (_callUp)
1251     _callUp->HypothesisModified();
1252
1253   SMESH_Algo *algo;
1254   const SMESH_HypoFilter* compatibleHypoKind;
1255   list <const SMESHDS_Hypothesis * > usedHyps;
1256
1257   // keep sub-meshes not to miss ones whose state can change due to notifying others
1258   vector< SMESH_subMesh* > smToNotify;
1259
1260   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1261   while ( smIt->more() )
1262   {
1263     SMESH_subMesh* aSubMesh = smIt->next();
1264
1265     // if aSubMesh meshing depends on hyp,
1266     // we call aSubMesh->AlgoStateEngine( MODIF_HYP, hyp ) that causes either
1267     // 1) clearing of already computed aSubMesh or
1268     // 2) changing algo_state from MISSING_HYP to HYP_OK when parameters of hyp becomes valid,
1269     // other possible changes are not interesting. (IPAL0052457 - assigning hyp performance pb)
1270     if ( aSubMesh->GetComputeState() != SMESH_subMesh::COMPUTE_OK &&
1271          aSubMesh->GetComputeState() != SMESH_subMesh::FAILED_TO_COMPUTE &&
1272          aSubMesh->GetAlgoState()    != SMESH_subMesh::MISSING_HYP &&
1273          !hyp->DataDependOnParams() )
1274       continue;
1275
1276     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
1277
1278     if (( aSubMesh->IsApplicableHypotesis( hyp )) &&
1279         ( algo = aSubMesh->GetAlgo() )            &&
1280         ( compatibleHypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() )) &&
1281         ( compatibleHypoKind->IsOk( hyp, aSubShape )))
1282     {
1283       // check if hyp is used by algo
1284       usedHyps.clear();
1285       if ( GetHypotheses( aSubMesh, *compatibleHypoKind, usedHyps, true ) &&
1286            find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() )
1287       {
1288         smToNotify.push_back( aSubMesh );
1289       }
1290     }
1291   }
1292
1293   for ( size_t i = 0; i < smToNotify.size(); ++i )
1294   {
1295     smToNotify[i]->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
1296                                    const_cast< SMESH_Hypothesis*>( hyp ));
1297   }
1298
1299   HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
1300   GetMeshDS()->Modified();
1301 }
1302
1303 //=============================================================================
1304 /*!
1305  *  Auto color functionality
1306  */
1307 //=============================================================================
1308 void SMESH_Mesh::SetAutoColor(bool theAutoColor) throw(SALOME_Exception)
1309 {
1310   Unexpect aCatch(SalomeException);
1311   _isAutoColor = theAutoColor;
1312 }
1313
1314 bool SMESH_Mesh::GetAutoColor() throw(SALOME_Exception)
1315 {
1316   Unexpect aCatch(SalomeException);
1317   return _isAutoColor;
1318 }
1319
1320 //=======================================================================
1321 //function : SetIsModified
1322 //purpose  : Set the flag meaning that the mesh has been edited "manually"
1323 //=======================================================================
1324
1325 void SMESH_Mesh::SetIsModified(bool isModified)
1326 {
1327   _isModified = isModified;
1328
1329   if ( _isModified )
1330     // check if mesh becomes empty as result of modification
1331     HasModificationsToDiscard();
1332 }
1333
1334 //=======================================================================
1335 //function : HasModificationsToDiscard
1336 //purpose  : Return true if the mesh has been edited since a total re-compute
1337 //           and those modifications may prevent successful partial re-compute.
1338 //           As a side effect reset _isModified flag if mesh is empty
1339 //issue    : 0020693
1340 //=======================================================================
1341
1342 bool SMESH_Mesh::HasModificationsToDiscard() const
1343 {
1344   if ( ! _isModified )
1345     return false;
1346
1347   // return true if the next Compute() will be partial and
1348   // existing but changed elements may prevent successful re-compute
1349   bool hasComputed = false, hasNotComputed = false;
1350   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1351   while ( smIt->more() )
1352   {
1353     const SMESH_subMesh* aSubMesh = smIt->next();
1354     switch ( aSubMesh->GetSubShape().ShapeType() )
1355     {
1356     case TopAbs_EDGE:
1357     case TopAbs_FACE:
1358     case TopAbs_SOLID:
1359       if ( aSubMesh->IsMeshComputed() )
1360         hasComputed = true;
1361       else
1362         hasNotComputed = true;
1363       if ( hasComputed && hasNotComputed)
1364         return true;
1365
1366     default:;
1367     }
1368   }
1369   if ( NbNodes() < 1 )
1370     const_cast<SMESH_Mesh*>(this)->_isModified = false;
1371
1372   return false;
1373 }
1374
1375 //================================================================================
1376 /*!
1377  * \brief Check if any groups of the same type have equal names
1378  */
1379 //================================================================================
1380
1381 bool SMESH_Mesh::HasDuplicatedGroupNamesMED()
1382 {
1383   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1384   map< SMDSAbs_ElementType, set<string> > aGroupNames;
1385   for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1386   {
1387     SMESH_Group* aGroup = it->second;
1388     SMDSAbs_ElementType aType = aGroup->GetGroupDS()->GetType();
1389     string aGroupName = aGroup->GetName();
1390     aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1391     if (!aGroupNames[aType].insert(aGroupName).second)
1392       return true;
1393   }
1394
1395   return false;
1396 }
1397
1398 //================================================================================
1399 /*!
1400  * \brief Export the mesh to a med file
1401  *  \param [in] file - name of the MED file
1402  *  \param [in] theMeshName - name of this mesh
1403  *  \param [in] theAutoGroups - boolean parameter for creating/not creating
1404  *              the groups Group_On_All_Nodes, Group_On_All_Faces, ... ;
1405  *              the typical use is auto_groups=false.
1406  *  \param [in] theVersion - defines the version of format of MED file, that will be created
1407  *  \param [in] meshPart - mesh data to export
1408  *  \param [in] theAutoDimension - if \c true, a space dimension of a MED mesh can be either
1409      *         - 1D if all mesh nodes lie on OX coordinate axis, or
1410      *         - 2D if all mesh nodes lie on XOY coordinate plane, or
1411      *         - 3D in the rest cases.
1412      *         If \a theAutoDimension is \c false, the space dimension is always 3.
1413  *  \return int - mesh index in the file
1414  */
1415 //================================================================================
1416
1417 void SMESH_Mesh::ExportMED(const char *        file, 
1418                            const char*         theMeshName, 
1419                            bool                theAutoGroups,
1420                            int                 theVersion,
1421                            const SMESHDS_Mesh* meshPart,
1422                            bool                theAutoDimension,
1423                            bool                theAddODOnVertices)
1424   throw(SALOME_Exception)
1425 {
1426   SMESH_TRY;
1427
1428   DriverMED_W_SMESHDS_Mesh myWriter;
1429   myWriter.SetFile         ( file, MED::EVersion(theVersion) );
1430   myWriter.SetMesh         ( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS   );
1431   myWriter.SetAutoDimension( theAutoDimension );
1432   myWriter.AddODOnVertices ( theAddODOnVertices );
1433   if ( !theMeshName ) 
1434     myWriter.SetMeshId     ( _id         );
1435   else {
1436     myWriter.SetMeshId     ( -1          );
1437     myWriter.SetMeshName   ( theMeshName );
1438   }
1439
1440   if ( theAutoGroups ) {
1441     myWriter.AddGroupOfNodes();
1442     myWriter.AddGroupOfEdges();
1443     myWriter.AddGroupOfFaces();
1444     myWriter.AddGroupOfVolumes();
1445   }
1446
1447   // Pass groups to writer. Provide unique group names.
1448   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1449   if ( !meshPart )
1450   {
1451     map< SMDSAbs_ElementType, set<string> > aGroupNames;
1452     char aString [256];
1453     int maxNbIter = 10000; // to guarantee cycle finish
1454     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1455       SMESH_Group*       aGroup   = it->second;
1456       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1457       if ( aGroupDS ) {
1458         SMDSAbs_ElementType aType = aGroupDS->GetType();
1459         string aGroupName0 = aGroup->GetName();
1460         aGroupName0.resize(MAX_MED_GROUP_NAME_LENGTH);
1461         string aGroupName = aGroupName0;
1462         for (int i = 1; !aGroupNames[aType].insert(aGroupName).second && i < maxNbIter; i++) {
1463           sprintf(&aString[0], "GR_%d_%s", i, aGroupName0.c_str());
1464           aGroupName = aString;
1465           aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1466         }
1467         aGroupDS->SetStoreName( aGroupName.c_str() );
1468         myWriter.AddGroup( aGroupDS );
1469       }
1470     }
1471   }
1472   // Perform export
1473   myWriter.Perform();
1474
1475   SMESH_CATCH( SMESH::throwSalomeEx );
1476 }
1477
1478 //================================================================================
1479 /*!
1480  * \brief Export the mesh to a SAUV file
1481  */
1482 //================================================================================
1483
1484 void SMESH_Mesh::ExportSAUV(const char *file, 
1485                             const char* theMeshName, 
1486                             bool theAutoGroups)
1487   throw(SALOME_Exception)
1488 {
1489   std::string medfilename(file);
1490   medfilename += ".med";
1491   std::string cmd;
1492 #ifdef WIN32
1493   cmd = "%PYTHONBIN% ";
1494 #else
1495   cmd = "python ";
1496 #endif
1497   cmd += "-c \"";
1498   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1499   cmd += "\"";
1500   system(cmd.c_str());
1501   ExportMED(medfilename.c_str(), theMeshName, theAutoGroups, 1);
1502 #ifdef WIN32
1503   cmd = "%PYTHONBIN% ";
1504 #else
1505   cmd = "python ";
1506 #endif
1507   cmd += "-c \"";
1508   cmd += "from medutilities import convert ; convert(r'" + medfilename + "', 'MED', 'GIBI', 1, r'" + file + "')";
1509   cmd += "\"";
1510   system(cmd.c_str());
1511 #ifdef WIN32
1512   cmd = "%PYTHONBIN% ";
1513 #else
1514   cmd = "python ";
1515 #endif
1516   cmd += "-c \"";
1517   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1518   cmd += "\"";
1519   system(cmd.c_str());
1520 }
1521
1522 //================================================================================
1523 /*!
1524  * \brief Export the mesh to a DAT file
1525  */
1526 //================================================================================
1527
1528 void SMESH_Mesh::ExportDAT(const char *        file,
1529                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1530 {
1531   Unexpect aCatch(SalomeException);
1532   DriverDAT_W_SMDS_Mesh myWriter;
1533   myWriter.SetFile( file );
1534   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1535   myWriter.SetMeshId(_id);
1536   myWriter.Perform();
1537 }
1538
1539 //================================================================================
1540 /*!
1541  * \brief Export the mesh to an UNV file
1542  */
1543 //================================================================================
1544
1545 void SMESH_Mesh::ExportUNV(const char *        file,
1546                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1547 {
1548   Unexpect aCatch(SalomeException);
1549   DriverUNV_W_SMDS_Mesh myWriter;
1550   myWriter.SetFile( file );
1551   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1552   myWriter.SetMeshId(_id);
1553   //  myWriter.SetGroups(_mapGroup);
1554
1555   if ( !meshPart )
1556   {
1557     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1558       SMESH_Group*       aGroup   = it->second;
1559       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1560       if ( aGroupDS ) {
1561         string aGroupName = aGroup->GetName();
1562         aGroupDS->SetStoreName( aGroupName.c_str() );
1563         myWriter.AddGroup( aGroupDS );
1564       }
1565     }
1566   }
1567   myWriter.Perform();
1568 }
1569
1570 //================================================================================
1571 /*!
1572  * \brief Export the mesh to an STL file
1573  */
1574 //================================================================================
1575
1576 void SMESH_Mesh::ExportSTL(const char *        file,
1577                            const bool          isascii,
1578                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1579 {
1580   Unexpect aCatch(SalomeException);
1581   DriverSTL_W_SMDS_Mesh myWriter;
1582   myWriter.SetFile( file );
1583   myWriter.SetIsAscii( isascii );
1584   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS);
1585   myWriter.SetMeshId(_id);
1586   myWriter.Perform();
1587 }
1588
1589 //================================================================================
1590 /*!
1591  * \brief Export the mesh to the CGNS file
1592  */
1593 //================================================================================
1594
1595 void SMESH_Mesh::ExportCGNS(const char *        file,
1596                             const SMESHDS_Mesh* meshDS,
1597                             const char *        meshName)
1598 {
1599   int res = Driver_Mesh::DRS_FAIL;
1600 #ifdef WITH_CGNS
1601   DriverCGNS_Write myWriter;
1602   myWriter.SetFile( file );
1603   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1604   myWriter.SetMeshName( SMESH_Comment("Mesh_") << meshDS->GetPersistentId());
1605   if ( meshName && meshName[0] )
1606     myWriter.SetMeshName( meshName );
1607   res = myWriter.Perform();
1608 #endif
1609   if ( res != Driver_Mesh::DRS_OK )
1610     throw SALOME_Exception("Export failed");
1611 }
1612
1613 //================================================================================
1614 /*!
1615  * \brief Export the mesh to a GMF file
1616  */
1617 //================================================================================
1618
1619 void SMESH_Mesh::ExportGMF(const char *        file,
1620                            const SMESHDS_Mesh* meshDS,
1621                            bool                withRequiredGroups)
1622 {
1623   DriverGMF_Write myWriter;
1624   myWriter.SetFile( file );
1625   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1626   myWriter.SetExportRequiredGroups( withRequiredGroups );
1627
1628   myWriter.Perform();
1629 }
1630
1631 //================================================================================
1632 /*!
1633  * \brief Return a ratio of "compute cost" of computed sub-meshes to the whole
1634  *        "compute cost".
1635  */
1636 //================================================================================
1637
1638 double SMESH_Mesh::GetComputeProgress() const
1639 {
1640   double totalCost = 1e-100, computedCost = 0;
1641   const SMESH_subMesh* curSM = _gen->GetCurrentSubMesh();
1642
1643   // get progress of a current algo
1644   TColStd_MapOfInteger currentSubIds; 
1645   if ( curSM )
1646     if ( SMESH_Algo* algo = curSM->GetAlgo() )
1647     {
1648       int algoNotDoneCost = 0, algoDoneCost = 0;
1649       const std::vector<SMESH_subMesh*>& smToCompute = algo->SubMeshesToCompute();
1650       for ( size_t i = 0; i < smToCompute.size(); ++i )
1651       {
1652         if ( smToCompute[i]->IsEmpty() )
1653           algoNotDoneCost += smToCompute[i]->GetComputeCost();
1654         else
1655           algoDoneCost += smToCompute[i]->GetComputeCost();
1656         currentSubIds.Add( smToCompute[i]->GetId() );
1657       }
1658       double rate = 0;
1659       try
1660       {
1661         OCC_CATCH_SIGNALS;
1662         rate = algo->GetProgress();
1663       }
1664       catch (...) {
1665 #ifdef _DEBUG_
1666         cerr << "Exception in " << algo->GetName() << "::GetProgress()" << endl;
1667 #endif
1668       }
1669       if ( 0. < rate && rate < 1.001 )
1670       {
1671         computedCost += rate * ( algoDoneCost + algoNotDoneCost );
1672       }
1673       else
1674       {
1675         rate = algo->GetProgressByTic();
1676         computedCost += algoDoneCost + rate * algoNotDoneCost;
1677       }
1678       // cout << "rate: "<<rate << " algoNotDoneCost: " << algoNotDoneCost << endl;
1679     }
1680
1681   // get cost of already treated sub-meshes
1682   if ( SMESH_subMesh* mainSM = GetSubMeshContaining( 1 ))
1683   {
1684     SMESH_subMeshIteratorPtr smIt = mainSM->getDependsOnIterator(/*includeSelf=*/true);
1685     while ( smIt->more() )
1686     {
1687       const SMESH_subMesh* sm = smIt->next();
1688       const int smCost = sm->GetComputeCost();
1689       totalCost += smCost;
1690       if ( !currentSubIds.Contains( sm->GetId() ) )
1691       {
1692         if (( !sm->IsEmpty() ) ||
1693             ( sm->GetComputeState() == SMESH_subMesh::FAILED_TO_COMPUTE &&
1694               !sm->DependsOn( curSM ) ))
1695           computedCost += smCost;
1696       }
1697     }
1698   }
1699   // cout << "Total: " << totalCost
1700   //      << " computed: " << computedCost << " progress: " << computedCost / totalCost
1701   //      << " nbElems: " << GetMeshDS()->GetMeshInfo().NbElements() << endl;
1702   return computedCost / totalCost;
1703 }
1704
1705 //================================================================================
1706 /*!
1707  * \brief Return number of nodes in the mesh
1708  */
1709 //================================================================================
1710
1711 int SMESH_Mesh::NbNodes() const throw(SALOME_Exception)
1712 {
1713   Unexpect aCatch(SalomeException);
1714   return _myMeshDS->NbNodes();
1715 }
1716
1717 //================================================================================
1718 /*!
1719  * \brief  Return number of edges of given order in the mesh
1720  */
1721 //================================================================================
1722
1723 int SMESH_Mesh::Nb0DElements() const throw(SALOME_Exception)
1724 {
1725   Unexpect aCatch(SalomeException);
1726   return _myMeshDS->GetMeshInfo().Nb0DElements();
1727 }
1728
1729 //================================================================================
1730 /*!
1731  * \brief  Return number of edges of given order in the mesh
1732  */
1733 //================================================================================
1734
1735 int SMESH_Mesh::NbEdges(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1736 {
1737   Unexpect aCatch(SalomeException);
1738   return _myMeshDS->GetMeshInfo().NbEdges(order);
1739 }
1740
1741 //================================================================================
1742 /*!
1743  * \brief Return number of faces of given order in the mesh
1744  */
1745 //================================================================================
1746
1747 int SMESH_Mesh::NbFaces(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1748 {
1749   Unexpect aCatch(SalomeException);
1750   return _myMeshDS->GetMeshInfo().NbFaces(order);
1751 }
1752
1753 //================================================================================
1754 /*!
1755  * \brief Return the number of faces in the mesh
1756  */
1757 //================================================================================
1758
1759 int SMESH_Mesh::NbTriangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1760 {
1761   Unexpect aCatch(SalomeException);
1762   return _myMeshDS->GetMeshInfo().NbTriangles(order);
1763 }
1764
1765 //================================================================================
1766 /*!
1767  * \brief Return number of biquadratic triangles in the mesh
1768  */
1769 //================================================================================
1770
1771 int SMESH_Mesh::NbBiQuadTriangles() const throw(SALOME_Exception)
1772 {
1773   Unexpect aCatch(SalomeException);
1774   return _myMeshDS->GetMeshInfo().NbBiQuadTriangles();
1775 }
1776
1777 //================================================================================
1778 /*!
1779  * \brief Return the number nodes faces in the mesh
1780  */
1781 //================================================================================
1782
1783 int SMESH_Mesh::NbQuadrangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1784 {
1785   Unexpect aCatch(SalomeException);
1786   return _myMeshDS->GetMeshInfo().NbQuadrangles(order);
1787 }
1788
1789 //================================================================================
1790 /*!
1791  * \brief Return number of biquadratic quadrangles in the mesh
1792  */
1793 //================================================================================
1794
1795 int SMESH_Mesh::NbBiQuadQuadrangles() const throw(SALOME_Exception)
1796 {
1797   Unexpect aCatch(SalomeException);
1798   return _myMeshDS->GetMeshInfo().NbBiQuadQuadrangles();
1799 }
1800
1801 //================================================================================
1802 /*!
1803  * \brief Return the number of polygonal faces in the mesh
1804  */
1805 //================================================================================
1806
1807 int SMESH_Mesh::NbPolygons(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1808 {
1809   Unexpect aCatch(SalomeException);
1810   return _myMeshDS->GetMeshInfo().NbPolygons(order);
1811 }
1812
1813 //================================================================================
1814 /*!
1815  * \brief Return number of volumes of given order in the mesh
1816  */
1817 //================================================================================
1818
1819 int SMESH_Mesh::NbVolumes(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1820 {
1821   Unexpect aCatch(SalomeException);
1822   return _myMeshDS->GetMeshInfo().NbVolumes(order);
1823 }
1824
1825 //================================================================================
1826 /*!
1827  * \brief  Return number of tetrahedrons of given order in the mesh
1828  */
1829 //================================================================================
1830
1831 int SMESH_Mesh::NbTetras(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1832 {
1833   Unexpect aCatch(SalomeException);
1834   return _myMeshDS->GetMeshInfo().NbTetras(order);
1835 }
1836
1837 //================================================================================
1838 /*!
1839  * \brief  Return number of hexahedrons of given order in the mesh
1840  */
1841 //================================================================================
1842
1843 int SMESH_Mesh::NbHexas(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1844 {
1845   Unexpect aCatch(SalomeException);
1846   return _myMeshDS->GetMeshInfo().NbHexas(order);
1847 }
1848
1849 //================================================================================
1850 /*!
1851  * \brief  Return number of triquadratic hexahedrons in the mesh
1852  */
1853 //================================================================================
1854
1855 int SMESH_Mesh::NbTriQuadraticHexas() const throw(SALOME_Exception)
1856 {
1857   Unexpect aCatch(SalomeException);
1858   return _myMeshDS->GetMeshInfo().NbTriQuadHexas();
1859 }
1860
1861 //================================================================================
1862 /*!
1863  * \brief  Return number of pyramids of given order in the mesh
1864  */
1865 //================================================================================
1866
1867 int SMESH_Mesh::NbPyramids(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1868 {
1869   Unexpect aCatch(SalomeException);
1870   return _myMeshDS->GetMeshInfo().NbPyramids(order);
1871 }
1872
1873 //================================================================================
1874 /*!
1875  * \brief  Return number of prisms (penthahedrons) of given order in the mesh
1876  */
1877 //================================================================================
1878
1879 int SMESH_Mesh::NbPrisms(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1880 {
1881   Unexpect aCatch(SalomeException);
1882   return _myMeshDS->GetMeshInfo().NbPrisms(order);
1883 }
1884
1885 //================================================================================
1886 /*!
1887  * \brief  Return number of hexagonal prisms in the mesh
1888  */
1889 //================================================================================
1890
1891 int SMESH_Mesh::NbHexagonalPrisms() const throw(SALOME_Exception)
1892 {
1893   Unexpect aCatch(SalomeException);
1894   return _myMeshDS->GetMeshInfo().NbHexPrisms();
1895 }
1896
1897 //================================================================================
1898 /*!
1899  * \brief  Return number of polyhedrons in the mesh
1900  */
1901 //================================================================================
1902
1903 int SMESH_Mesh::NbPolyhedrons() const throw(SALOME_Exception)
1904 {
1905   Unexpect aCatch(SalomeException);
1906   return _myMeshDS->GetMeshInfo().NbPolyhedrons();
1907 }
1908
1909 //================================================================================
1910 /*!
1911  * \brief  Return number of ball elements in the mesh
1912  */
1913 //================================================================================
1914
1915 int SMESH_Mesh::NbBalls() const throw(SALOME_Exception)
1916 {
1917   Unexpect aCatch(SalomeException);
1918   return _myMeshDS->GetMeshInfo().NbBalls();
1919 }
1920
1921 //================================================================================
1922 /*!
1923  * \brief  Return number of submeshes in the mesh
1924  */
1925 //================================================================================
1926
1927 int SMESH_Mesh::NbSubMesh() const throw(SALOME_Exception)
1928 {
1929   Unexpect aCatch(SalomeException);
1930   return _myMeshDS->NbSubMesh();
1931 }
1932
1933 //================================================================================
1934 /*!
1935  * \brief Returns number of meshes in the Study, that is supposed to be
1936  *        equal to SMESHDS_Document::NbMeshes()
1937  */
1938 //================================================================================
1939
1940 int SMESH_Mesh::NbMeshes() const // nb meshes in the Study
1941 {
1942   return _myDocument->NbMeshes();
1943 }
1944
1945 //=======================================================================
1946 //function : IsNotConformAllowed
1947 //purpose  : check if a hypothesis alowing notconform mesh is present
1948 //=======================================================================
1949
1950 bool SMESH_Mesh::IsNotConformAllowed() const
1951 {
1952   if(MYDEBUG) MESSAGE("SMESH_Mesh::IsNotConformAllowed");
1953
1954   static SMESH_HypoFilter filter( SMESH_HypoFilter::HasName( "NotConformAllowed" ));
1955   return GetHypothesis( _myMeshDS->ShapeToMesh(), filter, false );
1956 }
1957
1958 //=======================================================================
1959 //function : IsMainShape
1960 //purpose  : 
1961 //=======================================================================
1962
1963 bool SMESH_Mesh::IsMainShape(const TopoDS_Shape& theShape) const
1964 {
1965   return theShape.IsSame(_myMeshDS->ShapeToMesh() );
1966 }
1967
1968 //=============================================================================
1969 /*!
1970  *  
1971  */
1972 //=============================================================================
1973
1974 SMESH_Group* SMESH_Mesh::AddGroup (const SMDSAbs_ElementType theType,
1975                                    const char*               theName,
1976                                    int&                      theId,
1977                                    const TopoDS_Shape&       theShape,
1978                                    const SMESH_PredicatePtr& thePredicate)
1979 {
1980   if (_mapGroup.count(_groupId))
1981     return NULL;
1982   theId = _groupId;
1983   SMESH_Group* aGroup = new SMESH_Group (theId, this, theType, theName, theShape, thePredicate);
1984   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1985   _mapGroup[_groupId++] = aGroup;
1986   return aGroup;
1987 }
1988
1989 //================================================================================
1990 /*!
1991  * \brief Creates a group based on an existing SMESHDS group. Group ID should be unique
1992  */
1993 //================================================================================
1994
1995 SMESH_Group* SMESH_Mesh::AddGroup (SMESHDS_GroupBase* groupDS) throw(SALOME_Exception)
1996 {
1997   if ( !groupDS ) 
1998     throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup(): NULL SMESHDS_GroupBase"));
1999
2000   map <int, SMESH_Group*>::iterator i_g = _mapGroup.find( groupDS->GetID() );
2001   if ( i_g != _mapGroup.end() && i_g->second )
2002   {
2003     if ( i_g->second->GetGroupDS() == groupDS )
2004       return i_g->second;
2005     else
2006       throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup() wrong ID of SMESHDS_GroupBase"));
2007   }
2008   SMESH_Group* aGroup = new SMESH_Group (groupDS);
2009   _mapGroup[ groupDS->GetID() ] = aGroup;
2010   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
2011
2012   _groupId = 1 + _mapGroup.rbegin()->first;
2013
2014   return aGroup;
2015 }
2016
2017
2018 //================================================================================
2019 /*!
2020  * \brief Creates SMESH_Groups for not wrapped SMESHDS_Groups
2021  *  \retval bool - true if new SMESH_Groups have been created
2022  * 
2023  */
2024 //================================================================================
2025
2026 bool SMESH_Mesh::SynchronizeGroups()
2027 {
2028   size_t nbGroups = _mapGroup.size();
2029   const set<SMESHDS_GroupBase*>& groups = _myMeshDS->GetGroups();
2030   set<SMESHDS_GroupBase*>::const_iterator gIt = groups.begin();
2031   for ( ; gIt != groups.end(); ++gIt )
2032   {
2033     SMESHDS_GroupBase* groupDS = (SMESHDS_GroupBase*) *gIt;
2034     _groupId = groupDS->GetID();
2035     if ( !_mapGroup.count( _groupId ))
2036       _mapGroup[_groupId] = new SMESH_Group( groupDS );
2037   }
2038   if ( !_mapGroup.empty() )
2039     _groupId = _mapGroup.rbegin()->first + 1;
2040
2041   return nbGroups < _mapGroup.size();
2042 }
2043
2044 //================================================================================
2045 /*!
2046  * \brief Return iterator on all existing groups
2047  */
2048 //================================================================================
2049
2050 SMESH_Mesh::GroupIteratorPtr SMESH_Mesh::GetGroups() const
2051 {
2052   typedef map <int, SMESH_Group *> TMap;
2053   return GroupIteratorPtr( new SMDS_mapIterator<TMap>( _mapGroup ));
2054 }
2055
2056 //=============================================================================
2057 /*!
2058  * \brief Return a group by ID
2059  */
2060 //=============================================================================
2061
2062 SMESH_Group* SMESH_Mesh::GetGroup (const int theGroupID)
2063 {
2064   if (_mapGroup.find(theGroupID) == _mapGroup.end())
2065     return NULL;
2066   return _mapGroup[theGroupID];
2067 }
2068
2069
2070 //=============================================================================
2071 /*!
2072  * \brief Return IDs of all groups
2073  */
2074 //=============================================================================
2075
2076 list<int> SMESH_Mesh::GetGroupIds() const
2077 {
2078   list<int> anIds;
2079   for ( map<int, SMESH_Group*>::const_iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
2080     anIds.push_back( it->first );
2081   
2082   return anIds;
2083 }
2084
2085 //================================================================================
2086 /*!
2087  * \brief Set a caller of methods at level of CORBA API implementation.
2088  * The set upCaller will be deleted by SMESH_Mesh
2089  */
2090 //================================================================================
2091
2092 void SMESH_Mesh::SetCallUp( TCallUp* upCaller )
2093 {
2094   if ( _callUp ) delete _callUp;
2095   _callUp = upCaller;
2096 }
2097
2098 //=============================================================================
2099 /*!
2100  *  
2101  */
2102 //=============================================================================
2103
2104 bool SMESH_Mesh::RemoveGroup (const int theGroupID)
2105 {
2106   if (_mapGroup.find(theGroupID) == _mapGroup.end())
2107     return false;
2108   GetMeshDS()->RemoveGroup( _mapGroup[theGroupID]->GetGroupDS() );
2109   delete _mapGroup[theGroupID];
2110   _mapGroup.erase (theGroupID);
2111   if (_callUp)
2112     _callUp->RemoveGroup( theGroupID );
2113   return true;
2114 }
2115
2116 //=======================================================================
2117 //function : GetAncestors
2118 //purpose  : return list of ancestors of theSubShape in the order
2119 //           that lower dimention shapes come first.
2120 //=======================================================================
2121
2122 const TopTools_ListOfShape& SMESH_Mesh::GetAncestors(const TopoDS_Shape& theS) const
2123 {
2124   if ( _mapAncestors.Contains( theS ) )
2125     return _mapAncestors.FindFromKey( theS );
2126
2127   static TopTools_ListOfShape emptyList;
2128   return emptyList;
2129 }
2130
2131 //=======================================================================
2132 //function : Dump
2133 //purpose  : dumps contents of mesh to stream [ debug purposes ]
2134 //=======================================================================
2135
2136 ostream& SMESH_Mesh::Dump(ostream& save)
2137 {
2138   int clause = 0;
2139   save << "========================== Dump contents of mesh ==========================" << endl << endl;
2140   save << ++clause << ") Total number of nodes:      \t" << NbNodes() << endl;
2141   save << ++clause << ") Total number of edges:      \t" << NbEdges() << endl;
2142   save << ++clause << ") Total number of faces:      \t" << NbFaces() << endl;
2143   save << ++clause << ") Total number of polygons:   \t" << NbPolygons() << endl;
2144   save << ++clause << ") Total number of volumes:    \t" << NbVolumes() << endl;
2145   save << ++clause << ") Total number of polyhedrons:\t" << NbPolyhedrons() << endl << endl;
2146   for ( int isQuadratic = 0; isQuadratic < 2; ++isQuadratic )
2147   {
2148     string orderStr = isQuadratic ? "quadratic" : "linear";
2149     SMDSAbs_ElementOrder order  = isQuadratic ? ORDER_QUADRATIC : ORDER_LINEAR;
2150
2151     save << ++clause << ") Total number of " << orderStr << " edges:\t" << NbEdges(order) << endl;
2152     save << ++clause << ") Total number of " << orderStr << " faces:\t" << NbFaces(order) << endl;
2153     if ( NbFaces(order) > 0 ) {
2154       int nb3 = NbTriangles(order);
2155       int nb4 = NbQuadrangles(order);
2156       save << clause << ".1) Number of " << orderStr << " triangles:  \t" << nb3 << endl;
2157       save << clause << ".2) Number of " << orderStr << " quadrangles:\t" << nb4 << endl;
2158       if ( nb3 + nb4 !=  NbFaces(order) ) {
2159         map<int,int> myFaceMap;
2160         SMDS_FaceIteratorPtr itFaces=_myMeshDS->facesIterator();
2161         while( itFaces->more( ) ) {
2162           int nbNodes = itFaces->next()->NbNodes();
2163           if ( myFaceMap.find( nbNodes ) == myFaceMap.end() )
2164             myFaceMap[ nbNodes ] = 0;
2165           myFaceMap[ nbNodes ] = myFaceMap[ nbNodes ] + 1;
2166         }
2167         save << clause << ".3) Faces in detail: " << endl;
2168         map <int,int>::iterator itF;
2169         for (itF = myFaceMap.begin(); itF != myFaceMap.end(); itF++)
2170           save << "--> nb nodes: " << itF->first << " - nb elemens:\t" << itF->second << endl;
2171       }
2172     }
2173     save << ++clause << ") Total number of " << orderStr << " volumes:\t" << NbVolumes(order) << endl;
2174     if ( NbVolumes(order) > 0 ) {
2175       int nb8 = NbHexas(order);
2176       int nb4 = NbTetras(order);
2177       int nb5 = NbPyramids(order);
2178       int nb6 = NbPrisms(order);
2179       save << clause << ".1) Number of " << orderStr << " hexahedrons: \t" << nb8 << endl;
2180       save << clause << ".2) Number of " << orderStr << " tetrahedrons:\t" << nb4 << endl;
2181       save << clause << ".3) Number of " << orderStr << " prisms:      \t" << nb6 << endl;
2182       save << clause << ".4) Number of " << orderStr << " pyramids:    \t" << nb5 << endl;
2183       if ( nb8 + nb4 + nb5 + nb6 != NbVolumes(order) ) {
2184         map<int,int> myVolumesMap;
2185         SMDS_VolumeIteratorPtr itVolumes=_myMeshDS->volumesIterator();
2186         while( itVolumes->more( ) ) {
2187           int nbNodes = itVolumes->next()->NbNodes();
2188           if ( myVolumesMap.find( nbNodes ) == myVolumesMap.end() )
2189             myVolumesMap[ nbNodes ] = 0;
2190           myVolumesMap[ nbNodes ] = myVolumesMap[ nbNodes ] + 1;
2191         }
2192         save << clause << ".5) Volumes in detail: " << endl;
2193         map <int,int>::iterator itV;
2194         for (itV = myVolumesMap.begin(); itV != myVolumesMap.end(); itV++)
2195           save << "--> nb nodes: " << itV->first << " - nb elemens:\t" << itV->second << endl;
2196       }
2197     }
2198     save << endl;
2199   }
2200   save << "===========================================================================" << endl;
2201   return save;
2202 }
2203
2204 //=======================================================================
2205 //function : GetElementType
2206 //purpose  : Returns type of mesh element with certain id
2207 //=======================================================================
2208
2209 SMDSAbs_ElementType SMESH_Mesh::GetElementType( const int id, const bool iselem )
2210 {
2211   return _myMeshDS->GetElementType( id, iselem );
2212 }
2213
2214 //=============================================================================
2215 /*!
2216  *  \brief Convert group on geometry into standalone group
2217  */
2218 //=============================================================================
2219
2220 SMESH_Group* SMESH_Mesh::ConvertToStandalone ( int theGroupID )
2221 {
2222   SMESH_Group* aGroup = 0;
2223   map < int, SMESH_Group * >::iterator itg = _mapGroup.find( theGroupID );
2224   if ( itg == _mapGroup.end() )
2225     return aGroup;
2226
2227   SMESH_Group* anOldGrp = (*itg).second;
2228   if ( !anOldGrp || !anOldGrp->GetGroupDS() )
2229     return aGroup;
2230   SMESHDS_GroupBase* anOldGrpDS = anOldGrp->GetGroupDS();
2231
2232   // create new standalone group
2233   aGroup = new SMESH_Group (theGroupID, this, anOldGrpDS->GetType(), anOldGrp->GetName() );
2234   _mapGroup[theGroupID] = aGroup;
2235
2236   SMESHDS_Group* aNewGrpDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
2237   GetMeshDS()->RemoveGroup( anOldGrpDS );
2238   GetMeshDS()->AddGroup( aNewGrpDS );
2239
2240   // add elements (or nodes) into new created group
2241   SMDS_ElemIteratorPtr anItr = anOldGrpDS->GetElements();
2242   while ( anItr->more() )
2243     aNewGrpDS->Add( (anItr->next())->GetID() );
2244
2245   // set color
2246   aNewGrpDS->SetColor( anOldGrpDS->GetColor() );
2247
2248   // remove old group
2249   delete anOldGrp;
2250
2251   return aGroup;
2252 }
2253
2254 //=============================================================================
2255 /*!
2256  *  \brief remove submesh order  from Mesh
2257  */
2258 //=============================================================================
2259
2260 void SMESH_Mesh::ClearMeshOrder()
2261 {
2262   _mySubMeshOrder.clear();
2263 }
2264
2265 //=============================================================================
2266 /*!
2267  *  \brief remove submesh order  from Mesh
2268  */
2269 //=============================================================================
2270
2271 void SMESH_Mesh::SetMeshOrder(const TListOfListOfInt& theOrder )
2272 {
2273   _mySubMeshOrder = theOrder;
2274 }
2275
2276 //=============================================================================
2277 /*!
2278  *  \brief return submesh order if any
2279  */
2280 //=============================================================================
2281
2282 const TListOfListOfInt& SMESH_Mesh::GetMeshOrder() const
2283 {
2284   return _mySubMeshOrder;
2285 }
2286
2287 //=============================================================================
2288 /*!
2289  *  \brief fill _mapAncestors
2290  */
2291 //=============================================================================
2292
2293 void SMESH_Mesh::fillAncestorsMap(const TopoDS_Shape& theShape)
2294 {
2295   int desType, ancType;
2296   if ( !theShape.IsSame( GetShapeToMesh()) && theShape.ShapeType() == TopAbs_COMPOUND )
2297   {
2298     // a geom group is added. Insert it into lists of ancestors before
2299     // the first ancestor more complex than group members
2300     TopoDS_Iterator subIt( theShape );
2301     if ( !subIt.More() ) return;
2302     int memberType = subIt.Value().ShapeType();
2303     for ( desType = TopAbs_VERTEX; desType >= memberType; desType-- )
2304       for (TopExp_Explorer des( theShape, TopAbs_ShapeEnum( desType )); des.More(); des.Next())
2305       {
2306         if ( !_mapAncestors.Contains( des.Current() )) continue;// issue 0020982
2307         TopTools_ListOfShape& ancList = _mapAncestors.ChangeFromKey( des.Current() );
2308         TopTools_ListIteratorOfListOfShape ancIt (ancList);
2309         while ( ancIt.More() && ancIt.Value().ShapeType() >= memberType )
2310           ancIt.Next();
2311         if ( ancIt.More() ) ancList.InsertBefore( theShape, ancIt );
2312         else                ancList.Append( theShape );
2313       }
2314   }
2315   else // else added for 52457: Addition of hypotheses is 8 time longer than meshing
2316   {
2317     for ( desType = TopAbs_VERTEX; desType > TopAbs_COMPOUND; desType-- )
2318       for ( ancType = desType - 1; ancType >= TopAbs_COMPOUND; ancType-- )
2319         TopExp::MapShapesAndAncestors ( theShape,
2320                                         (TopAbs_ShapeEnum) desType,
2321                                         (TopAbs_ShapeEnum) ancType,
2322                                         _mapAncestors );
2323   }
2324   // visit COMPOUNDs inside a COMPOUND that are not reachable by TopExp_Explorer
2325   if ( theShape.ShapeType() == TopAbs_COMPOUND )
2326   {
2327     TopoDS_Iterator sIt(theShape);
2328     if ( sIt.More() && sIt.Value().ShapeType() == TopAbs_COMPOUND )
2329       for ( ; sIt.More(); sIt.Next() )
2330         if ( sIt.Value().ShapeType() == TopAbs_COMPOUND )
2331           fillAncestorsMap( sIt.Value() );
2332   }
2333 }
2334
2335 //=============================================================================
2336 /*!
2337  * \brief sort submeshes according to stored mesh order
2338  * \param theListToSort in out list to be sorted
2339  * \return FALSE if nothing sorted
2340  */
2341 //=============================================================================
2342
2343  bool SMESH_Mesh::SortByMeshOrder(std::vector<SMESH_subMesh*>& theListToSort) const
2344 {
2345   if ( !_mySubMeshOrder.size() || theListToSort.size() < 2)
2346     return true;
2347   
2348   bool res = false;
2349   vector<SMESH_subMesh*> onlyOrderedList;
2350   // collect all ordered submeshes in one list as pointers
2351   // and get their positions within theListToSort
2352   typedef vector<SMESH_subMesh*>::iterator TPosInList;
2353   map< int, TPosInList > sortedPos;
2354   TPosInList smBeg = theListToSort.begin(), smEnd = theListToSort.end();
2355   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2356   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2357   {
2358     const TListOfInt& listOfId = *listIdsIt;
2359     // convert sm ids to sm's
2360     vector<SMESH_subMesh*> smVec;
2361     TListOfInt::const_iterator idIt = listOfId.begin();
2362     for ( ; idIt != listOfId.end(); idIt++ ) {
2363       if ( SMESH_subMesh * sm = GetSubMeshContaining( *idIt )) {
2364         if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->IsComplexSubmesh() )
2365         {
2366           SMESHDS_SubMeshIteratorPtr smdsIt = sm->GetSubMeshDS()->GetSubMeshIterator();
2367           while ( smdsIt->more() )
2368           {
2369             const SMESHDS_SubMesh* smDS = smdsIt->next();
2370             if (( sm = GetSubMeshContaining( smDS->GetID() )))
2371               smVec.push_back( sm );
2372           }
2373         }
2374         else
2375         {
2376           smVec.push_back( sm );
2377         }
2378       }
2379     }
2380     // find smVec items in theListToSort
2381     for ( size_t i = 0; i < smVec.size(); ++i )
2382     {
2383       TPosInList smPos = find( smBeg, smEnd, smVec[i] );
2384       if ( smPos != smEnd ) {
2385         onlyOrderedList.push_back( smVec[i] );
2386         sortedPos[ distance( smBeg, smPos )] = smPos;
2387       }
2388     }
2389   }
2390   if (onlyOrderedList.size() < 2)
2391     return res;
2392   res = true;
2393
2394   vector<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
2395   vector<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
2396
2397   // iterate on ordered submeshes and insert them in detected positions
2398   map< int, TPosInList >::iterator i_pos = sortedPos.begin();
2399   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
2400     *(i_pos->second) = *onlyBIt;
2401
2402   return res;
2403 }
2404
2405 //================================================================================
2406 /*!
2407  * \brief Return true if given order of sub-meshes is OK
2408  */
2409 //================================================================================
2410
2411 bool SMESH_Mesh::IsOrderOK( const SMESH_subMesh* smBefore,
2412                             const SMESH_subMesh* smAfter ) const
2413 {
2414   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2415   TListOfInt::const_iterator idBef, idAft;
2416   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2417   {
2418     const TListOfInt& listOfId = *listIdsIt;
2419     idBef = std::find( listOfId.begin(), listOfId.end(), smBefore->GetId() );
2420     if ( idBef != listOfId.end() )
2421       idAft = std::find( listOfId.begin(), listOfId.end(), smAfter->GetId() );
2422     if ( idAft != listOfId.end () )
2423       return ( std::distance( listOfId.begin(), idBef ) <
2424                std::distance( listOfId.begin(), idAft )   );
2425   }
2426   return true; // no order imposed to given submeshes
2427
2428
2429 //=============================================================================
2430 /*!
2431  * \brief sort submeshes according to stored mesh order
2432  * \param theListToSort in out list to be sorted
2433  * \return FALSE if nothing sorted
2434  */
2435 //=============================================================================
2436
2437 void SMESH_Mesh::getAncestorsSubMeshes (const TopoDS_Shape&            theSubShape,
2438                                         std::vector< SMESH_subMesh* >& theSubMeshes) const
2439 {
2440   theSubMeshes.clear();
2441   TopTools_ListIteratorOfListOfShape it( GetAncestors( theSubShape ));
2442   for (; it.More(); it.Next() )
2443     if ( SMESH_subMesh* sm = GetSubMeshContaining( it.Value() ))
2444       theSubMeshes.push_back(sm);
2445
2446   // sort submeshes according to stored mesh order
2447   SortByMeshOrder( theSubMeshes );
2448 }