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