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