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