]> SALOME platform Git repositories - modules/smesh.git/blob - src/OBJECT/SMESH_Object.cxx
Salome HOME
PAL16842 (Genertion of groups when a mesh is transformed)
[modules/smesh.git] / src / OBJECT / SMESH_Object.cxx
1 //  SMESH OBJECT : interactive object for SMESH visualization
2 //
3 //  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
5 // 
6 //  This library is free software; you can redistribute it and/or 
7 //  modify it under the terms of the GNU Lesser General Public 
8 //  License as published by the Free Software Foundation; either 
9 //  version 2.1 of the License. 
10 // 
11 //  This library is distributed in the hope that it will be useful, 
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of 
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
14 //  Lesser General Public License for more details. 
15 // 
16 //  You should have received a copy of the GNU Lesser General Public 
17 //  License along with this library; if not, write to the Free Software 
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
19 // 
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //
23 //
24 //  File   : SMESH_Grid.cxx
25 //  Author : Nicolas REJNERI
26 //  Module : SMESH
27
28 #include "SMESH_ObjectDef.h"
29 #include "SMESH_ActorUtils.h"
30
31 #include "SMDS_Mesh.hxx"
32 #include "SMDS_PolyhedralVolumeOfNodes.hxx"
33 #include "SMESH_Actor.h"
34 #include "SMESH_ControlsDef.hxx"
35 #include "SalomeApp_Application.h"
36 #include "VTKViewer_ExtractUnstructuredGrid.h"
37 #include "VTKViewer_CellLocationsArray.h"
38
39 #include CORBA_SERVER_HEADER(SMESH_Gen)
40 #include CORBA_SERVER_HEADER(SALOME_Exception)
41
42 #include <vtkCell.h>
43 #include <vtkIdList.h>
44 #include <vtkCellArray.h>
45 #include <vtkUnsignedCharArray.h>
46
47 #include <vtkUnstructuredGrid.h>
48
49 #include <memory>
50 #include <sstream>      
51 #include <stdexcept>
52 #include <set>
53
54 #include "utilities.h"
55
56 using namespace std;
57
58 #ifndef EXCEPTION
59 #define EXCEPTION(TYPE, MSG) {\
60   std::ostringstream aStream;\
61   aStream<<__FILE__<<"["<<__LINE__<<"]::"<<MSG;\
62   throw TYPE(aStream.str());\
63 }
64 #endif
65
66 #ifdef _DEBUG_
67 static int MYDEBUG = 0;
68 static int MYDEBUGWITHFILES = 0;
69 #else
70 static int MYDEBUG = 0;
71 static int MYDEBUGWITHFILES = 0;
72 #endif
73
74
75 /*
76   Class       : SMESH_VisualObjDef
77   Description : Base class for all mesh objects to be visuilised
78 */
79
80 //=================================================================================
81 // function : getCellType
82 // purpose  : Get type of VTK cell
83 //=================================================================================
84 static inline vtkIdType getCellType( const SMDSAbs_ElementType theType,
85                                      const bool thePoly,
86                                      const int theNbNodes )
87 {
88   switch( theType )
89   {
90     case SMDSAbs_Edge: 
91       if( theNbNodes == 2 )         return VTK_LINE;
92       else if ( theNbNodes == 3 )   return VTK_QUADRATIC_EDGE;
93       else return VTK_EMPTY_CELL;
94
95     case SMDSAbs_Face  :
96       if (thePoly && theNbNodes>2 ) return VTK_POLYGON;
97       else if ( theNbNodes == 3 )   return VTK_TRIANGLE;
98       else if ( theNbNodes == 4 )   return VTK_QUAD;
99       else if ( theNbNodes == 6 )   return VTK_QUADRATIC_TRIANGLE;
100       else if ( theNbNodes == 8 )   return VTK_QUADRATIC_QUAD;
101       else return VTK_EMPTY_CELL;
102       
103     case SMDSAbs_Volume:
104       if (thePoly && theNbNodes>3 ) return VTK_CONVEX_POINT_SET;
105       else if ( theNbNodes == 4 )   return VTK_TETRA;
106       else if ( theNbNodes == 5 )   return VTK_PYRAMID;
107       else if ( theNbNodes == 6 )   return VTK_WEDGE;
108       else if ( theNbNodes == 8 )   return VTK_HEXAHEDRON;
109       else if ( theNbNodes == 10 )  {
110         return VTK_QUADRATIC_TETRA;
111       }
112       else if ( theNbNodes == 20 )  {
113         return VTK_QUADRATIC_HEXAHEDRON;
114       }
115       else if ( theNbNodes==13 || theNbNodes==15 )  {
116         return VTK_CONVEX_POINT_SET;
117       }
118       else return VTK_EMPTY_CELL;
119
120     default: return VTK_EMPTY_CELL;
121   }
122 }
123
124 //=================================================================================
125 // functions : SMESH_VisualObjDef
126 // purpose   : Constructor
127 //=================================================================================
128 SMESH_VisualObjDef::SMESH_VisualObjDef()
129 {
130   myGrid = vtkUnstructuredGrid::New();
131 }
132 SMESH_VisualObjDef::~SMESH_VisualObjDef()
133 {
134   if ( MYDEBUG )
135     MESSAGE( "~SMESH_MeshObj - myGrid->GetReferenceCount() = " << myGrid->GetReferenceCount() );
136   myGrid->Delete();
137 }
138
139 //=================================================================================
140 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
141 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
142 //=================================================================================
143 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
144 {
145   return myVTK2SMDSNodes.find(theVTKID) == myVTK2SMDSNodes.end() ? -1 : myVTK2SMDSNodes[theVTKID];
146 }
147
148 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
149 {
150   return mySMDS2VTKNodes.find(theObjID) == mySMDS2VTKNodes.end() ? -1 : mySMDS2VTKNodes[theObjID];
151 }
152
153 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
154 {
155   return myVTK2SMDSElems.find(theVTKID) == myVTK2SMDSElems.end() ? -1 : myVTK2SMDSElems[theVTKID];
156 }
157
158 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
159 {
160   return mySMDS2VTKElems.find(theObjID) == mySMDS2VTKElems.end() ? -1 : mySMDS2VTKElems[theObjID];
161 }
162
163 //=================================================================================
164 // function : SMESH_VisualObjDef::createPoints
165 // purpose  : Create points from nodes
166 //=================================================================================
167 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
168 {
169   if ( thePoints == 0 )
170     return;
171
172   TEntityList aNodes;
173   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
174   thePoints->SetNumberOfPoints( nbNodes );
175   
176   int nbPoints = 0;
177
178   TEntityList::const_iterator anIter;
179   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
180   {
181     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
182     if ( aNode != 0 )
183     {
184       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
185       int anId = aNode->GetID();
186       mySMDS2VTKNodes.insert( TMapOfIds::value_type( anId, nbPoints ) );
187       myVTK2SMDSNodes.insert( TMapOfIds::value_type( nbPoints, anId ) );
188       nbPoints++;
189     }
190   }
191
192   if ( nbPoints != nbNodes )
193     thePoints->SetNumberOfPoints( nbPoints );
194 }
195
196 //=================================================================================
197 // function : buildPrs
198 // purpose  : create VTK cells( fill unstructured grid )
199 //=================================================================================
200 void SMESH_VisualObjDef::buildPrs()
201 {
202   try 
203   {
204     mySMDS2VTKNodes.clear();
205     myVTK2SMDSNodes.clear();
206     mySMDS2VTKElems.clear();
207     myVTK2SMDSElems.clear();
208     
209     if ( IsNodePrs() )
210       buildNodePrs();
211     else
212       buildElemPrs();
213   }
214   catch(...)
215   {
216     // TODO: catch inside buildNodePrs and buildElemPrs and free
217     // temporarily allocated objects
218     mySMDS2VTKNodes.clear();
219     myVTK2SMDSNodes.clear();
220     mySMDS2VTKElems.clear();
221     myVTK2SMDSElems.clear();
222
223     myGrid->SetPoints( 0 );
224     myGrid->SetCells( 0, 0, 0 );
225     throw;
226   }
227   
228   if( MYDEBUG ) MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
229   if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
230 }
231
232 //=================================================================================
233 // function : buildNodePrs
234 // purpose  : create VTK cells for nodes
235 //=================================================================================
236 void SMESH_VisualObjDef::buildNodePrs()
237 {
238   // PAL16631: without swap, bad_alloc is not thrown but hung up and crash instead,
239   // so check remaining memory size for safety
240   SMDS_Mesh::CheckMemory(); // PAL16631
241   vtkPoints* aPoints = vtkPoints::New();
242   createPoints( aPoints );
243   SMDS_Mesh::CheckMemory();
244   myGrid->SetPoints( aPoints );
245   aPoints->Delete();
246
247   myGrid->SetCells( 0, 0, 0 );
248 }
249
250 //=================================================================================
251 // function : buildElemPrs
252 // purpose  : Create VTK cells for elements
253 //=================================================================================
254
255 namespace{
256   typedef std::vector<const SMDS_MeshElement*> TConnect;
257
258   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
259                  TConnect& theConnect)
260   {
261     theConnect.clear();
262     for(; theNodesIter->more();)
263       theConnect.push_back(theNodesIter->next());
264     return theConnect.size();
265   }
266   
267   inline 
268   void SetId(vtkIdList *theIdList, 
269              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
270              const TConnect& theConnect, 
271              int thePosition,
272              int theId)
273   {
274     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
275   }
276
277 }
278
279
280 void SMESH_VisualObjDef::buildElemPrs()
281 {
282   // Create points
283   
284   vtkPoints* aPoints = vtkPoints::New();
285   createPoints( aPoints );
286   myGrid->SetPoints( aPoints );
287   aPoints->Delete();
288     
289   if ( MYDEBUG )
290     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
291
292   // Calculate cells size
293
294   static SMDSAbs_ElementType aTypes[ 3 ] = { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
295
296   // get entity data
297   map<SMDSAbs_ElementType,int> nbEnts;
298   map<SMDSAbs_ElementType,TEntityList> anEnts;
299
300   for ( int i = 0; i <= 2; i++ )
301     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
302
303   // PAL16631: without swap, bad_alloc is not thrown but hung up and crash instead,
304   // so check remaining memory size for safety
305   SMDS_Mesh::CheckMemory(); // PAL16631
306
307   vtkIdType aCellsSize =  3 * nbEnts[ SMDSAbs_Edge ];
308
309   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
310   {
311     if ( nbEnts[ aTypes[ i ] ] )
312     {
313       const TEntityList& aList = anEnts[ aTypes[ i ] ];
314       TEntityList::const_iterator anIter;
315       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
316         aCellsSize += (*anIter)->NbNodes() + 1;
317     }
318   }
319
320   vtkIdType aNbCells = nbEnts[ SMDSAbs_Edge ] + nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
321   
322   if ( MYDEBUG )
323     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
324
325   // Create cells
326   
327   vtkCellArray* aConnectivity = vtkCellArray::New();
328   aConnectivity->Allocate( aCellsSize, 0 );
329   
330   SMDS_Mesh::CheckMemory(); // PAL16631
331
332   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
333   aCellTypesArray->SetNumberOfComponents( 1 );
334   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
335   
336   SMDS_Mesh::CheckMemory(); // PAL16631
337
338   vtkIdList *anIdList = vtkIdList::New();
339   vtkIdType iElem = 0;
340
341   TConnect aConnect;
342   aConnect.reserve(VTK_CELL_SIZE);
343
344   SMDS_Mesh::CheckMemory(); // PAL16631
345
346   for ( int i = 0; i <= 2; i++ ) // iterate through edges, faces and volumes
347   {
348     if( nbEnts[ aTypes[ i ] ] > 0 )
349     {
350       const SMDSAbs_ElementType& aType = aTypes[ i ];
351       const TEntityList& aList = anEnts[ aType ];
352       TEntityList::const_iterator anIter;
353       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
354       {
355         const SMDS_MeshElement* anElem = *anIter;
356         
357         vtkIdType aNbNodes = anElem->NbNodes();
358         anIdList->SetNumberOfIds( aNbNodes );
359
360         int anId = anElem->GetID();
361
362         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
363         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
364
365         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
366         switch(aType){
367         case SMDSAbs_Volume:{
368           aConnect.clear();
369           std::vector<int> aConnectivities;
370           // Convertions connectivities from SMDS to VTK
371           if (anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
372
373             if ( const SMDS_PolyhedralVolumeOfNodes* ph =
374                  dynamic_cast<const SMDS_PolyhedralVolumeOfNodes*> (anElem))
375             {
376               aNbNodes = GetConnect(ph->uniqueNodesIterator(),aConnect);
377               anIdList->SetNumberOfIds( aNbNodes );
378             }
379             for (int k = 0; k < aNbNodes; k++)
380               aConnectivities.push_back(k);
381
382           } else if (aNbNodes == 4) {
383             static int anIds[] = {0,2,1,3};
384             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
385
386           } else if (aNbNodes == 5) {
387             static int anIds[] = {0,3,2,1,4};
388             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
389
390           } else if (aNbNodes == 6) {
391             static int anIds[] = {0,1,2,3,4,5};
392             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
393
394           }
395           else if (aNbNodes == 8) {
396             static int anIds[] = {0,3,2,1,4,7,6,5};
397             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
398
399           }
400           else if (aNbNodes == 10) {
401             static int anIds[] = {0,2,1,3,6,5,4,7,9,8};
402             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
403           }
404           else if (aNbNodes == 13) {
405             static int anIds[] = {0,3,2,1,4,8,7,6,5,9,12,11,10};
406             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
407           }
408           else if (aNbNodes == 15) {
409             static int anIds[] = {0,2,1,3,5,4,8,7,6,11,10,9,12,14,13};
410             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
411             //for (int k = 0; k < aNbNodes; k++) {
412             //  int nn = aConnectivities[k];
413             //  const SMDS_MeshNode* N = static_cast<const SMDS_MeshNode*> (aConnect[nn]);
414             //  cout<<"k="<<k<<"  N("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
415             //}
416           }
417           else if (aNbNodes == 20) {
418             static int anIds[] = {0,3,2,1,4,7,6,5,11,10,9,8,15,14,13,12,16,19,18,17};
419             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
420           }
421           else {
422           }
423
424           if ( aConnect.empty() )
425             GetConnect(aNodesIter,aConnect);
426
427           if (aConnectivities.size() > 0) {
428             for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
429               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
430           }
431           break;
432         }
433         default:
434           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
435             const SMDS_MeshElement* aNode = aNodesIter->next();
436             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
437           }
438         }
439
440         aConnectivity->InsertNextCell( anIdList );
441         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
442
443         iElem++;
444       }
445     }
446     SMDS_Mesh::CheckMemory(); // PAL16631
447   }
448
449   // Insert cells in grid
450   
451   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
452   aCellLocationsArray->SetNumberOfComponents( 1 );
453   aCellLocationsArray->SetNumberOfTuples( aNbCells );
454   
455   SMDS_Mesh::CheckMemory(); // PAL16631
456
457   aConnectivity->InitTraversal();
458   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
459     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
460
461   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
462   
463   aCellLocationsArray->Delete();
464   aCellTypesArray->Delete();
465   aConnectivity->Delete();
466   anIdList->Delete();
467
468   SMDS_Mesh::CheckMemory(); // PAL16631
469 }
470
471 //=================================================================================
472 // function : GetEdgeNodes
473 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
474 //=================================================================================
475 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
476                                        const int theEdgeNum,
477                                        int&      theNodeId1,
478                                        int&      theNodeId2 ) const
479 {
480   const SMDS_Mesh* aMesh = GetMesh();
481   if ( aMesh == 0 )
482     return false;
483     
484   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
485   if ( anElem == 0 )
486     return false;
487     
488   int nbNodes = anElem->NbNodes();
489
490   if ( theEdgeNum < 0 || theEdgeNum > 3 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
491     return false;
492
493   int anIds[ nbNodes ];
494   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
495   int i = 0;
496   while( anIter->more() )
497     anIds[ i++ ] = anIter->next()->GetID();
498
499   if ( theEdgeNum < nbNodes - 1 )
500   {
501     theNodeId1 = anIds[ theEdgeNum ];
502     theNodeId2 = anIds[ theEdgeNum + 1 ];
503   }
504   else
505   {
506     theNodeId1 = anIds[ nbNodes - 1 ];
507     theNodeId2 = anIds[ 0 ];
508   }
509
510   return true;
511 }
512
513 /*
514   Class       : SMESH_MeshObj
515   Description : Class for visualisation of mesh
516 */
517
518 //=================================================================================
519 // function : SMESH_MeshObj
520 // purpose  : Constructor
521 //=================================================================================
522 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
523   myClient(SalomeApp_Application::orb(),theMesh)
524 {
525   if ( MYDEBUG ) 
526     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
527 }
528
529 //=================================================================================
530 // function : ~SMESH_MeshObj
531 // purpose  : Destructor
532 //=================================================================================
533 SMESH_MeshObj::~SMESH_MeshObj()
534 {
535   if ( MYDEBUG ) 
536     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
537 }
538
539 //=================================================================================
540 // function : Update
541 // purpose  : Update mesh and fill grid with new values if necessary 
542 //=================================================================================
543 bool SMESH_MeshObj::Update( int theIsClear )
544 {
545   // Update SMDS_Mesh on client part
546   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
547     buildPrs();  // Fill unstructured grid
548     return true;
549   }
550   return false;
551 }
552
553 //=================================================================================
554 // function : GetElemDimension
555 // purpose  : Get dimension of element
556 //=================================================================================
557 int SMESH_MeshObj::GetElemDimension( const int theObjId )
558 {
559   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
560   if ( anElem == 0 )
561     return 0;
562
563   int aType = anElem->GetType();
564   switch ( aType )
565   {
566     case SMDSAbs_Edge  : return 1;
567     case SMDSAbs_Face  : return 2;
568     case SMDSAbs_Volume: return 3;
569     default            : return 0;
570   }
571 }
572
573 //=================================================================================
574 // function : GetEntities
575 // purpose  : Get entities of specified type. Return number of entities
576 //=================================================================================
577 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
578 {
579   switch ( theType )
580   {
581     case SMDSAbs_Node:
582     {
583       return myClient->NbNodes();
584     }
585     break;
586     case SMDSAbs_Edge:
587     {
588       return myClient->NbEdges();
589     }
590     break;
591     case SMDSAbs_Face:
592     {
593       return myClient->NbFaces();
594     }
595     break;
596     case SMDSAbs_Volume:
597     {
598       return myClient->NbVolumes();
599     }
600     break;
601     default:
602       return 0;
603     break;
604   }
605 }
606
607 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
608 {
609   theObjs.clear();
610
611   switch ( theType )
612   {
613     case SMDSAbs_Node:
614     {
615       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
616       while ( anIter->more() ) theObjs.push_back( anIter->next() );
617     }
618     break;
619     case SMDSAbs_Edge:
620     {
621       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
622       while ( anIter->more() ) theObjs.push_back( anIter->next() );
623     }
624     break;
625     case SMDSAbs_Face:
626     {
627       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
628       while ( anIter->more() ) theObjs.push_back( anIter->next() );
629     }
630     break;
631     case SMDSAbs_Volume:
632     {
633       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
634       while ( anIter->more() ) theObjs.push_back( anIter->next() );
635     }
636     break;
637     default:
638     break;
639   }
640
641   return theObjs.size();
642 }
643
644 //=================================================================================
645 // function : UpdateFunctor
646 // purpose  : Update functor in accordance with current mesh
647 //=================================================================================
648 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
649 {
650   theFunctor->SetMesh( GetMesh() );
651 }
652
653 //=================================================================================
654 // function : IsNodePrs
655 // purpose  : Return true if node presentation is used
656 //=================================================================================
657 bool SMESH_MeshObj::IsNodePrs() const
658 {
659   return myClient->NbEdges() == 0 &&myClient->NbFaces() == 0 && myClient->NbVolumes() == 0 ;
660 }
661
662
663 /*
664   Class       : SMESH_SubMeshObj
665   Description : Base class for visualisation of submeshes and groups
666 */
667
668 //=================================================================================
669 // function : SMESH_SubMeshObj
670 // purpose  : Constructor
671 //=================================================================================
672 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
673 {
674   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
675   
676   myMeshObj = theMeshObj;
677 }
678
679 SMESH_SubMeshObj::~SMESH_SubMeshObj()
680 {
681 }
682
683 //=================================================================================
684 // function : GetElemDimension
685 // purpose  : Get dimension of element
686 //=================================================================================
687 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
688 {
689   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
690 }
691
692 //=================================================================================
693 // function : UpdateFunctor
694 // purpose  : Update functor in accordance with current mesh
695 //=================================================================================
696 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
697 {
698   theFunctor->SetMesh( myMeshObj->GetMesh() );
699 }
700
701 //=================================================================================
702 // function : Update
703 // purpose  : Update mesh object and fill grid with new values 
704 //=================================================================================
705 bool SMESH_SubMeshObj::Update( int theIsClear )
706 {
707   bool changed = myMeshObj->Update( theIsClear );
708   buildPrs();
709   return changed;
710 }
711
712
713 /*
714   Class       : SMESH_GroupObj
715   Description : Class for visualisation of groups
716 */
717
718 //=================================================================================
719 // function : SMESH_GroupObj
720 // purpose  : Constructor
721 //=================================================================================
722 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
723                                 SMESH_MeshObj*             theMeshObj )
724 : SMESH_SubMeshObj( theMeshObj ),
725   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
726 {
727   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
728   myGroupServer->Register();
729 }
730
731 SMESH_GroupObj::~SMESH_GroupObj()
732 {
733   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
734   myGroupServer->Destroy();
735 }
736
737 //=================================================================================
738 // function : IsNodePrs
739 // purpose  : Return true if node presentation is used
740 //=================================================================================
741 bool SMESH_GroupObj::IsNodePrs() const
742 {
743   return myGroupServer->GetType() == SMESH::NODE;
744 }
745
746 //=================================================================================
747 // function : getNodesFromElems
748 // purpose  : Retrieve nodes from elements
749 //=================================================================================
750 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
751                               const SMDS_Mesh*                    theMesh,
752                               std::list<const SMDS_MeshElement*>& theResList )
753 {
754   set<const SMDS_MeshElement*> aNodeSet;
755
756   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
757   {
758     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
759     if ( anElem != 0 )
760     {
761       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
762       while ( anIter->more() )
763       {
764         const SMDS_MeshElement* aNode = anIter->next();
765         if ( aNode != 0 )
766           aNodeSet.insert( aNode );
767       }
768     }
769   }
770
771   set<const SMDS_MeshElement*>::const_iterator anIter;
772   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
773     theResList.push_back( *anIter );
774
775   return theResList.size();    
776 }
777
778 //=================================================================================
779 // function : getPointers
780 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
781 //=================================================================================
782 static int getPointers( const SMDSAbs_ElementType            theRequestType,
783                         SMESH::long_array_var&              theElemIds,
784                         const SMDS_Mesh*                    theMesh,
785                         std::list<const SMDS_MeshElement*>& theResList )
786 {
787   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
788   {
789     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
790       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
791
792     if ( anElem != 0 )
793       theResList.push_back( anElem );
794   }
795
796   return theResList.size();
797 }
798
799
800 //=================================================================================
801 // function : GetEntities
802 // purpose  : Get entities of specified type. Return number of entities
803 //=================================================================================
804 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
805 {
806   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
807     return myGroupServer->Size();
808   }
809   return 0;
810 }
811
812 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
813 {
814   theResList.clear();
815   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
816   
817   if ( myGroupServer->Size() == 0 || aMesh == 0 )
818     return 0;
819
820   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
821   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
822
823   if ( aGrpType == theType )
824     return getPointers( theType, anIds, aMesh, theResList );
825   else if ( theType == SMDSAbs_Node )
826     return getNodesFromElems( anIds, aMesh, theResList );
827   else
828     return 0;
829 }    
830
831
832
833 /*
834   Class       : SMESH_subMeshObj
835   Description : Class for visualisation of submeshes
836 */
837
838 //=================================================================================
839 // function : SMESH_subMeshObj
840 // purpose  : Constructor
841 //=================================================================================
842 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
843                                     SMESH_MeshObj*           theMeshObj )
844 : SMESH_SubMeshObj( theMeshObj ),
845   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
846 {
847   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
848   
849   mySubMeshServer->Register();
850 }
851
852 SMESH_subMeshObj::~SMESH_subMeshObj()
853 {
854   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
855   mySubMeshServer->Destroy();
856 }
857
858 //=================================================================================
859 // function : GetEntities
860 // purpose  : Get entities of specified type. Return number of entities
861 //=================================================================================
862 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
863 {
864   switch ( theType )
865   {
866     case SMDSAbs_Node:
867     {
868       return mySubMeshServer->GetNumberOfNodes( false );
869     }
870     break;
871     case SMDSAbs_Edge:
872     case SMDSAbs_Face:
873     case SMDSAbs_Volume:
874     {
875       SMESH::long_array_var anIds = 
876         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
877       return anIds->length();
878     }
879     default:
880       return 0;
881     break;
882   }
883 }
884
885 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
886 {
887   theResList.clear();
888
889   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
890   if ( aMesh == 0 )
891     return 0;
892
893   bool isNodal = IsNodePrs();
894
895   if ( isNodal )
896   {
897     if ( theType == SMDSAbs_Node )
898     {
899       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
900       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
901     }
902   }
903   else
904   {
905     if ( theType == SMDSAbs_Node )
906     {
907       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
908       return getNodesFromElems( anIds, aMesh, theResList );
909     }
910     else
911     {
912       SMESH::long_array_var anIds = 
913         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
914       return getPointers( theType, anIds, aMesh, theResList );
915     }
916   }
917
918   return 0;
919 }
920
921 //=================================================================================
922 // function : IsNodePrs
923 // purpose  : Return true if node presentation is used
924 //=================================================================================
925 bool SMESH_subMeshObj::IsNodePrs() const
926 {
927   return mySubMeshServer->GetNumberOfElements() == 0;
928 }