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