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