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