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