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