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