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