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