Salome HOME
23081: [CEA 1496] Control merge nodes behaviour: set fixed nodes
[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_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 = 1;
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   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   MESSAGE("---------------------------------------------SMESH_VisualObjDef::~SMESH_VisualObjDef");
145   //if ( MYDEBUG )
146     MESSAGE( "~SMESH_MeshObj - myGrid->GetReferenceCount() = " << myGrid->GetReferenceCount() );
147   myGrid->Delete();
148 }
149
150 //=================================================================================
151 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
152 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
153 //=================================================================================
154 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
155 {
156   if (myLocalGrid)
157   {
158     TMapOfIds::const_iterator i = myVTK2SMDSNodes.find(theVTKID);
159     return i == myVTK2SMDSNodes.end() ? -1 : i->second;
160   }
161   const SMDS_MeshNode* aNode = 0;
162   if( this->GetMesh() )
163     aNode = this->GetMesh()->FindNodeVtk( theVTKID );
164
165   return aNode ? aNode->GetID() : -1;
166 }
167
168 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
169 {
170   if (myLocalGrid)
171   {
172     TMapOfIds::const_iterator i = mySMDS2VTKNodes.find(theObjID);
173     return i == mySMDS2VTKNodes.end() ? -1 : i->second;
174   }
175
176   const SMDS_MeshNode* aNode = 0;
177   if( this->GetMesh() ) {
178     aNode = this->GetMesh()->FindNode(theObjID);
179   }
180   return aNode ? aNode->getVtkId() : -1;
181 }
182
183 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
184 {
185   if (myLocalGrid)
186   {
187     TMapOfIds::const_iterator i = myVTK2SMDSElems.find(theVTKID);
188     return i == myVTK2SMDSElems.end() ? -1 : i->second;
189   }
190   return this->GetMesh()->fromVtkToSmds(theVTKID);
191 }
192
193 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
194 {
195   if (myLocalGrid)
196   {
197     TMapOfIds::const_iterator i = mySMDS2VTKElems.find(theObjID);
198     return i == mySMDS2VTKElems.end() ? -1 : i->second;
199   }
200
201   const SMDS_MeshElement* e = 0;
202   if ( this->GetMesh() )
203     e = this->GetMesh()->FindElement(theObjID);
204
205   return e ? e->getVtkId() : -1;
206 }
207
208 //=================================================================================
209 // function : SMESH_VisualObjDef::createPoints
210 // purpose  : Create points from nodes
211 //=================================================================================
212 /*! fills a vtkPoints structure for a submesh.
213  *  fills a std::list of SMDS_MeshElements*, then extract the points.
214  *  fills also conversion id maps between SMDS and VTK.
215  */
216 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
217 {
218   if ( thePoints == 0 )
219     return;
220
221   TEntityList aNodes;
222   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
223   thePoints->SetNumberOfPoints( nbNodes );
224
225   int nbPoints = 0;
226
227   TEntityList::const_iterator anIter;
228   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
229   {
230     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
231     if ( aNode != 0 )
232     {
233       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
234       int anId = aNode->GetID();
235       mySMDS2VTKNodes.insert( TMapOfIds::value_type( anId, nbPoints ) );
236       myVTK2SMDSNodes.insert( TMapOfIds::value_type( nbPoints, anId ) );
237       nbPoints++;
238     }
239   }
240
241   if ( nbPoints != nbNodes )
242     thePoints->SetNumberOfPoints( nbPoints );
243 }
244
245 //=================================================================================
246 // function : buildPrs
247 // purpose  : create VTK cells( fill unstructured grid )
248 //=================================================================================
249 void SMESH_VisualObjDef::buildPrs(bool buildGrid)
250 {
251   MESSAGE("----------------------------------------------------------SMESH_VisualObjDef::buildPrs " << buildGrid);
252   if (buildGrid)
253   {
254         myLocalGrid = true;
255         try
256         {
257                 mySMDS2VTKNodes.clear();
258                 myVTK2SMDSNodes.clear();
259                 mySMDS2VTKElems.clear();
260                 myVTK2SMDSElems.clear();
261
262                 if ( IsNodePrs() )
263                         buildNodePrs();
264                 else
265                         buildElemPrs();
266         }
267         catch(...)
268         {
269                 mySMDS2VTKNodes.clear();
270                 myVTK2SMDSNodes.clear();
271                 mySMDS2VTKElems.clear();
272                 myVTK2SMDSElems.clear();
273
274                 myGrid->SetPoints( 0 );
275                 myGrid->SetCells( 0, 0, 0, 0, 0 );
276                 throw;
277         }
278   }
279   else
280   {
281         myLocalGrid = false;
282         if (!GetMesh()->isCompacted())
283           {
284             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( TMapOfIds::value_type( anId, iElem ) );
458         myVTK2SMDSElems.insert( TMapOfIds::value_type( 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 1 )
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 || (nbNodes != 3 && nbNodes != 4) || theEdgeNum > nbNodes )
560     return false;
561
562   theNodeId1 = anElem->GetNode( theEdgeNum - 1       )->GetID();
563   theNodeId2 = anElem->GetNode( theEdgeNum % nbNodes )->GetID();
564
565   return true;
566 }
567
568 vtkUnstructuredGrid* SMESH_VisualObjDef::GetUnstructuredGrid()
569 {
570   if ( !myLocalGrid && !GetMesh()->isCompacted() )
571   {
572     GetMesh()->compactMesh();
573     updateEntitiesFlags();
574     vtkUnstructuredGrid *theGrid = GetMesh()->getGrid();
575     myGrid->ShallowCopy(theGrid);
576   }
577   return myGrid;
578 }
579
580
581 //=================================================================================
582 // function : IsValid
583 // purpose  : Return true if there are some entities
584 //=================================================================================
585 bool SMESH_VisualObjDef::IsValid() const
586 {
587         //MESSAGE("SMESH_VisualObjDef::IsValid");
588   return ( GetNbEntities(SMDSAbs_0DElement) > 0 || 
589            GetNbEntities(SMDSAbs_Ball     ) > 0 || 
590            GetNbEntities(SMDSAbs_Edge     ) > 0 || 
591            GetNbEntities(SMDSAbs_Face     ) > 0 ||
592            GetNbEntities(SMDSAbs_Volume   ) > 0 ||
593            GetNbEntities(SMDSAbs_Node     ) > 0 );
594 }
595
596 //=================================================================================
597 // function : updateEntitiesFlags
598 // purpose  : Update entities flags
599 //=================================================================================
600 void SMESH_VisualObjDef::updateEntitiesFlags() {
601
602         unsigned int tmp = myEntitiesState;
603         ClearEntitiesFlags();
604
605         map<SMDSAbs_ElementType,int> entities = SMESH::GetEntitiesFromObject(this);
606         
607
608         if( myEntitiesCache[SMDSAbs_0DElement] != 0 ||  myEntitiesCache[SMDSAbs_0DElement] >= entities[SMDSAbs_0DElement] )
609                 myEntitiesState &= ~SMESH_Actor::e0DElements;
610
611         if( myEntitiesCache[SMDSAbs_Ball] != 0 ||  myEntitiesCache[SMDSAbs_Ball] >= entities[SMDSAbs_Ball] )
612                 myEntitiesState &= ~SMESH_Actor::eBallElem;
613
614         if( myEntitiesCache[SMDSAbs_Edge] != 0 || myEntitiesCache[SMDSAbs_Edge] >= entities[SMDSAbs_Edge] )
615                 myEntitiesState &= ~SMESH_Actor::eEdges; 
616
617         if( myEntitiesCache[SMDSAbs_Face] != 0 || myEntitiesCache[SMDSAbs_Face] >= entities[SMDSAbs_Face] )
618                 myEntitiesState &= ~SMESH_Actor::eFaces; 
619
620         if( myEntitiesCache[SMDSAbs_Volume] != 0 || myEntitiesCache[SMDSAbs_Volume] >= entities[SMDSAbs_Volume] )
621                 myEntitiesState &= ~SMESH_Actor::eVolumes;
622
623         if( tmp != myEntitiesState ) {
624                 myEntitiesFlag = true;
625         }
626         
627         myEntitiesCache = entities;
628 }
629
630 //=================================================================================
631 // function : ClearEntitiesFlags
632 // purpose  : Clear the entities flags
633 //=================================================================================
634 void SMESH_VisualObjDef::ClearEntitiesFlags() {
635         myEntitiesState = SMESH_Actor::eAllEntity;
636         myEntitiesFlag = false;
637 }
638
639 //=================================================================================
640 // function : GetEntitiesFlag
641 // purpose  : Return the entities flag
642 //=================================================================================
643 bool SMESH_VisualObjDef::GetEntitiesFlag() {
644         return myEntitiesFlag;
645 }
646
647 //=================================================================================
648 // function : GetEntitiesState
649 // purpose  : Return the entities state
650 //=================================================================================
651 unsigned int SMESH_VisualObjDef::GetEntitiesState() {
652         return myEntitiesState;
653 }
654
655 /*
656   Class       : SMESH_MeshObj
657   Description : Class for visualisation of mesh
658 */
659
660 //=================================================================================
661 // function : SMESH_MeshObj
662 // purpose  : Constructor
663 //=================================================================================
664 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
665   myClient(SalomeApp_Application::orb(),theMesh)
666 {
667         myEmptyGrid = 0;
668   if ( MYDEBUG ) 
669     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
670 }
671
672 //=================================================================================
673 // function : ~SMESH_MeshObj
674 // purpose  : Destructor
675 //=================================================================================
676 SMESH_MeshObj::~SMESH_MeshObj()
677 {
678   if ( MYDEBUG ) 
679     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
680   if ( myEmptyGrid )
681     myEmptyGrid->Delete();
682 }
683
684 //=================================================================================
685 // function : Update
686 // purpose  : Update mesh and fill grid with new values if necessary 
687 //=================================================================================
688 bool SMESH_MeshObj::Update( int theIsClear )
689 {
690   // Update SMDS_Mesh on client part
691   MESSAGE("SMESH_MeshObj::Update " << this);
692   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
693     MESSAGE("buildPrs");
694     buildPrs();  // Fill unstructured grid
695     return true;
696   }
697   return false;
698 }
699
700 bool SMESH_MeshObj::NulData()
701 {
702         MESSAGE ("SMESH_MeshObj::NulData() ==================================================================================");
703         if (!myEmptyGrid)
704         {
705           myEmptyGrid = SMDS_UnstructuredGrid::New();
706           myEmptyGrid->Initialize();
707           myEmptyGrid->Allocate();
708           vtkPoints* points = vtkPoints::New();
709           points->SetNumberOfPoints(0);
710           myEmptyGrid->SetPoints( points );
711           points->Delete();
712           myEmptyGrid->BuildLinks();
713         }
714         myGrid->ShallowCopy(myEmptyGrid);
715         return true;
716 }
717 //=================================================================================
718 // function : GetElemDimension
719 // purpose  : Get dimension of element
720 //=================================================================================
721 int SMESH_MeshObj::GetElemDimension( const int theObjId )
722 {
723   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
724   if ( anElem == 0 )
725     return 0;
726
727   int aType = anElem->GetType();
728   switch ( aType )
729   {
730     case SMDSAbs_0DElement : return 0;
731     case SMDSAbs_Ball : return 0;
732     case SMDSAbs_Edge  : return 1;
733     case SMDSAbs_Face  : return 2;
734     case SMDSAbs_Volume: return 3;
735     default            : return 0;
736   }
737 }
738
739 //=================================================================================
740 // function : GetEntities
741 // purpose  : Get entities of specified type. Return number of entities
742 //=================================================================================
743 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
744 {
745   switch ( theType )
746   {
747     case SMDSAbs_Node:
748     {
749       return myClient->NbNodes();
750     }
751     break;
752     case SMDSAbs_0DElement:
753     {
754       return myClient->Nb0DElements();
755     }
756     case SMDSAbs_Ball:
757     {
758       return myClient->NbBalls();
759     }
760     break;
761     case SMDSAbs_Edge:
762     {
763       return myClient->NbEdges();
764     }
765     break;
766     case SMDSAbs_Face:
767     {
768       return myClient->NbFaces();
769     }
770     break;
771     case SMDSAbs_Volume:
772     {
773       return myClient->NbVolumes();
774     }
775     break;
776     default:
777       return 0;
778     break;
779   }
780 }
781
782 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
783 {
784   theObjs.clear();
785
786   switch ( theType )
787   {
788     case SMDSAbs_Node:
789     {
790       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
791       while ( anIter->more() ) theObjs.push_back( anIter->next() );
792     }
793     break;
794     case SMDSAbs_0DElement:
795     {
796       SMDS_ElemIteratorPtr anIter = myClient->elementsIterator(SMDSAbs_0DElement);
797       while ( anIter->more() ) theObjs.push_back( anIter->next() );
798     }
799     break;
800     case SMDSAbs_Ball:
801     {
802       SMDS_ElemIteratorPtr anIter = myClient->elementGeomIterator(SMDSGeom_BALL);
803       while ( anIter->more() ) theObjs.push_back( anIter->next() );
804     }
805     break;
806     case SMDSAbs_Edge:
807     {
808       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
809       while ( anIter->more() ) theObjs.push_back( anIter->next() );
810     }
811     break;
812     case SMDSAbs_Face:
813     {
814       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
815       while ( anIter->more() ) theObjs.push_back( anIter->next() );
816     }
817     break;
818     case SMDSAbs_Volume:
819     {
820       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
821       while ( anIter->more() ) theObjs.push_back( anIter->next() );
822     }
823     break;
824     default:
825     break;
826   }
827
828   return theObjs.size();
829 }
830
831 //=================================================================================
832 // function : UpdateFunctor
833 // purpose  : Update functor in accordance with current mesh
834 //=================================================================================
835 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
836 {
837   theFunctor->SetMesh( GetMesh() );
838 }
839
840 //=================================================================================
841 // function : IsNodePrs
842 // purpose  : Return true if node presentation is used
843 //=================================================================================
844 bool SMESH_MeshObj::IsNodePrs() const
845 {
846   return myClient->Nb0DElements() + myClient->NbEdges() + myClient->NbFaces() + myClient->NbVolumes() + myClient->NbBalls() == 0 ;
847 }
848
849
850 /*
851   Class       : SMESH_SubMeshObj
852   Description : Base class for visualisation of submeshes and groups
853 */
854
855 //=================================================================================
856 // function : SMESH_SubMeshObj
857 // purpose  : Constructor
858 //=================================================================================
859 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
860 {
861   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
862   
863   myMeshObj = theMeshObj;
864 }
865
866 SMESH_SubMeshObj::~SMESH_SubMeshObj()
867 {
868 }
869
870 //=================================================================================
871 // function : GetElemDimension
872 // purpose  : Get dimension of element
873 //=================================================================================
874 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
875 {
876   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
877 }
878
879 //=================================================================================
880 // function : UpdateFunctor
881 // purpose  : Update functor in accordance with current mesh
882 //=================================================================================
883
884 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
885 {
886   theFunctor->SetMesh( myMeshObj->GetMesh() );
887 }
888
889 //=================================================================================
890 // function : Update
891 // purpose  : Update mesh object and fill grid with new values 
892 //=================================================================================
893 bool SMESH_SubMeshObj::Update( int theIsClear )
894 {
895   MESSAGE("SMESH_SubMeshObj::Update " << this)
896   bool changed = myMeshObj->Update( theIsClear );
897   buildPrs(true);
898   return changed;
899 }
900
901
902 /*
903   Class       : SMESH_GroupObj
904   Description : Class for visualisation of groups
905 */
906
907 //=================================================================================
908 // function : SMESH_GroupObj
909 // purpose  : Constructor
910 //=================================================================================
911 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
912                                 SMESH_MeshObj*             theMeshObj )
913 : SMESH_SubMeshObj( theMeshObj ),
914   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
915 {
916   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
917   myGroupServer->Register();
918 }
919
920 SMESH_GroupObj::~SMESH_GroupObj()
921 {
922   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
923   myGroupServer->UnRegister();
924 }
925
926 //=================================================================================
927 // function : IsNodePrs
928 // purpose  : Return true if node presentation is used
929 //=================================================================================
930 bool SMESH_GroupObj::IsNodePrs() const
931 {
932   return myGroupServer->GetType() == SMESH::NODE;
933 }
934
935 //=================================================================================
936 // function : GetElementType
937 // purpose  : Return type of elements of group
938 //=================================================================================
939 SMDSAbs_ElementType SMESH_GroupObj::GetElementType() const
940 {
941   return SMDSAbs_ElementType(myGroupServer->GetType());
942 }
943
944 //=================================================================================
945 // function : getNodesFromElems
946 // purpose  : Retrieve nodes from elements
947 //=================================================================================
948 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
949                               const SMDS_Mesh*                    theMesh,
950                               std::list<const SMDS_MeshElement*>& theResList )
951 {
952   set<const SMDS_MeshElement*> aNodeSet;
953
954   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
955   {
956     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
957     if ( anElem != 0 )
958     {
959       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
960       while ( anIter->more() )
961       {
962         const SMDS_MeshElement* aNode = anIter->next();
963         if ( aNode != 0 )
964           aNodeSet.insert( aNode );
965       }
966     }
967   }
968
969   set<const SMDS_MeshElement*>::const_iterator anIter;
970   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
971     theResList.push_back( *anIter );
972
973   return theResList.size();    
974 }
975
976 //=================================================================================
977 // function : getPointers
978 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
979 //=================================================================================
980 static int getPointers( const SMDSAbs_ElementType            theRequestType,
981                         SMESH::long_array_var&              theElemIds,
982                         const SMDS_Mesh*                    theMesh,
983                         std::list<const SMDS_MeshElement*>& theResList )
984 {
985   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
986   {
987     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
988       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
989
990     if ( anElem != 0 )
991       theResList.push_back( anElem );
992   }
993
994   return theResList.size();
995 }
996
997
998 //=================================================================================
999 // function : GetEntities
1000 // purpose  : Get entities of specified type. Return number of entities
1001 //=================================================================================
1002 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1003 {
1004   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType) {
1005     return myGroupServer->Size();
1006   }
1007   if ( theType == SMDSAbs_Node ) {
1008     return myGroupServer->GetNumberOfNodes();
1009   }
1010   return 0;
1011 }
1012
1013 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1014 {
1015   theResList.clear();
1016   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1017   
1018   if ( aMesh == 0 )
1019     return 0;
1020
1021   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
1022   if ( aGrpType != theType && theType != SMDSAbs_Node )
1023     return 0;
1024
1025   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
1026   if ( anIds->length() == 0 )
1027     return 0;
1028
1029   if ( aGrpType == theType )
1030     return getPointers( theType, anIds, aMesh, theResList );
1031   else if ( theType == SMDSAbs_Node )
1032     return getNodesFromElems( anIds, aMesh, theResList );
1033   else
1034     return 0;
1035 }
1036
1037
1038
1039 /*
1040   Class       : SMESH_subMeshObj
1041   Description : Class for visualisation of submeshes
1042 */
1043
1044 //=================================================================================
1045 // function : SMESH_subMeshObj
1046 // purpose  : Constructor
1047 //=================================================================================
1048 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
1049                                     SMESH_MeshObj*           theMeshObj )
1050 : SMESH_SubMeshObj( theMeshObj ),
1051   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
1052 {
1053   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
1054   
1055   mySubMeshServer->Register();
1056 }
1057
1058 SMESH_subMeshObj::~SMESH_subMeshObj()
1059 {
1060   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1061   mySubMeshServer->UnRegister();
1062 }
1063
1064 //=================================================================================
1065 // function : GetEntities
1066 // purpose  : Get entities of specified type. Return number of entities
1067 //=================================================================================
1068 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1069 {
1070   switch ( theType )
1071   {
1072     case SMDSAbs_Node:
1073     {
1074       return mySubMeshServer->GetNumberOfNodes( /*all=*/true );
1075     }
1076     break;
1077     case SMDSAbs_Ball:
1078     case SMDSAbs_0DElement:
1079     case SMDSAbs_Edge:
1080     case SMDSAbs_Face:
1081     case SMDSAbs_Volume:
1082     {
1083       SMESH::long_array_var anIds = 
1084         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1085       return anIds->length();
1086     }
1087     default:
1088       return 0;
1089     break;
1090   }
1091 }
1092
1093 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1094 {
1095   theResList.clear();
1096
1097   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1098   if ( aMesh == 0 )
1099     return 0;
1100
1101   bool isNodal = IsNodePrs();
1102
1103   if ( isNodal )
1104   {
1105     if ( theType == SMDSAbs_Node )
1106     {
1107       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1108       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1109     }
1110   }
1111   else
1112   {
1113     if ( theType == SMDSAbs_Node )
1114     {
1115       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1116       return getNodesFromElems( anIds, aMesh, theResList );
1117     }
1118     else
1119     {
1120       SMESH::long_array_var anIds = 
1121         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1122       return getPointers( theType, anIds, aMesh, theResList );
1123     }
1124   }
1125
1126   return 0;
1127 }
1128
1129 //=================================================================================
1130 // function : IsNodePrs
1131 // purpose  : Return true if node presentation is used
1132 //=================================================================================
1133 bool SMESH_subMeshObj::IsNodePrs() const
1134 {
1135   return mySubMeshServer->GetNumberOfElements() == 0;
1136 }