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