Salome HOME
PAL16774,PAL16631(SALOME crash after a mesh computation that failed because of lack...
[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.salome-platform.org/ or email : webmaster.salome@opencascade.com
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 "SMDS_PolyhedralVolumeOfNodes.hxx"
33 #include "SMESH_Actor.h"
34 #include "SMESH_ControlsDef.hxx"
35 #include "SalomeApp_Application.h"
36 #include "VTKViewer_ExtractUnstructuredGrid.h"
37 #include "VTKViewer_CellLocationsArray.h"
38
39 #include CORBA_SERVER_HEADER(SMESH_Gen)
40 #include CORBA_SERVER_HEADER(SALOME_Exception)
41
42 #include <vtkCell.h>
43 #include <vtkIdList.h>
44 #include <vtkCellArray.h>
45 #include <vtkUnsignedCharArray.h>
46
47 #include <vtkUnstructuredGrid.h>
48
49 #include <memory>
50 #include <sstream>      
51 #include <stdexcept>
52 #include <set>
53
54 #include "utilities.h"
55
56 using namespace std;
57
58 #ifndef EXCEPTION
59 #define EXCEPTION(TYPE, MSG) {\
60   std::ostringstream aStream;\
61   aStream<<__FILE__<<"["<<__LINE__<<"]::"<<MSG;\
62   throw TYPE(aStream.str());\
63 }
64 #endif
65
66 #ifdef _DEBUG_
67 static int MYDEBUG = 0;
68 static int MYDEBUGWITHFILES = 0;
69 #else
70 static int MYDEBUG = 0;
71 static int MYDEBUGWITHFILES = 0;
72 #endif
73
74
75 /*
76   Class       : SMESH_VisualObjDef
77   Description : Base class for all mesh objects to be visuilised
78 */
79
80 //=================================================================================
81 // function : getCellType
82 // purpose  : Get type of VTK cell
83 //=================================================================================
84 static inline vtkIdType getCellType( const SMDSAbs_ElementType theType,
85                                      const bool thePoly,
86                                      const int theNbNodes )
87 {
88   switch( theType )
89   {
90     case SMDSAbs_Edge: 
91       if( theNbNodes == 2 )         return VTK_LINE;
92       else if ( theNbNodes == 3 )   return VTK_QUADRATIC_EDGE;
93       else return VTK_EMPTY_CELL;
94
95     case SMDSAbs_Face  :
96       if (thePoly && theNbNodes>2 ) return VTK_POLYGON;
97       else if ( theNbNodes == 3 )   return VTK_TRIANGLE;
98       else if ( theNbNodes == 4 )   return VTK_QUAD;
99       else if ( theNbNodes == 6 )   return VTK_QUADRATIC_TRIANGLE;
100       else if ( theNbNodes == 8 )   return VTK_QUADRATIC_QUAD;
101       else return VTK_EMPTY_CELL;
102       
103     case SMDSAbs_Volume:
104       if (thePoly && theNbNodes>3 ) return VTK_CONVEX_POINT_SET;
105       else if ( theNbNodes == 4 )   return VTK_TETRA;
106       else if ( theNbNodes == 5 )   return VTK_PYRAMID;
107       else if ( theNbNodes == 6 )   return VTK_WEDGE;
108       else if ( theNbNodes == 8 )   return VTK_HEXAHEDRON;
109       else if ( theNbNodes == 10 )  {
110         return VTK_QUADRATIC_TETRA;
111       }
112       else if ( theNbNodes == 20 )  {
113         return VTK_QUADRATIC_HEXAHEDRON;
114       }
115       else if ( theNbNodes==13 || theNbNodes==15 )  {
116         return VTK_CONVEX_POINT_SET;
117       }
118       else return VTK_EMPTY_CELL;
119
120     default: return VTK_EMPTY_CELL;
121   }
122 }
123
124 //=================================================================================
125 // functions : SMESH_VisualObjDef
126 // purpose   : Constructor
127 //=================================================================================
128 SMESH_VisualObjDef::SMESH_VisualObjDef()
129 {
130   myGrid = vtkUnstructuredGrid::New();
131 }
132 SMESH_VisualObjDef::~SMESH_VisualObjDef()
133 {
134   if ( MYDEBUG )
135     MESSAGE( "~SMESH_MeshObj - myGrid->GetReferenceCount() = " << myGrid->GetReferenceCount() );
136   myGrid->Delete();
137 }
138
139 //=================================================================================
140 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
141 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
142 //=================================================================================
143 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
144 {
145   return myVTK2SMDSNodes.find(theVTKID) == myVTK2SMDSNodes.end() ? -1 : myVTK2SMDSNodes[theVTKID];
146 }
147
148 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
149 {
150   return mySMDS2VTKNodes.find(theObjID) == mySMDS2VTKNodes.end() ? -1 : mySMDS2VTKNodes[theObjID];
151 }
152
153 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
154 {
155   return myVTK2SMDSElems.find(theVTKID) == myVTK2SMDSElems.end() ? -1 : myVTK2SMDSElems[theVTKID];
156 }
157
158 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
159 {
160   return mySMDS2VTKElems.find(theObjID) == mySMDS2VTKElems.end() ? -1 : mySMDS2VTKElems[theObjID];
161 }
162
163 //=================================================================================
164 // function : SMESH_VisualObjDef::createPoints
165 // purpose  : Create points from nodes
166 //=================================================================================
167 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
168 {
169   if ( thePoints == 0 )
170     return;
171
172   TEntityList aNodes;
173   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
174   thePoints->SetNumberOfPoints( nbNodes );
175   
176   int nbPoints = 0;
177
178   TEntityList::const_iterator anIter;
179   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
180   {
181     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
182     if ( aNode != 0 )
183     {
184       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
185       int anId = aNode->GetID();
186       mySMDS2VTKNodes.insert( TMapOfIds::value_type( anId, nbPoints ) );
187       myVTK2SMDSNodes.insert( TMapOfIds::value_type( nbPoints, anId ) );
188       nbPoints++;
189     }
190   }
191
192   if ( nbPoints != nbNodes )
193     thePoints->SetNumberOfPoints( nbPoints );
194 }
195
196 //=================================================================================
197 // function : buildPrs
198 // purpose  : create VTK cells( fill unstructured grid )
199 //=================================================================================
200 void SMESH_VisualObjDef::buildPrs()
201 {
202   // PAL16631(crash after a mesh computation that failed because of lack of memory):
203   // Catch exceptions upper by stack
204 //   try 
205 //   {
206     mySMDS2VTKNodes.clear();
207     myVTK2SMDSNodes.clear();
208     mySMDS2VTKElems.clear();
209     myVTK2SMDSElems.clear();
210     
211     if ( IsNodePrs() )
212       buildNodePrs();
213     else
214       buildElemPrs();
215 //   }
216 //   catch( const std::exception& exc )
217 //   {
218 //     INFOS("Follow exception was cought:\n\t"<<exc.what());
219 //   }
220 //   catch(...)
221 //   {
222 //     INFOS("Unknown exception was cought !!!");
223 //   }
224   
225   if( MYDEBUG ) MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
226   if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
227 }
228
229 //=================================================================================
230 // function : buildNodePrs
231 // purpose  : create VTK cells for nodes
232 //=================================================================================
233 void SMESH_VisualObjDef::buildNodePrs()
234 {
235   vtkPoints* aPoints = vtkPoints::New();
236   createPoints( aPoints );
237   myGrid->SetPoints( aPoints );
238   aPoints->Delete();
239
240   myGrid->SetCells( 0, 0, 0 );
241 }
242
243 //=================================================================================
244 // function : buildElemPrs
245 // purpose  : Create VTK cells for elements
246 //=================================================================================
247
248 namespace{
249   typedef std::vector<const SMDS_MeshElement*> TConnect;
250
251   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
252                  TConnect& theConnect)
253   {
254     theConnect.clear();
255     for(; theNodesIter->more();)
256       theConnect.push_back(theNodesIter->next());
257     return theConnect.size();
258   }
259   
260   inline 
261   void SetId(vtkIdList *theIdList, 
262              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
263              const TConnect& theConnect, 
264              int thePosition,
265              int theId)
266   {
267     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
268   }
269
270 }
271
272
273 void SMESH_VisualObjDef::buildElemPrs()
274 {
275   // Create points
276   
277   vtkPoints* aPoints = vtkPoints::New();
278   createPoints( aPoints );
279   myGrid->SetPoints( aPoints );
280   aPoints->Delete();
281     
282   if ( MYDEBUG )
283     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
284
285   // Calculate cells size
286
287   static SMDSAbs_ElementType aTypes[ 3 ] = { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
288
289   // get entity data
290   map<SMDSAbs_ElementType,int> nbEnts;
291   map<SMDSAbs_ElementType,TEntityList> anEnts;
292
293   for ( int i = 0; i <= 2; i++ )
294     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
295
296   vtkIdType aCellsSize =  3 * nbEnts[ SMDSAbs_Edge ];
297
298   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
299   {
300     if ( nbEnts[ aTypes[ i ] ] )
301     {
302       const TEntityList& aList = anEnts[ aTypes[ i ] ];
303       TEntityList::const_iterator anIter;
304       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
305         aCellsSize += (*anIter)->NbNodes() + 1;
306     }
307   }
308
309   vtkIdType aNbCells = nbEnts[ SMDSAbs_Edge ] + nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
310   
311   if ( MYDEBUG )
312     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
313
314   // Create cells
315   
316   vtkCellArray* aConnectivity = vtkCellArray::New();
317   aConnectivity->Allocate( aCellsSize, 0 );
318   
319   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
320   aCellTypesArray->SetNumberOfComponents( 1 );
321   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
322   
323   vtkIdList *anIdList = vtkIdList::New();
324   vtkIdType iElem = 0;
325
326   TConnect aConnect;
327   aConnect.reserve(VTK_CELL_SIZE);
328
329   for ( int i = 0; i <= 2; i++ ) // iterate through edges, faces and volumes
330   {
331     if( nbEnts[ aTypes[ i ] ] > 0 )
332     {
333       const SMDSAbs_ElementType& aType = aTypes[ i ];
334       const TEntityList& aList = anEnts[ aType ];
335       TEntityList::const_iterator anIter;
336       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
337       {
338         const SMDS_MeshElement* anElem = *anIter;
339         
340         vtkIdType aNbNodes = anElem->NbNodes();
341         anIdList->SetNumberOfIds( aNbNodes );
342
343         int anId = anElem->GetID();
344
345         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
346         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
347
348         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
349         switch(aType){
350         case SMDSAbs_Volume:{
351           aConnect.clear();
352           std::vector<int> aConnectivities;
353           // Convertions connectivities from SMDS to VTK
354           if (anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
355
356             if ( const SMDS_PolyhedralVolumeOfNodes* ph =
357                  dynamic_cast<const SMDS_PolyhedralVolumeOfNodes*> (anElem))
358             {
359               aNbNodes = GetConnect(ph->uniqueNodesIterator(),aConnect);
360               anIdList->SetNumberOfIds( aNbNodes );
361             }
362             for (int k = 0; k < aNbNodes; k++)
363               aConnectivities.push_back(k);
364
365           } else if (aNbNodes == 4) {
366             static int anIds[] = {0,2,1,3};
367             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
368
369           } else if (aNbNodes == 5) {
370             static int anIds[] = {0,3,2,1,4};
371             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
372
373           } else if (aNbNodes == 6) {
374             static int anIds[] = {0,1,2,3,4,5};
375             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
376
377           }
378           else if (aNbNodes == 8) {
379             static int anIds[] = {0,3,2,1,4,7,6,5};
380             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
381
382           }
383           else if (aNbNodes == 10) {
384             static int anIds[] = {0,2,1,3,6,5,4,7,9,8};
385             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
386           }
387           else if (aNbNodes == 13) {
388             static int anIds[] = {0,3,2,1,4,8,7,6,5,9,12,11,10};
389             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
390           }
391           else if (aNbNodes == 15) {
392             static int anIds[] = {0,2,1,3,5,4,8,7,6,11,10,9,12,14,13};
393             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
394             //for (int k = 0; k < aNbNodes; k++) {
395             //  int nn = aConnectivities[k];
396             //  const SMDS_MeshNode* N = static_cast<const SMDS_MeshNode*> (aConnect[nn]);
397             //  cout<<"k="<<k<<"  N("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
398             //}
399           }
400           else if (aNbNodes == 20) {
401             static int anIds[] = {0,3,2,1,4,7,6,5,11,10,9,8,15,14,13,12,16,19,18,17};
402             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
403           }
404           else {
405           }
406
407           if ( aConnect.empty() )
408             GetConnect(aNodesIter,aConnect);
409
410           if (aConnectivities.size() > 0) {
411             for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
412               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
413           }
414           break;
415         }
416         default:
417           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
418             const SMDS_MeshElement* aNode = aNodesIter->next();
419             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
420           }
421         }
422
423         aConnectivity->InsertNextCell( anIdList );
424         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
425
426         iElem++;
427       }
428     }
429   }
430
431   // Insert cells in grid
432   
433   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
434   aCellLocationsArray->SetNumberOfComponents( 1 );
435   aCellLocationsArray->SetNumberOfTuples( aNbCells );
436   
437   aConnectivity->InitTraversal();
438   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
439     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
440
441   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
442   
443   aCellLocationsArray->Delete();
444   aCellTypesArray->Delete();
445   aConnectivity->Delete();
446   anIdList->Delete();
447 }
448
449 //=================================================================================
450 // function : GetEdgeNodes
451 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
452 //=================================================================================
453 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
454                                        const int theEdgeNum,
455                                        int&      theNodeId1,
456                                        int&      theNodeId2 ) const
457 {
458   const SMDS_Mesh* aMesh = GetMesh();
459   if ( aMesh == 0 )
460     return false;
461     
462   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
463   if ( anElem == 0 )
464     return false;
465     
466   int nbNodes = anElem->NbNodes();
467
468   if ( theEdgeNum < 0 || theEdgeNum > 3 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
469     return false;
470
471   int anIds[ nbNodes ];
472   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
473   int i = 0;
474   while( anIter->more() )
475     anIds[ i++ ] = anIter->next()->GetID();
476
477   if ( theEdgeNum < nbNodes - 1 )
478   {
479     theNodeId1 = anIds[ theEdgeNum ];
480     theNodeId2 = anIds[ theEdgeNum + 1 ];
481   }
482   else
483   {
484     theNodeId1 = anIds[ nbNodes - 1 ];
485     theNodeId2 = anIds[ 0 ];
486   }
487
488   return true;
489 }
490
491 /*
492   Class       : SMESH_MeshObj
493   Description : Class for visualisation of mesh
494 */
495
496 //=================================================================================
497 // function : SMESH_MeshObj
498 // purpose  : Constructor
499 //=================================================================================
500 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
501   myClient(SalomeApp_Application::orb(),theMesh)
502 {
503   if ( MYDEBUG ) 
504     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
505 }
506
507 //=================================================================================
508 // function : ~SMESH_MeshObj
509 // purpose  : Destructor
510 //=================================================================================
511 SMESH_MeshObj::~SMESH_MeshObj()
512 {
513   if ( MYDEBUG ) 
514     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
515 }
516
517 //=================================================================================
518 // function : Update
519 // purpose  : Update mesh and fill grid with new values if necessary 
520 //=================================================================================
521 bool SMESH_MeshObj::Update( int theIsClear )
522 {
523   // Update SMDS_Mesh on client part
524   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
525     buildPrs();  // Fill unstructured grid
526     return true;
527   }
528   return false;
529 }
530
531 //=================================================================================
532 // function : GetElemDimension
533 // purpose  : Get dimension of element
534 //=================================================================================
535 int SMESH_MeshObj::GetElemDimension( const int theObjId )
536 {
537   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
538   if ( anElem == 0 )
539     return 0;
540
541   int aType = anElem->GetType();
542   switch ( aType )
543   {
544     case SMDSAbs_Edge  : return 1;
545     case SMDSAbs_Face  : return 2;
546     case SMDSAbs_Volume: return 3;
547     default            : return 0;
548   }
549 }
550
551 //=================================================================================
552 // function : GetEntities
553 // purpose  : Get entities of specified type. Return number of entities
554 //=================================================================================
555 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
556 {
557   switch ( theType )
558   {
559     case SMDSAbs_Node:
560     {
561       return myClient->NbNodes();
562     }
563     break;
564     case SMDSAbs_Edge:
565     {
566       return myClient->NbEdges();
567     }
568     break;
569     case SMDSAbs_Face:
570     {
571       return myClient->NbFaces();
572     }
573     break;
574     case SMDSAbs_Volume:
575     {
576       return myClient->NbVolumes();
577     }
578     break;
579     default:
580       return 0;
581     break;
582   }
583 }
584
585 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
586 {
587   theObjs.clear();
588
589   switch ( theType )
590   {
591     case SMDSAbs_Node:
592     {
593       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
594       while ( anIter->more() ) theObjs.push_back( anIter->next() );
595     }
596     break;
597     case SMDSAbs_Edge:
598     {
599       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
600       while ( anIter->more() ) theObjs.push_back( anIter->next() );
601     }
602     break;
603     case SMDSAbs_Face:
604     {
605       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
606       while ( anIter->more() ) theObjs.push_back( anIter->next() );
607     }
608     break;
609     case SMDSAbs_Volume:
610     {
611       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
612       while ( anIter->more() ) theObjs.push_back( anIter->next() );
613     }
614     break;
615     default:
616     break;
617   }
618
619   return theObjs.size();
620 }
621
622 //=================================================================================
623 // function : UpdateFunctor
624 // purpose  : Update functor in accordance with current mesh
625 //=================================================================================
626 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
627 {
628   theFunctor->SetMesh( GetMesh() );
629 }
630
631 //=================================================================================
632 // function : IsNodePrs
633 // purpose  : Return true if node presentation is used
634 //=================================================================================
635 bool SMESH_MeshObj::IsNodePrs() const
636 {
637   return myClient->NbEdges() == 0 &&myClient->NbFaces() == 0 && myClient->NbVolumes() == 0 ;
638 }
639
640
641 /*
642   Class       : SMESH_SubMeshObj
643   Description : Base class for visualisation of submeshes and groups
644 */
645
646 //=================================================================================
647 // function : SMESH_SubMeshObj
648 // purpose  : Constructor
649 //=================================================================================
650 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
651 {
652   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
653   
654   myMeshObj = theMeshObj;
655 }
656
657 SMESH_SubMeshObj::~SMESH_SubMeshObj()
658 {
659 }
660
661 //=================================================================================
662 // function : GetElemDimension
663 // purpose  : Get dimension of element
664 //=================================================================================
665 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
666 {
667   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
668 }
669
670 //=================================================================================
671 // function : UpdateFunctor
672 // purpose  : Update functor in accordance with current mesh
673 //=================================================================================
674 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
675 {
676   theFunctor->SetMesh( myMeshObj->GetMesh() );
677 }
678
679 //=================================================================================
680 // function : Update
681 // purpose  : Update mesh object and fill grid with new values 
682 //=================================================================================
683 bool SMESH_SubMeshObj::Update( int theIsClear )
684 {
685   bool changed = myMeshObj->Update( theIsClear );
686   buildPrs();
687   return changed;
688 }
689
690
691 /*
692   Class       : SMESH_GroupObj
693   Description : Class for visualisation of groups
694 */
695
696 //=================================================================================
697 // function : SMESH_GroupObj
698 // purpose  : Constructor
699 //=================================================================================
700 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
701                                 SMESH_MeshObj*             theMeshObj )
702 : SMESH_SubMeshObj( theMeshObj ),
703   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
704 {
705   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
706   myGroupServer->Register();
707 }
708
709 SMESH_GroupObj::~SMESH_GroupObj()
710 {
711   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
712   myGroupServer->Destroy();
713 }
714
715 //=================================================================================
716 // function : IsNodePrs
717 // purpose  : Return true if node presentation is used
718 //=================================================================================
719 bool SMESH_GroupObj::IsNodePrs() const
720 {
721   return myGroupServer->GetType() == SMESH::NODE;
722 }
723
724 //=================================================================================
725 // function : getNodesFromElems
726 // purpose  : Retrieve nodes from elements
727 //=================================================================================
728 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
729                               const SMDS_Mesh*                    theMesh,
730                               std::list<const SMDS_MeshElement*>& theResList )
731 {
732   set<const SMDS_MeshElement*> aNodeSet;
733
734   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
735   {
736     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
737     if ( anElem != 0 )
738     {
739       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
740       while ( anIter->more() )
741       {
742         const SMDS_MeshElement* aNode = anIter->next();
743         if ( aNode != 0 )
744           aNodeSet.insert( aNode );
745       }
746     }
747   }
748
749   set<const SMDS_MeshElement*>::const_iterator anIter;
750   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
751     theResList.push_back( *anIter );
752
753   return theResList.size();    
754 }
755
756 //=================================================================================
757 // function : getPointers
758 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
759 //=================================================================================
760 static int getPointers( const SMDSAbs_ElementType            theRequestType,
761                         SMESH::long_array_var&              theElemIds,
762                         const SMDS_Mesh*                    theMesh,
763                         std::list<const SMDS_MeshElement*>& theResList )
764 {
765   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
766   {
767     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
768       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
769
770     if ( anElem != 0 )
771       theResList.push_back( anElem );
772   }
773
774   return theResList.size();
775 }
776
777
778 //=================================================================================
779 // function : GetEntities
780 // purpose  : Get entities of specified type. Return number of entities
781 //=================================================================================
782 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
783 {
784   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
785     return myGroupServer->Size();
786   }
787   return 0;
788 }
789
790 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
791 {
792   theResList.clear();
793   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
794   
795   if ( myGroupServer->Size() == 0 || aMesh == 0 )
796     return 0;
797
798   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
799   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
800
801   if ( aGrpType == theType )
802     return getPointers( theType, anIds, aMesh, theResList );
803   else if ( theType == SMDSAbs_Node )
804     return getNodesFromElems( anIds, aMesh, theResList );
805   else
806     return 0;
807 }    
808
809
810
811 /*
812   Class       : SMESH_subMeshObj
813   Description : Class for visualisation of submeshes
814 */
815
816 //=================================================================================
817 // function : SMESH_subMeshObj
818 // purpose  : Constructor
819 //=================================================================================
820 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
821                                     SMESH_MeshObj*           theMeshObj )
822 : SMESH_SubMeshObj( theMeshObj ),
823   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
824 {
825   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
826   
827   mySubMeshServer->Register();
828 }
829
830 SMESH_subMeshObj::~SMESH_subMeshObj()
831 {
832   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
833   mySubMeshServer->Destroy();
834 }
835
836 //=================================================================================
837 // function : GetEntities
838 // purpose  : Get entities of specified type. Return number of entities
839 //=================================================================================
840 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
841 {
842   switch ( theType )
843   {
844     case SMDSAbs_Node:
845     {
846       return mySubMeshServer->GetNumberOfNodes( false );
847     }
848     break;
849     case SMDSAbs_Edge:
850     case SMDSAbs_Face:
851     case SMDSAbs_Volume:
852     {
853       SMESH::long_array_var anIds = 
854         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
855       return anIds->length();
856     }
857     default:
858       return 0;
859     break;
860   }
861 }
862
863 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
864 {
865   theResList.clear();
866
867   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
868   if ( aMesh == 0 )
869     return 0;
870
871   bool isNodal = IsNodePrs();
872
873   if ( isNodal )
874   {
875     if ( theType == SMDSAbs_Node )
876     {
877       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
878       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
879     }
880   }
881   else
882   {
883     if ( theType == SMDSAbs_Node )
884     {
885       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
886       return getNodesFromElems( anIds, aMesh, theResList );
887     }
888     else
889     {
890       SMESH::long_array_var anIds = 
891         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
892       return getPointers( theType, anIds, aMesh, theResList );
893     }
894   }
895
896   return 0;
897 }
898
899 //=================================================================================
900 // function : IsNodePrs
901 // purpose  : Return true if node presentation is used
902 //=================================================================================
903 bool SMESH_subMeshObj::IsNodePrs() const
904 {
905   return mySubMeshServer->GetNumberOfElements() == 0;
906 }