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