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