Salome HOME
61cfc4a94454bf03d7004a8fd63ce1455ecf086c
[modules/smesh.git] / src / SMESH / SMESH_Mesh.cxx
1 // Copyright (C) 2007-2011  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.
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 //  SMESH SMESH : implementaion of SMESH idl descriptions
24 //  File   : SMESH_Mesh.cxx
25 //  Author : Paul RASCLE, EDF
26 //  Module : SMESH
27 //
28 #include "SMESH_Mesh.hxx"
29 #include "SMESH_subMesh.hxx"
30 #include "SMESH_Gen.hxx"
31 #include "SMESH_Hypothesis.hxx"
32 #include "SMESH_Group.hxx"
33 #include "SMESH_HypoFilter.hxx"
34 #include "SMESHDS_Group.hxx"
35 #include "SMESHDS_Script.hxx"
36 #include "SMESHDS_GroupOnGeom.hxx"
37 #include "SMESHDS_Document.hxx"
38 #include "SMDS_MeshVolume.hxx"
39 #include "SMDS_SetIterator.hxx"
40
41 #include "utilities.h"
42
43 #include "DriverMED_W_SMESHDS_Mesh.h"
44 #include "DriverDAT_W_SMDS_Mesh.h"
45 #include "DriverUNV_W_SMDS_Mesh.h"
46 #include "DriverSTL_W_SMDS_Mesh.h"
47
48 #include "DriverMED_R_SMESHDS_Mesh.h"
49 #include "DriverUNV_R_SMDS_Mesh.h"
50 #include "DriverSTL_R_SMDS_Mesh.h"
51 #ifdef WITH_CGNS
52 #include "DriverCGNS_Read.hxx"
53 #include "DriverCGNS_Write.hxx"
54 #endif
55
56 #undef _Precision_HeaderFile
57 #include <BRepBndLib.hxx>
58 #include <BRepPrimAPI_MakeBox.hxx>
59 #include <Bnd_Box.hxx>
60 #include <TopExp.hxx>
61 #include <TopExp_Explorer.hxx>
62 #include <TopTools_ListIteratorOfListOfShape.hxx>
63 #include <TopTools_ListOfShape.hxx>
64 #include <TopTools_MapOfShape.hxx>
65 #include <TopoDS_Iterator.hxx>
66
67 #include "Utils_ExceptHandlers.hxx"
68
69 using namespace std;
70
71 // maximum stored group name length in MED file
72 #define MAX_MED_GROUP_NAME_LENGTH 80
73
74 #ifdef _DEBUG_
75 static int MYDEBUG = 0;
76 #else
77 static int MYDEBUG = 0;
78 #endif
79
80 #define cSMESH_Hyp(h) static_cast<const SMESH_Hypothesis*>(h)
81
82 typedef SMESH_HypoFilter THypType;
83
84 //=============================================================================
85 /*!
86  * 
87  */
88 //=============================================================================
89
90 SMESH_Mesh::SMESH_Mesh(int               theLocalId, 
91                        int               theStudyId, 
92                        SMESH_Gen*        theGen,
93                        bool              theIsEmbeddedMode,
94                        SMESHDS_Document* theDocument):
95   _groupId( 0 ), _nbSubShapes( 0 )
96 {
97   MESSAGE("SMESH_Mesh::SMESH_Mesh(int localId)");
98   _id            = theLocalId;
99   _studyId       = theStudyId;
100   _gen           = theGen;
101   _myDocument    = theDocument;
102   _idDoc         = theDocument->NewMesh(theIsEmbeddedMode);
103   _myMeshDS      = theDocument->GetMesh(_idDoc);
104   _isShapeToMesh = false;
105   _isAutoColor   = false;
106   _isModified    = false;
107   _shapeDiagonal = 0.0;
108   _rmGroupCallUp = 0;
109   _myMeshDS->ShapeToMesh( PseudoShape() );
110 }
111
112 //================================================================================
113 /*!
114  * \brief Constructor of SMESH_Mesh being a base of some descendant class
115  */
116 //================================================================================
117
118 SMESH_Mesh::SMESH_Mesh():
119   _groupId( 0 ), _nbSubShapes( 0 )
120 {
121   _myMeshDS      = 0;
122   _isShapeToMesh = false;
123   _isAutoColor   = false;
124   _isModified    = false;
125   _shapeDiagonal = 0.0;
126   _rmGroupCallUp = 0;
127 }
128
129 //=============================================================================
130 /*!
131  * 
132  */
133 //=============================================================================
134
135 SMESH_Mesh::~SMESH_Mesh()
136 {
137   MESSAGE("SMESH_Mesh::~SMESH_Mesh");
138
139   // issue 0020340: EDF 1022 SMESH : Crash with FindNodeClosestTo in a second new study
140   //   Notify event listeners at least that something happens
141   if ( SMESH_subMesh * sm = GetSubMeshContaining(1))
142     sm->ComputeStateEngine( SMESH_subMesh::MESH_ENTITY_REMOVED );
143
144   // delete groups
145   map < int, SMESH_Group * >::iterator itg;
146   for (itg = _mapGroup.begin(); itg != _mapGroup.end(); itg++) {
147     SMESH_Group *aGroup = (*itg).second;
148     delete aGroup;
149   }
150   _mapGroup.clear();
151
152   if ( _rmGroupCallUp) delete _rmGroupCallUp;
153   _rmGroupCallUp = 0;
154 }
155
156 //=============================================================================
157 /*!
158  * \brief Set geometry to be meshed
159  */
160 //=============================================================================
161
162 void SMESH_Mesh::ShapeToMesh(const TopoDS_Shape & aShape)
163 {
164   if(MYDEBUG) MESSAGE("SMESH_Mesh::ShapeToMesh");
165
166   if ( !aShape.IsNull() && _isShapeToMesh ) {
167     if ( aShape.ShapeType() != TopAbs_COMPOUND && // group contents is allowed to change
168          _myMeshDS->ShapeToMesh().ShapeType() != TopAbs_COMPOUND )
169       throw SALOME_Exception(LOCALIZED ("a shape to mesh has already been defined"));
170   }
171   // clear current data
172   if ( !_myMeshDS->ShapeToMesh().IsNull() )
173   {
174     // removal of a shape to mesh, delete objects referring to sub-shapes:
175     // - sub-meshes
176     map <int, SMESH_subMesh *>::iterator i_sm = _mapSubMesh.begin();
177     for ( ; i_sm != _mapSubMesh.end(); ++i_sm )
178       delete i_sm->second;
179     _mapSubMesh.clear();
180     //  - groups on geometry
181     map <int, SMESH_Group *>::iterator i_gr = _mapGroup.begin();
182     while ( i_gr != _mapGroup.end() ) {
183       if ( dynamic_cast<SMESHDS_GroupOnGeom*>( i_gr->second->GetGroupDS() )) {
184         _myMeshDS->RemoveGroup( i_gr->second->GetGroupDS() );
185         delete i_gr->second;
186         _mapGroup.erase( i_gr++ );
187       }
188       else
189         i_gr++;
190     }
191     _mapAncestors.Clear();
192
193     // clear SMESHDS
194     TopoDS_Shape aNullShape;
195     _myMeshDS->ShapeToMesh( aNullShape );
196
197     _shapeDiagonal = 0.0;
198   }
199
200   // set a new geometry
201   if ( !aShape.IsNull() )
202   {
203     _myMeshDS->ShapeToMesh(aShape);
204     _isShapeToMesh = true;
205     _nbSubShapes = _myMeshDS->MaxShapeIndex();
206
207     // fill map of ancestors
208     fillAncestorsMap(aShape);
209   }
210   else
211   {
212     _isShapeToMesh = false;
213     _shapeDiagonal = 0.0;
214     _myMeshDS->ShapeToMesh( PseudoShape() );
215   }
216   _isModified = false;
217 }
218
219 //=======================================================================
220 /*!
221  * \brief Return geometry to be meshed. (It may be a PseudoShape()!)
222  */
223 //=======================================================================
224
225 TopoDS_Shape SMESH_Mesh::GetShapeToMesh() const
226 {
227   return _myMeshDS->ShapeToMesh();
228 }
229
230 //=======================================================================
231 /*!
232  * \brief Return a solid which is returned by GetShapeToMesh() if
233  *        a real geometry to be meshed was not set
234  */
235 //=======================================================================
236
237 const TopoDS_Solid& SMESH_Mesh::PseudoShape()
238 {
239   static TopoDS_Solid aSolid;
240   if ( aSolid.IsNull() )
241   {
242     aSolid = BRepPrimAPI_MakeBox(1,1,1);
243   }
244   return aSolid;
245 }
246
247 //=======================================================================
248 /*!
249  * \brief Return diagonal size of bounding box of a shape
250  */
251 //=======================================================================
252
253 double SMESH_Mesh::GetShapeDiagonalSize(const TopoDS_Shape & aShape)
254 {
255   if ( !aShape.IsNull() ) {
256     Bnd_Box Box;
257     BRepBndLib::Add(aShape, Box);
258     return sqrt( Box.SquareExtent() );
259   }
260   return 0;
261 }
262
263 //=======================================================================
264 /*!
265  * \brief Return diagonal size of bounding box of shape to mesh
266  */
267 //=======================================================================
268
269 double SMESH_Mesh::GetShapeDiagonalSize() const
270 {
271   if ( _shapeDiagonal == 0. && _isShapeToMesh )
272     const_cast<SMESH_Mesh*>(this)->_shapeDiagonal = GetShapeDiagonalSize( GetShapeToMesh() );
273
274   return _shapeDiagonal;
275 }
276
277 //=======================================================================
278 /*!
279  * \brief Remove all nodes and elements
280  */
281 //=======================================================================
282
283 void SMESH_Mesh::Clear()
284 {
285   // clear mesh data
286   _myMeshDS->ClearMesh();
287
288   // update compute state of submeshes
289   if ( SMESH_subMesh *sm = GetSubMeshContaining( GetShapeToMesh() ) )
290   {
291     sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
292     sm->ComputeSubMeshStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
293     sm->ComputeStateEngine( SMESH_subMesh::CLEAN ); // for event listeners (issue 0020918)
294     sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
295   }
296   _isModified = false;
297 }
298
299 //=======================================================================
300 /*!
301  * \brief Remove all nodes and elements of indicated shape
302  */
303 //=======================================================================
304
305 void SMESH_Mesh::ClearSubMesh(const int theShapeId)
306 {
307   // clear sub-meshes; get ready to re-compute as a side-effect 
308   if ( SMESH_subMesh *sm = GetSubMeshContaining( theShapeId ) )
309   {
310     SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
311                                                              /*complexShapeFirst=*/false);
312     while ( smIt->more() )
313     {
314       sm = smIt->next();
315       TopAbs_ShapeEnum shapeType = sm->GetSubShape().ShapeType();      
316       if ( shapeType == TopAbs_VERTEX || shapeType < TopAbs_SOLID )
317         // all other shapes depends on vertices so they are already cleaned
318         sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
319       // to recompute even if failed
320       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
321     }
322   }
323 }
324
325 //=======================================================================
326 //function : UNVToMesh
327 //purpose  : 
328 //=======================================================================
329
330 int SMESH_Mesh::UNVToMesh(const char* theFileName)
331 {
332   if(MYDEBUG) MESSAGE("UNVToMesh - theFileName = "<<theFileName);
333   if(_isShapeToMesh)
334     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
335   _isShapeToMesh = false;
336   DriverUNV_R_SMDS_Mesh myReader;
337   myReader.SetMesh(_myMeshDS);
338   myReader.SetFile(theFileName);
339   myReader.SetMeshId(-1);
340   myReader.Perform();
341   if(MYDEBUG){
342     MESSAGE("UNVToMesh - _myMeshDS->NbNodes() = "<<_myMeshDS->NbNodes());
343     MESSAGE("UNVToMesh - _myMeshDS->NbEdges() = "<<_myMeshDS->NbEdges());
344     MESSAGE("UNVToMesh - _myMeshDS->NbFaces() = "<<_myMeshDS->NbFaces());
345     MESSAGE("UNVToMesh - _myMeshDS->NbVolumes() = "<<_myMeshDS->NbVolumes());
346   }
347   SMDS_MeshGroup* aGroup = (SMDS_MeshGroup*) myReader.GetGroup();
348   if (aGroup != 0) {
349     TGroupNamesMap aGroupNames = myReader.GetGroupNamesMap();
350     //const TGroupIdMap& aGroupId = myReader.GetGroupIdMap();
351     aGroup->InitSubGroupsIterator();
352     while (aGroup->MoreSubGroups()) {
353       SMDS_MeshGroup* aSubGroup = (SMDS_MeshGroup*) aGroup->NextSubGroup();
354       string aName = aGroupNames[aSubGroup];
355       int aId;
356
357       SMESH_Group* aSMESHGroup = AddGroup( aSubGroup->GetType(), aName.c_str(), aId );
358       if ( aSMESHGroup ) {
359         if(MYDEBUG) MESSAGE("UNVToMesh - group added: "<<aName);      
360         SMESHDS_Group* aGroupDS = dynamic_cast<SMESHDS_Group*>( aSMESHGroup->GetGroupDS() );
361         if ( aGroupDS ) {
362           aGroupDS->SetStoreName(aName.c_str());
363           aSubGroup->InitIterator();
364           const SMDS_MeshElement* aElement = 0;
365           while (aSubGroup->More()) {
366             aElement = aSubGroup->Next();
367             if (aElement) {
368               aGroupDS->SMDSGroup().Add(aElement);
369             }
370           }
371           if (aElement)
372             aGroupDS->SetType(aElement->GetType());
373         }
374       }
375     }
376   }
377   return 1;
378 }
379
380 //=======================================================================
381 //function : MEDToMesh
382 //purpose  : 
383 //=======================================================================
384
385 int SMESH_Mesh::MEDToMesh(const char* theFileName, const char* theMeshName)
386 {
387   if(MYDEBUG) MESSAGE("MEDToMesh - theFileName = "<<theFileName<<", mesh name = "<<theMeshName);
388   if(_isShapeToMesh)
389     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
390   _isShapeToMesh = false;
391   DriverMED_R_SMESHDS_Mesh myReader;
392   myReader.SetMesh(_myMeshDS);
393   myReader.SetMeshId(-1);
394   myReader.SetFile(theFileName);
395   myReader.SetMeshName(theMeshName);
396   Driver_Mesh::Status status = myReader.Perform();
397   if(MYDEBUG){
398     MESSAGE("MEDToMesh - _myMeshDS->NbNodes() = "<<_myMeshDS->NbNodes());
399     MESSAGE("MEDToMesh - _myMeshDS->NbEdges() = "<<_myMeshDS->NbEdges());
400     MESSAGE("MEDToMesh - _myMeshDS->NbFaces() = "<<_myMeshDS->NbFaces());
401     MESSAGE("MEDToMesh - _myMeshDS->NbVolumes() = "<<_myMeshDS->NbVolumes());
402   }
403
404   // Reading groups (sub-meshes are out of scope of MED import functionality)
405   list<TNameAndType> aGroupNames = myReader.GetGroupNamesAndTypes();
406   if(MYDEBUG) MESSAGE("MEDToMesh - Nb groups = "<<aGroupNames.size()); 
407   int anId;
408   list<TNameAndType>::iterator name_type = aGroupNames.begin();
409   for ( ; name_type != aGroupNames.end(); name_type++ ) {
410     SMESH_Group* aGroup = AddGroup( name_type->second, name_type->first.c_str(), anId );
411     if ( aGroup ) {
412       if(MYDEBUG) MESSAGE("MEDToMesh - group added: "<<name_type->first.c_str());      
413       SMESHDS_Group* aGroupDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
414       if ( aGroupDS ) {
415         aGroupDS->SetStoreName( name_type->first.c_str() );
416         myReader.GetGroup( aGroupDS );
417       }
418     }
419   }
420   return (int) status;
421 }
422
423 //=======================================================================
424 //function : STLToMesh
425 //purpose  : 
426 //=======================================================================
427
428 int SMESH_Mesh::STLToMesh(const char* theFileName)
429 {
430   if(MYDEBUG) MESSAGE("STLToMesh - theFileName = "<<theFileName);
431   if(_isShapeToMesh)
432     throw SALOME_Exception(LOCALIZED("a shape to mesh has already been defined"));
433   _isShapeToMesh = false;
434   DriverSTL_R_SMDS_Mesh myReader;
435   myReader.SetMesh(_myMeshDS);
436   myReader.SetFile(theFileName);
437   myReader.SetMeshId(-1);
438   myReader.Perform();
439   if(MYDEBUG){
440     MESSAGE("STLToMesh - _myMeshDS->NbNodes() = "<<_myMeshDS->NbNodes());
441     MESSAGE("STLToMesh - _myMeshDS->NbEdges() = "<<_myMeshDS->NbEdges());
442     MESSAGE("STLToMesh - _myMeshDS->NbFaces() = "<<_myMeshDS->NbFaces());
443     MESSAGE("STLToMesh - _myMeshDS->NbVolumes() = "<<_myMeshDS->NbVolumes());
444   }
445   return 1;
446 }
447
448 //================================================================================
449 /*!
450  * \brief Reads the given mesh from the CGNS file
451  *  \param theFileName - name of the file
452  *  \retval int - Driver_Mesh::Status
453  */
454 //================================================================================
455
456 int SMESH_Mesh::CGNSToMesh(const char*  theFileName,
457                            const int    theMeshIndex,
458                            std::string& theMeshName)
459 {
460   int res = Driver_Mesh::DRS_FAIL;
461 #ifdef WITH_CGNS
462
463   DriverCGNS_Read myReader;
464   myReader.SetMesh(_myMeshDS);
465   myReader.SetFile(theFileName);
466   myReader.SetMeshId(theMeshIndex);
467   res = myReader.Perform();
468   theMeshName = myReader.GetMeshName();
469
470   // create groups
471   SynchronizeGroups();
472
473 #endif
474   return res;
475 }
476
477 //=============================================================================
478 /*!
479  * 
480  */
481 //=============================================================================
482
483 SMESH_Hypothesis::Hypothesis_Status
484   SMESH_Mesh::AddHypothesis(const TopoDS_Shape & aSubShape,
485                             int                  anHypId  ) throw(SALOME_Exception)
486 {
487   Unexpect aCatch(SalomeException);
488   if(MYDEBUG) MESSAGE("SMESH_Mesh::AddHypothesis");
489
490   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
491   if ( !subMesh || !subMesh->GetId())
492     return SMESH_Hypothesis::HYP_BAD_SUBSHAPE;
493
494   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
495   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
496   {
497     if(MYDEBUG) MESSAGE("Hypothesis ID does not give an hypothesis");
498     if(MYDEBUG) {
499       SCRUTE(_studyId);
500       SCRUTE(anHypId);
501     }
502     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
503   }
504
505   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
506   MESSAGE( "SMESH_Mesh::AddHypothesis " << anHyp->GetName() );
507
508   bool isGlobalHyp = IsMainShape( aSubShape );
509
510   // NotConformAllowed can be only global
511   if ( !isGlobalHyp )
512   {
513     // NOTE: this is not a correct way to check a name of hypothesis,
514     // there should be an attribute of hypothesis saying that it can/can't
515     // be global/local
516     string hypName = anHyp->GetName();
517     if ( hypName == "NotConformAllowed" )
518     {
519       if(MYDEBUG) MESSAGE( "Hypotesis <NotConformAllowed> can be only global" );
520       return SMESH_Hypothesis::HYP_INCOMPATIBLE;
521     }
522   }
523
524   // shape 
525
526   bool isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
527   int event = isAlgo ? SMESH_subMesh::ADD_ALGO : SMESH_subMesh::ADD_HYP;
528
529   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
530
531   // subShapes
532   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
533       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is added on father
534   {
535     event = isAlgo ? SMESH_subMesh::ADD_FATHER_ALGO : SMESH_subMesh::ADD_FATHER_HYP;
536
537     SMESH_Hypothesis::Hypothesis_Status ret2 =
538       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
539     if (ret2 > ret)
540       ret = ret2;
541
542     // check concurent hypotheses on ancestors
543     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !isGlobalHyp )
544     {
545       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
546       while ( smIt->more() ) {
547         SMESH_subMesh* sm = smIt->next();
548         if ( sm->IsApplicableHypotesis( anHyp )) {
549           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
550           if (ret2 > ret) {
551             ret = ret2;
552             break;
553           }
554         }
555       }
556     }
557   }
558   HasModificationsToDiscard(); // to reset _isModified flag if a mesh becomes empty
559
560   GetMeshDS()->Modified();
561
562   if(MYDEBUG) subMesh->DumpAlgoState(true);
563   if(MYDEBUG) SCRUTE(ret);
564   return ret;
565 }
566
567 //=============================================================================
568 /*!
569  * 
570  */
571 //=============================================================================
572
573 SMESH_Hypothesis::Hypothesis_Status
574   SMESH_Mesh::RemoveHypothesis(const TopoDS_Shape & aSubShape,
575                                int anHypId)throw(SALOME_Exception)
576 {
577   Unexpect aCatch(SalomeException);
578   if(MYDEBUG) MESSAGE("SMESH_Mesh::RemoveHypothesis");
579   
580   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
581   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
582     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
583   
584   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
585   if(MYDEBUG) {
586     int hypType = anHyp->GetType();
587     SCRUTE(hypType);
588   }
589   
590   // shape 
591   
592   bool isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
593   int event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
594
595   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
596
597   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
598
599   // there may appear concurrent hyps that were covered by the removed hyp
600   if (ret < SMESH_Hypothesis::HYP_CONCURENT &&
601       subMesh->IsApplicableHypotesis( anHyp ) &&
602       subMesh->CheckConcurentHypothesis( anHyp->GetType() ) != SMESH_Hypothesis::HYP_OK)
603     ret = SMESH_Hypothesis::HYP_CONCURENT;
604
605   // subShapes
606   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
607       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is removed from father
608   {
609     event = isAlgo ? SMESH_subMesh::REMOVE_FATHER_ALGO : SMESH_subMesh::REMOVE_FATHER_HYP;
610
611     SMESH_Hypothesis::Hypothesis_Status ret2 =
612       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
613     if (ret2 > ret) // more severe
614       ret = ret2;
615
616     // check concurent hypotheses on ancestors
617     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !IsMainShape( aSubShape ) )
618     {
619       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
620       while ( smIt->more() ) {
621         SMESH_subMesh* sm = smIt->next();
622         if ( sm->IsApplicableHypotesis( anHyp )) {
623           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
624           if (ret2 > ret) {
625             ret = ret2;
626             break;
627           }
628         }
629       }
630     }
631   }
632
633   HasModificationsToDiscard(); // to reset _isModified flag if mesh become empty
634
635   GetMeshDS()->Modified();
636
637   if(MYDEBUG) subMesh->DumpAlgoState(true);
638   if(MYDEBUG) SCRUTE(ret);
639   return ret;
640 }
641
642 //=============================================================================
643 /*!
644  * 
645  */
646 //=============================================================================
647
648 const list<const SMESHDS_Hypothesis*>&
649 SMESH_Mesh::GetHypothesisList(const TopoDS_Shape & aSubShape) const
650   throw(SALOME_Exception)
651 {
652   Unexpect aCatch(SalomeException);
653   return _myMeshDS->GetHypothesis(aSubShape);
654 }
655
656 //=======================================================================
657 /*!
658  * \brief Return the hypothesis assigned to the shape
659  *  \param aSubShape    - the shape to check
660  *  \param aFilter      - the hypothesis filter
661  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
662  *  \param assignedTo   - to return the shape the found hypo is assigned to
663  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
664  */
665 //=======================================================================
666
667 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const TopoDS_Shape &    aSubShape,
668                                                    const SMESH_HypoFilter& aFilter,
669                                                    const bool              andAncestors,
670                                                    TopoDS_Shape*           assignedTo) const
671 {
672   {
673     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
674     list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
675     for ( ; hyp != hypList.end(); hyp++ ) {
676       const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
677       if ( aFilter.IsOk( h, aSubShape)) {
678         if ( assignedTo ) *assignedTo = aSubShape;
679         return h;
680       }
681     }
682   }
683   if ( andAncestors )
684   {
685     // user sorted submeshes of ancestors, according to stored submesh priority
686     const list<SMESH_subMesh*> smList = getAncestorsSubMeshes( aSubShape );
687     list<SMESH_subMesh*>::const_iterator smIt = smList.begin(); 
688     for ( ; smIt != smList.end(); smIt++ )
689     {
690       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
691       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
692       list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
693       for ( ; hyp != hypList.end(); hyp++ ) {
694         const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
695         if (aFilter.IsOk( h, curSh )) {
696           if ( assignedTo ) *assignedTo = curSh;
697           return h;
698         }
699       }
700     }
701   }
702   return 0;
703 }
704
705 //================================================================================
706 /*!
707  * \brief Return hypothesis assigned to the shape
708   * \param aSubShape - the shape to check
709   * \param aFilter - the hypothesis filter
710   * \param aHypList - the list of the found hypotheses
711   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
712   * \retval int - number of unique hypos in aHypList
713  */
714 //================================================================================
715
716 int SMESH_Mesh::GetHypotheses(const TopoDS_Shape &                aSubShape,
717                               const SMESH_HypoFilter&             aFilter,
718                               list <const SMESHDS_Hypothesis * >& aHypList,
719                               const bool                          andAncestors) const
720 {
721   set<string> hypTypes; // to exclude same type hypos from the result list
722   int nbHyps = 0;
723
724   // only one main hypothesis is allowed
725   bool mainHypFound = false;
726
727   // fill in hypTypes
728   list<const SMESHDS_Hypothesis*>::const_iterator hyp;
729   for ( hyp = aHypList.begin(); hyp != aHypList.end(); hyp++ ) {
730     if ( hypTypes.insert( (*hyp)->GetName() ).second )
731       nbHyps++;
732     if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
733       mainHypFound = true;
734   }
735
736   // get hypos from aSubShape
737   {
738     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
739     for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
740       if ( aFilter.IsOk (cSMESH_Hyp( *hyp ), aSubShape) &&
741            ( cSMESH_Hyp(*hyp)->IsAuxiliary() || !mainHypFound ) &&
742            hypTypes.insert( (*hyp)->GetName() ).second )
743       {
744         aHypList.push_back( *hyp );
745         nbHyps++;
746         if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
747           mainHypFound = true;
748       }
749   }
750
751   // get hypos from ancestors of aSubShape
752   if ( andAncestors )
753   {
754     TopTools_MapOfShape map;
755
756     // user sorted submeshes of ancestors, according to stored submesh priority
757     const list<SMESH_subMesh*> smList = getAncestorsSubMeshes( aSubShape );
758     list<SMESH_subMesh*>::const_iterator smIt = smList.begin(); 
759     for ( ; smIt != smList.end(); smIt++ )
760     {
761       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
762      if ( !map.Add( curSh ))
763         continue;
764       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
765       for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
766         if (aFilter.IsOk( cSMESH_Hyp( *hyp ), curSh ) &&
767             ( cSMESH_Hyp(*hyp)->IsAuxiliary() || !mainHypFound ) &&
768             hypTypes.insert( (*hyp)->GetName() ).second )
769         {
770           aHypList.push_back( *hyp );
771           nbHyps++;
772           if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
773             mainHypFound = true;
774         }
775     }
776   }
777   return nbHyps;
778 }
779
780 //=============================================================================
781 /*!
782  * 
783  */
784 //=============================================================================
785
786 const list<SMESHDS_Command*> & SMESH_Mesh::GetLog() throw(SALOME_Exception)
787 {
788   Unexpect aCatch(SalomeException);
789   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetLog");
790   return _myMeshDS->GetScript()->GetCommands();
791 }
792
793 //=============================================================================
794 /*!
795  * 
796  */
797 //=============================================================================
798 void SMESH_Mesh::ClearLog() throw(SALOME_Exception)
799 {
800   Unexpect aCatch(SalomeException);
801   if(MYDEBUG) MESSAGE("SMESH_Mesh::ClearLog");
802   _myMeshDS->GetScript()->Clear();
803 }
804
805 //=============================================================================
806 /*!
807  * Get or Create the SMESH_subMesh object implementation
808  */
809 //=============================================================================
810
811 SMESH_subMesh *SMESH_Mesh::GetSubMesh(const TopoDS_Shape & aSubShape)
812   throw(SALOME_Exception)
813 {
814   Unexpect aCatch(SalomeException);
815   SMESH_subMesh *aSubMesh;
816   int index = _myMeshDS->ShapeToIndex(aSubShape);
817
818   // for submeshes on GEOM Group
819   if (( !index || index > _nbSubShapes ) && aSubShape.ShapeType() == TopAbs_COMPOUND ) {
820     TopoDS_Iterator it( aSubShape );
821     if ( it.More() )
822     {
823       index = _myMeshDS->AddCompoundSubmesh( aSubShape, it.Value().ShapeType() );
824       if ( index > _nbSubShapes ) _nbSubShapes = index; // not to create sm for this group again
825
826       // fill map of Ancestors
827       fillAncestorsMap(aSubShape);
828     }
829   }
830 //   if ( !index )
831 //     return NULL; // neither sub-shape nor a group
832
833   map <int, SMESH_subMesh *>::iterator i_sm = _mapSubMesh.find(index);
834   if ( i_sm != _mapSubMesh.end())
835   {
836     aSubMesh = i_sm->second;
837   }
838   else
839   {
840     aSubMesh = new SMESH_subMesh(index, this, _myMeshDS, aSubShape);
841     _mapSubMesh[index] = aSubMesh;
842     ClearMeshOrder();
843   }
844   return aSubMesh;
845 }
846
847 //=============================================================================
848 /*!
849  * Get the SMESH_subMesh object implementation. Dont create it, return null
850  * if it does not exist.
851  */
852 //=============================================================================
853
854 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const TopoDS_Shape & aSubShape) const
855   throw(SALOME_Exception)
856 {
857   Unexpect aCatch(SalomeException);
858   SMESH_subMesh *aSubMesh = NULL;
859   
860   int index = _myMeshDS->ShapeToIndex(aSubShape);
861
862   map <int, SMESH_subMesh *>::const_iterator i_sm = _mapSubMesh.find(index);
863   if ( i_sm != _mapSubMesh.end())
864     aSubMesh = i_sm->second;
865
866   return aSubMesh;
867 }
868 //=============================================================================
869 /*!
870  * Get the SMESH_subMesh object implementation. Dont create it, return null
871  * if it does not exist.
872  */
873 //=============================================================================
874
875 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const int aShapeID) const
876 throw(SALOME_Exception)
877 {
878   Unexpect aCatch(SalomeException);
879   
880   map <int, SMESH_subMesh *>::const_iterator i_sm = _mapSubMesh.find(aShapeID);
881   if (i_sm == _mapSubMesh.end())
882     return NULL;
883   return i_sm->second;
884 }
885 //================================================================================
886 /*!
887  * \brief Return submeshes of groups containing the given subshape
888  */
889 //================================================================================
890
891 list<SMESH_subMesh*>
892 SMESH_Mesh::GetGroupSubMeshesContaining(const TopoDS_Shape & aSubShape) const
893   throw(SALOME_Exception)
894 {
895   Unexpect aCatch(SalomeException);
896   list<SMESH_subMesh*> found;
897
898   SMESH_subMesh * subMesh = GetSubMeshContaining(aSubShape);
899   if ( !subMesh )
900     return found;
901
902   // submeshes of groups have max IDs, so search from the map end
903   map<int, SMESH_subMesh *>::const_reverse_iterator i_sm;
904   for ( i_sm = _mapSubMesh.rbegin(); i_sm != _mapSubMesh.rend(); ++i_sm) {
905     SMESHDS_SubMesh * ds = i_sm->second->GetSubMeshDS();
906     if ( ds && ds->IsComplexSubmesh() ) {
907       TopExp_Explorer exp( i_sm->second->GetSubShape(), aSubShape.ShapeType() );
908       for ( ; exp.More(); exp.Next() ) {
909         if ( aSubShape.IsSame( exp.Current() )) {
910           found.push_back( i_sm->second );
911           break;
912         }
913       }
914     } else {
915       break;
916     }
917   }
918   return found;
919 }
920 //=======================================================================
921 //function : IsUsedHypothesis
922 //purpose  : Return True if anHyp is used to mesh aSubShape
923 //=======================================================================
924
925 bool SMESH_Mesh::IsUsedHypothesis(SMESHDS_Hypothesis * anHyp,
926                                   const SMESH_subMesh* aSubMesh)
927 {
928   SMESH_Hypothesis* hyp = static_cast<SMESH_Hypothesis*>(anHyp);
929
930   // check if anHyp can be used to mesh aSubMesh
931   if ( !aSubMesh || !aSubMesh->IsApplicableHypotesis( hyp ))
932     return false;
933
934   const TopoDS_Shape & aSubShape = const_cast<SMESH_subMesh*>( aSubMesh )->GetSubShape();
935
936   SMESH_Algo *algo = _gen->GetAlgo(*this, aSubShape );
937
938   // algorithm
939   if (anHyp->GetType() > SMESHDS_Hypothesis::PARAM_ALGO)
940     return ( anHyp == algo );
941
942   // algorithm parameter
943   if (algo)
944   {
945     // look trough hypotheses used by algo
946     SMESH_HypoFilter hypoKind;
947     if ( algo->InitCompatibleHypoFilter( hypoKind, !hyp->IsAuxiliary() )) {
948       list <const SMESHDS_Hypothesis * > usedHyps;
949       if ( GetHypotheses( aSubShape, hypoKind, usedHyps, true ))
950         return ( find( usedHyps.begin(), usedHyps.end(), anHyp ) != usedHyps.end() );
951     }
952   }
953
954   // look through all assigned hypotheses
955   //SMESH_HypoFilter filter( SMESH_HypoFilter::Is( hyp ));
956   return false; //GetHypothesis( aSubShape, filter, true );
957 }
958
959 //=============================================================================
960 /*!
961  *
962  */
963 //=============================================================================
964
965 const list < SMESH_subMesh * >&
966 SMESH_Mesh::GetSubMeshUsingHypothesis(SMESHDS_Hypothesis * anHyp)
967   throw(SALOME_Exception)
968 {
969   Unexpect aCatch(SalomeException);
970   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetSubMeshUsingHypothesis");
971   map < int, SMESH_subMesh * >::iterator itsm;
972   _subMeshesUsingHypothesisList.clear();
973   for (itsm = _mapSubMesh.begin(); itsm != _mapSubMesh.end(); itsm++)
974   {
975     SMESH_subMesh *aSubMesh = (*itsm).second;
976     if ( IsUsedHypothesis ( anHyp, aSubMesh ))
977       _subMeshesUsingHypothesisList.push_back(aSubMesh);
978   }
979   return _subMeshesUsingHypothesisList;
980 }
981
982 //=======================================================================
983 //function : NotifySubMeshesHypothesisModification
984 //purpose  : Say all submeshes using theChangedHyp that it has been modified
985 //=======================================================================
986
987 void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* hyp)
988 {
989   Unexpect aCatch(SalomeException);
990
991   const SMESH_Algo *foundAlgo = 0;
992   SMESH_HypoFilter algoKind, compatibleHypoKind;
993   list <const SMESHDS_Hypothesis * > usedHyps;
994
995
996   map < int, SMESH_subMesh * >::iterator itsm;
997   for (itsm = _mapSubMesh.begin(); itsm != _mapSubMesh.end(); itsm++)
998   {
999     SMESH_subMesh *aSubMesh = (*itsm).second;
1000     if ( aSubMesh->IsApplicableHypotesis( hyp ))
1001     {
1002       const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
1003
1004       if ( !foundAlgo ) // init filter for algo search
1005         algoKind.Init( THypType::IsAlgo() ).And( THypType::IsApplicableTo( aSubShape ));
1006       
1007       const SMESH_Algo *algo = static_cast<const SMESH_Algo*>
1008         ( GetHypothesis( aSubShape, algoKind, true ));
1009
1010       if ( algo )
1011       {
1012         bool sameAlgo = ( algo == foundAlgo );
1013         if ( !sameAlgo && foundAlgo )
1014           sameAlgo = ( strcmp( algo->GetName(), foundAlgo->GetName() ) == 0);
1015
1016         if ( !sameAlgo ) { // init filter for used hypos search
1017           if ( !algo->InitCompatibleHypoFilter( compatibleHypoKind, !hyp->IsAuxiliary() ))
1018             continue; // algo does not use any hypothesis
1019           foundAlgo = algo;
1020         }
1021
1022         // check if hyp is used by algo
1023         usedHyps.clear();
1024         if ( GetHypotheses( aSubShape, compatibleHypoKind, usedHyps, true ) &&
1025              find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() )
1026         {
1027           aSubMesh->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
1028                                     const_cast< SMESH_Hypothesis*>( hyp ));
1029         }
1030       }
1031     }
1032   }
1033   HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
1034   GetMeshDS()->Modified();
1035 }
1036
1037 //=============================================================================
1038 /*!
1039  *  Auto color functionality
1040  */
1041 //=============================================================================
1042 void SMESH_Mesh::SetAutoColor(bool theAutoColor) throw(SALOME_Exception)
1043 {
1044   Unexpect aCatch(SalomeException);
1045   _isAutoColor = theAutoColor;
1046 }
1047
1048 bool SMESH_Mesh::GetAutoColor() throw(SALOME_Exception)
1049 {
1050   Unexpect aCatch(SalomeException);
1051   return _isAutoColor;
1052 }
1053
1054 //=======================================================================
1055 //function : SetIsModified
1056 //purpose  : Set the flag meaning that the mesh has been edited "manually"
1057 //=======================================================================
1058
1059 void SMESH_Mesh::SetIsModified(bool isModified)
1060 {
1061   _isModified = isModified;
1062
1063   if ( _isModified )
1064     // check if mesh becomes empty as result of modification
1065     HasModificationsToDiscard();
1066 }
1067
1068 //=======================================================================
1069 //function : HasModificationsToDiscard
1070 //purpose  : Return true if the mesh has been edited since a total re-compute
1071 //           and those modifications may prevent successful partial re-compute.
1072 //           As a side effect reset _isModified flag if mesh is empty
1073 //issue    : 0020693
1074 //=======================================================================
1075
1076 bool SMESH_Mesh::HasModificationsToDiscard() const
1077 {
1078   if ( ! _isModified )
1079     return false;
1080
1081   // return true if the next Compute() will be partial and
1082   // existing but changed elements may prevent successful re-compute
1083   bool hasComputed = false, hasNotComputed = false;
1084   map <int, SMESH_subMesh*>::const_iterator i_sm = _mapSubMesh.begin();
1085   for ( ; i_sm != _mapSubMesh.end() ; ++i_sm )
1086     switch ( i_sm->second->GetSubShape().ShapeType() )
1087     {
1088     case TopAbs_EDGE:
1089     case TopAbs_FACE:
1090     case TopAbs_SOLID:
1091       if ( i_sm->second->IsMeshComputed() )
1092         hasComputed = true;
1093       else
1094         hasNotComputed = true;
1095       if ( hasComputed && hasNotComputed)
1096         return true;
1097     }
1098
1099   if ( NbNodes() < 1 )
1100     const_cast<SMESH_Mesh*>(this)->_isModified = false;
1101
1102   return false;
1103 }
1104
1105 //================================================================================
1106 /*!
1107  * \brief Check if any groups of the same type have equal names
1108  */
1109 //================================================================================
1110
1111 bool SMESH_Mesh::HasDuplicatedGroupNamesMED()
1112 {
1113   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1114   map< SMDSAbs_ElementType, set<string> > aGroupNames;
1115   for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1116   {
1117     SMESH_Group* aGroup = it->second;
1118     SMDSAbs_ElementType aType = aGroup->GetGroupDS()->GetType();
1119     string aGroupName = aGroup->GetName();
1120     aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1121     if (!aGroupNames[aType].insert(aGroupName).second)
1122       return true;
1123   }
1124
1125   return false;
1126 }
1127
1128 //================================================================================
1129 /*!
1130  * \brief Export the mesh to a med file
1131  */
1132 //================================================================================
1133
1134 void SMESH_Mesh::ExportMED(const char *        file, 
1135                            const char*         theMeshName, 
1136                            bool                theAutoGroups,
1137                            int                 theVersion,
1138                            const SMESHDS_Mesh* meshPart) 
1139   throw(SALOME_Exception)
1140 {
1141   Unexpect aCatch(SalomeException);
1142
1143   DriverMED_W_SMESHDS_Mesh myWriter;
1144   myWriter.SetFile    ( file, MED::EVersion(theVersion) );
1145   myWriter.SetMesh    ( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS   );
1146   if ( !theMeshName ) 
1147     myWriter.SetMeshId  ( _idDoc      );
1148   else {
1149     myWriter.SetMeshId  ( -1          );
1150     myWriter.SetMeshName( theMeshName );
1151   }
1152
1153   if ( theAutoGroups ) {
1154     myWriter.AddGroupOfNodes();
1155     myWriter.AddGroupOfEdges();
1156     myWriter.AddGroupOfFaces();
1157     myWriter.AddGroupOfVolumes();
1158   }
1159
1160   // Pass groups to writer. Provide unique group names.
1161   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1162   if ( !meshPart )
1163   {
1164     map< SMDSAbs_ElementType, set<string> > aGroupNames;
1165     char aString [256];
1166     int maxNbIter = 10000; // to guarantee cycle finish
1167     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1168       SMESH_Group*       aGroup   = it->second;
1169       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1170       if ( aGroupDS ) {
1171         SMDSAbs_ElementType aType = aGroupDS->GetType();
1172         string aGroupName0 = aGroup->GetName();
1173         aGroupName0.resize(MAX_MED_GROUP_NAME_LENGTH);
1174         string aGroupName = aGroupName0;
1175         for (int i = 1; !aGroupNames[aType].insert(aGroupName).second && i < maxNbIter; i++) {
1176           sprintf(&aString[0], "GR_%d_%s", i, aGroupName0.c_str());
1177           aGroupName = aString;
1178           aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1179         }
1180         aGroupDS->SetStoreName( aGroupName.c_str() );
1181         myWriter.AddGroup( aGroupDS );
1182       }
1183     }
1184   }
1185   // Perform export
1186   myWriter.Perform();
1187 }
1188
1189 void SMESH_Mesh::ExportSAUV(const char *file, 
1190                             const char* theMeshName, 
1191                             bool theAutoGroups)
1192   throw(SALOME_Exception)
1193 {
1194   std::string medfilename(file);
1195   medfilename += ".med";
1196   std::string cmd;
1197 #ifdef WNT
1198   cmd = "%PYTHONBIN% ";
1199 #else
1200   cmd = "python ";
1201 #endif
1202   cmd += "-c \"";
1203   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1204   cmd += "\"";
1205   system(cmd.c_str());
1206   ExportMED(medfilename.c_str(), theMeshName, theAutoGroups, 1);
1207 #ifdef WNT
1208   cmd = "%PYTHONBIN% ";
1209 #else
1210   cmd = "python ";
1211 #endif
1212   cmd += "-c \"";
1213   cmd += "from medutilities import convert ; convert(r'" + medfilename + "', 'MED', 'GIBI', 1, r'" + file + "')";
1214   cmd += "\"";
1215   system(cmd.c_str());
1216 #ifdef WNT
1217   cmd = "%PYTHONBIN% ";
1218 #else
1219   cmd = "python ";
1220 #endif
1221   cmd += "-c \"";
1222   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1223   cmd += "\"";
1224   system(cmd.c_str());
1225 }
1226
1227 //================================================================================
1228 /*!
1229  * \brief Export the mesh to a DAT file
1230  */
1231 //================================================================================
1232
1233 void SMESH_Mesh::ExportDAT(const char *        file,
1234                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1235 {
1236   Unexpect aCatch(SalomeException);
1237   DriverDAT_W_SMDS_Mesh myWriter;
1238   myWriter.SetFile( file );
1239   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1240   myWriter.SetMeshId(_idDoc);
1241   myWriter.Perform();
1242 }
1243
1244 //================================================================================
1245 /*!
1246  * \brief Export the mesh to an UNV file
1247  */
1248 //================================================================================
1249
1250 void SMESH_Mesh::ExportUNV(const char *        file,
1251                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1252 {
1253   Unexpect aCatch(SalomeException);
1254   DriverUNV_W_SMDS_Mesh myWriter;
1255   myWriter.SetFile( file );
1256   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1257   myWriter.SetMeshId(_idDoc);
1258   //  myWriter.SetGroups(_mapGroup);
1259
1260   if ( !meshPart )
1261   {
1262     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1263       SMESH_Group*       aGroup   = it->second;
1264       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1265       if ( aGroupDS ) {
1266         string aGroupName = aGroup->GetName();
1267         aGroupDS->SetStoreName( aGroupName.c_str() );
1268         myWriter.AddGroup( aGroupDS );
1269       }
1270     }
1271   }
1272   myWriter.Perform();
1273 }
1274
1275 //================================================================================
1276 /*!
1277  * \brief Export the mesh to an STL file
1278  */
1279 //================================================================================
1280
1281 void SMESH_Mesh::ExportSTL(const char *        file,
1282                            const bool          isascii,
1283                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1284 {
1285   Unexpect aCatch(SalomeException);
1286   DriverSTL_W_SMDS_Mesh myWriter;
1287   myWriter.SetFile( file );
1288   myWriter.SetIsAscii( isascii );
1289   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS);
1290   myWriter.SetMeshId(_idDoc);
1291   myWriter.Perform();
1292 }
1293
1294 //================================================================================
1295 /*!
1296  * \brief Export the mesh to the CGNS file
1297  */
1298 //================================================================================
1299
1300 void SMESH_Mesh::ExportCGNS(const char *        file,
1301                             const SMESHDS_Mesh* meshDS)
1302 {
1303   int res = Driver_Mesh::DRS_FAIL;
1304 #ifdef WITH_CGNS
1305   DriverCGNS_Write myWriter;
1306   myWriter.SetFile( file );
1307   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1308   myWriter.SetMeshName( SMESH_Comment("Mesh_") << meshDS->GetPersistentId());
1309   res = myWriter.Perform();
1310 #endif
1311   if ( res != Driver_Mesh::DRS_OK )
1312     throw SALOME_Exception("Export failed");
1313 }
1314
1315 //================================================================================
1316 /*!
1317  * \brief Return number of nodes in the mesh
1318  */
1319 //================================================================================
1320
1321 int SMESH_Mesh::NbNodes() const throw(SALOME_Exception)
1322 {
1323   Unexpect aCatch(SalomeException);
1324   return _myMeshDS->NbNodes();
1325 }
1326
1327 //================================================================================
1328 /*!
1329  * \brief  Return number of edges of given order in the mesh
1330  */
1331 //================================================================================
1332
1333 int SMESH_Mesh::Nb0DElements() const throw(SALOME_Exception)
1334 {
1335   Unexpect aCatch(SalomeException);
1336   return _myMeshDS->GetMeshInfo().Nb0DElements();
1337 }
1338
1339 //================================================================================
1340 /*!
1341  * \brief  Return number of edges of given order in the mesh
1342  */
1343 //================================================================================
1344
1345 int SMESH_Mesh::NbEdges(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1346 {
1347   Unexpect aCatch(SalomeException);
1348   return _myMeshDS->GetMeshInfo().NbEdges(order);
1349 }
1350
1351 //================================================================================
1352 /*!
1353  * \brief Return number of faces of given order in the mesh
1354  */
1355 //================================================================================
1356
1357 int SMESH_Mesh::NbFaces(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1358 {
1359   Unexpect aCatch(SalomeException);
1360   return _myMeshDS->GetMeshInfo().NbFaces(order);
1361 }
1362
1363 //================================================================================
1364 /*!
1365  * \brief Return the number of faces in the mesh
1366  */
1367 //================================================================================
1368
1369 int SMESH_Mesh::NbTriangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1370 {
1371   Unexpect aCatch(SalomeException);
1372   return _myMeshDS->GetMeshInfo().NbTriangles(order);
1373 }
1374
1375 //================================================================================
1376 /*!
1377  * \brief Return the number nodes faces in the mesh
1378  */
1379 //================================================================================
1380
1381 int SMESH_Mesh::NbQuadrangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1382 {
1383   Unexpect aCatch(SalomeException);
1384   return _myMeshDS->GetMeshInfo().NbQuadrangles(order);
1385 }
1386
1387 //================================================================================
1388 /*!
1389  * \brief Return the number of polygonal faces in the mesh
1390  */
1391 //================================================================================
1392
1393 int SMESH_Mesh::NbPolygons() const throw(SALOME_Exception)
1394 {
1395   Unexpect aCatch(SalomeException);
1396   return _myMeshDS->GetMeshInfo().NbPolygons();
1397 }
1398
1399 //================================================================================
1400 /*!
1401  * \brief Return number of volumes of given order in the mesh
1402  */
1403 //================================================================================
1404
1405 int SMESH_Mesh::NbVolumes(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1406 {
1407   Unexpect aCatch(SalomeException);
1408   return _myMeshDS->GetMeshInfo().NbVolumes(order);
1409 }
1410
1411 //================================================================================
1412 /*!
1413  * \brief  Return number of tetrahedrons of given order in the mesh
1414  */
1415 //================================================================================
1416
1417 int SMESH_Mesh::NbTetras(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1418 {
1419   Unexpect aCatch(SalomeException);
1420   return _myMeshDS->GetMeshInfo().NbTetras(order);
1421 }
1422
1423 //================================================================================
1424 /*!
1425  * \brief  Return number of hexahedrons of given order in the mesh
1426  */
1427 //================================================================================
1428
1429 int SMESH_Mesh::NbHexas(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1430 {
1431   Unexpect aCatch(SalomeException);
1432   return _myMeshDS->GetMeshInfo().NbHexas(order);
1433 }
1434
1435 //================================================================================
1436 /*!
1437  * \brief  Return number of pyramids of given order in the mesh
1438  */
1439 //================================================================================
1440
1441 int SMESH_Mesh::NbPyramids(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1442 {
1443   Unexpect aCatch(SalomeException);
1444   return _myMeshDS->GetMeshInfo().NbPyramids(order);
1445 }
1446
1447 //================================================================================
1448 /*!
1449  * \brief  Return number of prisms (penthahedrons) of given order in the mesh
1450  */
1451 //================================================================================
1452
1453 int SMESH_Mesh::NbPrisms(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1454 {
1455   Unexpect aCatch(SalomeException);
1456   return _myMeshDS->GetMeshInfo().NbPrisms(order);
1457 }
1458
1459 //================================================================================
1460 /*!
1461  * \brief  Return number of polyhedrons in the mesh
1462  */
1463 //================================================================================
1464
1465 int SMESH_Mesh::NbPolyhedrons() const throw(SALOME_Exception)
1466 {
1467   Unexpect aCatch(SalomeException);
1468   return _myMeshDS->GetMeshInfo().NbPolyhedrons();
1469 }
1470
1471 //================================================================================
1472 /*!
1473  * \brief  Return number of submeshes in the mesh
1474  */
1475 //================================================================================
1476
1477 int SMESH_Mesh::NbSubMesh() const throw(SALOME_Exception)
1478 {
1479   Unexpect aCatch(SalomeException);
1480   return _myMeshDS->NbSubMesh();
1481 }
1482
1483 //=======================================================================
1484 //function : IsNotConformAllowed
1485 //purpose  : check if a hypothesis alowing notconform mesh is present
1486 //=======================================================================
1487
1488 bool SMESH_Mesh::IsNotConformAllowed() const
1489 {
1490   if(MYDEBUG) MESSAGE("SMESH_Mesh::IsNotConformAllowed");
1491
1492   static SMESH_HypoFilter filter( SMESH_HypoFilter::HasName( "NotConformAllowed" ));
1493   return GetHypothesis( _myMeshDS->ShapeToMesh(), filter, false );
1494 }
1495
1496 //=======================================================================
1497 //function : IsMainShape
1498 //purpose  : 
1499 //=======================================================================
1500
1501 bool SMESH_Mesh::IsMainShape(const TopoDS_Shape& theShape) const
1502 {
1503   return theShape.IsSame(_myMeshDS->ShapeToMesh() );
1504 }
1505
1506 //=============================================================================
1507 /*!
1508  *  
1509  */
1510 //=============================================================================
1511
1512 SMESH_Group* SMESH_Mesh::AddGroup (const SMDSAbs_ElementType theType,
1513                                    const char*               theName,
1514                                    int&                      theId,
1515                                    const TopoDS_Shape&       theShape,
1516                                    const SMESH_PredicatePtr& thePredicate)
1517 {
1518   if (_mapGroup.count(_groupId))
1519     return NULL;
1520   theId = _groupId;
1521   SMESH_Group* aGroup = new SMESH_Group (theId, this, theType, theName, theShape, thePredicate);
1522   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1523   _mapGroup[_groupId++] = aGroup;
1524   return aGroup;
1525 }
1526
1527 //================================================================================
1528 /*!
1529  * \brief Creates SMESH_Groups for not wrapped SMESHDS_Groups
1530  *  \retval bool - true if new SMESH_Groups have been created
1531  * 
1532  */
1533 //================================================================================
1534
1535 bool SMESH_Mesh::SynchronizeGroups()
1536 {
1537   int nbGroups = _mapGroup.size();
1538   const set<SMESHDS_GroupBase*>& groups = _myMeshDS->GetGroups();
1539   set<SMESHDS_GroupBase*>::const_iterator gIt = groups.begin();
1540   for ( ; gIt != groups.end(); ++gIt )
1541   {
1542     SMESHDS_GroupBase* groupDS = (SMESHDS_GroupBase*) *gIt;
1543     _groupId = groupDS->GetID();
1544     if ( !_mapGroup.count( _groupId ))
1545       _mapGroup[_groupId] = new SMESH_Group( groupDS );
1546   }
1547   if ( !_mapGroup.empty() )
1548     _groupId = _mapGroup.rbegin()->first + 1;
1549
1550   return nbGroups < _mapGroup.size();
1551 }
1552
1553 //================================================================================
1554 /*!
1555  * \brief Return iterator on all existing groups
1556  */
1557 //================================================================================
1558
1559 SMESH_Mesh::GroupIteratorPtr SMESH_Mesh::GetGroups() const
1560 {
1561   typedef map <int, SMESH_Group *> TMap;
1562   return GroupIteratorPtr( new SMDS_mapIterator<TMap>( _mapGroup ));
1563 }
1564
1565 //=============================================================================
1566 /*!
1567  * \brief Return a group by ID
1568  */
1569 //=============================================================================
1570
1571 SMESH_Group* SMESH_Mesh::GetGroup (const int theGroupID)
1572 {
1573   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1574     return NULL;
1575   return _mapGroup[theGroupID];
1576 }
1577
1578
1579 //=============================================================================
1580 /*!
1581  * \brief Return IDs of all groups
1582  */
1583 //=============================================================================
1584
1585 list<int> SMESH_Mesh::GetGroupIds() const
1586 {
1587   list<int> anIds;
1588   for ( map<int, SMESH_Group*>::const_iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1589     anIds.push_back( it->first );
1590   
1591   return anIds;
1592 }
1593
1594 //================================================================================
1595 /*!
1596  * \brief Set a caller of RemoveGroup() at level of CORBA API implementation.
1597  * The set upCaller will be deleted by SMESH_Mesh
1598  */
1599 //================================================================================
1600
1601 void SMESH_Mesh::SetRemoveGroupCallUp( TRmGroupCallUp* upCaller )
1602 {
1603   if ( _rmGroupCallUp ) delete _rmGroupCallUp;
1604   _rmGroupCallUp = upCaller;
1605 }
1606
1607 //=============================================================================
1608 /*!
1609  *  
1610  */
1611 //=============================================================================
1612
1613 bool SMESH_Mesh::RemoveGroup (const int theGroupID)
1614 {
1615   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1616     return false;
1617   GetMeshDS()->RemoveGroup( _mapGroup[theGroupID]->GetGroupDS() );
1618   delete _mapGroup[theGroupID];
1619   _mapGroup.erase (theGroupID);
1620   if (_rmGroupCallUp)
1621     _rmGroupCallUp->RemoveGroup( theGroupID );
1622   return true;
1623 }
1624
1625 //=======================================================================
1626 //function : GetAncestors
1627 //purpose  : return list of ancestors of theSubShape in the order
1628 //           that lower dimention shapes come first.
1629 //=======================================================================
1630
1631 const TopTools_ListOfShape& SMESH_Mesh::GetAncestors(const TopoDS_Shape& theS) const
1632 {
1633   if ( _mapAncestors.Contains( theS ) )
1634     return _mapAncestors.FindFromKey( theS );
1635
1636   static TopTools_ListOfShape emptyList;
1637   return emptyList;
1638 }
1639
1640 //=======================================================================
1641 //function : Dump
1642 //purpose  : dumps contents of mesh to stream [ debug purposes ]
1643 //=======================================================================
1644
1645 ostream& SMESH_Mesh::Dump(ostream& save)
1646 {
1647   int clause = 0;
1648   save << "========================== Dump contents of mesh ==========================" << endl << endl;
1649   save << ++clause << ") Total number of nodes:   \t"    << NbNodes() << endl;
1650   save << ++clause << ") Total number of edges:   \t"    << NbEdges() << endl;
1651   save << ++clause << ") Total number of faces:   \t"    << NbFaces() << endl;
1652   save << ++clause << ") Total number of polygons:\t"    << NbPolygons() << endl;
1653   save << ++clause << ") Total number of volumes:\t"     << NbVolumes() << endl;
1654   save << ++clause << ") Total number of polyhedrons:\t" << NbPolyhedrons() << endl << endl;
1655   for ( int isQuadratic = 0; isQuadratic < 2; ++isQuadratic )
1656   {
1657     string orderStr = isQuadratic ? "quadratic" : "linear";
1658     SMDSAbs_ElementOrder order  = isQuadratic ? ORDER_QUADRATIC : ORDER_LINEAR;
1659
1660     save << ++clause << ") Total number of " << orderStr << " edges:\t" << NbEdges(order) << endl;
1661     save << ++clause << ") Total number of " << orderStr << " faces:\t" << NbFaces(order) << endl;
1662     if ( NbFaces(order) > 0 ) {
1663       int nb3 = NbTriangles(order);
1664       int nb4 = NbQuadrangles(order);
1665       save << clause << ".1) Number of " << orderStr << " triangles:  \t" << nb3 << endl;
1666       save << clause << ".2) Number of " << orderStr << " quadrangles:\t" << nb4 << endl;
1667       if ( nb3 + nb4 !=  NbFaces(order) ) {
1668         map<int,int> myFaceMap;
1669         SMDS_FaceIteratorPtr itFaces=_myMeshDS->facesIterator();
1670         while( itFaces->more( ) ) {
1671           int nbNodes = itFaces->next()->NbNodes();
1672           if ( myFaceMap.find( nbNodes ) == myFaceMap.end() )
1673             myFaceMap[ nbNodes ] = 0;
1674           myFaceMap[ nbNodes ] = myFaceMap[ nbNodes ] + 1;
1675         }
1676         save << clause << ".3) Faces in detail: " << endl;
1677         map <int,int>::iterator itF;
1678         for (itF = myFaceMap.begin(); itF != myFaceMap.end(); itF++)
1679           save << "--> nb nodes: " << itF->first << " - nb elemens:\t" << itF->second << endl;
1680       }
1681     }
1682     save << ++clause << ") Total number of " << orderStr << " volumes:\t" << NbVolumes(order) << endl;
1683     if ( NbVolumes(order) > 0 ) {
1684       int nb8 = NbHexas(order);
1685       int nb4 = NbTetras(order);
1686       int nb5 = NbPyramids(order);
1687       int nb6 = NbPrisms(order);
1688       save << clause << ".1) Number of " << orderStr << " hexahedrons:\t" << nb8 << endl;
1689       save << clause << ".2) Number of " << orderStr << " tetrahedrons:\t" << nb4 << endl;
1690       save << clause << ".3) Number of " << orderStr << " prisms:      \t" << nb6 << endl;
1691       save << clause << ".4) Number of " << orderStr << " pyramids:\t" << nb5 << endl;
1692       if ( nb8 + nb4 + nb5 + nb6 != NbVolumes(order) ) {
1693         map<int,int> myVolumesMap;
1694         SMDS_VolumeIteratorPtr itVolumes=_myMeshDS->volumesIterator();
1695         while( itVolumes->more( ) ) {
1696           int nbNodes = itVolumes->next()->NbNodes();
1697           if ( myVolumesMap.find( nbNodes ) == myVolumesMap.end() )
1698             myVolumesMap[ nbNodes ] = 0;
1699           myVolumesMap[ nbNodes ] = myVolumesMap[ nbNodes ] + 1;
1700         }
1701         save << clause << ".5) Volumes in detail: " << endl;
1702         map <int,int>::iterator itV;
1703         for (itV = myVolumesMap.begin(); itV != myVolumesMap.end(); itV++)
1704           save << "--> nb nodes: " << itV->first << " - nb elemens:\t" << itV->second << endl;
1705       }
1706     }
1707     save << endl;
1708   }
1709   save << "===========================================================================" << endl;
1710   return save;
1711 }
1712
1713 //=======================================================================
1714 //function : GetElementType
1715 //purpose  : Returns type of mesh element with certain id
1716 //=======================================================================
1717
1718 SMDSAbs_ElementType SMESH_Mesh::GetElementType( const int id, const bool iselem )
1719 {
1720   return _myMeshDS->GetElementType( id, iselem );
1721 }
1722
1723 //=============================================================================
1724 /*!
1725  *  \brief Convert group on geometry into standalone group
1726  */
1727 //=============================================================================
1728
1729 SMESH_Group* SMESH_Mesh::ConvertToStandalone ( int theGroupID )
1730 {
1731   SMESH_Group* aGroup = 0;
1732   map < int, SMESH_Group * >::iterator itg = _mapGroup.find( theGroupID );
1733   if ( itg == _mapGroup.end() )
1734     return aGroup;
1735
1736   SMESH_Group* anOldGrp = (*itg).second;
1737   SMESHDS_GroupBase* anOldGrpDS = anOldGrp->GetGroupDS();
1738   if ( !anOldGrp || !anOldGrpDS )
1739     return aGroup;
1740
1741   // create new standalone group
1742   aGroup = new SMESH_Group (theGroupID, this, anOldGrpDS->GetType(), anOldGrp->GetName() );
1743   _mapGroup[theGroupID] = aGroup;
1744
1745   SMESHDS_Group* aNewGrpDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1746   GetMeshDS()->RemoveGroup( anOldGrpDS );
1747   GetMeshDS()->AddGroup( aNewGrpDS );
1748
1749   // add elements (or nodes) into new created group
1750   SMDS_ElemIteratorPtr anItr = anOldGrpDS->GetElements();
1751   while ( anItr->more() )
1752     aNewGrpDS->Add( (anItr->next())->GetID() );
1753
1754   // remove old group
1755   delete anOldGrp;
1756
1757   return aGroup;
1758 }
1759
1760 //=============================================================================
1761 /*!
1762  *  \brief remove submesh order  from Mesh
1763  */
1764 //=============================================================================
1765
1766 void SMESH_Mesh::ClearMeshOrder()
1767 {
1768   _mySubMeshOrder.clear();
1769 }
1770
1771 //=============================================================================
1772 /*!
1773  *  \brief remove submesh order  from Mesh
1774  */
1775 //=============================================================================
1776
1777 void SMESH_Mesh::SetMeshOrder(const TListOfListOfInt& theOrder )
1778 {
1779   _mySubMeshOrder = theOrder;
1780 }
1781
1782 //=============================================================================
1783 /*!
1784  *  \brief return submesh order if any
1785  */
1786 //=============================================================================
1787
1788 const TListOfListOfInt& SMESH_Mesh::GetMeshOrder() const
1789 {
1790   return _mySubMeshOrder;
1791 }
1792
1793 //=============================================================================
1794 /*!
1795  *  \brief fill _mapAncestors
1796  */
1797 //=============================================================================
1798
1799 void SMESH_Mesh::fillAncestorsMap(const TopoDS_Shape& theShape)
1800 {
1801
1802   int desType, ancType;
1803   if ( !theShape.IsSame( GetShapeToMesh()) && theShape.ShapeType() == TopAbs_COMPOUND )
1804   {
1805     // a geom group is added. Insert it into lists of ancestors before
1806     // the first ancestor more complex than group members
1807     int memberType = TopoDS_Iterator( theShape ).Value().ShapeType();
1808     for ( desType = TopAbs_VERTEX; desType >= memberType; desType-- )
1809       for (TopExp_Explorer des( theShape, TopAbs_ShapeEnum( desType )); des.More(); des.Next())
1810       {
1811         if ( !_mapAncestors.Contains( des.Current() )) continue;// issue 0020982
1812         TopTools_ListOfShape& ancList = _mapAncestors.ChangeFromKey( des.Current() );
1813         TopTools_ListIteratorOfListOfShape ancIt (ancList);
1814         while ( ancIt.More() && ancIt.Value().ShapeType() >= memberType )
1815           ancIt.Next();
1816         if ( ancIt.More() )
1817           ancList.InsertBefore( theShape, ancIt );
1818       }
1819   }
1820   {
1821     for ( desType = TopAbs_VERTEX; desType > TopAbs_COMPOUND; desType-- )
1822       for ( ancType = desType - 1; ancType >= TopAbs_COMPOUND; ancType-- )
1823         TopExp::MapShapesAndAncestors ( theShape,
1824                                         (TopAbs_ShapeEnum) desType,
1825                                         (TopAbs_ShapeEnum) ancType,
1826                                         _mapAncestors );
1827   }
1828 }
1829
1830 //=============================================================================
1831 /*!
1832  * \brief sort submeshes according to stored mesh order
1833  * \param theListToSort in out list to be sorted
1834  * \return FALSE if nothing sorted
1835  */
1836 //=============================================================================
1837
1838 bool SMESH_Mesh::SortByMeshOrder(list<SMESH_subMesh*>& theListToSort) const
1839 {
1840   if ( !_mySubMeshOrder.size() || theListToSort.size() < 2)
1841     return true;
1842   
1843   bool res = false;
1844   list<SMESH_subMesh*> onlyOrderedList;
1845   // collect all ordered submeshes in one list as pointers
1846   // and get their positions within theListToSort
1847   typedef list<SMESH_subMesh*>::iterator TPosInList;
1848   map< int, TPosInList > sortedPos;
1849   TPosInList smBeg = theListToSort.begin(), smEnd = theListToSort.end();
1850   TListOfListOfInt::const_iterator listIddIt = _mySubMeshOrder.begin();
1851   for( ; listIddIt != _mySubMeshOrder.end(); listIddIt++) {
1852     const TListOfInt& listOfId = *listIddIt;
1853     TListOfInt::const_iterator idIt = listOfId.begin();
1854     for ( ; idIt != listOfId.end(); idIt++ ) {
1855       if ( SMESH_subMesh * sm = GetSubMeshContaining( *idIt )) {
1856         TPosInList smPos = find( smBeg, smEnd, sm );
1857         if ( smPos != smEnd ) {
1858           onlyOrderedList.push_back( sm );
1859           sortedPos[ distance( smBeg, smPos )] = smPos;
1860         }
1861       }
1862     }
1863   }
1864   if (onlyOrderedList.size() < 2)
1865     return res;
1866   res = true;
1867
1868   list<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
1869   list<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
1870
1871   // iterate on ordered submeshes and insert them in detected positions
1872   map< int, TPosInList >::iterator i_pos = sortedPos.begin();
1873   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
1874     *(i_pos->second) = *onlyBIt;
1875
1876   return res;
1877 }
1878
1879 //=============================================================================
1880 /*!
1881  * \brief sort submeshes according to stored mesh order
1882  * \param theListToSort in out list to be sorted
1883  * \return FALSE if nothing sorted
1884  */
1885 //=============================================================================
1886
1887 list<SMESH_subMesh*> SMESH_Mesh::getAncestorsSubMeshes
1888   (const TopoDS_Shape& theSubShape) const
1889 {
1890   list<SMESH_subMesh*> listOfSubMesh;
1891   TopTools_ListIteratorOfListOfShape it( GetAncestors( theSubShape ));
1892   for (; it.More(); it.Next() )
1893     if ( SMESH_subMesh* sm = GetSubMeshContaining( it.Value() ))
1894       listOfSubMesh.push_back(sm);
1895
1896   // sort submeshes according to stored mesh order
1897   SortByMeshOrder( listOfSubMesh );
1898
1899   return listOfSubMesh;
1900 }