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