Salome HOME
fix connectivity, PAL15152 (hexa mesh is built wrong)
[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   try
203   {
204     mySMDS2VTKNodes.clear();
205     myVTK2SMDSNodes.clear();
206     mySMDS2VTKElems.clear();
207     myVTK2SMDSElems.clear();
208     
209     if ( IsNodePrs() )
210       buildNodePrs();
211     else
212       buildElemPrs();
213   }
214   catch( const std::exception& exc )
215   {
216     INFOS("Follow exception was cought:\n\t"<<exc.what());
217   }
218   catch(...)
219   {
220     INFOS("Unknown exception was cought !!!");
221   }
222   
223   if( MYDEBUG ) MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
224   if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
225 }
226
227 //=================================================================================
228 // function : buildNodePrs
229 // purpose  : create VTK cells for nodes
230 //=================================================================================
231 void SMESH_VisualObjDef::buildNodePrs()
232 {
233   vtkPoints* aPoints = vtkPoints::New();
234   createPoints( aPoints );
235   myGrid->SetPoints( aPoints );
236   aPoints->Delete();
237
238   myGrid->SetCells( 0, 0, 0 );
239 }
240
241 //=================================================================================
242 // function : buildElemPrs
243 // purpose  : Create VTK cells for elements
244 //=================================================================================
245
246 namespace{
247   typedef std::vector<const SMDS_MeshElement*> TConnect;
248
249   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
250                  TConnect& theConnect)
251   {
252     theConnect.clear();
253     for(; theNodesIter->more();)
254       theConnect.push_back(theNodesIter->next());
255     return theConnect.size();
256   }
257   
258   inline 
259   void SetId(vtkIdList *theIdList, 
260              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
261              const TConnect& theConnect, 
262              int thePosition,
263              int theId)
264   {
265     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
266   }
267
268 }
269
270
271 void SMESH_VisualObjDef::buildElemPrs()
272 {
273   // Create points
274   
275   vtkPoints* aPoints = vtkPoints::New();
276   createPoints( aPoints );
277   myGrid->SetPoints( aPoints );
278   aPoints->Delete();
279     
280   if ( MYDEBUG )
281     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
282
283   // Calculate cells size
284
285   static SMDSAbs_ElementType aTypes[ 3 ] = { SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
286
287   // get entity data
288   map<SMDSAbs_ElementType,int> nbEnts;
289   map<SMDSAbs_ElementType,TEntityList> anEnts;
290
291   for ( int i = 0; i <= 2; i++ )
292     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
293
294   vtkIdType aCellsSize =  3 * nbEnts[ SMDSAbs_Edge ];
295
296   for ( int i = 1; i <= 2; i++ ) // iterate through faces and volumes
297   {
298     if ( nbEnts[ aTypes[ i ] ] )
299     {
300       const TEntityList& aList = anEnts[ aTypes[ i ] ];
301       TEntityList::const_iterator anIter;
302       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
303         aCellsSize += (*anIter)->NbNodes() + 1;
304     }
305   }
306
307   vtkIdType aNbCells = nbEnts[ SMDSAbs_Edge ] + nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
308   
309   if ( MYDEBUG )
310     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
311
312   // Create cells
313   
314   vtkCellArray* aConnectivity = vtkCellArray::New();
315   aConnectivity->Allocate( aCellsSize, 0 );
316   
317   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
318   aCellTypesArray->SetNumberOfComponents( 1 );
319   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
320   
321   vtkIdList *anIdList = vtkIdList::New();
322   vtkIdType iElem = 0;
323
324   TConnect aConnect;
325   aConnect.reserve(VTK_CELL_SIZE);
326
327   for ( int i = 0; i <= 2; i++ ) // iterate through edges, faces and volumes
328   {
329     if( nbEnts[ aTypes[ i ] ] > 0 )
330     {
331       const SMDSAbs_ElementType& aType = aTypes[ i ];
332       const TEntityList& aList = anEnts[ aType ];
333       TEntityList::const_iterator anIter;
334       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
335       {
336         const SMDS_MeshElement* anElem = *anIter;
337         
338         vtkIdType aNbNodes = anElem->NbNodes();
339         anIdList->SetNumberOfIds( aNbNodes );
340
341         int anId = anElem->GetID();
342
343         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
344         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
345
346         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
347         switch(aType){
348         case SMDSAbs_Volume:{
349           aConnect.clear();
350           std::vector<int> aConnectivities;
351           // Convertions connectivities from SMDS to VTK
352           if (anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
353
354             if ( const SMDS_PolyhedralVolumeOfNodes* ph =
355                  dynamic_cast<const SMDS_PolyhedralVolumeOfNodes*> (anElem))
356             {
357               aNbNodes = GetConnect(ph->uniqueNodesIterator(),aConnect);
358               anIdList->SetNumberOfIds( aNbNodes );
359             }
360             for (int k = 0; k < aNbNodes; k++)
361               aConnectivities.push_back(k);
362
363           } else if (aNbNodes == 4) {
364             static int anIds[] = {0,2,1,3};
365             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
366
367           } else if (aNbNodes == 5) {
368             static int anIds[] = {0,3,2,1,4};
369             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
370
371           } else if (aNbNodes == 6) {
372             static int anIds[] = {0,1,2,3,4,5};
373             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
374
375           }
376           else if (aNbNodes == 8) {
377             static int anIds[] = {0,3,2,1,4,7,6,5};
378             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
379
380           }
381           else if (aNbNodes == 10) {
382             static int anIds[] = {0,2,1,3,6,5,4,7,9,8};
383             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
384           }
385           else if (aNbNodes == 13) {
386             static int anIds[] = {0,3,2,1,4,8,7,6,5,9,12,11,10};
387             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
388           }
389           else if (aNbNodes == 15) {
390             static int anIds[] = {0,2,1,3,5,4,8,7,6,11,10,9,12,14,13};
391             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
392             //for (int k = 0; k < aNbNodes; k++) {
393             //  int nn = aConnectivities[k];
394             //  const SMDS_MeshNode* N = static_cast<const SMDS_MeshNode*> (aConnect[nn]);
395             //  cout<<"k="<<k<<"  N("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
396             //}
397           }
398           else if (aNbNodes == 20) {
399             static int anIds[] = {0,3,2,1,4,7,6,5,11,10,9,8,15,14,13,12,16,19,18,17};
400             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
401           }
402           else {
403           }
404
405           if ( aConnect.empty() )
406             GetConnect(aNodesIter,aConnect);
407
408           if (aConnectivities.size() > 0) {
409             for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
410               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
411           }
412           break;
413         }
414         default:
415           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
416             const SMDS_MeshElement* aNode = aNodesIter->next();
417             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
418           }
419         }
420
421         aConnectivity->InsertNextCell( anIdList );
422         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
423
424         iElem++;
425       }
426     }
427   }
428
429   // Insert cells in grid
430   
431   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
432   aCellLocationsArray->SetNumberOfComponents( 1 );
433   aCellLocationsArray->SetNumberOfTuples( aNbCells );
434   
435   aConnectivity->InitTraversal();
436   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
437     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
438
439   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
440   
441   aCellLocationsArray->Delete();
442   aCellTypesArray->Delete();
443   aConnectivity->Delete();
444   anIdList->Delete();
445 }
446
447 //=================================================================================
448 // function : GetEdgeNodes
449 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
450 //=================================================================================
451 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
452                                        const int theEdgeNum,
453                                        int&      theNodeId1,
454                                        int&      theNodeId2 ) const
455 {
456   const SMDS_Mesh* aMesh = GetMesh();
457   if ( aMesh == 0 )
458     return false;
459     
460   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
461   if ( anElem == 0 )
462     return false;
463     
464   int nbNodes = anElem->NbNodes();
465
466   if ( theEdgeNum < 0 || theEdgeNum > 3 || nbNodes != 3 && nbNodes != 4 || theEdgeNum > nbNodes )
467     return false;
468
469   int anIds[ nbNodes ];
470   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
471   int i = 0;
472   while( anIter->more() )
473     anIds[ i++ ] = anIter->next()->GetID();
474
475   if ( theEdgeNum < nbNodes - 1 )
476   {
477     theNodeId1 = anIds[ theEdgeNum ];
478     theNodeId2 = anIds[ theEdgeNum + 1 ];
479   }
480   else
481   {
482     theNodeId1 = anIds[ nbNodes - 1 ];
483     theNodeId2 = anIds[ 0 ];
484   }
485
486   return true;
487 }
488
489 /*
490   Class       : SMESH_MeshObj
491   Description : Class for visualisation of mesh
492 */
493
494 //=================================================================================
495 // function : SMESH_MeshObj
496 // purpose  : Constructor
497 //=================================================================================
498 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
499   myClient(SalomeApp_Application::orb(),theMesh)
500 {
501   if ( MYDEBUG ) 
502     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
503 }
504
505 //=================================================================================
506 // function : ~SMESH_MeshObj
507 // purpose  : Destructor
508 //=================================================================================
509 SMESH_MeshObj::~SMESH_MeshObj()
510 {
511   if ( MYDEBUG ) 
512     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
513 }
514
515 //=================================================================================
516 // function : Update
517 // purpose  : Update mesh and fill grid with new values if necessary 
518 //=================================================================================
519 void SMESH_MeshObj::Update( int theIsClear )
520 {
521   // Update SMDS_Mesh on client part
522   if ( myClient.Update(theIsClear) )
523     buildPrs();  // Fill unstructured grid
524 }
525
526 //=================================================================================
527 // function : GetElemDimension
528 // purpose  : Get dimension of element
529 //=================================================================================
530 int SMESH_MeshObj::GetElemDimension( const int theObjId )
531 {
532   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
533   if ( anElem == 0 )
534     return 0;
535
536   int aType = anElem->GetType();
537   switch ( aType )
538   {
539     case SMDSAbs_Edge  : return 1;
540     case SMDSAbs_Face  : return 2;
541     case SMDSAbs_Volume: return 3;
542     default            : return 0;
543   }
544 }
545
546 //=================================================================================
547 // function : GetEntities
548 // purpose  : Get entities of specified type. Return number of entities
549 //=================================================================================
550 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
551 {
552   switch ( theType )
553   {
554     case SMDSAbs_Node:
555     {
556       return myClient->NbNodes();
557     }
558     break;
559     case SMDSAbs_Edge:
560     {
561       return myClient->NbEdges();
562     }
563     break;
564     case SMDSAbs_Face:
565     {
566       return myClient->NbFaces();
567     }
568     break;
569     case SMDSAbs_Volume:
570     {
571       return myClient->NbVolumes();
572     }
573     break;
574     default:
575       return 0;
576     break;
577   }
578 }
579
580 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
581 {
582   theObjs.clear();
583
584   switch ( theType )
585   {
586     case SMDSAbs_Node:
587     {
588       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
589       while ( anIter->more() ) theObjs.push_back( anIter->next() );
590     }
591     break;
592     case SMDSAbs_Edge:
593     {
594       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
595       while ( anIter->more() ) theObjs.push_back( anIter->next() );
596     }
597     break;
598     case SMDSAbs_Face:
599     {
600       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
601       while ( anIter->more() ) theObjs.push_back( anIter->next() );
602     }
603     break;
604     case SMDSAbs_Volume:
605     {
606       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
607       while ( anIter->more() ) theObjs.push_back( anIter->next() );
608     }
609     break;
610     default:
611     break;
612   }
613
614   return theObjs.size();
615 }
616
617 //=================================================================================
618 // function : UpdateFunctor
619 // purpose  : Update functor in accordance with current mesh
620 //=================================================================================
621 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
622 {
623   theFunctor->SetMesh( GetMesh() );
624 }
625
626 //=================================================================================
627 // function : IsNodePrs
628 // purpose  : Return true if node presentation is used
629 //=================================================================================
630 bool SMESH_MeshObj::IsNodePrs() const
631 {
632   return myClient->NbEdges() == 0 &&myClient->NbFaces() == 0 && myClient->NbVolumes() == 0 ;
633 }
634
635
636 /*
637   Class       : SMESH_SubMeshObj
638   Description : Base class for visualisation of submeshes and groups
639 */
640
641 //=================================================================================
642 // function : SMESH_SubMeshObj
643 // purpose  : Constructor
644 //=================================================================================
645 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
646 {
647   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
648   
649   myMeshObj = theMeshObj;
650 }
651
652 SMESH_SubMeshObj::~SMESH_SubMeshObj()
653 {
654 }
655
656 //=================================================================================
657 // function : GetElemDimension
658 // purpose  : Get dimension of element
659 //=================================================================================
660 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
661 {
662   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
663 }
664
665 //=================================================================================
666 // function : UpdateFunctor
667 // purpose  : Update functor in accordance with current mesh
668 //=================================================================================
669 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
670 {
671   theFunctor->SetMesh( myMeshObj->GetMesh() );
672 }
673
674 //=================================================================================
675 // function : Update
676 // purpose  : Update mesh object and fill grid with new values 
677 //=================================================================================
678 void SMESH_SubMeshObj::Update( int theIsClear )
679 {
680   myMeshObj->Update( theIsClear );
681   buildPrs();
682 }
683
684
685 /*
686   Class       : SMESH_GroupObj
687   Description : Class for visualisation of groups
688 */
689
690 //=================================================================================
691 // function : SMESH_GroupObj
692 // purpose  : Constructor
693 //=================================================================================
694 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
695                                 SMESH_MeshObj*             theMeshObj )
696 : SMESH_SubMeshObj( theMeshObj ),
697   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
698 {
699   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
700   myGroupServer->Register();
701 }
702
703 SMESH_GroupObj::~SMESH_GroupObj()
704 {
705   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
706   myGroupServer->Destroy();
707 }
708
709 //=================================================================================
710 // function : IsNodePrs
711 // purpose  : Return true if node presentation is used
712 //=================================================================================
713 bool SMESH_GroupObj::IsNodePrs() const
714 {
715   return myGroupServer->GetType() == SMESH::NODE;
716 }
717
718 //=================================================================================
719 // function : getNodesFromElems
720 // purpose  : Retrieve nodes from elements
721 //=================================================================================
722 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
723                               const SMDS_Mesh*                    theMesh,
724                               std::list<const SMDS_MeshElement*>& theResList )
725 {
726   set<const SMDS_MeshElement*> aNodeSet;
727
728   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
729   {
730     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
731     if ( anElem != 0 )
732     {
733       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
734       while ( anIter->more() )
735       {
736         const SMDS_MeshElement* aNode = anIter->next();
737         if ( aNode != 0 )
738           aNodeSet.insert( aNode );
739       }
740     }
741   }
742
743   set<const SMDS_MeshElement*>::const_iterator anIter;
744   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
745     theResList.push_back( *anIter );
746
747   return theResList.size();    
748 }
749
750 //=================================================================================
751 // function : getPointers
752 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
753 //=================================================================================
754 static int getPointers( const SMDSAbs_ElementType            theRequestType,
755                         SMESH::long_array_var&              theElemIds,
756                         const SMDS_Mesh*                    theMesh,
757                         std::list<const SMDS_MeshElement*>& theResList )
758 {
759   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
760   {
761     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
762       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
763
764     if ( anElem != 0 )
765       theResList.push_back( anElem );
766   }
767
768   return theResList.size();
769 }
770
771
772 //=================================================================================
773 // function : GetEntities
774 // purpose  : Get entities of specified type. Return number of entities
775 //=================================================================================
776 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
777 {
778   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
779     return myGroupServer->Size();
780   }
781   return 0;
782 }
783
784 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
785 {
786   theResList.clear();
787   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
788   
789   if ( myGroupServer->Size() == 0 || aMesh == 0 )
790     return 0;
791
792   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
793   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
794
795   if ( aGrpType == theType )
796     return getPointers( theType, anIds, aMesh, theResList );
797   else if ( theType == SMDSAbs_Node )
798     return getNodesFromElems( anIds, aMesh, theResList );
799   else
800     return 0;
801 }    
802
803
804
805 /*
806   Class       : SMESH_subMeshObj
807   Description : Class for visualisation of submeshes
808 */
809
810 //=================================================================================
811 // function : SMESH_subMeshObj
812 // purpose  : Constructor
813 //=================================================================================
814 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
815                                     SMESH_MeshObj*           theMeshObj )
816 : SMESH_SubMeshObj( theMeshObj ),
817   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
818 {
819   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
820   
821   mySubMeshServer->Register();
822 }
823
824 SMESH_subMeshObj::~SMESH_subMeshObj()
825 {
826   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
827   mySubMeshServer->Destroy();
828 }
829
830 //=================================================================================
831 // function : GetEntities
832 // purpose  : Get entities of specified type. Return number of entities
833 //=================================================================================
834 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
835 {
836   switch ( theType )
837   {
838     case SMDSAbs_Node:
839     {
840       return mySubMeshServer->GetNumberOfNodes( false );
841     }
842     break;
843     case SMDSAbs_Edge:
844     case SMDSAbs_Face:
845     case SMDSAbs_Volume:
846     {
847       SMESH::long_array_var anIds = 
848         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
849       return anIds->length();
850     }
851     default:
852       return 0;
853     break;
854   }
855 }
856
857 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
858 {
859   theResList.clear();
860
861   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
862   if ( aMesh == 0 )
863     return 0;
864
865   bool isNodal = IsNodePrs();
866
867   if ( isNodal )
868   {
869     if ( theType == SMDSAbs_Node )
870     {
871       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
872       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
873     }
874   }
875   else
876   {
877     if ( theType == SMDSAbs_Node )
878     {
879       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
880       return getNodesFromElems( anIds, aMesh, theResList );
881     }
882     else
883     {
884       SMESH::long_array_var anIds = 
885         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
886       return getPointers( theType, anIds, aMesh, theResList );
887     }
888   }
889
890   return 0;
891 }
892
893 //=================================================================================
894 // function : IsNodePrs
895 // purpose  : Return true if node presentation is used
896 //=================================================================================
897 bool SMESH_subMeshObj::IsNodePrs() const
898 {
899   return mySubMeshServer->GetNumberOfElements() == 0;
900 }