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