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