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