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