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