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