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