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