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