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