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