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