Salome HOME
23418: [OCC] Mesh: Minimization of memory usage of SMESH
[modules/smesh.git] / src / OBJECT / SMESH_Object.cxx
1 // Copyright (C) 2007-2016  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, or (at your option) any later version.
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_BallElement.hxx"
32 #include "SMDS_Mesh.hxx"
33 #include "SMDS_MeshCell.hxx"
34 #include "SMESHDS_Mesh.hxx"
35 #include "SMESHDS_Script.hxx"
36 #include "SMESH_Actor.h"
37 #include "SMESH_ControlsDef.hxx"
38
39 #include <SalomeApp_Application.h>
40 #include <VTKViewer_ExtractUnstructuredGrid.h>
41 #include <VTKViewer_CellLocationsArray.h>
42
43 #include CORBA_SERVER_HEADER(SMESH_Gen)
44 #include CORBA_SERVER_HEADER(SALOME_Exception)
45
46 #include <vtkCell.h>
47 #include <vtkIdList.h>
48 #include <vtkCellArray.h>
49 #include <vtkUnsignedCharArray.h>
50 #include <vtkCellData.h>
51 #include <vtkUnstructuredGrid.h>
52
53 #include <memory>
54 #include <sstream>      
55 #include <stdexcept>
56 #include <set>
57
58 #include <utilities.h>
59
60 using namespace std;
61
62 #ifndef EXCEPTION
63 #define EXCEPTION(TYPE, MSG) {\
64   std::ostringstream aStream;\
65   aStream<<__FILE__<<"["<<__LINE__<<"]::"<<MSG;\
66   throw TYPE(aStream.str());\
67 }
68 #endif
69
70 #ifdef _DEBUG_
71 static int MYDEBUG = 0;
72 static int MYDEBUGWITHFILES = 0;//1;
73 #else
74 static int MYDEBUG = 0;
75 static int MYDEBUGWITHFILES = 0;
76 #endif
77
78
79 /*
80   Class       : SMESH_VisualObjDef
81   Description : Base class for all mesh objects to be visuilised
82 */
83 //=================================================================================
84 // functions : SMESH_VisualObjDef
85 // purpose   : Constructor
86 //=================================================================================
87 SMESH_VisualObjDef::SMESH_VisualObjDef()
88 {
89   if ( MYDEBUG ) MESSAGE("-------------------------------SMESH_VisualObjDef::SMESH_VisualObjDef");
90   myGrid = vtkUnstructuredGrid::New();
91   myLocalGrid = false;
92   ClearEntitiesFlags();
93   SMESH::GetEntitiesFromObject(NULL);
94 }
95 SMESH_VisualObjDef::~SMESH_VisualObjDef()
96 {
97   if ( MYDEBUG ) MESSAGE("--------------------------------SMESH_VisualObjDef::~SMESH_VisualObjDef");
98   if ( MYDEBUG ) MESSAGE( "myGrid->GetReferenceCount() = " << myGrid->GetReferenceCount() );
99   myGrid->Delete();
100 }
101
102 //=================================================================================
103 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
104 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
105 //=================================================================================
106 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
107 {
108   if (myLocalGrid)
109   {
110     TMapOfIds::const_iterator i = myVTK2SMDSNodes.find(theVTKID);
111     return i == myVTK2SMDSNodes.end() ? -1 : i->second;
112   }
113   const SMDS_MeshNode* aNode = 0;
114   if( this->GetMesh() )
115     aNode = this->GetMesh()->FindNodeVtk( theVTKID );
116
117   return aNode ? aNode->GetID() : -1;
118 }
119
120 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
121 {
122   if (myLocalGrid)
123   {
124     TMapOfIds::const_iterator i = mySMDS2VTKNodes.find(theObjID);
125     return i == mySMDS2VTKNodes.end() ? -1 : i->second;
126   }
127
128   const SMDS_MeshNode* aNode = 0;
129   if( this->GetMesh() ) {
130     aNode = this->GetMesh()->FindNode(theObjID);
131   }
132   return aNode ? aNode->GetVtkID() : -1;
133 }
134
135 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
136 {
137   if (myLocalGrid)
138   {
139     TMapOfIds::const_iterator i = myVTK2SMDSElems.find(theVTKID);
140     return i == myVTK2SMDSElems.end() ? -1 : i->second;
141   }
142   return this->GetMesh()->FromVtkToSmds(theVTKID);
143 }
144
145 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
146 {
147   if (myLocalGrid)
148   {
149     TMapOfIds::const_iterator i = mySMDS2VTKElems.find(theObjID);
150     return i == mySMDS2VTKElems.end() ? -1 : i->second;
151   }
152
153   const SMDS_MeshElement* e = 0;
154   if ( this->GetMesh() )
155     e = this->GetMesh()->FindElement(theObjID);
156
157   return e ? e->GetVtkID() : -1;
158 }
159
160 //=================================================================================
161 // function : SMESH_VisualObjDef::createPoints
162 // purpose  : Create points from nodes
163 //=================================================================================
164 /*! fills a vtkPoints structure for a submesh.
165  *  fills a std::list of SMDS_MeshElements*, then extract the points.
166  *  fills also conversion id maps between SMDS and VTK.
167  */
168 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
169 {
170   if ( thePoints == 0 )
171     return;
172
173   TEntityList aNodes;
174   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
175   thePoints->SetNumberOfPoints( nbNodes );
176
177   int nbPoints = 0;
178
179   TEntityList::const_iterator anIter;
180   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
181   {
182     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
183     if ( aNode != 0 )
184     {
185       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
186       int anId = aNode->GetID();
187       mySMDS2VTKNodes.insert( mySMDS2VTKNodes.end(), std::make_pair( anId, nbPoints ));
188       myVTK2SMDSNodes.insert( myVTK2SMDSNodes.end(), std::make_pair( nbPoints, anId ));
189       nbPoints++;
190     }
191   }
192
193   if ( nbPoints != nbNodes )
194     thePoints->SetNumberOfPoints( nbPoints );
195 }
196
197 //=================================================================================
198 // function : buildPrs
199 // purpose  : create VTK cells( fill unstructured grid )
200 //=================================================================================
201 void SMESH_VisualObjDef::buildPrs(bool buildGrid)
202 {
203   if ( MYDEBUG ) MESSAGE("---------------------------SMESH_VisualObjDef::buildPrs " << buildGrid);
204   if (buildGrid)
205   {
206     myLocalGrid = true;
207     try
208     {
209       mySMDS2VTKNodes.clear();
210       myVTK2SMDSNodes.clear();
211       mySMDS2VTKElems.clear();
212       myVTK2SMDSElems.clear();
213
214       if ( IsNodePrs() )
215         buildNodePrs();
216       else
217         buildElemPrs();
218     }
219     catch(...)
220     {
221       mySMDS2VTKNodes.clear();
222       myVTK2SMDSNodes.clear();
223       mySMDS2VTKElems.clear();
224       myVTK2SMDSElems.clear();
225
226       myGrid->SetPoints( 0 );
227       myGrid->SetCells( 0, 0, 0, 0, 0 );
228       throw;
229     }
230   }
231   else
232   {
233     myLocalGrid = false;
234     if (!GetMesh()->IsCompacted())
235     {
236       NulData(); // detach from the SMDS grid to allow immediate memory de-allocation in compactMesh()
237       if ( MYDEBUG ) MESSAGE("*** buildPrs ==> compactMesh!");
238       GetMesh()->CompactMesh();
239       if ( SMESHDS_Mesh* m = dynamic_cast<SMESHDS_Mesh*>( GetMesh() )) // IPAL53915
240         m->GetScript()->SetModified(false); // drop IsModified set in compactMesh()
241     }
242     vtkUnstructuredGrid *theGrid = GetMesh()->GetGrid();
243     updateEntitiesFlags();
244     myGrid->ShallowCopy(theGrid);
245     //MESSAGE(myGrid->GetReferenceCount());
246     //MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
247     //MESSAGE( "Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints() );
248     if( MYDEBUGWITHFILES ) {
249       SMESH::WriteUnstructuredGrid( myGrid,"myPrs.vtu" );
250     }
251   }
252 }
253
254 //=================================================================================
255 // function : buildNodePrs
256 // purpose  : create VTK cells for nodes
257 //=================================================================================
258
259 void SMESH_VisualObjDef::buildNodePrs()
260 {
261   // PAL16631: without swap, bad_alloc is not thrown but hung up and crash instead,
262   // so check remaining memory size for safety
263   SMDS_Mesh::CheckMemory(); // PAL16631
264   vtkPoints* aPoints = vtkPoints::New();
265   createPoints( aPoints );
266   SMDS_Mesh::CheckMemory();
267   myGrid->SetPoints( aPoints );
268   aPoints->Delete();
269
270   myGrid->SetCells( 0, 0, 0, 0, 0 );
271 }
272
273 //=================================================================================
274 // function : buildElemPrs
275 // purpose  : Create VTK cells for elements
276 //=================================================================================
277
278 namespace{
279   typedef std::vector<const SMDS_MeshElement*> TConnect;
280
281   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
282                  TConnect& theConnect)
283   {
284     theConnect.clear();
285     for(; theNodesIter->more();)
286       theConnect.push_back(theNodesIter->next());
287     return theConnect.size();
288   }
289   
290   inline 
291   void SetId(vtkIdList *theIdList, 
292              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
293              const TConnect& theConnect, 
294              int thePosition,
295              int theId)
296   {
297     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
298   }
299
300 }
301
302
303 void SMESH_VisualObjDef::buildElemPrs()
304 {
305   // Create points
306
307   vtkPoints* aPoints = vtkPoints::New();
308   createPoints( aPoints );
309   myGrid->SetPoints( aPoints );
310   aPoints->Delete();
311
312   if ( MYDEBUG )
313     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
314
315   // Calculate cells size
316
317   const int nbTypes = 5;
318   static SMDSAbs_ElementType aTypes[ nbTypes ] =
319     { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume, SMDSAbs_Ball, SMDSAbs_0DElement };
320
321   // get entity data
322   map<SMDSAbs_ElementType,int> nbEnts;
323   map<SMDSAbs_ElementType,TEntityList> anEnts;
324
325   vtkIdType aNbCells = 0;
326
327   for ( int i = 0; i < nbTypes; i++ )
328   {
329     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
330     aNbCells += nbEnts[ aTypes [ i ]];
331   }
332   // PAL16631: without swap, bad_alloc is not thrown but hung up and crash instead,
333   // so check remaining memory size for safety
334   SMDS_Mesh::CheckMemory(); // PAL16631
335
336   vtkIdType aCellsSize =  2 * nbEnts[ SMDSAbs_0DElement ] + 3 * nbEnts[ SMDSAbs_Edge ];
337   aCellsSize += 2 * nbEnts[ SMDSAbs_Ball ];
338
339   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
340   {
341     if ( nbEnts[ aTypes[ i ] ] )
342     {
343       const TEntityList& aList = anEnts[ aTypes[ i ] ];
344       TEntityList::const_iterator anIter;
345       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter ) {
346         if((*anIter)->GetEntityType() != SMDSEntity_Polyhedra &&
347            (*anIter)->GetEntityType() != SMDSEntity_Quad_Polyhedra) {
348           aCellsSize += (*anIter)->NbNodes() + 1;
349         }
350         // Special case for the VTK_POLYHEDRON:
351         // itsinput cellArray is of special format.
352         //  [nCellFaces, nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...]
353         else {
354           if ( const SMDS_MeshVolume* ph = SMDS_Mesh::DownCast<SMDS_MeshVolume>( *anIter )) {
355             int nbFaces = ph->NbFaces();
356             aCellsSize += (1 + ph->NbFaces());
357             for( int i = 1; i <= nbFaces; i++ )
358               aCellsSize += ph->NbFaceNodes(i);
359           }
360         }
361       }
362     }
363   }
364   if ( MYDEBUG )
365     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
366
367   // Create cells
368
369   vtkCellArray* aConnectivity = vtkCellArray::New();
370   aConnectivity->Allocate( aCellsSize, 0 );
371
372   SMDS_Mesh::CheckMemory(); // PAL16631
373
374   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
375   aCellTypesArray->SetNumberOfComponents( 1 );
376   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
377
378   SMDS_Mesh::CheckMemory(); // PAL16631
379
380   vtkIdList *anIdList = vtkIdList::New();
381   vtkIdType iElem = 0;
382
383   TConnect aConnect;
384   aConnect.reserve(VTK_CELL_SIZE);
385
386   SMDS_Mesh::CheckMemory(); // PAL16631
387   bool hasBalls = nbEnts[ SMDSAbs_Ball ] > 0;
388   vtkDataArray* aScalars = 0;
389   if(hasBalls) {
390     aScalars = vtkDataArray::CreateDataArray(VTK_DOUBLE);
391     aScalars->SetNumberOfComponents(1);
392     aScalars->SetNumberOfTuples(aNbCells);
393   }
394   for ( int i = 0; i < nbTypes; i++ ) // iterate through all types of elements
395   {
396     if ( nbEnts[ aTypes[ i ] ] > 0 ) {
397
398       const SMDSAbs_ElementType& aType = aTypes[ i ];
399       const TEntityList& aList = anEnts[ aType ];
400       TEntityList::const_iterator anIter;
401       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
402       {
403         const SMDS_MeshElement* anElem = *anIter;
404
405         vtkIdType aNbNodes = anElem->NbNodes();
406         anIdList->SetNumberOfIds( aNbNodes );
407         const vtkIdType vtkElemType = SMDS_MeshCell::toVtkType( anElem->GetEntityType() );
408
409         int anId = anElem->GetID();
410
411         mySMDS2VTKElems.insert( mySMDS2VTKElems.end(), std::make_pair( anId, iElem ));
412         myVTK2SMDSElems.insert( myVTK2SMDSElems.end(), std::make_pair( iElem, anId ));
413
414         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
415         {
416           // Convert connectivities from SMDS to VTK
417
418           if (aType == SMDSAbs_Volume && anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
419             anIdList->Reset();
420             if ( const SMDS_MeshVolume* ph = SMDS_Mesh::DownCast<SMDS_MeshVolume>( anElem )) {
421               int nbFaces = ph->NbFaces();
422               anIdList->InsertNextId(nbFaces);
423               for( int i = 1; i <= nbFaces; i++ ) {
424                 anIdList->InsertNextId(ph->NbFaceNodes(i));
425                 for(int j = 1; j <= ph->NbFaceNodes(i); j++) {
426                   if ( const SMDS_MeshNode* n = ph->GetFaceNode( i, j ))
427                     anIdList->InsertNextId( mySMDS2VTKNodes[ n->GetID() ]);
428                 }
429               }
430             }
431           }
432           else {
433             const std::vector<int>& aConnectivities =
434               SMDS_MeshCell::toVtkOrder( VTKCellType( vtkElemType ));
435             if (aConnectivities.size() > 0) {
436               aConnect.clear();
437               GetConnect(aNodesIter,aConnect);
438               for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
439                 SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
440             }
441             else {
442               for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
443                 const SMDS_MeshElement* aNode = aNodesIter->next();
444                 anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
445               }
446             }
447           }
448         }
449         vtkIdType aCurId = aConnectivity->InsertNextCell( anIdList );
450         aCellTypesArray->InsertNextValue( vtkElemType );
451         
452         //Store diameters of the balls
453         if(aScalars) {
454           double aDiam = 0;
455           if (const SMDS_BallElement* ball = SMDS_Mesh::DownCast<SMDS_BallElement>(anElem) )
456             aDiam = ball->GetDiameter();
457           aScalars->SetTuple(aCurId,&aDiam);
458         }
459
460         iElem++;
461       }
462     }
463     SMDS_Mesh::CheckMemory(); // PAL16631
464   }
465
466   // Insert cells in grid
467
468   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
469   aCellLocationsArray->SetNumberOfComponents( 1 );
470   aCellLocationsArray->SetNumberOfTuples( aNbCells );
471
472   SMDS_Mesh::CheckMemory(); // PAL16631
473
474   aConnectivity->InitTraversal();
475   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
476     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
477
478   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
479   myGrid->GetCellData()->SetScalars(aScalars);
480
481   aCellLocationsArray->Delete();
482   aCellTypesArray->Delete();
483   aConnectivity->Delete();
484   anIdList->Delete();
485
486   SMDS_Mesh::CheckMemory(); // PAL16631
487 }
488
489 //=================================================================================
490 // function : GetEdgeNodes
491 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 0 )
492 //=================================================================================
493 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
494                                        const int theEdgeNum,
495                                        int&      theNodeId1,
496                                        int&      theNodeId2 ) const
497 {
498   const SMDS_Mesh* aMesh = GetMesh();
499   if ( aMesh == 0 )
500     return false;
501     
502   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
503   if ( anElem == 0 )
504     return false;
505     
506   int nbNodes = anElem->NbCornerNodes();
507
508   if (( theEdgeNum < 0 || theEdgeNum > 3 ) ||
509       ( nbNodes != 3 && nbNodes != 4 ) ||
510       ( theEdgeNum >= nbNodes ))
511     return false;
512
513   theNodeId1 = anElem->GetNode(  theEdgeNum                 )->GetID();
514   theNodeId2 = anElem->GetNode(( theEdgeNum + 1 ) % nbNodes )->GetID();
515
516   return true;
517 }
518
519 vtkUnstructuredGrid* SMESH_VisualObjDef::GetUnstructuredGrid()
520 {
521   if ( !myLocalGrid && !GetMesh()->IsCompacted() )
522   {
523     NulData(); // detach from the SMDS grid to allow immediate memory de-allocation in CompactMesh()
524     GetMesh()->CompactMesh();
525     if ( SMESHDS_Mesh* m = dynamic_cast<SMESHDS_Mesh*>( GetMesh() )) // IPAL53915
526       m->GetScript()->SetModified(false); // drop IsModified set in CompactMesh()
527     updateEntitiesFlags();
528     vtkUnstructuredGrid *theGrid = GetMesh()->GetGrid();
529     myGrid->ShallowCopy(theGrid);
530   }
531   return myGrid;
532 }
533
534
535 //=================================================================================
536 // function : IsValid
537 // purpose  : Return true if there are some entities
538 //=================================================================================
539 bool SMESH_VisualObjDef::IsValid() const
540 {
541   return ( GetNbEntities(SMDSAbs_0DElement) > 0 ||
542            GetNbEntities(SMDSAbs_Ball     ) > 0 ||
543            GetNbEntities(SMDSAbs_Edge     ) > 0 ||
544            GetNbEntities(SMDSAbs_Face     ) > 0 ||
545            GetNbEntities(SMDSAbs_Volume   ) > 0 ||
546            GetNbEntities(SMDSAbs_Node     ) > 0 );
547 }
548
549 //=================================================================================
550 // function : updateEntitiesFlags
551 // purpose  : Update entities flags
552 //=================================================================================
553 void SMESH_VisualObjDef::updateEntitiesFlags()
554 {
555   unsigned int tmp = myEntitiesState;
556   ClearEntitiesFlags();
557
558   map<SMDSAbs_ElementType,int> entities = SMESH::GetEntitiesFromObject(this);
559
560
561   if( myEntitiesCache[SMDSAbs_0DElement] != 0 ||
562       myEntitiesCache[SMDSAbs_0DElement] >= entities[SMDSAbs_0DElement] )
563     myEntitiesState &= ~SMESH_Actor::e0DElements;
564
565   if( myEntitiesCache[SMDSAbs_Ball] != 0 ||
566       myEntitiesCache[SMDSAbs_Ball] >= entities[SMDSAbs_Ball] )
567     myEntitiesState &= ~SMESH_Actor::eBallElem;
568
569   if( myEntitiesCache[SMDSAbs_Edge] != 0 ||
570       myEntitiesCache[SMDSAbs_Edge] >= entities[SMDSAbs_Edge] )
571     myEntitiesState &= ~SMESH_Actor::eEdges;
572
573   if( myEntitiesCache[SMDSAbs_Face] != 0 ||
574       myEntitiesCache[SMDSAbs_Face] >= entities[SMDSAbs_Face] )
575     myEntitiesState &= ~SMESH_Actor::eFaces;
576
577   if( myEntitiesCache[SMDSAbs_Volume] != 0 ||
578       myEntitiesCache[SMDSAbs_Volume] >= entities[SMDSAbs_Volume] )
579     myEntitiesState &= ~SMESH_Actor::eVolumes;
580
581   if( tmp != myEntitiesState ) {
582     myEntitiesFlag = true;
583   }
584
585   myEntitiesCache = entities;
586 }
587
588 //=================================================================================
589 // function : ClearEntitiesFlags
590 // purpose  : Clear the entities flags
591 //=================================================================================
592 void SMESH_VisualObjDef::ClearEntitiesFlags()
593 {
594   myEntitiesState = SMESH_Actor::eAllEntity;
595   myEntitiesFlag = false;
596 }
597
598 //=================================================================================
599 // function : GetEntitiesFlag
600 // purpose  : Return the entities flag
601 //=================================================================================
602 bool SMESH_VisualObjDef::GetEntitiesFlag()
603 {
604   return myEntitiesFlag;
605 }
606
607 //=================================================================================
608 // function : GetEntitiesState
609 // purpose  : Return the entities state
610 //=================================================================================
611 unsigned int SMESH_VisualObjDef::GetEntitiesState()
612 {
613   return myEntitiesState;
614 }
615
616 /*
617   Class       : SMESH_MeshObj
618   Description : Class for visualisation of mesh
619 */
620
621 //=================================================================================
622 // function : SMESH_MeshObj
623 // purpose  : Constructor
624 //=================================================================================
625 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
626   myClient(SalomeApp_Application::orb(),theMesh)
627 {
628         myEmptyGrid = 0;
629   if ( MYDEBUG ) 
630     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
631 }
632
633 //=================================================================================
634 // function : ~SMESH_MeshObj
635 // purpose  : Destructor
636 //=================================================================================
637 SMESH_MeshObj::~SMESH_MeshObj()
638 {
639   if ( MYDEBUG ) 
640     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
641   if ( myEmptyGrid )
642     myEmptyGrid->Delete();
643 }
644
645 //=================================================================================
646 // function : Update
647 // purpose  : Update mesh and fill grid with new values if necessary 
648 //=================================================================================
649 bool SMESH_MeshObj::Update( int theIsClear )
650 {
651   // Update SMDS_Mesh on client part
652   if ( MYDEBUG ) MESSAGE("SMESH_MeshObj::Update " << this);
653   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
654     if ( MYDEBUG ) MESSAGE("buildPrs");
655     buildPrs();  // Fill unstructured grid
656     return true;
657   }
658   return false;
659 }
660
661 bool SMESH_MeshObj::NulData()
662 {
663   if ( MYDEBUG ) MESSAGE ("SMESH_MeshObj::NulData() =============================================");
664   if (!myEmptyGrid)
665   {
666     myEmptyGrid = SMDS_UnstructuredGrid::New();
667     myEmptyGrid->Initialize();
668     myEmptyGrid->Allocate();
669     vtkPoints* points = vtkPoints::New();
670     points->SetNumberOfPoints(0);
671     myEmptyGrid->SetPoints( points );
672     points->Delete();
673     //myEmptyGrid->BuildLinks();
674   }
675   myGrid->ShallowCopy(myEmptyGrid);
676   return true;
677 }
678 //=================================================================================
679 // function : GetElemDimension
680 // purpose  : Get dimension of element
681 //=================================================================================
682 int SMESH_MeshObj::GetElemDimension( const int theObjId )
683 {
684   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
685   if ( anElem == 0 )
686     return 0;
687
688   int aType = anElem->GetType();
689   switch ( aType )
690   {
691     case SMDSAbs_0DElement : return 0;
692     case SMDSAbs_Ball : return 0;
693     case SMDSAbs_Edge  : return 1;
694     case SMDSAbs_Face  : return 2;
695     case SMDSAbs_Volume: return 3;
696     default            : return 0;
697   }
698 }
699
700 //=================================================================================
701 // function : GetEntities
702 // purpose  : Get entities of specified type. Return number of entities
703 //=================================================================================
704 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
705 {
706   switch ( theType )
707   {
708     case SMDSAbs_Node:
709     {
710       return myClient->NbNodes();
711     }
712     break;
713     case SMDSAbs_0DElement:
714     {
715       return myClient->Nb0DElements();
716     }
717     case SMDSAbs_Ball:
718     {
719       return myClient->NbBalls();
720     }
721     break;
722     case SMDSAbs_Edge:
723     {
724       return myClient->NbEdges();
725     }
726     break;
727     case SMDSAbs_Face:
728     {
729       return myClient->NbFaces();
730     }
731     break;
732     case SMDSAbs_Volume:
733     {
734       return myClient->NbVolumes();
735     }
736     break;
737     default:
738       return 0;
739     break;
740   }
741 }
742
743 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
744 {
745   theObjs.clear();
746
747   switch ( theType )
748   {
749     case SMDSAbs_Node:
750     {
751       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
752       while ( anIter->more() ) theObjs.push_back( anIter->next() );
753     }
754     break;
755     case SMDSAbs_0DElement:
756     {
757       SMDS_ElemIteratorPtr anIter = myClient->elementsIterator(SMDSAbs_0DElement);
758       while ( anIter->more() ) theObjs.push_back( anIter->next() );
759     }
760     break;
761     case SMDSAbs_Ball:
762     {
763       SMDS_ElemIteratorPtr anIter = myClient->elementGeomIterator(SMDSGeom_BALL);
764       while ( anIter->more() ) theObjs.push_back( anIter->next() );
765     }
766     break;
767     case SMDSAbs_Edge:
768     {
769       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
770       while ( anIter->more() ) theObjs.push_back( anIter->next() );
771     }
772     break;
773     case SMDSAbs_Face:
774     {
775       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
776       while ( anIter->more() ) theObjs.push_back( anIter->next() );
777     }
778     break;
779     case SMDSAbs_Volume:
780     {
781       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
782       while ( anIter->more() ) theObjs.push_back( anIter->next() );
783     }
784     break;
785     default:
786     break;
787   }
788
789   return theObjs.size();
790 }
791
792 //=================================================================================
793 // function : UpdateFunctor
794 // purpose  : Update functor in accordance with current mesh
795 //=================================================================================
796 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
797 {
798   theFunctor->SetMesh( GetMesh() );
799 }
800
801 //=================================================================================
802 // function : IsNodePrs
803 // purpose  : Return true if node presentation is used
804 //=================================================================================
805 bool SMESH_MeshObj::IsNodePrs() const
806 {
807   return myClient->Nb0DElements() + myClient->NbEdges() + myClient->NbFaces() + myClient->NbVolumes() + myClient->NbBalls() == 0 ;
808 }
809
810
811 /*
812   Class       : SMESH_SubMeshObj
813   Description : Base class for visualisation of submeshes and groups
814 */
815
816 //=================================================================================
817 // function : SMESH_SubMeshObj
818 // purpose  : Constructor
819 //=================================================================================
820 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
821 {
822   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
823   
824   myMeshObj = theMeshObj;
825 }
826
827 SMESH_SubMeshObj::~SMESH_SubMeshObj()
828 {
829 }
830
831 //=================================================================================
832 // function : GetElemDimension
833 // purpose  : Get dimension of element
834 //=================================================================================
835 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
836 {
837   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
838 }
839
840 //=================================================================================
841 // function : UpdateFunctor
842 // purpose  : Update functor in accordance with current mesh
843 //=================================================================================
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   if ( MYDEBUG ) 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   if ( theType == SMDSAbs_Node ) {
969     return myGroupServer->GetNumberOfNodes();
970   }
971   return 0;
972 }
973
974 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
975 {
976   theResList.clear();
977   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
978   
979   if ( aMesh == 0 )
980     return 0;
981
982   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
983   if ( aGrpType != theType && theType != SMDSAbs_Node )
984     return 0;
985
986   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
987   if ( anIds->length() == 0 )
988     return 0;
989
990   if ( aGrpType == theType )
991     return getPointers( theType, anIds, aMesh, theResList );
992   else if ( theType == SMDSAbs_Node )
993     return getNodesFromElems( anIds, aMesh, theResList );
994   else
995     return 0;
996 }
997
998
999
1000 /*
1001   Class       : SMESH_subMeshObj
1002   Description : Class for visualisation of submeshes
1003 */
1004
1005 //=================================================================================
1006 // function : SMESH_subMeshObj
1007 // purpose  : Constructor
1008 //=================================================================================
1009 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
1010                                     SMESH_MeshObj*           theMeshObj )
1011 : SMESH_SubMeshObj( theMeshObj ),
1012   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
1013 {
1014   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
1015   
1016   mySubMeshServer->Register();
1017 }
1018
1019 SMESH_subMeshObj::~SMESH_subMeshObj()
1020 {
1021   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1022   mySubMeshServer->UnRegister();
1023 }
1024
1025 //=================================================================================
1026 // function : GetEntities
1027 // purpose  : Get entities of specified type. Return number of entities
1028 //=================================================================================
1029 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1030 {
1031   switch ( theType )
1032   {
1033     case SMDSAbs_Node:
1034     {
1035       return mySubMeshServer->GetNumberOfNodes( /*all=*/true );
1036     }
1037     break;
1038     case SMDSAbs_Ball:
1039     case SMDSAbs_0DElement:
1040     case SMDSAbs_Edge:
1041     case SMDSAbs_Face:
1042     case SMDSAbs_Volume:
1043     {
1044       SMESH::long_array_var anIds = 
1045         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1046       return anIds->length();
1047     }
1048     default:
1049       return 0;
1050     break;
1051   }
1052 }
1053
1054 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1055 {
1056   theResList.clear();
1057
1058   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1059   if ( aMesh == 0 )
1060     return 0;
1061
1062   bool isNodal = IsNodePrs();
1063
1064   if ( isNodal )
1065   {
1066     if ( theType == SMDSAbs_Node )
1067     {
1068       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1069       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1070     }
1071   }
1072   else
1073   {
1074     if ( theType == SMDSAbs_Node )
1075     {
1076       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1077       return getNodesFromElems( anIds, aMesh, theResList );
1078     }
1079     else
1080     {
1081       SMESH::long_array_var anIds = 
1082         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1083       return getPointers( theType, anIds, aMesh, theResList );
1084     }
1085   }
1086
1087   return 0;
1088 }
1089
1090 //=================================================================================
1091 // function : IsNodePrs
1092 // purpose  : Return true if node presentation is used
1093 //=================================================================================
1094 bool SMESH_subMeshObj::IsNodePrs() const
1095 {
1096   return mySubMeshServer->GetNumberOfElements() == 0;
1097 }