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