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