Salome HOME
Do not calc Controls when clipping is changed.
[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   myGrid = vtkUnstructuredGrid::New();
291 }
292 SMESH_VisualObjDef::~SMESH_VisualObjDef()
293 {
294   if ( MYDEBUG )
295     MESSAGE( "~SMESH_MeshObj - myGrid->GetReferenceCount() = " << myGrid->GetReferenceCount() );
296   myGrid->Delete();
297 }
298
299 //=================================================================================
300 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
301 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
302 //=================================================================================
303 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
304 {
305   return myVTK2SMDSNodes.find(theVTKID) == myVTK2SMDSNodes.end() ? -1 : myVTK2SMDSNodes[theVTKID];
306 }
307
308 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
309 {
310   return mySMDS2VTKNodes.find(theObjID) == mySMDS2VTKNodes.end() ? -1 : mySMDS2VTKNodes[theObjID];
311 }
312
313 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
314 {
315   return myVTK2SMDSElems.find(theVTKID) == myVTK2SMDSElems.end() ? -1 : myVTK2SMDSElems[theVTKID];
316 }
317
318 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
319 {
320   return mySMDS2VTKElems.find(theObjID) == mySMDS2VTKElems.end() ? -1 : mySMDS2VTKElems[theObjID];
321 }
322
323 //=================================================================================
324 // function : SMESH_VisualObjDef::createPoints
325 // purpose  : Create points from nodes
326 //=================================================================================
327 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
328 {
329   if ( thePoints == 0 )
330     return;
331
332   TEntityList aNodes;
333   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
334   thePoints->SetNumberOfPoints( nbNodes );
335   
336   int nbPoints = 0;
337
338   TEntityList::const_iterator anIter;
339   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
340   {
341     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
342     if ( aNode != 0 )
343     {
344       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
345       int anId = aNode->GetID();
346       mySMDS2VTKNodes.insert( TMapOfIds::value_type( anId, nbPoints ) );
347       myVTK2SMDSNodes.insert( TMapOfIds::value_type( nbPoints, anId ) );
348       nbPoints++;
349     }
350   }
351
352   if ( nbPoints != nbNodes )
353     thePoints->SetNumberOfPoints( nbPoints );
354 }
355
356 //=================================================================================
357 // function : buildPrs
358 // purpose  : create VTK cells( fill unstructured grid )
359 //=================================================================================
360 void SMESH_VisualObjDef::buildPrs()
361 {
362   try
363   {
364     mySMDS2VTKNodes.clear();
365     myVTK2SMDSNodes.clear();
366     mySMDS2VTKElems.clear();
367     myVTK2SMDSElems.clear();
368     
369     if ( IsNodePrs() )
370       buildNodePrs();
371     else
372       buildElemPrs();
373   }
374   catch( const std::exception& exc )
375   {
376     INFOS("Follow exception was cought:\n\t"<<exc.what());
377   }
378   catch(...)
379   {
380     INFOS("Unknown exception was cought !!!");
381   }
382   
383   if( MYDEBUG ) MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
384   if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
385 }
386
387 //=================================================================================
388 // function : buildNodePrs
389 // purpose  : create VTK cells for nodes
390 //=================================================================================
391 void SMESH_VisualObjDef::buildNodePrs()
392 {
393   vtkPoints* aPoints = vtkPoints::New();
394   createPoints( aPoints );
395   myGrid->SetPoints( aPoints );
396   aPoints->Delete();
397
398   myGrid->SetCells( 0, 0, 0 );
399
400   // Create cells
401   /*
402   int nbPoints = aPoints->GetNumberOfPoints();
403   vtkIdList *anIdList = vtkIdList::New();
404   anIdList->SetNumberOfIds( 1 );
405
406   vtkCellArray *aCells = vtkCellArray::New();
407   aCells->Allocate( 2 * nbPoints, 0 );
408
409   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
410   aCellTypesArray->SetNumberOfComponents( 1 );
411   aCellTypesArray->Allocate( nbPoints );
412
413   for( vtkIdType aCellId = 0; aCellId < nbPoints; aCellId++ )
414   {
415     anIdList->SetId( 0, aCellId );
416     aCells->InsertNextCell( anIdList );
417     aCellTypesArray->InsertNextValue( VTK_VERTEX );
418   }
419
420   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
421   aCellLocationsArray->SetNumberOfComponents( 1 );
422   aCellLocationsArray->SetNumberOfTuples( nbPoints );
423
424   aCells->InitTraversal();
425   for( vtkIdType i = 0, *pts, npts; aCells->GetNextCell( npts, pts ); i++ )
426     aCellLocationsArray->SetValue( i, aCells->GetTraversalLocation( npts ) );
427
428   myGrid->SetCells( aCellTypesArray, aCellLocationsArray, aCells );
429
430   aCellLocationsArray->Delete();
431   aCellTypesArray->Delete();
432   aCells->Delete();
433   anIdList->Delete(); 
434   */
435 }
436
437 //=================================================================================
438 // function : buildElemPrs
439 // purpose  : Create VTK cells for elements
440 //=================================================================================
441
442 namespace{
443   typedef std::vector<const SMDS_MeshElement*> TConnect;
444
445   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
446                  TConnect& theConnect)
447   {
448     theConnect.clear();
449     for(; theNodesIter->more();)
450       theConnect.push_back(theNodesIter->next());
451     return theConnect.size();
452   }
453   
454   inline 
455   void SetId(vtkIdList *theIdList, 
456              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
457              const TConnect& theConnect, 
458              int thePosition,
459              int theId)
460   {
461     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
462   }
463
464 }
465
466
467 void SMESH_VisualObjDef::buildElemPrs()
468 {
469   // Create points
470   
471   vtkPoints* aPoints = vtkPoints::New();
472   createPoints( aPoints );
473   myGrid->SetPoints( aPoints );
474   aPoints->Delete();
475     
476   if ( MYDEBUG )
477     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
478
479   // Calculate cells size
480
481   static SMDSAbs_ElementType aTypes[ 3 ] = { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
482
483   // get entity data
484   map<SMDSAbs_ElementType,int> nbEnts;
485   map<SMDSAbs_ElementType,TEntityList> anEnts;
486
487   for ( int i = 0; i <= 2; i++ )
488     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
489
490   vtkIdType aCellsSize =  3 * nbEnts[ SMDSAbs_Edge ];
491
492   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
493   {
494     if ( nbEnts[ aTypes[ i ] ] )
495     {
496       const TEntityList& aList = anEnts[ aTypes[ i ] ];
497       TEntityList::const_iterator anIter;
498       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
499         aCellsSize += (*anIter)->NbNodes() + 1;
500     }
501   }
502
503   vtkIdType aNbCells = nbEnts[ SMDSAbs_Edge ] + nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
504   
505   if ( MYDEBUG )
506     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
507
508   // Create cells
509   
510   vtkCellArray* aConnectivity = vtkCellArray::New();
511   aConnectivity->Allocate( aCellsSize, 0 );
512   
513   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
514   aCellTypesArray->SetNumberOfComponents( 1 );
515   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
516   
517   vtkIdList *anIdList = vtkIdList::New();
518   vtkIdType iElem = 0;
519
520   TConnect aConnect;
521   aConnect.reserve(VTK_CELL_SIZE);
522
523   for ( int i = 0; i <= 2; i++ ) // iterate through edges, faces and volumes
524   {
525     if( nbEnts[ aTypes[ i ] ] > 0 )
526     {
527       const SMDSAbs_ElementType& aType = aTypes[ i ];
528       const TEntityList& aList = anEnts[ aType ];
529       TEntityList::const_iterator anIter;
530       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
531       {
532         const SMDS_MeshElement* anElem = *anIter;
533         
534         vtkIdType aNbNodes = anElem->NbNodes();
535         anIdList->SetNumberOfIds( aNbNodes );
536
537         int anId = anElem->GetID();
538
539         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
540         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
541
542         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
543         switch(aType){
544         case SMDSAbs_Volume:{
545           int* aConnectivities = NULL;
546           GetConnect(aNodesIter,aConnect);
547           // Convertions connectivities from SMDS to VTK
548           switch(aNbNodes){
549           case 4:{
550             static int anIds[] = {0,2,1,3};
551             aConnectivities = anIds;
552             break;
553           }
554           case 5:{
555             static int anIds[] = {0,3,2,1,4};
556             aConnectivities = anIds;
557             break;
558           }
559           case 6:{
560             static int anIds[] = {0,1,2,3,4,5};
561             aConnectivities = anIds;
562             break;
563           }
564           case 8:{
565             static int anIds[] = {0,3,2,1,4,7,6,5};
566             aConnectivities = anIds;
567             break;
568           }}
569
570           if(aConnectivities)
571             for( vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++ )
572               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
573           
574           break;
575         }
576         default:
577           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
578             const SMDS_MeshElement* aNode = aNodesIter->next();
579             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
580           }
581         }
582
583         aConnectivity->InsertNextCell( anIdList );
584         aCellTypesArray->InsertNextValue( getCellType( aType, aNbNodes ) );
585
586         iElem++;
587       }
588     }
589   }
590
591   // Insert cells in grid
592   
593   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
594   aCellLocationsArray->SetNumberOfComponents( 1 );
595   aCellLocationsArray->SetNumberOfTuples( aNbCells );
596   
597   aConnectivity->InitTraversal();
598   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
599     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
600
601   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
602   
603   aCellLocationsArray->Delete();
604   aCellTypesArray->Delete();
605   aConnectivity->Delete();
606   anIdList->Delete();
607 }
608
609 //=================================================================================
610 // function : GetEdgeNodes
611 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
612 //=================================================================================
613 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
614                                        const int theEdgeNum,
615                                        int&      theNodeId1,
616                                        int&      theNodeId2 ) const
617 {
618   const SMDS_Mesh* aMesh = GetMesh();
619   if ( aMesh == 0 )
620     return false;
621     
622   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
623   if ( anElem == 0 )
624     return false;
625     
626   int nbNodes = anElem->NbNodes();
627
628   if ( theEdgeNum < 1 || theEdgeNum > 4 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
629     return false;
630
631   int anIds[ nbNodes ];
632   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
633   int i = 0;
634   while( anIter->more() )
635     anIds[ i++ ] = anIter->next()->GetID();
636
637   if ( nbNodes != theEdgeNum )
638   {
639     theNodeId1 = anIds[ theEdgeNum - 1 ];
640     theNodeId2 = anIds[ theEdgeNum ];
641   }
642   else
643   {
644     theNodeId1 = anIds[ nbNodes - 1 ];
645     theNodeId2 = anIds[ 0 ];
646   }
647
648   return true;
649 }
650
651 /*
652   Class       : SMESH_MeshObj
653   Description : Class for visualisation of mesh
654 */
655
656 //=================================================================================
657 // function : SMESH_MeshObj
658 // purpose  : Constructor
659 //=================================================================================
660 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh)
661 {
662   if ( MYDEBUG ) 
663     MESSAGE("SMESH_MeshObj - theMesh->_is_nil() = "<<theMesh->_is_nil());
664     
665   myMeshServer = SMESH::SMESH_Mesh::_duplicate( theMesh );
666   myMeshServer->Register();
667   myMesh = new SMDS_Mesh();
668 }
669
670 //=================================================================================
671 // function : ~SMESH_MeshObj
672 // purpose  : Destructor
673 //=================================================================================
674 SMESH_MeshObj::~SMESH_MeshObj()
675 {
676   myMeshServer->Destroy();
677   delete myMesh;
678 }
679
680 //=================================================================================
681 // function : Update
682 // purpose  : Update mesh and fill grid with new values if necessary 
683 //=================================================================================
684 void SMESH_MeshObj::Update( int theIsClear )
685 {
686   // Update SMDS_Mesh on client part
687   
688   try
689   {
690     SMESH::log_array_var aSeq = myMeshServer->GetLog( theIsClear );
691     CORBA::Long aLength = aSeq->length();
692     
693     if( MYDEBUG ) MESSAGE( "Update: length of the script is "<<aLength );
694     
695     if( !aLength )
696       return;
697       
698     for ( CORBA::Long anId = 0; anId < aLength; anId++)
699     {
700       const SMESH::double_array& aCoords = aSeq[anId].coords;
701       const SMESH::long_array& anIndexes = aSeq[anId].indexes;
702       CORBA::Long anElemId = 0, aNbElems = aSeq[anId].number;
703       CORBA::Long aCommand = aSeq[anId].commandType;
704       
705       switch(aCommand)
706       {
707         case SMESH::ADD_NODE       : AddNodesWithID   ( myMesh, aSeq, anId ); break;
708         case SMESH::ADD_EDGE       : AddEdgesWithID   ( myMesh, aSeq, anId ); break;
709         case SMESH::ADD_TRIANGLE   : AddTriasWithID   ( myMesh, aSeq, anId ); break;
710         case SMESH::ADD_QUADRANGLE : AddQuadsWithID   ( myMesh, aSeq, anId ); break;
711         case SMESH::ADD_TETRAHEDRON: AddTetrasWithID  ( myMesh, aSeq, anId ); break;
712         case SMESH::ADD_PYRAMID    : AddPiramidsWithID( myMesh, aSeq, anId ); break;
713         case SMESH::ADD_PRISM      : AddPrismsWithID  ( myMesh, aSeq, anId ); break;
714         case SMESH::ADD_HEXAHEDRON : AddHexasWithID   ( myMesh, aSeq, anId ); break;
715         
716         case SMESH::REMOVE_NODE:
717           for( ; anElemId < aNbElems; anElemId++ )
718             myMesh->RemoveNode( FindNode( myMesh, anIndexes[anElemId] ) );
719         break;
720         
721         case SMESH::REMOVE_ELEMENT:
722           for( ; anElemId < aNbElems; anElemId++ )
723             myMesh->RemoveElement( FindElement( myMesh, anIndexes[anElemId] ) );
724         break;
725
726         case SMESH::MOVE_NODE:
727           for(CORBA::Long aCoordId=0; anElemId < aNbElems; anElemId++, aCoordId+=3)
728           {
729             SMDS_MeshNode* node =
730               const_cast<SMDS_MeshNode*>( FindNode( myMesh, anIndexes[anElemId] ));
731             node->setXYZ( aCoords[aCoordId], aCoords[aCoordId+1], aCoords[aCoordId+2] );
732           }
733         break;
734
735         case SMESH::CHANGE_ELEMENT_NODES:
736           for ( CORBA::Long i = 0; anElemId < aNbElems; anElemId++ )
737           {
738             // find element
739             const SMDS_MeshElement* elem = FindElement( myMesh, anIndexes[i++] );
740             // nb nodes
741             int nbNodes = anIndexes[i++];
742             // nodes
743             ASSERT( nbNodes < 9 );
744             const SMDS_MeshNode* aNodes[ 8 ];
745             for ( int iNode = 0; iNode < nbNodes; iNode++ )
746               aNodes[ iNode ] = FindNode( myMesh, anIndexes[i++] );
747             // change
748             myMesh->ChangeElementNodes( elem, aNodes, nbNodes );
749           }
750           break;
751
752         case SMESH::RENUMBER:
753           for(CORBA::Long i=0; anElemId < aNbElems; anElemId++, i+=3)
754           {
755             myMesh->Renumber( anIndexes[i], anIndexes[i+1], anIndexes[i+2] );
756           }
757           break;
758           
759         default:;
760       }
761     }
762   }
763   catch ( SALOME::SALOME_Exception& exc )
764   {
765     INFOS("Following exception was cought:\n\t"<<exc.details.text);
766   }
767   catch( const std::exception& exc)
768   {
769     INFOS("Following exception was cought:\n\t"<<exc.what());
770   }
771   catch(...)
772   {
773     INFOS("Unknown exception was cought !!!");
774   }
775   
776   if ( MYDEBUG )
777   {
778     MESSAGE("Update - myMesh->NbNodes() = "<<myMesh->NbNodes());
779     MESSAGE("Update - myMesh->NbEdges() = "<<myMesh->NbEdges());
780     MESSAGE("Update - myMesh->NbFaces() = "<<myMesh->NbFaces());
781     MESSAGE("Update - myMesh->NbVolumes() = "<<myMesh->NbVolumes());
782   }
783
784   // Fill unstructured grid
785   buildPrs();
786 }
787
788 //=================================================================================
789 // function : GetElemDimension
790 // purpose  : Get dimension of element
791 //=================================================================================
792 int SMESH_MeshObj::GetElemDimension( const int theObjId )
793 {
794   const SMDS_MeshElement* anElem = myMesh->FindElement( theObjId );
795   if ( anElem == 0 )
796     return 0;
797
798   int aType = anElem->GetType();
799   switch ( aType )
800   {
801     case SMDSAbs_Edge  : return 1;
802     case SMDSAbs_Face  : return 2;
803     case SMDSAbs_Volume: return 3;
804     default            : return 0;
805   }
806 }
807
808 //=================================================================================
809 // function : GetEntities
810 // purpose  : Get entities of specified type. Return number of entities
811 //=================================================================================
812 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
813 {
814   switch ( theType )
815   {
816     case SMDSAbs_Node:
817     {
818       return myMesh->NbNodes();
819     }
820     break;
821     case SMDSAbs_Edge:
822     {
823       return myMesh->NbEdges();
824     }
825     break;
826     case SMDSAbs_Face:
827     {
828       return myMesh->NbFaces();
829     }
830     break;
831     case SMDSAbs_Volume:
832     {
833       return myMesh->NbVolumes();
834     }
835     break;
836     default:
837       return 0;
838     break;
839   }
840 }
841
842 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
843 {
844   theObjs.clear();
845
846   switch ( theType )
847   {
848     case SMDSAbs_Node:
849     {
850       SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
851       while ( anIter->more() ) theObjs.push_back( anIter->next() );
852     }
853     break;
854     case SMDSAbs_Edge:
855     {
856       SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
857       while ( anIter->more() ) theObjs.push_back( anIter->next() );
858     }
859     break;
860     case SMDSAbs_Face:
861     {
862       SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
863       while ( anIter->more() ) theObjs.push_back( anIter->next() );
864     }
865     break;
866     case SMDSAbs_Volume:
867     {
868       SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
869       while ( anIter->more() ) theObjs.push_back( anIter->next() );
870     }
871     break;
872     default:
873     break;
874   }
875
876   return theObjs.size();
877 }
878
879 //=================================================================================
880 // function : UpdateFunctor
881 // purpose  : Update functor in accordance with current mesh
882 //=================================================================================
883 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
884 {
885   theFunctor->SetMesh( GetMesh() );
886 }
887
888 //=================================================================================
889 // function : IsNodePrs
890 // purpose  : Return true if node presentation is used
891 //=================================================================================
892 bool SMESH_MeshObj::IsNodePrs() const
893 {
894   return myMesh->NbEdges() == 0 &&myMesh->NbFaces() == 0 &&myMesh->NbVolumes() == 0 ;
895 }
896
897
898 /*
899   Class       : SMESH_SubMeshObj
900   Description : Base class for visualisation of submeshes and groups
901 */
902
903 //=================================================================================
904 // function : SMESH_SubMeshObj
905 // purpose  : Constructor
906 //=================================================================================
907 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
908 {
909   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
910   
911   myMeshObj = theMeshObj;
912 }
913
914 SMESH_SubMeshObj::~SMESH_SubMeshObj()
915 {
916 }
917
918 //=================================================================================
919 // function : GetElemDimension
920 // purpose  : Get dimension of element
921 //=================================================================================
922 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
923 {
924   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
925 }
926
927 //=================================================================================
928 // function : UpdateFunctor
929 // purpose  : Update functor in accordance with current mesh
930 //=================================================================================
931 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
932 {
933   theFunctor->SetMesh( myMeshObj->GetMesh() );
934 }
935
936 //=================================================================================
937 // function : Update
938 // purpose  : Update mesh object and fill grid with new values 
939 //=================================================================================
940 void SMESH_SubMeshObj::Update( int theIsClear )
941 {
942   myMeshObj->Update( theIsClear );
943   buildPrs();
944 }
945
946
947 /*
948   Class       : SMESH_GroupObj
949   Description : Class for visualisation of groups
950 */
951
952 //=================================================================================
953 // function : SMESH_GroupObj
954 // purpose  : Constructor
955 //=================================================================================
956 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
957                                 SMESH_MeshObj*             theMeshObj )
958 : SMESH_SubMeshObj( theMeshObj ),
959   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
960 {
961   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
962   myGroupServer->Register();
963 }
964
965 SMESH_GroupObj::~SMESH_GroupObj()
966 {
967   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
968   myGroupServer->Destroy();
969 }
970
971 //=================================================================================
972 // function : IsNodePrs
973 // purpose  : Return true if node presentation is used
974 //=================================================================================
975 bool SMESH_GroupObj::IsNodePrs() const
976 {
977   return myGroupServer->GetType() == SMESH::NODE;
978 }
979
980 //=================================================================================
981 // function : getNodesFromElems
982 // purpose  : Retrieve nodes from elements
983 //=================================================================================
984 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
985                               const SMDS_Mesh*                    theMesh,
986                               std::list<const SMDS_MeshElement*>& theResList )
987 {
988   set<const SMDS_MeshElement*> aNodeSet;
989
990   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
991   {
992     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
993     if ( anElem != 0 )
994     {
995       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
996       while ( anIter->more() )
997       {
998         const SMDS_MeshElement* aNode = anIter->next();
999         if ( aNode != 0 )
1000           aNodeSet.insert( aNode );
1001       }
1002     }
1003   }
1004
1005   set<const SMDS_MeshElement*>::const_iterator anIter;
1006   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
1007     theResList.push_back( *anIter );
1008
1009   return theResList.size();    
1010 }
1011
1012 //=================================================================================
1013 // function : getPointers
1014 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
1015 //=================================================================================
1016 static int getPointers( const SMDSAbs_ElementType            theRequestType,
1017                         SMESH::long_array_var&              theElemIds,
1018                         const SMDS_Mesh*                    theMesh,
1019                         std::list<const SMDS_MeshElement*>& theResList )
1020 {
1021   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
1022   {
1023     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
1024       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
1025
1026     if ( anElem != 0 )
1027       theResList.push_back( anElem );
1028   }
1029
1030   return theResList.size();
1031 }
1032
1033
1034 //=================================================================================
1035 // function : GetEntities
1036 // purpose  : Get entities of specified type. Return number of entities
1037 //=================================================================================
1038 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1039 {
1040   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
1041     return myGroupServer->Size();
1042   }
1043   return 0;
1044 }
1045
1046 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1047 {
1048   theResList.clear();
1049   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1050   
1051   if ( myGroupServer->Size() == 0 || aMesh == 0 )
1052     return 0;
1053
1054   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
1055   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
1056
1057   if ( aGrpType == theType )
1058     return getPointers( theType, anIds, aMesh, theResList );
1059   else if ( theType == SMDSAbs_Node )
1060     return getNodesFromElems( anIds, aMesh, theResList );
1061   else
1062     return 0;
1063 }    
1064
1065
1066
1067 /*
1068   Class       : SMESH_subMeshObj
1069   Description : Class for visualisation of submeshes
1070 */
1071
1072 //=================================================================================
1073 // function : SMESH_subMeshObj
1074 // purpose  : Constructor
1075 //=================================================================================
1076 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
1077                                     SMESH_MeshObj*           theMeshObj )
1078 : SMESH_SubMeshObj( theMeshObj ),
1079   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
1080 {
1081   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
1082   
1083   mySubMeshServer->Register();
1084 }
1085
1086 SMESH_subMeshObj::~SMESH_subMeshObj()
1087 {
1088   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1089   mySubMeshServer->Destroy();
1090 }
1091
1092 //=================================================================================
1093 // function : GetEntities
1094 // purpose  : Get entities of specified type. Return number of entities
1095 //=================================================================================
1096 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1097 {
1098   switch ( theType )
1099   {
1100     case SMDSAbs_Node:
1101     {
1102       return mySubMeshServer->GetNumberOfNodes( false );
1103     }
1104     break;
1105     case SMDSAbs_Edge:
1106     case SMDSAbs_Face:
1107     case SMDSAbs_Volume:
1108     {
1109       SMESH::long_array_var anIds = 
1110         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1111       return anIds->length();
1112     }
1113     default:
1114       return 0;
1115     break;
1116   }
1117 }
1118
1119 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1120 {
1121   theResList.clear();
1122
1123   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1124   if ( aMesh == 0 )
1125     return 0;
1126
1127   bool isNodal = IsNodePrs();
1128
1129   if ( isNodal )
1130   {
1131     if ( theType == SMDSAbs_Node )
1132     {
1133       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1134       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1135     }
1136   }
1137   else
1138   {
1139     if ( theType == SMDSAbs_Node )
1140     {
1141       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1142       return getNodesFromElems( anIds, aMesh, theResList );
1143     }
1144     else
1145     {
1146       SMESH::long_array_var anIds = 
1147         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1148       return getPointers( theType, anIds, aMesh, theResList );
1149     }
1150   }
1151
1152   return 0;
1153 }
1154
1155 //=================================================================================
1156 // function : IsNodePrs
1157 // purpose  : Return true if node presentation is used
1158 //=================================================================================
1159 bool SMESH_subMeshObj::IsNodePrs() const
1160 {
1161   return mySubMeshServer->GetNumberOfElements() == 0;
1162 }
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174