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