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