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