Salome HOME
IPAL54382: TC8.5.0: SMESH: Regression in mesh visualization (body fitting)
[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() || // !IsCompacted() is needed ???
522                          GetMesh()->GetGrid()->GetMTime() > myGrid->GetMTime() ))
523   {
524     NulData(); // detach from the SMDS grid to allow immediate memory de-allocation in CompactMesh()
525     GetMesh()->CompactMesh();
526     if ( SMESHDS_Mesh* m = dynamic_cast<SMESHDS_Mesh*>( GetMesh() )) // IPAL53915
527       m->GetScript()->SetModified(false); // drop IsModified set in CompactMesh()
528     updateEntitiesFlags();
529     vtkUnstructuredGrid *theGrid = GetMesh()->GetGrid();
530     myGrid->ShallowCopy(theGrid);
531   }
532   return myGrid;
533 }
534
535
536 //=================================================================================
537 // function : IsValid
538 // purpose  : Return true if there are some entities
539 //=================================================================================
540 bool SMESH_VisualObjDef::IsValid() const
541 {
542   return ( GetNbEntities(SMDSAbs_0DElement) > 0 ||
543            GetNbEntities(SMDSAbs_Ball     ) > 0 ||
544            GetNbEntities(SMDSAbs_Edge     ) > 0 ||
545            GetNbEntities(SMDSAbs_Face     ) > 0 ||
546            GetNbEntities(SMDSAbs_Volume   ) > 0 ||
547            GetNbEntities(SMDSAbs_Node     ) > 0 );
548 }
549
550 //=================================================================================
551 // function : updateEntitiesFlags
552 // purpose  : Update entities flags
553 //=================================================================================
554 void SMESH_VisualObjDef::updateEntitiesFlags()
555 {
556   unsigned int tmp = myEntitiesState;
557   ClearEntitiesFlags();
558
559   map<SMDSAbs_ElementType,int> entities = SMESH::GetEntitiesFromObject(this);
560
561
562   if( myEntitiesCache[SMDSAbs_0DElement] != 0 ||
563       myEntitiesCache[SMDSAbs_0DElement] >= entities[SMDSAbs_0DElement] )
564     myEntitiesState &= ~SMESH_Actor::e0DElements;
565
566   if( myEntitiesCache[SMDSAbs_Ball] != 0 ||
567       myEntitiesCache[SMDSAbs_Ball] >= entities[SMDSAbs_Ball] )
568     myEntitiesState &= ~SMESH_Actor::eBallElem;
569
570   if( myEntitiesCache[SMDSAbs_Edge] != 0 ||
571       myEntitiesCache[SMDSAbs_Edge] >= entities[SMDSAbs_Edge] )
572     myEntitiesState &= ~SMESH_Actor::eEdges;
573
574   if( myEntitiesCache[SMDSAbs_Face] != 0 ||
575       myEntitiesCache[SMDSAbs_Face] >= entities[SMDSAbs_Face] )
576     myEntitiesState &= ~SMESH_Actor::eFaces;
577
578   if( myEntitiesCache[SMDSAbs_Volume] != 0 ||
579       myEntitiesCache[SMDSAbs_Volume] >= entities[SMDSAbs_Volume] )
580     myEntitiesState &= ~SMESH_Actor::eVolumes;
581
582   if( tmp != myEntitiesState ) {
583     myEntitiesFlag = true;
584   }
585
586   myEntitiesCache = entities;
587 }
588
589 //=================================================================================
590 // function : ClearEntitiesFlags
591 // purpose  : Clear the entities flags
592 //=================================================================================
593 void SMESH_VisualObjDef::ClearEntitiesFlags()
594 {
595   myEntitiesState = SMESH_Actor::eAllEntity;
596   myEntitiesFlag = false;
597 }
598
599 //=================================================================================
600 // function : GetEntitiesFlag
601 // purpose  : Return the entities flag
602 //=================================================================================
603 bool SMESH_VisualObjDef::GetEntitiesFlag()
604 {
605   return myEntitiesFlag;
606 }
607
608 //=================================================================================
609 // function : GetEntitiesState
610 // purpose  : Return the entities state
611 //=================================================================================
612 unsigned int SMESH_VisualObjDef::GetEntitiesState()
613 {
614   return myEntitiesState;
615 }
616
617 /*
618   Class       : SMESH_MeshObj
619   Description : Class for visualisation of mesh
620 */
621
622 //=================================================================================
623 // function : SMESH_MeshObj
624 // purpose  : Constructor
625 //=================================================================================
626 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
627   myClient(SalomeApp_Application::orb(),theMesh)
628 {
629         myEmptyGrid = 0;
630   if ( MYDEBUG ) 
631     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
632 }
633
634 //=================================================================================
635 // function : ~SMESH_MeshObj
636 // purpose  : Destructor
637 //=================================================================================
638 SMESH_MeshObj::~SMESH_MeshObj()
639 {
640   if ( MYDEBUG ) 
641     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
642   if ( myEmptyGrid )
643     myEmptyGrid->Delete();
644 }
645
646 //=================================================================================
647 // function : Update
648 // purpose  : Update mesh and fill grid with new values if necessary 
649 //=================================================================================
650 bool SMESH_MeshObj::Update( int theIsClear )
651 {
652   // Update SMDS_Mesh on client part
653   if ( MYDEBUG ) MESSAGE("SMESH_MeshObj::Update " << this);
654   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
655     if ( MYDEBUG ) MESSAGE("buildPrs");
656     buildPrs();  // Fill unstructured grid
657     return true;
658   }
659   return false;
660 }
661
662 bool SMESH_MeshObj::NulData()
663 {
664   if ( MYDEBUG ) MESSAGE ("SMESH_MeshObj::NulData() =============================================");
665   if (!myEmptyGrid)
666   {
667     myEmptyGrid = SMDS_UnstructuredGrid::New();
668     myEmptyGrid->Initialize();
669     myEmptyGrid->Allocate();
670     vtkPoints* points = vtkPoints::New();
671     points->SetNumberOfPoints(0);
672     myEmptyGrid->SetPoints( points );
673     points->Delete();
674     //myEmptyGrid->BuildLinks();
675   }
676   myGrid->ShallowCopy(myEmptyGrid);
677   return true;
678 }
679 //=================================================================================
680 // function : GetElemDimension
681 // purpose  : Get dimension of element
682 //=================================================================================
683 int SMESH_MeshObj::GetElemDimension( const int theObjId )
684 {
685   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
686   if ( anElem == 0 )
687     return 0;
688
689   int aType = anElem->GetType();
690   switch ( aType )
691   {
692     case SMDSAbs_0DElement : return 0;
693     case SMDSAbs_Ball : return 0;
694     case SMDSAbs_Edge  : return 1;
695     case SMDSAbs_Face  : return 2;
696     case SMDSAbs_Volume: return 3;
697     default            : return 0;
698   }
699 }
700
701 //=================================================================================
702 // function : GetEntities
703 // purpose  : Get entities of specified type. Return number of entities
704 //=================================================================================
705 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
706 {
707   switch ( theType )
708   {
709     case SMDSAbs_Node:
710     {
711       return myClient->NbNodes();
712     }
713     break;
714     case SMDSAbs_0DElement:
715     {
716       return myClient->Nb0DElements();
717     }
718     case SMDSAbs_Ball:
719     {
720       return myClient->NbBalls();
721     }
722     break;
723     case SMDSAbs_Edge:
724     {
725       return myClient->NbEdges();
726     }
727     break;
728     case SMDSAbs_Face:
729     {
730       return myClient->NbFaces();
731     }
732     break;
733     case SMDSAbs_Volume:
734     {
735       return myClient->NbVolumes();
736     }
737     break;
738     default:
739       return 0;
740     break;
741   }
742 }
743
744 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
745 {
746   theObjs.clear();
747
748   switch ( theType )
749   {
750     case SMDSAbs_Node:
751     {
752       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
753       while ( anIter->more() ) theObjs.push_back( anIter->next() );
754     }
755     break;
756     case SMDSAbs_0DElement:
757     {
758       SMDS_ElemIteratorPtr anIter = myClient->elementsIterator(SMDSAbs_0DElement);
759       while ( anIter->more() ) theObjs.push_back( anIter->next() );
760     }
761     break;
762     case SMDSAbs_Ball:
763     {
764       SMDS_ElemIteratorPtr anIter = myClient->elementGeomIterator(SMDSGeom_BALL);
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() + myClient->NbEdges() + myClient->NbFaces() + myClient->NbVolumes() + myClient->NbBalls() == 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
846 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
847 {
848   theFunctor->SetMesh( myMeshObj->GetMesh() );
849 }
850
851 //=================================================================================
852 // function : Update
853 // purpose  : Update mesh object and fill grid with new values 
854 //=================================================================================
855 bool SMESH_SubMeshObj::Update( int theIsClear )
856 {
857   if ( MYDEBUG ) MESSAGE("SMESH_SubMeshObj::Update " << this)
858   bool changed = myMeshObj->Update( theIsClear );
859   buildPrs(true);
860   return changed;
861 }
862
863
864 /*
865   Class       : SMESH_GroupObj
866   Description : Class for visualisation of groups
867 */
868
869 //=================================================================================
870 // function : SMESH_GroupObj
871 // purpose  : Constructor
872 //=================================================================================
873 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
874                                 SMESH_MeshObj*             theMeshObj )
875 : SMESH_SubMeshObj( theMeshObj ),
876   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
877 {
878   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
879   myGroupServer->Register();
880 }
881
882 SMESH_GroupObj::~SMESH_GroupObj()
883 {
884   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
885   myGroupServer->UnRegister();
886 }
887
888 //=================================================================================
889 // function : IsNodePrs
890 // purpose  : Return true if node presentation is used
891 //=================================================================================
892 bool SMESH_GroupObj::IsNodePrs() const
893 {
894   return myGroupServer->GetType() == SMESH::NODE;
895 }
896
897 //=================================================================================
898 // function : GetElementType
899 // purpose  : Return type of elements of group
900 //=================================================================================
901 SMDSAbs_ElementType SMESH_GroupObj::GetElementType() const
902 {
903   return SMDSAbs_ElementType(myGroupServer->GetType());
904 }
905
906 //=================================================================================
907 // function : getNodesFromElems
908 // purpose  : Retrieve nodes from elements
909 //=================================================================================
910 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
911                               const SMDS_Mesh*                    theMesh,
912                               std::list<const SMDS_MeshElement*>& theResList )
913 {
914   set<const SMDS_MeshElement*> aNodeSet;
915
916   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
917   {
918     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
919     if ( anElem != 0 )
920     {
921       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
922       while ( anIter->more() )
923       {
924         const SMDS_MeshElement* aNode = anIter->next();
925         if ( aNode != 0 )
926           aNodeSet.insert( aNode );
927       }
928     }
929   }
930
931   set<const SMDS_MeshElement*>::const_iterator anIter;
932   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
933     theResList.push_back( *anIter );
934
935   return theResList.size();    
936 }
937
938 //=================================================================================
939 // function : getPointers
940 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
941 //=================================================================================
942 static int getPointers( const SMDSAbs_ElementType           theRequestType,
943                         SMESH::long_array_var&              theElemIds,
944                         const SMDS_Mesh*                    theMesh,
945                         std::list<const SMDS_MeshElement*>& theResList )
946 {
947   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
948   {
949     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
950       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
951
952     if ( anElem != 0 )
953       theResList.push_back( anElem );
954   }
955
956   return theResList.size();
957 }
958
959
960 //=================================================================================
961 // function : GetEntities
962 // purpose  : Get entities of specified type. Return number of entities
963 //=================================================================================
964 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
965 {
966   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType) {
967     return myGroupServer->Size();
968   }
969   if ( theType == SMDSAbs_Node ) {
970     return myGroupServer->GetNumberOfNodes();
971   }
972   return 0;
973 }
974
975 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
976 {
977   theResList.clear();
978   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
979   
980   if ( aMesh == 0 )
981     return 0;
982
983   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
984   if ( aGrpType != theType && theType != SMDSAbs_Node )
985     return 0;
986
987   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
988   if ( anIds->length() == 0 )
989     return 0;
990
991   if ( aGrpType == theType )
992     return getPointers( theType, anIds, aMesh, theResList );
993   else if ( theType == SMDSAbs_Node )
994     return getNodesFromElems( anIds, aMesh, theResList );
995   else
996     return 0;
997 }
998
999
1000
1001 /*
1002   Class       : SMESH_subMeshObj
1003   Description : Class for visualisation of submeshes
1004 */
1005
1006 //=================================================================================
1007 // function : SMESH_subMeshObj
1008 // purpose  : Constructor
1009 //=================================================================================
1010 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
1011                                     SMESH_MeshObj*           theMeshObj )
1012 : SMESH_SubMeshObj( theMeshObj ),
1013   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
1014 {
1015   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
1016   
1017   mySubMeshServer->Register();
1018 }
1019
1020 SMESH_subMeshObj::~SMESH_subMeshObj()
1021 {
1022   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1023   mySubMeshServer->UnRegister();
1024 }
1025
1026 //=================================================================================
1027 // function : GetEntities
1028 // purpose  : Get entities of specified type. Return number of entities
1029 //=================================================================================
1030 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1031 {
1032   switch ( theType )
1033   {
1034     case SMDSAbs_Node:
1035     {
1036       return mySubMeshServer->GetNumberOfNodes( /*all=*/true );
1037     }
1038     break;
1039     case SMDSAbs_Ball:
1040     case SMDSAbs_0DElement:
1041     case SMDSAbs_Edge:
1042     case SMDSAbs_Face:
1043     case SMDSAbs_Volume:
1044     {
1045       SMESH::long_array_var anIds = 
1046         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1047       return anIds->length();
1048     }
1049     default:
1050       return 0;
1051     break;
1052   }
1053 }
1054
1055 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1056 {
1057   theResList.clear();
1058
1059   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1060   if ( aMesh == 0 )
1061     return 0;
1062
1063   bool isNodal = IsNodePrs();
1064
1065   if ( isNodal )
1066   {
1067     if ( theType == SMDSAbs_Node )
1068     {
1069       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1070       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1071     }
1072   }
1073   else
1074   {
1075     if ( theType == SMDSAbs_Node )
1076     {
1077       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1078       return getNodesFromElems( anIds, aMesh, theResList );
1079     }
1080     else
1081     {
1082       SMESH::long_array_var anIds = 
1083         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1084       return getPointers( theType, anIds, aMesh, theResList );
1085     }
1086   }
1087
1088   return 0;
1089 }
1090
1091 //=================================================================================
1092 // function : IsNodePrs
1093 // purpose  : Return true if node presentation is used
1094 //=================================================================================
1095 bool SMESH_subMeshObj::IsNodePrs() const
1096 {
1097   return mySubMeshServer->GetNumberOfElements() == 0;
1098 }