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