Salome HOME
Merge from V6_main (04/10/2012)
[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 //  File   : SMESH_Mesh.cxx
24 //  Author : Paul RASCLE, EDF
25 //  Module : SMESH
26 //
27 #include "SMESH_Mesh.hxx"
28 #include "SMESH_MesherHelper.hxx"
29 #include "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 "DriverDAT_W_SMDS_Mesh.h"
44 #include "DriverGMF_Read.hxx"
45 #include "DriverGMF_Write.hxx"
46 #include "DriverMED_R_SMESHDS_Mesh.h"
47 #include "DriverMED_W_SMESHDS_Mesh.h"
48 #include "DriverSTL_R_SMDS_Mesh.h"
49 #include "DriverSTL_W_SMDS_Mesh.h"
50 #include "DriverUNV_R_SMDS_Mesh.h"
51 #include "DriverUNV_W_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  * \brief Fill its data by reading a GMF file
545  */
546 //================================================================================
547
548 SMESH_ComputeErrorPtr SMESH_Mesh::GMFToMesh(const char* theFileName)
549 {
550   DriverGMF_Read myReader;
551   myReader.SetMesh(_myMeshDS);
552   myReader.SetFile(theFileName);
553   myReader.Perform();
554   //theMeshName = myReader.GetMeshName();
555
556   // create groups
557   SynchronizeGroups();
558
559   return myReader.GetError();
560 }
561
562 //=============================================================================
563 /*!
564  * 
565  */
566 //=============================================================================
567
568 SMESH_Hypothesis::Hypothesis_Status
569   SMESH_Mesh::AddHypothesis(const TopoDS_Shape & aSubShape,
570                             int                  anHypId  ) throw(SALOME_Exception)
571 {
572   Unexpect aCatch(SalomeException);
573   if(MYDEBUG) MESSAGE("SMESH_Mesh::AddHypothesis");
574
575   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
576   if ( !subMesh || !subMesh->GetId())
577     return SMESH_Hypothesis::HYP_BAD_SUBSHAPE;
578
579   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
580   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
581   {
582     if(MYDEBUG) MESSAGE("Hypothesis ID does not give an hypothesis");
583     if(MYDEBUG) {
584       SCRUTE(_studyId);
585       SCRUTE(anHypId);
586     }
587     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
588   }
589
590   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
591   MESSAGE( "SMESH_Mesh::AddHypothesis " << anHyp->GetName() );
592
593   bool isGlobalHyp = IsMainShape( aSubShape );
594
595   // NotConformAllowed can be only global
596   if ( !isGlobalHyp )
597   {
598     // NOTE: this is not a correct way to check a name of hypothesis,
599     // there should be an attribute of hypothesis saying that it can/can't
600     // be global/local
601     string hypName = anHyp->GetName();
602     if ( hypName == "NotConformAllowed" )
603     {
604       if(MYDEBUG) MESSAGE( "Hypotesis <NotConformAllowed> can be only global" );
605       return SMESH_Hypothesis::HYP_INCOMPATIBLE;
606     }
607   }
608
609   // shape 
610
611   bool isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
612   int event = isAlgo ? SMESH_subMesh::ADD_ALGO : SMESH_subMesh::ADD_HYP;
613
614   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
615
616   // sub-shapes
617   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
618       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is added on father
619   {
620     event = isAlgo ? SMESH_subMesh::ADD_FATHER_ALGO : SMESH_subMesh::ADD_FATHER_HYP;
621
622     SMESH_Hypothesis::Hypothesis_Status ret2 =
623       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
624     if (ret2 > ret)
625       ret = ret2;
626
627     // check concurent hypotheses on ancestors
628     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !isGlobalHyp )
629     {
630       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
631       while ( smIt->more() ) {
632         SMESH_subMesh* sm = smIt->next();
633         if ( sm->IsApplicableHypotesis( anHyp )) {
634           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
635           if (ret2 > ret) {
636             ret = ret2;
637             break;
638           }
639         }
640       }
641     }
642   }
643   HasModificationsToDiscard(); // to reset _isModified flag if a mesh becomes empty
644
645   GetMeshDS()->Modified();
646
647   if(MYDEBUG) subMesh->DumpAlgoState(true);
648   if(MYDEBUG) SCRUTE(ret);
649   return ret;
650 }
651
652 //=============================================================================
653 /*!
654  * 
655  */
656 //=============================================================================
657
658 SMESH_Hypothesis::Hypothesis_Status
659   SMESH_Mesh::RemoveHypothesis(const TopoDS_Shape & aSubShape,
660                                int anHypId)throw(SALOME_Exception)
661 {
662   Unexpect aCatch(SalomeException);
663   if(MYDEBUG) MESSAGE("SMESH_Mesh::RemoveHypothesis");
664   
665   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
666   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
667     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
668   
669   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
670   if(MYDEBUG) {
671     int hypType = anHyp->GetType();
672     SCRUTE(hypType);
673   }
674   
675   // shape 
676   
677   bool isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
678   int event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
679
680   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
681
682   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
683
684   // there may appear concurrent hyps that were covered by the removed hyp
685   if (ret < SMESH_Hypothesis::HYP_CONCURENT &&
686       subMesh->IsApplicableHypotesis( anHyp ) &&
687       subMesh->CheckConcurentHypothesis( anHyp->GetType() ) != SMESH_Hypothesis::HYP_OK)
688     ret = SMESH_Hypothesis::HYP_CONCURENT;
689
690   // sub-shapes
691   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
692       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is removed from father
693   {
694     event = isAlgo ? SMESH_subMesh::REMOVE_FATHER_ALGO : SMESH_subMesh::REMOVE_FATHER_HYP;
695
696     SMESH_Hypothesis::Hypothesis_Status ret2 =
697       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
698     if (ret2 > ret) // more severe
699       ret = ret2;
700
701     // check concurent hypotheses on ancestors
702     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !IsMainShape( aSubShape ) )
703     {
704       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
705       while ( smIt->more() ) {
706         SMESH_subMesh* sm = smIt->next();
707         if ( sm->IsApplicableHypotesis( anHyp )) {
708           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
709           if (ret2 > ret) {
710             ret = ret2;
711             break;
712           }
713         }
714       }
715     }
716   }
717
718   HasModificationsToDiscard(); // to reset _isModified flag if mesh become empty
719
720   GetMeshDS()->Modified();
721
722   if(MYDEBUG) subMesh->DumpAlgoState(true);
723   if(MYDEBUG) SCRUTE(ret);
724   return ret;
725 }
726
727 //=============================================================================
728 /*!
729  * 
730  */
731 //=============================================================================
732
733 const list<const SMESHDS_Hypothesis*>&
734 SMESH_Mesh::GetHypothesisList(const TopoDS_Shape & aSubShape) const
735   throw(SALOME_Exception)
736 {
737   Unexpect aCatch(SalomeException);
738   return _myMeshDS->GetHypothesis(aSubShape);
739 }
740
741 //=======================================================================
742 /*!
743  * \brief Return the hypothesis assigned to the shape
744  *  \param aSubShape    - the shape to check
745  *  \param aFilter      - the hypothesis filter
746  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
747  *  \param assignedTo   - to return the shape the found hypo is assigned to
748  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
749  */
750 //=======================================================================
751
752 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const TopoDS_Shape &    aSubShape,
753                                                    const SMESH_HypoFilter& aFilter,
754                                                    const bool              andAncestors,
755                                                    TopoDS_Shape*           assignedTo) const
756 {
757   {
758     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
759     list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
760     for ( ; hyp != hypList.end(); hyp++ ) {
761       const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
762       if ( aFilter.IsOk( h, aSubShape)) {
763         if ( assignedTo ) *assignedTo = aSubShape;
764         return h;
765       }
766     }
767   }
768   if ( andAncestors )
769   {
770     // user sorted submeshes of ancestors, according to stored submesh priority
771     const list<SMESH_subMesh*> smList = getAncestorsSubMeshes( aSubShape );
772     list<SMESH_subMesh*>::const_iterator smIt = smList.begin(); 
773     for ( ; smIt != smList.end(); smIt++ )
774     {
775       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
776       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
777       list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
778       for ( ; hyp != hypList.end(); hyp++ ) {
779         const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
780         if (aFilter.IsOk( h, curSh )) {
781           if ( assignedTo ) *assignedTo = curSh;
782           return h;
783         }
784       }
785     }
786   }
787   return 0;
788 }
789
790 //================================================================================
791 /*!
792  * \brief Return hypothesis assigned to the shape
793   * \param aSubShape - the shape to check
794   * \param aFilter - the hypothesis filter
795   * \param aHypList - the list of the found hypotheses
796   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
797   * \retval int - number of unique hypos in aHypList
798  */
799 //================================================================================
800
801 int SMESH_Mesh::GetHypotheses(const TopoDS_Shape &                aSubShape,
802                               const SMESH_HypoFilter&             aFilter,
803                               list <const SMESHDS_Hypothesis * >& aHypList,
804                               const bool                          andAncestors) const
805 {
806   set<string> hypTypes; // to exclude same type hypos from the result list
807   int nbHyps = 0;
808
809   // only one main hypothesis is allowed
810   bool mainHypFound = false;
811
812   // fill in hypTypes
813   list<const SMESHDS_Hypothesis*>::const_iterator hyp;
814   for ( hyp = aHypList.begin(); hyp != aHypList.end(); hyp++ ) {
815     if ( hypTypes.insert( (*hyp)->GetName() ).second )
816       nbHyps++;
817     if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
818       mainHypFound = true;
819   }
820
821   // get hypos from aSubShape
822   {
823     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
824     for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
825       if ( aFilter.IsOk (cSMESH_Hyp( *hyp ), aSubShape) &&
826            ( cSMESH_Hyp(*hyp)->IsAuxiliary() || !mainHypFound ) &&
827            hypTypes.insert( (*hyp)->GetName() ).second )
828       {
829         aHypList.push_back( *hyp );
830         nbHyps++;
831         if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
832           mainHypFound = true;
833       }
834   }
835
836   // get hypos from ancestors of aSubShape
837   if ( andAncestors )
838   {
839     TopTools_MapOfShape map;
840
841     // user sorted submeshes of ancestors, according to stored submesh priority
842     const list<SMESH_subMesh*> smList = getAncestorsSubMeshes( aSubShape );
843     list<SMESH_subMesh*>::const_iterator smIt = smList.begin(); 
844     for ( ; smIt != smList.end(); smIt++ )
845     {
846       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
847      if ( !map.Add( curSh ))
848         continue;
849       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
850       for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
851         if (aFilter.IsOk( cSMESH_Hyp( *hyp ), curSh ) &&
852             ( cSMESH_Hyp(*hyp)->IsAuxiliary() || !mainHypFound ) &&
853             hypTypes.insert( (*hyp)->GetName() ).second )
854         {
855           aHypList.push_back( *hyp );
856           nbHyps++;
857           if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
858             mainHypFound = true;
859         }
860     }
861   }
862   return nbHyps;
863 }
864
865 //=============================================================================
866 /*!
867  * 
868  */
869 //=============================================================================
870
871 const list<SMESHDS_Command*> & SMESH_Mesh::GetLog() throw(SALOME_Exception)
872 {
873   Unexpect aCatch(SalomeException);
874   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetLog");
875   return _myMeshDS->GetScript()->GetCommands();
876 }
877
878 //=============================================================================
879 /*!
880  * 
881  */
882 //=============================================================================
883 void SMESH_Mesh::ClearLog() throw(SALOME_Exception)
884 {
885   Unexpect aCatch(SalomeException);
886   if(MYDEBUG) MESSAGE("SMESH_Mesh::ClearLog");
887   _myMeshDS->GetScript()->Clear();
888 }
889
890 //=============================================================================
891 /*!
892  * Get or Create the SMESH_subMesh object implementation
893  */
894 //=============================================================================
895
896 SMESH_subMesh *SMESH_Mesh::GetSubMesh(const TopoDS_Shape & aSubShape)
897   throw(SALOME_Exception)
898 {
899   Unexpect aCatch(SalomeException);
900   SMESH_subMesh *aSubMesh;
901   int index = _myMeshDS->ShapeToIndex(aSubShape);
902
903   // for submeshes on GEOM Group
904   if (( !index || index > _nbSubShapes ) && aSubShape.ShapeType() == TopAbs_COMPOUND ) {
905     TopoDS_Iterator it( aSubShape );
906     if ( it.More() )
907     {
908       index = _myMeshDS->AddCompoundSubmesh( aSubShape, it.Value().ShapeType() );
909       if ( index > _nbSubShapes ) _nbSubShapes = index; // not to create sm for this group again
910
911       // fill map of Ancestors
912       fillAncestorsMap(aSubShape);
913     }
914   }
915 //   if ( !index )
916 //     return NULL; // neither sub-shape nor a group
917
918   map <int, SMESH_subMesh *>::iterator i_sm = _mapSubMesh.find(index);
919   if ( i_sm != _mapSubMesh.end())
920   {
921     aSubMesh = i_sm->second;
922   }
923   else
924   {
925     aSubMesh = new SMESH_subMesh(index, this, _myMeshDS, aSubShape);
926     _mapSubMesh[index] = aSubMesh;
927   }
928   return aSubMesh;
929 }
930
931 //=============================================================================
932 /*!
933  * Get the SMESH_subMesh object implementation. Dont create it, return null
934  * if it does not exist.
935  */
936 //=============================================================================
937
938 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const TopoDS_Shape & aSubShape) const
939   throw(SALOME_Exception)
940 {
941   Unexpect aCatch(SalomeException);
942   SMESH_subMesh *aSubMesh = NULL;
943   
944   int index = _myMeshDS->ShapeToIndex(aSubShape);
945
946   map <int, SMESH_subMesh *>::const_iterator i_sm = _mapSubMesh.find(index);
947   if ( i_sm != _mapSubMesh.end())
948     aSubMesh = i_sm->second;
949
950   return aSubMesh;
951 }
952 //=============================================================================
953 /*!
954  * Get the SMESH_subMesh object implementation. Dont create it, return null
955  * if it does not exist.
956  */
957 //=============================================================================
958
959 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const int aShapeID) const
960 throw(SALOME_Exception)
961 {
962   Unexpect aCatch(SalomeException);
963   
964   map <int, SMESH_subMesh *>::const_iterator i_sm = _mapSubMesh.find(aShapeID);
965   if (i_sm == _mapSubMesh.end())
966     return NULL;
967   return i_sm->second;
968 }
969 //================================================================================
970 /*!
971  * \brief Return submeshes of groups containing the given sub-shape
972  */
973 //================================================================================
974
975 list<SMESH_subMesh*>
976 SMESH_Mesh::GetGroupSubMeshesContaining(const TopoDS_Shape & aSubShape) const
977   throw(SALOME_Exception)
978 {
979   Unexpect aCatch(SalomeException);
980   list<SMESH_subMesh*> found;
981
982   SMESH_subMesh * subMesh = GetSubMeshContaining(aSubShape);
983   if ( !subMesh )
984     return found;
985
986   // submeshes of groups have max IDs, so search from the map end
987   map<int, SMESH_subMesh *>::const_reverse_iterator i_sm;
988   for ( i_sm = _mapSubMesh.rbegin(); i_sm != _mapSubMesh.rend(); ++i_sm) {
989     SMESHDS_SubMesh * ds = i_sm->second->GetSubMeshDS();
990     if ( ds && ds->IsComplexSubmesh() ) {
991       if ( SMESH_MesherHelper::IsSubShape( aSubShape, i_sm->second->GetSubShape() ))
992       {
993         found.push_back( i_sm->second );
994         //break;
995       }
996     } else {
997       break; // the rest sub-meshes are not those of groups
998     }
999   }
1000
1001   if ( found.empty() ) // maybe the main shape is a COMPOUND (issue 0021530)
1002   {
1003     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1))
1004       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1005       {
1006         TopoDS_Iterator it( mainSM->GetSubShape() );
1007         if ( it.Value().ShapeType() == aSubShape.ShapeType() &&
1008              SMESH_MesherHelper::IsSubShape( aSubShape, mainSM->GetSubShape() ))
1009           found.push_back( mainSM );
1010       }
1011   }
1012   return found;
1013 }
1014 //=======================================================================
1015 //function : IsUsedHypothesis
1016 //purpose  : Return True if anHyp is used to mesh aSubShape
1017 //=======================================================================
1018
1019 bool SMESH_Mesh::IsUsedHypothesis(SMESHDS_Hypothesis * anHyp,
1020                                   const SMESH_subMesh* aSubMesh)
1021 {
1022   SMESH_Hypothesis* hyp = static_cast<SMESH_Hypothesis*>(anHyp);
1023
1024   // check if anHyp can be used to mesh aSubMesh
1025   if ( !aSubMesh || !aSubMesh->IsApplicableHypotesis( hyp ))
1026     return false;
1027
1028   const TopoDS_Shape & aSubShape = const_cast<SMESH_subMesh*>( aSubMesh )->GetSubShape();
1029
1030   SMESH_Algo *algo = _gen->GetAlgo(*this, aSubShape );
1031
1032   // algorithm
1033   if (anHyp->GetType() > SMESHDS_Hypothesis::PARAM_ALGO)
1034     return ( anHyp == algo );
1035
1036   // algorithm parameter
1037   if (algo)
1038   {
1039     // look trough hypotheses used by algo
1040     SMESH_HypoFilter hypoKind;
1041     if ( algo->InitCompatibleHypoFilter( hypoKind, !hyp->IsAuxiliary() )) {
1042       list <const SMESHDS_Hypothesis * > usedHyps;
1043       if ( GetHypotheses( aSubShape, hypoKind, usedHyps, true ))
1044         return ( find( usedHyps.begin(), usedHyps.end(), anHyp ) != usedHyps.end() );
1045     }
1046   }
1047
1048   // look through all assigned hypotheses
1049   //SMESH_HypoFilter filter( SMESH_HypoFilter::Is( hyp ));
1050   return false; //GetHypothesis( aSubShape, filter, true );
1051 }
1052
1053 //=============================================================================
1054 /*!
1055  *
1056  */
1057 //=============================================================================
1058
1059 const list < SMESH_subMesh * >&
1060 SMESH_Mesh::GetSubMeshUsingHypothesis(SMESHDS_Hypothesis * anHyp)
1061   throw(SALOME_Exception)
1062 {
1063   Unexpect aCatch(SalomeException);
1064   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetSubMeshUsingHypothesis");
1065   map < int, SMESH_subMesh * >::iterator itsm;
1066   _subMeshesUsingHypothesisList.clear();
1067   for (itsm = _mapSubMesh.begin(); itsm != _mapSubMesh.end(); itsm++)
1068   {
1069     SMESH_subMesh *aSubMesh = (*itsm).second;
1070     if ( IsUsedHypothesis ( anHyp, aSubMesh ))
1071       _subMeshesUsingHypothesisList.push_back(aSubMesh);
1072   }
1073   return _subMeshesUsingHypothesisList;
1074 }
1075
1076 //=======================================================================
1077 //function : NotifySubMeshesHypothesisModification
1078 //purpose  : Say all submeshes using theChangedHyp that it has been modified
1079 //=======================================================================
1080
1081 void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* hyp)
1082 {
1083   Unexpect aCatch(SalomeException);
1084
1085   if ( !GetMeshDS()->IsUsedHypothesis( hyp ))
1086     return;
1087
1088   if (_callUp)
1089     _callUp->HypothesisModified();
1090
1091   const SMESH_Algo *foundAlgo = 0;
1092   SMESH_HypoFilter algoKind, compatibleHypoKind;
1093   list <const SMESHDS_Hypothesis * > usedHyps;
1094
1095
1096   map < int, SMESH_subMesh * >::iterator itsm;
1097   for (itsm = _mapSubMesh.begin(); itsm != _mapSubMesh.end(); itsm++)
1098   {
1099     SMESH_subMesh *aSubMesh = (*itsm).second;
1100     if ( aSubMesh->IsApplicableHypotesis( hyp ))
1101     {
1102       const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
1103
1104       if ( !foundAlgo ) // init filter for algo search
1105         algoKind.Init( THypType::IsAlgo() ).And( THypType::IsApplicableTo( aSubShape ));
1106       
1107       const SMESH_Algo *algo = static_cast<const SMESH_Algo*>
1108         ( GetHypothesis( aSubShape, algoKind, true ));
1109
1110       if ( algo )
1111       {
1112         bool sameAlgo = ( algo == foundAlgo );
1113         if ( !sameAlgo && foundAlgo )
1114           sameAlgo = ( strcmp( algo->GetName(), foundAlgo->GetName() ) == 0);
1115
1116         if ( !sameAlgo ) { // init filter for used hypos search
1117           if ( !algo->InitCompatibleHypoFilter( compatibleHypoKind, !hyp->IsAuxiliary() ))
1118             continue; // algo does not use any hypothesis
1119           foundAlgo = algo;
1120         }
1121
1122         // check if hyp is used by algo
1123         usedHyps.clear();
1124         if ( GetHypotheses( aSubShape, compatibleHypoKind, usedHyps, true ) &&
1125              find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() )
1126         {
1127           aSubMesh->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
1128                                     const_cast< SMESH_Hypothesis*>( hyp ));
1129         }
1130       }
1131     }
1132   }
1133   HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
1134   GetMeshDS()->Modified();
1135 }
1136
1137 //=============================================================================
1138 /*!
1139  *  Auto color functionality
1140  */
1141 //=============================================================================
1142 void SMESH_Mesh::SetAutoColor(bool theAutoColor) throw(SALOME_Exception)
1143 {
1144   Unexpect aCatch(SalomeException);
1145   _isAutoColor = theAutoColor;
1146 }
1147
1148 bool SMESH_Mesh::GetAutoColor() throw(SALOME_Exception)
1149 {
1150   Unexpect aCatch(SalomeException);
1151   return _isAutoColor;
1152 }
1153
1154 //=======================================================================
1155 //function : SetIsModified
1156 //purpose  : Set the flag meaning that the mesh has been edited "manually"
1157 //=======================================================================
1158
1159 void SMESH_Mesh::SetIsModified(bool isModified)
1160 {
1161   _isModified = isModified;
1162
1163   if ( _isModified )
1164     // check if mesh becomes empty as result of modification
1165     HasModificationsToDiscard();
1166 }
1167
1168 //=======================================================================
1169 //function : HasModificationsToDiscard
1170 //purpose  : Return true if the mesh has been edited since a total re-compute
1171 //           and those modifications may prevent successful partial re-compute.
1172 //           As a side effect reset _isModified flag if mesh is empty
1173 //issue    : 0020693
1174 //=======================================================================
1175
1176 bool SMESH_Mesh::HasModificationsToDiscard() const
1177 {
1178   if ( ! _isModified )
1179     return false;
1180
1181   // return true if the next Compute() will be partial and
1182   // existing but changed elements may prevent successful re-compute
1183   bool hasComputed = false, hasNotComputed = false;
1184   map <int, SMESH_subMesh*>::const_iterator i_sm = _mapSubMesh.begin();
1185   for ( ; i_sm != _mapSubMesh.end() ; ++i_sm )
1186     switch ( i_sm->second->GetSubShape().ShapeType() )
1187     {
1188     case TopAbs_EDGE:
1189     case TopAbs_FACE:
1190     case TopAbs_SOLID:
1191       if ( i_sm->second->IsMeshComputed() )
1192         hasComputed = true;
1193       else
1194         hasNotComputed = true;
1195       if ( hasComputed && hasNotComputed)
1196         return true;
1197     }
1198
1199   if ( NbNodes() < 1 )
1200     const_cast<SMESH_Mesh*>(this)->_isModified = false;
1201
1202   return false;
1203 }
1204
1205 //================================================================================
1206 /*!
1207  * \brief Check if any groups of the same type have equal names
1208  */
1209 //================================================================================
1210
1211 bool SMESH_Mesh::HasDuplicatedGroupNamesMED()
1212 {
1213   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1214   map< SMDSAbs_ElementType, set<string> > aGroupNames;
1215   for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1216   {
1217     SMESH_Group* aGroup = it->second;
1218     SMDSAbs_ElementType aType = aGroup->GetGroupDS()->GetType();
1219     string aGroupName = aGroup->GetName();
1220     aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1221     if (!aGroupNames[aType].insert(aGroupName).second)
1222       return true;
1223   }
1224
1225   return false;
1226 }
1227
1228 //================================================================================
1229 /*!
1230  * \brief Export the mesh to a med file
1231  */
1232 //================================================================================
1233
1234 void SMESH_Mesh::ExportMED(const char *        file, 
1235                            const char*         theMeshName, 
1236                            bool                theAutoGroups,
1237                            int                 theVersion,
1238                            const SMESHDS_Mesh* meshPart) 
1239   throw(SALOME_Exception)
1240 {
1241   Unexpect aCatch(SalomeException);
1242
1243   DriverMED_W_SMESHDS_Mesh myWriter;
1244   myWriter.SetFile    ( file, MED::EVersion(theVersion) );
1245   myWriter.SetMesh    ( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS   );
1246   if ( !theMeshName ) 
1247     myWriter.SetMeshId  ( _idDoc      );
1248   else {
1249     myWriter.SetMeshId  ( -1          );
1250     myWriter.SetMeshName( theMeshName );
1251   }
1252
1253   if ( theAutoGroups ) {
1254     myWriter.AddGroupOfNodes();
1255     myWriter.AddGroupOfEdges();
1256     myWriter.AddGroupOfFaces();
1257     myWriter.AddGroupOfVolumes();
1258   }
1259
1260   // Pass groups to writer. Provide unique group names.
1261   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1262   if ( !meshPart )
1263   {
1264     map< SMDSAbs_ElementType, set<string> > aGroupNames;
1265     char aString [256];
1266     int maxNbIter = 10000; // to guarantee cycle finish
1267     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1268       SMESH_Group*       aGroup   = it->second;
1269       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1270       if ( aGroupDS ) {
1271         SMDSAbs_ElementType aType = aGroupDS->GetType();
1272         string aGroupName0 = aGroup->GetName();
1273         aGroupName0.resize(MAX_MED_GROUP_NAME_LENGTH);
1274         string aGroupName = aGroupName0;
1275         for (int i = 1; !aGroupNames[aType].insert(aGroupName).second && i < maxNbIter; i++) {
1276           sprintf(&aString[0], "GR_%d_%s", i, aGroupName0.c_str());
1277           aGroupName = aString;
1278           aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1279         }
1280         aGroupDS->SetStoreName( aGroupName.c_str() );
1281         myWriter.AddGroup( aGroupDS );
1282       }
1283     }
1284   }
1285   // Perform export
1286   myWriter.Perform();
1287 }
1288
1289 void SMESH_Mesh::ExportSAUV(const char *file, 
1290                             const char* theMeshName, 
1291                             bool theAutoGroups)
1292   throw(SALOME_Exception)
1293 {
1294   std::string medfilename(file);
1295   medfilename += ".med";
1296   std::string cmd;
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   ExportMED(medfilename.c_str(), theMeshName, theAutoGroups, 1);
1307 #ifdef WNT
1308   cmd = "%PYTHONBIN% ";
1309 #else
1310   cmd = "python ";
1311 #endif
1312   cmd += "-c \"";
1313   cmd += "from medutilities import convert ; convert(r'" + medfilename + "', 'MED', 'GIBI', 1, r'" + file + "')";
1314   cmd += "\"";
1315   system(cmd.c_str());
1316 #ifdef WNT
1317   cmd = "%PYTHONBIN% ";
1318 #else
1319   cmd = "python ";
1320 #endif
1321   cmd += "-c \"";
1322   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1323   cmd += "\"";
1324   system(cmd.c_str());
1325 }
1326
1327 //================================================================================
1328 /*!
1329  * \brief Export the mesh to a DAT file
1330  */
1331 //================================================================================
1332
1333 void SMESH_Mesh::ExportDAT(const char *        file,
1334                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1335 {
1336   Unexpect aCatch(SalomeException);
1337   DriverDAT_W_SMDS_Mesh myWriter;
1338   myWriter.SetFile( file );
1339   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1340   myWriter.SetMeshId(_idDoc);
1341   myWriter.Perform();
1342 }
1343
1344 //================================================================================
1345 /*!
1346  * \brief Export the mesh to an UNV file
1347  */
1348 //================================================================================
1349
1350 void SMESH_Mesh::ExportUNV(const char *        file,
1351                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1352 {
1353   Unexpect aCatch(SalomeException);
1354   DriverUNV_W_SMDS_Mesh myWriter;
1355   myWriter.SetFile( file );
1356   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1357   myWriter.SetMeshId(_idDoc);
1358   //  myWriter.SetGroups(_mapGroup);
1359
1360   if ( !meshPart )
1361   {
1362     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1363       SMESH_Group*       aGroup   = it->second;
1364       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1365       if ( aGroupDS ) {
1366         string aGroupName = aGroup->GetName();
1367         aGroupDS->SetStoreName( aGroupName.c_str() );
1368         myWriter.AddGroup( aGroupDS );
1369       }
1370     }
1371   }
1372   myWriter.Perform();
1373 }
1374
1375 //================================================================================
1376 /*!
1377  * \brief Export the mesh to an STL file
1378  */
1379 //================================================================================
1380
1381 void SMESH_Mesh::ExportSTL(const char *        file,
1382                            const bool          isascii,
1383                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1384 {
1385   Unexpect aCatch(SalomeException);
1386   DriverSTL_W_SMDS_Mesh myWriter;
1387   myWriter.SetFile( file );
1388   myWriter.SetIsAscii( isascii );
1389   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS);
1390   myWriter.SetMeshId(_idDoc);
1391   myWriter.Perform();
1392 }
1393
1394 //================================================================================
1395 /*!
1396  * \brief Export the mesh to the CGNS file
1397  */
1398 //================================================================================
1399
1400 void SMESH_Mesh::ExportCGNS(const char *        file,
1401                             const SMESHDS_Mesh* meshDS)
1402 {
1403   int res = Driver_Mesh::DRS_FAIL;
1404 #ifdef WITH_CGNS
1405   DriverCGNS_Write myWriter;
1406   myWriter.SetFile( file );
1407   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1408   myWriter.SetMeshName( SMESH_Comment("Mesh_") << meshDS->GetPersistentId());
1409   res = myWriter.Perform();
1410 #endif
1411   if ( res != Driver_Mesh::DRS_OK )
1412     throw SALOME_Exception("Export failed");
1413 }
1414
1415 //================================================================================
1416 /*!
1417  * \brief Export the mesh to a GMF file
1418  */
1419 //================================================================================
1420
1421 void SMESH_Mesh::ExportGMF(const char *        file,
1422                            const SMESHDS_Mesh* meshDS)
1423 {
1424   DriverGMF_Write myWriter;
1425   myWriter.SetFile( file );
1426   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1427   myWriter.Perform();
1428 }
1429
1430 //================================================================================
1431 /*!
1432  * \brief Return number of nodes in the mesh
1433  */
1434 //================================================================================
1435
1436 int SMESH_Mesh::NbNodes() const throw(SALOME_Exception)
1437 {
1438   Unexpect aCatch(SalomeException);
1439   return _myMeshDS->NbNodes();
1440 }
1441
1442 //================================================================================
1443 /*!
1444  * \brief  Return number of edges of given order in the mesh
1445  */
1446 //================================================================================
1447
1448 int SMESH_Mesh::Nb0DElements() const throw(SALOME_Exception)
1449 {
1450   Unexpect aCatch(SalomeException);
1451   return _myMeshDS->GetMeshInfo().Nb0DElements();
1452 }
1453
1454 //================================================================================
1455 /*!
1456  * \brief  Return number of edges of given order in the mesh
1457  */
1458 //================================================================================
1459
1460 int SMESH_Mesh::NbEdges(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1461 {
1462   Unexpect aCatch(SalomeException);
1463   return _myMeshDS->GetMeshInfo().NbEdges(order);
1464 }
1465
1466 //================================================================================
1467 /*!
1468  * \brief Return number of faces of given order in the mesh
1469  */
1470 //================================================================================
1471
1472 int SMESH_Mesh::NbFaces(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1473 {
1474   Unexpect aCatch(SalomeException);
1475   return _myMeshDS->GetMeshInfo().NbFaces(order);
1476 }
1477
1478 //================================================================================
1479 /*!
1480  * \brief Return the number of faces in the mesh
1481  */
1482 //================================================================================
1483
1484 int SMESH_Mesh::NbTriangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1485 {
1486   Unexpect aCatch(SalomeException);
1487   return _myMeshDS->GetMeshInfo().NbTriangles(order);
1488 }
1489
1490 //================================================================================
1491 /*!
1492  * \brief Return the number nodes faces in the mesh
1493  */
1494 //================================================================================
1495
1496 int SMESH_Mesh::NbQuadrangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1497 {
1498   Unexpect aCatch(SalomeException);
1499   return _myMeshDS->GetMeshInfo().NbQuadrangles(order);
1500 }
1501
1502 //================================================================================
1503 /*!
1504  * \brief Return number of biquadratic quadrangles in the mesh
1505  */
1506 //================================================================================
1507
1508 int SMESH_Mesh::NbBiQuadQuadrangles() const throw(SALOME_Exception)
1509 {
1510   Unexpect aCatch(SalomeException);
1511   return _myMeshDS->GetMeshInfo().NbBiQuadQuadrangles();
1512 }
1513
1514 //================================================================================
1515 /*!
1516  * \brief Return the number of polygonal faces in the mesh
1517  */
1518 //================================================================================
1519
1520 int SMESH_Mesh::NbPolygons() const throw(SALOME_Exception)
1521 {
1522   Unexpect aCatch(SalomeException);
1523   return _myMeshDS->GetMeshInfo().NbPolygons();
1524 }
1525
1526 //================================================================================
1527 /*!
1528  * \brief Return number of volumes of given order in the mesh
1529  */
1530 //================================================================================
1531
1532 int SMESH_Mesh::NbVolumes(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1533 {
1534   Unexpect aCatch(SalomeException);
1535   return _myMeshDS->GetMeshInfo().NbVolumes(order);
1536 }
1537
1538 //================================================================================
1539 /*!
1540  * \brief  Return number of tetrahedrons of given order in the mesh
1541  */
1542 //================================================================================
1543
1544 int SMESH_Mesh::NbTetras(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1545 {
1546   Unexpect aCatch(SalomeException);
1547   return _myMeshDS->GetMeshInfo().NbTetras(order);
1548 }
1549
1550 //================================================================================
1551 /*!
1552  * \brief  Return number of hexahedrons of given order in the mesh
1553  */
1554 //================================================================================
1555
1556 int SMESH_Mesh::NbHexas(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1557 {
1558   Unexpect aCatch(SalomeException);
1559   return _myMeshDS->GetMeshInfo().NbHexas(order);
1560 }
1561
1562 //================================================================================
1563 /*!
1564  * \brief  Return number of triquadratic hexahedrons in the mesh
1565  */
1566 //================================================================================
1567
1568 int SMESH_Mesh::NbTriQuadraticHexas() const throw(SALOME_Exception)
1569 {
1570   Unexpect aCatch(SalomeException);
1571   return _myMeshDS->GetMeshInfo().NbTriQuadHexas();
1572 }
1573
1574 //================================================================================
1575 /*!
1576  * \brief  Return number of pyramids of given order in the mesh
1577  */
1578 //================================================================================
1579
1580 int SMESH_Mesh::NbPyramids(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1581 {
1582   Unexpect aCatch(SalomeException);
1583   return _myMeshDS->GetMeshInfo().NbPyramids(order);
1584 }
1585
1586 //================================================================================
1587 /*!
1588  * \brief  Return number of prisms (penthahedrons) of given order in the mesh
1589  */
1590 //================================================================================
1591
1592 int SMESH_Mesh::NbPrisms(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1593 {
1594   Unexpect aCatch(SalomeException);
1595   return _myMeshDS->GetMeshInfo().NbPrisms(order);
1596 }
1597
1598 //================================================================================
1599 /*!
1600  * \brief  Return number of hexagonal prisms in the mesh
1601  */
1602 //================================================================================
1603
1604 int SMESH_Mesh::NbHexagonalPrisms() const throw(SALOME_Exception)
1605 {
1606   Unexpect aCatch(SalomeException);
1607   return _myMeshDS->GetMeshInfo().NbHexPrisms();
1608 }
1609
1610 //================================================================================
1611 /*!
1612  * \brief  Return number of polyhedrons in the mesh
1613  */
1614 //================================================================================
1615
1616 int SMESH_Mesh::NbPolyhedrons() const throw(SALOME_Exception)
1617 {
1618   Unexpect aCatch(SalomeException);
1619   return _myMeshDS->GetMeshInfo().NbPolyhedrons();
1620 }
1621
1622 //================================================================================
1623 /*!
1624  * \brief  Return number of ball elements in the mesh
1625  */
1626 //================================================================================
1627
1628 int SMESH_Mesh::NbBalls() const throw(SALOME_Exception)
1629 {
1630   Unexpect aCatch(SalomeException);
1631   return _myMeshDS->GetMeshInfo().NbBalls();
1632 }
1633
1634 //================================================================================
1635 /*!
1636  * \brief  Return number of submeshes in the mesh
1637  */
1638 //================================================================================
1639
1640 int SMESH_Mesh::NbSubMesh() const throw(SALOME_Exception)
1641 {
1642   Unexpect aCatch(SalomeException);
1643   return _myMeshDS->NbSubMesh();
1644 }
1645
1646 //=======================================================================
1647 //function : IsNotConformAllowed
1648 //purpose  : check if a hypothesis alowing notconform mesh is present
1649 //=======================================================================
1650
1651 bool SMESH_Mesh::IsNotConformAllowed() const
1652 {
1653   if(MYDEBUG) MESSAGE("SMESH_Mesh::IsNotConformAllowed");
1654
1655   static SMESH_HypoFilter filter( SMESH_HypoFilter::HasName( "NotConformAllowed" ));
1656   return GetHypothesis( _myMeshDS->ShapeToMesh(), filter, false );
1657 }
1658
1659 //=======================================================================
1660 //function : IsMainShape
1661 //purpose  : 
1662 //=======================================================================
1663
1664 bool SMESH_Mesh::IsMainShape(const TopoDS_Shape& theShape) const
1665 {
1666   return theShape.IsSame(_myMeshDS->ShapeToMesh() );
1667 }
1668
1669 //=============================================================================
1670 /*!
1671  *  
1672  */
1673 //=============================================================================
1674
1675 SMESH_Group* SMESH_Mesh::AddGroup (const SMDSAbs_ElementType theType,
1676                                    const char*               theName,
1677                                    int&                      theId,
1678                                    const TopoDS_Shape&       theShape,
1679                                    const SMESH_PredicatePtr& thePredicate)
1680 {
1681   if (_mapGroup.count(_groupId))
1682     return NULL;
1683   theId = _groupId;
1684   SMESH_Group* aGroup = new SMESH_Group (theId, this, theType, theName, theShape, thePredicate);
1685   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1686   _mapGroup[_groupId++] = aGroup;
1687   return aGroup;
1688 }
1689
1690 //================================================================================
1691 /*!
1692  * \brief Creates a group based on an existing SMESHDS group. Group ID should be unique
1693  */
1694 //================================================================================
1695
1696 SMESH_Group* SMESH_Mesh::AddGroup (SMESHDS_GroupBase* groupDS) throw(SALOME_Exception)
1697 {
1698   if ( !groupDS ) 
1699     throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup(): NULL SMESHDS_GroupBase"));
1700
1701   map <int, SMESH_Group*>::iterator i_g = _mapGroup.find( groupDS->GetID() );
1702   if ( i_g != _mapGroup.end() && i_g->second )
1703   {
1704     if ( i_g->second->GetGroupDS() == groupDS )
1705       return i_g->second;
1706     else
1707       throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup() wrong ID of SMESHDS_GroupBase"));
1708   }
1709   SMESH_Group* aGroup = new SMESH_Group (groupDS);
1710   _mapGroup[ groupDS->GetID() ] = aGroup;
1711   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1712
1713   _groupId = 1 + _mapGroup.rbegin()->first;
1714
1715   return aGroup;
1716 }
1717
1718
1719 //================================================================================
1720 /*!
1721  * \brief Creates SMESH_Groups for not wrapped SMESHDS_Groups
1722  *  \retval bool - true if new SMESH_Groups have been created
1723  * 
1724  */
1725 //================================================================================
1726
1727 bool SMESH_Mesh::SynchronizeGroups()
1728 {
1729   int nbGroups = _mapGroup.size();
1730   const set<SMESHDS_GroupBase*>& groups = _myMeshDS->GetGroups();
1731   set<SMESHDS_GroupBase*>::const_iterator gIt = groups.begin();
1732   for ( ; gIt != groups.end(); ++gIt )
1733   {
1734     SMESHDS_GroupBase* groupDS = (SMESHDS_GroupBase*) *gIt;
1735     _groupId = groupDS->GetID();
1736     if ( !_mapGroup.count( _groupId ))
1737       _mapGroup[_groupId] = new SMESH_Group( groupDS );
1738   }
1739   if ( !_mapGroup.empty() )
1740     _groupId = _mapGroup.rbegin()->first + 1;
1741
1742   return nbGroups < _mapGroup.size();
1743 }
1744
1745 //================================================================================
1746 /*!
1747  * \brief Return iterator on all existing groups
1748  */
1749 //================================================================================
1750
1751 SMESH_Mesh::GroupIteratorPtr SMESH_Mesh::GetGroups() const
1752 {
1753   typedef map <int, SMESH_Group *> TMap;
1754   return GroupIteratorPtr( new SMDS_mapIterator<TMap>( _mapGroup ));
1755 }
1756
1757 //=============================================================================
1758 /*!
1759  * \brief Return a group by ID
1760  */
1761 //=============================================================================
1762
1763 SMESH_Group* SMESH_Mesh::GetGroup (const int theGroupID)
1764 {
1765   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1766     return NULL;
1767   return _mapGroup[theGroupID];
1768 }
1769
1770
1771 //=============================================================================
1772 /*!
1773  * \brief Return IDs of all groups
1774  */
1775 //=============================================================================
1776
1777 list<int> SMESH_Mesh::GetGroupIds() const
1778 {
1779   list<int> anIds;
1780   for ( map<int, SMESH_Group*>::const_iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1781     anIds.push_back( it->first );
1782   
1783   return anIds;
1784 }
1785
1786 //================================================================================
1787 /*!
1788  * \brief Set a caller of methods at level of CORBA API implementation.
1789  * The set upCaller will be deleted by SMESH_Mesh
1790  */
1791 //================================================================================
1792
1793 void SMESH_Mesh::SetCallUp( TCallUp* upCaller )
1794 {
1795   if ( _callUp ) delete _callUp;
1796   _callUp = upCaller;
1797 }
1798
1799 //=============================================================================
1800 /*!
1801  *  
1802  */
1803 //=============================================================================
1804
1805 bool SMESH_Mesh::RemoveGroup (const int theGroupID)
1806 {
1807   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1808     return false;
1809   GetMeshDS()->RemoveGroup( _mapGroup[theGroupID]->GetGroupDS() );
1810   delete _mapGroup[theGroupID];
1811   _mapGroup.erase (theGroupID);
1812   if (_callUp)
1813     _callUp->RemoveGroup( theGroupID );
1814   return true;
1815 }
1816
1817 //=======================================================================
1818 //function : GetAncestors
1819 //purpose  : return list of ancestors of theSubShape in the order
1820 //           that lower dimention shapes come first.
1821 //=======================================================================
1822
1823 const TopTools_ListOfShape& SMESH_Mesh::GetAncestors(const TopoDS_Shape& theS) const
1824 {
1825   if ( _mapAncestors.Contains( theS ) )
1826     return _mapAncestors.FindFromKey( theS );
1827
1828   static TopTools_ListOfShape emptyList;
1829   return emptyList;
1830 }
1831
1832 //=======================================================================
1833 //function : Dump
1834 //purpose  : dumps contents of mesh to stream [ debug purposes ]
1835 //=======================================================================
1836
1837 ostream& SMESH_Mesh::Dump(ostream& save)
1838 {
1839   int clause = 0;
1840   save << "========================== Dump contents of mesh ==========================" << endl << endl;
1841   save << ++clause << ") Total number of nodes:   \t"    << NbNodes() << endl;
1842   save << ++clause << ") Total number of edges:   \t"    << NbEdges() << endl;
1843   save << ++clause << ") Total number of faces:   \t"    << NbFaces() << endl;
1844   save << ++clause << ") Total number of polygons:\t"    << NbPolygons() << endl;
1845   save << ++clause << ") Total number of volumes:\t"     << NbVolumes() << endl;
1846   save << ++clause << ") Total number of polyhedrons:\t" << NbPolyhedrons() << endl << endl;
1847   for ( int isQuadratic = 0; isQuadratic < 2; ++isQuadratic )
1848   {
1849     string orderStr = isQuadratic ? "quadratic" : "linear";
1850     SMDSAbs_ElementOrder order  = isQuadratic ? ORDER_QUADRATIC : ORDER_LINEAR;
1851
1852     save << ++clause << ") Total number of " << orderStr << " edges:\t" << NbEdges(order) << endl;
1853     save << ++clause << ") Total number of " << orderStr << " faces:\t" << NbFaces(order) << endl;
1854     if ( NbFaces(order) > 0 ) {
1855       int nb3 = NbTriangles(order);
1856       int nb4 = NbQuadrangles(order);
1857       save << clause << ".1) Number of " << orderStr << " triangles:  \t" << nb3 << endl;
1858       save << clause << ".2) Number of " << orderStr << " quadrangles:\t" << nb4 << endl;
1859       if ( nb3 + nb4 !=  NbFaces(order) ) {
1860         map<int,int> myFaceMap;
1861         SMDS_FaceIteratorPtr itFaces=_myMeshDS->facesIterator();
1862         while( itFaces->more( ) ) {
1863           int nbNodes = itFaces->next()->NbNodes();
1864           if ( myFaceMap.find( nbNodes ) == myFaceMap.end() )
1865             myFaceMap[ nbNodes ] = 0;
1866           myFaceMap[ nbNodes ] = myFaceMap[ nbNodes ] + 1;
1867         }
1868         save << clause << ".3) Faces in detail: " << endl;
1869         map <int,int>::iterator itF;
1870         for (itF = myFaceMap.begin(); itF != myFaceMap.end(); itF++)
1871           save << "--> nb nodes: " << itF->first << " - nb elemens:\t" << itF->second << endl;
1872       }
1873     }
1874     save << ++clause << ") Total number of " << orderStr << " volumes:\t" << NbVolumes(order) << endl;
1875     if ( NbVolumes(order) > 0 ) {
1876       int nb8 = NbHexas(order);
1877       int nb4 = NbTetras(order);
1878       int nb5 = NbPyramids(order);
1879       int nb6 = NbPrisms(order);
1880       save << clause << ".1) Number of " << orderStr << " hexahedrons:\t" << nb8 << endl;
1881       save << clause << ".2) Number of " << orderStr << " tetrahedrons:\t" << nb4 << endl;
1882       save << clause << ".3) Number of " << orderStr << " prisms:      \t" << nb6 << endl;
1883       save << clause << ".4) Number of " << orderStr << " pyramids:\t" << nb5 << endl;
1884       if ( nb8 + nb4 + nb5 + nb6 != NbVolumes(order) ) {
1885         map<int,int> myVolumesMap;
1886         SMDS_VolumeIteratorPtr itVolumes=_myMeshDS->volumesIterator();
1887         while( itVolumes->more( ) ) {
1888           int nbNodes = itVolumes->next()->NbNodes();
1889           if ( myVolumesMap.find( nbNodes ) == myVolumesMap.end() )
1890             myVolumesMap[ nbNodes ] = 0;
1891           myVolumesMap[ nbNodes ] = myVolumesMap[ nbNodes ] + 1;
1892         }
1893         save << clause << ".5) Volumes in detail: " << endl;
1894         map <int,int>::iterator itV;
1895         for (itV = myVolumesMap.begin(); itV != myVolumesMap.end(); itV++)
1896           save << "--> nb nodes: " << itV->first << " - nb elemens:\t" << itV->second << endl;
1897       }
1898     }
1899     save << endl;
1900   }
1901   save << "===========================================================================" << endl;
1902   return save;
1903 }
1904
1905 //=======================================================================
1906 //function : GetElementType
1907 //purpose  : Returns type of mesh element with certain id
1908 //=======================================================================
1909
1910 SMDSAbs_ElementType SMESH_Mesh::GetElementType( const int id, const bool iselem )
1911 {
1912   return _myMeshDS->GetElementType( id, iselem );
1913 }
1914
1915 //=============================================================================
1916 /*!
1917  *  \brief Convert group on geometry into standalone group
1918  */
1919 //=============================================================================
1920
1921 SMESH_Group* SMESH_Mesh::ConvertToStandalone ( int theGroupID )
1922 {
1923   SMESH_Group* aGroup = 0;
1924   map < int, SMESH_Group * >::iterator itg = _mapGroup.find( theGroupID );
1925   if ( itg == _mapGroup.end() )
1926     return aGroup;
1927
1928   SMESH_Group* anOldGrp = (*itg).second;
1929   SMESHDS_GroupBase* anOldGrpDS = anOldGrp->GetGroupDS();
1930   if ( !anOldGrp || !anOldGrpDS )
1931     return aGroup;
1932
1933   // create new standalone group
1934   aGroup = new SMESH_Group (theGroupID, this, anOldGrpDS->GetType(), anOldGrp->GetName() );
1935   _mapGroup[theGroupID] = aGroup;
1936
1937   SMESHDS_Group* aNewGrpDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1938   GetMeshDS()->RemoveGroup( anOldGrpDS );
1939   GetMeshDS()->AddGroup( aNewGrpDS );
1940
1941   // add elements (or nodes) into new created group
1942   SMDS_ElemIteratorPtr anItr = anOldGrpDS->GetElements();
1943   while ( anItr->more() )
1944     aNewGrpDS->Add( (anItr->next())->GetID() );
1945
1946   // set color
1947   aNewGrpDS->SetColor( anOldGrpDS->GetColor() );
1948
1949   // remove old group
1950   delete anOldGrp;
1951
1952   return aGroup;
1953 }
1954
1955 //=============================================================================
1956 /*!
1957  *  \brief remove submesh order  from Mesh
1958  */
1959 //=============================================================================
1960
1961 void SMESH_Mesh::ClearMeshOrder()
1962 {
1963   _mySubMeshOrder.clear();
1964 }
1965
1966 //=============================================================================
1967 /*!
1968  *  \brief remove submesh order  from Mesh
1969  */
1970 //=============================================================================
1971
1972 void SMESH_Mesh::SetMeshOrder(const TListOfListOfInt& theOrder )
1973 {
1974   _mySubMeshOrder = theOrder;
1975 }
1976
1977 //=============================================================================
1978 /*!
1979  *  \brief return submesh order if any
1980  */
1981 //=============================================================================
1982
1983 const TListOfListOfInt& SMESH_Mesh::GetMeshOrder() const
1984 {
1985   return _mySubMeshOrder;
1986 }
1987
1988 //=============================================================================
1989 /*!
1990  *  \brief fill _mapAncestors
1991  */
1992 //=============================================================================
1993
1994 void SMESH_Mesh::fillAncestorsMap(const TopoDS_Shape& theShape)
1995 {
1996
1997   int desType, ancType;
1998   if ( !theShape.IsSame( GetShapeToMesh()) && theShape.ShapeType() == TopAbs_COMPOUND )
1999   {
2000     // a geom group is added. Insert it into lists of ancestors before
2001     // the first ancestor more complex than group members
2002     int memberType = TopoDS_Iterator( theShape ).Value().ShapeType();
2003     for ( desType = TopAbs_VERTEX; desType >= memberType; desType-- )
2004       for (TopExp_Explorer des( theShape, TopAbs_ShapeEnum( desType )); des.More(); des.Next())
2005       {
2006         if ( !_mapAncestors.Contains( des.Current() )) continue;// issue 0020982
2007         TopTools_ListOfShape& ancList = _mapAncestors.ChangeFromKey( des.Current() );
2008         TopTools_ListIteratorOfListOfShape ancIt (ancList);
2009         while ( ancIt.More() && ancIt.Value().ShapeType() >= memberType )
2010           ancIt.Next();
2011         if ( ancIt.More() )
2012           ancList.InsertBefore( theShape, ancIt );
2013       }
2014   }
2015   {
2016     for ( desType = TopAbs_VERTEX; desType > TopAbs_COMPOUND; desType-- )
2017       for ( ancType = desType - 1; ancType >= TopAbs_COMPOUND; ancType-- )
2018         TopExp::MapShapesAndAncestors ( theShape,
2019                                         (TopAbs_ShapeEnum) desType,
2020                                         (TopAbs_ShapeEnum) ancType,
2021                                         _mapAncestors );
2022   }
2023   // visit COMPOUNDs inside a COMPOUND that are not reachable by TopExp_Explorer
2024   if ( theShape.ShapeType() == TopAbs_COMPOUND )
2025   {
2026     for ( TopoDS_Iterator sIt(theShape); sIt.More(); sIt.Next() )
2027       if ( sIt.Value().ShapeType() == TopAbs_COMPOUND )
2028         fillAncestorsMap( sIt.Value() );
2029   }
2030 }
2031
2032 //=============================================================================
2033 /*!
2034  * \brief sort submeshes according to stored mesh order
2035  * \param theListToSort in out list to be sorted
2036  * \return FALSE if nothing sorted
2037  */
2038 //=============================================================================
2039
2040 bool SMESH_Mesh::SortByMeshOrder(list<SMESH_subMesh*>& theListToSort) const
2041 {
2042   if ( !_mySubMeshOrder.size() || theListToSort.size() < 2)
2043     return true;
2044   
2045   bool res = false;
2046   list<SMESH_subMesh*> onlyOrderedList;
2047   // collect all ordered submeshes in one list as pointers
2048   // and get their positions within theListToSort
2049   typedef list<SMESH_subMesh*>::iterator TPosInList;
2050   map< int, TPosInList > sortedPos;
2051   TPosInList smBeg = theListToSort.begin(), smEnd = theListToSort.end();
2052   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2053   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++) {
2054     const TListOfInt& listOfId = *listIdsIt;
2055     TListOfInt::const_iterator idIt = listOfId.begin();
2056     for ( ; idIt != listOfId.end(); idIt++ ) {
2057       if ( SMESH_subMesh * sm = GetSubMeshContaining( *idIt )) {
2058         TPosInList smPos = find( smBeg, smEnd, sm );
2059         if ( smPos != smEnd ) {
2060           onlyOrderedList.push_back( sm );
2061           sortedPos[ distance( smBeg, smPos )] = smPos;
2062         }
2063       }
2064     }
2065   }
2066   if (onlyOrderedList.size() < 2)
2067     return res;
2068   res = true;
2069
2070   list<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
2071   list<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
2072
2073   // iterate on ordered submeshes and insert them in detected positions
2074   map< int, TPosInList >::iterator i_pos = sortedPos.begin();
2075   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
2076     *(i_pos->second) = *onlyBIt;
2077
2078   return res;
2079 }
2080
2081 //================================================================================
2082 /*!
2083  * \brief Return true if given order of sub-meshes is OK
2084  */
2085 //================================================================================
2086
2087 bool SMESH_Mesh::IsOrderOK( const SMESH_subMesh* smBefore,
2088                             const SMESH_subMesh* smAfter ) const
2089 {
2090   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2091   TListOfInt::const_iterator idBef, idAft;
2092   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2093   {
2094     const TListOfInt& listOfId = *listIdsIt;
2095     idBef = std::find( listOfId.begin(), listOfId.end(), smBefore->GetId() );
2096     if ( idBef != listOfId.end() )
2097       idAft = std::find( listOfId.begin(), listOfId.end(), smAfter->GetId() );
2098     if ( idAft != listOfId.end () )
2099       return ( std::distance( listOfId.begin(), idBef ) <
2100                std::distance( listOfId.begin(), idAft )   );
2101   }
2102   return true; // no order imposed to given submeshes
2103
2104
2105 //=============================================================================
2106 /*!
2107  * \brief sort submeshes according to stored mesh order
2108  * \param theListToSort in out list to be sorted
2109  * \return FALSE if nothing sorted
2110  */
2111 //=============================================================================
2112
2113 list<SMESH_subMesh*>
2114 SMESH_Mesh::getAncestorsSubMeshes (const TopoDS_Shape& theSubShape) const
2115 {
2116   list<SMESH_subMesh*> listOfSubMesh;
2117   TopTools_ListIteratorOfListOfShape it( GetAncestors( theSubShape ));
2118   for (; it.More(); it.Next() )
2119     if ( SMESH_subMesh* sm = GetSubMeshContaining( it.Value() ))
2120       listOfSubMesh.push_back(sm);
2121
2122   // sort submeshes according to stored mesh order
2123   SortByMeshOrder( listOfSubMesh );
2124
2125   return listOfSubMesh;
2126 }