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