Salome HOME
Update version number: 3.1.0a2
[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 <VTKViewer_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             for (int k = 0; k < aNbNodes; k++) {
639               aConnectivities.push_back(k);
640             }
641
642           } else if (aNbNodes == 4) {
643             static int anIds[] = {0,2,1,3};
644             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
645
646           } else if (aNbNodes == 5) {
647             static int anIds[] = {0,3,2,1,4};
648             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
649
650           } else if (aNbNodes == 6) {
651             static int anIds[] = {0,1,2,3,4,5};
652             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
653
654           } else if (aNbNodes == 8) {
655             static int anIds[] = {0,3,2,1,4,7,6,5};
656             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
657
658           } else {
659           }
660
661           if (aConnectivities.size() > 0) {
662             for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
663               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
664           }
665           break;
666         }
667         default:
668           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
669             const SMDS_MeshElement* aNode = aNodesIter->next();
670             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
671           }
672         }
673
674         aConnectivity->InsertNextCell( anIdList );
675         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
676
677         iElem++;
678       }
679     }
680   }
681
682   // Insert cells in grid
683   
684   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
685   aCellLocationsArray->SetNumberOfComponents( 1 );
686   aCellLocationsArray->SetNumberOfTuples( aNbCells );
687   
688   aConnectivity->InitTraversal();
689   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
690     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
691
692   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
693   
694   aCellLocationsArray->Delete();
695   aCellTypesArray->Delete();
696   aConnectivity->Delete();
697   anIdList->Delete();
698 }
699
700 //=================================================================================
701 // function : GetEdgeNodes
702 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
703 //=================================================================================
704 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
705                                        const int theEdgeNum,
706                                        int&      theNodeId1,
707                                        int&      theNodeId2 ) const
708 {
709   const SMDS_Mesh* aMesh = GetMesh();
710   if ( aMesh == 0 )
711     return false;
712     
713   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
714   if ( anElem == 0 )
715     return false;
716     
717   int nbNodes = anElem->NbNodes();
718
719   if ( theEdgeNum < 0 || theEdgeNum > 3 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
720     return false;
721
722   int anIds[ nbNodes ];
723   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
724   int i = 0;
725   while( anIter->more() )
726     anIds[ i++ ] = anIter->next()->GetID();
727
728   if ( theEdgeNum < nbNodes - 1 )
729   {
730     theNodeId1 = anIds[ theEdgeNum ];
731     theNodeId2 = anIds[ theEdgeNum + 1 ];
732   }
733   else
734   {
735     theNodeId1 = anIds[ nbNodes - 1 ];
736     theNodeId2 = anIds[ 0 ];
737   }
738
739   return true;
740 }
741
742 /*
743   Class       : SMESH_MeshObj
744   Description : Class for visualisation of mesh
745 */
746
747 //=================================================================================
748 // function : SMESH_MeshObj
749 // purpose  : Constructor
750 //=================================================================================
751 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh)
752 {
753   if ( MYDEBUG ) 
754     MESSAGE("SMESH_MeshObj - theMesh->_is_nil() = "<<theMesh->_is_nil());
755     
756   myMeshServer = SMESH::SMESH_Mesh::_duplicate( theMesh );
757   myMeshServer->Register();
758   myMesh = new SMDS_Mesh();
759 }
760
761 //=================================================================================
762 // function : ~SMESH_MeshObj
763 // purpose  : Destructor
764 //=================================================================================
765 SMESH_MeshObj::~SMESH_MeshObj()
766 {
767   myMeshServer->Destroy();
768   delete myMesh;
769 }
770
771 //=================================================================================
772 // function : Update
773 // purpose  : Update mesh and fill grid with new values if necessary 
774 //=================================================================================
775 void SMESH_MeshObj::Update( int theIsClear )
776 {
777   // Update SMDS_Mesh on client part
778   
779   try
780   {
781     SMESH::log_array_var aSeq = myMeshServer->GetLog( theIsClear );
782     CORBA::Long aLength = aSeq->length();
783     
784     if( MYDEBUG ) MESSAGE( "Update: length of the script is "<<aLength );
785     
786     if( !aLength )
787       return;
788
789     for ( CORBA::Long anId = 0; anId < aLength; anId++)
790     {
791       const SMESH::double_array& aCoords = aSeq[anId].coords;
792       const SMESH::long_array& anIndexes = aSeq[anId].indexes;
793       CORBA::Long anElemId = 0, aNbElems = aSeq[anId].number;
794       CORBA::Long aCommand = aSeq[anId].commandType;
795
796       switch(aCommand)
797       {
798         case SMESH::ADD_NODE       : AddNodesWithID      ( myMesh, aSeq, anId ); break;
799         case SMESH::ADD_EDGE       : AddEdgesWithID      ( myMesh, aSeq, anId ); break;
800         case SMESH::ADD_TRIANGLE   : AddTriasWithID      ( myMesh, aSeq, anId ); break;
801         case SMESH::ADD_QUADRANGLE : AddQuadsWithID      ( myMesh, aSeq, anId ); break;
802         case SMESH::ADD_POLYGON    : AddPolygonsWithID   ( myMesh, aSeq, anId ); break;
803         case SMESH::ADD_TETRAHEDRON: AddTetrasWithID     ( myMesh, aSeq, anId ); break;
804         case SMESH::ADD_PYRAMID    : AddPiramidsWithID   ( myMesh, aSeq, anId ); break;
805         case SMESH::ADD_PRISM      : AddPrismsWithID     ( myMesh, aSeq, anId ); break;
806         case SMESH::ADD_HEXAHEDRON : AddHexasWithID      ( myMesh, aSeq, anId ); break;
807         case SMESH::ADD_POLYHEDRON : AddPolyhedronsWithID( myMesh, aSeq, anId ); break;
808
809         case SMESH::REMOVE_NODE:
810           for( ; anElemId < aNbElems; anElemId++ )
811             myMesh->RemoveNode( FindNode( myMesh, anIndexes[anElemId] ) );
812         break;
813         
814         case SMESH::REMOVE_ELEMENT:
815           for( ; anElemId < aNbElems; anElemId++ )
816             myMesh->RemoveElement( FindElement( myMesh, anIndexes[anElemId] ) );
817         break;
818
819         case SMESH::MOVE_NODE:
820           for(CORBA::Long aCoordId=0; anElemId < aNbElems; anElemId++, aCoordId+=3)
821           {
822             SMDS_MeshNode* node =
823               const_cast<SMDS_MeshNode*>( FindNode( myMesh, anIndexes[anElemId] ));
824             node->setXYZ( aCoords[aCoordId], aCoords[aCoordId+1], aCoords[aCoordId+2] );
825           }
826         break;
827
828         case SMESH::CHANGE_ELEMENT_NODES:
829           for ( CORBA::Long i = 0; anElemId < aNbElems; anElemId++ )
830           {
831             // find element
832             const SMDS_MeshElement* elem = FindElement( myMesh, anIndexes[i++] );
833             // nb nodes
834             int nbNodes = anIndexes[i++];
835             // nodes
836             //ASSERT( nbNodes < 9 );
837             const SMDS_MeshNode* aNodes[ nbNodes ];
838             for ( int iNode = 0; iNode < nbNodes; iNode++ )
839               aNodes[ iNode ] = FindNode( myMesh, anIndexes[i++] );
840             // change
841             myMesh->ChangeElementNodes( elem, aNodes, nbNodes );
842           }
843           break;
844
845         case SMESH::CHANGE_POLYHEDRON_NODES:
846           ChangePolyhedronNodes(myMesh, aSeq, anId);
847           break;
848         case SMESH::RENUMBER:
849           for(CORBA::Long i=0; anElemId < aNbElems; anElemId++, i+=3)
850           {
851             myMesh->Renumber( anIndexes[i], anIndexes[i+1], anIndexes[i+2] );
852           }
853           break;
854           
855         default:;
856       }
857     }
858   }
859   catch ( SALOME::SALOME_Exception& exc )
860   {
861     INFOS("Following exception was cought:\n\t"<<exc.details.text);
862   }
863   catch( const std::exception& exc)
864   {
865     INFOS("Following exception was cought:\n\t"<<exc.what());
866   }
867   catch(...)
868   {
869     INFOS("Unknown exception was cought !!!");
870   }
871   
872   if ( MYDEBUG )
873   {
874     MESSAGE("Update - myMesh->NbNodes() = "<<myMesh->NbNodes());
875     MESSAGE("Update - myMesh->NbEdges() = "<<myMesh->NbEdges());
876     MESSAGE("Update - myMesh->NbFaces() = "<<myMesh->NbFaces());
877     MESSAGE("Update - myMesh->NbVolumes() = "<<myMesh->NbVolumes());
878   }
879
880   // Fill unstructured grid
881   buildPrs();
882 }
883
884 //=================================================================================
885 // function : GetElemDimension
886 // purpose  : Get dimension of element
887 //=================================================================================
888 int SMESH_MeshObj::GetElemDimension( const int theObjId )
889 {
890   const SMDS_MeshElement* anElem = myMesh->FindElement( theObjId );
891   if ( anElem == 0 )
892     return 0;
893
894   int aType = anElem->GetType();
895   switch ( aType )
896   {
897     case SMDSAbs_Edge  : return 1;
898     case SMDSAbs_Face  : return 2;
899     case SMDSAbs_Volume: return 3;
900     default            : return 0;
901   }
902 }
903
904 //=================================================================================
905 // function : GetEntities
906 // purpose  : Get entities of specified type. Return number of entities
907 //=================================================================================
908 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
909 {
910   switch ( theType )
911   {
912     case SMDSAbs_Node:
913     {
914       return myMesh->NbNodes();
915     }
916     break;
917     case SMDSAbs_Edge:
918     {
919       return myMesh->NbEdges();
920     }
921     break;
922     case SMDSAbs_Face:
923     {
924       return myMesh->NbFaces();
925     }
926     break;
927     case SMDSAbs_Volume:
928     {
929       return myMesh->NbVolumes();
930     }
931     break;
932     default:
933       return 0;
934     break;
935   }
936 }
937
938 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
939 {
940   theObjs.clear();
941
942   switch ( theType )
943   {
944     case SMDSAbs_Node:
945     {
946       SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
947       while ( anIter->more() ) theObjs.push_back( anIter->next() );
948     }
949     break;
950     case SMDSAbs_Edge:
951     {
952       SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
953       while ( anIter->more() ) theObjs.push_back( anIter->next() );
954     }
955     break;
956     case SMDSAbs_Face:
957     {
958       SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
959       while ( anIter->more() ) theObjs.push_back( anIter->next() );
960     }
961     break;
962     case SMDSAbs_Volume:
963     {
964       SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
965       while ( anIter->more() ) theObjs.push_back( anIter->next() );
966     }
967     break;
968     default:
969     break;
970   }
971
972   return theObjs.size();
973 }
974
975 //=================================================================================
976 // function : UpdateFunctor
977 // purpose  : Update functor in accordance with current mesh
978 //=================================================================================
979 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
980 {
981   theFunctor->SetMesh( GetMesh() );
982 }
983
984 //=================================================================================
985 // function : IsNodePrs
986 // purpose  : Return true if node presentation is used
987 //=================================================================================
988 bool SMESH_MeshObj::IsNodePrs() const
989 {
990   return myMesh->NbEdges() == 0 &&myMesh->NbFaces() == 0 &&myMesh->NbVolumes() == 0 ;
991 }
992
993
994 /*
995   Class       : SMESH_SubMeshObj
996   Description : Base class for visualisation of submeshes and groups
997 */
998
999 //=================================================================================
1000 // function : SMESH_SubMeshObj
1001 // purpose  : Constructor
1002 //=================================================================================
1003 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
1004 {
1005   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
1006   
1007   myMeshObj = theMeshObj;
1008 }
1009
1010 SMESH_SubMeshObj::~SMESH_SubMeshObj()
1011 {
1012 }
1013
1014 //=================================================================================
1015 // function : GetElemDimension
1016 // purpose  : Get dimension of element
1017 //=================================================================================
1018 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
1019 {
1020   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
1021 }
1022
1023 //=================================================================================
1024 // function : UpdateFunctor
1025 // purpose  : Update functor in accordance with current mesh
1026 //=================================================================================
1027 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
1028 {
1029   theFunctor->SetMesh( myMeshObj->GetMesh() );
1030 }
1031
1032 //=================================================================================
1033 // function : Update
1034 // purpose  : Update mesh object and fill grid with new values 
1035 //=================================================================================
1036 void SMESH_SubMeshObj::Update( int theIsClear )
1037 {
1038   myMeshObj->Update( theIsClear );
1039   buildPrs();
1040 }
1041
1042
1043 /*
1044   Class       : SMESH_GroupObj
1045   Description : Class for visualisation of groups
1046 */
1047
1048 //=================================================================================
1049 // function : SMESH_GroupObj
1050 // purpose  : Constructor
1051 //=================================================================================
1052 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
1053                                 SMESH_MeshObj*             theMeshObj )
1054 : SMESH_SubMeshObj( theMeshObj ),
1055   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
1056 {
1057   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
1058   myGroupServer->Register();
1059 }
1060
1061 SMESH_GroupObj::~SMESH_GroupObj()
1062 {
1063   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
1064   myGroupServer->Destroy();
1065 }
1066
1067 //=================================================================================
1068 // function : IsNodePrs
1069 // purpose  : Return true if node presentation is used
1070 //=================================================================================
1071 bool SMESH_GroupObj::IsNodePrs() const
1072 {
1073   return myGroupServer->GetType() == SMESH::NODE;
1074 }
1075
1076 //=================================================================================
1077 // function : getNodesFromElems
1078 // purpose  : Retrieve nodes from elements
1079 //=================================================================================
1080 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
1081                               const SMDS_Mesh*                    theMesh,
1082                               std::list<const SMDS_MeshElement*>& theResList )
1083 {
1084   set<const SMDS_MeshElement*> aNodeSet;
1085
1086   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
1087   {
1088     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
1089     if ( anElem != 0 )
1090     {
1091       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
1092       while ( anIter->more() )
1093       {
1094         const SMDS_MeshElement* aNode = anIter->next();
1095         if ( aNode != 0 )
1096           aNodeSet.insert( aNode );
1097       }
1098     }
1099   }
1100
1101   set<const SMDS_MeshElement*>::const_iterator anIter;
1102   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
1103     theResList.push_back( *anIter );
1104
1105   return theResList.size();    
1106 }
1107
1108 //=================================================================================
1109 // function : getPointers
1110 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
1111 //=================================================================================
1112 static int getPointers( const SMDSAbs_ElementType            theRequestType,
1113                         SMESH::long_array_var&              theElemIds,
1114                         const SMDS_Mesh*                    theMesh,
1115                         std::list<const SMDS_MeshElement*>& theResList )
1116 {
1117   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
1118   {
1119     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
1120       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
1121
1122     if ( anElem != 0 )
1123       theResList.push_back( anElem );
1124   }
1125
1126   return theResList.size();
1127 }
1128
1129
1130 //=================================================================================
1131 // function : GetEntities
1132 // purpose  : Get entities of specified type. Return number of entities
1133 //=================================================================================
1134 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1135 {
1136   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
1137     return myGroupServer->Size();
1138   }
1139   return 0;
1140 }
1141
1142 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1143 {
1144   theResList.clear();
1145   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1146   
1147   if ( myGroupServer->Size() == 0 || aMesh == 0 )
1148     return 0;
1149
1150   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
1151   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
1152
1153   if ( aGrpType == theType )
1154     return getPointers( theType, anIds, aMesh, theResList );
1155   else if ( theType == SMDSAbs_Node )
1156     return getNodesFromElems( anIds, aMesh, theResList );
1157   else
1158     return 0;
1159 }    
1160
1161
1162
1163 /*
1164   Class       : SMESH_subMeshObj
1165   Description : Class for visualisation of submeshes
1166 */
1167
1168 //=================================================================================
1169 // function : SMESH_subMeshObj
1170 // purpose  : Constructor
1171 //=================================================================================
1172 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
1173                                     SMESH_MeshObj*           theMeshObj )
1174 : SMESH_SubMeshObj( theMeshObj ),
1175   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
1176 {
1177   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
1178   
1179   mySubMeshServer->Register();
1180 }
1181
1182 SMESH_subMeshObj::~SMESH_subMeshObj()
1183 {
1184   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1185   mySubMeshServer->Destroy();
1186 }
1187
1188 //=================================================================================
1189 // function : GetEntities
1190 // purpose  : Get entities of specified type. Return number of entities
1191 //=================================================================================
1192 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1193 {
1194   switch ( theType )
1195   {
1196     case SMDSAbs_Node:
1197     {
1198       return mySubMeshServer->GetNumberOfNodes( false );
1199     }
1200     break;
1201     case SMDSAbs_Edge:
1202     case SMDSAbs_Face:
1203     case SMDSAbs_Volume:
1204     {
1205       SMESH::long_array_var anIds = 
1206         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1207       return anIds->length();
1208     }
1209     default:
1210       return 0;
1211     break;
1212   }
1213 }
1214
1215 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1216 {
1217   theResList.clear();
1218
1219   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1220   if ( aMesh == 0 )
1221     return 0;
1222
1223   bool isNodal = IsNodePrs();
1224
1225   if ( isNodal )
1226   {
1227     if ( theType == SMDSAbs_Node )
1228     {
1229       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1230       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1231     }
1232   }
1233   else
1234   {
1235     if ( theType == SMDSAbs_Node )
1236     {
1237       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1238       return getNodesFromElems( anIds, aMesh, theResList );
1239     }
1240     else
1241     {
1242       SMESH::long_array_var anIds = 
1243         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1244       return getPointers( theType, anIds, aMesh, theResList );
1245     }
1246   }
1247
1248   return 0;
1249 }
1250
1251 //=================================================================================
1252 // function : IsNodePrs
1253 // purpose  : Return true if node presentation is used
1254 //=================================================================================
1255 bool SMESH_subMeshObj::IsNodePrs() const
1256 {
1257   return mySubMeshServer->GetNumberOfElements() == 0;
1258 }
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270