Salome HOME
31efebf587e13dd4e861969349f617f6a688b7ef
[modules/smesh.git] / src / SMESH / SMESH_Mesh.cxx
1 // Copyright (C) 2007-2019  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   std::list<TNameAndType>::iterator name_type = aGroupNames.begin();
525   for ( ; name_type != aGroupNames.end(); name_type++ )
526   {
527     SMESH_Group* aGroup = AddGroup( name_type->second, name_type->first.c_str() );
528     if ( aGroup ) {
529       SMESHDS_Group* aGroupDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
530       if ( aGroupDS ) {
531         aGroupDS->SetStoreName( name_type->first.c_str() );
532         myReader.GetGroup( aGroupDS );
533       }
534     }
535   }
536
537   _myMeshDS->Modified();
538   _myMeshDS->CompactMesh();
539
540   return (int) status;
541 }
542
543 //=======================================================================
544 //function : STLToMesh
545 //purpose  :
546 //=======================================================================
547
548 std::string SMESH_Mesh::STLToMesh(const char* theFileName)
549 {
550   if(_isShapeToMesh)
551     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
552   _isShapeToMesh = false;
553
554   DriverSTL_R_SMDS_Mesh myReader;
555   myReader.SetMesh(_myMeshDS);
556   myReader.SetFile(theFileName);
557   myReader.SetMeshId(-1);
558   myReader.Perform();
559
560   return myReader.GetName();
561 }
562
563 //================================================================================
564 /*!
565  * \brief Reads the given mesh from the CGNS file
566  *  \param theFileName - name of the file
567  *  \retval int - Driver_Mesh::Status
568  */
569 //================================================================================
570
571 int SMESH_Mesh::CGNSToMesh(const char*  theFileName,
572                            const int    theMeshIndex,
573                            std::string& theMeshName)
574 {
575   int res = Driver_Mesh::DRS_FAIL;
576 #ifdef WITH_CGNS
577
578   DriverCGNS_Read myReader;
579   myReader.SetMesh(_myMeshDS);
580   myReader.SetFile(theFileName);
581   myReader.SetMeshId(theMeshIndex);
582   res = myReader.Perform();
583   theMeshName = myReader.GetMeshName();
584
585   // create groups
586   SynchronizeGroups();
587
588 #endif
589   return res;
590 }
591
592 //================================================================================
593 /*!
594  * \brief Fill its data by reading a GMF file
595  */
596 //================================================================================
597
598 SMESH_ComputeErrorPtr SMESH_Mesh::GMFToMesh(const char* theFileName,
599                                             bool        theMakeRequiredGroups)
600 {
601   DriverGMF_Read myReader;
602   myReader.SetMesh(_myMeshDS);
603   myReader.SetFile(theFileName);
604   myReader.SetMakeRequiredGroups( theMakeRequiredGroups );
605   myReader.Perform();
606   //theMeshName = myReader.GetMeshName();
607
608   // create groups
609   SynchronizeGroups();
610
611   return myReader.GetError();
612 }
613
614 //=============================================================================
615 /*!
616  * 
617  */
618 //=============================================================================
619
620 SMESH_Hypothesis::Hypothesis_Status
621 SMESH_Mesh::AddHypothesis(const TopoDS_Shape & aSubShape,
622                           int                  anHypId,
623                           std::string*         anError  ) throw(SALOME_Exception)
624 {
625   Unexpect aCatch(SalomeException);
626   if(MYDEBUG) MESSAGE("SMESH_Mesh::AddHypothesis");
627
628   if ( anError )
629     anError->clear();
630
631   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
632   if ( !subMesh || !subMesh->GetId())
633     return SMESH_Hypothesis::HYP_BAD_SUBSHAPE;
634
635   SMESH_Hypothesis *anHyp = GetHypothesis( anHypId );
636   if ( !anHyp )
637     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
638
639   bool isGlobalHyp = IsMainShape( aSubShape );
640
641   // NotConformAllowed can be only global
642   if ( !isGlobalHyp )
643   {
644     // NOTE: this is not a correct way to check a name of hypothesis,
645     // there should be an attribute of hypothesis saying that it can/can't
646     // be global/local
647     std::string hypName = anHyp->GetName();
648     if ( hypName == "NotConformAllowed" )
649     {
650       if(MYDEBUG) MESSAGE( "Hypothesis <NotConformAllowed> can be only global" );
651       return SMESH_Hypothesis::HYP_INCOMPATIBLE;
652     }
653   }
654
655   // shape
656
657   bool                     isAlgo = ( anHyp->GetType() != SMESHDS_Hypothesis::PARAM_ALGO );
658   SMESH_subMesh::algo_event event = isAlgo ? SMESH_subMesh::ADD_ALGO : SMESH_subMesh::ADD_HYP;
659
660   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
661
662   if ( anError && SMESH_Hypothesis::IsStatusFatal(ret) && subMesh->GetComputeError() )
663     *anError = subMesh->GetComputeError()->myComment;
664
665   // sub-shapes
666   if ( !SMESH_Hypothesis::IsStatusFatal(ret) &&
667        anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is added on father
668   {
669     event = isAlgo ? SMESH_subMesh::ADD_FATHER_ALGO : SMESH_subMesh::ADD_FATHER_HYP;
670
671     SMESH_Hypothesis::Hypothesis_Status ret2 =
672       subMesh->SubMeshesAlgoStateEngine(event, anHyp, /*exitOnFatal=*/true);
673     if (ret2 > ret)
674     {
675       ret = ret2;
676       if ( SMESH_Hypothesis::IsStatusFatal( ret ))
677       {
678         if ( anError && subMesh->GetComputeError() )
679           *anError = subMesh->GetComputeError()->myComment;
680         // remove anHyp
681         event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
682         subMesh->AlgoStateEngine(event, anHyp);
683       }
684     }
685
686     // check concurrent hypotheses on ancestors
687     if (ret < SMESH_Hypothesis::HYP_CONCURRENT && !isGlobalHyp )
688     {
689       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
690       while ( smIt->more() ) {
691         SMESH_subMesh* sm = smIt->next();
692         if ( sm->IsApplicableHypothesis( anHyp )) {
693           ret2 = sm->CheckConcurrentHypothesis( anHyp->GetType() );
694           if (ret2 > ret) {
695             ret = ret2;
696             break;
697           }
698         }
699       }
700     }
701   }
702   HasModificationsToDiscard(); // to reset _isModified flag if a mesh becomes empty
703   GetMeshDS()->Modified();
704
705   if(MYDEBUG) subMesh->DumpAlgoState(true);
706   if(MYDEBUG) SCRUTE(ret);
707   return ret;
708 }
709
710 //=============================================================================
711 /*!
712  * 
713  */
714 //=============================================================================
715
716 SMESH_Hypothesis::Hypothesis_Status
717 SMESH_Mesh::RemoveHypothesis(const TopoDS_Shape & aSubShape,
718                              int                    anHypId) throw( SALOME_Exception )
719 {
720   Unexpect aCatch(SalomeException);
721   if(MYDEBUG) MESSAGE("SMESH_Mesh::RemoveHypothesis");
722
723   StudyContextStruct *sc = _gen->GetStudyContext();
724   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
725     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
726
727   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
728   if(MYDEBUG) { SCRUTE(anHyp->GetType()); }
729
730   // shape 
731
732   bool                     isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
733   SMESH_subMesh::algo_event event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
734
735   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
736
737   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
738
739   // there may appear concurrent hyps that were covered by the removed hyp
740   if (ret < SMESH_Hypothesis::HYP_CONCURRENT &&
741       subMesh->IsApplicableHypothesis( anHyp ) &&
742       subMesh->CheckConcurrentHypothesis( anHyp->GetType() ) != SMESH_Hypothesis::HYP_OK)
743     ret = SMESH_Hypothesis::HYP_CONCURRENT;
744
745   // sub-shapes
746   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
747       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is removed from father
748   {
749     event = isAlgo ? SMESH_subMesh::REMOVE_FATHER_ALGO : SMESH_subMesh::REMOVE_FATHER_HYP;
750
751     SMESH_Hypothesis::Hypothesis_Status ret2 =
752       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
753     if (ret2 > ret) // more severe
754       ret = ret2;
755
756     // check concurrent hypotheses on ancestors
757     if (ret < SMESH_Hypothesis::HYP_CONCURRENT && !IsMainShape( aSubShape ) )
758     {
759       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
760       while ( smIt->more() ) {
761         SMESH_subMesh* sm = smIt->next();
762         if ( sm->IsApplicableHypothesis( anHyp )) {
763           ret2 = sm->CheckConcurrentHypothesis( anHyp->GetType() );
764           if (ret2 > ret) {
765             ret = ret2;
766             break;
767           }
768         }
769       }
770     }
771   }
772
773   HasModificationsToDiscard(); // to reset _isModified flag if mesh become empty
774   GetMeshDS()->Modified();
775
776   if(MYDEBUG) subMesh->DumpAlgoState(true);
777   if(MYDEBUG) SCRUTE(ret);
778   return ret;
779 }
780
781 //=============================================================================
782 /*!
783  * 
784  */
785 //=============================================================================
786
787 const std::list<const SMESHDS_Hypothesis*>&
788 SMESH_Mesh::GetHypothesisList(const TopoDS_Shape & aSubShape) const
789   throw(SALOME_Exception)
790 {
791   return _myMeshDS->GetHypothesis(aSubShape);
792 }
793
794 //=======================================================================
795 /*!
796  * \brief Return the hypothesis assigned to the shape
797  *  \param aSubShape    - the shape to check
798  *  \param aFilter      - the hypothesis filter
799  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
800  *  \param assignedTo   - to return the shape the found hypo is assigned to
801  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
802  */
803 //=======================================================================
804
805 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const TopoDS_Shape &    aSubShape,
806                                                    const SMESH_HypoFilter& aFilter,
807                                                    const bool              andAncestors,
808                                                    TopoDS_Shape*           assignedTo) const
809 {
810   return GetHypothesis( const_cast< SMESH_Mesh* >(this)->GetSubMesh( aSubShape ),
811                         aFilter, andAncestors, assignedTo );
812 }
813
814 //=======================================================================
815 /*!
816  * \brief Return the hypothesis assigned to the shape of a sub-mesh
817  *  \param aSubMesh     - the sub-mesh to check
818  *  \param aFilter      - the hypothesis filter
819  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
820  *  \param assignedTo   - to return the shape the found hypo is assigned to
821  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
822  */
823 //=======================================================================
824
825 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const SMESH_subMesh *   aSubMesh,
826                                                    const SMESH_HypoFilter& aFilter,
827                                                    const bool              andAncestors,
828                                                    TopoDS_Shape*           assignedTo) const
829 {
830   if ( !aSubMesh ) return 0;
831
832   {
833     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
834     const std::list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
835     std::list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
836     for ( ; hyp != hypList.end(); hyp++ ) {
837       const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
838       if ( aFilter.IsOk( h, aSubShape)) {
839         if ( assignedTo ) *assignedTo = aSubShape;
840         return h;
841       }
842     }
843   }
844   if ( andAncestors )
845   {
846     // user sorted submeshes of ancestors, according to stored submesh priority
847     std::vector< SMESH_subMesh * > & ancestors =
848       const_cast< std::vector< SMESH_subMesh * > & > ( aSubMesh->GetAncestors() );
849     SortByMeshOrder( ancestors );
850
851     std::vector<SMESH_subMesh*>::const_iterator smIt = ancestors.begin(); 
852     for ( ; smIt != ancestors.end(); smIt++ )
853     {
854       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
855       const std::list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
856       std::list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
857       for ( ; hyp != hypList.end(); hyp++ ) {
858         const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
859         if (aFilter.IsOk( h, curSh )) {
860           if ( assignedTo ) *assignedTo = curSh;
861           return h;
862         }
863       }
864     }
865   }
866   return 0;
867 }
868
869 //================================================================================
870 /*!
871  * \brief Return hypotheses assigned to the shape
872   * \param aSubShape - the shape to check
873   * \param aFilter - the hypothesis filter
874   * \param aHypList - the list of the found hypotheses
875   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
876   * \retval int - number of unique hypos in aHypList
877  */
878 //================================================================================
879
880 int SMESH_Mesh::GetHypotheses(const TopoDS_Shape &                     aSubShape,
881                               const SMESH_HypoFilter&                  aFilter,
882                               std::list <const SMESHDS_Hypothesis * >& aHypList,
883                               const bool                               andAncestors,
884                               std::list< TopoDS_Shape > *              assignedTo/*=0*/) const
885 {
886   return GetHypotheses( const_cast< SMESH_Mesh* >(this)->GetSubMesh( aSubShape ),
887                         aFilter, aHypList, andAncestors, assignedTo );
888 }
889
890 //================================================================================
891 /*!
892  * \brief Return hypotheses assigned to the shape of a sub-mesh
893   * \param aSubShape - the sub-mesh to check
894   * \param aFilter - the hypothesis filter
895   * \param aHypList - the list of the found hypotheses
896   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
897   * \retval int - number of unique hypos in aHypList
898  */
899 //================================================================================
900
901 int SMESH_Mesh::GetHypotheses(const SMESH_subMesh *                    aSubMesh,
902                               const SMESH_HypoFilter&                  aFilter,
903                               std::list <const SMESHDS_Hypothesis * >& aHypList,
904                               const bool                               andAncestors,
905                               std::list< TopoDS_Shape > *              assignedTo/*=0*/) const
906 {
907   if ( !aSubMesh ) return 0;
908
909   std::set< std::string > hypTypes; // to exclude same type hypos from the result list
910   int nbHyps = 0;
911
912   // only one main hypothesis is allowed
913   bool mainHypFound = false;
914
915   // fill in hypTypes
916   std::list<const SMESHDS_Hypothesis*>::const_iterator hyp;
917   for ( hyp = aHypList.begin(); hyp != aHypList.end(); hyp++ ) {
918     if ( hypTypes.insert( (*hyp)->GetName() ).second )
919       nbHyps++;
920     if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
921       mainHypFound = true;
922   }
923
924   // get hypos from aSubShape
925   {
926     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
927     const std::list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
928     for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
929     {
930       const SMESH_Hypothesis* h = cSMESH_Hyp( *hyp );
931       if (( aFilter.IsOk( h, aSubShape )) &&
932           ( h->IsAuxiliary() || !mainHypFound ) &&
933           ( h->IsAuxiliary() || hypTypes.insert( h->GetName() ).second ))
934       {
935         aHypList.push_back( *hyp );
936         nbHyps++;
937         if ( !h->IsAuxiliary() )
938           mainHypFound = true;
939         if ( assignedTo ) assignedTo->push_back( aSubShape );
940       }
941     }
942   }
943
944   // get hypos from ancestors of aSubShape
945   if ( andAncestors )
946   {
947     // user sorted submeshes of ancestors, according to stored submesh priority
948     std::vector< SMESH_subMesh * > & ancestors =
949       const_cast< std::vector< SMESH_subMesh * > & > ( aSubMesh->GetAncestors() );
950     SortByMeshOrder( ancestors );
951
952     std::vector<SMESH_subMesh*>::const_iterator smIt = ancestors.begin();
953     for ( ; smIt != ancestors.end(); smIt++ )
954     {
955       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
956       const std::list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
957       for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
958       {
959         const SMESH_Hypothesis* h = cSMESH_Hyp( *hyp );
960         if (( aFilter.IsOk( h, curSh )) &&
961             ( h->IsAuxiliary() || !mainHypFound ) &&
962             ( h->IsAuxiliary() || hypTypes.insert( h->GetName() ).second ))
963         {
964           aHypList.push_back( *hyp );
965           nbHyps++;
966           if ( !h->IsAuxiliary() )
967             mainHypFound = true;
968           if ( assignedTo ) assignedTo->push_back( curSh );
969         }
970       }
971     }
972   }
973   return nbHyps;
974 }
975
976 //================================================================================
977 /*!
978  * \brief Return a hypothesis by its ID
979  */
980 //================================================================================
981
982 SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const int anHypId) const
983 {
984   StudyContextStruct *sc = _gen->GetStudyContext();
985   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
986     return NULL;
987
988   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
989   return anHyp;
990 }
991
992 //=============================================================================
993 /*!
994  * 
995  */
996 //=============================================================================
997
998 const std::list<SMESHDS_Command*> & SMESH_Mesh::GetLog() throw(SALOME_Exception)
999 {
1000   Unexpect aCatch(SalomeException);
1001   return _myMeshDS->GetScript()->GetCommands();
1002 }
1003
1004 //=============================================================================
1005 /*!
1006  * 
1007  */
1008 //=============================================================================
1009 void SMESH_Mesh::ClearLog() throw(SALOME_Exception)
1010 {
1011   Unexpect aCatch(SalomeException);
1012   _myMeshDS->GetScript()->Clear();
1013 }
1014
1015 //=============================================================================
1016 /*!
1017  * Get or Create the SMESH_subMesh object implementation
1018  */
1019 //=============================================================================
1020
1021 SMESH_subMesh * SMESH_Mesh::GetSubMesh(const TopoDS_Shape & aSubShape)
1022   throw(SALOME_Exception)
1023 {
1024   int index = _myMeshDS->ShapeToIndex(aSubShape);
1025   if ( !index && aSubShape.IsNull() )
1026     return 0;
1027
1028   // for submeshes on GEOM Group
1029   if (( !index || index > _nbSubShapes ) && aSubShape.ShapeType() == TopAbs_COMPOUND )
1030   {
1031     TopoDS_Iterator it( aSubShape );
1032     if ( it.More() )
1033     {
1034       index = _myMeshDS->AddCompoundSubmesh( aSubShape, it.Value().ShapeType() );
1035       // fill map of Ancestors
1036       while ( _nbSubShapes < index )
1037         fillAncestorsMap( _myMeshDS->IndexToShape( ++_nbSubShapes ));
1038     }
1039   }
1040   // if ( !index )
1041   //   return NULL; // neither sub-shape nor a group
1042
1043   SMESH_subMesh* aSubMesh = _subMeshHolder->Get( index );
1044   if ( !aSubMesh )
1045   {
1046     aSubMesh = new SMESH_subMesh(index, this, _myMeshDS, aSubShape);
1047     _subMeshHolder->Add( index, aSubMesh );
1048
1049     // include non-computable sub-meshes in SMESH_subMesh::_ancestors of sub-submeshes
1050     switch ( aSubShape.ShapeType() ) {
1051     case TopAbs_COMPOUND:
1052     case TopAbs_WIRE:
1053     case TopAbs_SHELL:
1054       for ( TopoDS_Iterator subIt( aSubShape ); subIt.More(); subIt.Next() )
1055       {
1056         SMESH_subMesh* sm = GetSubMesh( subIt.Value() );
1057         SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*inclideSelf=*/true);
1058         while ( smIt->more() )
1059           smIt->next()->ClearAncestors();
1060       }
1061     default:;
1062     }
1063   }
1064   return aSubMesh;
1065 }
1066
1067 //=============================================================================
1068 /*!
1069  * Get the SMESH_subMesh object implementation. Don't create it, return null
1070  * if it does not exist.
1071  */
1072 //=============================================================================
1073
1074 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const TopoDS_Shape & aSubShape) const
1075   throw(SALOME_Exception)
1076 {
1077   int index = _myMeshDS->ShapeToIndex(aSubShape);
1078   return GetSubMeshContaining( index );
1079 }
1080
1081 //=============================================================================
1082 /*!
1083  * Get the SMESH_subMesh object implementation. Don't create it, return null
1084  * if it does not exist.
1085  */
1086 //=============================================================================
1087
1088 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const int aShapeID) const
1089 throw(SALOME_Exception)
1090 {
1091   SMESH_subMesh *aSubMesh = _subMeshHolder->Get( aShapeID );
1092
1093   return aSubMesh;
1094 }
1095
1096 //================================================================================
1097 /*!
1098  * \brief Return sub-meshes of groups containing the given sub-shape
1099  */
1100 //================================================================================
1101
1102 std::list<SMESH_subMesh*>
1103 SMESH_Mesh::GetGroupSubMeshesContaining(const TopoDS_Shape & aSubShape) const
1104   throw(SALOME_Exception)
1105 {
1106   std::list<SMESH_subMesh*> found;
1107
1108   SMESH_subMesh * subMesh = GetSubMeshContaining(aSubShape);
1109   if ( !subMesh )
1110     return found;
1111
1112   // sub-meshes of groups have max IDs, so search from the map end
1113   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator( /*reverse=*/true ) );
1114   while ( smIt->more() ) {
1115     SMESH_subMesh*    sm = smIt->next();
1116     SMESHDS_SubMesh * ds = sm->GetSubMeshDS();
1117     if ( ds && ds->IsComplexSubmesh() ) {
1118       if ( SMESH_MesherHelper::IsSubShape( aSubShape, sm->GetSubShape() ))
1119       {
1120         found.push_back( sm );
1121         //break;
1122       }
1123     } else {
1124       break; // the rest sub-meshes are not those of groups
1125     }
1126   }
1127
1128   if ( found.empty() ) // maybe the main shape is a COMPOUND (issue 0021530)
1129   {
1130     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1) )
1131       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1132       {
1133         TopoDS_Iterator it( mainSM->GetSubShape() );
1134         if ( it.Value().ShapeType() == aSubShape.ShapeType() &&
1135              SMESH_MesherHelper::IsSubShape( aSubShape, mainSM->GetSubShape() ))
1136           found.push_back( mainSM );
1137       }
1138   }
1139   else // issue 0023068
1140   {
1141     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1) )
1142       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1143         found.push_back( mainSM );
1144   }
1145   return found;
1146 }
1147 //=======================================================================
1148 //function : IsUsedHypothesis
1149 //purpose  : Return True if anHyp is used to mesh aSubShape
1150 //=======================================================================
1151
1152 bool SMESH_Mesh::IsUsedHypothesis(SMESHDS_Hypothesis * anHyp,
1153                                   const SMESH_subMesh* aSubMesh)
1154 {
1155   SMESH_Hypothesis* hyp = static_cast<SMESH_Hypothesis*>(anHyp);
1156
1157   // check if anHyp can be used to mesh aSubMesh
1158   if ( !aSubMesh || !aSubMesh->IsApplicableHypothesis( hyp ))
1159     return false;
1160
1161   SMESH_Algo *algo = aSubMesh->GetAlgo();
1162
1163   // algorithm
1164   if (anHyp->GetType() > SMESHDS_Hypothesis::PARAM_ALGO)
1165     return ( anHyp == algo );
1166
1167   // algorithm parameter
1168   if (algo)
1169   {
1170     // look trough hypotheses used by algo
1171     const SMESH_HypoFilter* hypoKind;
1172     if (( hypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() ))) {
1173       std::list <const SMESHDS_Hypothesis * > usedHyps;
1174       if ( GetHypotheses( aSubMesh, *hypoKind, usedHyps, true ))
1175         return ( find( usedHyps.begin(), usedHyps.end(), anHyp ) != usedHyps.end() );
1176     }
1177   }
1178
1179   return false;
1180 }
1181
1182 //=======================================================================
1183 //function : NotifySubMeshesHypothesisModification
1184 //purpose  : Say all submeshes using theChangedHyp that it has been modified
1185 //=======================================================================
1186
1187 void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* hyp)
1188 {
1189   Unexpect aCatch(SalomeException);
1190
1191   if ( !GetMeshDS()->IsUsedHypothesis( hyp ))
1192     return;
1193
1194   if (_callUp && hyp)
1195     _callUp->HypothesisModified( hyp->GetID() );
1196
1197   SMESH_Algo *algo;
1198   const SMESH_HypoFilter* compatibleHypoKind;
1199   std::list <const SMESHDS_Hypothesis * > usedHyps;
1200   std::vector< SMESH_subMesh* > smToNotify;
1201   bool allMeshedEdgesNotified = true;
1202
1203   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1204   while ( smIt->more() )
1205   {
1206     SMESH_subMesh* aSubMesh = smIt->next();
1207     bool toNotify = false;
1208
1209     // if aSubMesh meshing depends on hyp,
1210     // we call aSubMesh->AlgoStateEngine( MODIF_HYP, hyp ) that causes either
1211     // 1) clearing already computed aSubMesh or
1212     // 2) changing algo_state from MISSING_HYP to HYP_OK when parameters of hyp becomes valid,
1213     // other possible changes are not interesting. (IPAL0052457 - assigning hyp performance pb)
1214     if ( aSubMesh->GetComputeState() == SMESH_subMesh::COMPUTE_OK ||
1215          aSubMesh->GetComputeState() == SMESH_subMesh::FAILED_TO_COMPUTE ||
1216          aSubMesh->GetAlgoState()    == SMESH_subMesh::MISSING_HYP ||
1217          hyp->DataDependOnParams() )
1218     {
1219       const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
1220
1221       if (( aSubMesh->IsApplicableHypothesis( hyp )) &&
1222           ( algo = aSubMesh->GetAlgo() )            &&
1223           ( compatibleHypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() )) &&
1224           ( compatibleHypoKind->IsOk( hyp, aSubShape )))
1225       {
1226         // check if hyp is used by algo
1227         usedHyps.clear();
1228         toNotify = ( GetHypotheses( aSubMesh, *compatibleHypoKind, usedHyps, true ) &&
1229                      std::find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() );
1230       }
1231     }
1232     if ( toNotify )
1233     {
1234       smToNotify.push_back( aSubMesh );
1235       if ( aSubMesh->GetAlgoState() == SMESH_subMesh::MISSING_HYP )
1236         allMeshedEdgesNotified = false; //  update of algo state needed, not mesh clearing
1237     }
1238     else
1239     {
1240       if ( !aSubMesh->IsEmpty() &&
1241            aSubMesh->GetSubShape().ShapeType() == TopAbs_EDGE )
1242         allMeshedEdgesNotified = false;
1243     }
1244   }
1245   if ( smToNotify.empty() )
1246     return;
1247
1248   // if all meshed EDGEs will be notified then the notification is equivalent
1249   // to the whole mesh clearing, which is usually faster
1250   if ( allMeshedEdgesNotified && NbNodes() > 0 )
1251   {
1252     Clear();
1253   }
1254   else
1255   {
1256     // notify in reverse order to avoid filling the pool of IDs
1257     for ( int i = smToNotify.size()-1; i >= 0; --i )
1258       smToNotify[i]->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
1259                                      const_cast< SMESH_Hypothesis*>( hyp ));
1260   }
1261   HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
1262   GetMeshDS()->Modified();
1263 }
1264
1265 //=============================================================================
1266 /*!
1267  *  Auto color functionality
1268  */
1269 //=============================================================================
1270 void SMESH_Mesh::SetAutoColor(bool theAutoColor) throw(SALOME_Exception)
1271 {
1272   Unexpect aCatch(SalomeException);
1273   _isAutoColor = theAutoColor;
1274 }
1275
1276 bool SMESH_Mesh::GetAutoColor() throw(SALOME_Exception)
1277 {
1278   Unexpect aCatch(SalomeException);
1279   return _isAutoColor;
1280 }
1281
1282 //=======================================================================
1283 //function : SetIsModified
1284 //purpose  : Set the flag meaning that the mesh has been edited "manually"
1285 //=======================================================================
1286
1287 void SMESH_Mesh::SetIsModified(bool isModified)
1288 {
1289   _isModified = isModified;
1290
1291   if ( _isModified )
1292     // check if mesh becomes empty as result of modification
1293     HasModificationsToDiscard();
1294 }
1295
1296 //=======================================================================
1297 //function : HasModificationsToDiscard
1298 //purpose  : Return true if the mesh has been edited since a total re-compute
1299 //           and those modifications may prevent successful partial re-compute.
1300 //           As a side effect reset _isModified flag if mesh is empty
1301 //issue    : 0020693
1302 //=======================================================================
1303
1304 bool SMESH_Mesh::HasModificationsToDiscard() const
1305 {
1306   if ( ! _isModified )
1307     return false;
1308
1309   // return true if the next Compute() will be partial and
1310   // existing but changed elements may prevent successful re-compute
1311   bool hasComputed = false, hasNotComputed = false;
1312   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1313   while ( smIt->more() )
1314   {
1315     const SMESH_subMesh* aSubMesh = smIt->next();
1316     switch ( aSubMesh->GetSubShape().ShapeType() )
1317     {
1318     case TopAbs_EDGE:
1319     case TopAbs_FACE:
1320     case TopAbs_SOLID:
1321       if ( aSubMesh->IsMeshComputed() )
1322         hasComputed = true;
1323       else
1324         hasNotComputed = true;
1325       if ( hasComputed && hasNotComputed)
1326         return true;
1327
1328     default:;
1329     }
1330   }
1331   if ( NbNodes() < 1 )
1332     const_cast<SMESH_Mesh*>(this)->_isModified = false;
1333
1334   return false;
1335 }
1336
1337 //=============================================================================
1338 /*!
1339  * \brief Return true if all sub-meshes are computed OK - to update an icon
1340  */
1341 //=============================================================================
1342
1343 bool SMESH_Mesh::IsComputedOK()
1344 {
1345   if ( NbNodes() == 0 )
1346     return false;
1347
1348   // if ( !HasShapeToMesh() )
1349   //   return true;
1350
1351   if ( SMESH_subMesh* mainSM = GetSubMeshContaining( 1 ))
1352   {
1353     SMESH_subMeshIteratorPtr smIt = mainSM->getDependsOnIterator(/*includeSelf=*/true);
1354     while ( smIt->more() )
1355     {
1356       const SMESH_subMesh* sm = smIt->next();
1357       if ( !sm->IsAlwaysComputed() )
1358         switch ( sm->GetComputeState() )
1359         {
1360         case SMESH_subMesh::NOT_READY:
1361         case SMESH_subMesh::COMPUTE_OK:
1362           continue; // ok
1363         case SMESH_subMesh::FAILED_TO_COMPUTE:
1364         case SMESH_subMesh::READY_TO_COMPUTE:
1365           return false;
1366         }
1367     }
1368   }
1369   return true;
1370 }
1371
1372 //================================================================================
1373 /*!
1374  * \brief Check if any groups of the same type have equal names
1375  */
1376 //================================================================================
1377
1378 bool SMESH_Mesh::HasDuplicatedGroupNamesMED()
1379 {
1380   // Corrected for Mantis issue 0020028
1381   std::map< SMDSAbs_ElementType, std::set< std::string > > aGroupNames;
1382   for ( std::map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1383   {
1384     SMESH_Group*       aGroup = it->second;
1385     SMDSAbs_ElementType aType = aGroup->GetGroupDS()->GetType();
1386     std::string    aGroupName = aGroup->GetName();
1387     aGroupName.resize( MAX_MED_GROUP_NAME_LENGTH );
1388     if ( !aGroupNames[aType].insert(aGroupName).second )
1389       return true;
1390   }
1391
1392   return false;
1393 }
1394
1395 //================================================================================
1396 /*!
1397  * \brief Export the mesh to a med file
1398  *  \param [in] file - name of the MED file
1399  *  \param [in] theMeshName - name of this mesh
1400  *  \param [in] theAutoGroups - boolean parameter for creating/not creating
1401  *              the groups Group_On_All_Nodes, Group_On_All_Faces, ... ;
1402  *              the typical use is auto_groups=false.
1403  *  \param [in] theMinor - define the minor version (y, where version is x.y.z) of MED file format.
1404  *              The theMinor must be between 0 and the current minor version of MED file library.
1405  *              If theMinor is equal to -1, the minor version is not changed (default).
1406  *              The major version (x, where version is x.y.z) cannot be changed.
1407  *  \param [in] meshPart - mesh data to export
1408  *  \param [in] theAutoDimension - if \c true, a space dimension of a MED mesh can be either
1409  *              - 1D if all mesh nodes lie on OX coordinate axis, or
1410  *              - 2D if all mesh nodes lie on XOY coordinate plane, or
1411  *              - 3D in the rest cases.
1412  *              If \a theAutoDimension is \c false, the space dimension is always 3.
1413  *  \param [in] theAddODOnVertices - to create 0D elements on all vertices
1414  *  \param [in] theAllElemsToGroup - to make every element to belong to any group (PAL23413)
1415  *  \param [in] ZTolerance - tolerance in Z direction. If Z coordinate of a node is close to zero
1416  *              within a given tolerance, the coordinate is set to zero.
1417  *              If \a ZTolerance is negative, the node coordinates are kept as is.
1418  *  \return int - mesh index in the file
1419  */
1420 //================================================================================
1421
1422 void SMESH_Mesh::ExportMED(const char *        file,
1423                            const char*         theMeshName,
1424                            bool                theAutoGroups,
1425                            int                 theVersion,
1426                            const SMESHDS_Mesh* meshPart,
1427                            bool                theAutoDimension,
1428                            bool                theAddODOnVertices,
1429                            double              theZTolerance,
1430                            bool                theAllElemsToGroup)
1431 throw(SALOME_Exception)
1432 {
1433   MESSAGE("MED_VERSION:"<< theVersion);
1434   SMESH_TRY;
1435
1436   DriverMED_W_SMESHDS_Mesh myWriter;
1437   myWriter.SetFile         ( file , theVersion);
1438   myWriter.SetMesh         ( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS   );
1439   myWriter.SetAutoDimension( theAutoDimension );
1440   myWriter.AddODOnVertices ( theAddODOnVertices );
1441   myWriter.SetZTolerance   ( theZTolerance );
1442   if ( !theMeshName )
1443     myWriter.SetMeshId     ( _id         );
1444   else {
1445     myWriter.SetMeshId     ( -1          );
1446     myWriter.SetMeshName   ( theMeshName );
1447   }
1448
1449   if ( theAutoGroups ) {
1450     myWriter.AddGroupOfNodes();
1451     myWriter.AddGroupOfEdges();
1452     myWriter.AddGroupOfFaces();
1453     myWriter.AddGroupOfVolumes();
1454     myWriter.AddGroupOf0DElems();
1455     myWriter.AddGroupOfBalls();
1456   }
1457   if ( theAllElemsToGroup )
1458     myWriter.AddAllToGroup();
1459
1460   // Pass groups to writer. Provide unique group names.
1461   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1462   if ( !meshPart )
1463   {
1464     std::map< SMDSAbs_ElementType, std::set<std::string> > aGroupNames;
1465     char aString [256];
1466     int maxNbIter = 10000; // to guarantee cycle finish
1467     for ( std::map<int, SMESH_Group*>::iterator it = _mapGroup.begin();
1468           it != _mapGroup.end();
1469           it++ ) {
1470       SMESH_Group*       aGroup   = it->second;
1471       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1472       if ( aGroupDS ) {
1473         SMDSAbs_ElementType aType = aGroupDS->GetType();
1474         std::string aGroupName0 = aGroup->GetName();
1475         aGroupName0.resize(MAX_MED_GROUP_NAME_LENGTH);
1476         std::string aGroupName = aGroupName0;
1477         for (int i = 1; !aGroupNames[aType].insert(aGroupName).second && i < maxNbIter; i++) {
1478           sprintf(&aString[0], "GR_%d_%s", i, aGroupName0.c_str());
1479           aGroupName = aString;
1480           aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1481         }
1482         aGroupDS->SetStoreName( aGroupName.c_str() );
1483         myWriter.AddGroup( aGroupDS );
1484       }
1485     }
1486   }
1487   // Perform export
1488   myWriter.Perform();
1489
1490   SMESH_CATCH( SMESH::throwSalomeEx );
1491 }
1492
1493 //================================================================================
1494 /*!
1495  * \brief Export the mesh to a SAUV file
1496  */
1497 //================================================================================
1498
1499 void SMESH_Mesh::ExportSAUV(const char *file, 
1500                             const char* theMeshName, 
1501                             bool theAutoGroups)
1502   throw(SALOME_Exception)
1503 {
1504   std::string medfilename(file);
1505   medfilename += ".med";
1506   std::string cmd;
1507 #ifdef WIN32
1508   cmd = "%PYTHONBIN% ";
1509 #else
1510   cmd = "python3 ";
1511 #endif
1512   cmd += "-c \"";
1513   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1514   cmd += "\"";
1515   system(cmd.c_str());
1516   ExportMED(medfilename.c_str(), theMeshName, theAutoGroups, /*minor=*/-1,
1517             /*meshPart=*/NULL, /*theAutoDimension=*/false, /*theAddODOnVertices=*/false,
1518             /*zTol=*/-1, /*theAllElemsToGroup=*/true ); // theAllElemsToGroup is for PAL0023413
1519 #ifdef WIN32
1520   cmd = "%PYTHONBIN% ";
1521 #else
1522   cmd = "python3 ";
1523 #endif
1524   cmd += "-c \"";
1525   cmd += "from medutilities import convert ; convert(r'" + medfilename + "', 'MED', 'GIBI', 1, r'" + file + "')";
1526   cmd += "\"";
1527   system(cmd.c_str());
1528 #ifdef WIN32
1529   cmd = "%PYTHONBIN% ";
1530 #else
1531   cmd = "python3 ";
1532 #endif
1533   cmd += "-c \"";
1534   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1535   cmd += "\"";
1536   system(cmd.c_str());
1537 }
1538
1539 //================================================================================
1540 /*!
1541  * \brief Export the mesh to a DAT file
1542  */
1543 //================================================================================
1544
1545 void SMESH_Mesh::ExportDAT(const char *        file,
1546                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1547 {
1548   Unexpect aCatch(SalomeException);
1549   DriverDAT_W_SMDS_Mesh myWriter;
1550   myWriter.SetFile( file );
1551   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1552   myWriter.SetMeshId(_id);
1553   myWriter.Perform();
1554 }
1555
1556 //================================================================================
1557 /*!
1558  * \brief Export the mesh to an UNV file
1559  */
1560 //================================================================================
1561
1562 void SMESH_Mesh::ExportUNV(const char *        file,
1563                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1564 {
1565   Unexpect aCatch(SalomeException);
1566   DriverUNV_W_SMDS_Mesh myWriter;
1567   myWriter.SetFile( file );
1568   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1569   myWriter.SetMeshId(_id);
1570   //  myWriter.SetGroups(_mapGroup);
1571
1572   // pass group names to SMESHDS
1573   if ( !meshPart )
1574   {
1575     std::map<int, SMESH_Group*>::iterator it = _mapGroup.begin();
1576     for ( ; it != _mapGroup.end(); it++ ) {
1577       SMESH_Group*       aGroup   = it->second;
1578       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1579       if ( aGroupDS ) {
1580         std::string aGroupName = aGroup->GetName();
1581         aGroupDS->SetStoreName( aGroupName.c_str() );
1582         myWriter.AddGroup( aGroupDS );
1583       }
1584     }
1585   }
1586   myWriter.Perform();
1587 }
1588
1589 //================================================================================
1590 /*!
1591  * \brief Export the mesh to an STL file
1592  */
1593 //================================================================================
1594
1595 void SMESH_Mesh::ExportSTL(const char *        file,
1596                            const bool          isascii,
1597                            const char *        name,
1598                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1599 {
1600   Unexpect aCatch(SalomeException);
1601   DriverSTL_W_SMDS_Mesh myWriter;
1602   myWriter.SetFile( file );
1603   myWriter.SetIsAscii( isascii );
1604   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS);
1605   myWriter.SetMeshId(_id);
1606   if ( name ) myWriter.SetName( name );
1607   myWriter.Perform();
1608 }
1609
1610 //================================================================================
1611 /*!
1612  * \brief Export the mesh to the CGNS file
1613  */
1614 //================================================================================
1615
1616 void SMESH_Mesh::ExportCGNS(const char *        file,
1617                             const SMESHDS_Mesh* meshDS,
1618                             const char *        meshName,
1619                             const bool          groupElemsByType)
1620 {
1621   int res = Driver_Mesh::DRS_FAIL;
1622
1623   // pass group names to SMESHDS
1624   std::map<int, SMESH_Group*>::iterator it = _mapGroup.begin();
1625   for ( ; it != _mapGroup.end(); it++ ) {
1626     SMESH_Group*       group   = it->second;
1627     SMESHDS_GroupBase* groupDS = group->GetGroupDS();
1628     if ( groupDS ) {
1629       std::string groupName = group->GetName();
1630       groupDS->SetStoreName( groupName.c_str() );
1631     }
1632   }
1633 #ifdef WITH_CGNS
1634
1635   DriverCGNS_Write myWriter;
1636   myWriter.SetFile( file );
1637   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1638   myWriter.SetMeshName( SMESH_Comment("Mesh_") << meshDS->GetPersistentId());
1639   if ( meshName && meshName[0] )
1640     myWriter.SetMeshName( meshName );
1641   myWriter.SetElementsByType( groupElemsByType );
1642   res = myWriter.Perform();
1643   if ( res != Driver_Mesh::DRS_OK )
1644   {
1645     SMESH_ComputeErrorPtr err = myWriter.GetError();
1646     if ( err && !err->IsOK() && !err->myComment.empty() )
1647       throw SALOME_Exception(("Export failed: " + err->myComment ).c_str() );
1648   }
1649
1650 #endif
1651   if ( res != Driver_Mesh::DRS_OK )
1652     throw SALOME_Exception("Export failed");
1653 }
1654
1655 //================================================================================
1656 /*!
1657  * \brief Export the mesh to a GMF file
1658  */
1659 //================================================================================
1660
1661 void SMESH_Mesh::ExportGMF(const char *        file,
1662                            const SMESHDS_Mesh* meshDS,
1663                            bool                withRequiredGroups)
1664 {
1665   DriverGMF_Write myWriter;
1666   myWriter.SetFile( file );
1667   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1668   myWriter.SetExportRequiredGroups( withRequiredGroups );
1669
1670   myWriter.Perform();
1671 }
1672
1673 //================================================================================
1674 /*!
1675  * \brief Return a ratio of "compute cost" of computed sub-meshes to the whole
1676  *        "compute cost".
1677  */
1678 //================================================================================
1679
1680 double SMESH_Mesh::GetComputeProgress() const
1681 {
1682   double totalCost = 1e-100, computedCost = 0;
1683   const SMESH_subMesh* curSM = _gen->GetCurrentSubMesh();
1684
1685   // get progress of a current algo
1686   TColStd_MapOfInteger currentSubIds; 
1687   if ( curSM )
1688     if ( SMESH_Algo* algo = curSM->GetAlgo() )
1689     {
1690       int algoNotDoneCost = 0, algoDoneCost = 0;
1691       const std::vector<SMESH_subMesh*>& smToCompute = algo->SubMeshesToCompute();
1692       for ( size_t i = 0; i < smToCompute.size(); ++i )
1693       {
1694         if ( smToCompute[i]->IsEmpty() || smToCompute.size() == 1 )
1695           algoNotDoneCost += smToCompute[i]->GetComputeCost();
1696         else
1697           algoDoneCost += smToCompute[i]->GetComputeCost();
1698         currentSubIds.Add( smToCompute[i]->GetId() );
1699       }
1700       double rate = 0;
1701       try
1702       {
1703         OCC_CATCH_SIGNALS;
1704         rate = algo->GetProgress();
1705       }
1706       catch (...) {
1707 #ifdef _DEBUG_
1708         std::cerr << "Exception in " << algo->GetName() << "::GetProgress()" << std::endl;
1709 #endif
1710       }
1711       if ( 0. < rate && rate < 1.001 )
1712       {
1713         computedCost += rate * ( algoDoneCost + algoNotDoneCost );
1714       }
1715       else
1716       {
1717         rate = algo->GetProgressByTic();
1718         computedCost += algoDoneCost + rate * algoNotDoneCost;
1719       }
1720       // cout << "rate: "<<rate << " algoNotDoneCost: " << algoNotDoneCost << endl;
1721     }
1722
1723   // get cost of already treated sub-meshes
1724   if ( SMESH_subMesh* mainSM = GetSubMeshContaining( 1 ))
1725   {
1726     SMESH_subMeshIteratorPtr smIt = mainSM->getDependsOnIterator(/*includeSelf=*/true);
1727     while ( smIt->more() )
1728     {
1729       const SMESH_subMesh* sm = smIt->next();
1730       const int smCost = sm->GetComputeCost();
1731       totalCost += smCost;
1732       if ( !currentSubIds.Contains( sm->GetId() ) )
1733       {
1734         if (( !sm->IsEmpty() ) ||
1735             ( sm->GetComputeState() == SMESH_subMesh::FAILED_TO_COMPUTE &&
1736               !sm->DependsOn( curSM ) ))
1737           computedCost += smCost;
1738       }
1739     }
1740   }
1741   // cout << "Total: " << totalCost
1742   //      << " computed: " << computedCost << " progress: " << computedCost / totalCost
1743   //      << " nbElems: " << GetMeshDS()->GetMeshInfo().NbElements() << endl;
1744   return computedCost / totalCost;
1745 }
1746
1747 //================================================================================
1748 /*!
1749  * \brief Return number of nodes in the mesh
1750  */
1751 //================================================================================
1752
1753 int SMESH_Mesh::NbNodes() const throw(SALOME_Exception)
1754 {
1755   Unexpect aCatch(SalomeException);
1756   return _myMeshDS->NbNodes();
1757 }
1758
1759 //================================================================================
1760 /*!
1761  * \brief  Return number of edges of given order in the mesh
1762  */
1763 //================================================================================
1764
1765 int SMESH_Mesh::Nb0DElements() const throw(SALOME_Exception)
1766 {
1767   Unexpect aCatch(SalomeException);
1768   return _myMeshDS->GetMeshInfo().Nb0DElements();
1769 }
1770
1771 //================================================================================
1772 /*!
1773  * \brief  Return number of edges of given order in the mesh
1774  */
1775 //================================================================================
1776
1777 int SMESH_Mesh::NbEdges(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1778 {
1779   Unexpect aCatch(SalomeException);
1780   return _myMeshDS->GetMeshInfo().NbEdges(order);
1781 }
1782
1783 //================================================================================
1784 /*!
1785  * \brief Return number of faces of given order in the mesh
1786  */
1787 //================================================================================
1788
1789 int SMESH_Mesh::NbFaces(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1790 {
1791   Unexpect aCatch(SalomeException);
1792   return _myMeshDS->GetMeshInfo().NbFaces(order);
1793 }
1794
1795 //================================================================================
1796 /*!
1797  * \brief Return the number of faces in the mesh
1798  */
1799 //================================================================================
1800
1801 int SMESH_Mesh::NbTriangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1802 {
1803   Unexpect aCatch(SalomeException);
1804   return _myMeshDS->GetMeshInfo().NbTriangles(order);
1805 }
1806
1807 //================================================================================
1808 /*!
1809  * \brief Return number of biquadratic triangles in the mesh
1810  */
1811 //================================================================================
1812
1813 int SMESH_Mesh::NbBiQuadTriangles() const throw(SALOME_Exception)
1814 {
1815   Unexpect aCatch(SalomeException);
1816   return _myMeshDS->GetMeshInfo().NbBiQuadTriangles();
1817 }
1818
1819 //================================================================================
1820 /*!
1821  * \brief Return the number nodes faces in the mesh
1822  */
1823 //================================================================================
1824
1825 int SMESH_Mesh::NbQuadrangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1826 {
1827   Unexpect aCatch(SalomeException);
1828   return _myMeshDS->GetMeshInfo().NbQuadrangles(order);
1829 }
1830
1831 //================================================================================
1832 /*!
1833  * \brief Return number of biquadratic quadrangles in the mesh
1834  */
1835 //================================================================================
1836
1837 int SMESH_Mesh::NbBiQuadQuadrangles() const throw(SALOME_Exception)
1838 {
1839   Unexpect aCatch(SalomeException);
1840   return _myMeshDS->GetMeshInfo().NbBiQuadQuadrangles();
1841 }
1842
1843 //================================================================================
1844 /*!
1845  * \brief Return the number of polygonal faces in the mesh
1846  */
1847 //================================================================================
1848
1849 int SMESH_Mesh::NbPolygons(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1850 {
1851   Unexpect aCatch(SalomeException);
1852   return _myMeshDS->GetMeshInfo().NbPolygons(order);
1853 }
1854
1855 //================================================================================
1856 /*!
1857  * \brief Return number of volumes of given order in the mesh
1858  */
1859 //================================================================================
1860
1861 int SMESH_Mesh::NbVolumes(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1862 {
1863   Unexpect aCatch(SalomeException);
1864   return _myMeshDS->GetMeshInfo().NbVolumes(order);
1865 }
1866
1867 //================================================================================
1868 /*!
1869  * \brief  Return number of tetrahedrons of given order in the mesh
1870  */
1871 //================================================================================
1872
1873 int SMESH_Mesh::NbTetras(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1874 {
1875   Unexpect aCatch(SalomeException);
1876   return _myMeshDS->GetMeshInfo().NbTetras(order);
1877 }
1878
1879 //================================================================================
1880 /*!
1881  * \brief  Return number of hexahedrons of given order in the mesh
1882  */
1883 //================================================================================
1884
1885 int SMESH_Mesh::NbHexas(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1886 {
1887   Unexpect aCatch(SalomeException);
1888   return _myMeshDS->GetMeshInfo().NbHexas(order);
1889 }
1890
1891 //================================================================================
1892 /*!
1893  * \brief  Return number of triquadratic hexahedrons in the mesh
1894  */
1895 //================================================================================
1896
1897 int SMESH_Mesh::NbTriQuadraticHexas() const throw(SALOME_Exception)
1898 {
1899   Unexpect aCatch(SalomeException);
1900   return _myMeshDS->GetMeshInfo().NbTriQuadHexas();
1901 }
1902
1903 //================================================================================
1904 /*!
1905  * \brief  Return number of pyramids of given order in the mesh
1906  */
1907 //================================================================================
1908
1909 int SMESH_Mesh::NbPyramids(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1910 {
1911   Unexpect aCatch(SalomeException);
1912   return _myMeshDS->GetMeshInfo().NbPyramids(order);
1913 }
1914
1915 //================================================================================
1916 /*!
1917  * \brief  Return number of prisms (penthahedrons) of given order in the mesh
1918  */
1919 //================================================================================
1920
1921 int SMESH_Mesh::NbPrisms(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1922 {
1923   Unexpect aCatch(SalomeException);
1924   return _myMeshDS->GetMeshInfo().NbPrisms(order);
1925 }
1926
1927 int SMESH_Mesh::NbQuadPrisms() const throw (SALOME_Exception)
1928 {
1929   Unexpect aCatch(SalomeException);
1930   return _myMeshDS->GetMeshInfo().NbQuadPrisms();
1931 }
1932
1933 int SMESH_Mesh::NbBiQuadPrisms() const throw (SALOME_Exception)
1934 {
1935   Unexpect aCatch(SalomeException);
1936   return _myMeshDS->GetMeshInfo().NbBiQuadPrisms();
1937 }
1938
1939
1940 //================================================================================
1941 /*!
1942  * \brief  Return number of hexagonal prisms in the mesh
1943  */
1944 //================================================================================
1945
1946 int SMESH_Mesh::NbHexagonalPrisms() const throw(SALOME_Exception)
1947 {
1948   Unexpect aCatch(SalomeException);
1949   return _myMeshDS->GetMeshInfo().NbHexPrisms();
1950 }
1951
1952 //================================================================================
1953 /*!
1954  * \brief  Return number of polyhedrons in the mesh
1955  */
1956 //================================================================================
1957
1958 int SMESH_Mesh::NbPolyhedrons() const throw(SALOME_Exception)
1959 {
1960   Unexpect aCatch(SalomeException);
1961   return _myMeshDS->GetMeshInfo().NbPolyhedrons();
1962 }
1963
1964 //================================================================================
1965 /*!
1966  * \brief  Return number of ball elements in the mesh
1967  */
1968 //================================================================================
1969
1970 int SMESH_Mesh::NbBalls() const throw(SALOME_Exception)
1971 {
1972   Unexpect aCatch(SalomeException);
1973   return _myMeshDS->GetMeshInfo().NbBalls();
1974 }
1975
1976 //================================================================================
1977 /*!
1978  * \brief  Return number of submeshes in the mesh
1979  */
1980 //================================================================================
1981
1982 int SMESH_Mesh::NbSubMesh() const throw(SALOME_Exception)
1983 {
1984   Unexpect aCatch(SalomeException);
1985   return _myMeshDS->NbSubMesh();
1986 }
1987
1988 //================================================================================
1989 /*!
1990  * \brief Returns number of meshes in the Study, that is supposed to be
1991  *        equal to SMESHDS_Document::NbMeshes()
1992  */
1993 //================================================================================
1994
1995 int SMESH_Mesh::NbMeshes() const // nb meshes in the Study
1996 {
1997   return _myDocument->NbMeshes();
1998 }
1999
2000 //=======================================================================
2001 //function : IsNotConformAllowed
2002 //purpose  : check if a hypothesis allowing notconform mesh is present
2003 //=======================================================================
2004
2005 bool SMESH_Mesh::IsNotConformAllowed() const
2006 {
2007   if(MYDEBUG) MESSAGE("SMESH_Mesh::IsNotConformAllowed");
2008
2009   static SMESH_HypoFilter filter( SMESH_HypoFilter::HasName( "NotConformAllowed" ));
2010   return GetHypothesis( _myMeshDS->ShapeToMesh(), filter, false );
2011 }
2012
2013 //=======================================================================
2014 //function : IsMainShape
2015 //purpose  : 
2016 //=======================================================================
2017
2018 bool SMESH_Mesh::IsMainShape(const TopoDS_Shape& theShape) const
2019 {
2020   return theShape.IsSame(_myMeshDS->ShapeToMesh() );
2021 }
2022
2023 //=============================================================================
2024 /*!
2025  *  
2026  */
2027 //=============================================================================
2028
2029 SMESH_Group* SMESH_Mesh::AddGroup (const SMDSAbs_ElementType theType,
2030                                    const char*               theName,
2031                                    const int                 theId,
2032                                    const TopoDS_Shape&       theShape,
2033                                    const SMESH_PredicatePtr& thePredicate)
2034 {
2035   if ( _mapGroup.count( theId ))
2036     return NULL;
2037   int id = ( theId < 0 ) ? _groupId : theId;
2038   SMESH_Group* aGroup = new SMESH_Group ( id, this, theType, theName, theShape, thePredicate );
2039   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
2040   _mapGroup[ id ] = aGroup;
2041   _groupId = 1 + _mapGroup.rbegin()->first;
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 = 1 + _mapGroup.rbegin()->first;
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) const
2119 {
2120   std::map <int, SMESH_Group*>::const_iterator id_grp = _mapGroup.find( theGroupID );
2121   if ( id_grp == _mapGroup.end() )
2122     return NULL;
2123   return id_grp->second;
2124 }
2125
2126
2127 //=============================================================================
2128 /*!
2129  * \brief Return IDs of all groups
2130  */
2131 //=============================================================================
2132
2133 std::list<int> SMESH_Mesh::GetGroupIds() const
2134 {
2135   std::list<int> anIds;
2136   std::map<int, SMESH_Group*>::const_iterator it = _mapGroup.begin();
2137   for ( ; it != _mapGroup.end(); it++ )
2138     anIds.push_back( it->first );
2139   
2140   return anIds;
2141 }
2142
2143 //================================================================================
2144 /*!
2145  * \brief Set a caller of methods at level of CORBA API implementation.
2146  * The set upCaller will be deleted by SMESH_Mesh
2147  */
2148 //================================================================================
2149
2150 void SMESH_Mesh::SetCallUp( TCallUp* upCaller )
2151 {
2152   if ( _callUp ) delete _callUp;
2153   _callUp = upCaller;
2154 }
2155
2156 //=============================================================================
2157 /*!
2158  *  
2159  */
2160 //=============================================================================
2161
2162 bool SMESH_Mesh::RemoveGroup( const int theGroupID )
2163 {
2164   if (_mapGroup.find(theGroupID) == _mapGroup.end())
2165     return false;
2166   GetMeshDS()->RemoveGroup( _mapGroup[theGroupID]->GetGroupDS() );
2167   delete _mapGroup[theGroupID];
2168   _mapGroup.erase (theGroupID);
2169   if (_callUp)
2170     _callUp->RemoveGroup( theGroupID );
2171   return true;
2172 }
2173
2174 //=======================================================================
2175 //function : GetAncestors
2176 //purpose  : return list of ancestors of theSubShape in the order
2177 //           that lower dimension shapes come first.
2178 //=======================================================================
2179
2180 const TopTools_ListOfShape& SMESH_Mesh::GetAncestors(const TopoDS_Shape& theS) const
2181 {
2182   if ( _mapAncestors.Contains( theS ) )
2183     return _mapAncestors.FindFromKey( theS );
2184
2185   static TopTools_ListOfShape emptyList;
2186   return emptyList;
2187 }
2188
2189 //=======================================================================
2190 //function : Dump
2191 //purpose  : dumps contents of mesh to stream [ debug purposes ]
2192 //=======================================================================
2193
2194 ostream& SMESH_Mesh::Dump(ostream& save)
2195 {
2196   int clause = 0;
2197   save << "========================== Dump contents of mesh ==========================" << endl << endl;
2198   save << ++clause << ") Total number of nodes:      \t" << NbNodes() << endl;
2199   save << ++clause << ") Total number of edges:      \t" << NbEdges() << endl;
2200   save << ++clause << ") Total number of faces:      \t" << NbFaces() << endl;
2201   save << ++clause << ") Total number of polygons:   \t" << NbPolygons() << endl;
2202   save << ++clause << ") Total number of volumes:    \t" << NbVolumes() << endl;
2203   save << ++clause << ") Total number of polyhedrons:\t" << NbPolyhedrons() << endl << endl;
2204   for ( int isQuadratic = 0; isQuadratic < 2; ++isQuadratic )
2205   {
2206     std::string orderStr = isQuadratic ? "quadratic" : "linear";
2207     SMDSAbs_ElementOrder order  = isQuadratic ? ORDER_QUADRATIC : ORDER_LINEAR;
2208
2209     save << ++clause << ") Total number of " << orderStr << " edges:\t" << NbEdges(order) << endl;
2210     save << ++clause << ") Total number of " << orderStr << " faces:\t" << NbFaces(order) << endl;
2211     if ( NbFaces(order) > 0 ) {
2212       int nb3 = NbTriangles(order);
2213       int nb4 = NbQuadrangles(order);
2214       save << clause << ".1) Number of " << orderStr << " triangles:  \t" << nb3 << endl;
2215       save << clause << ".2) Number of " << orderStr << " quadrangles:\t" << nb4 << endl;
2216       if ( nb3 + nb4 !=  NbFaces(order) ) {
2217         std::map<int,int> myFaceMap;
2218         SMDS_FaceIteratorPtr itFaces=_myMeshDS->facesIterator();
2219         while( itFaces->more( ) ) {
2220           int nbNodes = itFaces->next()->NbNodes();
2221           if ( myFaceMap.find( nbNodes ) == myFaceMap.end() )
2222             myFaceMap[ nbNodes ] = 0;
2223           myFaceMap[ nbNodes ] = myFaceMap[ nbNodes ] + 1;
2224         }
2225         save << clause << ".3) Faces in detail: " << endl;
2226         std::map <int,int>::iterator itF;
2227         for (itF = myFaceMap.begin(); itF != myFaceMap.end(); itF++)
2228           save << "--> nb nodes: " << itF->first << " - nb elements:\t" << itF->second << endl;
2229       }
2230     }
2231     save << ++clause << ") Total number of " << orderStr << " volumes:\t" << NbVolumes(order) << endl;
2232     if ( NbVolumes(order) > 0 ) {
2233       int nb8 = NbHexas(order);
2234       int nb4 = NbTetras(order);
2235       int nb5 = NbPyramids(order);
2236       int nb6 = NbPrisms(order);
2237       save << clause << ".1) Number of " << orderStr << " hexahedrons: \t" << nb8 << endl;
2238       save << clause << ".2) Number of " << orderStr << " tetrahedrons:\t" << nb4 << endl;
2239       save << clause << ".3) Number of " << orderStr << " prisms:      \t" << nb6 << endl;
2240       save << clause << ".4) Number of " << orderStr << " pyramids:    \t" << nb5 << endl;
2241       if ( nb8 + nb4 + nb5 + nb6 != NbVolumes(order) ) {
2242         std::map<int,int> myVolumesMap;
2243         SMDS_VolumeIteratorPtr itVolumes=_myMeshDS->volumesIterator();
2244         while( itVolumes->more( ) ) {
2245           int nbNodes = itVolumes->next()->NbNodes();
2246           if ( myVolumesMap.find( nbNodes ) == myVolumesMap.end() )
2247             myVolumesMap[ nbNodes ] = 0;
2248           myVolumesMap[ nbNodes ] = myVolumesMap[ nbNodes ] + 1;
2249         }
2250         save << clause << ".5) Volumes in detail: " << endl;
2251         std::map <int,int>::iterator itV;
2252         for (itV = myVolumesMap.begin(); itV != myVolumesMap.end(); itV++)
2253           save << "--> nb nodes: " << itV->first << " - nb elements:\t" << itV->second << endl;
2254       }
2255     }
2256     save << endl;
2257   }
2258   save << "===========================================================================" << endl;
2259   return save;
2260 }
2261
2262 //=======================================================================
2263 //function : GetElementType
2264 //purpose  : Returns type of mesh element with certain id
2265 //=======================================================================
2266
2267 SMDSAbs_ElementType SMESH_Mesh::GetElementType( const int id, const bool iselem )
2268 {
2269   return _myMeshDS->GetElementType( id, iselem );
2270 }
2271
2272 //=============================================================================
2273 /*!
2274  *  \brief Convert group on geometry into standalone group
2275  */
2276 //=============================================================================
2277
2278 SMESH_Group* SMESH_Mesh::ConvertToStandalone ( int theGroupID )
2279 {
2280   SMESH_Group* aGroup = 0;
2281   std::map < int, SMESH_Group * >::iterator itg = _mapGroup.find( theGroupID );
2282   if ( itg == _mapGroup.end() )
2283     return aGroup;
2284
2285   SMESH_Group* anOldGrp = (*itg).second;
2286   if ( !anOldGrp || !anOldGrp->GetGroupDS() )
2287     return aGroup;
2288   SMESHDS_GroupBase* anOldGrpDS = anOldGrp->GetGroupDS();
2289
2290   // create new standalone group
2291   aGroup = new SMESH_Group (theGroupID, this, anOldGrpDS->GetType(), anOldGrp->GetName() );
2292   _mapGroup[theGroupID] = aGroup;
2293
2294   SMESHDS_Group* aNewGrpDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
2295   GetMeshDS()->RemoveGroup( anOldGrpDS );
2296   GetMeshDS()->AddGroup( aNewGrpDS );
2297
2298   // add elements (or nodes) into new created group
2299   SMDS_ElemIteratorPtr anItr = anOldGrpDS->GetElements();
2300   while ( anItr->more() )
2301     aNewGrpDS->Add( (anItr->next())->GetID() );
2302
2303   // set color
2304   aNewGrpDS->SetColor( anOldGrpDS->GetColor() );
2305
2306   // remove old group
2307   delete anOldGrp;
2308
2309   return aGroup;
2310 }
2311
2312 //=============================================================================
2313 /*!
2314  *  \brief remove submesh order  from Mesh
2315  */
2316 //=============================================================================
2317
2318 void SMESH_Mesh::ClearMeshOrder()
2319 {
2320   _mySubMeshOrder.clear();
2321 }
2322
2323 //=============================================================================
2324 /*!
2325  *  \brief remove submesh order  from Mesh
2326  */
2327 //=============================================================================
2328
2329 void SMESH_Mesh::SetMeshOrder(const TListOfListOfInt& theOrder )
2330 {
2331   _mySubMeshOrder = theOrder;
2332 }
2333
2334 //=============================================================================
2335 /*!
2336  *  \brief return submesh order if any
2337  */
2338 //=============================================================================
2339
2340 const TListOfListOfInt& SMESH_Mesh::GetMeshOrder() const
2341 {
2342   return _mySubMeshOrder;
2343 }
2344
2345 //=============================================================================
2346 /*!
2347  *  \brief fill _mapAncestors
2348  */
2349 //=============================================================================
2350
2351 void SMESH_Mesh::fillAncestorsMap(const TopoDS_Shape& theShape)
2352 {
2353   int desType, ancType;
2354   if ( !theShape.IsSame( GetShapeToMesh()) && theShape.ShapeType() == TopAbs_COMPOUND )
2355   {
2356     // a geom group is added. Insert it into lists of ancestors before
2357     // the first ancestor more complex than group members
2358     TopoDS_Iterator subIt( theShape );
2359     if ( !subIt.More() ) return;
2360     int memberType = subIt.Value().ShapeType();
2361     for ( desType = TopAbs_VERTEX; desType >= memberType; desType-- )
2362       for (TopExp_Explorer des( theShape, TopAbs_ShapeEnum( desType )); des.More(); des.Next())
2363       {
2364         if ( !_mapAncestors.Contains( des.Current() )) continue;// issue 0020982
2365         TopTools_ListOfShape& ancList = _mapAncestors.ChangeFromKey( des.Current() );
2366         TopTools_ListIteratorOfListOfShape ancIt (ancList);
2367         while ( ancIt.More() && ancIt.Value().ShapeType() >= memberType )
2368           ancIt.Next();
2369         if ( ancIt.More() ) ancList.InsertBefore( theShape, ancIt );
2370         else                ancList.Append( theShape );
2371       }
2372   }
2373   else // else added for 52457: Addition of hypotheses is 8 time longer than meshing
2374   {
2375     for ( desType = TopAbs_VERTEX; desType > TopAbs_COMPOUND; desType-- )
2376       for ( ancType = desType - 1; ancType >= TopAbs_COMPOUND; ancType-- )
2377         TopExp::MapShapesAndAncestors ( theShape,
2378                                         (TopAbs_ShapeEnum) desType,
2379                                         (TopAbs_ShapeEnum) ancType,
2380                                         _mapAncestors );
2381   }
2382   // visit COMPOUNDs inside a COMPOUND that are not reachable by TopExp_Explorer
2383   if ( theShape.ShapeType() == TopAbs_COMPOUND )
2384   {
2385     TopoDS_Iterator sIt(theShape);
2386     if ( sIt.More() && sIt.Value().ShapeType() == TopAbs_COMPOUND )
2387       for ( ; sIt.More(); sIt.Next() )
2388         if ( sIt.Value().ShapeType() == TopAbs_COMPOUND )
2389           fillAncestorsMap( sIt.Value() );
2390   }
2391 }
2392
2393 //=============================================================================
2394 /*!
2395  * \brief sort submeshes according to stored mesh order
2396  * \param theListToSort in out list to be sorted
2397  * \return FALSE if nothing sorted
2398  */
2399 //=============================================================================
2400
2401 bool SMESH_Mesh::SortByMeshOrder(std::vector<SMESH_subMesh*>& theListToSort) const
2402 {
2403   if ( !_mySubMeshOrder.size() || theListToSort.size() < 2)
2404     return true;
2405   
2406   bool res = false;
2407   std::vector<SMESH_subMesh*> onlyOrderedList, smVec;
2408
2409   // collect all ordered submeshes in one list as pointers
2410   // and get their positions within theListToSort
2411   typedef std::vector<SMESH_subMesh*>::iterator TPosInList;
2412   std::map< int, TPosInList > sortedPos;
2413   TPosInList smBeg = theListToSort.begin(), smEnd = theListToSort.end();
2414   TListOfListOfInt::const_iterator      listIdsIt = _mySubMeshOrder.begin();
2415   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2416   {
2417     const TListOfInt& listOfId = *listIdsIt;
2418     // convert sm ids to sm's
2419     smVec.clear();
2420     TListOfInt::const_iterator idIt = listOfId.begin();
2421     for ( ; idIt != listOfId.end(); idIt++ )
2422     {
2423       if ( SMESH_subMesh * sm = GetSubMeshContaining( *idIt ))
2424       {
2425         smVec.push_back( sm );
2426         if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->IsComplexSubmesh() )
2427         {
2428           SMESHDS_SubMeshIteratorPtr smdsIt = sm->GetSubMeshDS()->GetSubMeshIterator();
2429           while ( smdsIt->more() )
2430           {
2431             const SMESHDS_SubMesh* smDS = smdsIt->next();
2432             if (( sm = GetSubMeshContaining( smDS->GetID() )))
2433               smVec.push_back( sm );
2434           }
2435         }
2436       }
2437     }
2438     // find smVec items in theListToSort
2439     for ( size_t i = 0; i < smVec.size(); ++i )
2440     {
2441       TPosInList smPos = find( smBeg, smEnd, smVec[i] );
2442       if ( smPos != smEnd ) {
2443         sortedPos[ std::distance( smBeg, smPos )] = smPos;
2444         if ( sortedPos.size() > onlyOrderedList.size() )
2445           onlyOrderedList.push_back( smVec[i] );
2446       }
2447     }
2448   }
2449   if (onlyOrderedList.size() < 2)
2450     return res;
2451   res = true;
2452
2453   std::vector<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
2454   std::vector<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
2455
2456   // iterate on ordered sub-meshes and insert them in detected positions
2457   std::map< int, TPosInList >::iterator i_pos = sortedPos.begin();
2458   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
2459     *(i_pos->second) = *onlyBIt;
2460
2461   return res;
2462 }
2463
2464 //================================================================================
2465 /*!
2466  * \brief Return true if given order of sub-meshes is OK
2467  */
2468 //================================================================================
2469
2470 bool SMESH_Mesh::IsOrderOK( const SMESH_subMesh* smBefore,
2471                             const SMESH_subMesh* smAfter ) const
2472 {
2473   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2474   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2475   {
2476     const TListOfInt& listOfId = *listIdsIt;
2477     int iB = -1, iA = -1, i = 0;
2478     for ( TListOfInt::const_iterator id = listOfId.begin(); id != listOfId.end(); ++id, ++i )
2479     {
2480       if ( *id == smBefore->GetId() )
2481       {
2482         iB = i;
2483         if ( iA > -1 )
2484           return iB < iA;
2485       }
2486       else if ( *id == smAfter->GetId() )
2487       {
2488         iA = i;
2489         if ( iB > -1 )
2490           return iB < iA;
2491       }
2492     }
2493   }
2494   return true; // no order imposed to given sub-meshes
2495
2496
2497 //=============================================================================
2498 /*!
2499  * \brief sort submeshes according to stored mesh order
2500  * \param theListToSort in out list to be sorted
2501  * \return FALSE if nothing sorted
2502  */
2503 //=============================================================================
2504
2505 void SMESH_Mesh::getAncestorsSubMeshes (const TopoDS_Shape&            theSubShape,
2506                                         std::vector< SMESH_subMesh* >& theSubMeshes) const
2507 {
2508   theSubMeshes.clear();
2509   TopTools_ListIteratorOfListOfShape it( GetAncestors( theSubShape ));
2510   for (; it.More(); it.Next() )
2511     if ( SMESH_subMesh* sm = GetSubMeshContaining( it.Value() ))
2512       theSubMeshes.push_back(sm);
2513
2514   // sort submeshes according to stored mesh order
2515   SortByMeshOrder( theSubMeshes );
2516 }