Salome HOME
PAL13460 (PAL EDF 301 force the mesh to go through a point)
[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( const std::exception& exc )
215   {
216     INFOS("Follow exception was cought:\n\t"<<exc.what());
217   }
218   catch(...)
219   {
220     INFOS("Unknown exception was cought !!!");
221   }
222   
223   if( MYDEBUG ) MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
224   if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
225 }
226
227 //=================================================================================
228 // function : buildNodePrs
229 // purpose  : create VTK cells for nodes
230 //=================================================================================
231 void SMESH_VisualObjDef::buildNodePrs()
232 {
233   vtkPoints* aPoints = vtkPoints::New();
234   createPoints( aPoints );
235   myGrid->SetPoints( aPoints );
236   aPoints->Delete();
237
238   myGrid->SetCells( 0, 0, 0 );
239 }
240
241 //=================================================================================
242 // function : buildElemPrs
243 // purpose  : Create VTK cells for elements
244 //=================================================================================
245
246 namespace{
247   typedef std::vector<const SMDS_MeshElement*> TConnect;
248
249   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
250                  TConnect& theConnect)
251   {
252     theConnect.clear();
253     for(; theNodesIter->more();)
254       theConnect.push_back(theNodesIter->next());
255     return theConnect.size();
256   }
257   
258   inline 
259   void SetId(vtkIdList *theIdList, 
260              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
261              const TConnect& theConnect, 
262              int thePosition,
263              int theId)
264   {
265     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
266   }
267
268 }
269
270
271 void SMESH_VisualObjDef::buildElemPrs()
272 {
273   // Create points
274   
275   vtkPoints* aPoints = vtkPoints::New();
276   createPoints( aPoints );
277   myGrid->SetPoints( aPoints );
278   aPoints->Delete();
279     
280   if ( MYDEBUG )
281     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
282
283   // Calculate cells size
284
285   static SMDSAbs_ElementType aTypes[ 3 ] = { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
286
287   // get entity data
288   map<SMDSAbs_ElementType,int> nbEnts;
289   map<SMDSAbs_ElementType,TEntityList> anEnts;
290
291   for ( int i = 0; i <= 2; i++ )
292     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
293
294   vtkIdType aCellsSize =  3 * nbEnts[ SMDSAbs_Edge ];
295
296   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
297   {
298     if ( nbEnts[ aTypes[ i ] ] )
299     {
300       const TEntityList& aList = anEnts[ aTypes[ i ] ];
301       TEntityList::const_iterator anIter;
302       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
303         aCellsSize += (*anIter)->NbNodes() + 1;
304     }
305   }
306
307   vtkIdType aNbCells = nbEnts[ SMDSAbs_Edge ] + nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
308   
309   if ( MYDEBUG )
310     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
311
312   // Create cells
313   
314   vtkCellArray* aConnectivity = vtkCellArray::New();
315   aConnectivity->Allocate( aCellsSize, 0 );
316   
317   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
318   aCellTypesArray->SetNumberOfComponents( 1 );
319   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
320   
321   vtkIdList *anIdList = vtkIdList::New();
322   vtkIdType iElem = 0;
323
324   TConnect aConnect;
325   aConnect.reserve(VTK_CELL_SIZE);
326
327   for ( int i = 0; i <= 2; i++ ) // iterate through edges, faces and volumes
328   {
329     if( nbEnts[ aTypes[ i ] ] > 0 )
330     {
331       const SMDSAbs_ElementType& aType = aTypes[ i ];
332       const TEntityList& aList = anEnts[ aType ];
333       TEntityList::const_iterator anIter;
334       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
335       {
336         const SMDS_MeshElement* anElem = *anIter;
337         
338         vtkIdType aNbNodes = anElem->NbNodes();
339         anIdList->SetNumberOfIds( aNbNodes );
340
341         int anId = anElem->GetID();
342
343         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
344         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
345
346         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
347         switch(aType){
348         case SMDSAbs_Volume:{
349           std::vector<int> aConnectivities;
350           // Convertions connectivities from SMDS to VTK
351           if (anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
352
353             if ( const SMDS_PolyhedralVolumeOfNodes* ph =
354                  dynamic_cast<const SMDS_PolyhedralVolumeOfNodes*> (anElem))
355             {
356               aNbNodes = GetConnect(ph->uniqueNodesIterator(),aConnect);
357               anIdList->SetNumberOfIds( aNbNodes );
358             }
359             for (int k = 0; k < aNbNodes; k++)
360               aConnectivities.push_back(k);
361
362           } else if (aNbNodes == 4) {
363             static int anIds[] = {0,2,1,3};
364             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
365
366           } else if (aNbNodes == 5) {
367             static int anIds[] = {0,3,2,1,4};
368             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
369
370           } else if (aNbNodes == 6) {
371             static int anIds[] = {0,1,2,3,4,5};
372             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
373
374           }
375           else if (aNbNodes == 8) {
376             static int anIds[] = {0,3,2,1,4,7,6,5};
377             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
378
379           }
380           else if (aNbNodes == 10) {
381             static int anIds[] = {0,2,1,3,6,5,4,7,9,8};
382             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
383           }
384           else if (aNbNodes == 13) {
385             static int anIds[] = {0,3,2,1,4,8,7,6,5,9,12,11,10};
386             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
387           }
388           else if (aNbNodes == 15) {
389             static int anIds[] = {0,2,1,3,5,4,8,7,6,11,10,9,12,14,13};
390             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
391             //for (int k = 0; k < aNbNodes; k++) {
392             //  int nn = aConnectivities[k];
393             //  const SMDS_MeshNode* N = static_cast<const SMDS_MeshNode*> (aConnect[nn]);
394             //  cout<<"k="<<k<<"  N("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
395             //}
396           }
397           else if (aNbNodes == 20) {
398             static int anIds[] = {0,3,2,1,4,7,6,5,11,10,9,8,15,14,13,12,16,19,18,17};
399             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
400           }
401           else {
402           }
403
404           if ( aConnect.empty() )
405             GetConnect(aNodesIter,aConnect);
406
407           if (aConnectivities.size() > 0) {
408             for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
409               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
410           }
411           break;
412         }
413         default:
414           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
415             const SMDS_MeshElement* aNode = aNodesIter->next();
416             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
417           }
418         }
419
420         aConnectivity->InsertNextCell( anIdList );
421         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
422
423         iElem++;
424       }
425     }
426   }
427
428   // Insert cells in grid
429   
430   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
431   aCellLocationsArray->SetNumberOfComponents( 1 );
432   aCellLocationsArray->SetNumberOfTuples( aNbCells );
433   
434   aConnectivity->InitTraversal();
435   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
436     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
437
438   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
439   
440   aCellLocationsArray->Delete();
441   aCellTypesArray->Delete();
442   aConnectivity->Delete();
443   anIdList->Delete();
444 }
445
446 //=================================================================================
447 // function : GetEdgeNodes
448 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
449 //=================================================================================
450 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
451                                        const int theEdgeNum,
452                                        int&      theNodeId1,
453                                        int&      theNodeId2 ) const
454 {
455   const SMDS_Mesh* aMesh = GetMesh();
456   if ( aMesh == 0 )
457     return false;
458     
459   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
460   if ( anElem == 0 )
461     return false;
462     
463   int nbNodes = anElem->NbNodes();
464
465   if ( theEdgeNum < 0 || theEdgeNum > 3 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
466     return false;
467
468   int anIds[ nbNodes ];
469   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
470   int i = 0;
471   while( anIter->more() )
472     anIds[ i++ ] = anIter->next()->GetID();
473
474   if ( theEdgeNum < nbNodes - 1 )
475   {
476     theNodeId1 = anIds[ theEdgeNum ];
477     theNodeId2 = anIds[ theEdgeNum + 1 ];
478   }
479   else
480   {
481     theNodeId1 = anIds[ nbNodes - 1 ];
482     theNodeId2 = anIds[ 0 ];
483   }
484
485   return true;
486 }
487
488 /*
489   Class       : SMESH_MeshObj
490   Description : Class for visualisation of mesh
491 */
492
493 //=================================================================================
494 // function : SMESH_MeshObj
495 // purpose  : Constructor
496 //=================================================================================
497 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
498   myClient(SalomeApp_Application::orb(),theMesh)
499 {
500   if ( MYDEBUG ) 
501     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
502 }
503
504 //=================================================================================
505 // function : ~SMESH_MeshObj
506 // purpose  : Destructor
507 //=================================================================================
508 SMESH_MeshObj::~SMESH_MeshObj()
509 {
510   if ( MYDEBUG ) 
511     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
512 }
513
514 //=================================================================================
515 // function : Update
516 // purpose  : Update mesh and fill grid with new values if necessary 
517 //=================================================================================
518 void SMESH_MeshObj::Update( int theIsClear )
519 {
520   // Update SMDS_Mesh on client part
521   if ( myClient.Update(theIsClear) )
522     buildPrs();  // Fill unstructured grid
523 }
524
525 //=================================================================================
526 // function : GetElemDimension
527 // purpose  : Get dimension of element
528 //=================================================================================
529 int SMESH_MeshObj::GetElemDimension( const int theObjId )
530 {
531   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
532   if ( anElem == 0 )
533     return 0;
534
535   int aType = anElem->GetType();
536   switch ( aType )
537   {
538     case SMDSAbs_Edge  : return 1;
539     case SMDSAbs_Face  : return 2;
540     case SMDSAbs_Volume: return 3;
541     default            : return 0;
542   }
543 }
544
545 //=================================================================================
546 // function : GetEntities
547 // purpose  : Get entities of specified type. Return number of entities
548 //=================================================================================
549 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
550 {
551   switch ( theType )
552   {
553     case SMDSAbs_Node:
554     {
555       return myClient->NbNodes();
556     }
557     break;
558     case SMDSAbs_Edge:
559     {
560       return myClient->NbEdges();
561     }
562     break;
563     case SMDSAbs_Face:
564     {
565       return myClient->NbFaces();
566     }
567     break;
568     case SMDSAbs_Volume:
569     {
570       return myClient->NbVolumes();
571     }
572     break;
573     default:
574       return 0;
575     break;
576   }
577 }
578
579 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
580 {
581   theObjs.clear();
582
583   switch ( theType )
584   {
585     case SMDSAbs_Node:
586     {
587       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
588       while ( anIter->more() ) theObjs.push_back( anIter->next() );
589     }
590     break;
591     case SMDSAbs_Edge:
592     {
593       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
594       while ( anIter->more() ) theObjs.push_back( anIter->next() );
595     }
596     break;
597     case SMDSAbs_Face:
598     {
599       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
600       while ( anIter->more() ) theObjs.push_back( anIter->next() );
601     }
602     break;
603     case SMDSAbs_Volume:
604     {
605       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
606       while ( anIter->more() ) theObjs.push_back( anIter->next() );
607     }
608     break;
609     default:
610     break;
611   }
612
613   return theObjs.size();
614 }
615
616 //=================================================================================
617 // function : UpdateFunctor
618 // purpose  : Update functor in accordance with current mesh
619 //=================================================================================
620 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
621 {
622   theFunctor->SetMesh( GetMesh() );
623 }
624
625 //=================================================================================
626 // function : IsNodePrs
627 // purpose  : Return true if node presentation is used
628 //=================================================================================
629 bool SMESH_MeshObj::IsNodePrs() const
630 {
631   return myClient->NbEdges() == 0 &&myClient->NbFaces() == 0 && myClient->NbVolumes() == 0 ;
632 }
633
634
635 /*
636   Class       : SMESH_SubMeshObj
637   Description : Base class for visualisation of submeshes and groups
638 */
639
640 //=================================================================================
641 // function : SMESH_SubMeshObj
642 // purpose  : Constructor
643 //=================================================================================
644 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
645 {
646   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
647   
648   myMeshObj = theMeshObj;
649 }
650
651 SMESH_SubMeshObj::~SMESH_SubMeshObj()
652 {
653 }
654
655 //=================================================================================
656 // function : GetElemDimension
657 // purpose  : Get dimension of element
658 //=================================================================================
659 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
660 {
661   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
662 }
663
664 //=================================================================================
665 // function : UpdateFunctor
666 // purpose  : Update functor in accordance with current mesh
667 //=================================================================================
668 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
669 {
670   theFunctor->SetMesh( myMeshObj->GetMesh() );
671 }
672
673 //=================================================================================
674 // function : Update
675 // purpose  : Update mesh object and fill grid with new values 
676 //=================================================================================
677 void SMESH_SubMeshObj::Update( int theIsClear )
678 {
679   myMeshObj->Update( theIsClear );
680   buildPrs();
681 }
682
683
684 /*
685   Class       : SMESH_GroupObj
686   Description : Class for visualisation of groups
687 */
688
689 //=================================================================================
690 // function : SMESH_GroupObj
691 // purpose  : Constructor
692 //=================================================================================
693 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
694                                 SMESH_MeshObj*             theMeshObj )
695 : SMESH_SubMeshObj( theMeshObj ),
696   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
697 {
698   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
699   myGroupServer->Register();
700 }
701
702 SMESH_GroupObj::~SMESH_GroupObj()
703 {
704   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
705   myGroupServer->Destroy();
706 }
707
708 //=================================================================================
709 // function : IsNodePrs
710 // purpose  : Return true if node presentation is used
711 //=================================================================================
712 bool SMESH_GroupObj::IsNodePrs() const
713 {
714   return myGroupServer->GetType() == SMESH::NODE;
715 }
716
717 //=================================================================================
718 // function : getNodesFromElems
719 // purpose  : Retrieve nodes from elements
720 //=================================================================================
721 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
722                               const SMDS_Mesh*                    theMesh,
723                               std::list<const SMDS_MeshElement*>& theResList )
724 {
725   set<const SMDS_MeshElement*> aNodeSet;
726
727   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
728   {
729     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
730     if ( anElem != 0 )
731     {
732       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
733       while ( anIter->more() )
734       {
735         const SMDS_MeshElement* aNode = anIter->next();
736         if ( aNode != 0 )
737           aNodeSet.insert( aNode );
738       }
739     }
740   }
741
742   set<const SMDS_MeshElement*>::const_iterator anIter;
743   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
744     theResList.push_back( *anIter );
745
746   return theResList.size();    
747 }
748
749 //=================================================================================
750 // function : getPointers
751 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
752 //=================================================================================
753 static int getPointers( const SMDSAbs_ElementType            theRequestType,
754                         SMESH::long_array_var&              theElemIds,
755                         const SMDS_Mesh*                    theMesh,
756                         std::list<const SMDS_MeshElement*>& theResList )
757 {
758   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
759   {
760     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
761       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
762
763     if ( anElem != 0 )
764       theResList.push_back( anElem );
765   }
766
767   return theResList.size();
768 }
769
770
771 //=================================================================================
772 // function : GetEntities
773 // purpose  : Get entities of specified type. Return number of entities
774 //=================================================================================
775 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
776 {
777   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
778     return myGroupServer->Size();
779   }
780   return 0;
781 }
782
783 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
784 {
785   theResList.clear();
786   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
787   
788   if ( myGroupServer->Size() == 0 || aMesh == 0 )
789     return 0;
790
791   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
792   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
793
794   if ( aGrpType == theType )
795     return getPointers( theType, anIds, aMesh, theResList );
796   else if ( theType == SMDSAbs_Node )
797     return getNodesFromElems( anIds, aMesh, theResList );
798   else
799     return 0;
800 }    
801
802
803
804 /*
805   Class       : SMESH_subMeshObj
806   Description : Class for visualisation of submeshes
807 */
808
809 //=================================================================================
810 // function : SMESH_subMeshObj
811 // purpose  : Constructor
812 //=================================================================================
813 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
814                                     SMESH_MeshObj*           theMeshObj )
815 : SMESH_SubMeshObj( theMeshObj ),
816   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
817 {
818   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
819   
820   mySubMeshServer->Register();
821 }
822
823 SMESH_subMeshObj::~SMESH_subMeshObj()
824 {
825   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
826   mySubMeshServer->Destroy();
827 }
828
829 //=================================================================================
830 // function : GetEntities
831 // purpose  : Get entities of specified type. Return number of entities
832 //=================================================================================
833 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
834 {
835   switch ( theType )
836   {
837     case SMDSAbs_Node:
838     {
839       return mySubMeshServer->GetNumberOfNodes( false );
840     }
841     break;
842     case SMDSAbs_Edge:
843     case SMDSAbs_Face:
844     case SMDSAbs_Volume:
845     {
846       SMESH::long_array_var anIds = 
847         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
848       return anIds->length();
849     }
850     default:
851       return 0;
852     break;
853   }
854 }
855
856 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
857 {
858   theResList.clear();
859
860   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
861   if ( aMesh == 0 )
862     return 0;
863
864   bool isNodal = IsNodePrs();
865
866   if ( isNodal )
867   {
868     if ( theType == SMDSAbs_Node )
869     {
870       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
871       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
872     }
873   }
874   else
875   {
876     if ( theType == SMDSAbs_Node )
877     {
878       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
879       return getNodesFromElems( anIds, aMesh, theResList );
880     }
881     else
882     {
883       SMESH::long_array_var anIds = 
884         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
885       return getPointers( theType, anIds, aMesh, theResList );
886     }
887   }
888
889   return 0;
890 }
891
892 //=================================================================================
893 // function : IsNodePrs
894 // purpose  : Return true if node presentation is used
895 //=================================================================================
896 bool SMESH_subMeshObj::IsNodePrs() const
897 {
898   return mySubMeshServer->GetNumberOfElements() == 0;
899 }