Salome HOME
Add missing files for make dist step
[modules/smesh.git] / src / OBJECT / SMESH_Object.cxx
1 // Copyright (C) 2007-2011  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         if((*anIter)->GetEntityType() != SMDSEntity_Polyhedra &&
373            (*anIter)->GetEntityType() != SMDSEntity_Quad_Polyhedra) {
374           aCellsSize += (*anIter)->NbNodes() + 1;
375         } 
376         // Special case for the VTK_POLYHEDRON:
377         // itsinput cellArray is of special format.
378         //  [nCellFaces, nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...]   
379         else {
380           if( const SMDS_VtkVolume* ph = dynamic_cast<const SMDS_VtkVolume*>(*anIter) ) {
381             int nbFaces = ph->NbFaces();
382             aCellsSize += (1 + ph->NbFaces());
383             for( int i = 1; i <= nbFaces; i++ ) {
384               aCellsSize += ph->NbFaceNodes(i);
385             }
386           }
387         }
388       }
389     }
390   }
391   
392   vtkIdType aNbCells = nbEnts[ SMDSAbs_0DElement ] + nbEnts[ SMDSAbs_Edge ] +
393                        nbEnts[ SMDSAbs_Face ] + nbEnts[ SMDSAbs_Volume ];
394
395   if ( MYDEBUG )
396     MESSAGE( "Update - aNbCells = "<<aNbCells<<"; aCellsSize = "<<aCellsSize );
397
398   // Create cells
399
400   vtkCellArray* aConnectivity = vtkCellArray::New();
401   aConnectivity->Allocate( aCellsSize, 0 );
402
403   SMDS_Mesh::CheckMemory(); // PAL16631
404
405   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
406   aCellTypesArray->SetNumberOfComponents( 1 );
407   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
408
409   SMDS_Mesh::CheckMemory(); // PAL16631
410
411   vtkIdList *anIdList = vtkIdList::New();
412   vtkIdType iElem = 0;
413
414   TConnect aConnect;
415   aConnect.reserve(VTK_CELL_SIZE);
416
417   SMDS_Mesh::CheckMemory(); // PAL16631
418
419   for ( int i = 0; i <= 3; i++ ) // iterate through 0d elements, edges, faces and volumes
420   {
421     if ( nbEnts[ aTypes[ i ] ] > 0 ) {
422       
423       const SMDSAbs_ElementType& aType = aTypes[ i ];
424       const TEntityList& aList = anEnts[ aType ];
425       TEntityList::const_iterator anIter;
426       for ( anIter = aList.begin(); anIter != aList.end(); ++anIter )
427       {
428         const SMDS_MeshElement* anElem = *anIter;
429         
430         vtkIdType aNbNodes = anElem->NbNodes();
431         anIdList->SetNumberOfIds( aNbNodes );
432         
433         int anId = anElem->GetID();
434         
435         mySMDS2VTKElems.insert( TMapOfIds::value_type( anId, iElem ) );
436         myVTK2SMDSElems.insert( TMapOfIds::value_type( iElem, anId ) );
437         
438         SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
439         switch (aType) {
440         case SMDSAbs_Volume: {
441           aConnect.clear();
442           std::vector<int> aConnectivities;
443           // Convertions connectivities from SMDS to VTK
444           
445           if (anElem->IsPoly() && aNbNodes > 3) { // POLYEDRE
446             anIdList->Reset();
447             if ( const SMDS_VtkVolume* ph = dynamic_cast<const SMDS_VtkVolume*>(anElem) ) {
448               int nbFaces = ph->NbFaces();
449               anIdList->InsertNextId(nbFaces);
450               for( int i = 1; i <= nbFaces; i++ ) {
451                 anIdList->InsertNextId(ph->NbFaceNodes(i));
452                 for(int j = 1; j <= ph->NbFaceNodes(i); j++) {
453                   const SMDS_MeshNode* n = ph->GetFaceNode(i,j);
454                   if(n) {
455                     anIdList->InsertNextId(mySMDS2VTKNodes[n->GetID()]);
456                   }
457                 }
458               }
459             }
460             
461           } else if (aNbNodes == 4) {
462             static int anIds[] = {0,2,1,3};
463             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
464             
465           } else if (aNbNodes == 5) {
466             static int anIds[] = {0,3,2,1,4};
467             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
468
469           } else if (aNbNodes == 6) {
470             static int anIds[] = {0,1,2,3,4,5};
471             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
472
473           }
474           else if (aNbNodes == 8) {
475             static int anIds[] = {0,3,2,1,4,7,6,5};
476             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
477
478           }
479           else if (aNbNodes == 10) {
480             static int anIds[] = {0,2,1,3,6,5,4,7,9,8};
481             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
482           }
483           else if (aNbNodes == 13) {
484             static int anIds[] = {0,3,2,1,4,8,7,6,5,9,12,11,10};
485             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
486           }
487           else if (aNbNodes == 15) {
488             //static int anIds[] = {0,2,1,3,5,4,8,7,6,11,10,9,12,14,13};
489             static int anIds[] = {0,1,2,3,4,5,6,7,8,9,10,11,12,13,14};
490             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
491             //for (int k = 0; k < aNbNodes; k++) {
492             //  int nn = aConnectivities[k];
493             //  const SMDS_MeshNode* N = static_cast<const SMDS_MeshNode*> (aConnect[nn]);
494             //  cout<<"k="<<k<<"  N("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
495             //}
496           }
497           else if (aNbNodes == 20) {
498             static int anIds[] = {0,3,2,1,4,7,6,5,11,10,9,8,15,14,13,12,16,19,18,17};
499             for (int k = 0; k < aNbNodes; k++) aConnectivities.push_back(anIds[k]);
500           }
501           else {
502           }
503
504           if (!(anElem->IsPoly() && aNbNodes > 3)) {
505             if ( aConnect.empty() )
506               GetConnect(aNodesIter,aConnect);
507
508             if (aConnectivities.size() > 0) {
509               for (vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++)
510                 SetId(anIdList,mySMDS2VTKNodes,aConnect,aNodeId,aConnectivities[aNodeId]);
511             }
512           }
513           break;
514         }
515         default:
516           for( vtkIdType aNodeId = 0; aNodesIter->more(); aNodeId++ ){
517             const SMDS_MeshElement* aNode = aNodesIter->next();
518             anIdList->SetId( aNodeId, mySMDS2VTKNodes[aNode->GetID()] );
519           }
520         }
521
522         
523         aConnectivity->InsertNextCell( anIdList );
524         aCellTypesArray->InsertNextValue( getCellType( aType, anElem->IsPoly(), aNbNodes ) );
525
526         iElem++;
527       }
528     }
529     SMDS_Mesh::CheckMemory(); // PAL16631
530   }
531
532   // Insert cells in grid
533
534   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
535   aCellLocationsArray->SetNumberOfComponents( 1 );
536   aCellLocationsArray->SetNumberOfTuples( aNbCells );
537
538   SMDS_Mesh::CheckMemory(); // PAL16631
539
540   aConnectivity->InitTraversal();
541   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
542     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
543
544   myGrid->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
545
546   aCellLocationsArray->Delete();
547   aCellTypesArray->Delete();
548   aConnectivity->Delete();
549   anIdList->Delete();
550
551   SMDS_Mesh::CheckMemory(); // PAL16631
552 }
553
554 //=================================================================================
555 // function : GetEdgeNodes
556 // purpose  : Retrieve ids of nodes from edge of elements ( edge is numbered from 1 )
557 //=================================================================================
558 bool SMESH_VisualObjDef::GetEdgeNodes( const int theElemId,
559                                        const int theEdgeNum,
560                                        int&      theNodeId1,
561                                        int&      theNodeId2 ) const
562 {
563   const SMDS_Mesh* aMesh = GetMesh();
564   if ( aMesh == 0 )
565     return false;
566     
567   const SMDS_MeshElement* anElem = aMesh->FindElement( theElemId );
568   if ( anElem == 0 )
569     return false;
570     
571   int nbNodes = anElem->NbNodes();
572
573   if ( theEdgeNum < 0 || theEdgeNum > 3 || (nbNodes != 3 && nbNodes != 4) || theEdgeNum > nbNodes )
574     return false;
575
576   vector<int> anIds( nbNodes );
577   SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
578   int i = 0;
579   while( anIter->more() )
580     anIds[ i++ ] = anIter->next()->GetID();
581
582   if ( theEdgeNum < nbNodes - 1 )
583   {
584     theNodeId1 = anIds[ theEdgeNum ];
585     theNodeId2 = anIds[ theEdgeNum + 1 ];
586   }
587   else
588   {
589     theNodeId1 = anIds[ nbNodes - 1 ];
590     theNodeId2 = anIds[ 0 ];
591   }
592
593   return true;
594 }
595
596 vtkUnstructuredGrid* SMESH_VisualObjDef::GetUnstructuredGrid()
597 {
598         //MESSAGE("SMESH_VisualObjDef::GetUnstructuredGrid " << myGrid);
599         return myGrid;
600 }
601
602
603 //=================================================================================
604 // function : IsValid
605 // purpose  : Return true if there are some entities
606 //=================================================================================
607 bool SMESH_VisualObjDef::IsValid() const
608 {
609         //MESSAGE("SMESH_VisualObjDef::IsValid");
610   return GetNbEntities(SMDSAbs_Node) > 0      || 
611          GetNbEntities(SMDSAbs_0DElement) > 0 || 
612          GetNbEntities(SMDSAbs_Edge) > 0      || 
613          GetNbEntities(SMDSAbs_Face) > 0      ||
614          GetNbEntities(SMDSAbs_Volume) > 0 ;
615 }
616
617 /*
618   Class       : SMESH_MeshObj
619   Description : Class for visualisation of mesh
620 */
621
622 //=================================================================================
623 // function : SMESH_MeshObj
624 // purpose  : Constructor
625 //=================================================================================
626 SMESH_MeshObj::SMESH_MeshObj(SMESH::SMESH_Mesh_ptr theMesh):
627   myClient(SalomeApp_Application::orb(),theMesh)
628 {
629         myEmptyGrid = 0;
630   if ( MYDEBUG ) 
631     MESSAGE("SMESH_MeshObj - this = "<<this<<"; theMesh->_is_nil() = "<<theMesh->_is_nil());
632 }
633
634 //=================================================================================
635 // function : ~SMESH_MeshObj
636 // purpose  : Destructor
637 //=================================================================================
638 SMESH_MeshObj::~SMESH_MeshObj()
639 {
640   if ( MYDEBUG ) 
641     MESSAGE("SMESH_MeshObj - this = "<<this<<"\n");
642 }
643
644 //=================================================================================
645 // function : Update
646 // purpose  : Update mesh and fill grid with new values if necessary 
647 //=================================================================================
648 bool SMESH_MeshObj::Update( int theIsClear )
649 {
650   // Update SMDS_Mesh on client part
651   MESSAGE("SMESH_MeshObj::Update " << this);
652   if ( myClient.Update(theIsClear) || GetUnstructuredGrid()->GetNumberOfPoints()==0) {
653     MESSAGE("buildPrs");
654     buildPrs();  // Fill unstructured grid
655     return true;
656   }
657   return false;
658 }
659
660 bool SMESH_MeshObj::NulData()
661 {
662         MESSAGE ("SMESH_MeshObj::NulData() ==================================================================================");
663         if (!myEmptyGrid)
664         {
665           myEmptyGrid = SMDS_UnstructuredGrid::New();
666           myEmptyGrid->Initialize();
667           myEmptyGrid->Allocate();
668           vtkPoints* points = vtkPoints::New();
669           points->SetNumberOfPoints(0);
670           myEmptyGrid->SetPoints( points );
671           points->Delete();
672           myEmptyGrid->BuildLinks();
673         }
674         myGrid->ShallowCopy(myEmptyGrid);
675         return true;
676 }
677 //=================================================================================
678 // function : GetElemDimension
679 // purpose  : Get dimension of element
680 //=================================================================================
681 int SMESH_MeshObj::GetElemDimension( const int theObjId )
682 {
683   const SMDS_MeshElement* anElem = myClient->FindElement( theObjId );
684   if ( anElem == 0 )
685     return 0;
686
687   int aType = anElem->GetType();
688   switch ( aType )
689   {
690     case SMDSAbs_0DElement : return 0;
691     case SMDSAbs_Edge  : return 1;
692     case SMDSAbs_Face  : return 2;
693     case SMDSAbs_Volume: return 3;
694     default            : return 0;
695   }
696 }
697
698 //=================================================================================
699 // function : GetEntities
700 // purpose  : Get entities of specified type. Return number of entities
701 //=================================================================================
702 int SMESH_MeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
703 {
704   switch ( theType )
705   {
706     case SMDSAbs_Node:
707     {
708       return myClient->NbNodes();
709     }
710     break;
711     case SMDSAbs_0DElement:
712     {
713       return myClient->Nb0DElements();
714     }
715     break;
716     case SMDSAbs_Edge:
717     {
718       return myClient->NbEdges();
719     }
720     break;
721     case SMDSAbs_Face:
722     {
723       return myClient->NbFaces();
724     }
725     break;
726     case SMDSAbs_Volume:
727     {
728       return myClient->NbVolumes();
729     }
730     break;
731     default:
732       return 0;
733     break;
734   }
735 }
736
737 int SMESH_MeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theObjs ) const
738 {
739   theObjs.clear();
740
741   switch ( theType )
742   {
743     case SMDSAbs_Node:
744     {
745       SMDS_NodeIteratorPtr anIter = myClient->nodesIterator();
746       while ( anIter->more() ) theObjs.push_back( anIter->next() );
747     }
748     break;
749     case SMDSAbs_0DElement:
750     {
751       SMDS_0DElementIteratorPtr anIter = myClient->elements0dIterator();
752       while ( anIter->more() ) theObjs.push_back( anIter->next() );
753     }
754     break;
755     case SMDSAbs_Edge:
756     {
757       SMDS_EdgeIteratorPtr anIter = myClient->edgesIterator();
758       while ( anIter->more() ) theObjs.push_back( anIter->next() );
759     }
760     break;
761     case SMDSAbs_Face:
762     {
763       SMDS_FaceIteratorPtr anIter = myClient->facesIterator();
764       while ( anIter->more() ) theObjs.push_back( anIter->next() );
765     }
766     break;
767     case SMDSAbs_Volume:
768     {
769       SMDS_VolumeIteratorPtr anIter = myClient->volumesIterator();
770       while ( anIter->more() ) theObjs.push_back( anIter->next() );
771     }
772     break;
773     default:
774     break;
775   }
776
777   return theObjs.size();
778 }
779
780 //=================================================================================
781 // function : UpdateFunctor
782 // purpose  : Update functor in accordance with current mesh
783 //=================================================================================
784 void SMESH_MeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
785 {
786   theFunctor->SetMesh( GetMesh() );
787 }
788
789 //=================================================================================
790 // function : IsNodePrs
791 // purpose  : Return true if node presentation is used
792 //=================================================================================
793 bool SMESH_MeshObj::IsNodePrs() const
794 {
795   return myClient->Nb0DElements() == 0 && myClient->NbEdges() == 0 && myClient->NbFaces() == 0 && myClient->NbVolumes() == 0 ;
796 }
797
798
799 /*
800   Class       : SMESH_SubMeshObj
801   Description : Base class for visualisation of submeshes and groups
802 */
803
804 //=================================================================================
805 // function : SMESH_SubMeshObj
806 // purpose  : Constructor
807 //=================================================================================
808 SMESH_SubMeshObj::SMESH_SubMeshObj( SMESH_MeshObj* theMeshObj )
809 {
810   if ( MYDEBUG ) MESSAGE( "SMESH_SubMeshObj - theMeshObj = " << theMeshObj );
811   
812   myMeshObj = theMeshObj;
813 }
814
815 SMESH_SubMeshObj::~SMESH_SubMeshObj()
816 {
817 }
818
819 //=================================================================================
820 // function : GetElemDimension
821 // purpose  : Get dimension of element
822 //=================================================================================
823 int SMESH_SubMeshObj::GetElemDimension( const int theObjId )
824 {
825   return myMeshObj == 0 ? 0 : myMeshObj->GetElemDimension( theObjId );
826 }
827
828 //=================================================================================
829 // function : UpdateFunctor
830 // purpose  : Update functor in accordance with current mesh
831 //=================================================================================
832 void SMESH_SubMeshObj::UpdateFunctor( const SMESH::Controls::FunctorPtr& theFunctor )
833 {
834   theFunctor->SetMesh( myMeshObj->GetMesh() );
835 }
836
837 //=================================================================================
838 // function : Update
839 // purpose  : Update mesh object and fill grid with new values 
840 //=================================================================================
841 bool SMESH_SubMeshObj::Update( int theIsClear )
842 {
843         MESSAGE("SMESH_SubMeshObj::Update " << this)
844   bool changed = myMeshObj->Update( theIsClear );
845   buildPrs(true);
846   return changed;
847 }
848
849
850 /*
851   Class       : SMESH_GroupObj
852   Description : Class for visualisation of groups
853 */
854
855 //=================================================================================
856 // function : SMESH_GroupObj
857 // purpose  : Constructor
858 //=================================================================================
859 SMESH_GroupObj::SMESH_GroupObj( SMESH::SMESH_GroupBase_ptr theGroup, 
860                                 SMESH_MeshObj*             theMeshObj )
861 : SMESH_SubMeshObj( theMeshObj ),
862   myGroupServer( SMESH::SMESH_GroupBase::_duplicate(theGroup) )
863 {
864   if ( MYDEBUG ) MESSAGE("SMESH_GroupObj - theGroup->_is_nil() = "<<theGroup->_is_nil());
865   myGroupServer->Register();
866 }
867
868 SMESH_GroupObj::~SMESH_GroupObj()
869 {
870   if ( MYDEBUG ) MESSAGE("~SMESH_GroupObj");
871   myGroupServer->UnRegister();
872 }
873
874 //=================================================================================
875 // function : IsNodePrs
876 // purpose  : Return true if node presentation is used
877 //=================================================================================
878 bool SMESH_GroupObj::IsNodePrs() const
879 {
880   return myGroupServer->GetType() == SMESH::NODE;
881 }
882
883 //=================================================================================
884 // function : GetElementType
885 // purpose  : Return type of elements of group
886 //=================================================================================
887 SMDSAbs_ElementType SMESH_GroupObj::GetElementType() const
888 {
889   return SMDSAbs_ElementType(myGroupServer->GetType());
890 }
891
892 //=================================================================================
893 // function : getNodesFromElems
894 // purpose  : Retrieve nodes from elements
895 //=================================================================================
896 static int getNodesFromElems( SMESH::long_array_var&              theElemIds,
897                               const SMDS_Mesh*                    theMesh,
898                               std::list<const SMDS_MeshElement*>& theResList )
899 {
900   set<const SMDS_MeshElement*> aNodeSet;
901
902   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
903   {
904     const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
905     if ( anElem != 0 )
906     {
907       SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
908       while ( anIter->more() )
909       {
910         const SMDS_MeshElement* aNode = anIter->next();
911         if ( aNode != 0 )
912           aNodeSet.insert( aNode );
913       }
914     }
915   }
916
917   set<const SMDS_MeshElement*>::const_iterator anIter;
918   for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
919     theResList.push_back( *anIter );
920
921   return theResList.size();    
922 }
923
924 //=================================================================================
925 // function : getPointers
926 // purpose  : Get std::list<const SMDS_MeshElement*> from list of IDs
927 //=================================================================================
928 static int getPointers( const SMDSAbs_ElementType            theRequestType,
929                         SMESH::long_array_var&              theElemIds,
930                         const SMDS_Mesh*                    theMesh,
931                         std::list<const SMDS_MeshElement*>& theResList )
932 {
933   for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
934   {
935     const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
936       ? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
937
938     if ( anElem != 0 )
939       theResList.push_back( anElem );
940   }
941
942   return theResList.size();
943 }
944
945
946 //=================================================================================
947 // function : GetEntities
948 // purpose  : Get entities of specified type. Return number of entities
949 //=================================================================================
950 int SMESH_GroupObj::GetNbEntities( const SMDSAbs_ElementType theType) const
951 {
952   if(SMDSAbs_ElementType(myGroupServer->GetType()) == theType){
953     return myGroupServer->Size();
954   }
955   return 0;
956 }
957
958 int SMESH_GroupObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
959 {
960   theResList.clear();
961   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
962   
963   if ( myGroupServer->Size() == 0 || aMesh == 0 )
964     return 0;
965
966   SMDSAbs_ElementType aGrpType = SMDSAbs_ElementType(myGroupServer->GetType());
967   SMESH::long_array_var anIds = myGroupServer->GetListOfID();
968
969   if ( aGrpType == theType )
970     return getPointers( theType, anIds, aMesh, theResList );
971   else if ( theType == SMDSAbs_Node )
972     return getNodesFromElems( anIds, aMesh, theResList );
973   else
974     return 0;
975 }    
976
977
978
979 /*
980   Class       : SMESH_subMeshObj
981   Description : Class for visualisation of submeshes
982 */
983
984 //=================================================================================
985 // function : SMESH_subMeshObj
986 // purpose  : Constructor
987 //=================================================================================
988 SMESH_subMeshObj::SMESH_subMeshObj( SMESH::SMESH_subMesh_ptr theSubMesh,
989                                     SMESH_MeshObj*           theMeshObj )
990 : SMESH_SubMeshObj( theMeshObj ),
991   mySubMeshServer( SMESH::SMESH_subMesh::_duplicate( theSubMesh ) )
992 {
993   if ( MYDEBUG ) MESSAGE( "SMESH_subMeshObj - theSubMesh->_is_nil() = " << theSubMesh->_is_nil() );
994   
995   mySubMeshServer->Register();
996 }
997
998 SMESH_subMeshObj::~SMESH_subMeshObj()
999 {
1000   if ( MYDEBUG ) MESSAGE( "~SMESH_subMeshObj" );
1001   mySubMeshServer->UnRegister();
1002 }
1003
1004 //=================================================================================
1005 // function : GetEntities
1006 // purpose  : Get entities of specified type. Return number of entities
1007 //=================================================================================
1008 int SMESH_subMeshObj::GetNbEntities( const SMDSAbs_ElementType theType) const
1009 {
1010   switch ( theType )
1011   {
1012     case SMDSAbs_Node:
1013     {
1014       return mySubMeshServer->GetNumberOfNodes( false );
1015     }
1016     break;
1017     case SMDSAbs_0DElement:
1018     case SMDSAbs_Edge:
1019     case SMDSAbs_Face:
1020     case SMDSAbs_Volume:
1021     {
1022       SMESH::long_array_var anIds = 
1023         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1024       return anIds->length();
1025     }
1026     default:
1027       return 0;
1028     break;
1029   }
1030 }
1031
1032 int SMESH_subMeshObj::GetEntities( const SMDSAbs_ElementType theType, TEntityList& theResList ) const
1033 {
1034   theResList.clear();
1035
1036   SMDS_Mesh* aMesh = myMeshObj->GetMesh();
1037   if ( aMesh == 0 )
1038     return 0;
1039
1040   bool isNodal = IsNodePrs();
1041
1042   if ( isNodal )
1043   {
1044     if ( theType == SMDSAbs_Node )
1045     {
1046       SMESH::long_array_var anIds = mySubMeshServer->GetNodesId();
1047       return getPointers( SMDSAbs_Node, anIds, aMesh, theResList );
1048     }
1049   }
1050   else
1051   {
1052     if ( theType == SMDSAbs_Node )
1053     {
1054       SMESH::long_array_var anIds = mySubMeshServer->GetElementsId();
1055       return getNodesFromElems( anIds, aMesh, theResList );
1056     }
1057     else
1058     {
1059       SMESH::long_array_var anIds = 
1060         mySubMeshServer->GetElementsByType( SMESH::ElementType(theType) );
1061       return getPointers( theType, anIds, aMesh, theResList );
1062     }
1063   }
1064
1065   return 0;
1066 }
1067
1068 //=================================================================================
1069 // function : IsNodePrs
1070 // purpose  : Return true if node presentation is used
1071 //=================================================================================
1072 bool SMESH_subMeshObj::IsNodePrs() const
1073 {
1074   return mySubMeshServer->GetNumberOfElements() == 0;
1075 }