Salome HOME
0022100: EDF 2413 SMESH: Take into account TRIA7
[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) const
817 {
818   set<string> hypTypes; // to exclude same type hypos from the result list
819   int nbHyps = 0;
820
821   // only one main hypothesis is allowed
822   bool mainHypFound = false;
823
824   // fill in hypTypes
825   list<const SMESHDS_Hypothesis*>::const_iterator hyp;
826   for ( hyp = aHypList.begin(); hyp != aHypList.end(); hyp++ ) {
827     if ( hypTypes.insert( (*hyp)->GetName() ).second )
828       nbHyps++;
829     if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
830       mainHypFound = true;
831   }
832
833   // get hypos from aSubShape
834   {
835     const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(aSubShape);
836     for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
837       if ( aFilter.IsOk (cSMESH_Hyp( *hyp ), aSubShape) &&
838            ( cSMESH_Hyp(*hyp)->IsAuxiliary() || !mainHypFound ) &&
839            hypTypes.insert( (*hyp)->GetName() ).second )
840       {
841         aHypList.push_back( *hyp );
842         nbHyps++;
843         if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
844           mainHypFound = true;
845       }
846   }
847
848   // get hypos from ancestors of aSubShape
849   if ( andAncestors )
850   {
851     TopTools_MapOfShape map;
852
853     // user sorted submeshes of ancestors, according to stored submesh priority
854     const list<SMESH_subMesh*> smList = getAncestorsSubMeshes( aSubShape );
855     list<SMESH_subMesh*>::const_iterator smIt = smList.begin(); 
856     for ( ; smIt != smList.end(); smIt++ )
857     {
858       const TopoDS_Shape& curSh = (*smIt)->GetSubShape();
859      if ( !map.Add( curSh ))
860         continue;
861       const list<const SMESHDS_Hypothesis*>& hypList = _myMeshDS->GetHypothesis(curSh);
862       for ( hyp = hypList.begin(); hyp != hypList.end(); hyp++ )
863         if (aFilter.IsOk( cSMESH_Hyp( *hyp ), curSh ) &&
864             ( cSMESH_Hyp(*hyp)->IsAuxiliary() || !mainHypFound ) &&
865             hypTypes.insert( (*hyp)->GetName() ).second )
866         {
867           aHypList.push_back( *hyp );
868           nbHyps++;
869           if ( !cSMESH_Hyp(*hyp)->IsAuxiliary() )
870             mainHypFound = true;
871         }
872     }
873   }
874   return nbHyps;
875 }
876
877 //=============================================================================
878 /*!
879  * 
880  */
881 //=============================================================================
882
883 const list<SMESHDS_Command*> & SMESH_Mesh::GetLog() throw(SALOME_Exception)
884 {
885   Unexpect aCatch(SalomeException);
886   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetLog");
887   return _myMeshDS->GetScript()->GetCommands();
888 }
889
890 //=============================================================================
891 /*!
892  * 
893  */
894 //=============================================================================
895 void SMESH_Mesh::ClearLog() throw(SALOME_Exception)
896 {
897   Unexpect aCatch(SalomeException);
898   if(MYDEBUG) MESSAGE("SMESH_Mesh::ClearLog");
899   _myMeshDS->GetScript()->Clear();
900 }
901
902 //=============================================================================
903 /*!
904  * Get or Create the SMESH_subMesh object implementation
905  */
906 //=============================================================================
907
908 SMESH_subMesh *SMESH_Mesh::GetSubMesh(const TopoDS_Shape & aSubShape)
909   throw(SALOME_Exception)
910 {
911   Unexpect aCatch(SalomeException);
912   SMESH_subMesh *aSubMesh;
913   int index = _myMeshDS->ShapeToIndex(aSubShape);
914
915   // for submeshes on GEOM Group
916   if (( !index || index > _nbSubShapes ) && aSubShape.ShapeType() == TopAbs_COMPOUND ) {
917     TopoDS_Iterator it( aSubShape );
918     if ( it.More() )
919     {
920       index = _myMeshDS->AddCompoundSubmesh( aSubShape, it.Value().ShapeType() );
921       if ( index > _nbSubShapes ) _nbSubShapes = index; // not to create sm for this group again
922
923       // fill map of Ancestors
924       fillAncestorsMap(aSubShape);
925     }
926   }
927 //   if ( !index )
928 //     return NULL; // neither sub-shape nor a group
929
930   map <int, SMESH_subMesh *>::iterator i_sm = _mapSubMesh.find(index);
931   if ( i_sm != _mapSubMesh.end())
932   {
933     aSubMesh = i_sm->second;
934   }
935   else
936   {
937     aSubMesh = new SMESH_subMesh(index, this, _myMeshDS, aSubShape);
938     _mapSubMesh[index] = aSubMesh;
939   }
940   return aSubMesh;
941 }
942
943 //=============================================================================
944 /*!
945  * Get the SMESH_subMesh object implementation. Dont create it, return null
946  * if it does not exist.
947  */
948 //=============================================================================
949
950 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const TopoDS_Shape & aSubShape) const
951   throw(SALOME_Exception)
952 {
953   Unexpect aCatch(SalomeException);
954   SMESH_subMesh *aSubMesh = NULL;
955   
956   int index = _myMeshDS->ShapeToIndex(aSubShape);
957
958   map <int, SMESH_subMesh *>::const_iterator i_sm = _mapSubMesh.find(index);
959   if ( i_sm != _mapSubMesh.end())
960     aSubMesh = i_sm->second;
961
962   return aSubMesh;
963 }
964 //=============================================================================
965 /*!
966  * Get the SMESH_subMesh object implementation. Dont create it, return null
967  * if it does not exist.
968  */
969 //=============================================================================
970
971 SMESH_subMesh *SMESH_Mesh::GetSubMeshContaining(const int aShapeID) const
972 throw(SALOME_Exception)
973 {
974   Unexpect aCatch(SalomeException);
975   
976   map <int, SMESH_subMesh *>::const_iterator i_sm = _mapSubMesh.find(aShapeID);
977   if (i_sm == _mapSubMesh.end())
978     return NULL;
979   return i_sm->second;
980 }
981 //================================================================================
982 /*!
983  * \brief Return submeshes of groups containing the given sub-shape
984  */
985 //================================================================================
986
987 list<SMESH_subMesh*>
988 SMESH_Mesh::GetGroupSubMeshesContaining(const TopoDS_Shape & aSubShape) const
989   throw(SALOME_Exception)
990 {
991   Unexpect aCatch(SalomeException);
992   list<SMESH_subMesh*> found;
993
994   SMESH_subMesh * subMesh = GetSubMeshContaining(aSubShape);
995   if ( !subMesh )
996     return found;
997
998   // submeshes of groups have max IDs, so search from the map end
999   map<int, SMESH_subMesh *>::const_reverse_iterator i_sm;
1000   for ( i_sm = _mapSubMesh.rbegin(); i_sm != _mapSubMesh.rend(); ++i_sm) {
1001     SMESHDS_SubMesh * ds = i_sm->second->GetSubMeshDS();
1002     if ( ds && ds->IsComplexSubmesh() ) {
1003       if ( SMESH_MesherHelper::IsSubShape( aSubShape, i_sm->second->GetSubShape() ))
1004       {
1005         found.push_back( i_sm->second );
1006         //break;
1007       }
1008     } else {
1009       break; // the rest sub-meshes are not those of groups
1010     }
1011   }
1012
1013   if ( found.empty() ) // maybe the main shape is a COMPOUND (issue 0021530)
1014   {
1015     if ( SMESH_subMesh * mainSM = GetSubMeshContaining(1))
1016       if ( mainSM->GetSubShape().ShapeType() == TopAbs_COMPOUND )
1017       {
1018         TopoDS_Iterator it( mainSM->GetSubShape() );
1019         if ( it.Value().ShapeType() == aSubShape.ShapeType() &&
1020              SMESH_MesherHelper::IsSubShape( aSubShape, mainSM->GetSubShape() ))
1021           found.push_back( mainSM );
1022       }
1023   }
1024   return found;
1025 }
1026 //=======================================================================
1027 //function : IsUsedHypothesis
1028 //purpose  : Return True if anHyp is used to mesh aSubShape
1029 //=======================================================================
1030
1031 bool SMESH_Mesh::IsUsedHypothesis(SMESHDS_Hypothesis * anHyp,
1032                                   const SMESH_subMesh* aSubMesh)
1033 {
1034   SMESH_Hypothesis* hyp = static_cast<SMESH_Hypothesis*>(anHyp);
1035
1036   // check if anHyp can be used to mesh aSubMesh
1037   if ( !aSubMesh || !aSubMesh->IsApplicableHypotesis( hyp ))
1038     return false;
1039
1040   const TopoDS_Shape & aSubShape = const_cast<SMESH_subMesh*>( aSubMesh )->GetSubShape();
1041
1042   SMESH_Algo *algo = _gen->GetAlgo(*this, aSubShape );
1043
1044   // algorithm
1045   if (anHyp->GetType() > SMESHDS_Hypothesis::PARAM_ALGO)
1046     return ( anHyp == algo );
1047
1048   // algorithm parameter
1049   if (algo)
1050   {
1051     // look trough hypotheses used by algo
1052     SMESH_HypoFilter hypoKind;
1053     if ( algo->InitCompatibleHypoFilter( hypoKind, !hyp->IsAuxiliary() )) {
1054       list <const SMESHDS_Hypothesis * > usedHyps;
1055       if ( GetHypotheses( aSubShape, hypoKind, usedHyps, true ))
1056         return ( find( usedHyps.begin(), usedHyps.end(), anHyp ) != usedHyps.end() );
1057     }
1058   }
1059
1060   // look through all assigned hypotheses
1061   //SMESH_HypoFilter filter( SMESH_HypoFilter::Is( hyp ));
1062   return false; //GetHypothesis( aSubShape, filter, true );
1063 }
1064
1065 //=============================================================================
1066 /*!
1067  *
1068  */
1069 //=============================================================================
1070
1071 const list < SMESH_subMesh * >&
1072 SMESH_Mesh::GetSubMeshUsingHypothesis(SMESHDS_Hypothesis * anHyp)
1073   throw(SALOME_Exception)
1074 {
1075   Unexpect aCatch(SalomeException);
1076   if(MYDEBUG) MESSAGE("SMESH_Mesh::GetSubMeshUsingHypothesis");
1077   map < int, SMESH_subMesh * >::iterator itsm;
1078   _subMeshesUsingHypothesisList.clear();
1079   for (itsm = _mapSubMesh.begin(); itsm != _mapSubMesh.end(); itsm++)
1080   {
1081     SMESH_subMesh *aSubMesh = (*itsm).second;
1082     if ( IsUsedHypothesis ( anHyp, aSubMesh ))
1083       _subMeshesUsingHypothesisList.push_back(aSubMesh);
1084   }
1085   return _subMeshesUsingHypothesisList;
1086 }
1087
1088 //=======================================================================
1089 //function : NotifySubMeshesHypothesisModification
1090 //purpose  : Say all submeshes using theChangedHyp that it has been modified
1091 //=======================================================================
1092
1093 void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* hyp)
1094 {
1095   Unexpect aCatch(SalomeException);
1096
1097   if ( !GetMeshDS()->IsUsedHypothesis( hyp ))
1098     return;
1099
1100   if (_callUp)
1101     _callUp->HypothesisModified();
1102
1103   const SMESH_Algo *foundAlgo = 0;
1104   SMESH_HypoFilter algoKind, compatibleHypoKind;
1105   list <const SMESHDS_Hypothesis * > usedHyps;
1106
1107
1108   map < int, SMESH_subMesh * >::iterator itsm;
1109   for (itsm = _mapSubMesh.begin(); itsm != _mapSubMesh.end(); itsm++)
1110   {
1111     SMESH_subMesh *aSubMesh = (*itsm).second;
1112     if ( aSubMesh->IsApplicableHypotesis( hyp ))
1113     {
1114       const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
1115
1116       if ( !foundAlgo ) // init filter for algo search
1117         algoKind.Init( THypType::IsAlgo() ).And( THypType::IsApplicableTo( aSubShape ));
1118       
1119       const SMESH_Algo *algo = static_cast<const SMESH_Algo*>
1120         ( GetHypothesis( aSubShape, algoKind, true ));
1121
1122       if ( algo )
1123       {
1124         bool sameAlgo = ( algo == foundAlgo );
1125         if ( !sameAlgo && foundAlgo )
1126           sameAlgo = ( strcmp( algo->GetName(), foundAlgo->GetName() ) == 0);
1127
1128         if ( !sameAlgo ) { // init filter for used hypos search
1129           if ( !algo->InitCompatibleHypoFilter( compatibleHypoKind, !hyp->IsAuxiliary() ))
1130             continue; // algo does not use any hypothesis
1131           foundAlgo = algo;
1132         }
1133
1134         // check if hyp is used by algo
1135         usedHyps.clear();
1136         if ( GetHypotheses( aSubShape, compatibleHypoKind, usedHyps, true ) &&
1137              find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() )
1138         {
1139           aSubMesh->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
1140                                     const_cast< SMESH_Hypothesis*>( hyp ));
1141         }
1142       }
1143     }
1144   }
1145   HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
1146   GetMeshDS()->Modified();
1147 }
1148
1149 //=============================================================================
1150 /*!
1151  *  Auto color functionality
1152  */
1153 //=============================================================================
1154 void SMESH_Mesh::SetAutoColor(bool theAutoColor) throw(SALOME_Exception)
1155 {
1156   Unexpect aCatch(SalomeException);
1157   _isAutoColor = theAutoColor;
1158 }
1159
1160 bool SMESH_Mesh::GetAutoColor() throw(SALOME_Exception)
1161 {
1162   Unexpect aCatch(SalomeException);
1163   return _isAutoColor;
1164 }
1165
1166 //=======================================================================
1167 //function : SetIsModified
1168 //purpose  : Set the flag meaning that the mesh has been edited "manually"
1169 //=======================================================================
1170
1171 void SMESH_Mesh::SetIsModified(bool isModified)
1172 {
1173   _isModified = isModified;
1174
1175   if ( _isModified )
1176     // check if mesh becomes empty as result of modification
1177     HasModificationsToDiscard();
1178 }
1179
1180 //=======================================================================
1181 //function : HasModificationsToDiscard
1182 //purpose  : Return true if the mesh has been edited since a total re-compute
1183 //           and those modifications may prevent successful partial re-compute.
1184 //           As a side effect reset _isModified flag if mesh is empty
1185 //issue    : 0020693
1186 //=======================================================================
1187
1188 bool SMESH_Mesh::HasModificationsToDiscard() const
1189 {
1190   if ( ! _isModified )
1191     return false;
1192
1193   // return true if the next Compute() will be partial and
1194   // existing but changed elements may prevent successful re-compute
1195   bool hasComputed = false, hasNotComputed = false;
1196   map <int, SMESH_subMesh*>::const_iterator i_sm = _mapSubMesh.begin();
1197   for ( ; i_sm != _mapSubMesh.end() ; ++i_sm )
1198     switch ( i_sm->second->GetSubShape().ShapeType() )
1199     {
1200     case TopAbs_EDGE:
1201     case TopAbs_FACE:
1202     case TopAbs_SOLID:
1203       if ( i_sm->second->IsMeshComputed() )
1204         hasComputed = true;
1205       else
1206         hasNotComputed = true;
1207       if ( hasComputed && hasNotComputed)
1208         return true;
1209     }
1210
1211   if ( NbNodes() < 1 )
1212     const_cast<SMESH_Mesh*>(this)->_isModified = false;
1213
1214   return false;
1215 }
1216
1217 //================================================================================
1218 /*!
1219  * \brief Check if any groups of the same type have equal names
1220  */
1221 //================================================================================
1222
1223 bool SMESH_Mesh::HasDuplicatedGroupNamesMED()
1224 {
1225   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1226   map< SMDSAbs_ElementType, set<string> > aGroupNames;
1227   for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1228   {
1229     SMESH_Group* aGroup = it->second;
1230     SMDSAbs_ElementType aType = aGroup->GetGroupDS()->GetType();
1231     string aGroupName = aGroup->GetName();
1232     aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1233     if (!aGroupNames[aType].insert(aGroupName).second)
1234       return true;
1235   }
1236
1237   return false;
1238 }
1239
1240 //================================================================================
1241 /*!
1242  * \brief Export the mesh to a med file
1243  */
1244 //================================================================================
1245
1246 void SMESH_Mesh::ExportMED(const char *        file, 
1247                            const char*         theMeshName, 
1248                            bool                theAutoGroups,
1249                            int                 theVersion,
1250                            const SMESHDS_Mesh* meshPart) 
1251   throw(SALOME_Exception)
1252 {
1253   Unexpect aCatch(SalomeException);
1254
1255   DriverMED_W_SMESHDS_Mesh myWriter;
1256   myWriter.SetFile    ( file, MED::EVersion(theVersion) );
1257   myWriter.SetMesh    ( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS   );
1258   if ( !theMeshName ) 
1259     myWriter.SetMeshId  ( _id         );
1260   else {
1261     myWriter.SetMeshId  ( -1          );
1262     myWriter.SetMeshName( theMeshName );
1263   }
1264
1265   if ( theAutoGroups ) {
1266     myWriter.AddGroupOfNodes();
1267     myWriter.AddGroupOfEdges();
1268     myWriter.AddGroupOfFaces();
1269     myWriter.AddGroupOfVolumes();
1270   }
1271
1272   // Pass groups to writer. Provide unique group names.
1273   //set<string> aGroupNames; // Corrected for Mantis issue 0020028
1274   if ( !meshPart )
1275   {
1276     map< SMDSAbs_ElementType, set<string> > aGroupNames;
1277     char aString [256];
1278     int maxNbIter = 10000; // to guarantee cycle finish
1279     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1280       SMESH_Group*       aGroup   = it->second;
1281       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1282       if ( aGroupDS ) {
1283         SMDSAbs_ElementType aType = aGroupDS->GetType();
1284         string aGroupName0 = aGroup->GetName();
1285         aGroupName0.resize(MAX_MED_GROUP_NAME_LENGTH);
1286         string aGroupName = aGroupName0;
1287         for (int i = 1; !aGroupNames[aType].insert(aGroupName).second && i < maxNbIter; i++) {
1288           sprintf(&aString[0], "GR_%d_%s", i, aGroupName0.c_str());
1289           aGroupName = aString;
1290           aGroupName.resize(MAX_MED_GROUP_NAME_LENGTH);
1291         }
1292         aGroupDS->SetStoreName( aGroupName.c_str() );
1293         myWriter.AddGroup( aGroupDS );
1294       }
1295     }
1296   }
1297   // Perform export
1298   myWriter.Perform();
1299 }
1300
1301 void SMESH_Mesh::ExportSAUV(const char *file, 
1302                             const char* theMeshName, 
1303                             bool theAutoGroups)
1304   throw(SALOME_Exception)
1305 {
1306   std::string medfilename(file);
1307   medfilename += ".med";
1308   std::string cmd;
1309 #ifdef WNT
1310   cmd = "%PYTHONBIN% ";
1311 #else
1312   cmd = "python ";
1313 #endif
1314   cmd += "-c \"";
1315   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1316   cmd += "\"";
1317   system(cmd.c_str());
1318   ExportMED(medfilename.c_str(), theMeshName, theAutoGroups, 1);
1319 #ifdef WNT
1320   cmd = "%PYTHONBIN% ";
1321 #else
1322   cmd = "python ";
1323 #endif
1324   cmd += "-c \"";
1325   cmd += "from medutilities import convert ; convert(r'" + medfilename + "', 'MED', 'GIBI', 1, r'" + file + "')";
1326   cmd += "\"";
1327   system(cmd.c_str());
1328 #ifdef WNT
1329   cmd = "%PYTHONBIN% ";
1330 #else
1331   cmd = "python ";
1332 #endif
1333   cmd += "-c \"";
1334   cmd += "from medutilities import my_remove ; my_remove(r'" + medfilename + "')";
1335   cmd += "\"";
1336   system(cmd.c_str());
1337 }
1338
1339 //================================================================================
1340 /*!
1341  * \brief Export the mesh to a DAT file
1342  */
1343 //================================================================================
1344
1345 void SMESH_Mesh::ExportDAT(const char *        file,
1346                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1347 {
1348   Unexpect aCatch(SalomeException);
1349   DriverDAT_W_SMDS_Mesh myWriter;
1350   myWriter.SetFile( file );
1351   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1352   myWriter.SetMeshId(_id);
1353   myWriter.Perform();
1354 }
1355
1356 //================================================================================
1357 /*!
1358  * \brief Export the mesh to an UNV file
1359  */
1360 //================================================================================
1361
1362 void SMESH_Mesh::ExportUNV(const char *        file,
1363                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1364 {
1365   Unexpect aCatch(SalomeException);
1366   DriverUNV_W_SMDS_Mesh myWriter;
1367   myWriter.SetFile( file );
1368   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS );
1369   myWriter.SetMeshId(_id);
1370   //  myWriter.SetGroups(_mapGroup);
1371
1372   if ( !meshPart )
1373   {
1374     for ( map<int, SMESH_Group*>::iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ ) {
1375       SMESH_Group*       aGroup   = it->second;
1376       SMESHDS_GroupBase* aGroupDS = aGroup->GetGroupDS();
1377       if ( aGroupDS ) {
1378         string aGroupName = aGroup->GetName();
1379         aGroupDS->SetStoreName( aGroupName.c_str() );
1380         myWriter.AddGroup( aGroupDS );
1381       }
1382     }
1383   }
1384   myWriter.Perform();
1385 }
1386
1387 //================================================================================
1388 /*!
1389  * \brief Export the mesh to an STL file
1390  */
1391 //================================================================================
1392
1393 void SMESH_Mesh::ExportSTL(const char *        file,
1394                            const bool          isascii,
1395                            const SMESHDS_Mesh* meshPart) throw(SALOME_Exception)
1396 {
1397   Unexpect aCatch(SalomeException);
1398   DriverSTL_W_SMDS_Mesh myWriter;
1399   myWriter.SetFile( file );
1400   myWriter.SetIsAscii( isascii );
1401   myWriter.SetMesh( meshPart ? (SMESHDS_Mesh*) meshPart : _myMeshDS);
1402   myWriter.SetMeshId(_id);
1403   myWriter.Perform();
1404 }
1405
1406 //================================================================================
1407 /*!
1408  * \brief Export the mesh to the CGNS file
1409  */
1410 //================================================================================
1411
1412 void SMESH_Mesh::ExportCGNS(const char *        file,
1413                             const SMESHDS_Mesh* meshDS)
1414 {
1415   int res = Driver_Mesh::DRS_FAIL;
1416 #ifdef WITH_CGNS
1417   DriverCGNS_Write myWriter;
1418   myWriter.SetFile( file );
1419   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1420   myWriter.SetMeshName( SMESH_Comment("Mesh_") << meshDS->GetPersistentId());
1421   res = myWriter.Perform();
1422 #endif
1423   if ( res != Driver_Mesh::DRS_OK )
1424     throw SALOME_Exception("Export failed");
1425 }
1426
1427 //================================================================================
1428 /*!
1429  * \brief Export the mesh to a GMF file
1430  */
1431 //================================================================================
1432
1433 void SMESH_Mesh::ExportGMF(const char *        file,
1434                            const SMESHDS_Mesh* meshDS,
1435                            bool                withRequiredGroups)
1436 {
1437   DriverGMF_Write myWriter;
1438   myWriter.SetFile( file );
1439   myWriter.SetMesh( const_cast<SMESHDS_Mesh*>( meshDS ));
1440   myWriter.SetExportRequiredGroups( withRequiredGroups );
1441
1442   myWriter.Perform();
1443 }
1444
1445 //================================================================================
1446 /*!
1447  * \brief Return number of nodes in the mesh
1448  */
1449 //================================================================================
1450
1451 int SMESH_Mesh::NbNodes() const throw(SALOME_Exception)
1452 {
1453   Unexpect aCatch(SalomeException);
1454   return _myMeshDS->NbNodes();
1455 }
1456
1457 //================================================================================
1458 /*!
1459  * \brief  Return number of edges of given order in the mesh
1460  */
1461 //================================================================================
1462
1463 int SMESH_Mesh::Nb0DElements() const throw(SALOME_Exception)
1464 {
1465   Unexpect aCatch(SalomeException);
1466   return _myMeshDS->GetMeshInfo().Nb0DElements();
1467 }
1468
1469 //================================================================================
1470 /*!
1471  * \brief  Return number of edges of given order in the mesh
1472  */
1473 //================================================================================
1474
1475 int SMESH_Mesh::NbEdges(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1476 {
1477   Unexpect aCatch(SalomeException);
1478   return _myMeshDS->GetMeshInfo().NbEdges(order);
1479 }
1480
1481 //================================================================================
1482 /*!
1483  * \brief Return number of faces of given order in the mesh
1484  */
1485 //================================================================================
1486
1487 int SMESH_Mesh::NbFaces(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1488 {
1489   Unexpect aCatch(SalomeException);
1490   return _myMeshDS->GetMeshInfo().NbFaces(order);
1491 }
1492
1493 //================================================================================
1494 /*!
1495  * \brief Return the number of faces in the mesh
1496  */
1497 //================================================================================
1498
1499 int SMESH_Mesh::NbTriangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1500 {
1501   Unexpect aCatch(SalomeException);
1502   return _myMeshDS->GetMeshInfo().NbTriangles(order);
1503 }
1504
1505 //================================================================================
1506 /*!
1507  * \brief Return number of biquadratic triangles in the mesh
1508  */
1509 //================================================================================
1510
1511 int SMESH_Mesh::NbBiQuadTriangles() const throw(SALOME_Exception)
1512 {
1513   Unexpect aCatch(SalomeException);
1514   return _myMeshDS->GetMeshInfo().NbBiQuadTriangles();
1515 }
1516
1517 //================================================================================
1518 /*!
1519  * \brief Return the number nodes faces in the mesh
1520  */
1521 //================================================================================
1522
1523 int SMESH_Mesh::NbQuadrangles(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1524 {
1525   Unexpect aCatch(SalomeException);
1526   return _myMeshDS->GetMeshInfo().NbQuadrangles(order);
1527 }
1528
1529 //================================================================================
1530 /*!
1531  * \brief Return number of biquadratic quadrangles in the mesh
1532  */
1533 //================================================================================
1534
1535 int SMESH_Mesh::NbBiQuadQuadrangles() const throw(SALOME_Exception)
1536 {
1537   Unexpect aCatch(SalomeException);
1538   return _myMeshDS->GetMeshInfo().NbBiQuadQuadrangles();
1539 }
1540
1541 //================================================================================
1542 /*!
1543  * \brief Return the number of polygonal faces in the mesh
1544  */
1545 //================================================================================
1546
1547 int SMESH_Mesh::NbPolygons() const throw(SALOME_Exception)
1548 {
1549   Unexpect aCatch(SalomeException);
1550   return _myMeshDS->GetMeshInfo().NbPolygons();
1551 }
1552
1553 //================================================================================
1554 /*!
1555  * \brief Return number of volumes of given order in the mesh
1556  */
1557 //================================================================================
1558
1559 int SMESH_Mesh::NbVolumes(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1560 {
1561   Unexpect aCatch(SalomeException);
1562   return _myMeshDS->GetMeshInfo().NbVolumes(order);
1563 }
1564
1565 //================================================================================
1566 /*!
1567  * \brief  Return number of tetrahedrons of given order in the mesh
1568  */
1569 //================================================================================
1570
1571 int SMESH_Mesh::NbTetras(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1572 {
1573   Unexpect aCatch(SalomeException);
1574   return _myMeshDS->GetMeshInfo().NbTetras(order);
1575 }
1576
1577 //================================================================================
1578 /*!
1579  * \brief  Return number of hexahedrons of given order in the mesh
1580  */
1581 //================================================================================
1582
1583 int SMESH_Mesh::NbHexas(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1584 {
1585   Unexpect aCatch(SalomeException);
1586   return _myMeshDS->GetMeshInfo().NbHexas(order);
1587 }
1588
1589 //================================================================================
1590 /*!
1591  * \brief  Return number of triquadratic hexahedrons in the mesh
1592  */
1593 //================================================================================
1594
1595 int SMESH_Mesh::NbTriQuadraticHexas() const throw(SALOME_Exception)
1596 {
1597   Unexpect aCatch(SalomeException);
1598   return _myMeshDS->GetMeshInfo().NbTriQuadHexas();
1599 }
1600
1601 //================================================================================
1602 /*!
1603  * \brief  Return number of pyramids of given order in the mesh
1604  */
1605 //================================================================================
1606
1607 int SMESH_Mesh::NbPyramids(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1608 {
1609   Unexpect aCatch(SalomeException);
1610   return _myMeshDS->GetMeshInfo().NbPyramids(order);
1611 }
1612
1613 //================================================================================
1614 /*!
1615  * \brief  Return number of prisms (penthahedrons) of given order in the mesh
1616  */
1617 //================================================================================
1618
1619 int SMESH_Mesh::NbPrisms(SMDSAbs_ElementOrder order) const throw(SALOME_Exception)
1620 {
1621   Unexpect aCatch(SalomeException);
1622   return _myMeshDS->GetMeshInfo().NbPrisms(order);
1623 }
1624
1625 //================================================================================
1626 /*!
1627  * \brief  Return number of hexagonal prisms in the mesh
1628  */
1629 //================================================================================
1630
1631 int SMESH_Mesh::NbHexagonalPrisms() const throw(SALOME_Exception)
1632 {
1633   Unexpect aCatch(SalomeException);
1634   return _myMeshDS->GetMeshInfo().NbHexPrisms();
1635 }
1636
1637 //================================================================================
1638 /*!
1639  * \brief  Return number of polyhedrons in the mesh
1640  */
1641 //================================================================================
1642
1643 int SMESH_Mesh::NbPolyhedrons() const throw(SALOME_Exception)
1644 {
1645   Unexpect aCatch(SalomeException);
1646   return _myMeshDS->GetMeshInfo().NbPolyhedrons();
1647 }
1648
1649 //================================================================================
1650 /*!
1651  * \brief  Return number of ball elements in the mesh
1652  */
1653 //================================================================================
1654
1655 int SMESH_Mesh::NbBalls() const throw(SALOME_Exception)
1656 {
1657   Unexpect aCatch(SalomeException);
1658   return _myMeshDS->GetMeshInfo().NbBalls();
1659 }
1660
1661 //================================================================================
1662 /*!
1663  * \brief  Return number of submeshes in the mesh
1664  */
1665 //================================================================================
1666
1667 int SMESH_Mesh::NbSubMesh() const throw(SALOME_Exception)
1668 {
1669   Unexpect aCatch(SalomeException);
1670   return _myMeshDS->NbSubMesh();
1671 }
1672
1673 //=======================================================================
1674 //function : IsNotConformAllowed
1675 //purpose  : check if a hypothesis alowing notconform mesh is present
1676 //=======================================================================
1677
1678 bool SMESH_Mesh::IsNotConformAllowed() const
1679 {
1680   if(MYDEBUG) MESSAGE("SMESH_Mesh::IsNotConformAllowed");
1681
1682   static SMESH_HypoFilter filter( SMESH_HypoFilter::HasName( "NotConformAllowed" ));
1683   return GetHypothesis( _myMeshDS->ShapeToMesh(), filter, false );
1684 }
1685
1686 //=======================================================================
1687 //function : IsMainShape
1688 //purpose  : 
1689 //=======================================================================
1690
1691 bool SMESH_Mesh::IsMainShape(const TopoDS_Shape& theShape) const
1692 {
1693   return theShape.IsSame(_myMeshDS->ShapeToMesh() );
1694 }
1695
1696 //=============================================================================
1697 /*!
1698  *  
1699  */
1700 //=============================================================================
1701
1702 SMESH_Group* SMESH_Mesh::AddGroup (const SMDSAbs_ElementType theType,
1703                                    const char*               theName,
1704                                    int&                      theId,
1705                                    const TopoDS_Shape&       theShape,
1706                                    const SMESH_PredicatePtr& thePredicate)
1707 {
1708   if (_mapGroup.count(_groupId))
1709     return NULL;
1710   theId = _groupId;
1711   SMESH_Group* aGroup = new SMESH_Group (theId, this, theType, theName, theShape, thePredicate);
1712   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1713   _mapGroup[_groupId++] = aGroup;
1714   return aGroup;
1715 }
1716
1717 //================================================================================
1718 /*!
1719  * \brief Creates a group based on an existing SMESHDS group. Group ID should be unique
1720  */
1721 //================================================================================
1722
1723 SMESH_Group* SMESH_Mesh::AddGroup (SMESHDS_GroupBase* groupDS) throw(SALOME_Exception)
1724 {
1725   if ( !groupDS ) 
1726     throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup(): NULL SMESHDS_GroupBase"));
1727
1728   map <int, SMESH_Group*>::iterator i_g = _mapGroup.find( groupDS->GetID() );
1729   if ( i_g != _mapGroup.end() && i_g->second )
1730   {
1731     if ( i_g->second->GetGroupDS() == groupDS )
1732       return i_g->second;
1733     else
1734       throw SALOME_Exception(LOCALIZED ("SMESH_Mesh::AddGroup() wrong ID of SMESHDS_GroupBase"));
1735   }
1736   SMESH_Group* aGroup = new SMESH_Group (groupDS);
1737   _mapGroup[ groupDS->GetID() ] = aGroup;
1738   GetMeshDS()->AddGroup( aGroup->GetGroupDS() );
1739
1740   _groupId = 1 + _mapGroup.rbegin()->first;
1741
1742   return aGroup;
1743 }
1744
1745
1746 //================================================================================
1747 /*!
1748  * \brief Creates SMESH_Groups for not wrapped SMESHDS_Groups
1749  *  \retval bool - true if new SMESH_Groups have been created
1750  * 
1751  */
1752 //================================================================================
1753
1754 bool SMESH_Mesh::SynchronizeGroups()
1755 {
1756   int nbGroups = _mapGroup.size();
1757   const set<SMESHDS_GroupBase*>& groups = _myMeshDS->GetGroups();
1758   set<SMESHDS_GroupBase*>::const_iterator gIt = groups.begin();
1759   for ( ; gIt != groups.end(); ++gIt )
1760   {
1761     SMESHDS_GroupBase* groupDS = (SMESHDS_GroupBase*) *gIt;
1762     _groupId = groupDS->GetID();
1763     if ( !_mapGroup.count( _groupId ))
1764       _mapGroup[_groupId] = new SMESH_Group( groupDS );
1765   }
1766   if ( !_mapGroup.empty() )
1767     _groupId = _mapGroup.rbegin()->first + 1;
1768
1769   return nbGroups < _mapGroup.size();
1770 }
1771
1772 //================================================================================
1773 /*!
1774  * \brief Return iterator on all existing groups
1775  */
1776 //================================================================================
1777
1778 SMESH_Mesh::GroupIteratorPtr SMESH_Mesh::GetGroups() const
1779 {
1780   typedef map <int, SMESH_Group *> TMap;
1781   return GroupIteratorPtr( new SMDS_mapIterator<TMap>( _mapGroup ));
1782 }
1783
1784 //=============================================================================
1785 /*!
1786  * \brief Return a group by ID
1787  */
1788 //=============================================================================
1789
1790 SMESH_Group* SMESH_Mesh::GetGroup (const int theGroupID)
1791 {
1792   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1793     return NULL;
1794   return _mapGroup[theGroupID];
1795 }
1796
1797
1798 //=============================================================================
1799 /*!
1800  * \brief Return IDs of all groups
1801  */
1802 //=============================================================================
1803
1804 list<int> SMESH_Mesh::GetGroupIds() const
1805 {
1806   list<int> anIds;
1807   for ( map<int, SMESH_Group*>::const_iterator it = _mapGroup.begin(); it != _mapGroup.end(); it++ )
1808     anIds.push_back( it->first );
1809   
1810   return anIds;
1811 }
1812
1813 //================================================================================
1814 /*!
1815  * \brief Set a caller of methods at level of CORBA API implementation.
1816  * The set upCaller will be deleted by SMESH_Mesh
1817  */
1818 //================================================================================
1819
1820 void SMESH_Mesh::SetCallUp( TCallUp* upCaller )
1821 {
1822   if ( _callUp ) delete _callUp;
1823   _callUp = upCaller;
1824 }
1825
1826 //=============================================================================
1827 /*!
1828  *  
1829  */
1830 //=============================================================================
1831
1832 bool SMESH_Mesh::RemoveGroup (const int theGroupID)
1833 {
1834   if (_mapGroup.find(theGroupID) == _mapGroup.end())
1835     return false;
1836   GetMeshDS()->RemoveGroup( _mapGroup[theGroupID]->GetGroupDS() );
1837   delete _mapGroup[theGroupID];
1838   _mapGroup.erase (theGroupID);
1839   if (_callUp)
1840     _callUp->RemoveGroup( theGroupID );
1841   return true;
1842 }
1843
1844 //=======================================================================
1845 //function : GetAncestors
1846 //purpose  : return list of ancestors of theSubShape in the order
1847 //           that lower dimention shapes come first.
1848 //=======================================================================
1849
1850 const TopTools_ListOfShape& SMESH_Mesh::GetAncestors(const TopoDS_Shape& theS) const
1851 {
1852   if ( _mapAncestors.Contains( theS ) )
1853     return _mapAncestors.FindFromKey( theS );
1854
1855   static TopTools_ListOfShape emptyList;
1856   return emptyList;
1857 }
1858
1859 //=======================================================================
1860 //function : Dump
1861 //purpose  : dumps contents of mesh to stream [ debug purposes ]
1862 //=======================================================================
1863
1864 ostream& SMESH_Mesh::Dump(ostream& save)
1865 {
1866   int clause = 0;
1867   save << "========================== Dump contents of mesh ==========================" << endl << endl;
1868   save << ++clause << ") Total number of nodes:   \t"    << NbNodes() << endl;
1869   save << ++clause << ") Total number of edges:   \t"    << NbEdges() << endl;
1870   save << ++clause << ") Total number of faces:   \t"    << NbFaces() << endl;
1871   save << ++clause << ") Total number of polygons:\t"    << NbPolygons() << endl;
1872   save << ++clause << ") Total number of volumes:\t"     << NbVolumes() << endl;
1873   save << ++clause << ") Total number of polyhedrons:\t" << NbPolyhedrons() << endl << endl;
1874   for ( int isQuadratic = 0; isQuadratic < 2; ++isQuadratic )
1875   {
1876     string orderStr = isQuadratic ? "quadratic" : "linear";
1877     SMDSAbs_ElementOrder order  = isQuadratic ? ORDER_QUADRATIC : ORDER_LINEAR;
1878
1879     save << ++clause << ") Total number of " << orderStr << " edges:\t" << NbEdges(order) << endl;
1880     save << ++clause << ") Total number of " << orderStr << " faces:\t" << NbFaces(order) << endl;
1881     if ( NbFaces(order) > 0 ) {
1882       int nb3 = NbTriangles(order);
1883       int nb4 = NbQuadrangles(order);
1884       save << clause << ".1) Number of " << orderStr << " triangles:  \t" << nb3 << endl;
1885       save << clause << ".2) Number of " << orderStr << " quadrangles:\t" << nb4 << endl;
1886       if ( nb3 + nb4 !=  NbFaces(order) ) {
1887         map<int,int> myFaceMap;
1888         SMDS_FaceIteratorPtr itFaces=_myMeshDS->facesIterator();
1889         while( itFaces->more( ) ) {
1890           int nbNodes = itFaces->next()->NbNodes();
1891           if ( myFaceMap.find( nbNodes ) == myFaceMap.end() )
1892             myFaceMap[ nbNodes ] = 0;
1893           myFaceMap[ nbNodes ] = myFaceMap[ nbNodes ] + 1;
1894         }
1895         save << clause << ".3) Faces in detail: " << endl;
1896         map <int,int>::iterator itF;
1897         for (itF = myFaceMap.begin(); itF != myFaceMap.end(); itF++)
1898           save << "--> nb nodes: " << itF->first << " - nb elemens:\t" << itF->second << endl;
1899       }
1900     }
1901     save << ++clause << ") Total number of " << orderStr << " volumes:\t" << NbVolumes(order) << endl;
1902     if ( NbVolumes(order) > 0 ) {
1903       int nb8 = NbHexas(order);
1904       int nb4 = NbTetras(order);
1905       int nb5 = NbPyramids(order);
1906       int nb6 = NbPrisms(order);
1907       save << clause << ".1) Number of " << orderStr << " hexahedrons:\t" << nb8 << endl;
1908       save << clause << ".2) Number of " << orderStr << " tetrahedrons:\t" << nb4 << endl;
1909       save << clause << ".3) Number of " << orderStr << " prisms:      \t" << nb6 << endl;
1910       save << clause << ".4) Number of " << orderStr << " pyramids:\t" << nb5 << endl;
1911       if ( nb8 + nb4 + nb5 + nb6 != NbVolumes(order) ) {
1912         map<int,int> myVolumesMap;
1913         SMDS_VolumeIteratorPtr itVolumes=_myMeshDS->volumesIterator();
1914         while( itVolumes->more( ) ) {
1915           int nbNodes = itVolumes->next()->NbNodes();
1916           if ( myVolumesMap.find( nbNodes ) == myVolumesMap.end() )
1917             myVolumesMap[ nbNodes ] = 0;
1918           myVolumesMap[ nbNodes ] = myVolumesMap[ nbNodes ] + 1;
1919         }
1920         save << clause << ".5) Volumes in detail: " << endl;
1921         map <int,int>::iterator itV;
1922         for (itV = myVolumesMap.begin(); itV != myVolumesMap.end(); itV++)
1923           save << "--> nb nodes: " << itV->first << " - nb elemens:\t" << itV->second << endl;
1924       }
1925     }
1926     save << endl;
1927   }
1928   save << "===========================================================================" << endl;
1929   return save;
1930 }
1931
1932 //=======================================================================
1933 //function : GetElementType
1934 //purpose  : Returns type of mesh element with certain id
1935 //=======================================================================
1936
1937 SMDSAbs_ElementType SMESH_Mesh::GetElementType( const int id, const bool iselem )
1938 {
1939   return _myMeshDS->GetElementType( id, iselem );
1940 }
1941
1942 //=============================================================================
1943 /*!
1944  *  \brief Convert group on geometry into standalone group
1945  */
1946 //=============================================================================
1947
1948 SMESH_Group* SMESH_Mesh::ConvertToStandalone ( int theGroupID )
1949 {
1950   SMESH_Group* aGroup = 0;
1951   map < int, SMESH_Group * >::iterator itg = _mapGroup.find( theGroupID );
1952   if ( itg == _mapGroup.end() )
1953     return aGroup;
1954
1955   SMESH_Group* anOldGrp = (*itg).second;
1956   SMESHDS_GroupBase* anOldGrpDS = anOldGrp->GetGroupDS();
1957   if ( !anOldGrp || !anOldGrpDS )
1958     return aGroup;
1959
1960   // create new standalone group
1961   aGroup = new SMESH_Group (theGroupID, this, anOldGrpDS->GetType(), anOldGrp->GetName() );
1962   _mapGroup[theGroupID] = aGroup;
1963
1964   SMESHDS_Group* aNewGrpDS = dynamic_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1965   GetMeshDS()->RemoveGroup( anOldGrpDS );
1966   GetMeshDS()->AddGroup( aNewGrpDS );
1967
1968   // add elements (or nodes) into new created group
1969   SMDS_ElemIteratorPtr anItr = anOldGrpDS->GetElements();
1970   while ( anItr->more() )
1971     aNewGrpDS->Add( (anItr->next())->GetID() );
1972
1973   // set color
1974   aNewGrpDS->SetColor( anOldGrpDS->GetColor() );
1975
1976   // remove old group
1977   delete anOldGrp;
1978
1979   return aGroup;
1980 }
1981
1982 //=============================================================================
1983 /*!
1984  *  \brief remove submesh order  from Mesh
1985  */
1986 //=============================================================================
1987
1988 void SMESH_Mesh::ClearMeshOrder()
1989 {
1990   _mySubMeshOrder.clear();
1991 }
1992
1993 //=============================================================================
1994 /*!
1995  *  \brief remove submesh order  from Mesh
1996  */
1997 //=============================================================================
1998
1999 void SMESH_Mesh::SetMeshOrder(const TListOfListOfInt& theOrder )
2000 {
2001   _mySubMeshOrder = theOrder;
2002 }
2003
2004 //=============================================================================
2005 /*!
2006  *  \brief return submesh order if any
2007  */
2008 //=============================================================================
2009
2010 const TListOfListOfInt& SMESH_Mesh::GetMeshOrder() const
2011 {
2012   return _mySubMeshOrder;
2013 }
2014
2015 //=============================================================================
2016 /*!
2017  *  \brief fill _mapAncestors
2018  */
2019 //=============================================================================
2020
2021 void SMESH_Mesh::fillAncestorsMap(const TopoDS_Shape& theShape)
2022 {
2023
2024   int desType, ancType;
2025   if ( !theShape.IsSame( GetShapeToMesh()) && theShape.ShapeType() == TopAbs_COMPOUND )
2026   {
2027     // a geom group is added. Insert it into lists of ancestors before
2028     // the first ancestor more complex than group members
2029     TopoDS_Iterator subIt( theShape );
2030     if ( !subIt.More() ) return;
2031     int memberType = subIt.Value().ShapeType();
2032     for ( desType = TopAbs_VERTEX; desType >= memberType; desType-- )
2033       for (TopExp_Explorer des( theShape, TopAbs_ShapeEnum( desType )); des.More(); des.Next())
2034       {
2035         if ( !_mapAncestors.Contains( des.Current() )) continue;// issue 0020982
2036         TopTools_ListOfShape& ancList = _mapAncestors.ChangeFromKey( des.Current() );
2037         TopTools_ListIteratorOfListOfShape ancIt (ancList);
2038         while ( ancIt.More() && ancIt.Value().ShapeType() >= memberType )
2039           ancIt.Next();
2040         if ( ancIt.More() )
2041           ancList.InsertBefore( theShape, ancIt );
2042       }
2043   }
2044   {
2045     for ( desType = TopAbs_VERTEX; desType > TopAbs_COMPOUND; desType-- )
2046       for ( ancType = desType - 1; ancType >= TopAbs_COMPOUND; ancType-- )
2047         TopExp::MapShapesAndAncestors ( theShape,
2048                                         (TopAbs_ShapeEnum) desType,
2049                                         (TopAbs_ShapeEnum) ancType,
2050                                         _mapAncestors );
2051   }
2052   // visit COMPOUNDs inside a COMPOUND that are not reachable by TopExp_Explorer
2053   if ( theShape.ShapeType() == TopAbs_COMPOUND )
2054   {
2055     for ( TopoDS_Iterator sIt(theShape); sIt.More(); sIt.Next() )
2056       if ( sIt.Value().ShapeType() == TopAbs_COMPOUND )
2057         fillAncestorsMap( sIt.Value() );
2058   }
2059 }
2060
2061 //=============================================================================
2062 /*!
2063  * \brief sort submeshes according to stored mesh order
2064  * \param theListToSort in out list to be sorted
2065  * \return FALSE if nothing sorted
2066  */
2067 //=============================================================================
2068
2069 bool SMESH_Mesh::SortByMeshOrder(list<SMESH_subMesh*>& theListToSort) const
2070 {
2071   if ( !_mySubMeshOrder.size() || theListToSort.size() < 2)
2072     return true;
2073   
2074   bool res = false;
2075   list<SMESH_subMesh*> onlyOrderedList;
2076   // collect all ordered submeshes in one list as pointers
2077   // and get their positions within theListToSort
2078   typedef list<SMESH_subMesh*>::iterator TPosInList;
2079   map< int, TPosInList > sortedPos;
2080   TPosInList smBeg = theListToSort.begin(), smEnd = theListToSort.end();
2081   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2082   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++) {
2083     const TListOfInt& listOfId = *listIdsIt;
2084     TListOfInt::const_iterator idIt = listOfId.begin();
2085     for ( ; idIt != listOfId.end(); idIt++ ) {
2086       if ( SMESH_subMesh * sm = GetSubMeshContaining( *idIt )) {
2087         TPosInList smPos = find( smBeg, smEnd, sm );
2088         if ( smPos != smEnd ) {
2089           onlyOrderedList.push_back( sm );
2090           sortedPos[ distance( smBeg, smPos )] = smPos;
2091         }
2092       }
2093     }
2094   }
2095   if (onlyOrderedList.size() < 2)
2096     return res;
2097   res = true;
2098
2099   list<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
2100   list<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
2101
2102   // iterate on ordered submeshes and insert them in detected positions
2103   map< int, TPosInList >::iterator i_pos = sortedPos.begin();
2104   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
2105     *(i_pos->second) = *onlyBIt;
2106
2107   return res;
2108 }
2109
2110 //================================================================================
2111 /*!
2112  * \brief Return true if given order of sub-meshes is OK
2113  */
2114 //================================================================================
2115
2116 bool SMESH_Mesh::IsOrderOK( const SMESH_subMesh* smBefore,
2117                             const SMESH_subMesh* smAfter ) const
2118 {
2119   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
2120   TListOfInt::const_iterator idBef, idAft;
2121   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
2122   {
2123     const TListOfInt& listOfId = *listIdsIt;
2124     idBef = std::find( listOfId.begin(), listOfId.end(), smBefore->GetId() );
2125     if ( idBef != listOfId.end() )
2126       idAft = std::find( listOfId.begin(), listOfId.end(), smAfter->GetId() );
2127     if ( idAft != listOfId.end () )
2128       return ( std::distance( listOfId.begin(), idBef ) <
2129                std::distance( listOfId.begin(), idAft )   );
2130   }
2131   return true; // no order imposed to given submeshes
2132
2133
2134 //=============================================================================
2135 /*!
2136  * \brief sort submeshes according to stored mesh order
2137  * \param theListToSort in out list to be sorted
2138  * \return FALSE if nothing sorted
2139  */
2140 //=============================================================================
2141
2142 list<SMESH_subMesh*>
2143 SMESH_Mesh::getAncestorsSubMeshes (const TopoDS_Shape& theSubShape) const
2144 {
2145   list<SMESH_subMesh*> listOfSubMesh;
2146   TopTools_ListIteratorOfListOfShape it( GetAncestors( theSubShape ));
2147   for (; it.More(); it.Next() )
2148     if ( SMESH_subMesh* sm = GetSubMeshContaining( it.Value() ))
2149       listOfSubMesh.push_back(sm);
2150
2151   // sort submeshes according to stored mesh order
2152   SortByMeshOrder( listOfSubMesh );
2153
2154   return listOfSubMesh;
2155 }