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