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