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