Salome HOME
Replace oe by ?
[modules/smesh.git] / src / OBJECT / SMESH_Object.cxx
1 //  Copyright (C) 2007-2010  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  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 //  SMESH OBJECT : interactive object for SMESH visualization
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 = 1;
68 static int MYDEBUGWITHFILES = 0;//1;
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_0DElement: 
91       return VTK_VERTEX;
92
93     case SMDSAbs_Edge: 
94       if( theNbNodes == 2 )         return VTK_LINE;
95       else if ( theNbNodes == 3 )   return VTK_QUADRATIC_EDGE;
96       else return VTK_EMPTY_CELL;
97
98     case SMDSAbs_Face  :
99       if (thePoly && theNbNodes>2 ) return VTK_POLYGON;
100       else if ( theNbNodes == 3 )   return VTK_TRIANGLE;
101       else if ( theNbNodes == 4 )   return VTK_QUAD;
102       else if ( theNbNodes == 6 )   return VTK_QUADRATIC_TRIANGLE;
103       else if ( theNbNodes == 8 )   return VTK_QUADRATIC_QUAD;
104       else return VTK_EMPTY_CELL;
105       
106     case SMDSAbs_Volume:
107       if (thePoly && theNbNodes>3 ) return VTK_POLYHEDRON; //VTK_CONVEX_POINT_SET;
108       else if ( theNbNodes == 4 )   return VTK_TETRA;
109       else if ( theNbNodes == 5 )   return VTK_PYRAMID;
110       else if ( theNbNodes == 6 )   return VTK_WEDGE;
111       else if ( theNbNodes == 8 )   return VTK_HEXAHEDRON;
112       else if ( theNbNodes == 10 )  {
113         return VTK_QUADRATIC_TETRA;
114       }
115       else if ( theNbNodes == 20 )  {
116         return VTK_QUADRATIC_HEXAHEDRON;
117       }
118       else if ( theNbNodes == 15 )  {
119         return VTK_QUADRATIC_WEDGE;
120       }
121       else if ( theNbNodes==13 )  {
122         return VTK_QUADRATIC_PYRAMID; //VTK_CONVEX_POINT_SET;
123       }
124       else return VTK_EMPTY_CELL;
125
126     default: return VTK_EMPTY_CELL;
127   }
128 }
129
130 //=================================================================================
131 // functions : SMESH_VisualObjDef
132 // purpose   : Constructor
133 //=================================================================================
134 SMESH_VisualObjDef::SMESH_VisualObjDef()
135 {
136   MESSAGE("---------------------------------------------SMESH_VisualObjDef::SMESH_VisualObjDef");
137   myGrid = vtkUnstructuredGrid::New();
138   myLocalGrid = false;
139 }
140 SMESH_VisualObjDef::~SMESH_VisualObjDef()
141 {
142   MESSAGE("---------------------------------------------SMESH_VisualObjDef::~SMESH_VisualObjDef");
143   //if ( MYDEBUG )
144     MESSAGE( "~SMESH_MeshObj - myGrid->GetReferenceCount() = " << myGrid->GetReferenceCount() );
145   myGrid->Delete();
146 }
147
148 //=================================================================================
149 // functions : GetNodeObjId, GetNodeVTKId, GetElemObjId, GetElemVTKId
150 // purpose   : Methods for retrieving VTK IDs by SMDS IDs and  vice versa
151 //=================================================================================
152 vtkIdType SMESH_VisualObjDef::GetNodeObjId( int theVTKID )
153 {
154         if (myLocalGrid)
155         {
156                 TMapOfIds::const_iterator i = myVTK2SMDSNodes.find(theVTKID);
157                 return i == myVTK2SMDSNodes.end() ? -1 : i->second;
158         }
159   return this->GetMesh()->FindNodeVtk(theVTKID)->GetID();
160 }
161
162 vtkIdType SMESH_VisualObjDef::GetNodeVTKId( int theObjID )
163 {
164         if (myLocalGrid)
165         {
166                 TMapOfIds::const_iterator i = mySMDS2VTKNodes.find(theObjID);
167     return i == mySMDS2VTKNodes.end() ? -1 : i->second;
168         }
169
170         const SMDS_MeshNode* aNode = 0;
171         if( this->GetMesh() ) {
172           aNode = this->GetMesh()->FindNode(theObjID);
173         }
174         return aNode ? aNode->getVtkId() : -1;
175 }
176
177 vtkIdType SMESH_VisualObjDef::GetElemObjId( int theVTKID )
178 {
179         if (myLocalGrid)
180         {
181                 TMapOfIds::const_iterator i = myVTK2SMDSElems.find(theVTKID);
182                 return i == myVTK2SMDSElems.end() ? -1 : i->second;
183         }
184   return this->GetMesh()->fromVtkToSmds(theVTKID);
185 }
186
187 vtkIdType SMESH_VisualObjDef::GetElemVTKId( int theObjID )
188 {
189         if (myLocalGrid)
190         {
191                 TMapOfIds::const_iterator i = mySMDS2VTKElems.find(theObjID);
192                 return i == mySMDS2VTKElems.end() ? -1 : i->second;
193         }
194   return this->GetMesh()->FindElement(theObjID)->getVtkId();
195   //return this->GetMesh()->fromSmdsToVtk(theObjID);
196 }
197
198 //=================================================================================
199 // function : SMESH_VisualObjDef::createPoints
200 // purpose  : Create points from nodes
201 //=================================================================================
202 /*! fills a vtkPoints structure for a submesh.
203  *  fills a std::list of SMDS_MeshElements*, then extract the points.
204  *  fills also conversion id maps between SMDS and VTK.
205  */
206 void SMESH_VisualObjDef::createPoints( vtkPoints* thePoints )
207 {
208   if ( thePoints == 0 )
209     return;
210
211   TEntityList aNodes;
212   vtkIdType nbNodes = GetEntities( SMDSAbs_Node, aNodes );
213   thePoints->SetNumberOfPoints( nbNodes );
214
215   int nbPoints = 0;
216
217   TEntityList::const_iterator anIter;
218   for ( anIter = aNodes.begin(); anIter != aNodes.end(); ++anIter )
219   {
220     const SMDS_MeshNode* aNode = ( const SMDS_MeshNode* )(*anIter);
221     if ( aNode != 0 )
222     {
223       thePoints->SetPoint( nbPoints, aNode->X(), aNode->Y(), aNode->Z() );
224       int anId = aNode->GetID();
225       mySMDS2VTKNodes.insert( TMapOfIds::value_type( anId, nbPoints ) );
226       myVTK2SMDSNodes.insert( TMapOfIds::value_type( nbPoints, anId ) );
227       nbPoints++;
228     }
229   }
230
231   if ( nbPoints != nbNodes )
232     thePoints->SetNumberOfPoints( nbPoints );
233 }
234
235 //=================================================================================
236 // function : buildPrs
237 // purpose  : create VTK cells( fill unstructured grid )
238 //=================================================================================
239 void SMESH_VisualObjDef::buildPrs(bool buildGrid)
240 {
241   MESSAGE("----------------------------------------------------------SMESH_VisualObjDef::buildPrs " << buildGrid);
242   if (buildGrid)
243   {
244         myLocalGrid = true;
245         try
246         {
247                 mySMDS2VTKNodes.clear();
248                 myVTK2SMDSNodes.clear();
249                 mySMDS2VTKElems.clear();
250                 myVTK2SMDSElems.clear();
251
252                 if ( IsNodePrs() )
253                         buildNodePrs();
254                 else
255                         buildElemPrs();
256         }
257         catch(...)
258         {
259                 mySMDS2VTKNodes.clear();
260                 myVTK2SMDSNodes.clear();
261                 mySMDS2VTKElems.clear();
262                 myVTK2SMDSElems.clear();
263
264                 myGrid->SetPoints( 0 );
265                 myGrid->SetCells( 0, 0, 0, 0, 0 );
266                 throw;
267         }
268   }
269   else
270   {
271         myLocalGrid = false;
272         if (!GetMesh()->isCompacted())
273           {
274             MESSAGE("*** buildPrs ==> compactMesh!");
275             GetMesh()->compactMesh();
276           }
277         vtkUnstructuredGrid *theGrid = GetMesh()->getGrid();
278         myGrid->ShallowCopy(theGrid);
279         //MESSAGE(myGrid->GetReferenceCount());
280         //MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
281         //MESSAGE( "Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints() );
282         if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"buildPrs.vtu" );
283   }
284 }
285
286 //=================================================================================
287 // function : buildNodePrs
288 // purpose  : create VTK cells for nodes
289 //=================================================================================
290
291 void SMESH_VisualObjDef::buildNodePrs()
292 {
293   // PAL16631: without swap, bad_alloc is not thrown but hung up and crash instead,
294   // so check remaining memory size for safety
295   // SMDS_Mesh::CheckMemory(); // PAL16631
296   vtkPoints* aPoints = vtkPoints::New();
297   createPoints( aPoints );
298   // SMDS_Mesh::CheckMemory();
299   myGrid->SetPoints( aPoints );
300   aPoints->Delete();
301
302   myGrid->SetCells( 0, 0, 0, 0, 0 );
303 }
304
305 //=================================================================================
306 // function : buildElemPrs
307 // purpose  : Create VTK cells for elements
308 //=================================================================================
309
310 namespace{
311   typedef std::vector<const SMDS_MeshElement*> TConnect;
312
313   int GetConnect(const SMDS_ElemIteratorPtr& theNodesIter, 
314                  TConnect& theConnect)
315   {
316     theConnect.clear();
317     for(; theNodesIter->more();)
318       theConnect.push_back(theNodesIter->next());
319     return theConnect.size();
320   }
321   
322   inline 
323   void SetId(vtkIdList *theIdList, 
324              const SMESH_VisualObjDef::TMapOfIds& theSMDS2VTKNodes, 
325              const TConnect& theConnect, 
326              int thePosition,
327              int theId)
328   {
329     theIdList->SetId(thePosition,theSMDS2VTKNodes.find(theConnect[theId]->GetID())->second);
330   }
331
332 }
333
334
335 void SMESH_VisualObjDef::buildElemPrs()
336 {
337   // Create points
338
339   vtkPoints* aPoints = vtkPoints::New();
340   createPoints( aPoints );
341   myGrid->SetPoints( aPoints );
342   aPoints->Delete();
343
344   if ( MYDEBUG )
345     MESSAGE("Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints());
346
347   // Calculate cells size
348
349   static SMDSAbs_ElementType aTypes[ 4 ] =
350     { SMDSAbs_0DElement, SMDSAbs_Edge, SMDSAbs_Face, SMDSAbs_Volume };
351
352   // get entity data
353   map<SMDSAbs_ElementType,int> nbEnts;
354   map<SMDSAbs_ElementType,TEntityList> anEnts;
355
356   for ( int i = 0; i <= 3; i++ )
357     nbEnts[ aTypes[ i ] ] = GetEntities( aTypes[ i ], anEnts[ aTypes[ i ] ] );
358
359   // PAL16631: without swap, bad_alloc is not thrown but hung up and crash instead,
360   // so check remaining memory size for safety
361   // SMDS_Mesh::CheckMemory(); // PAL16631
362
363   vtkIdType aCellsSize =  2 * nbEnts[ SMDSAbs_0DElement ] + 3 * nbEnts[ SMDSAbs_Edge ];
364
365   for ( int i = 2; i <= 3; i++ ) // iterate through faces and volumes
366   {
367     if ( nbEnts[ aTypes[ i ] ] )
368     {
369       const TEntityList& aList = anEnts[ aTypes[ i ] ];
370       TEntityList::const_iterator anIter;
371       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
372         aCellsSize += (*anIter)->NbNodes() + 1;
373     }
374   }
375
376   vtkIdType aNbCells = nbEnts[ SMDSAbs_0DElement ] + nbEnts[ SMDSAbs_Edge ] +
377                        nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
378
379   if ( MYDEBUG )
380     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
381
382   // Create cells
383
384   vtkCellArray* aConnectivity = vtkCellArray::New();
385   aConnectivity->Allocate( aCellsSize, 0 );
386
387   // SMDS_Mesh::CheckMemory(); // PAL16631
388
389   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
390   aCellTypesArray->SetNumberOfComponents( 1 );
391   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
392
393   // SMDS_Mesh::CheckMemory(); // PAL16631
394
395   vtkIdList *anIdList = vtkIdList::New();
396   vtkIdType iElem = 0;
397
398   TConnect aConnect;
399   aConnect.reserve(VTK_CELL_SIZE);
400
401   // SMDS_Mesh::CheckMemory(); // PAL16631
402
403   for ( int i = 0; i <= 3; i++ ) // iterate through 0d elements, edges, faces and volumes
404   {
405     if ( nbEnts[ aTypes[ i ] ] > 0 )
406     {
407       const SMDSAbs_ElementType& aType = aTypes[ i ];
408       const TEntityList& aList = anEnts[ aType ];
409       TEntityList::const_iterator anIter;
410       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
411       {
412         const SMDS_MeshElement* anElem = *anIter;
413
414         vtkIdType aNbNodes = anElem->NbNodes();
415         anIdList->SetNumberOfIds( aNbNodes );
416
417         int anId = anElem->GetID();
418
419         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
420         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
421
422         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
423         switch (aType) {
424         case SMDSAbs_Volume:{
425           aConnect.clear();
426           std::vector<int> aConnectivities;
427           // Convertions connectivities from SMDS to VTK
428           if (anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
429
430             if ( const SMDS_VtkVolume* ph =
431                  dynamic_cast<const SMDS_VtkVolume*> (anElem))
432             {
433               aNbNodes = GetConnect(ph->uniqueNodesIterator(),aConnect);
434               anIdList->SetNumberOfIds( aNbNodes );
435             }
436             for (int k = 0; k < aNbNodes; k++)
437               aConnectivities.push_back(k);
438
439           } else if (aNbNodes == 4) {
440             static int anIds[] = {0,2,1,3};
441             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
442
443           } else if (aNbNodes == 5) {
444             static int anIds[] = {0,3,2,1,4};
445             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
446
447           } else if (aNbNodes == 6) {
448             static int anIds[] = {0,1,2,3,4,5};
449             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
450
451           }
452           else if (aNbNodes == 8) {
453             static int anIds[] = {0,3,2,1,4,7,6,5};
454             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
455
456           }
457           else if (aNbNodes == 10) {
458             static int anIds[] = {0,2,1,3,6,5,4,7,9,8};
459             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
460           }
461           else if (aNbNodes == 13) {
462             static int anIds[] = {0,3,2,1,4,8,7,6,5,9,12,11,10};
463             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
464           }
465           else if (aNbNodes == 15) {
466             //static int anIds[] = {0,2,1,3,5,4,8,7,6,11,10,9,12,14,13};
467             static int anIds[] = {0,1,2,3,4,5,6,7,8,9,10,11,12,13,14};
468             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
469             //for (int k = 0; k < aNbNodes; k++) {
470             //  int nn = aConnectivities[k];
471             //  const SMDS_MeshNode* N = static_cast<const SMDS_MeshNode*> (aConnect[nn]);
472             //  cout<<"k="<<k<<"  N("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
473             //}
474           }
475           else if (aNbNodes == 20) {
476             static int anIds[] = {0,3,2,1,4,7,6,5,11,10,9,8,15,14,13,12,16,19,18,17};
477             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
478           }
479           else {
480           }
481
482           if ( aConnect.empty() )
483             GetConnect(aNodesIter,aConnect);
484
485           if (aConnectivities.size() > 0) {
486             for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
487               SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
488           }
489           break;
490         }
491         default:
492           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
493             const SMDS_MeshElement* aNode = aNodesIter->next();
494             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
495           }
496         }
497
498         aConnectivity->InsertNextCell( anIdList );
499         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
500
501         iElem++;
502       }
503     }
504     // SMDS_Mesh::CheckMemory(); // PAL16631
505   }
506
507   // Insert cells in grid
508
509   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
510   aCellLocationsArray->SetNumberOfComponents( 1 );
511   aCellLocationsArray->SetNumberOfTuples( aNbCells );
512
513   // SMDS_Mesh::CheckMemory(); // PAL16631
514
515   aConnectivity->InitTraversal();
516   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
517     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
518
519   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
520
521   aCellLocationsArray->Delete();
522   aCellTypesArray->Delete();
523   aConnectivity->Delete();
524   anIdList->Delete();
525
526   // SMDS_Mesh::CheckMemory(); // PAL16631
527 }
528
529 //=================================================================================
530 // function : GetEdgeNodes
531 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
532 //=================================================================================
533 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
534                                        const int theEdgeNum,
535                                        int&      theNodeId1,
536                                        int&      theNodeId2 ) const
537 {
538   const SMDS_Mesh* aMesh = GetMesh();
539   if ( aMesh == 0 )
540     return false;
541     
542   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
543   if ( anElem == 0 )
544     return false;
545     
546   int nbNodes = anElem->NbNodes();
547
548   if ( theEdgeNum < 0 || theEdgeNum > 3 || (nbNodes != 3 && nbNodes != 4) || theEdgeNum > nbNodes )
549     return false;
550
551   vector<int> anIds( nbNodes );
552   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
553   int i = 0;
554   while( anIter->more() )
555     anIds[ i++ ] = anIter->next()->GetID();
556
557   if ( theEdgeNum < nbNodes - 1 )
558   {
559     theNodeId1 = anIds[ theEdgeNum ];
560     theNodeId2 = anIds[ theEdgeNum + 1 ];
561   }
562   else
563   {
564     theNodeId1 = anIds[ nbNodes - 1 ];
565     theNodeId2 = anIds[ 0 ];
566   }
567
568   return true;
569 }
570
571 vtkUnstructuredGrid* SMESH_VisualObjDef::GetUnstructuredGrid()
572 {
573         //MESSAGE("SMESH_VisualObjDef::GetUnstructuredGrid " << myGrid);
574         return myGrid;
575 }
576
577
578 //=================================================================================
579 // function : IsValid
580 // purpose  : Return true if there are some entities
581 //=================================================================================
582 bool SMESH_VisualObjDef::IsValid() const
583 {
584         //MESSAGE("SMESH_VisualObjDef::IsValid");
585   return GetNbEntities(SMDSAbs_Node) > 0      || 
586          GetNbEntities(SMDSAbs_0DElement) > 0 || 
587          GetNbEntities(SMDSAbs_Edge) > 0      || 
588          GetNbEntities(SMDSAbs_Face) > 0      ||
589          GetNbEntities(SMDSAbs_Volume) > 0 ;
590 }
591
592 /*
593   Class       : SMESH_MeshObj
594   Description : Class for visualisation of mesh
595 */
596
597 //=================================================================================
598 // function : SMESH_MeshObj
599 // purpose  : Constructor
600 //=================================================================================
601 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
602   myClient(SalomeApp_Application::orb(),theMesh)
603 {
604         myEmptyGrid = 0;
605   if ( MYDEBUG ) 
606     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
607 }
608
609 //=================================================================================
610 // function : ~SMESH_MeshObj
611 // purpose  : Destructor
612 //=================================================================================
613 SMESH_MeshObj::~SMESH_MeshObj()
614 {
615   if ( MYDEBUG ) 
616     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
617 }
618
619 //=================================================================================
620 // function : Update
621 // purpose  : Update mesh and fill grid with new values if necessary 
622 //=================================================================================
623 bool SMESH_MeshObj::Update( int theIsClear )
624 {
625   // Update SMDS_Mesh on client part
626   MESSAGE("SMESH_MeshObj::Update " << this);
627   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
628     MESSAGE("buildPrs");
629     buildPrs();  // Fill unstructured grid
630     return true;
631   }
632   return false;
633 }
634
635 bool SMESH_MeshObj::NulData()
636 {
637         MESSAGE ("SMESH_MeshObj::NulData() ==================================================================================");
638         if (!myEmptyGrid)
639         {
640           myEmptyGrid = SMDS_UnstructuredGrid::New();
641           myEmptyGrid->Initialize();
642           myEmptyGrid->Allocate();
643           vtkPoints* points = vtkPoints::New();
644           points->SetNumberOfPoints(0);
645           myEmptyGrid->SetPoints( points );
646           points->Delete();
647           myEmptyGrid->BuildLinks();
648         }
649         myGrid->ShallowCopy(myEmptyGrid);
650         return true;
651 }
652 //=================================================================================
653 // function : GetElemDimension
654 // purpose  : Get dimension of element
655 //=================================================================================
656 int SMESH_MeshObj::GetElemDimension( const int theObjId )
657 {
658   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
659   if ( anElem == 0 )
660     return 0;
661
662   int aType = anElem->GetType();
663   switch ( aType )
664   {
665     case SMDSAbs_0DElement : return 0;
666     case SMDSAbs_Edge  : return 1;
667     case SMDSAbs_Face  : return 2;
668     case SMDSAbs_Volume: return 3;
669     default            : return 0;
670   }
671 }
672
673 //=================================================================================
674 // function : GetEntities
675 // purpose  : Get entities of specified type. Return number of entities
676 //=================================================================================
677 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
678 {
679   switch ( theType )
680   {
681     case SMDSAbs_Node:
682     {
683       return myClient->NbNodes();
684     }
685     break;
686     case SMDSAbs_0DElement:
687     {
688       return myClient->Nb0DElements();
689     }
690     break;
691     case SMDSAbs_Edge:
692     {
693       return myClient->NbEdges();
694     }
695     break;
696     case SMDSAbs_Face:
697     {
698       return myClient->NbFaces();
699     }
700     break;
701     case SMDSAbs_Volume:
702     {
703       return myClient->NbVolumes();
704     }
705     break;
706     default:
707       return 0;
708     break;
709   }
710 }
711
712 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
713 {
714   theObjs.clear();
715
716   switch ( theType )
717   {
718     case SMDSAbs_Node:
719     {
720       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
721       while ( anIter->more() ) theObjs.push_back( anIter->next() );
722     }
723     break;
724     case SMDSAbs_0DElement:
725     {
726       SMDS_0DElementIteratorPtr anIter = myClient->elements0dIterator();
727       while ( anIter->more() ) theObjs.push_back( anIter->next() );
728     }
729     break;
730     case SMDSAbs_Edge:
731     {
732       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
733       while ( anIter->more() ) theObjs.push_back( anIter->next() );
734     }
735     break;
736     case SMDSAbs_Face:
737     {
738       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
739       while ( anIter->more() ) theObjs.push_back( anIter->next() );
740     }
741     break;
742     case SMDSAbs_Volume:
743     {
744       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
745       while ( anIter->more() ) theObjs.push_back( anIter->next() );
746     }
747     break;
748     default:
749     break;
750   }
751
752   return theObjs.size();
753 }
754
755 //=================================================================================
756 // function : UpdateFunctor
757 // purpose  : Update functor in accordance with current mesh
758 //=================================================================================
759 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
760 {
761   theFunctor->SetMesh( GetMesh() );
762 }
763
764 //=================================================================================
765 // function : IsNodePrs
766 // purpose  : Return true if node presentation is used
767 //=================================================================================
768 bool SMESH_MeshObj::IsNodePrs() const
769 {
770   return myClient->Nb0DElements() == 0 && myClient->NbEdges() == 0 && myClient->NbFaces() == 0 && myClient->NbVolumes() == 0 ;
771 }
772
773
774 /*
775   Class       : SMESH_SubMeshObj
776   Description : Base class for visualisation of submeshes and groups
777 */
778
779 //=================================================================================
780 // function : SMESH_SubMeshObj
781 // purpose  : Constructor
782 //=================================================================================
783 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
784 {
785   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
786   
787   myMeshObj = theMeshObj;
788 }
789
790 SMESH_SubMeshObj::~SMESH_SubMeshObj()
791 {
792 }
793
794 //=================================================================================
795 // function : GetElemDimension
796 // purpose  : Get dimension of element
797 //=================================================================================
798 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
799 {
800   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
801 }
802
803 //=================================================================================
804 // function : UpdateFunctor
805 // purpose  : Update functor in accordance with current mesh
806 //=================================================================================
807 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
808 {
809   theFunctor->SetMesh( myMeshObj->GetMesh() );
810 }
811
812 //=================================================================================
813 // function : Update
814 // purpose  : Update mesh object and fill grid with new values 
815 //=================================================================================
816 bool SMESH_SubMeshObj::Update( int theIsClear )
817 {
818         MESSAGE("SMESH_SubMeshObj::Update " << this)
819   bool changed = myMeshObj->Update( theIsClear );
820   buildPrs(true);
821   return changed;
822 }
823
824
825 /*
826   Class       : SMESH_GroupObj
827   Description : Class for visualisation of groups
828 */
829
830 //=================================================================================
831 // function : SMESH_GroupObj
832 // purpose  : Constructor
833 //=================================================================================
834 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
835                                 SMESH_MeshObj*             theMeshObj )
836 : SMESH_SubMeshObj( theMeshObj ),
837   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
838 {
839   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
840   myGroupServer->Register();
841 }
842
843 SMESH_GroupObj::~SMESH_GroupObj()
844 {
845   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
846   myGroupServer->UnRegister();
847 }
848
849 //=================================================================================
850 // function : IsNodePrs
851 // purpose  : Return true if node presentation is used
852 //=================================================================================
853 bool SMESH_GroupObj::IsNodePrs() const
854 {
855   return myGroupServer->GetType() == SMESH::NODE;
856 }
857
858 //=================================================================================
859 // function : GetElementType
860 // purpose  : Return type of elements of group
861 //=================================================================================
862 SMDSAbs_ElementType SMESH_GroupObj::GetElementType() const
863 {
864   return SMDSAbs_ElementType(myGroupServer->GetType());
865 }
866
867 //=================================================================================
868 // function : getNodesFromElems
869 // purpose  : Retrieve nodes from elements
870 //=================================================================================
871 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
872                               const SMDS_Mesh*                    theMesh,
873                               std::list<const SMDS_MeshElement*>& theResList )
874 {
875   set<const SMDS_MeshElement*> aNodeSet;
876
877   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
878   {
879     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
880     if ( anElem != 0 )
881     {
882       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
883       while ( anIter->more() )
884       {
885         const SMDS_MeshElement* aNode = anIter->next();
886         if ( aNode != 0 )
887           aNodeSet.insert( aNode );
888       }
889     }
890   }
891
892   set<const SMDS_MeshElement*>::const_iterator anIter;
893   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
894     theResList.push_back( *anIter );
895
896   return theResList.size();    
897 }
898
899 //=================================================================================
900 // function : getPointers
901 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
902 //=================================================================================
903 static int getPointers( const SMDSAbs_ElementType            theRequestType,
904                         SMESH::long_array_var&              theElemIds,
905                         const SMDS_Mesh*                    theMesh,
906                         std::list<const SMDS_MeshElement*>& theResList )
907 {
908   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
909   {
910     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
911       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
912
913     if ( anElem != 0 )
914       theResList.push_back( anElem );
915   }
916
917   return theResList.size();
918 }
919
920
921 //=================================================================================
922 // function : GetEntities
923 // purpose  : Get entities of specified type. Return number of entities
924 //=================================================================================
925 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
926 {
927   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
928     return myGroupServer->Size();
929   }
930   return 0;
931 }
932
933 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
934 {
935   theResList.clear();
936   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
937   
938   if ( myGroupServer->Size() == 0 || aMesh == 0 )
939     return 0;
940
941   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
942   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
943
944   if ( aGrpType == theType )
945     return getPointers( theType, anIds, aMesh, theResList );
946   else if ( theType == SMDSAbs_Node )
947     return getNodesFromElems( anIds, aMesh, theResList );
948   else
949     return 0;
950 }    
951
952
953
954 /*
955   Class       : SMESH_subMeshObj
956   Description : Class for visualisation of submeshes
957 */
958
959 //=================================================================================
960 // function : SMESH_subMeshObj
961 // purpose  : Constructor
962 //=================================================================================
963 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
964                                     SMESH_MeshObj*           theMeshObj )
965 : SMESH_SubMeshObj( theMeshObj ),
966   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
967 {
968   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
969   
970   mySubMeshServer->Register();
971 }
972
973 SMESH_subMeshObj::~SMESH_subMeshObj()
974 {
975   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
976   mySubMeshServer->UnRegister();
977 }
978
979 //=================================================================================
980 // function : GetEntities
981 // purpose  : Get entities of specified type. Return number of entities
982 //=================================================================================
983 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
984 {
985   switch ( theType )
986   {
987     case SMDSAbs_Node:
988     {
989       return mySubMeshServer->GetNumberOfNodes( false );
990     }
991     break;
992     case SMDSAbs_0DElement:
993     case SMDSAbs_Edge:
994     case SMDSAbs_Face:
995     case SMDSAbs_Volume:
996     {
997       SMESH::long_array_var anIds = 
998         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
999       return anIds->length();
1000     }
1001     default:
1002       return 0;
1003     break;
1004   }
1005 }
1006
1007 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1008 {
1009   theResList.clear();
1010
1011   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1012   if ( aMesh == 0 )
1013     return 0;
1014
1015   bool isNodal = IsNodePrs();
1016
1017   if ( isNodal )
1018   {
1019     if ( theType == SMDSAbs_Node )
1020     {
1021       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1022       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1023     }
1024   }
1025   else
1026   {
1027     if ( theType == SMDSAbs_Node )
1028     {
1029       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1030       return getNodesFromElems( anIds, aMesh, theResList );
1031     }
1032     else
1033     {
1034       SMESH::long_array_var anIds = 
1035         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1036       return getPointers( theType, anIds, aMesh, theResList );
1037     }
1038   }
1039
1040   return 0;
1041 }
1042
1043 //=================================================================================
1044 // function : IsNodePrs
1045 // purpose  : Return true if node presentation is used
1046 //=================================================================================
1047 bool SMESH_subMeshObj::IsNodePrs() const
1048 {
1049   return mySubMeshServer->GetNumberOfElements() == 0;
1050 }