Salome HOME
52546: Viscous Layers construction fails
[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 ? bool( _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,
606                           std::string*         anError  ) throw(SALOME_Exception)
607 {
608   Unexpect aCatch(SalomeException);
609   if(MYDEBUG) MESSAGE("SMESH_Mesh::AddHypothesis");
610
611   if ( anError )
612     anError->clear();
613
614   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
615   if ( !subMesh || !subMesh->GetId())
616     return SMESH_Hypothesis::HYP_BAD_SUBSHAPE;
617
618   SMESH_Hypothesis *anHyp = GetHypothesis( anHypId );
619   if ( !anHyp )
620     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
621
622   bool isGlobalHyp = IsMainShape( aSubShape );
623
624   // NotConformAllowed can be only global
625   if ( !isGlobalHyp )
626   {
627     // NOTE: this is not a correct way to check a name of hypothesis,
628     // there should be an attribute of hypothesis saying that it can/can't
629     // be global/local
630     string hypName = anHyp->GetName();
631     if ( hypName == "NotConformAllowed" )
632     {
633       if(MYDEBUG) MESSAGE( "Hypotesis <NotConformAllowed> can be only global" );
634       return SMESH_Hypothesis::HYP_INCOMPATIBLE;
635     }
636   }
637
638   // shape
639
640   bool isAlgo = ( anHyp->GetType() != SMESHDS_Hypothesis::PARAM_ALGO );
641   int   event = isAlgo ? SMESH_subMesh::ADD_ALGO : SMESH_subMesh::ADD_HYP;
642
643   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
644
645   if ( anError && SMESH_Hypothesis::IsStatusFatal(ret) && subMesh->GetComputeError() )
646     *anError = subMesh->GetComputeError()->myComment;
647
648   // sub-shapes
649   if ( !SMESH_Hypothesis::IsStatusFatal(ret) &&
650        anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is added on father
651   {
652     event = isAlgo ? SMESH_subMesh::ADD_FATHER_ALGO : SMESH_subMesh::ADD_FATHER_HYP;
653
654     SMESH_Hypothesis::Hypothesis_Status ret2 =
655       subMesh->SubMeshesAlgoStateEngine(event, anHyp, /*exitOnFatal=*/true);
656     if (ret2 > ret)
657     {
658       ret = ret2;
659       if ( SMESH_Hypothesis::IsStatusFatal( ret ))
660       {
661         if ( anError && subMesh->GetComputeError() )
662           *anError = subMesh->GetComputeError()->myComment;
663         // remove anHyp
664         event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
665         subMesh->AlgoStateEngine(event, anHyp);
666       }
667     }
668
669     // check concurent hypotheses on ancestors
670     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !isGlobalHyp )
671     {
672       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
673       while ( smIt->more() ) {
674         SMESH_subMesh* sm = smIt->next();
675         if ( sm->IsApplicableHypotesis( anHyp )) {
676           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
677           if (ret2 > ret) {
678             ret = ret2;
679             break;
680           }
681         }
682       }
683     }
684   }
685   HasModificationsToDiscard(); // to reset _isModified flag if a mesh becomes empty
686
687   if(MYDEBUG) subMesh->DumpAlgoState(true);
688   if(MYDEBUG) SCRUTE(ret);
689   return ret;
690 }
691
692 //=============================================================================
693 /*!
694  * 
695  */
696 //=============================================================================
697
698 SMESH_Hypothesis::Hypothesis_Status
699   SMESH_Mesh::RemoveHypothesis(const TopoDS_Shape & aSubShape,
700                                int anHypId)throw(SALOME_Exception)
701 {
702   Unexpect aCatch(SalomeException);
703   if(MYDEBUG) MESSAGE("SMESH_Mesh::RemoveHypothesis");
704   
705   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
706   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
707     throw SALOME_Exception(LOCALIZED("hypothesis does not exist"));
708   
709   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
710   if(MYDEBUG) {
711     SCRUTE(anHyp->GetType());
712   }
713   
714   // shape 
715   
716   bool isAlgo = ( !anHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO );
717   int event = isAlgo ? SMESH_subMesh::REMOVE_ALGO : SMESH_subMesh::REMOVE_HYP;
718
719   SMESH_subMesh *subMesh = GetSubMesh(aSubShape);
720
721   SMESH_Hypothesis::Hypothesis_Status ret = subMesh->AlgoStateEngine(event, anHyp);
722
723   // there may appear concurrent hyps that were covered by the removed hyp
724   if (ret < SMESH_Hypothesis::HYP_CONCURENT &&
725       subMesh->IsApplicableHypotesis( anHyp ) &&
726       subMesh->CheckConcurentHypothesis( anHyp->GetType() ) != SMESH_Hypothesis::HYP_OK)
727     ret = SMESH_Hypothesis::HYP_CONCURENT;
728
729   // sub-shapes
730   if (!SMESH_Hypothesis::IsStatusFatal(ret) &&
731       anHyp->GetDim() <= SMESH_Gen::GetShapeDim(aSubShape)) // is removed from father
732   {
733     event = isAlgo ? SMESH_subMesh::REMOVE_FATHER_ALGO : SMESH_subMesh::REMOVE_FATHER_HYP;
734
735     SMESH_Hypothesis::Hypothesis_Status ret2 =
736       subMesh->SubMeshesAlgoStateEngine(event, anHyp);
737     if (ret2 > ret) // more severe
738       ret = ret2;
739
740     // check concurent hypotheses on ancestors
741     if (ret < SMESH_Hypothesis::HYP_CONCURENT && !IsMainShape( aSubShape ) )
742     {
743       SMESH_subMeshIteratorPtr smIt = subMesh->getDependsOnIterator(false,false);
744       while ( smIt->more() ) {
745         SMESH_subMesh* sm = smIt->next();
746         if ( sm->IsApplicableHypotesis( anHyp )) {
747           ret2 = sm->CheckConcurentHypothesis( anHyp->GetType() );
748           if (ret2 > ret) {
749             ret = ret2;
750             break;
751           }
752         }
753       }
754     }
755   }
756
757   HasModificationsToDiscard(); // to reset _isModified flag if mesh become empty
758
759   if(MYDEBUG) subMesh->DumpAlgoState(true);
760   if(MYDEBUG) SCRUTE(ret);
761   return ret;
762 }
763
764 //=============================================================================
765 /*!
766  * 
767  */
768 //=============================================================================
769
770 const list<const SMESHDS_Hypothesis*>&
771 SMESH_Mesh::GetHypothesisList(const TopoDS_Shape & aSubShape) const
772   throw(SALOME_Exception)
773 {
774   return _myMeshDS->GetHypothesis(aSubShape);
775 }
776
777 //=======================================================================
778 /*!
779  * \brief Return the hypothesis assigned to the shape
780  *  \param aSubShape    - the shape to check
781  *  \param aFilter      - the hypothesis filter
782  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
783  *  \param assignedTo   - to return the shape the found hypo is assigned to
784  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
785  */
786 //=======================================================================
787
788 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const TopoDS_Shape &    aSubShape,
789                                                    const SMESH_HypoFilter& aFilter,
790                                                    const bool              andAncestors,
791                                                    TopoDS_Shape*           assignedTo) const
792 {
793   return GetHypothesis( const_cast< SMESH_Mesh* >(this)->GetSubMesh( aSubShape ),
794                         aFilter, andAncestors, assignedTo );
795 }
796
797 //=======================================================================
798 /*!
799  * \brief Return the hypothesis assigned to the shape of a sub-mesh
800  *  \param aSubMesh     - the sub-mesh to check
801  *  \param aFilter      - the hypothesis filter
802  *  \param andAncestors - flag to check hypos assigned to ancestors of the shape
803  *  \param assignedTo   - to return the shape the found hypo is assigned to
804  *  \retval SMESH_Hypothesis* - the first hypo passed through aFilter
805  */
806 //=======================================================================
807
808 const SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const SMESH_subMesh *   aSubMesh,
809                                                    const SMESH_HypoFilter& aFilter,
810                                                    const bool              andAncestors,
811                                                    TopoDS_Shape*           assignedTo) const
812 {
813   if ( !aSubMesh ) return 0;
814
815   {
816     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
817     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
818     list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
819     for ( ; hyp != hypList.end(); hyp++ ) {
820       const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
821       if ( aFilter.IsOk( h, aSubShape)) {
822         if ( assignedTo ) *assignedTo = aSubShape;
823         return h;
824       }
825     }
826   }
827   if ( andAncestors )
828   {
829     // user sorted submeshes of ancestors, according to stored submesh priority
830     std::vector< SMESH_subMesh * > & ancestors =
831       const_cast< std::vector< SMESH_subMesh * > & > ( aSubMesh->GetAncestors() );
832     SortByMeshOrder( ancestors );
833
834     vector<SMESH_subMesh*>::const_iterator smIt = ancestors.begin(); 
835     for ( ; smIt != ancestors.end(); smIt++ )
836     {
837       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
838       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
839       list<const SMESHDS_Hypothesis*>::const_iterator hyp = hypList.begin();
840       for ( ; hyp != hypList.end(); hyp++ ) {
841         const SMESH_Hypothesis * h = cSMESH_Hyp( *hyp );
842         if (aFilter.IsOk( h, curSh )) {
843           if ( assignedTo ) *assignedTo = curSh;
844           return h;
845         }
846       }
847     }
848   }
849   return 0;
850 }
851
852 //================================================================================
853 /*!
854  * \brief Return hypotheses assigned to the shape
855   * \param aSubShape - the shape to check
856   * \param aFilter - the hypothesis filter
857   * \param aHypList - the list of the found hypotheses
858   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
859   * \retval int - number of unique hypos in aHypList
860  */
861 //================================================================================
862
863 int SMESH_Mesh::GetHypotheses(const TopoDS_Shape &                aSubShape,
864                               const SMESH_HypoFilter&             aFilter,
865                               list <const SMESHDS_Hypothesis * >& aHypList,
866                               const bool                          andAncestors,
867                               list< TopoDS_Shape > *              assignedTo/*=0*/) const
868 {
869   return GetHypotheses( const_cast< SMESH_Mesh* >(this)->GetSubMesh( aSubShape ),
870                         aFilter, aHypList, andAncestors, assignedTo );
871 }
872
873 //================================================================================
874 /*!
875  * \brief Return hypotheses assigned to the shape of a sub-mesh
876   * \param aSubShape - the sub-mesh to check
877   * \param aFilter - the hypothesis filter
878   * \param aHypList - the list of the found hypotheses
879   * \param andAncestors - flag to check hypos assigned to ancestors of the shape
880   * \retval int - number of unique hypos in aHypList
881  */
882 //================================================================================
883
884 int SMESH_Mesh::GetHypotheses(const SMESH_subMesh *               aSubMesh,
885                               const SMESH_HypoFilter&             aFilter,
886                               list <const SMESHDS_Hypothesis * >& aHypList,
887                               const bool                          andAncestors,
888                               list< TopoDS_Shape > *              assignedTo/*=0*/) const
889 {
890   if ( !aSubMesh ) return 0;
891
892   set<string> hypTypes; // to exclude same type hypos from the result list
893   int nbHyps = 0;
894
895   // only one main hypothesis is allowed
896   bool mainHypFound = false;
897
898   // fill in hypTypes
899   list<const SMESHDS_Hypothesis*>::const_iterator hyp;
900   for ( hyp = aHypList.begin(); hyp != aHypList.end(); hyp++ ) {
901     if ( hypTypes.insert( (*hyp)->GetName() ).second )
902       nbHyps++;
903     if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
904       mainHypFound = true;
905   }
906
907   // get hypos from aSubShape
908   {
909     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
910     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
911     for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
912     {
913       const SMESH_Hypothesis* h = cSMESH_Hyp( *hyp );
914       if (( aFilter.IsOk( h, aSubShape )) &&
915           ( h->IsAuxiliary() || !mainHypFound ) &&
916           ( h->IsAuxiliary() || hypTypes.insert( h->GetName() ).second ))
917       {
918         aHypList.push_back( *hyp );
919         nbHyps++;
920         if ( !h->IsAuxiliary() )
921           mainHypFound = true;
922         if ( assignedTo ) assignedTo->push_back( aSubShape );
923       }
924     }
925   }
926
927   // get hypos from ancestors of aSubShape
928   if ( andAncestors )
929   {
930     // user sorted submeshes of ancestors, according to stored submesh priority
931     std::vector< SMESH_subMesh * > & ancestors =
932       const_cast< std::vector< SMESH_subMesh * > & > ( aSubMesh->GetAncestors() );
933     SortByMeshOrder( ancestors );
934
935     vector<SMESH_subMesh*>::const_iterator smIt = ancestors.begin();
936     for ( ; smIt != ancestors.end(); smIt++ )
937     {
938       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
939       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
940       for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
941       {
942         const SMESH_Hypothesis* h = cSMESH_Hyp( *hyp );
943         if (( aFilter.IsOk( h, curSh )) &&
944             ( h->IsAuxiliary() || !mainHypFound ) &&
945             ( h->IsAuxiliary() || hypTypes.insert( h->GetName() ).second ))
946         {
947           aHypList.push_back( *hyp );
948           nbHyps++;
949           if ( !h->IsAuxiliary() )
950             mainHypFound = true;
951           if ( assignedTo ) assignedTo->push_back( curSh );
952         }
953       }
954     }
955   }
956   return nbHyps;
957 }
958
959 //================================================================================
960 /*!
961  * \brief Return a hypothesis by its ID
962  */
963 //================================================================================
964
965 SMESH_Hypothesis * SMESH_Mesh::GetHypothesis(const int anHypId) const
966 {
967   StudyContextStruct *sc = _gen->GetStudyContext(_studyId);
968   if (sc->mapHypothesis.find(anHypId) == sc->mapHypothesis.end())
969     return NULL;
970
971   SMESH_Hypothesis *anHyp = sc->mapHypothesis[anHypId];
972   return anHyp;
973 }
974
975 //=============================================================================
976 /*!
977  * 
978  */
979 //=============================================================================
980
981 const list<SMESHDS_Command*> & SMESH_Mesh::GetLog() throw(SALOME_Exception)
982 {
983   Unexpect aCatch(SalomeException);
984   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetLog");
985   return _myMeshDS->GetScript()->GetCommands();
986 }
987
988 //=============================================================================
989 /*!
990  * 
991  */
992 //=============================================================================
993 void SMESH_Mesh::ClearLog() throw(SALOME_Exception)
994 {
995   Unexpect aCatch(SalomeException);
996   if(MYDEBUG) MESSAGE("SMESH_Mesh::ClearLog");
997   _myMeshDS->GetScript()->Clear();
998 }
999
1000 //=============================================================================
1001 /*!
1002  * Get or Create the SMESH_subMesh object implementation
1003  */
1004 //=============================================================================
1005
1006 SMESH_subMesh *SMESH_Mesh::GetSubMesh(const TopoDS_Shape & aSubShape)
1007   throw(SALOME_Exception)
1008 {
1009   int index = _myMeshDS->ShapeToIndex(aSubShape);
1010   if ( !index && aSubShape.IsNull() )
1011     return 0;
1012
1013   // for submeshes on GEOM Group
1014   if (( !index || index > _nbSubShapes ) && aSubShape.ShapeType() == TopAbs_COMPOUND ) {
1015     TopoDS_Iterator it( aSubShape );
1016     if ( it.More() )
1017     {
1018       index = _myMeshDS->AddCompoundSubmesh( aSubShape, it.Value().ShapeType() );
1019       // fill map of Ancestors
1020       while ( _nbSubShapes < index )
1021         fillAncestorsMap( _myMeshDS->IndexToShape( ++_nbSubShapes ));
1022     }
1023   }
1024   // if ( !index )
1025   //   return NULL; // neither sub-shape nor a group
1026
1027   SMESH_subMesh* aSubMesh = _subMeshHolder->Get( index );
1028   if ( !aSubMesh )
1029   {
1030     aSubMesh = new SMESH_subMesh(index, this, _myMeshDS, aSubShape);
1031     _subMeshHolder->Add( index, aSubMesh );
1032
1033     // include non-computable sub-meshes in SMESH_subMesh::_ancestors of sub-submeshes
1034     switch ( aSubShape.ShapeType() ) {
1035     case TopAbs_COMPOUND:
1036     case TopAbs_WIRE:
1037     case TopAbs_SHELL:
1038       for ( TopoDS_Iterator subIt( aSubShape ); subIt.More(); subIt.Next() )
1039       {
1040         SMESH_subMesh* sm = GetSubMesh( subIt.Value() );
1041         SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*inclideSelf=*/true);
1042         while ( smIt->more() )
1043           smIt->next()->ClearAncestors();
1044       }
1045     default:;
1046     }
1047   }
1048   return aSubMesh;
1049 }
1050
1051 //=============================================================================
1052 /*!
1053  * Get the SMESH_subMesh object implementation. Dont create it, return null
1054  * if it does not exist.
1055  */
1056 //=============================================================================
1057
1058 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const TopoDS_Shape & aSubShape) const
1059   throw(SALOME_Exception)
1060 {
1061   int index = _myMeshDS->ShapeToIndex(aSubShape);
1062   return GetSubMeshContaining( index );
1063 }
1064
1065 //=============================================================================
1066 /*!
1067  * Get the SMESH_subMesh object implementation. Dont create it, return null
1068  * if it does not exist.
1069  */
1070 //=============================================================================
1071
1072 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const int aShapeID) const
1073 throw(SALOME_Exception)
1074 {
1075   SMESH_subMesh *aSubMesh = _subMeshHolder->Get( aShapeID );
1076
1077   return aSubMesh;
1078 }
1079
1080 //================================================================================
1081 /*!
1082  * \brief Return submeshes of groups containing the given sub-shape
1083  */
1084 //================================================================================
1085
1086 list<SMESH_subMesh*>
1087 SMESH_Mesh::GetGroupSubMeshesContaining(const TopoDS_Shape & aSubShape) const
1088   throw(SALOME_Exception)
1089 {
1090   list<SMESH_subMesh*> found;
1091
1092   SMESH_subMesh * subMesh = GetSubMeshContaining(aSubShape);
1093   if ( !subMesh )
1094     return found;
1095
1096   // submeshes of groups have max IDs, so search from the map end
1097 SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator( /*reverse=*/true ) );
1098   while ( smIt->more() ) {
1099     SMESH_subMesh*    sm = smIt->next();
1100     SMESHDS_SubMesh * ds = sm->GetSubMeshDS();
1101     if ( ds && ds->IsComplexSubmesh() ) {
1102       if ( SMESH_MesherHelper::IsSubShape( aSubShape, sm->GetSubShape() ))
1103       {
1104         found.push_back( sm );
1105         //break;
1106       }
1107     } else {
1108       break; // the rest sub-meshes are not those of groups
1109     }
1110   }
1111
1112   if ( found.empty() ) // maybe the main shape is a COMPOUND (issue 0021530)
1113   {
1114     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1) )
1115       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1116       {
1117         TopoDS_Iterator it( mainSM->GetSubShape() );
1118         if ( it.Value().ShapeType() == aSubShape.ShapeType() &&
1119              SMESH_MesherHelper::IsSubShape( aSubShape, mainSM->GetSubShape() ))
1120           found.push_back( mainSM );
1121       }
1122   }
1123   return found;
1124 }
1125 //=======================================================================
1126 //function : IsUsedHypothesis
1127 //purpose  : Return True if anHyp is used to mesh aSubShape
1128 //=======================================================================
1129
1130 bool SMESH_Mesh::IsUsedHypothesis(SMESHDS_Hypothesis * anHyp,
1131                                   const SMESH_subMesh* aSubMesh)
1132 {
1133   SMESH_Hypothesis* hyp = static_cast<SMESH_Hypothesis*>(anHyp);
1134
1135   // check if anHyp can be used to mesh aSubMesh
1136   if ( !aSubMesh || !aSubMesh->IsApplicableHypotesis( hyp ))
1137     return false;
1138
1139   SMESH_Algo *algo = aSubMesh->GetAlgo();
1140
1141   // algorithm
1142   if (anHyp->GetType() > SMESHDS_Hypothesis::PARAM_ALGO)
1143     return ( anHyp == algo );
1144
1145   // algorithm parameter
1146   if (algo)
1147   {
1148     // look trough hypotheses used by algo
1149     const SMESH_HypoFilter* hypoKind;
1150     if (( hypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() ))) {
1151       list <const SMESHDS_Hypothesis * > usedHyps;
1152       if ( GetHypotheses( aSubMesh, *hypoKind, usedHyps, true ))
1153         return ( find( usedHyps.begin(), usedHyps.end(), anHyp ) != usedHyps.end() );
1154     }
1155   }
1156
1157   return false;
1158 }
1159
1160 //=============================================================================
1161 /*!
1162  *
1163  */
1164 //=============================================================================
1165
1166 // const list < SMESH_subMesh * >&
1167 // SMESH_Mesh::GetSubMeshUsingHypothesis(SMESHDS_Hypothesis * anHyp)
1168 //   throw(SALOME_Exception)
1169 // {
1170 //   _subMeshesUsingHypothesisList.clear();
1171 //   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1172 //   while ( smIt->more() )
1173 //   {
1174 //     SMESH_subMesh* aSubMesh = smIt->next();
1175 //     if ( IsUsedHypothesis ( anHyp, aSubMesh ))
1176 //       _subMeshesUsingHypothesisList.push_back( aSubMesh );
1177 //   }
1178 //   return _subMeshesUsingHypothesisList;
1179 // }
1180
1181 //=======================================================================
1182 //function : NotifySubMeshesHypothesisModification
1183 //purpose  : Say all submeshes using theChangedHyp that it has been modified
1184 //=======================================================================
1185
1186 void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* hyp)
1187 {
1188   Unexpect aCatch(SalomeException);
1189
1190   if ( !GetMeshDS()->IsUsedHypothesis( hyp ))
1191     return;
1192
1193   if (_callUp)
1194     _callUp->HypothesisModified();
1195
1196   SMESH_Algo *algo;
1197   const SMESH_HypoFilter* compatibleHypoKind;
1198   list <const SMESHDS_Hypothesis * > usedHyps;
1199
1200   // keep sub-meshes not to miss ones whose state can change due to notifying others
1201   vector< SMESH_subMesh* > smToNotify;
1202
1203   SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1204   while ( smIt->more() )
1205   {
1206     SMESH_subMesh* aSubMesh = smIt->next();
1207
1208     // if aSubMesh meshing depends on hyp,
1209     // we call aSubMesh->AlgoStateEngine( MODIF_HYP, hyp ) that causes either
1210     // 1) clearing of already computed aSubMesh or
1211     // 2) changing algo_state from MISSING_HYP to HYP_OK when parameters of hyp becomes valid,
1212     // other possible changes are not interesting. (IPAL0052457 - assigning hyp performance pb)
1213     if ( aSubMesh->GetComputeState() != SMESH_subMesh::COMPUTE_OK &&
1214          aSubMesh->GetComputeState() != SMESH_subMesh::FAILED_TO_COMPUTE &&
1215          aSubMesh->GetAlgoState()    != SMESH_subMesh::MISSING_HYP )
1216       continue;
1217
1218     const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
1219
1220     if (( aSubMesh->IsApplicableHypotesis( hyp )) &&
1221         ( algo = aSubMesh->GetAlgo() )            &&
1222         ( compatibleHypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() )) &&
1223         ( compatibleHypoKind->IsOk( hyp, aSubShape )))
1224     {
1225       // check if hyp is used by algo
1226       usedHyps.clear();
1227       if ( GetHypotheses( aSubMesh, *compatibleHypoKind, usedHyps, true ) &&
1228            find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() )
1229       {
1230         smToNotify.push_back( aSubMesh );
1231       }
1232     }
1233   }
1234
1235   for ( size_t i = 0; i < smToNotify.size(); ++i )
1236   {
1237     smToNotify[i]->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
1238                                    const_cast< SMESH_Hypothesis*>( hyp ));
1239   }
1240
1241   HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
1242   GetMeshDS()->Modified();
1243 }
1244
1245 //=============================================================================
1246 /*!
1247  *  Auto color functionality
1248  */
1249 //=============================================================================
1250 void SMESH_Mesh::SetAutoColor(bool theAutoColor) throw(SALOME_Exception)
1251 {
1252   Unexpect aCatch(SalomeException);
1253   _isAutoColor = theAutoColor;
1254 }
1255
1256 bool SMESH_Mesh::GetAutoColor() throw(SALOME_Exception)
1257 {
1258   Unexpect aCatch(SalomeException);
1259   return _isAutoColor;
1260 }
1261
1262 //=======================================================================
1263 //function : SetIsModified
1264 //purpose  : Set the flag meaning that the mesh has been edited "manually"
1265 //=======================================================================
1266
1267 void SMESH_Mesh::SetIsModified(bool isModified)
1268 {
1269   _isModified = isModified;
1270
1271   if ( _isModified )
1272     // check if mesh becomes empty as result of modification
1273     HasModificationsToDiscard();
1274 }
1275
1276 //=======================================================================
1277 //function : HasModificationsToDiscard
1278 //purpose  : Return true if the mesh has been edited since a total re-compute
1279 //           and those modifications may prevent successful partial re-compute.
1280 //           As a side effect reset _isModified flag if mesh is empty
1281 //issue    : 0020693
1282 //=======================================================================
1283
1284 bool SMESH_Mesh::HasModificationsToDiscard() const
1285 {
1286   if ( ! _isModified )
1287     return false;
1288
1289   // return true if the next Compute() will be partial and
1290   // existing but changed elements may prevent successful re-compute
1291   bool hasComputed = false, hasNotComputed = false;
1292 SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
1293   while ( smIt->more() )
1294   {
1295     const SMESH_subMesh* aSubMesh = smIt->next();
1296     switch ( aSubMesh->GetSubShape().ShapeType() )
1297     {
1298     case TopAbs_EDGE:
1299     case TopAbs_FACE:
1300     case TopAbs_SOLID:
1301       if ( aSubMesh->IsMeshComputed() )
1302         hasComputed = true;
1303       else
1304         hasNotComputed = true;
1305       if ( hasComputed && hasNotComputed)
1306         return true;
1307     }
1308   }
1309   if ( NbNodes() < 1 )
1310     const_cast<SMESH_Mesh*>(this)->_isModified = false;
1311
1312   return false;
1313 }
1314
1315 //================================================================================
1316 /*!
1317  * \brief Check if any groups of the same type have equal names
1318  */
1319 //================================================================================
1320
1321 bool SMESH_Mesh::HasDuplicatedGroupNamesMED()
1322 {
1323   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1324   map< SMDSAbs_ElementType, set<string> > aGroupNames;
1325   for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1326   {
1327     SMESH_Group* aGroup = it->second;
1328     SMDSAbs_ElementType aType = aGroup->GetGroupDS()->GetType();
1329     string aGroupName = aGroup->GetName();
1330     aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1331     if (!aGroupNames[aType].insert(aGroupName).second)
1332       return true;
1333   }
1334
1335   return false;
1336 }
1337
1338 //================================================================================
1339 /*!
1340  * \brief Export the mesh to a med file
1341  *  \param [in] file - name of the MED file
1342  *  \param [in] theMeshName - name of this mesh
1343  *  \param [in] theAutoGroups - boolean parameter for creating/not creating
1344  *              the groups Group_On_All_Nodes, Group_On_All_Faces, ... ;
1345  *              the typical use is auto_groups=false.
1346  *  \param [in] theVersion - defines the version of format of MED file, that will be created
1347  *  \param [in] meshPart - mesh data to export
1348  *  \param [in] theAutoDimension - if \c true, a space dimension of a MED mesh can be either
1349      *         - 1D if all mesh nodes lie on OX coordinate axis, or
1350      *         - 2D if all mesh nodes lie on XOY coordinate plane, or
1351      *         - 3D in the rest cases.
1352      *         If \a theAutoDimension is \c false, the space dimension is always 3.
1353  *  \return int - mesh index in the file
1354  */
1355 //================================================================================
1356
1357 void SMESH_Mesh::ExportMED(const char *        file, 
1358                            const char*         theMeshName, 
1359                            bool                theAutoGroups,
1360                            int                 theVersion,
1361                            const SMESHDS_Mesh* meshPart,
1362                            bool                theAutoDimension,
1363                            bool                theAddODOnVertices)
1364   throw(SALOME_Exception)
1365 {
1366   SMESH_TRY;
1367
1368   DriverMED_W_SMESHDS_Mesh myWriter;
1369   myWriter.SetFile         ( file, MED::EVersion(theVersion) );
1370   myWriter.SetMesh         ( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS   );
1371   myWriter.SetAutoDimension( theAutoDimension );
1372   myWriter.AddODOnVertices ( theAddODOnVertices );
1373   if ( !theMeshName ) 
1374     myWriter.SetMeshId     ( _id         );
1375   else {
1376     myWriter.SetMeshId     ( -1          );
1377     myWriter.SetMeshName   ( theMeshName );
1378   }
1379
1380   if ( theAutoGroups ) {
1381     myWriter.AddGroupOfNodes();
1382     myWriter.AddGroupOfEdges();
1383     myWriter.AddGroupOfFaces();
1384     myWriter.AddGroupOfVolumes();
1385   }
1386
1387   // Pass groups to writer. Provide unique group names.
1388   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1389   if ( !meshPart )
1390   {
1391     map< SMDSAbs_ElementType, set<string> > aGroupNames;
1392     char aString [256];
1393     int maxNbIter = 10000; // to guarantee cycle finish
1394     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1395       SMESH_Group*       aGroup   = it->second;
1396       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1397       if ( aGroupDS ) {
1398         SMDSAbs_ElementType aType = aGroupDS->GetType();
1399         string aGroupName0 = aGroup->GetName();
1400         aGroupName0.resize(MAX_MED_GROUP_NAME_LENGTH);
1401         string aGroupName = aGroupName0;
1402         for (int i = 1; !aGroupNames[aType].insert(aGroupName).second && i < maxNbIter; i++) {
1403           sprintf(&aString[0], "GR_%d_%s", i, aGroupName0.c_str());
1404           aGroupName = aString;
1405           aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1406         }
1407         aGroupDS->SetStoreName( aGroupName.c_str() );
1408         myWriter.AddGroup( aGroupDS );
1409       }
1410     }
1411   }
1412   // Perform export
1413   myWriter.Perform();
1414
1415   SMESH_CATCH( SMESH::throwSalomeEx );
1416 }
1417
1418 //================================================================================
1419 /*!
1420  * \brief Export the mesh to a SAUV file
1421  */
1422 //================================================================================
1423
1424 void SMESH_Mesh::ExportSAUV(const char *file, 
1425                             const char* theMeshName, 
1426                             bool theAutoGroups)
1427   throw(SALOME_Exception)
1428 {
1429   std::string medfilename(file);
1430   medfilename += ".med";
1431   std::string cmd;
1432 #ifdef WIN32
1433   cmd = "%PYTHONBIN% ";
1434 #else
1435   cmd = "python ";
1436 #endif
1437   cmd += "-c \"";
1438   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1439   cmd += "\"";
1440   system(cmd.c_str());
1441   ExportMED(medfilename.c_str(), theMeshName, theAutoGroups, 1);
1442 #ifdef WIN32
1443   cmd = "%PYTHONBIN% ";
1444 #else
1445   cmd = "python ";
1446 #endif
1447   cmd += "-c \"";
1448   cmd += "from medutilities import convert ; convert(r'" + medfilename + "', 'MED', 'GIBI', 1, r'" + file + "')";
1449   cmd += "\"";
1450   system(cmd.c_str());
1451 #ifdef WIN32
1452   cmd = "%PYTHONBIN% ";
1453 #else
1454   cmd = "python ";
1455 #endif
1456   cmd += "-c \"";
1457   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1458   cmd += "\"";
1459   system(cmd.c_str());
1460 }
1461
1462 //================================================================================
1463 /*!
1464  * \brief Export the mesh to a DAT file
1465  */
1466 //================================================================================
1467
1468 void SMESH_Mesh::ExportDAT(const char *        file,
1469                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1470 {
1471   Unexpect aCatch(SalomeException);
1472   DriverDAT_W_SMDS_Mesh myWriter;
1473   myWriter.SetFile( file );
1474   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1475   myWriter.SetMeshId(_id);
1476   myWriter.Perform();
1477 }
1478
1479 //================================================================================
1480 /*!
1481  * \brief Export the mesh to an UNV file
1482  */
1483 //================================================================================
1484
1485 void SMESH_Mesh::ExportUNV(const char *        file,
1486                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1487 {
1488   Unexpect aCatch(SalomeException);
1489   DriverUNV_W_SMDS_Mesh myWriter;
1490   myWriter.SetFile( file );
1491   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1492   myWriter.SetMeshId(_id);
1493   //  myWriter.SetGroups(_mapGroup);
1494
1495   if ( !meshPart )
1496   {
1497     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1498       SMESH_Group*       aGroup   = it->second;
1499       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1500       if ( aGroupDS ) {
1501         string aGroupName = aGroup->GetName();
1502         aGroupDS->SetStoreName( aGroupName.c_str() );
1503         myWriter.AddGroup( aGroupDS );
1504       }
1505     }
1506   }
1507   myWriter.Perform();
1508 }
1509
1510 //================================================================================
1511 /*!
1512  * \brief Export the mesh to an STL file
1513  */
1514 //================================================================================
1515
1516 void SMESH_Mesh::ExportSTL(const char *        file,
1517                            const bool          isascii,
1518                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1519 {
1520   Unexpect aCatch(SalomeException);
1521   DriverSTL_W_SMDS_Mesh myWriter;
1522   myWriter.SetFile( file );
1523   myWriter.SetIsAscii( isascii );
1524   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS);
1525   myWriter.SetMeshId(_id);
1526   myWriter.Perform();
1527 }
1528
1529 //================================================================================
1530 /*!
1531  * \brief Export the mesh to the CGNS file
1532  */
1533 //================================================================================
1534
1535 void SMESH_Mesh::ExportCGNS(const char *        file,
1536                             const SMESHDS_Mesh* meshDS)
1537 {
1538   int res = Driver_Mesh::DRS_FAIL;
1539 #ifdef WITH_CGNS
1540   DriverCGNS_Write myWriter;
1541   myWriter.SetFile( file );
1542   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1543   myWriter.SetMeshName( SMESH_Comment("Mesh_") << meshDS->GetPersistentId());
1544   res = myWriter.Perform();
1545 #endif
1546   if ( res != Driver_Mesh::DRS_OK )
1547     throw SALOME_Exception("Export failed");
1548 }
1549
1550 //================================================================================
1551 /*!
1552  * \brief Export the mesh to a GMF file
1553  */
1554 //================================================================================
1555
1556 void SMESH_Mesh::ExportGMF(const char *        file,
1557                            const SMESHDS_Mesh* meshDS,
1558                            bool                withRequiredGroups)
1559 {
1560   DriverGMF_Write myWriter;
1561   myWriter.SetFile( file );
1562   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1563   myWriter.SetExportRequiredGroups( withRequiredGroups );
1564
1565   myWriter.Perform();
1566 }
1567
1568 //================================================================================
1569 /*!
1570  * \brief Return a ratio of "compute cost" of computed sub-meshes to the whole
1571  *        "compute cost".
1572  */
1573 //================================================================================
1574
1575 double SMESH_Mesh::GetComputeProgress() const
1576 {
1577   double totalCost = 1e-100, computedCost = 0;
1578   const SMESH_subMesh* curSM = _gen->GetCurrentSubMesh();
1579
1580   // get progress of a current algo
1581   TColStd_MapOfInteger currentSubIds; 
1582   if ( curSM )
1583     if ( SMESH_Algo* algo = curSM->GetAlgo() )
1584     {
1585       int algoNotDoneCost = 0, algoDoneCost = 0;
1586       const std::vector<SMESH_subMesh*>& smToCompute = algo->SubMeshesToCompute();
1587       for ( size_t i = 0; i < smToCompute.size(); ++i )
1588       {
1589         if ( smToCompute[i]->IsEmpty() )
1590           algoNotDoneCost += smToCompute[i]->GetComputeCost();
1591         else
1592           algoDoneCost += smToCompute[i]->GetComputeCost();
1593         currentSubIds.Add( smToCompute[i]->GetId() );
1594       }
1595       double rate = 0;
1596       try
1597       {
1598         OCC_CATCH_SIGNALS;
1599         rate = algo->GetProgress();
1600       }
1601       catch (...) {
1602 #ifdef _DEBUG_
1603         cerr << "Exception in " << algo->GetName() << "::GetProgress()" << endl;
1604 #endif
1605       }
1606       if ( 0. < rate && rate < 1.001 )
1607       {
1608         computedCost += rate * ( algoDoneCost + algoNotDoneCost );
1609       }
1610       else
1611       {
1612         rate = algo->GetProgressByTic();
1613         computedCost += algoDoneCost + rate * algoNotDoneCost;
1614       }
1615       // cout << "rate: "<<rate << " algoNotDoneCost: " << algoNotDoneCost << endl;
1616     }
1617
1618   // get cost of already treated sub-meshes
1619   if ( SMESH_subMesh* mainSM = GetSubMeshContaining( 1 ))
1620   {
1621     SMESH_subMeshIteratorPtr smIt = mainSM->getDependsOnIterator(/*includeSelf=*/true);
1622     while ( smIt->more() )
1623     {
1624       const SMESH_subMesh* sm = smIt->next();
1625       const int smCost = sm->GetComputeCost();
1626       totalCost += smCost;
1627       if ( !currentSubIds.Contains( sm->GetId() ) )
1628       {
1629         if (( !sm->IsEmpty() ) ||
1630             ( sm->GetComputeState() == SMESH_subMesh::FAILED_TO_COMPUTE &&
1631               !sm->DependsOn( curSM ) ))
1632           computedCost += smCost;
1633       }
1634     }
1635   }
1636   // cout << "Total: " << totalCost
1637   //      << " computed: " << computedCost << " progress: " << computedCost / totalCost
1638   //      << " nbElems: " << GetMeshDS()->GetMeshInfo().NbElements() << endl;
1639   return computedCost / totalCost;
1640 }
1641
1642 //================================================================================
1643 /*!
1644  * \brief Return number of nodes in the mesh
1645  */
1646 //================================================================================
1647
1648 int SMESH_Mesh::NbNodes() const throw(SALOME_Exception)
1649 {
1650   Unexpect aCatch(SalomeException);
1651   return _myMeshDS->NbNodes();
1652 }
1653
1654 //================================================================================
1655 /*!
1656  * \brief  Return number of edges of given order in the mesh
1657  */
1658 //================================================================================
1659
1660 int SMESH_Mesh::Nb0DElements() const throw(SALOME_Exception)
1661 {
1662   Unexpect aCatch(SalomeException);
1663   return _myMeshDS->GetMeshInfo().Nb0DElements();
1664 }
1665
1666 //================================================================================
1667 /*!
1668  * \brief  Return number of edges of given order in the mesh
1669  */
1670 //================================================================================
1671
1672 int SMESH_Mesh::NbEdges(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1673 {
1674   Unexpect aCatch(SalomeException);
1675   return _myMeshDS->GetMeshInfo().NbEdges(order);
1676 }
1677
1678 //================================================================================
1679 /*!
1680  * \brief Return number of faces of given order in the mesh
1681  */
1682 //================================================================================
1683
1684 int SMESH_Mesh::NbFaces(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1685 {
1686   Unexpect aCatch(SalomeException);
1687   return _myMeshDS->GetMeshInfo().NbFaces(order);
1688 }
1689
1690 //================================================================================
1691 /*!
1692  * \brief Return the number of faces in the mesh
1693  */
1694 //================================================================================
1695
1696 int SMESH_Mesh::NbTriangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1697 {
1698   Unexpect aCatch(SalomeException);
1699   return _myMeshDS->GetMeshInfo().NbTriangles(order);
1700 }
1701
1702 //================================================================================
1703 /*!
1704  * \brief Return number of biquadratic triangles in the mesh
1705  */
1706 //================================================================================
1707
1708 int SMESH_Mesh::NbBiQuadTriangles() const throw(SALOME_Exception)
1709 {
1710   Unexpect aCatch(SalomeException);
1711   return _myMeshDS->GetMeshInfo().NbBiQuadTriangles();
1712 }
1713
1714 //================================================================================
1715 /*!
1716  * \brief Return the number nodes faces in the mesh
1717  */
1718 //================================================================================
1719
1720 int SMESH_Mesh::NbQuadrangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1721 {
1722   Unexpect aCatch(SalomeException);
1723   return _myMeshDS->GetMeshInfo().NbQuadrangles(order);
1724 }
1725
1726 //================================================================================
1727 /*!
1728  * \brief Return number of biquadratic quadrangles in the mesh
1729  */
1730 //================================================================================
1731
1732 int SMESH_Mesh::NbBiQuadQuadrangles() const throw(SALOME_Exception)
1733 {
1734   Unexpect aCatch(SalomeException);
1735   return _myMeshDS->GetMeshInfo().NbBiQuadQuadrangles();
1736 }
1737
1738 //================================================================================
1739 /*!
1740  * \brief Return the number of polygonal faces in the mesh
1741  */
1742 //================================================================================
1743
1744 int SMESH_Mesh::NbPolygons() const throw(SALOME_Exception)
1745 {
1746   Unexpect aCatch(SalomeException);
1747   return _myMeshDS->GetMeshInfo().NbPolygons();
1748 }
1749
1750 //================================================================================
1751 /*!
1752  * \brief Return number of volumes of given order in the mesh
1753  */
1754 //================================================================================
1755
1756 int SMESH_Mesh::NbVolumes(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1757 {
1758   Unexpect aCatch(SalomeException);
1759   return _myMeshDS->GetMeshInfo().NbVolumes(order);
1760 }
1761
1762 //================================================================================
1763 /*!
1764  * \brief  Return number of tetrahedrons of given order in the mesh
1765  */
1766 //================================================================================
1767
1768 int SMESH_Mesh::NbTetras(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1769 {
1770   Unexpect aCatch(SalomeException);
1771   return _myMeshDS->GetMeshInfo().NbTetras(order);
1772 }
1773
1774 //================================================================================
1775 /*!
1776  * \brief  Return number of hexahedrons of given order in the mesh
1777  */
1778 //================================================================================
1779
1780 int SMESH_Mesh::NbHexas(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1781 {
1782   Unexpect aCatch(SalomeException);
1783   return _myMeshDS->GetMeshInfo().NbHexas(order);
1784 }
1785
1786 //================================================================================
1787 /*!
1788  * \brief  Return number of triquadratic hexahedrons in the mesh
1789  */
1790 //================================================================================
1791
1792 int SMESH_Mesh::NbTriQuadraticHexas() const throw(SALOME_Exception)
1793 {
1794   Unexpect aCatch(SalomeException);
1795   return _myMeshDS->GetMeshInfo().NbTriQuadHexas();
1796 }
1797
1798 //================================================================================
1799 /*!
1800  * \brief  Return number of pyramids of given order in the mesh
1801  */
1802 //================================================================================
1803
1804 int SMESH_Mesh::NbPyramids(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1805 {
1806   Unexpect aCatch(SalomeException);
1807   return _myMeshDS->GetMeshInfo().NbPyramids(order);
1808 }
1809
1810 //================================================================================
1811 /*!
1812  * \brief  Return number of prisms (penthahedrons) of given order in the mesh
1813  */
1814 //================================================================================
1815
1816 int SMESH_Mesh::NbPrisms(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1817 {
1818   Unexpect aCatch(SalomeException);
1819   return _myMeshDS->GetMeshInfo().NbPrisms(order);
1820 }
1821
1822 //================================================================================
1823 /*!
1824  * \brief  Return number of hexagonal prisms in the mesh
1825  */
1826 //================================================================================
1827
1828 int SMESH_Mesh::NbHexagonalPrisms() const throw(SALOME_Exception)
1829 {
1830   Unexpect aCatch(SalomeException);
1831   return _myMeshDS->GetMeshInfo().NbHexPrisms();
1832 }
1833
1834 //================================================================================
1835 /*!
1836  * \brief  Return number of polyhedrons in the mesh
1837  */
1838 //================================================================================
1839
1840 int SMESH_Mesh::NbPolyhedrons() const throw(SALOME_Exception)
1841 {
1842   Unexpect aCatch(SalomeException);
1843   return _myMeshDS->GetMeshInfo().NbPolyhedrons();
1844 }
1845
1846 //================================================================================
1847 /*!
1848  * \brief  Return number of ball elements in the mesh
1849  */
1850 //================================================================================
1851
1852 int SMESH_Mesh::NbBalls() const throw(SALOME_Exception)
1853 {
1854   Unexpect aCatch(SalomeException);
1855   return _myMeshDS->GetMeshInfo().NbBalls();
1856 }
1857
1858 //================================================================================
1859 /*!
1860  * \brief  Return number of submeshes in the mesh
1861  */
1862 //================================================================================
1863
1864 int SMESH_Mesh::NbSubMesh() const throw(SALOME_Exception)
1865 {
1866   Unexpect aCatch(SalomeException);
1867   return _myMeshDS->NbSubMesh();
1868 }
1869
1870 //=======================================================================
1871 //function : IsNotConformAllowed
1872 //purpose  : check if a hypothesis alowing notconform mesh is present
1873 //=======================================================================
1874
1875 bool SMESH_Mesh::IsNotConformAllowed() const
1876 {
1877   if(MYDEBUG) MESSAGE("SMESH_Mesh::IsNotConformAllowed");
1878
1879   static SMESH_HypoFilter filter( SMESH_HypoFilter::HasName( "NotConformAllowed" ));
1880   return GetHypothesis( _myMeshDS->ShapeToMesh(), filter, false );
1881 }
1882
1883 //=======================================================================
1884 //function : IsMainShape
1885 //purpose  : 
1886 //=======================================================================
1887
1888 bool SMESH_Mesh::IsMainShape(const TopoDS_Shape& theShape) const
1889 {
1890   return theShape.IsSame(_myMeshDS->ShapeToMesh() );
1891 }
1892
1893 //=============================================================================
1894 /*!
1895  *  
1896  */
1897 //=============================================================================
1898
1899 SMESH_Group* SMESH_Mesh::AddGroup (const SMDSAbs_ElementType theType,
1900                                    const char*               theName,
1901                                    int&                      theId,
1902                                    const TopoDS_Shape&       theShape,
1903                                    const SMESH_PredicatePtr& thePredicate)
1904 {
1905   if (_mapGroup.count(_groupId))
1906     return NULL;
1907   theId = _groupId;
1908   SMESH_Group* aGroup = new SMESH_Group (theId, this, theType, theName, theShape, thePredicate);
1909   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1910   _mapGroup[_groupId++] = aGroup;
1911   return aGroup;
1912 }
1913
1914 //================================================================================
1915 /*!
1916  * \brief Creates a group based on an existing SMESHDS group. Group ID should be unique
1917  */
1918 //================================================================================
1919
1920 SMESH_Group* SMESH_Mesh::AddGroup (SMESHDS_GroupBase* groupDS) throw(SALOME_Exception)
1921 {
1922   if ( !groupDS ) 
1923     throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup(): NULL SMESHDS_GroupBase"));
1924
1925   map <int, SMESH_Group*>::iterator i_g = _mapGroup.find( groupDS->GetID() );
1926   if ( i_g != _mapGroup.end() && i_g->second )
1927   {
1928     if ( i_g->second->GetGroupDS() == groupDS )
1929       return i_g->second;
1930     else
1931       throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup() wrong ID of SMESHDS_GroupBase"));
1932   }
1933   SMESH_Group* aGroup = new SMESH_Group (groupDS);
1934   _mapGroup[ groupDS->GetID() ] = aGroup;
1935   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1936
1937   _groupId = 1 + _mapGroup.rbegin()->first;
1938
1939   return aGroup;
1940 }
1941
1942
1943 //================================================================================
1944 /*!
1945  * \brief Creates SMESH_Groups for not wrapped SMESHDS_Groups
1946  *  \retval bool - true if new SMESH_Groups have been created
1947  * 
1948  */
1949 //================================================================================
1950
1951 bool SMESH_Mesh::SynchronizeGroups()
1952 {
1953   int nbGroups = _mapGroup.size();
1954   const set<SMESHDS_GroupBase*>& groups = _myMeshDS->GetGroups();
1955   set<SMESHDS_GroupBase*>::const_iterator gIt = groups.begin();
1956   for ( ; gIt != groups.end(); ++gIt )
1957   {
1958     SMESHDS_GroupBase* groupDS = (SMESHDS_GroupBase*) *gIt;
1959     _groupId = groupDS->GetID();
1960     if ( !_mapGroup.count( _groupId ))
1961       _mapGroup[_groupId] = new SMESH_Group( groupDS );
1962   }
1963   if ( !_mapGroup.empty() )
1964     _groupId = _mapGroup.rbegin()->first + 1;
1965
1966   return nbGroups < _mapGroup.size();
1967 }
1968
1969 //================================================================================
1970 /*!
1971  * \brief Return iterator on all existing groups
1972  */
1973 //================================================================================
1974
1975 SMESH_Mesh::GroupIteratorPtr SMESH_Mesh::GetGroups() const
1976 {
1977   typedef map <int, SMESH_Group *> TMap;
1978   return GroupIteratorPtr( new SMDS_mapIterator<TMap>( _mapGroup ));
1979 }
1980
1981 //=============================================================================
1982 /*!
1983  * \brief Return a group by ID
1984  */
1985 //=============================================================================
1986
1987 SMESH_Group* SMESH_Mesh::GetGroup (const int theGroupID)
1988 {
1989   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1990     return NULL;
1991   return _mapGroup[theGroupID];
1992 }
1993
1994
1995 //=============================================================================
1996 /*!
1997  * \brief Return IDs of all groups
1998  */
1999 //=============================================================================
2000
2001 list<int> SMESH_Mesh::GetGroupIds() const
2002 {
2003   list<int> anIds;
2004   for ( map<int, SMESH_Group*>::const_iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
2005     anIds.push_back( it->first );
2006   
2007   return anIds;
2008 }
2009
2010 //================================================================================
2011 /*!
2012  * \brief Set a caller of methods at level of CORBA API implementation.
2013  * The set upCaller will be deleted by SMESH_Mesh
2014  */
2015 //================================================================================
2016
2017 void SMESH_Mesh::SetCallUp( TCallUp* upCaller )
2018 {
2019   if ( _callUp ) delete _callUp;
2020   _callUp = upCaller;
2021 }
2022
2023 //=============================================================================
2024 /*!
2025  *  
2026  */
2027 //=============================================================================
2028
2029 bool SMESH_Mesh::RemoveGroup (const int theGroupID)
2030 {
2031   if (_mapGroup.find(theGroupID) == _mapGroup.end())
2032     return false;
2033   GetMeshDS()->RemoveGroup( _mapGroup[theGroupID]->GetGroupDS() );
2034   delete _mapGroup[theGroupID];
2035   _mapGroup.erase (theGroupID);
2036   if (_callUp)
2037     _callUp->RemoveGroup( theGroupID );
2038   return true;
2039 }
2040
2041 //=======================================================================
2042 //function : GetAncestors
2043 //purpose  : return list of ancestors of theSubShape in the order
2044 //           that lower dimention shapes come first.
2045 //=======================================================================
2046
2047 const TopTools_ListOfShape& SMESH_Mesh::GetAncestors(const TopoDS_Shape& theS) const
2048 {
2049   if ( _mapAncestors.Contains( theS ) )
2050     return _mapAncestors.FindFromKey( theS );
2051
2052   static TopTools_ListOfShape emptyList;
2053   return emptyList;
2054 }
2055
2056 //=======================================================================
2057 //function : Dump
2058 //purpose  : dumps contents of mesh to stream [ debug purposes ]
2059 //=======================================================================
2060
2061 ostream& SMESH_Mesh::Dump(ostream& save)
2062 {
2063   int clause = 0;
2064   save << "========================== Dump contents of mesh ==========================" << endl << endl;
2065   save << ++clause << ") Total number of nodes:   \t"    << NbNodes() << endl;
2066   save << ++clause << ") Total number of edges:   \t"    << NbEdges() << endl;
2067   save << ++clause << ") Total number of faces:   \t"    << NbFaces() << endl;
2068   save << ++clause << ") Total number of polygons:\t"    << NbPolygons() << endl;
2069   save << ++clause << ") Total number of volumes:\t"     << NbVolumes() << endl;
2070   save << ++clause << ") Total number of polyhedrons:\t" << NbPolyhedrons() << endl << endl;
2071   for ( int isQuadratic = 0; isQuadratic < 2; ++isQuadratic )
2072   {
2073     string orderStr = isQuadratic ? "quadratic" : "linear";
2074     SMDSAbs_ElementOrder order  = isQuadratic ? ORDER_QUADRATIC : ORDER_LINEAR;
2075
2076     save << ++clause << ") Total number of " << orderStr << " edges:\t" << NbEdges(order) << endl;
2077     save << ++clause << ") Total number of " << orderStr << " faces:\t" << NbFaces(order) << endl;
2078     if ( NbFaces(order) > 0 ) {
2079       int nb3 = NbTriangles(order);
2080       int nb4 = NbQuadrangles(order);
2081       save << clause << ".1) Number of " << orderStr << " triangles:  \t" << nb3 << endl;
2082       save << clause << ".2) Number of " << orderStr << " quadrangles:\t" << nb4 << endl;
2083       if ( nb3 + nb4 !=  NbFaces(order) ) {
2084         map<int,int> myFaceMap;
2085         SMDS_FaceIteratorPtr itFaces=_myMeshDS->facesIterator();
2086         while( itFaces->more( ) ) {
2087           int nbNodes = itFaces->next()->NbNodes();
2088           if ( myFaceMap.find( nbNodes ) == myFaceMap.end() )
2089             myFaceMap[ nbNodes ] = 0;
2090           myFaceMap[ nbNodes ] = myFaceMap[ nbNodes ] + 1;
2091         }
2092         save << clause << ".3) Faces in detail: " << endl;
2093         map <int,int>::iterator itF;
2094         for (itF = myFaceMap.begin(); itF != myFaceMap.end(); itF++)
2095           save << "--> nb nodes: " << itF->first << " - nb elemens:\t" << itF->second << endl;
2096       }
2097     }
2098     save << ++clause << ") Total number of " << orderStr << " volumes:\t" << NbVolumes(order) << endl;
2099     if ( NbVolumes(order) > 0 ) {
2100       int nb8 = NbHexas(order);
2101       int nb4 = NbTetras(order);
2102       int nb5 = NbPyramids(order);
2103       int nb6 = NbPrisms(order);
2104       save << clause << ".1) Number of " << orderStr << " hexahedrons:\t" << nb8 << endl;
2105       save << clause << ".2) Number of " << orderStr << " tetrahedrons:\t" << nb4 << endl;
2106       save << clause << ".3) Number of " << orderStr << " prisms:      \t" << nb6 << endl;
2107       save << clause << ".4) Number of " << orderStr << " pyramids:\t" << nb5 << endl;
2108       if ( nb8 + nb4 + nb5 + nb6 != NbVolumes(order) ) {
2109         map<int,int> myVolumesMap;
2110         SMDS_VolumeIteratorPtr itVolumes=_myMeshDS->volumesIterator();
2111         while( itVolumes->more( ) ) {
2112           int nbNodes = itVolumes->next()->NbNodes();
2113           if ( myVolumesMap.find( nbNodes ) == myVolumesMap.end() )
2114             myVolumesMap[ nbNodes ] = 0;
2115           myVolumesMap[ nbNodes ] = myVolumesMap[ nbNodes ] + 1;
2116         }
2117         save << clause << ".5) Volumes in detail: " << endl;
2118         map <int,int>::iterator itV;
2119         for (itV = myVolumesMap.begin(); itV != myVolumesMap.end(); itV++)
2120           save << "--> nb nodes: " << itV->first << " - nb elemens:\t" << itV->second << endl;
2121       }
2122     }
2123     save << endl;
2124   }
2125   save << "===========================================================================" << endl;
2126   return save;
2127 }
2128
2129 //=======================================================================
2130 //function : GetElementType
2131 //purpose  : Returns type of mesh element with certain id
2132 //=======================================================================
2133
2134 SMDSAbs_ElementType SMESH_Mesh::GetElementType( const int id, const bool iselem )
2135 {
2136   return _myMeshDS->GetElementType( id, iselem );
2137 }
2138
2139 //=============================================================================
2140 /*!
2141  *  \brief Convert group on geometry into standalone group
2142  */
2143 //=============================================================================
2144
2145 SMESH_Group* SMESH_Mesh::ConvertToStandalone ( int theGroupID )
2146 {
2147   SMESH_Group* aGroup = 0;
2148   map < int, SMESH_Group * >::iterator itg = _mapGroup.find( theGroupID );
2149   if ( itg == _mapGroup.end() )
2150     return aGroup;
2151
2152   SMESH_Group* anOldGrp = (*itg).second;
2153   SMESHDS_GroupBase* anOldGrpDS = anOldGrp->GetGroupDS();
2154   if ( !anOldGrp || !anOldGrpDS )
2155     return aGroup;
2156
2157   // create new standalone group
2158   aGroup = new SMESH_Group (theGroupID, this, anOldGrpDS->GetType(), anOldGrp->GetName() );
2159   _mapGroup[theGroupID] = aGroup;
2160
2161   SMESHDS_Group* aNewGrpDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
2162   GetMeshDS()->RemoveGroup( anOldGrpDS );
2163   GetMeshDS()->AddGroup( aNewGrpDS );
2164
2165   // add elements (or nodes) into new created group
2166   SMDS_ElemIteratorPtr anItr = anOldGrpDS->GetElements();
2167   while ( anItr->more() )
2168     aNewGrpDS->Add( (anItr->next())->GetID() );
2169
2170   // set color
2171   aNewGrpDS->SetColor( anOldGrpDS->GetColor() );
2172
2173   // remove old group
2174   delete anOldGrp;
2175
2176   return aGroup;
2177 }
2178
2179 //=============================================================================
2180 /*!
2181  *  \brief remove submesh order  from Mesh
2182  */
2183 //=============================================================================
2184
2185 void SMESH_Mesh::ClearMeshOrder()
2186 {
2187   _mySubMeshOrder.clear();
2188 }
2189
2190 //=============================================================================
2191 /*!
2192  *  \brief remove submesh order  from Mesh
2193  */
2194 //=============================================================================
2195
2196 void SMESH_Mesh::SetMeshOrder(const TListOfListOfInt& theOrder )
2197 {
2198   _mySubMeshOrder = theOrder;
2199 }
2200
2201 //=============================================================================
2202 /*!
2203  *  \brief return submesh order if any
2204  */
2205 //=============================================================================
2206
2207 const TListOfListOfInt& SMESH_Mesh::GetMeshOrder() const
2208 {
2209   return _mySubMeshOrder;
2210 }
2211
2212 //=============================================================================
2213 /*!
2214  *  \brief fill _mapAncestors
2215  */
2216 //=============================================================================
2217
2218 void SMESH_Mesh::fillAncestorsMap(const TopoDS_Shape& theShape)
2219 {
2220   int desType, ancType;
2221   if ( !theShape.IsSame( GetShapeToMesh()) && theShape.ShapeType() == TopAbs_COMPOUND )
2222   {
2223     // a geom group is added. Insert it into lists of ancestors before
2224     // the first ancestor more complex than group members
2225     TopoDS_Iterator subIt( theShape );
2226     if ( !subIt.More() ) return;
2227     int memberType = subIt.Value().ShapeType();
2228     for ( desType = TopAbs_VERTEX; desType >= memberType; desType-- )
2229       for (TopExp_Explorer des( theShape, TopAbs_ShapeEnum( desType )); des.More(); des.Next())
2230       {
2231         if ( !_mapAncestors.Contains( des.Current() )) continue;// issue 0020982
2232         TopTools_ListOfShape& ancList = _mapAncestors.ChangeFromKey( des.Current() );
2233         TopTools_ListIteratorOfListOfShape ancIt (ancList);
2234         while ( ancIt.More() && ancIt.Value().ShapeType() >= memberType )
2235           ancIt.Next();
2236         if ( ancIt.More() )
2237           ancList.InsertBefore( theShape, ancIt );
2238       }
2239   }
2240   else // else added for 52457: Addition of hypotheses is 8 time longer than meshing
2241   {
2242     for ( desType = TopAbs_VERTEX; desType > TopAbs_COMPOUND; desType-- )
2243       for ( ancType = desType - 1; ancType >= TopAbs_COMPOUND; ancType-- )
2244         TopExp::MapShapesAndAncestors ( theShape,
2245                                         (TopAbs_ShapeEnum) desType,
2246                                         (TopAbs_ShapeEnum) ancType,
2247                                         _mapAncestors );
2248   }
2249   // visit COMPOUNDs inside a COMPOUND that are not reachable by TopExp_Explorer
2250   if ( theShape.ShapeType() == TopAbs_COMPOUND )
2251   {
2252     TopoDS_Iterator sIt(theShape);
2253     if ( sIt.More() && sIt.Value().ShapeType() == TopAbs_COMPOUND )
2254       for ( ; sIt.More(); sIt.Next() )
2255         if ( sIt.Value().ShapeType() == TopAbs_COMPOUND )
2256           fillAncestorsMap( sIt.Value() );
2257   }
2258 }
2259
2260 //=============================================================================
2261 /*!
2262  * \brief sort submeshes according to stored mesh order
2263  * \param theListToSort in out list to be sorted
2264  * \return FALSE if nothing sorted
2265  */
2266 //=============================================================================
2267
2268  bool SMESH_Mesh::SortByMeshOrder(std::vector<SMESH_subMesh*>& theListToSort) const
2269 {
2270   if ( !_mySubMeshOrder.size() || theListToSort.size() < 2)
2271     return true;
2272   
2273   bool res = false;
2274   vector<SMESH_subMesh*> onlyOrderedList;
2275   // collect all ordered submeshes in one list as pointers
2276   // and get their positions within theListToSort
2277   typedef vector<SMESH_subMesh*>::iterator TPosInList;
2278   map< int, TPosInList > sortedPos;
2279   TPosInList smBeg = theListToSort.begin(), smEnd = theListToSort.end();
2280   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2281   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2282   {
2283     const TListOfInt& listOfId = *listIdsIt;
2284     // convert sm ids to sm's
2285     vector<SMESH_subMesh*> smVec;
2286     TListOfInt::const_iterator idIt = listOfId.begin();
2287     for ( ; idIt != listOfId.end(); idIt++ ) {
2288       if ( SMESH_subMesh * sm = GetSubMeshContaining( *idIt )) {
2289         if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->IsComplexSubmesh() )
2290         {
2291           SMESHDS_SubMeshIteratorPtr smdsIt = sm->GetSubMeshDS()->GetSubMeshIterator();
2292           while ( smdsIt->more() )
2293           {
2294             const SMESHDS_SubMesh* smDS = smdsIt->next();
2295             if (( sm = GetSubMeshContaining( smDS->GetID() )))
2296               smVec.push_back( sm );
2297           }
2298         }
2299         else
2300         {
2301           smVec.push_back( sm );
2302         }
2303       }
2304     }
2305     // find smVec items in theListToSort
2306     for ( size_t i = 0; i < smVec.size(); ++i )
2307     {
2308       TPosInList smPos = find( smBeg, smEnd, smVec[i] );
2309       if ( smPos != smEnd ) {
2310         onlyOrderedList.push_back( smVec[i] );
2311         sortedPos[ distance( smBeg, smPos )] = smPos;
2312       }
2313     }
2314   }
2315   if (onlyOrderedList.size() < 2)
2316     return res;
2317   res = true;
2318
2319   vector<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
2320   vector<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
2321
2322   // iterate on ordered submeshes and insert them in detected positions
2323   map< int, TPosInList >::iterator i_pos = sortedPos.begin();
2324   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
2325     *(i_pos->second) = *onlyBIt;
2326
2327   return res;
2328 }
2329
2330 //================================================================================
2331 /*!
2332  * \brief Return true if given order of sub-meshes is OK
2333  */
2334 //================================================================================
2335
2336 bool SMESH_Mesh::IsOrderOK( const SMESH_subMesh* smBefore,
2337                             const SMESH_subMesh* smAfter ) const
2338 {
2339   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2340   TListOfInt::const_iterator idBef, idAft;
2341   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2342   {
2343     const TListOfInt& listOfId = *listIdsIt;
2344     idBef = std::find( listOfId.begin(), listOfId.end(), smBefore->GetId() );
2345     if ( idBef != listOfId.end() )
2346       idAft = std::find( listOfId.begin(), listOfId.end(), smAfter->GetId() );
2347     if ( idAft != listOfId.end () )
2348       return ( std::distance( listOfId.begin(), idBef ) <
2349                std::distance( listOfId.begin(), idAft )   );
2350   }
2351   return true; // no order imposed to given submeshes
2352
2353
2354 //=============================================================================
2355 /*!
2356  * \brief sort submeshes according to stored mesh order
2357  * \param theListToSort in out list to be sorted
2358  * \return FALSE if nothing sorted
2359  */
2360 //=============================================================================
2361
2362 void SMESH_Mesh::getAncestorsSubMeshes (const TopoDS_Shape&            theSubShape,
2363                                         std::vector< SMESH_subMesh* >& theSubMeshes) const
2364 {
2365   theSubMeshes.clear();
2366   TopTools_ListIteratorOfListOfShape it( GetAncestors( theSubShape ));
2367   for (; it.More(); it.Next() )
2368     if ( SMESH_subMesh* sm = GetSubMeshContaining( it.Value() ))
2369       theSubMeshes.push_back(sm);
2370
2371   // sort submeshes according to stored mesh order
2372   SortByMeshOrder( theSubMeshes );
2373 }