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