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