Salome HOME
Bug PAL7444 - display mesh takes a lot of more memory in 2.1.0 than in 2.0.0.
[modules/smesh.git] / src / OBJECT / SMESH_Object.cxx
1 //  SMESH OBJECT : interactive object for SMESH visualization
2 //
3 //  Copyright (C) 2003  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.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
21 //
22 //
23 //
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 "SMESH_Actor.h"
33 #include "SMESH_ControlsDef.hxx"
34 #include "SALOME_ExtractUnstructuredGrid.h"
35
36 #include CORBA_SERVER_HEADER(SALOME_Exception)
37
38 #include <vtkCell.h>
39 #include <vtkIdList.h>
40 #include <vtkIntArray.h>
41 #include <vtkCellArray.h>
42 #include <vtkUnsignedCharArray.h>
43
44 #include <vtkUnstructuredGrid.h>
45
46 #include <memory>
47 #include <sstream>      
48 #include <stdexcept>
49 #include <set>
50
51 #include "utilities.h"
52
53 using namespace std;
54
55 #ifndef EXCEPTION
56 #define EXCEPTION(TYPE, MSG) {\
57   std::ostringstream aStream;\
58   aStream<<__FILE__<<"["<<__LINE__<<"]::"<<MSG;\
59   throw TYPE(aStream.str());\
60 }
61 #endif
62
63 #ifdef _DEBUG_
64 static int MYDEBUG = 0;
65 static int MYDEBUGWITHFILES = 0;
66 #else
67 static int MYDEBUG = 0;
68 static int MYDEBUGWITHFILES = 0;
69 #endif
70
71
72 namespace{
73
74   inline const SMDS_MeshNode* FindNode(const SMDS_Mesh* theMesh, int theId){
75     if(const SMDS_MeshNode* anElem = theMesh->FindNode(theId)) return anElem;
76     EXCEPTION(runtime_error,"SMDS_Mesh::FindNode - cannot find a SMDS_MeshNode for ID = "<<theId);
77   }
78
79
80   inline const SMDS_MeshElement* FindElement(const SMDS_Mesh* theMesh, int theId){
81     if(const SMDS_MeshElement* anElem = theMesh->FindElement(theId)) return anElem;
82     EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot find a SMDS_MeshElement for ID = "<<theId);
83   }
84
85
86   inline void AddNodesWithID(SMDS_Mesh* theMesh, 
87                              SMESH::log_array_var& theSeq,
88                              CORBA::Long theId)
89   {
90     const SMESH::double_array& aCoords = theSeq[theId].coords;
91     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
92     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
93     if(3*aNbElems != aCoords.length())
94       EXCEPTION(runtime_error,"AddNodesWithID - 3*aNbElems != aCoords.length()");
95     for(CORBA::Long aCoordId = 0; anElemId < aNbElems; anElemId++, aCoordId+=3){
96       SMDS_MeshElement* anElem = theMesh->AddNodeWithID(aCoords[aCoordId],
97                                                         aCoords[aCoordId+1],
98                                                         aCoords[aCoordId+2],
99                                                         anIndexes[anElemId]);
100       if(!anElem)
101         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddNodeWithID for ID = "<<anElemId);
102     }
103   }
104
105
106   inline void AddEdgesWithID(SMDS_Mesh* theMesh, 
107                              SMESH::log_array_var& theSeq,
108                              CORBA::Long theId)
109   {
110     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
111     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
112     if(3*aNbElems != anIndexes.length())
113       EXCEPTION(runtime_error,"AddEdgeWithID - 3*aNbElems != aCoords.length()");
114     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=3){
115       SMDS_MeshElement* anElem = theMesh->AddEdgeWithID(anIndexes[anIndexId+1],
116                                                         anIndexes[anIndexId+2],
117                                                         anIndexes[anIndexId]);
118       if(!anElem)
119         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddEdgeWithID for ID = "<<anElemId);
120     }
121   }
122
123
124   inline void AddTriasWithID(SMDS_Mesh* theMesh, 
125                              SMESH::log_array_var& theSeq,
126                              CORBA::Long theId)
127   {
128     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
129     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
130     if(4*aNbElems != anIndexes.length())
131       EXCEPTION(runtime_error,"AddEdgeWithID - 4*aNbElems != anIndexes.length()");
132     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=4){
133       SMDS_MeshElement* anElem = theMesh->AddFaceWithID(anIndexes[anIndexId+1],
134                                                         anIndexes[anIndexId+2],
135                                                         anIndexes[anIndexId+3],
136                                                         anIndexes[anIndexId]);
137       if(!anElem)
138         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddFaceWithID for ID = "<<anElemId);
139     }
140   }
141
142
143   inline void AddQuadsWithID(SMDS_Mesh* theMesh, 
144                              SMESH::log_array_var theSeq,
145                              CORBA::Long theId)
146   {
147     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
148     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
149     if(5*aNbElems != anIndexes.length())
150       EXCEPTION(runtime_error,"AddEdgeWithID - 4*aNbElems != anIndexes.length()");
151     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=5){
152       SMDS_MeshElement* anElem = theMesh->AddFaceWithID(anIndexes[anIndexId+1],
153                                                         anIndexes[anIndexId+2],
154                                                         anIndexes[anIndexId+3],
155                                                         anIndexes[anIndexId+4],
156                                                         anIndexes[anIndexId]);
157       if(!anElem)
158         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddFaceWithID for ID = "<<anElemId);
159     }
160   }
161
162
163   inline void AddTetrasWithID(SMDS_Mesh* theMesh, 
164                               SMESH::log_array_var& theSeq,
165                               CORBA::Long theId)
166   {
167     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
168     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
169     if(5*aNbElems != anIndexes.length())
170       EXCEPTION(runtime_error,"AddEdgeWithID - 5*aNbElems != anIndexes.length()");
171     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=5){
172       SMDS_MeshElement* anElem = theMesh->AddVolumeWithID(anIndexes[anIndexId+1],
173                                                           anIndexes[anIndexId+2],
174                                                           anIndexes[anIndexId+3],
175                                                           anIndexes[anIndexId+4],
176                                                           anIndexes[anIndexId]);
177       if(!anElem)
178         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddFaceWithID for ID = "<<anElemId);
179     }
180   }
181
182
183   inline void AddPiramidsWithID(SMDS_Mesh* theMesh, 
184                                 SMESH::log_array_var& theSeq,
185                                 CORBA::Long theId)
186   {
187     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
188     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
189     if(6*aNbElems != anIndexes.length())
190       EXCEPTION(runtime_error,"AddEdgeWithID - 6*aNbElems != anIndexes.length()");
191     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=6){
192       SMDS_MeshElement* anElem = theMesh->AddVolumeWithID(anIndexes[anIndexId+1],
193                                                           anIndexes[anIndexId+2],
194                                                           anIndexes[anIndexId+3],
195                                                           anIndexes[anIndexId+4],
196                                                           anIndexes[anIndexId+5],
197                                                           anIndexes[anIndexId]);
198       if(!anElem)
199         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddFaceWithID for ID = "<<anElemId);
200     }
201   }
202
203
204   inline void AddPrismsWithID(SMDS_Mesh* theMesh, 
205                               SMESH::log_array_var& theSeq,
206                               CORBA::Long theId)
207   {
208     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
209     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
210     if(7*aNbElems != anIndexes.length())
211       EXCEPTION(runtime_error,"AddEdgeWithID - 7*aNbElems != anIndexes.length()");
212     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=7){
213       SMDS_MeshElement* anElem = theMesh->AddVolumeWithID(anIndexes[anIndexId+1],
214                                                           anIndexes[anIndexId+2],
215                                                           anIndexes[anIndexId+3],
216                                                           anIndexes[anIndexId+4],
217                                                           anIndexes[anIndexId+5],
218                                                           anIndexes[anIndexId+6],
219                                                           anIndexes[anIndexId]);
220       if(!anElem)
221         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddFaceWithID for ID = "<<anElemId);
222     }
223   }
224
225
226   inline void AddHexasWithID(SMDS_Mesh* theMesh, 
227                              SMESH::log_array_var& theSeq,
228                              CORBA::Long theId)
229   {
230     const SMESH::long_array& anIndexes = theSeq[theId].indexes;
231     CORBA::Long anElemId = 0, aNbElems = theSeq[theId].number;
232     if(9*aNbElems != anIndexes.length())
233       EXCEPTION(runtime_error,"AddEdgeWithID - 9*aNbElems != anIndexes.length()");
234     for(CORBA::Long anIndexId = 0; anElemId < aNbElems; anElemId++, anIndexId+=9){
235       SMDS_MeshElement* anElem = theMesh->AddVolumeWithID(anIndexes[anIndexId+1],
236                                                           anIndexes[anIndexId+2],
237                                                           anIndexes[anIndexId+3],
238                                                           anIndexes[anIndexId+4],
239                                                           anIndexes[anIndexId+5],
240                                                           anIndexes[anIndexId+6],
241                                                           anIndexes[anIndexId+7],
242                                                           anIndexes[anIndexId+8],
243                                                           anIndexes[anIndexId]);
244       if(!anElem)
245         EXCEPTION(runtime_error,"SMDS_Mesh::FindElement - cannot AddFaceWithID for ID = "<<anElemId);
246     }
247   }
248
249
250 }
251 /*
252   Class       : SMESH_VisualObjDef
253   Description : Base class for all mesh objects to be visuilised
254 */
255
256 //=================================================================================
257 // function : getCellType
258 // purpose  : Get type of VTK cell
259 //=================================================================================
260 static inline vtkIdType getCellType( const SMDSAbs_ElementType theType,
261                                      const int theNbNodes )
262 {
263   switch( theType )
264   {
265     case SMDSAbs_Edge: 
266       return theNbNodes == 2 ? VTK_LINE : VTK_EMPTY_CELL;
267
268     case SMDSAbs_Face  : 
269       if      ( theNbNodes == 3 ) return VTK_TRIANGLE;
270       else if ( theNbNodes == 4 ) return VTK_QUAD;
271       else                        return VTK_EMPTY_CELL;
272       
273     case SMDSAbs_Volume: 
274       if      ( theNbNodes == 4 ) return VTK_TETRA;
275       else if ( theNbNodes == 5 ) return VTK_PYRAMID;
276       else if ( theNbNodes == 6 ) return VTK_WEDGE;
277       else if ( theNbNodes == 8 ) return VTK_HEXAHEDRON;
278       else                        return VTK_EMPTY_CELL;
279
280     default: return VTK_EMPTY_CELL;
281   }
282 }
283
284 //=================================================================================
285 // functions : SMESH_VisualObjDef
286 // purpose   : Constructor
287 //=================================================================================
288 SMESH_VisualObjDef::SMESH_VisualObjDef()
289 {
290   if(MYDEBUG) MESSAGE("SMESH_MeshObj - "<<this);
291   myGrid = vtkUnstructuredGrid::New();
292 }
293 SMESH_VisualObjDef::~SMESH_VisualObjDef()
294 {
295   if(MYDEBUG) {
296     MESSAGE("~SMESH_MeshObj - "<<this);
297     myGrid->DebugOn();
298   }
299   myGrid->Delete();
300 }
301
302 //=================================================================================
303 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
304 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
305 //=================================================================================
306 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
307 {
308   return myVTK2SMDSNodes.find(theVTKID) == myVTK2SMDSNodes.end() ? -1 : myVTK2SMDSNodes[theVTKID];
309 }
310
311 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
312 {
313   return mySMDS2VTKNodes.find(theObjID) == mySMDS2VTKNodes.end() ? -1 : mySMDS2VTKNodes[theObjID];
314 }
315
316 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
317 {
318   return myVTK2SMDSElems.find(theVTKID) == myVTK2SMDSElems.end() ? -1 : myVTK2SMDSElems[theVTKID];
319 }
320
321 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
322 {
323   return mySMDS2VTKElems.find(theObjID) == mySMDS2VTKElems.end() ? -1 : mySMDS2VTKElems[theObjID];
324 }
325
326 //=================================================================================
327 // function : SMESH_VisualObjDef::createPoints
328 // purpose  : Create points from nodes
329 //=================================================================================
330 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
331 {
332   if ( thePoints == 0 )
333     return;
334
335   TEntityList aNodes;
336   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
337   thePoints->SetNumberOfPoints( nbNodes );
338   
339   int nbPoints = 0;
340
341   TEntityList::const_iterator anIter;
342   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
343   {
344     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
345     if ( aNode != 0 )
346     {
347       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
348       int anId = aNode->GetID();
349       mySMDS2VTKNodes.insert( TMapOfIds::value_type( anId, nbPoints ) );
350       myVTK2SMDSNodes.insert( TMapOfIds::value_type( nbPoints, anId ) );
351       nbPoints++;
352     }
353   }
354
355   if ( nbPoints != nbNodes )
356     thePoints->SetNumberOfPoints( nbPoints );
357 }
358
359 //=================================================================================
360 // function : buildPrs
361 // purpose  : create VTK cells( fill unstructured grid )
362 //=================================================================================
363 void SMESH_VisualObjDef::buildPrs()
364 {
365   try
366   {
367     mySMDS2VTKNodes.clear();
368     myVTK2SMDSNodes.clear();
369     mySMDS2VTKElems.clear();
370     myVTK2SMDSElems.clear();
371     
372     if ( IsNodePrs() )
373       buildNodePrs();
374     else
375       buildElemPrs();
376   }
377   catch( const std::exception& exc )
378   {
379     INFOS("Follow exception was cought:\n\t"<<exc.what());
380   }
381   catch(...)
382   {
383     INFOS("Unknown exception was cought !!!");
384   }
385   
386   if( MYDEBUG ) MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
387   if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
388 }
389
390 //=================================================================================
391 // function : buildNodePrs
392 // purpose  : create VTK cells for nodes
393 //=================================================================================
394 void SMESH_VisualObjDef::buildNodePrs()
395 {
396   vtkPoints* aPoints = vtkPoints::New();
397   createPoints( aPoints );
398   myGrid->SetPoints( aPoints );
399   aPoints->Delete();
400
401   myGrid->SetCells( 0, 0, 0 );
402
403   // Create cells
404   /*
405   int nbPoints = aPoints->GetNumberOfPoints();
406   vtkIdList *anIdList = vtkIdList::New();
407   anIdList->SetNumberOfIds( 1 );
408
409   vtkCellArray *aCells = vtkCellArray::New();
410   aCells->Allocate( 2 * nbPoints, 0 );
411
412   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
413   aCellTypesArray->SetNumberOfComponents( 1 );
414   aCellTypesArray->Allocate( nbPoints );
415
416   for( vtkIdType aCellId = 0; aCellId < nbPoints; aCellId++ )
417   {
418     anIdList->SetId( 0, aCellId );
419     aCells->InsertNextCell( anIdList );
420     aCellTypesArray->InsertNextValue( VTK_VERTEX );
421   }
422
423   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
424   aCellLocationsArray->SetNumberOfComponents( 1 );
425   aCellLocationsArray->SetNumberOfTuples( nbPoints );
426
427   aCells->InitTraversal();
428   for( vtkIdType i = 0, *pts, npts; aCells->GetNextCell( npts, pts ); i++ )
429     aCellLocationsArray->SetValue( i, aCells->GetTraversalLocation( npts ) );
430
431   myGrid->SetCells( aCellTypesArray, aCellLocationsArray, aCells );
432
433   aCellLocationsArray->Delete();
434   aCellTypesArray->Delete();
435   aCells->Delete();
436   anIdList->Delete(); 
437   */
438 }
439
440 //=================================================================================
441 // function : buildElemPrs
442 // purpose  : Create VTK cells for elements
443 //=================================================================================
444
445 namespace{
446   typedef std::vector<const SMDS_MeshElement*> TConnect;
447
448   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
449                  TConnect& theConnect)
450   {
451     theConnect.clear();
452     for(; theNodesIter->more();)
453       theConnect.push_back(theNodesIter->next());
454     return theConnect.size();
455   }
456   
457   inline 
458   void SetId(vtkIdList *theIdList, 
459              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
460              const TConnect& theConnect, 
461              int thePosition,
462              int theId)
463   {
464     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
465   }
466
467 }
468
469
470 void SMESH_VisualObjDef::buildElemPrs()
471 {
472   // Create points
473   
474   vtkPoints* aPoints = vtkPoints::New();
475   createPoints( aPoints );
476   myGrid->SetPoints( aPoints );
477   aPoints->Delete();
478     
479   if ( MYDEBUG )
480     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
481
482   // Calculate cells size
483
484   static SMDSAbs_ElementType aTypes[ 3 ] = { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
485
486   // get entity data
487   map<SMDSAbs_ElementType,int> nbEnts;
488   map<SMDSAbs_ElementType,TEntityList> anEnts;
489
490   for ( int i = 0; i <= 2; i++ )
491     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
492
493   vtkIdType aCellsSize =  3 * nbEnts[ SMDSAbs_Edge ];
494
495   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
496   {
497     if ( nbEnts[ aTypes[ i ] ] )
498     {
499       const TEntityList& aList = anEnts[ aTypes[ i ] ];
500       TEntityList::const_iterator anIter;
501       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
502         aCellsSize += (*anIter)->NbNodes() + 1;
503     }
504   }
505
506   vtkIdType aNbCells = nbEnts[ SMDSAbs_Edge ] + nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
507   
508   if ( MYDEBUG )
509     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
510
511   // Create cells
512   
513   vtkCellArray* aConnectivity = vtkCellArray::New();
514   aConnectivity->Allocate( aCellsSize, 0 );
515   
516   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
517   aCellTypesArray->SetNumberOfComponents( 1 );
518   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
519   
520   vtkIdList *anIdList = vtkIdList::New();
521   vtkIdType iElem = 0;
522
523   TConnect aConnect;
524   aConnect.reserve(VTK_CELL_SIZE);
525
526   for ( int i = 0; i <= 2; i++ ) // iterate through edges, faces and volumes
527   {
528     if( nbEnts[ aTypes[ i ] ] > 0 )
529     {
530       const SMDSAbs_ElementType& aType = aTypes[ i ];
531       const TEntityList& aList = anEnts[ aType ];
532       TEntityList::const_iterator anIter;
533       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
534       {
535         const SMDS_MeshElement* anElem = *anIter;
536         
537         vtkIdType aNbNodes = anElem->NbNodes();
538         anIdList->SetNumberOfIds( aNbNodes );
539
540         int anId = anElem->GetID();
541
542         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
543         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
544
545         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
546         switch(aType){
547         case SMDSAbs_Volume:{
548           int* aConnectivities = NULL;
549           GetConnect(aNodesIter,aConnect);
550           // Convertions connectivities from SMDS to VTK
551           switch(aNbNodes){
552           case 4:{
553             static int anIds[] = {0,2,1,3};
554             aConnectivities = anIds;
555             break;
556           }
557           case 5:{
558             static int anIds[] = {0,3,2,1,4};
559             aConnectivities = anIds;
560             break;
561           }
562           case 6:{
563             static int anIds[] = {0,1,2,3,4,5};
564             aConnectivities = anIds;
565             break;
566           }
567           case 8:{
568             static int anIds[] = {0,3,2,1,4,7,6,5};
569             aConnectivities = anIds;
570             break;
571           }}
572
573           if(aConnectivities)
574             for( vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++ )
575               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
576           
577           break;
578         }
579         default:
580           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
581             const SMDS_MeshElement* aNode = aNodesIter->next();
582             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
583           }
584         }
585
586         aConnectivity->InsertNextCell( anIdList );
587         aCellTypesArray->InsertNextValue( getCellType( aType, aNbNodes ) );
588
589         iElem++;
590       }
591     }
592   }
593
594   // Insert cells in grid
595   
596   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
597   aCellLocationsArray->SetNumberOfComponents( 1 );
598   aCellLocationsArray->SetNumberOfTuples( aNbCells );
599   
600   aConnectivity->InitTraversal();
601   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
602     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
603
604   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
605   
606   aCellLocationsArray->Delete();
607   aCellTypesArray->Delete();
608   aConnectivity->Delete();
609   anIdList->Delete();
610 }
611
612 //=================================================================================
613 // function : GetEdgeNodes
614 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
615 //=================================================================================
616 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
617                                        const int theEdgeNum,
618                                        int&      theNodeId1,
619                                        int&      theNodeId2 ) const
620 {
621   const SMDS_Mesh* aMesh = GetMesh();
622   if ( aMesh == 0 )
623     return false;
624     
625   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
626   if ( anElem == 0 )
627     return false;
628     
629   int nbNodes = anElem->NbNodes();
630
631   if ( theEdgeNum < 1 || theEdgeNum > 4 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
632     return false;
633
634   int anIds[ nbNodes ];
635   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
636   int i = 0;
637   while( anIter->more() )
638     anIds[ i++ ] = anIter->next()->GetID();
639
640   if ( nbNodes != theEdgeNum )
641   {
642     theNodeId1 = anIds[ theEdgeNum - 1 ];
643     theNodeId2 = anIds[ theEdgeNum ];
644   }
645   else
646   {
647     theNodeId1 = anIds[ nbNodes - 1 ];
648     theNodeId2 = anIds[ 0 ];
649   }
650
651   return true;
652 }
653
654 /*
655   Class       : SMESH_MeshObj
656   Description : Class for visualisation of mesh
657 */
658
659 //=================================================================================
660 // function : SMESH_MeshObj
661 // purpose  : Constructor
662 //=================================================================================
663 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh)
664 {
665   if ( MYDEBUG ) 
666     MESSAGE("SMESH_MeshObj - theMesh->_is_nil() = "<<theMesh->_is_nil());
667     
668   myMeshServer = SMESH::SMESH_Mesh::_duplicate( theMesh );
669   myMeshServer->Register();
670   myMesh = new SMDS_Mesh();
671 }
672
673 //=================================================================================
674 // function : ~SMESH_MeshObj
675 // purpose  : Destructor
676 //=================================================================================
677 SMESH_MeshObj::~SMESH_MeshObj()
678 {
679   myMeshServer->Destroy();
680   delete myMesh;
681 }
682
683 //=================================================================================
684 // function : Update
685 // purpose  : Update mesh and fill grid with new values if necessary 
686 //=================================================================================
687 void SMESH_MeshObj::Update( int theIsClear )
688 {
689   // Update SMDS_Mesh on client part
690   
691   try
692   {
693     SMESH::log_array_var aSeq = myMeshServer->GetLog( theIsClear );
694     CORBA::Long aLength = aSeq->length();
695     
696     if( MYDEBUG ) MESSAGE( "Update: length of the script is "<<aLength );
697     
698     if( !aLength )
699       return;
700       
701     for ( CORBA::Long anId = 0; anId < aLength; anId++)
702     {
703       const SMESH::double_array& aCoords = aSeq[anId].coords;
704       const SMESH::long_array& anIndexes = aSeq[anId].indexes;
705       CORBA::Long anElemId = 0, aNbElems = aSeq[anId].number;
706       CORBA::Long aCommand = aSeq[anId].commandType;
707       
708       switch(aCommand)
709       {
710         case SMESH::ADD_NODE       : AddNodesWithID   ( myMesh, aSeq, anId ); break;
711         case SMESH::ADD_EDGE       : AddEdgesWithID   ( myMesh, aSeq, anId ); break;
712         case SMESH::ADD_TRIANGLE   : AddTriasWithID   ( myMesh, aSeq, anId ); break;
713         case SMESH::ADD_QUADRANGLE : AddQuadsWithID   ( myMesh, aSeq, anId ); break;
714         case SMESH::ADD_TETRAHEDRON: AddTetrasWithID  ( myMesh, aSeq, anId ); break;
715         case SMESH::ADD_PYRAMID    : AddPiramidsWithID( myMesh, aSeq, anId ); break;
716         case SMESH::ADD_PRISM      : AddPrismsWithID  ( myMesh, aSeq, anId ); break;
717         case SMESH::ADD_HEXAHEDRON : AddHexasWithID   ( myMesh, aSeq, anId ); break;
718         
719         case SMESH::REMOVE_NODE:
720           for( ; anElemId < aNbElems; anElemId++ )
721             myMesh->RemoveNode( FindNode( myMesh, anIndexes[anElemId] ) );
722         break;
723         
724         case SMESH::REMOVE_ELEMENT:
725           for( ; anElemId < aNbElems; anElemId++ )
726             myMesh->RemoveElement( FindElement( myMesh, anIndexes[anElemId] ) );
727         break;
728
729         case SMESH::MOVE_NODE:
730           for(CORBA::Long aCoordId=0; anElemId < aNbElems; anElemId++, aCoordId+=3)
731           {
732             SMDS_MeshNode* node =
733               const_cast<SMDS_MeshNode*>( FindNode( myMesh, anIndexes[anElemId] ));
734             node->setXYZ( aCoords[aCoordId], aCoords[aCoordId+1], aCoords[aCoordId+2] );
735           }
736         break;
737
738         case SMESH::CHANGE_ELEMENT_NODES:
739           for ( CORBA::Long i = 0; anElemId < aNbElems; anElemId++ )
740           {
741             // find element
742             const SMDS_MeshElement* elem = FindElement( myMesh, anIndexes[i++] );
743             // nb nodes
744             int nbNodes = anIndexes[i++];
745             // nodes
746             ASSERT( nbNodes < 9 );
747             const SMDS_MeshNode* aNodes[ 8 ];
748             for ( int iNode = 0; iNode < nbNodes; iNode++ )
749               aNodes[ iNode ] = FindNode( myMesh, anIndexes[i++] );
750             // change
751             myMesh->ChangeElementNodes( elem, aNodes, nbNodes );
752           }
753           break;
754
755         case SMESH::RENUMBER:
756           for(CORBA::Long i=0; anElemId < aNbElems; anElemId++, i+=3)
757           {
758             myMesh->Renumber( anIndexes[i], anIndexes[i+1], anIndexes[i+2] );
759           }
760           break;
761           
762         default:;
763       }
764     }
765   }
766   catch ( SALOME::SALOME_Exception& exc )
767   {
768     INFOS("Following exception was cought:\n\t"<<exc.details.text);
769   }
770   catch( const std::exception& exc)
771   {
772     INFOS("Following exception was cought:\n\t"<<exc.what());
773   }
774   catch(...)
775   {
776     INFOS("Unknown exception was cought !!!");
777   }
778   
779   if ( MYDEBUG )
780   {
781     MESSAGE("Update - myMesh->NbNodes() = "<<myMesh->NbNodes());
782     MESSAGE("Update - myMesh->NbEdges() = "<<myMesh->NbEdges());
783     MESSAGE("Update - myMesh->NbFaces() = "<<myMesh->NbFaces());
784     MESSAGE("Update - myMesh->NbVolumes() = "<<myMesh->NbVolumes());
785   }
786
787   // Fill unstructured grid
788   buildPrs();
789 }
790
791 //=================================================================================
792 // function : GetElemDimension
793 // purpose  : Get dimension of element
794 //=================================================================================
795 int SMESH_MeshObj::GetElemDimension( const int theObjId )
796 {
797   const SMDS_MeshElement* anElem = myMesh->FindElement( theObjId );
798   if ( anElem == 0 )
799     return 0;
800
801   int aType = anElem->GetType();
802   switch ( aType )
803   {
804     case SMDSAbs_Edge  : return 1;
805     case SMDSAbs_Face  : return 2;
806     case SMDSAbs_Volume: return 3;
807     default            : return 0;
808   }
809 }
810
811 //=================================================================================
812 // function : GetEntities
813 // purpose  : Get entities of specified type. Return number of entities
814 //=================================================================================
815 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
816 {
817   switch ( theType )
818   {
819     case SMDSAbs_Node:
820     {
821       return myMesh->NbNodes();
822     }
823     break;
824     case SMDSAbs_Edge:
825     {
826       return myMesh->NbEdges();
827     }
828     break;
829     case SMDSAbs_Face:
830     {
831       return myMesh->NbFaces();
832     }
833     break;
834     case SMDSAbs_Volume:
835     {
836       return myMesh->NbVolumes();
837     }
838     break;
839     default:
840       return 0;
841     break;
842   }
843 }
844
845 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
846 {
847   theObjs.clear();
848
849   switch ( theType )
850   {
851     case SMDSAbs_Node:
852     {
853       SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
854       while ( anIter->more() ) theObjs.push_back( anIter->next() );
855     }
856     break;
857     case SMDSAbs_Edge:
858     {
859       SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
860       while ( anIter->more() ) theObjs.push_back( anIter->next() );
861     }
862     break;
863     case SMDSAbs_Face:
864     {
865       SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
866       while ( anIter->more() ) theObjs.push_back( anIter->next() );
867     }
868     break;
869     case SMDSAbs_Volume:
870     {
871       SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
872       while ( anIter->more() ) theObjs.push_back( anIter->next() );
873     }
874     break;
875     default:
876     break;
877   }
878
879   return theObjs.size();
880 }
881
882 //=================================================================================
883 // function : UpdateFunctor
884 // purpose  : Update functor in accordance with current mesh
885 //=================================================================================
886 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
887 {
888   theFunctor->SetMesh( GetMesh() );
889 }
890
891 //=================================================================================
892 // function : IsNodePrs
893 // purpose  : Return true if node presentation is used
894 //=================================================================================
895 bool SMESH_MeshObj::IsNodePrs() const
896 {
897   return myMesh->NbEdges() == 0 &&myMesh->NbFaces() == 0 &&myMesh->NbVolumes() == 0 ;
898 }
899
900
901 /*
902   Class       : SMESH_SubMeshObj
903   Description : Base class for visualisation of submeshes and groups
904 */
905
906 //=================================================================================
907 // function : SMESH_SubMeshObj
908 // purpose  : Constructor
909 //=================================================================================
910 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
911 {
912   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
913   
914   myMeshObj = theMeshObj;
915 }
916
917 SMESH_SubMeshObj::~SMESH_SubMeshObj()
918 {
919 }
920
921 //=================================================================================
922 // function : GetElemDimension
923 // purpose  : Get dimension of element
924 //=================================================================================
925 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
926 {
927   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
928 }
929
930 //=================================================================================
931 // function : UpdateFunctor
932 // purpose  : Update functor in accordance with current mesh
933 //=================================================================================
934 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
935 {
936   theFunctor->SetMesh( myMeshObj->GetMesh() );
937 }
938
939 //=================================================================================
940 // function : Update
941 // purpose  : Update mesh object and fill grid with new values 
942 //=================================================================================
943 void SMESH_SubMeshObj::Update( int theIsClear )
944 {
945   myMeshObj->Update( theIsClear );
946   buildPrs();
947 }
948
949
950 /*
951   Class       : SMESH_GroupObj
952   Description : Class for visualisation of groups
953 */
954
955 //=================================================================================
956 // function : SMESH_GroupObj
957 // purpose  : Constructor
958 //=================================================================================
959 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
960                                 SMESH_MeshObj*             theMeshObj )
961 : SMESH_SubMeshObj( theMeshObj ),
962   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
963 {
964   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
965   myGroupServer->Register();
966 }
967
968 SMESH_GroupObj::~SMESH_GroupObj()
969 {
970   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
971   myGroupServer->Destroy();
972 }
973
974 //=================================================================================
975 // function : IsNodePrs
976 // purpose  : Return true if node presentation is used
977 //=================================================================================
978 bool SMESH_GroupObj::IsNodePrs() const
979 {
980   return myGroupServer->GetType() == SMESH::NODE;
981 }
982
983 //=================================================================================
984 // function : getNodesFromElems
985 // purpose  : Retrieve nodes from elements
986 //=================================================================================
987 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
988                               const SMDS_Mesh*                    theMesh,
989                               std::list<const SMDS_MeshElement*>& theResList )
990 {
991   set<const SMDS_MeshElement*> aNodeSet;
992
993   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
994   {
995     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
996     if ( anElem != 0 )
997     {
998       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
999       while ( anIter->more() )
1000       {
1001         const SMDS_MeshElement* aNode = anIter->next();
1002         if ( aNode != 0 )
1003           aNodeSet.insert( aNode );
1004       }
1005     }
1006   }
1007
1008   set<const SMDS_MeshElement*>::const_iterator anIter;
1009   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
1010     theResList.push_back( *anIter );
1011
1012   return theResList.size();    
1013 }
1014
1015 //=================================================================================
1016 // function : getPointers
1017 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
1018 //=================================================================================
1019 static int getPointers( const SMDSAbs_ElementType            theRequestType,
1020                         SMESH::long_array_var&              theElemIds,
1021                         const SMDS_Mesh*                    theMesh,
1022                         std::list<const SMDS_MeshElement*>& theResList )
1023 {
1024   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
1025   {
1026     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
1027       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
1028
1029     if ( anElem != 0 )
1030       theResList.push_back( anElem );
1031   }
1032
1033   return theResList.size();
1034 }
1035
1036
1037 //=================================================================================
1038 // function : GetEntities
1039 // purpose  : Get entities of specified type. Return number of entities
1040 //=================================================================================
1041 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1042 {
1043   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
1044     return myGroupServer->Size();
1045   }
1046   return 0;
1047 }
1048
1049 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1050 {
1051   theResList.clear();
1052   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1053   
1054   if ( myGroupServer->Size() == 0 || aMesh == 0 )
1055     return 0;
1056
1057   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
1058   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
1059
1060   if ( aGrpType == theType )
1061     return getPointers( theType, anIds, aMesh, theResList );
1062   else if ( theType == SMDSAbs_Node )
1063     return getNodesFromElems( anIds, aMesh, theResList );
1064   else
1065     return 0;
1066 }    
1067
1068
1069
1070 /*
1071   Class       : SMESH_subMeshObj
1072   Description : Class for visualisation of submeshes
1073 */
1074
1075 //=================================================================================
1076 // function : SMESH_subMeshObj
1077 // purpose  : Constructor
1078 //=================================================================================
1079 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
1080                                     SMESH_MeshObj*           theMeshObj )
1081 : SMESH_SubMeshObj( theMeshObj ),
1082   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
1083 {
1084   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
1085   
1086   mySubMeshServer->Register();
1087 }
1088
1089 SMESH_subMeshObj::~SMESH_subMeshObj()
1090 {
1091   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1092   mySubMeshServer->Destroy();
1093 }
1094
1095 //=================================================================================
1096 // function : GetEntities
1097 // purpose  : Get entities of specified type. Return number of entities
1098 //=================================================================================
1099 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1100 {
1101   switch ( theType )
1102   {
1103     case SMDSAbs_Node:
1104     {
1105       return mySubMeshServer->GetNumberOfNodes( false );
1106     }
1107     break;
1108     case SMDSAbs_Edge:
1109     case SMDSAbs_Face:
1110     case SMDSAbs_Volume:
1111     {
1112       SMESH::long_array_var anIds = 
1113         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1114       return anIds->length();
1115     }
1116     default:
1117       return 0;
1118     break;
1119   }
1120 }
1121
1122 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1123 {
1124   theResList.clear();
1125
1126   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1127   if ( aMesh == 0 )
1128     return 0;
1129
1130   bool isNodal = IsNodePrs();
1131
1132   if ( isNodal )
1133   {
1134     if ( theType == SMDSAbs_Node )
1135     {
1136       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1137       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1138     }
1139   }
1140   else
1141   {
1142     if ( theType == SMDSAbs_Node )
1143     {
1144       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1145       return getNodesFromElems( anIds, aMesh, theResList );
1146     }
1147     else
1148     {
1149       SMESH::long_array_var anIds = 
1150         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1151       return getPointers( theType, anIds, aMesh, theResList );
1152     }
1153   }
1154
1155   return 0;
1156 }
1157
1158 //=================================================================================
1159 // function : IsNodePrs
1160 // purpose  : Return true if node presentation is used
1161 //=================================================================================
1162 bool SMESH_subMeshObj::IsNodePrs() const
1163 {
1164   return mySubMeshServer->GetNumberOfElements() == 0;
1165 }
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177