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