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