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