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