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