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