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