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