]> SALOME platform Git repositories - modules/smesh.git/blob - src/OBJECT/SMESH_DeviceActor.cxx
Salome HOME
0021825: Error in the example of "Projection Algorithms" in the user's guide
[modules/smesh.git] / src / OBJECT / SMESH_DeviceActor.cxx
1 // Copyright (C) 2007-2012  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_DeviceActor.cxx
25 //  Author : 
26 //  Module : SMESH
27 //
28 #include "SMESH_DeviceActor.h"
29 #include "SMESH_ScalarBarActor.h"
30 #include "SMESH_ExtractGeometry.h"
31 #include "SMESH_ControlsDef.hxx"
32 #include "SMESH_ActorUtils.h"
33 #include "SMESH_FaceOrientationFilter.h"
34 #include "VTKViewer_CellLocationsArray.h"
35 #include "VTKViewer_PolyDataMapper.h"
36
37 #include <VTKViewer_Transform.h>
38 #include <VTKViewer_TransformFilter.h>
39 #include <VTKViewer_ExtractUnstructuredGrid.h>
40
41 // VTK Includes
42 #include <vtkObjectFactory.h>
43 #include <vtkShrinkFilter.h>
44 #include <vtkShrinkPolyData.h>
45
46 #include <vtkProperty.h>
47 #include <vtkPolyData.h>
48 #include <vtkMergeFilter.h>
49 #include <vtkPolyDataMapper.h>
50 #include <vtkUnstructuredGrid.h>
51
52 #include <vtkLookupTable.h>
53 #include <vtkDoubleArray.h>
54 #include <vtkCellData.h>
55
56 #include <vtkCell.h>
57 #include <vtkIdList.h>
58 #include <vtkCellArray.h>
59 #include <vtkUnsignedCharArray.h>
60
61 #include <vtkImplicitBoolean.h>
62 #include <vtkPassThroughFilter.h>
63
64 #include <vtkRenderer.h>
65
66 #include "utilities.h"
67
68 #ifdef _DEBUG_
69 static int MYDEBUG = 0;
70 #else
71 static int MYDEBUG = 0;
72 #endif
73
74 using namespace std;
75
76
77 vtkStandardNewMacro(SMESH_DeviceActor);
78
79
80 SMESH_DeviceActor
81 ::SMESH_DeviceActor()
82 {
83   if(MYDEBUG) MESSAGE("SMESH_DeviceActor - "<<this);
84
85   myIsShrinkable = false;
86   myIsShrunk = false;
87   myIsHighlited = false;
88
89   myRepresentation = eSurface;
90
91   myProperty = vtkProperty::New();
92   myMapper = VTKViewer_PolyDataMapper::New();
93
94   vtkMapper::GetResolveCoincidentTopologyPolygonOffsetParameters(myPolygonOffsetFactor,
95                                                                  myPolygonOffsetUnits);
96
97   myMapper->UseLookupTableScalarRangeOn();
98   myMapper->SetColorModeToMapScalars();
99
100   myShrinkFilter = vtkShrinkFilter::New();
101
102   myStoreClippingMapping = false;
103
104   myExtractGeometry = SMESH_ExtractGeometry::New();
105   myExtractGeometry->SetReleaseDataFlag(true);
106   myIsImplicitFunctionUsed = false;
107
108   myExtractUnstructuredGrid = VTKViewer_ExtractUnstructuredGrid::New();
109     
110   myMergeFilter = vtkMergeFilter::New();
111
112   myGeomFilter = VTKViewer_GeometryFilter::New();
113
114   myTransformFilter = VTKViewer_TransformFilter::New();
115
116   for(int i = 0; i < 6; i++)
117     myPassFilter.push_back(vtkPassThroughFilter::New());
118
119   // Orientation of faces
120   myIsFacesOriented = false;
121
122   vtkFloatingPointType anRGB[3] = { 1, 1, 1 };
123   SMESH::GetColor( "SMESH", "orientation_color", anRGB[0], anRGB[1], anRGB[2], QColor( 255, 255, 255 ) );
124
125   myFaceOrientationFilter = SMESH_FaceOrientationFilter::New();
126
127   myFaceOrientationDataMapper = vtkPolyDataMapper::New();
128   myFaceOrientationDataMapper->SetInput(myFaceOrientationFilter->GetOutput());
129
130   myFaceOrientation = vtkActor::New();
131   myFaceOrientation->SetMapper(myFaceOrientationDataMapper);
132   myFaceOrientation->GetProperty()->SetColor(anRGB[0], anRGB[1], anRGB[2]);
133 }
134
135
136 SMESH_DeviceActor
137 ::~SMESH_DeviceActor()
138 {
139   if(MYDEBUG) MESSAGE("~SMESH_DeviceActor - "<<this);
140
141   myProperty->Delete();
142
143   myMapper->Delete();
144
145   myShrinkFilter->Delete();
146
147   myExtractUnstructuredGrid->Delete();
148
149   myMergeFilter->Delete();
150
151   myGeomFilter->Delete();
152
153   myExtractGeometry->Delete();
154
155   myTransformFilter->Delete();
156
157   for(int i = 0, iEnd = myPassFilter.size(); i < iEnd; i++){
158     myPassFilter[i]->Delete();
159   }
160
161   // Orientation of faces
162   myFaceOrientationFilter->Delete();
163
164   myFaceOrientationDataMapper->RemoveAllInputs();
165   myFaceOrientationDataMapper->Delete();
166
167   myFaceOrientation->Delete();
168 }
169
170
171 void
172 SMESH_DeviceActor
173 ::SetStoreGemetryMapping(bool theStoreMapping)
174 {
175   myGeomFilter->SetStoreMapping(theStoreMapping);
176   SetStoreClippingMapping(theStoreMapping);
177 }
178
179
180 void
181 SMESH_DeviceActor
182 ::SetStoreClippingMapping(bool theStoreMapping)
183 {
184   myStoreClippingMapping = theStoreMapping;
185   myExtractGeometry->SetStoreMapping(theStoreMapping && myIsImplicitFunctionUsed);
186   SetStoreIDMapping(theStoreMapping);
187 }
188
189
190 void
191 SMESH_DeviceActor
192 ::SetStoreIDMapping(bool theStoreMapping)
193 {
194   myExtractUnstructuredGrid->SetStoreMapping(theStoreMapping);
195 }
196
197
198 void 
199 SMESH_DeviceActor
200 ::Init(TVisualObjPtr theVisualObj, 
201        vtkImplicitBoolean* theImplicitBoolean)
202 {
203   myVisualObj = theVisualObj;
204   myExtractGeometry->SetImplicitFunction(theImplicitBoolean);
205   SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
206 }
207
208
209 void
210 SMESH_DeviceActor
211 ::SetImplicitFunctionUsed(bool theIsImplicitFunctionUsed)
212 {
213   int anId = 0;
214   if(theIsImplicitFunctionUsed)
215     myPassFilter[ anId ]->SetInput( myExtractGeometry->GetOutput() );
216   else
217     myPassFilter[ anId ]->SetInput( myMergeFilter->GetOutput() );
218     
219   myIsImplicitFunctionUsed = theIsImplicitFunctionUsed;
220   SetStoreClippingMapping(myStoreClippingMapping);
221 }
222
223
224 void
225 SMESH_DeviceActor
226 ::SetUnstructuredGrid(vtkUnstructuredGrid* theGrid)
227 {
228   if(theGrid){
229     //myIsShrinkable = theGrid->GetNumberOfCells() > 10;
230     myIsShrinkable = true;
231
232     myExtractUnstructuredGrid->SetInput(theGrid);
233
234     myMergeFilter->SetGeometry(myExtractUnstructuredGrid->GetOutput());
235
236     myExtractGeometry->SetInput(myMergeFilter->GetOutput());
237
238     int anId = 0;
239     SetImplicitFunctionUsed(myIsImplicitFunctionUsed);
240     myPassFilter[ anId + 1]->SetInput( myPassFilter[ anId ]->GetOutput() );
241     
242     anId++; // 1
243     myTransformFilter->SetInput( myPassFilter[ anId ]->GetOutput() );
244
245     anId++; // 2
246     myPassFilter[ anId ]->SetInput( myTransformFilter->GetOutput() );
247     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
248
249     anId++; // 3
250     myGeomFilter->SetInput( myPassFilter[ anId ]->GetOutput() );
251
252     anId++; // 4
253     myPassFilter[ anId ]->SetInput( myGeomFilter->GetOutput() ); 
254     myPassFilter[ anId + 1 ]->SetInput( myPassFilter[ anId ]->GetOutput() );
255
256     anId++; // 5
257     myMapper->SetInput( myPassFilter[ anId ]->GetPolyDataOutput() );
258
259     vtkLODActor::SetMapper( myMapper );
260     Modified();
261   }
262 }
263
264
265 VTKViewer_ExtractUnstructuredGrid* 
266 SMESH_DeviceActor
267 ::GetExtractUnstructuredGrid()
268 {
269   return myExtractUnstructuredGrid;
270 }
271
272
273 vtkUnstructuredGrid* 
274 SMESH_DeviceActor
275 ::GetUnstructuredGrid()
276 {
277   myExtractUnstructuredGrid->Update();
278   return myExtractUnstructuredGrid->GetOutput();
279 }
280
281
282 void
283 SMESH_DeviceActor
284 ::SetControlMode(SMESH::Controls::FunctorPtr theFunctor,
285                  SMESH_ScalarBarActor* theScalarBarActor,
286                  vtkLookupTable* theLookupTable)
287 {
288   bool anIsInitialized = theFunctor;
289   if(anIsInitialized){
290     vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
291
292     SetStoreIDMapping(true);
293     myExtractUnstructuredGrid->Update();
294     vtkUnstructuredGrid* aGrid = myExtractUnstructuredGrid->GetOutput();
295
296     aDataSet->ShallowCopy(aGrid);
297     
298     vtkDoubleArray *aScalars = vtkDoubleArray::New();
299     vtkIdType aNbCells = aGrid->GetNumberOfCells();
300     aScalars->SetNumberOfComponents(1);
301     aScalars->SetNumberOfTuples(aNbCells);
302     
303     myVisualObj->UpdateFunctor(theFunctor);
304
305     using namespace SMESH::Controls;
306     if(NumericalFunctor* aNumericalFunctor = dynamic_cast<NumericalFunctor*>(theFunctor.get())){
307       for(vtkIdType i = 0; i < aNbCells; i++){
308         vtkIdType anId = myExtractUnstructuredGrid->GetInputId(i);
309         vtkIdType anObjId = myVisualObj->GetElemObjId(anId);
310         double aValue = aNumericalFunctor->GetValue(anObjId);
311         aScalars->SetValue(i,aValue);
312       }
313     }else if(Predicate* aPredicate = dynamic_cast<Predicate*>(theFunctor.get())){
314       for(vtkIdType i = 0; i < aNbCells; i++){
315         vtkIdType anId = myExtractUnstructuredGrid->GetInputId(i);
316         vtkIdType anObjId = myVisualObj->GetElemObjId(anId);
317         bool aValue = aPredicate->IsSatisfy(anObjId);
318         aScalars->SetValue(i,aValue);
319       }
320     }
321
322     aDataSet->GetCellData()->SetScalars(aScalars);
323     aScalars->Delete();
324         
325     theLookupTable->SetRange(aScalars->GetRange());
326     theLookupTable->SetNumberOfTableValues(theScalarBarActor->GetMaximumNumberOfColors());
327     theLookupTable->Build();
328     
329     myMergeFilter->SetScalars(aDataSet);
330     aDataSet->Delete();
331   }
332   GetMapper()->SetScalarVisibility(anIsInitialized);
333   theScalarBarActor->SetVisibility(anIsInitialized);
334 }
335
336 void
337 SMESH_DeviceActor
338 ::SetExtControlMode(SMESH::Controls::FunctorPtr theFunctor,
339                     SMESH_ScalarBarActor* theScalarBarActor,
340                     vtkLookupTable* theLookupTable)
341 {
342   bool anIsInitialized = theFunctor;
343   myExtractUnstructuredGrid->ClearRegisteredCells();
344   myExtractUnstructuredGrid->ClearRegisteredCellsWithType();
345   myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::ePassAll);
346   myVisualObj->UpdateFunctor(theFunctor);
347
348   using namespace SMESH::Controls;
349   if (anIsInitialized){
350     if (Length2D* aLength2D = dynamic_cast<Length2D*>(theFunctor.get())){
351       SMESH::Controls::Length2D::TValues aValues;
352
353       aLength2D->GetValues(aValues);
354       vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
355       vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
356
357       aDataSet->SetPoints(aGrid->GetPoints());
358       
359       vtkIdType aNbCells = aValues.size();
360       
361       vtkDoubleArray *aScalars = vtkDoubleArray::New();
362       aScalars->SetNumberOfComponents(1);
363       aScalars->SetNumberOfTuples(aNbCells);
364
365       vtkIdType aCellsSize = 3*aNbCells;
366       vtkCellArray* aConnectivity = vtkCellArray::New();
367       aConnectivity->Allocate( aCellsSize, 0 );
368       
369       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
370       aCellTypesArray->SetNumberOfComponents( 1 );
371       aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
372       
373       vtkIdList *anIdList = vtkIdList::New();
374       anIdList->SetNumberOfIds(2);
375       
376       Length2D::TValues::const_iterator anIter = aValues.begin();
377       for(vtkIdType aVtkId = 0; anIter != aValues.end(); anIter++,aVtkId++){
378         const Length2D::Value& aValue = *anIter;
379         int aNode[2] = {
380           myVisualObj->GetNodeVTKId(aValue.myPntId[0]),
381           myVisualObj->GetNodeVTKId(aValue.myPntId[1])
382         };
383         if(aNode[0] >= 0 && aNode[1] >= 0){
384           anIdList->SetId( 0, aNode[0] );
385           anIdList->SetId( 1, aNode[1] );
386           aConnectivity->InsertNextCell( anIdList );
387           aCellTypesArray->InsertNextValue( VTK_LINE );
388           aScalars->SetValue(aVtkId,aValue.myLength);
389         }
390       }
391       
392       VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
393       aCellLocationsArray->SetNumberOfComponents( 1 );
394       aCellLocationsArray->SetNumberOfTuples( aNbCells );
395       
396       aConnectivity->InitTraversal();
397       for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
398         aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
399       
400       aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
401       SetUnstructuredGrid(aDataSet);
402
403       aDataSet->GetCellData()->SetScalars(aScalars);
404       aScalars->Delete();
405       
406       theLookupTable->SetRange(aScalars->GetRange());
407       theLookupTable->Build();
408       
409       myMergeFilter->SetScalars(aDataSet);
410       aDataSet->Delete();
411     }
412     else if (MultiConnection2D* aMultiConnection2D = dynamic_cast<MultiConnection2D*>(theFunctor.get())){
413       SMESH::Controls::MultiConnection2D::MValues aValues;
414
415       aMultiConnection2D->GetValues(aValues);
416       vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
417       vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
418       aDataSet->SetPoints(aGrid->GetPoints());
419       
420       vtkIdType aNbCells = aValues.size();
421       vtkDoubleArray *aScalars = vtkDoubleArray::New();
422       aScalars->SetNumberOfComponents(1);
423       aScalars->SetNumberOfTuples(aNbCells);
424
425       vtkIdType aCellsSize = 3*aNbCells;
426       vtkCellArray* aConnectivity = vtkCellArray::New();
427       aConnectivity->Allocate( aCellsSize, 0 );
428       
429       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
430       aCellTypesArray->SetNumberOfComponents( 1 );
431       aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
432       
433       vtkIdList *anIdList = vtkIdList::New();
434       anIdList->SetNumberOfIds(2);
435       
436       MultiConnection2D::MValues::const_iterator anIter = aValues.begin();
437       for(vtkIdType aVtkId = 0; anIter != aValues.end(); anIter++,aVtkId++){
438         const MultiConnection2D::Value& aValue = (*anIter).first;
439         int aNode[2] = {
440           myVisualObj->GetNodeVTKId(aValue.myPntId[0]),
441           myVisualObj->GetNodeVTKId(aValue.myPntId[1])
442         };
443         if(aNode[0] >= 0 && aNode[1] >= 0){
444           anIdList->SetId( 0, aNode[0] );
445           anIdList->SetId( 1, aNode[1] );
446           aConnectivity->InsertNextCell( anIdList );
447           aCellTypesArray->InsertNextValue( VTK_LINE );
448           aScalars->SetValue(aVtkId,(*anIter).second);
449         }
450       }
451       
452       VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
453       aCellLocationsArray->SetNumberOfComponents( 1 );
454       aCellLocationsArray->SetNumberOfTuples( aNbCells );
455       
456       aConnectivity->InitTraversal();
457       for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
458         aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
459       
460       aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
461       SetUnstructuredGrid(aDataSet);
462
463       aDataSet->GetCellData()->SetScalars(aScalars);
464       aScalars->Delete();
465       
466       theLookupTable->SetRange(aScalars->GetRange());
467       theLookupTable->Build();
468       
469       myMergeFilter->SetScalars(aDataSet);
470       aDataSet->Delete();
471     }
472   }
473   GetMapper()->SetScalarVisibility(anIsInitialized);
474   theScalarBarActor->SetVisibility(anIsInitialized);
475 }
476
477 void
478 SMESH_DeviceActor
479 ::SetExtControlMode(SMESH::Controls::FunctorPtr theFunctor)
480 {
481   myExtractUnstructuredGrid->ClearRegisteredCells();
482   myExtractUnstructuredGrid->ClearRegisteredCellsWithType();
483   myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::ePassAll);
484   myVisualObj->UpdateFunctor(theFunctor);
485
486   using namespace SMESH::Controls;
487   if ( dynamic_cast<FreeBorders          *>(theFunctor.get()) ||
488        dynamic_cast<FreeFaces            *>(theFunctor.get()) ||
489        dynamic_cast<BareBorderVolume     *>(theFunctor.get()) ||
490        dynamic_cast<BareBorderFace       *>(theFunctor.get()) ||
491        dynamic_cast<OverConstrainedVolume*>(theFunctor.get()) ||
492        dynamic_cast<CoincidentElements1D *>(theFunctor.get()) ||
493        dynamic_cast<CoincidentElements2D *>(theFunctor.get()) ||
494        dynamic_cast<CoincidentElements3D *>(theFunctor.get()) ||
495        dynamic_cast<OverConstrainedFace  *>(theFunctor.get()))
496   {
497     Predicate* aPredicate = dynamic_cast<Predicate*>(theFunctor.get());
498     myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::eAdding);
499     vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
500     vtkIdType aNbCells = aGrid->GetNumberOfCells();
501     for( vtkIdType i = 0; i < aNbCells; i++ ){
502       vtkIdType anObjId = myVisualObj->GetElemObjId(i);
503       if(aPredicate->IsSatisfy(anObjId))
504         myExtractUnstructuredGrid->RegisterCell(i);
505     }
506     if(!myExtractUnstructuredGrid->IsCellsRegistered())
507       myExtractUnstructuredGrid->RegisterCell(-1);
508     SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
509   }
510   else if(FreeEdges* aFreeEdges = dynamic_cast<FreeEdges*>(theFunctor.get()))
511   {
512     SMESH::Controls::FreeEdges::TBorders aBorders;
513     aFreeEdges->GetBoreders(aBorders);
514     vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
515     vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
516     aDataSet->SetPoints(aGrid->GetPoints());
517
518     vtkIdType aNbCells = aBorders.size();
519     vtkIdType aCellsSize = 3*aNbCells;
520     vtkCellArray* aConnectivity = vtkCellArray::New();
521     aConnectivity->Allocate( aCellsSize, 0 );
522     
523     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
524     aCellTypesArray->SetNumberOfComponents( 1 );
525     aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
526     
527     vtkIdList *anIdList = vtkIdList::New();
528     anIdList->SetNumberOfIds(2);
529     
530     FreeEdges::TBorders::const_iterator anIter = aBorders.begin();
531     for(; anIter != aBorders.end(); anIter++){
532       const FreeEdges::Border& aBorder = *anIter;
533       int aNode[2] = {
534         myVisualObj->GetNodeVTKId(aBorder.myPntId[0]),
535         myVisualObj->GetNodeVTKId(aBorder.myPntId[1])
536       };
537       //cout<<"aNode = "<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<endl;
538       if(aNode[0] >= 0 && aNode[1] >= 0){
539         anIdList->SetId( 0, aNode[0] );
540         anIdList->SetId( 1, aNode[1] );
541         aConnectivity->InsertNextCell( anIdList );
542         aCellTypesArray->InsertNextValue( VTK_LINE );
543       }
544     }
545     
546     VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
547     aCellLocationsArray->SetNumberOfComponents( 1 );
548     aCellLocationsArray->SetNumberOfTuples( aNbCells );
549     
550     aConnectivity->InitTraversal();
551     for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
552       aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
553     
554     aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
555
556     SetUnstructuredGrid(aDataSet);
557     aDataSet->Delete();
558   }
559   else if(dynamic_cast<FreeNodes      *>(theFunctor.get()) ||
560           dynamic_cast<CoincidentNodes*>(theFunctor.get()))
561   {
562     Predicate* aPredicate = dynamic_cast<Predicate*>(theFunctor.get());
563     myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::eAdding);
564     vtkIdType aNbNodes = myVisualObj->GetNbEntities(SMDSAbs_Node);
565     for( vtkIdType i = 0; i < aNbNodes; i++ ){
566       vtkIdType anObjId = myVisualObj->GetNodeObjId(i);
567       if(aPredicate->IsSatisfy(anObjId))
568         myExtractUnstructuredGrid->RegisterCell(i);
569     }
570     if(!myExtractUnstructuredGrid->IsCellsRegistered())
571       myExtractUnstructuredGrid->RegisterCell(-1);
572     SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
573   }
574 }
575
576
577
578
579 unsigned long int 
580 SMESH_DeviceActor
581 ::GetMTime()
582 {
583   unsigned long mTime = this->Superclass::GetMTime();
584   mTime = max(mTime,myExtractGeometry->GetMTime());
585   mTime = max(mTime,myExtractUnstructuredGrid->GetMTime());
586   mTime = max(mTime,myMergeFilter->GetMTime());
587   mTime = max(mTime,myGeomFilter->GetMTime());
588   mTime = max(mTime,myTransformFilter->GetMTime());
589   mTime = max(mTime,myFaceOrientationFilter->GetMTime());
590   return mTime;
591 }
592
593
594 void
595 SMESH_DeviceActor
596 ::SetTransform(VTKViewer_Transform* theTransform)
597 {
598   myTransformFilter->SetTransform(theTransform);
599 }
600
601
602 void
603 SMESH_DeviceActor
604 ::SetShrink() 
605 {
606   if ( !myIsShrinkable ) return;
607   if ( vtkDataSet* aDataSet = myPassFilter[ 0 ]->GetOutput() )
608   {
609     myShrinkFilter->SetInput( aDataSet );
610     myPassFilter[ 1 ]->SetInput( myShrinkFilter->GetOutput() );
611     myIsShrunk = true;
612   }
613 }
614
615 void
616 SMESH_DeviceActor
617 ::UnShrink() 
618 {
619   if ( !myIsShrunk ) return;
620   if ( vtkDataSet* aDataSet = myPassFilter[ 0 ]->GetOutput() )
621   {    
622     myPassFilter[ 1 ]->SetInput( aDataSet );
623     myPassFilter[ 1 ]->Modified();
624     myIsShrunk = false;
625     Modified();
626   }
627 }
628
629
630 void
631 SMESH_DeviceActor
632 ::SetFacesOriented(bool theIsFacesOriented) 
633 {
634   if ( vtkDataSet* aDataSet = myTransformFilter->GetOutput() )
635   {
636     myIsFacesOriented = theIsFacesOriented;
637     if( theIsFacesOriented )
638       myFaceOrientationFilter->SetInput( aDataSet );
639     UpdateFaceOrientation();
640   }
641 }
642
643 void
644 SMESH_DeviceActor
645 ::SetFacesOrientationColor(vtkFloatingPointType theColor[3])
646 {
647   myFaceOrientation->GetProperty()->SetColor( theColor );
648 }
649
650 void
651 SMESH_DeviceActor
652 ::GetFacesOrientationColor(vtkFloatingPointType theColor[3])
653 {
654   myFaceOrientation->GetProperty()->GetColor( theColor );
655 }
656
657 void
658 SMESH_DeviceActor
659 ::SetFacesOrientationScale(vtkFloatingPointType theScale)
660 {
661   myFaceOrientationFilter->SetOrientationScale( theScale );
662 }
663
664 vtkFloatingPointType
665 SMESH_DeviceActor
666 ::GetFacesOrientationScale()
667 {
668   return myFaceOrientationFilter->GetOrientationScale();
669 }
670
671 void
672 SMESH_DeviceActor
673 ::SetFacesOrientation3DVectors(bool theState)
674 {
675   myFaceOrientationFilter->Set3dVectors( theState );
676 }
677
678 bool
679 SMESH_DeviceActor
680 ::GetFacesOrientation3DVectors()
681 {
682   return myFaceOrientationFilter->Get3dVectors();
683 }
684
685 void
686 SMESH_DeviceActor
687 ::UpdateFaceOrientation()
688 {
689   bool aShowFaceOrientation = myIsFacesOriented;
690   aShowFaceOrientation &= GetVisibility();
691   aShowFaceOrientation &= myRepresentation == eSurface;
692   myFaceOrientation->SetVisibility(aShowFaceOrientation);
693 }
694
695
696 void
697 SMESH_DeviceActor
698 ::SetRepresentation(EReperesent theMode)
699 {
700   switch(theMode){
701   case ePoint:
702     myGeomFilter->SetInside(true);
703     myGeomFilter->SetWireframeMode(false);
704     GetProperty()->SetRepresentation(0);
705     break;
706   case eWireframe:
707     myGeomFilter->SetInside(false);
708     myGeomFilter->SetWireframeMode(true);
709     GetProperty()->SetRepresentation(theMode);
710     break;
711   case eInsideframe:
712     myGeomFilter->SetInside(true);
713     myGeomFilter->SetWireframeMode(true);
714     GetProperty()->SetRepresentation(1);
715     break;
716   case eSurface:
717     myGeomFilter->SetInside(false);
718     myGeomFilter->SetWireframeMode(false);
719     GetProperty()->SetRepresentation(theMode);
720   }
721   SetMarkerEnabled(theMode == ePoint);
722   myRepresentation = theMode;
723   UpdateFaceOrientation();
724   GetProperty()->Modified();
725   myMapper->Modified();
726   Modified();
727 }
728
729
730 void
731 SMESH_DeviceActor
732 ::SetVisibility(int theMode)
733 {
734   if(!myExtractUnstructuredGrid->GetInput() || 
735      GetUnstructuredGrid()->GetNumberOfCells())
736   {
737     vtkLODActor::SetVisibility(theMode);
738   }else{
739     vtkLODActor::SetVisibility(false);
740   }
741   UpdateFaceOrientation();
742 }
743
744
745 int
746 SMESH_DeviceActor
747 ::GetVisibility()
748 {
749   if(!GetUnstructuredGrid()->GetNumberOfCells()){
750     vtkLODActor::SetVisibility(false);
751   }
752   return vtkLODActor::GetVisibility();
753 }
754
755
756 void
757 SMESH_DeviceActor
758 ::AddToRender(vtkRenderer* theRenderer)
759 {
760   theRenderer->AddActor(this);
761   theRenderer->AddActor(myFaceOrientation);
762 }
763
764 void
765 SMESH_DeviceActor
766 ::RemoveFromRender(vtkRenderer* theRenderer)
767 {
768   theRenderer->RemoveActor(this);
769   theRenderer->RemoveActor(myFaceOrientation);
770 }
771
772
773 int
774 SMESH_DeviceActor
775 ::GetNodeObjId(int theVtkID)
776 {
777   vtkIdType anID = theVtkID;
778
779   if(IsImplicitFunctionUsed())
780     anID = myExtractGeometry->GetNodeObjId(theVtkID);
781
782   vtkIdType aRetID = myVisualObj->GetNodeObjId(anID);
783   if(MYDEBUG) MESSAGE("GetNodeObjId - theVtkID = "<<theVtkID<<"; anID = "<<anID<<"; aRetID = "<<aRetID);
784   return aRetID;
785 }
786
787 vtkFloatingPointType* 
788 SMESH_DeviceActor
789 ::GetNodeCoord(int theObjID)
790 {
791   vtkDataSet* aDataSet = myMergeFilter->GetOutput();
792   vtkIdType anID = myVisualObj->GetNodeVTKId(theObjID);
793   vtkFloatingPointType* aCoord = (anID >=0) ? aDataSet->GetPoint(anID) : NULL;
794   if(MYDEBUG) MESSAGE("GetNodeCoord - theObjID = "<<theObjID<<"; anID = "<<anID);
795   return aCoord;
796 }
797
798
799 int
800 SMESH_DeviceActor
801 ::GetElemObjId(int theVtkID)
802 {
803   vtkIdType anId = myGeomFilter->GetElemObjId(theVtkID);
804   if(anId < 0) 
805     return -1;
806
807   vtkIdType anId2 = anId;
808   if(IsImplicitFunctionUsed())
809     anId2 = myExtractGeometry->GetElemObjId(anId);
810   if(anId2 < 0) 
811     return -1;
812
813   vtkIdType anId3 = myExtractUnstructuredGrid->GetInputId(anId2);
814   if(anId3 < 0) 
815     return -1;
816
817   vtkIdType aRetID = myVisualObj->GetElemObjId(anId3);
818   if(MYDEBUG) 
819      MESSAGE("GetElemObjId - theVtkID = "<<theVtkID<<"; anId2 = "<<anId2<<"; anId3 = "<<anId3<<"; aRetID = "<<aRetID);
820   return aRetID;
821 }
822
823 vtkCell* 
824 SMESH_DeviceActor
825 ::GetElemCell(int theObjID)
826 {
827   vtkDataSet* aDataSet = myVisualObj->GetUnstructuredGrid();
828   vtkIdType aGridID = myVisualObj->GetElemVTKId(theObjID);
829   vtkCell* aCell = (aGridID >= 0 ) ? aDataSet->GetCell(aGridID) : NULL;
830   if(MYDEBUG) 
831     MESSAGE("GetElemCell - theObjID = "<<theObjID<<"; aGridID = "<<aGridID);
832   return aCell;
833 }
834
835
836 vtkFloatingPointType 
837 SMESH_DeviceActor
838 ::GetShrinkFactor()
839 {
840   return myShrinkFilter->GetShrinkFactor();
841 }
842
843 void
844 SMESH_DeviceActor
845 ::SetShrinkFactor(vtkFloatingPointType theValue)
846 {
847   theValue = theValue > 0.1? theValue: 0.8;
848   myShrinkFilter->SetShrinkFactor(theValue);
849   Modified();
850 }
851
852
853 void
854 SMESH_DeviceActor
855 ::SetHighlited(bool theIsHighlited)
856 {
857   if ( myIsHighlited == theIsHighlited )
858     return;
859   myIsHighlited = theIsHighlited;
860   Modified();
861 }
862
863 void
864 SMESH_DeviceActor
865 ::Render(vtkRenderer *ren, vtkMapper* m)
866 {
867   int aResolveCoincidentTopology = vtkMapper::GetResolveCoincidentTopology();
868   vtkFloatingPointType aStoredFactor, aStoredUnit; 
869   vtkMapper::GetResolveCoincidentTopologyPolygonOffsetParameters(aStoredFactor,aStoredUnit);
870
871   vtkMapper::SetResolveCoincidentTopologyToPolygonOffset();
872   vtkFloatingPointType aFactor = myPolygonOffsetFactor, aUnits = myPolygonOffsetUnits;
873   if(myIsHighlited){
874     static vtkFloatingPointType EPS = .01;
875     aUnits *= (1.0-EPS);
876   }
877   vtkMapper::SetResolveCoincidentTopologyPolygonOffsetParameters(aFactor,aUnits);
878   vtkLODActor::Render(ren,m);
879
880   vtkMapper::SetResolveCoincidentTopologyPolygonOffsetParameters(aStoredFactor,aStoredUnit);
881   vtkMapper::SetResolveCoincidentTopology(aResolveCoincidentTopology);
882 }
883
884
885 void
886 SMESH_DeviceActor
887 ::SetPolygonOffsetParameters(vtkFloatingPointType factor, 
888                              vtkFloatingPointType units)
889 {
890   myPolygonOffsetFactor = factor;
891   myPolygonOffsetUnits = units;
892 }
893
894 /*!
895  * On/Off representation 2D quadratic element as arked polygon
896  */
897 void SMESH_DeviceActor::SetQuadraticArcMode(bool theFlag){
898   myGeomFilter->SetQuadraticArcMode(theFlag);
899 }
900
901 /*!
902  * Return true if 2D quadratic element displayed as arked polygon
903  */
904 bool SMESH_DeviceActor::GetQuadraticArcMode(){
905   return myGeomFilter->GetQuadraticArcMode();
906 }
907 /*!
908  * Set Max angle for representation 2D quadratic element as arked polygon
909  */
910 void SMESH_DeviceActor::SetQuadraticArcAngle(vtkFloatingPointType theMaxAngle){
911   myGeomFilter->SetQuadraticArcAngle(theMaxAngle);
912 }
913
914 /*!
915  * Return Max angle of the representation 2D quadratic element as arked polygon
916  */
917 vtkFloatingPointType SMESH_DeviceActor::GetQuadraticArcAngle(){
918   return myGeomFilter->GetQuadraticArcAngle();
919 }
920
921 /*!
922  * Set point marker enabled
923  * \param theMarkerEnabled flag to enable/disable point marker
924  */
925 void SMESH_DeviceActor::SetMarkerEnabled( bool theMarkerEnabled )
926 {
927   myMapper->SetMarkerEnabled( theMarkerEnabled );
928 }
929
930 /*!
931  * Set standard point marker
932  * \param theMarkerType type of the marker
933  */
934 void SMESH_DeviceActor::SetMarkerStd( VTK::MarkerType theMarkerType, VTK::MarkerScale theMarkerScale )
935 {
936   myMapper->SetMarkerStd( theMarkerType, theMarkerScale );
937 }
938
939 /*!
940  * Set custom point marker
941  * \param theMarkerId id of the marker texture
942  * \param theMarkerTexture marker texture
943  */
944 void SMESH_DeviceActor::SetMarkerTexture( int theMarkerId, VTK::MarkerTexture theMarkerTexture )
945 {
946   myMapper->SetMarkerTexture( theMarkerId, theMarkerTexture );
947 }
948
949 /*!
950  * Get type of the point marker
951  * \return type of the point marker
952  */
953 VTK::MarkerType SMESH_DeviceActor::GetMarkerType()
954 {
955   return myMapper->GetMarkerType();
956 }
957
958 /*!
959   Get scale of the point marker
960   \return scale of the point marker
961 */
962 VTK::MarkerScale SMESH_DeviceActor::GetMarkerScale()
963 {
964   return myMapper->GetMarkerScale();
965 }
966
967 /*!
968  * Get texture identifier of the point marker
969  * \return texture identifier of the point marker
970  */
971 int SMESH_DeviceActor::GetMarkerTexture()
972 {
973   return myMapper->GetMarkerTexture();
974 }
975
976 void SMESH_DeviceActor::SetCoincident3DAllowed(bool theFlag) {
977   myGeomFilter->SetAppendCoincident3D(theFlag);
978 }
979
980 bool SMESH_DeviceActor::IsCoincident3DAllowed() const {
981   return myGeomFilter->GetAppendCoincident3D();
982 }