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