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