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