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