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