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