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