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