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