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