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