Salome HOME
[Bug PAL7444] display mesh takes a lot of more memory in 2.1.0 than in 2.0.0.
[modules/smesh.git] / src / OBJECT / SMESH_DeviceActor.cxx
1 //  SMESH OBJECT : interactive object for SMESH visualization
2 //
3 //  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
5 // 
6 //  This library is free software; you can redistribute it and/or 
7 //  modify it under the terms of the GNU Lesser General Public 
8 //  License as published by the Free Software Foundation; either 
9 //  version 2.1 of the License. 
10 // 
11 //  This library is distributed in the hope that it will be useful, 
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of 
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
14 //  Lesser General Public License for more details. 
15 // 
16 //  You should have received a copy of the GNU Lesser General Public 
17 //  License along with this library; if not, write to the Free Software 
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
19 // 
20 //  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
21 //
22 //
23 //
24 //  File   : SMESH_DeviceActor.cxx
25 //  Author : 
26 //  Module : SMESH
27 //  $Header$
28
29
30 #include "SMESH_DeviceActor.h"
31 #include "SMESH_ExtractGeometry.h"
32 #include "SMESH_ControlsDef.hxx"
33 #include "SMESH_ActorUtils.h"
34
35 #include "SALOME_Transform.h"
36 #include "SALOME_TransformFilter.h"
37 #include "SALOME_PassThroughFilter.h"
38 #include "SALOME_ExtractUnstructuredGrid.h"
39
40 // VTK Includes
41 #include <vtkObjectFactory.h>
42 #include <vtkShrinkFilter.h>
43 #include <vtkShrinkPolyData.h>
44
45 #include <vtkProperty.h>
46 #include <vtkPolyData.h>
47 #include <vtkMergeFilter.h>
48 #include <vtkPolyDataMapper.h>
49 #include <vtkUnstructuredGrid.h>
50
51 #include <vtkScalarBarActor.h>
52 #include <vtkLookupTable.h>
53 #include <vtkDoubleArray.h>
54 #include <vtkCellData.h>
55
56 #include <vtkCell.h>
57 #include <vtkIdList.h>
58 #include <vtkIntArray.h>
59 #include <vtkCellArray.h>
60 #include <vtkUnsignedCharArray.h>
61
62 #include <vtkImplicitBoolean.h>
63
64 #include "utilities.h"
65
66 #ifdef _DEBUG_
67 static int MYDEBUG = 0;
68 #else
69 static int MYDEBUG = 0;
70 #endif
71
72 using namespace std;
73
74
75 vtkStandardNewMacro(SMESH_DeviceActor);
76
77
78 SMESH_DeviceActor::SMESH_DeviceActor(){
79   myIsShrunk = false;
80   myIsShrinkable = false;
81   myRepresentation = eSurface;
82
83   myProperty = vtkProperty::New();
84   myMapper = vtkPolyDataMapper::New();
85
86   vtkMapper::GetResolveCoincidentTopologyPolygonOffsetParameters(myPolygonOffsetFactor,
87                                                                  myPolygonOffsetUnits);
88
89   myMapper->UseLookupTableScalarRangeOn();
90   myMapper->SetColorModeToMapScalars();
91
92   myShrinkFilter = vtkShrinkFilter::New();
93
94   myExtractGeometry = SMESH_ExtractGeometry::New();
95   myExtractGeometry->SetReleaseDataFlag(true);
96   myExtractGeometry->SetStoreMapping(true);
97   myIsImplicitFunctionUsed = false;
98
99   myExtractUnstructuredGrid = SALOME_ExtractUnstructuredGrid::New();
100   myExtractUnstructuredGrid->SetStoreMapping(true);
101
102   myMergeFilter = vtkMergeFilter::New();
103
104   myStoreMapping = false;
105   myGeomFilter = SALOME_GeometryFilter::New();
106
107   myTransformFilter = SALOME_TransformFilter::New();
108
109   for(int i = 0; i < 6; i++)
110     myPassFilter.push_back(SALOME_PassThroughFilter::New());
111 }
112
113
114 SMESH_DeviceActor::~SMESH_DeviceActor(){
115   if(MYDEBUG) MESSAGE("~SMESH_DeviceActor");
116   myProperty->Delete();
117
118   myMapper->RemoveAllInputs();
119   myMapper->Delete();
120
121   myShrinkFilter->UnRegisterAllOutputs();
122   myShrinkFilter->Delete();
123
124   myExtractUnstructuredGrid->UnRegisterAllOutputs();
125   myExtractUnstructuredGrid->Delete();
126
127   myMergeFilter->UnRegisterAllOutputs();
128   myMergeFilter->Delete();
129
130   myGeomFilter->UnRegisterAllOutputs();
131   myGeomFilter->Delete();
132
133   myExtractGeometry->UnRegisterAllOutputs();
134   myExtractGeometry->Delete();
135
136   myTransformFilter->UnRegisterAllOutputs();
137   myTransformFilter->Delete();
138
139   for(int i = 0, iEnd = myPassFilter.size(); i < iEnd; i++){
140     myPassFilter[i]->UnRegisterAllOutputs(); 
141     myPassFilter[i]->Delete();
142   }
143 }
144
145
146 void SMESH_DeviceActor::SetStoreMapping(int theStoreMapping){
147   if (myStoreMapping == theStoreMapping)
148     return;
149   myStoreMapping = theStoreMapping;
150   myGeomFilter->SetStoreMapping( myStoreMapping );
151   Modified();
152 }
153
154
155 void SMESH_DeviceActor::Init(TVisualObjPtr theVisualObj, 
156                              vtkImplicitBoolean* theImplicitBoolean)
157 {
158   myVisualObj = theVisualObj;
159   myExtractGeometry->SetImplicitFunction(theImplicitBoolean);
160   SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
161 }
162
163
164 void
165 SMESH_DeviceActor::
166 SetImplicitFunctionUsed(bool theIsImplicitFunctionUsed)
167 {
168   if(myIsImplicitFunctionUsed == theIsImplicitFunctionUsed)
169     return;
170
171   int anId = 0;
172   if(theIsImplicitFunctionUsed)
173     myPassFilter[ anId ]->SetInput( myExtractGeometry->GetOutput() );
174   else
175     myPassFilter[ anId ]->SetInput( myMergeFilter->GetOutput() );
176     
177   myIsImplicitFunctionUsed = theIsImplicitFunctionUsed;
178 }
179
180
181 void SMESH_DeviceActor::SetUnstructuredGrid(vtkUnstructuredGrid* theGrid){
182   if(theGrid){
183     //myIsShrinkable = theGrid->GetNumberOfCells() > 10;
184     myIsShrinkable = true;
185
186     myExtractUnstructuredGrid->SetInput(theGrid);
187
188     myMergeFilter->SetGeometry(myExtractUnstructuredGrid->GetOutput());
189
190     myExtractGeometry->SetInput(myMergeFilter->GetOutput());
191
192     int anId = 0;
193     myPassFilter[ anId ]->SetInput( myMergeFilter->GetOutput() );
194     myPassFilter[ anId + 1]->SetInput( myPassFilter[ anId ]->GetOutput() );
195     
196     anId++; // 1
197     myGeomFilter->SetStoreMapping( myStoreMapping );
198     myGeomFilter->SetInput( myPassFilter[ anId ]->GetOutput() );
199
200     anId++; // 2
201     myPassFilter[ anId ]->SetInput( myGeomFilter->GetOutput() ); 
202     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
203
204     anId++; // 3
205     myTransformFilter->SetInput( myPassFilter[ anId ]->GetPolyDataOutput() );
206     myTransformFilter->SetInput( myPassFilter[ anId ]->GetPolyDataOutput() );
207
208     anId++; // 4
209     myPassFilter[ anId ]->SetInput( myTransformFilter->GetOutput() );
210     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
211     myPassFilter[ anId ]->SetInput( myTransformFilter->GetOutput() );
212     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
213
214     anId++; // 5
215     myMapper->SetInput( myPassFilter[ anId ]->GetPolyDataOutput() );
216     myMapper->SetInput( myPassFilter[ anId ]->GetPolyDataOutput() );
217
218     vtkLODActor::SetMapper( myMapper );
219     Modified();
220   }
221 }
222
223
224 SALOME_ExtractUnstructuredGrid* SMESH_DeviceActor::GetExtractUnstructuredGrid(){
225   return myExtractUnstructuredGrid;
226 }
227
228
229 vtkUnstructuredGrid* SMESH_DeviceActor::GetUnstructuredGrid(){
230   myExtractUnstructuredGrid->Update();
231   return myExtractUnstructuredGrid->GetOutput();
232 }
233
234
235 void SMESH_DeviceActor::SetControlMode(SMESH::Controls::FunctorPtr theFunctor,
236                                        vtkScalarBarActor* theScalarBarActor,
237                                        vtkLookupTable* theLookupTable)
238 {
239   bool anIsInitialized = theFunctor;
240   if(anIsInitialized){
241     vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
242     vtkUnstructuredGrid* aGrid = myExtractUnstructuredGrid->GetOutput();
243     aDataSet->ShallowCopy(aGrid);
244     
245     vtkDoubleArray *aScalars = vtkDoubleArray::New();
246     vtkIdType aNbCells = aGrid->GetNumberOfCells();
247     aScalars->SetNumberOfComponents(1);
248     aScalars->SetNumberOfTuples(aNbCells);
249     
250     myVisualObj->UpdateFunctor(theFunctor);
251
252     using namespace SMESH::Controls;
253     if(NumericalFunctor* aNumericalFunctor = dynamic_cast<NumericalFunctor*>(theFunctor.get())){
254       for(vtkIdType i = 0; i < aNbCells; i++){
255         vtkIdType anId = myExtractUnstructuredGrid->GetInputId(i);
256         vtkIdType anObjId = myVisualObj->GetElemObjId(anId);
257         double aValue = aNumericalFunctor->GetValue(anObjId);
258         aScalars->SetValue(i,aValue);
259       }
260     }else if(Predicate* aPredicate = dynamic_cast<Predicate*>(theFunctor.get())){
261       for(vtkIdType i = 0; i < aNbCells; i++){
262         vtkIdType anId = myExtractUnstructuredGrid->GetInputId(i);
263         vtkIdType anObjId = myVisualObj->GetElemObjId(anId);
264         bool aValue = aPredicate->IsSatisfy(anObjId);
265         aScalars->SetValue(i,aValue);
266       }
267     }
268
269     aDataSet->GetCellData()->SetScalars(aScalars);
270     aScalars->Delete();
271         
272     theLookupTable->SetRange(aScalars->GetRange());
273     theLookupTable->Build();
274     
275     myMergeFilter->SetScalars(aDataSet);
276     aDataSet->Delete();
277   }
278   GetMapper()->SetScalarVisibility(anIsInitialized);
279   theScalarBarActor->SetVisibility(anIsInitialized);
280 }
281
282 void SMESH_DeviceActor::SetExtControlMode(SMESH::Controls::FunctorPtr theFunctor,
283                                           SMESH_DeviceActor* theDeviceActor,
284                                           vtkScalarBarActor* theScalarBarActor,
285                                           vtkLookupTable* theLookupTable)
286 {
287   bool anIsInitialized = theFunctor;
288   myExtractUnstructuredGrid->ClearRegisteredCells();
289   myExtractUnstructuredGrid->ClearRegisteredCellsWithType();
290   myExtractUnstructuredGrid->SetModeOfChanging(SALOME_ExtractUnstructuredGrid::ePassAll);
291   myVisualObj->UpdateFunctor(theFunctor);
292
293   using namespace SMESH::Controls;
294   if (anIsInitialized){
295     if (Length2D* aLength2D = dynamic_cast<Length2D*>(theFunctor.get())){
296       SMESH::Controls::Length2D::TValues aValues;
297
298       aLength2D->GetValues(aValues);
299       vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
300       vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
301
302       aDataSet->SetPoints(aGrid->GetPoints());
303       
304       vtkIdType aNbCells = aValues.size();
305       
306       vtkDoubleArray *aScalars = vtkDoubleArray::New();
307       aScalars->SetNumberOfComponents(1);
308       aScalars->SetNumberOfTuples(aNbCells);
309
310       vtkIdType aCellsSize = 3*aNbCells;
311       vtkCellArray* aConnectivity = vtkCellArray::New();
312       aConnectivity->Allocate( aCellsSize, 0 );
313       
314       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
315       aCellTypesArray->SetNumberOfComponents( 1 );
316       aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
317       
318       vtkIdList *anIdList = vtkIdList::New();
319       anIdList->SetNumberOfIds(2);
320       
321       Length2D::TValues::const_iterator anIter = aValues.begin();
322       for(vtkIdType aVtkId = 0; anIter != aValues.end(); anIter++,aVtkId++){
323         const Length2D::Value& aValue = *anIter;
324         int aNode[2] = {
325           myVisualObj->GetNodeVTKId(aValue.myPntId[0]),
326           myVisualObj->GetNodeVTKId(aValue.myPntId[1])
327         };
328         if(aNode[0] >= 0 && aNode[1] >= 0){
329           anIdList->SetId( 0, aNode[0] );
330           anIdList->SetId( 1, aNode[1] );
331           aConnectivity->InsertNextCell( anIdList );
332           aCellTypesArray->InsertNextValue( VTK_LINE );
333           aScalars->SetValue(aVtkId,aValue.myLength);
334         }
335       }
336       
337       vtkIntArray* aCellLocationsArray = vtkIntArray::New();
338       aCellLocationsArray->SetNumberOfComponents( 1 );
339       aCellLocationsArray->SetNumberOfTuples( aNbCells );
340       
341       aConnectivity->InitTraversal();
342       for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
343         aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
344       
345       aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
346       SetUnstructuredGrid(aDataSet);
347
348       aDataSet->GetCellData()->SetScalars(aScalars);
349       aScalars->Delete();
350       
351       theLookupTable->SetRange(aScalars->GetRange());
352       theLookupTable->Build();
353       
354       myMergeFilter->SetScalars(aDataSet);
355       aDataSet->Delete();
356     }
357     else if (MultiConnection2D* aMultiConnection2D = dynamic_cast<MultiConnection2D*>(theFunctor.get())){
358       SMESH::Controls::MultiConnection2D::MValues aValues;
359
360       aMultiConnection2D->GetValues(aValues);
361       vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
362       vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
363       aDataSet->SetPoints(aGrid->GetPoints());
364       
365       vtkIdType aNbCells = aValues.size();
366       vtkDoubleArray *aScalars = vtkDoubleArray::New();
367       aScalars->SetNumberOfComponents(1);
368       aScalars->SetNumberOfTuples(aNbCells);
369
370       vtkIdType aCellsSize = 3*aNbCells;
371       vtkCellArray* aConnectivity = vtkCellArray::New();
372       aConnectivity->Allocate( aCellsSize, 0 );
373       
374       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
375       aCellTypesArray->SetNumberOfComponents( 1 );
376       aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
377       
378       vtkIdList *anIdList = vtkIdList::New();
379       anIdList->SetNumberOfIds(2);
380       
381       MultiConnection2D::MValues::const_iterator anIter = aValues.begin();
382       for(vtkIdType aVtkId = 0; anIter != aValues.end(); anIter++,aVtkId++){
383         const MultiConnection2D::Value& aValue = (*anIter).first;
384         int aNode[2] = {
385           myVisualObj->GetNodeVTKId(aValue.myPntId[0]),
386           myVisualObj->GetNodeVTKId(aValue.myPntId[1])
387         };
388         if(aNode[0] >= 0 && aNode[1] >= 0){
389           anIdList->SetId( 0, aNode[0] );
390           anIdList->SetId( 1, aNode[1] );
391           aConnectivity->InsertNextCell( anIdList );
392           aCellTypesArray->InsertNextValue( VTK_LINE );
393           aScalars->SetValue(aVtkId,(*anIter).second);
394         }
395       }
396       
397       vtkIntArray* aCellLocationsArray = vtkIntArray::New();
398       aCellLocationsArray->SetNumberOfComponents( 1 );
399       aCellLocationsArray->SetNumberOfTuples( aNbCells );
400       
401       aConnectivity->InitTraversal();
402       for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
403         aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
404       
405       aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
406       SetUnstructuredGrid(aDataSet);
407
408       aDataSet->GetCellData()->SetScalars(aScalars);
409       aScalars->Delete();
410       
411       theLookupTable->SetRange(aScalars->GetRange());
412       theLookupTable->Build();
413       
414       myMergeFilter->SetScalars(aDataSet);
415       aDataSet->Delete();
416     }
417   }
418   GetMapper()->SetScalarVisibility(anIsInitialized);
419   theScalarBarActor->SetVisibility(anIsInitialized);
420 }
421
422 void SMESH_DeviceActor::SetExtControlMode(SMESH::Controls::FunctorPtr theFunctor,
423                                           SMESH_DeviceActor* theDeviceActor)
424 {
425   myExtractUnstructuredGrid->ClearRegisteredCells();
426   myExtractUnstructuredGrid->ClearRegisteredCellsWithType();
427   myExtractUnstructuredGrid->SetModeOfChanging(SALOME_ExtractUnstructuredGrid::ePassAll);
428   myVisualObj->UpdateFunctor(theFunctor);
429
430   using namespace SMESH::Controls;
431   if(FreeBorders* aFreeBorders = dynamic_cast<FreeBorders*>(theFunctor.get())){
432     myExtractUnstructuredGrid->SetModeOfChanging(SALOME_ExtractUnstructuredGrid::eAdding);
433     vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
434     vtkIdType aNbCells = aGrid->GetNumberOfCells();
435     for( vtkIdType i = 0; i < aNbCells; i++ ){
436       vtkIdType anObjId = myVisualObj->GetElemObjId(i);
437       if(aFreeBorders->IsSatisfy(anObjId))
438         myExtractUnstructuredGrid->RegisterCell(i);
439     }
440     if(!myExtractUnstructuredGrid->IsCellsRegistered())
441       myExtractUnstructuredGrid->RegisterCell(-1);
442     SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
443   }else if(FreeEdges* aFreeEdges = dynamic_cast<FreeEdges*>(theFunctor.get())){
444     SMESH::Controls::FreeEdges::TBorders aBorders;
445     aFreeEdges->GetBoreders(aBorders);
446     vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
447     vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
448     aDataSet->SetPoints(aGrid->GetPoints());
449
450     vtkIdType aNbCells = aBorders.size();
451     vtkIdType aCellsSize = 3*aNbCells;
452     vtkCellArray* aConnectivity = vtkCellArray::New();
453     aConnectivity->Allocate( aCellsSize, 0 );
454     
455     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
456     aCellTypesArray->SetNumberOfComponents( 1 );
457     aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
458     
459     vtkIdList *anIdList = vtkIdList::New();
460     anIdList->SetNumberOfIds(2);
461     
462     FreeEdges::TBorders::const_iterator anIter = aBorders.begin();
463     for(; anIter != aBorders.end(); anIter++){
464       const FreeEdges::Border& aBorder = *anIter;
465       int aNode[2] = {
466         myVisualObj->GetNodeVTKId(aBorder.myPntId[0]),
467         myVisualObj->GetNodeVTKId(aBorder.myPntId[1])
468       };
469       //cout<<"aNode = "<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<endl;
470       if(aNode[0] >= 0 && aNode[1] >= 0){
471         anIdList->SetId( 0, aNode[0] );
472         anIdList->SetId( 1, aNode[1] );
473         aConnectivity->InsertNextCell( anIdList );
474         aCellTypesArray->InsertNextValue( VTK_LINE );
475       }
476     }
477     
478     vtkIntArray* aCellLocationsArray = vtkIntArray::New();
479     aCellLocationsArray->SetNumberOfComponents( 1 );
480     aCellLocationsArray->SetNumberOfTuples( aNbCells );
481     
482     aConnectivity->InitTraversal();
483     for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
484       aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
485     
486     aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
487
488     SetUnstructuredGrid(aDataSet);
489     aDataSet->Delete();
490   }
491 }
492
493
494
495
496 unsigned long int SMESH_DeviceActor::GetMTime(){
497   unsigned long mTime = this->Superclass::GetMTime();
498   mTime = max(mTime,myExtractGeometry->GetMTime());
499   mTime = max(mTime,myExtractUnstructuredGrid->GetMTime());
500   mTime = max(mTime,myMergeFilter->GetMTime());
501   mTime = max(mTime,myGeomFilter->GetMTime());
502   mTime = max(mTime,myTransformFilter->GetMTime());
503   return mTime;
504 }
505
506
507 void SMESH_DeviceActor::SetTransform(SALOME_Transform* theTransform){
508   myTransformFilter->SetTransform(theTransform);
509 }
510
511
512 void SMESH_DeviceActor::SetShrink() {
513   if ( !myIsShrinkable ) return;
514   if ( vtkDataSet* aDataSet = myPassFilter[ 0 ]->GetOutput() )
515   {
516     myShrinkFilter->SetInput( aDataSet );
517     myPassFilter[ 1 ]->SetInput( myShrinkFilter->GetOutput() );
518     myIsShrunk = true;
519   }
520 }
521
522 void SMESH_DeviceActor::UnShrink() {
523   if ( !myIsShrunk ) return;
524   if ( vtkDataSet* aDataSet = myPassFilter[ 0 ]->GetOutput() )
525   {    
526     myPassFilter[ 1 ]->SetInput( aDataSet );
527     myPassFilter[ 1 ]->Modified();
528     myIsShrunk = false;
529     Modified();
530   }
531 }
532
533
534 void SMESH_DeviceActor::SetRepresentation(EReperesent theMode){
535   switch(theMode){
536   case ePoint:
537     myGeomFilter->SetInside(true);
538     GetProperty()->SetRepresentation(0);
539     break;
540   case eInsideframe:
541     myGeomFilter->SetInside(true);
542     GetProperty()->SetRepresentation(1);
543     break;
544   default :
545     myGeomFilter->SetInside(false);
546     GetProperty()->SetRepresentation(theMode);
547   }
548   myRepresentation = theMode;
549   GetProperty()->Modified();
550   myMapper->Modified();
551   Modified();
552 }
553
554
555 void SMESH_DeviceActor::SetVisibility(int theMode){
556   if(!myExtractUnstructuredGrid->GetInput() || GetUnstructuredGrid()->GetNumberOfCells()){
557     vtkLODActor::SetVisibility(theMode);
558   }else{
559     vtkLODActor::SetVisibility(false);
560   }
561 }
562
563
564 int SMESH_DeviceActor::GetVisibility(){
565   if(!GetUnstructuredGrid()->GetNumberOfCells()){
566     vtkLODActor::SetVisibility(false);
567   }
568   return vtkLODActor::GetVisibility();
569 }
570
571
572 int SMESH_DeviceActor::GetNodeObjId(int theVtkID){
573   vtkIdType anID = theVtkID;
574
575   if(IsImplicitFunctionUsed())
576     anID = myExtractGeometry->GetNodeObjId(theVtkID);
577
578   vtkIdType aRetID = myVisualObj->GetNodeObjId(anID);
579   if(MYDEBUG) MESSAGE("GetNodeObjId - theVtkID = "<<theVtkID<<"; aRetID = "<<aRetID);
580   return aRetID;
581 }
582
583 float* SMESH_DeviceActor::GetNodeCoord(int theObjID){
584   vtkDataSet* aDataSet = myMergeFilter->GetOutput();
585   vtkIdType anID = myVisualObj->GetNodeVTKId(theObjID);
586   float* aCoord = aDataSet->GetPoint(anID);
587   if(MYDEBUG) MESSAGE("GetNodeCoord - theObjID = "<<theObjID<<"; anID = "<<anID);
588   return aCoord;
589 }
590
591
592 int SMESH_DeviceActor::GetElemObjId(int theVtkID){
593   vtkIdType anId = myGeomFilter->GetElemObjId(theVtkID);
594   if(anId < 0) 
595     return -1;
596
597   vtkIdType anId2 = anId;
598   if(IsImplicitFunctionUsed())
599     anId2 = myExtractGeometry->GetElemObjId(anId);
600   if(anId2 < 0) 
601     return -1;
602
603   vtkIdType anId3 = myExtractUnstructuredGrid->GetInputId(anId2);
604   if(anId3 < 0) 
605     return -1;
606
607   vtkIdType aRetID = myVisualObj->GetElemObjId(anId3);
608   if(MYDEBUG) 
609      MESSAGE("GetElemObjId - theVtkID = "<<theVtkID<<"; anId2 = "<<anId2<<"; anId3 = "<<anId3<<"; aRetID = "<<aRetID);
610   return aRetID;
611 }
612
613 vtkCell* SMESH_DeviceActor::GetElemCell(int theObjID){
614   vtkDataSet* aDataSet = myVisualObj->GetUnstructuredGrid();
615   vtkIdType aGridID = myVisualObj->GetElemVTKId(theObjID);
616   vtkCell* aCell = aDataSet->GetCell(aGridID);
617   if(MYDEBUG) 
618     MESSAGE("GetElemCell - theObjID = "<<theObjID<<"; aGridID = "<<aGridID);
619   return aCell;
620 }
621
622
623 float SMESH_DeviceActor::GetShrinkFactor(){
624   return myShrinkFilter->GetShrinkFactor();
625 }
626
627 void SMESH_DeviceActor::SetShrinkFactor(float theValue){
628   theValue = theValue > 0.1? theValue: 0.8;
629   myShrinkFilter->SetShrinkFactor(theValue);
630   Modified();
631 }
632
633
634 void SMESH_DeviceActor::SetHighlited(bool theIsHighlited){
635   myIsHighlited = theIsHighlited;
636   Modified();
637 }
638
639 void SMESH_DeviceActor::Render(vtkRenderer *ren, vtkMapper* m){
640   int aResolveCoincidentTopology = vtkMapper::GetResolveCoincidentTopology();
641   float aStoredFactor, aStoredUnit; 
642   vtkMapper::GetResolveCoincidentTopologyPolygonOffsetParameters(aStoredFactor,aStoredUnit);
643
644   vtkMapper::SetResolveCoincidentTopologyToPolygonOffset();
645   float aFactor = myPolygonOffsetFactor, aUnits = myPolygonOffsetUnits;
646   if(myIsHighlited){
647     static float EPS = .01;
648     aUnits *= (1.0-EPS);
649   }
650   vtkMapper::SetResolveCoincidentTopologyPolygonOffsetParameters(aFactor,aUnits);
651   vtkLODActor::Render(ren,m);
652
653   vtkMapper::SetResolveCoincidentTopologyPolygonOffsetParameters(aStoredFactor,aStoredUnit);
654   vtkMapper::SetResolveCoincidentTopology(aResolveCoincidentTopology);
655 }
656
657
658 void SMESH_DeviceActor::SetPolygonOffsetParameters(float factor, float units){
659   myPolygonOffsetFactor = factor;
660   myPolygonOffsetUnits = units;
661 }
662