Salome HOME
[bos #40653][CEA] New mesh import export formats with meshio.
[modules/smesh.git] / src / OBJECT / SMESH_DeviceActor.cxx
1 // Copyright (C) 2007-2024  CEA, EDF, 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     else if (Warping3D* aWarping3D = dynamic_cast<Warping3D*>(theFunctor.get())){
511
512       SMESH::Controls::Warping3D::WValues aValues;
513
514       aWarping3D->GetValues(aValues);
515       vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
516       vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
517
518       aDataSet->SetPoints(aGrid->GetPoints());
519
520       vtkIdType aNbCells = aValues.size();
521
522       vtkDoubleArray* aScalars = vtkDoubleArray::New();
523       aScalars->SetNumberOfComponents(1);
524       aScalars->SetNumberOfTuples(aNbCells);
525
526       vtkIdType aCellsSize = 3 * aNbCells;
527       vtkCellArray* aConnectivity = vtkCellArray::New();
528       aConnectivity->Allocate(aCellsSize, 0);
529
530       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
531       aCellTypesArray->SetNumberOfComponents(1);
532       aCellTypesArray->Allocate(aNbCells* aCellTypesArray->GetNumberOfComponents());
533
534       Warping3D::WValues::const_iterator anIter = aValues.begin();
535       aNbCells = 0;
536       for (; anIter != aValues.end(); anIter++) {
537
538         const Warping3D::Value& aValue = *anIter;
539         vtkIdList* anIdList = vtkIdList::New();
540         anIdList->SetNumberOfIds(aValue.myPntIds.size());
541         bool isExist = true;
542         for (int i = 0; i < aValue.myPntIds.size(); ++i)
543         {
544           int aVTKId = myVisualObj->GetNodeVTKId(aValue.myPntIds[i]);
545           if (aVTKId < 0)
546           {
547             isExist = false;
548             break;
549           }
550           anIdList->SetId(i, aVTKId);
551         }
552         if (isExist)
553         {
554           aConnectivity->InsertNextCell(anIdList);
555           aCellTypesArray->InsertNextValue(VTK_POLYGON);
556           aScalars->SetValue(aNbCells, aValue.myWarp);
557           aNbCells++;
558         }
559       }
560       aCellTypesArray->SetNumberOfTuples(aNbCells);
561       aScalars->SetNumberOfTuples(aNbCells);
562
563       vtkIdTypeArray* aCellLocationsArray = vtkIdTypeArray::New();
564       aCellLocationsArray->SetNumberOfComponents(1);
565       aCellLocationsArray->SetNumberOfTuples(aNbCells);
566
567       aConnectivity->InitTraversal();
568       vtkIdType const* pts(nullptr);
569       for (vtkIdType idType = 0, npts; aConnectivity->GetNextCell(npts, pts); idType++)
570         aCellLocationsArray->SetValue(idType, aConnectivity->GetTraversalLocation(npts));
571
572       aDataSet->SetCells(aCellTypesArray, aCellLocationsArray, aConnectivity);
573       SetUnstructuredGrid(aDataSet);
574
575       aDataSet->GetCellData()->SetScalars(aScalars);
576       aScalars->Delete();
577
578       theLookupTable->SetRange(aScalars->GetRange());
579       theLookupTable->Build();
580
581       myMergeFilter->SetScalarsData(aDataSet);
582       aDataSet->Delete();
583     }
584   }
585   GetMapper()->SetScalarVisibility(anIsInitialized);
586   theScalarBarActor->SetVisibility(anIsInitialized);
587 }
588
589 void
590 SMESH_DeviceActor
591 ::SetExtControlMode(SMESH::Controls::FunctorPtr theFunctor)
592 {
593   myExtractUnstructuredGrid->ClearRegisteredCells();
594   myExtractUnstructuredGrid->ClearRegisteredCellsWithType();
595   myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::ePassAll);
596   myVisualObj->UpdateFunctor(theFunctor);
597
598   using namespace SMESH::Controls;
599   Predicate* aPredicate = 0;
600   if (( aPredicate =  dynamic_cast<FreeBorders          *>(theFunctor.get())) ||
601       ( aPredicate =  dynamic_cast<FreeFaces            *>(theFunctor.get())) ||
602       ( aPredicate =  dynamic_cast<BareBorderVolume     *>(theFunctor.get())) ||
603       ( aPredicate =  dynamic_cast<BareBorderFace       *>(theFunctor.get())) ||
604       ( aPredicate =  dynamic_cast<OverConstrainedVolume*>(theFunctor.get())) ||
605       ( aPredicate =  dynamic_cast<CoincidentElements1D *>(theFunctor.get())) ||
606       ( aPredicate =  dynamic_cast<CoincidentElements2D *>(theFunctor.get())) ||
607       ( aPredicate =  dynamic_cast<CoincidentElements3D *>(theFunctor.get())) ||
608       ( aPredicate =  dynamic_cast<OverConstrainedFace  *>(theFunctor.get())))
609   {
610     myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::eAdding);
611     vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
612     vtkIdType aNbCells = aGrid->GetNumberOfCells();
613     for( vtkIdType i = 0; i < aNbCells; i++ ){
614       vtkIdType anObjId = myVisualObj->GetElemObjId(i);
615       if(aPredicate->IsSatisfy(anObjId))
616         myExtractUnstructuredGrid->RegisterCell(i);
617     }
618     if(!myExtractUnstructuredGrid->IsCellsRegistered())
619       myExtractUnstructuredGrid->RegisterCell(-1);
620     SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
621   }
622   else if(FreeEdges* aFreeEdges = dynamic_cast<FreeEdges*>(theFunctor.get()))
623   {
624     SMESH::Controls::FreeEdges::TBorders aBorders;
625     aFreeEdges->GetBoreders(aBorders);
626     vtkUnstructuredGrid* aDataSet = vtkUnstructuredGrid::New();
627     vtkUnstructuredGrid* aGrid = myVisualObj->GetUnstructuredGrid();
628     aDataSet->SetPoints(aGrid->GetPoints());
629
630     vtkIdType aNbCells = aBorders.size();
631     vtkIdType aCellsSize = 3*aNbCells;
632     vtkCellArray* aConnectivity = vtkCellArray::New();
633     aConnectivity->Allocate( aCellsSize, 0 );
634     
635     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
636     aCellTypesArray->SetNumberOfComponents( 1 );
637     aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
638     
639     vtkIdList *anIdList = vtkIdList::New();
640     anIdList->SetNumberOfIds(2);
641     
642     FreeEdges::TBorders::const_iterator anIter = aBorders.begin();
643     for(; anIter != aBorders.end(); anIter++){
644       const FreeEdges::Border& aBorder = *anIter;
645       vtkIdType aNode[2] = {
646         myVisualObj->GetNodeVTKId(aBorder.myPntId[0]),
647         myVisualObj->GetNodeVTKId(aBorder.myPntId[1])
648       };
649       //cout<<"aNode = "<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<endl;
650       if(aNode[0] >= 0 && aNode[1] >= 0){
651         anIdList->SetId( 0, aNode[0] );
652         anIdList->SetId( 1, aNode[1] );
653         aConnectivity->InsertNextCell( anIdList );
654         aCellTypesArray->InsertNextValue( VTK_LINE );
655       }
656     }
657     
658     vtkIdTypeArray* aCellLocationsArray = vtkIdTypeArray::New();
659     aCellLocationsArray->SetNumberOfComponents( 1 );
660     aCellLocationsArray->SetNumberOfTuples( aNbCells );
661     
662     aConnectivity->InitTraversal();
663     vtkIdType const *pts(nullptr);
664     for( vtkIdType idType = 0, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
665       aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
666     
667     aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity );
668
669     SetUnstructuredGrid(aDataSet);
670     aDataSet->Delete();
671   }
672   else if (( aPredicate = dynamic_cast<FreeNodes      *>(theFunctor.get())) ||
673            ( aPredicate = dynamic_cast<CoincidentNodes*>(theFunctor.get())))
674   {
675     myExtractUnstructuredGrid->SetModeOfChanging(VTKViewer_ExtractUnstructuredGrid::eAdding);
676     vtkIdType aNbNodes = FromSmIdType<vtkIdType>(myVisualObj->GetNbEntities(SMDSAbs_Node));
677     for( vtkIdType i = 0; i < aNbNodes; i++ ){
678       vtkIdType anObjId = myVisualObj->GetNodeObjId(i);
679       if(aPredicate->IsSatisfy(anObjId))
680         myExtractUnstructuredGrid->RegisterCell(i);
681     }
682     if(!myExtractUnstructuredGrid->IsCellsRegistered())
683       myExtractUnstructuredGrid->RegisterCell(-1);
684     SetUnstructuredGrid(myVisualObj->GetUnstructuredGrid());
685   }
686 }
687
688
689
690
691 vtkMTimeType
692 SMESH_DeviceActor
693 ::GetMTime()
694 {
695   // cout << "DA " << this
696   //      << " GF " << myGeomFilter;
697   // if ( this->Property )
698   //   cout << " P " << this->Property->GetMTime();
699   // if ( this->BackfaceProperty != NULL )
700   //   cout << " BP " << BackfaceProperty->GetMTime();
701   // if ( this->Texture != NULL )
702   //   cout << " T " << this->Texture->GetMTime();
703   // cout << " U " << this->GetUserTransformMatrixMTime()
704   //      << " M " << this->MTime.GetMTime() << endl;
705
706   // cout << "DA " << this
707   //      << " GF " << myGeomFilter
708   //      << " " << this->Superclass::GetMTime()
709   //      << " " << myExtractGeometry->GetMTime()
710   //      << " " << myExtractUnstructuredGrid->GetMTime()
711   //      << " " << myMergeFilter->GetMTime()
712   //      << " " << myGeomFilter->GetMTime()
713   //      << " " << myTransformFilter->GetMTime()
714   //      << " " << myFaceOrientationFilter->GetMTime() << endl;
715
716   vtkMTimeType mTime = this->Superclass::GetMTime();
717   mTime = max(mTime,myExtractGeometry->GetMTime());
718   mTime = max(mTime,myExtractUnstructuredGrid->GetMTime());
719   mTime = max(mTime,myMergeFilter->GetMTime());
720   mTime = max(mTime,myGeomFilter->GetMTime());
721   mTime = max(mTime,myTransformFilter->GetMTime());
722   mTime = max(mTime,myFaceOrientationFilter->GetMTime());
723   return mTime;
724 }
725
726
727 void
728 SMESH_DeviceActor
729 ::SetTransform(VTKViewer_Transform* theTransform)
730 {
731   myTransformFilter->SetTransform(theTransform);
732 }
733
734
735 void
736 SMESH_DeviceActor
737 ::SetShrink() 
738 {
739   if ( !myIsShrinkable ) return;
740   if ( vtkAlgorithmOutput* aDataSet = myPassFilter[ 0 ]->GetOutputPort() )
741   {
742     myShrinkFilter->SetInputConnection( aDataSet );
743     myPassFilter[ 1 ]->SetInputConnection( myShrinkFilter->GetOutputPort() );
744     myIsShrunk = true;
745   }
746 }
747
748 void
749 SMESH_DeviceActor
750 ::UnShrink() 
751 {
752   if ( !myIsShrunk ) return;
753   if ( vtkAlgorithmOutput* aDataSet = myPassFilter[ 0 ]->GetOutputPort() )
754   {    
755     myPassFilter[ 1 ]->SetInputConnection( aDataSet );
756     myPassFilter[ 1 ]->Modified();
757     myIsShrunk = false;
758     Modified();
759   }
760 }
761
762
763 void
764 SMESH_DeviceActor
765 ::SetFacesOriented(bool theIsFacesOriented) 
766 {
767   if ( vtkAlgorithmOutput* aDataSet = myTransformFilter->GetOutputPort() )
768   {
769     myIsFacesOriented = theIsFacesOriented;
770     if( theIsFacesOriented )
771       myFaceOrientationFilter->SetInputConnection( aDataSet );
772     UpdateFaceOrientation();
773   }
774 }
775
776 void
777 SMESH_DeviceActor
778 ::SetFacesOrientationColor(double r,double g,double b)
779 {
780   myFaceOrientation->GetProperty()->SetColor( r, g, b );
781 }
782
783 void
784 SMESH_DeviceActor
785 ::GetFacesOrientationColor(double& r,double& g,double& b)
786 {
787   myFaceOrientation->GetProperty()->GetColor( r, g, b );
788 }
789
790 void
791 SMESH_DeviceActor
792 ::SetFacesOrientationScale(double theScale)
793 {
794   myFaceOrientationFilter->SetOrientationScale( theScale );
795 }
796
797 double
798 SMESH_DeviceActor
799 ::GetFacesOrientationScale()
800 {
801   return myFaceOrientationFilter->GetOrientationScale();
802 }
803
804 void
805 SMESH_DeviceActor
806 ::SetFacesOrientation3DVectors(bool theState)
807 {
808   myFaceOrientationFilter->Set3dVectors( theState );
809 }
810
811 bool
812 SMESH_DeviceActor
813 ::GetFacesOrientation3DVectors()
814 {
815   return myFaceOrientationFilter->Get3dVectors();
816 }
817
818 void
819 SMESH_DeviceActor
820 ::UpdateFaceOrientation()
821 {
822   bool aShowFaceOrientation = myIsFacesOriented;
823   aShowFaceOrientation &= vtkLODActor::GetVisibility(); //GetVisibility(); -- avoid calling GetUnstructuredGrid()  
824   aShowFaceOrientation &= ( myRepresentation != ePoint );
825   myFaceOrientation->SetVisibility(aShowFaceOrientation);
826 }
827
828
829 void
830 SMESH_DeviceActor
831 ::SetRepresentation(EReperesent theMode)
832 {
833   if ( myRepresentation == theMode )
834     return;
835   switch(theMode){
836   case ePoint:
837     myGeomFilter->SetInside(true);
838     myGeomFilter->SetWireframeMode(false);
839     GetProperty()->SetRepresentation(0);
840     break;
841   case eWireframe:
842     myGeomFilter->SetInside(false);
843     myGeomFilter->SetWireframeMode(true);
844     GetProperty()->SetRepresentation(theMode);
845     break;
846   case eInsideframe:
847     myGeomFilter->SetInside(true);
848     myGeomFilter->SetWireframeMode(true);
849     GetProperty()->SetRepresentation(1);
850     break;
851   case eSurface:
852     myGeomFilter->SetInside(false);
853     myGeomFilter->SetWireframeMode(false);
854     GetProperty()->SetRepresentation(theMode);
855   case eNoneRepr:
856     return;
857   }
858   SetMarkerEnabled(theMode == ePoint);
859   myRepresentation = theMode;
860   UpdateFaceOrientation();
861   GetProperty()->Modified();
862   myMapper->Modified();
863   Modified();
864 }
865
866
867 void
868 SMESH_DeviceActor
869 ::SetVisibility(int theMode)
870 {
871   if(( theMode ) &&
872      ( !myExtractUnstructuredGrid->GetInput() || 
873        GetUnstructuredGrid()->GetNumberOfCells()))
874   {
875     vtkLODActor::SetVisibility(theMode);
876   }else{
877     vtkLODActor::SetVisibility(false);
878   }
879   UpdateFaceOrientation();
880 }
881
882
883 int
884 SMESH_DeviceActor
885 ::GetVisibility()
886 {
887   int visibi = vtkLODActor::GetVisibility();
888   if(visibi && !GetUnstructuredGrid()->GetNumberOfCells()){
889     vtkLODActor::SetVisibility(false);
890     visibi = 0;
891   }
892   return visibi;
893 }
894
895
896 void
897 SMESH_DeviceActor
898 ::AddToRender(vtkRenderer* theRenderer)
899 {
900   theRenderer->AddActor(this);
901   theRenderer->AddActor(myFaceOrientation);
902 }
903
904 void
905 SMESH_DeviceActor
906 ::RemoveFromRender(vtkRenderer* theRenderer)
907 {
908   theRenderer->RemoveActor(this);
909   theRenderer->RemoveActor(myFaceOrientation);
910 }
911
912
913 vtkIdType
914 SMESH_DeviceActor
915 ::GetNodeObjId(vtkIdType theVtkID)
916 {
917   vtkIdType anID = theVtkID;
918
919   if(IsImplicitFunctionUsed())
920     anID = myExtractGeometry->GetNodeObjId(theVtkID);
921
922   vtkIdType aRetID = myVisualObj->GetNodeObjId(anID);
923   MESSAGE("GetNodeObjId - theVtkID = "<<theVtkID<<"; anID = "<<anID<<"; aRetID = "<<aRetID);
924   return aRetID;
925 }
926
927 double* 
928 SMESH_DeviceActor
929 ::GetNodeCoord(vtkIdType theObjID)
930 {
931   vtkDataSet* aDataSet = myMergeFilter->GetOutput();
932   vtkIdType anID = myVisualObj->GetNodeVTKId(theObjID);
933   double* aCoord = (anID >=0 && anID < aDataSet->GetNumberOfPoints()) ? aDataSet->GetPoint(anID) : NULL;
934   MESSAGE("GetNodeCoord - theObjID = "<<theObjID<<"; anID = "<<anID);
935   return aCoord;
936 }
937
938 vtkIdType
939 SMESH_DeviceActor
940 ::GetNodeVtkId(vtkIdType theObjID) 
941 {
942   return myVisualObj->GetNodeVTKId(theObjID);
943 }
944
945 vtkIdType
946 SMESH_DeviceActor
947 ::GetElemObjId(vtkIdType theVtkID)
948 {
949   vtkIdType anId = myGeomFilter->GetElemObjId(theVtkID);
950   if(anId < 0) 
951     return -1;
952
953   vtkIdType anId2 = anId;
954   if(IsImplicitFunctionUsed())
955     anId2 = myExtractGeometry->GetElemObjId(anId);
956   if(anId2 < 0) 
957     return -1;
958
959   vtkIdType anId3 = myExtractUnstructuredGrid->GetInputId(anId2);
960   if(anId3 < 0) 
961     return -1;
962
963   vtkIdType aRetID = myVisualObj->GetElemObjId(anId3);
964
965   MESSAGE("GetElemObjId - theVtkID = "<<theVtkID<<"; anId2 = "<<anId2<<"; anId3 = "<<anId3<<"; aRetID = "<<aRetID);
966   return aRetID;
967 }
968
969 vtkCell* 
970 SMESH_DeviceActor
971 ::GetElemCell(vtkIdType theObjID)
972 {
973   vtkDataSet* aDataSet = myVisualObj->GetUnstructuredGrid();
974   vtkIdType aGridID = myVisualObj->GetElemVTKId(theObjID);
975   vtkCell* aCell = (aGridID >= 0 ) ? aDataSet->GetCell(aGridID) : NULL;
976
977   MESSAGE("GetElemCell - theObjID = "<<theObjID<<"; aGridID = "<<aGridID);
978   return aCell;
979 }
980
981
982 double 
983 SMESH_DeviceActor
984 ::GetShrinkFactor()
985 {
986   return myShrinkFilter->GetShrinkFactor();
987 }
988
989 void
990 SMESH_DeviceActor
991 ::SetShrinkFactor(double theValue)
992 {
993   theValue = theValue > 0.1? theValue: 0.8;
994   myShrinkFilter->SetShrinkFactor(theValue);
995   Modified();
996 }
997
998
999 void
1000 SMESH_DeviceActor
1001 ::SetHighlited(bool theIsHighlited)
1002 {
1003   if ( myIsHighlited == theIsHighlited )
1004     return;
1005   myIsHighlited = theIsHighlited;
1006   Modified();
1007 }
1008
1009 void
1010 SMESH_DeviceActor
1011 ::Render(vtkRenderer *ren, vtkMapper* m)
1012 {
1013   int aResolveCoincidentTopology = vtkMapper::GetResolveCoincidentTopology();
1014   double aStoredFactor, aStoredUnit; 
1015   vtkMapper::GetResolveCoincidentTopologyPolygonOffsetParameters(aStoredFactor,aStoredUnit);
1016
1017   vtkMapper::SetResolveCoincidentTopologyToPolygonOffset();
1018   double aFactor = myPolygonOffsetFactor, aUnits = myPolygonOffsetUnits;
1019   if(myIsHighlited){
1020     static double EPS = .01;
1021     aUnits *= (1.0-EPS);
1022   }
1023   vtkMapper::SetResolveCoincidentTopologyPolygonOffsetParameters(aFactor,aUnits);
1024   vtkLODActor::Render(ren,m);
1025
1026   vtkMapper::SetResolveCoincidentTopologyPolygonOffsetParameters(aStoredFactor,aStoredUnit);
1027   vtkMapper::SetResolveCoincidentTopology(aResolveCoincidentTopology);
1028 }
1029
1030
1031 void
1032 SMESH_DeviceActor
1033 ::SetPolygonOffsetParameters(double factor, 
1034                              double units)
1035 {
1036   myPolygonOffsetFactor = factor;
1037   myPolygonOffsetUnits = units;
1038 }
1039
1040 /*!
1041  * On/Off representation 2D quadratic element as arked polygon
1042  */
1043 void SMESH_DeviceActor::SetQuadraticArcMode(bool theFlag){
1044   myGeomFilter->SetQuadraticArcMode(theFlag);
1045 }
1046
1047 /*!
1048  * Return true if 2D quadratic element displayed as arked polygon
1049  */
1050 bool SMESH_DeviceActor::GetQuadraticArcMode(){
1051   return myGeomFilter->GetQuadraticArcMode();
1052 }
1053 /*!
1054  * Set Max angle for representation 2D quadratic element as arked polygon
1055  */
1056 void SMESH_DeviceActor::SetQuadraticArcAngle(double theMaxAngle){
1057   myGeomFilter->SetQuadraticArcAngle(theMaxAngle);
1058 }
1059
1060 /*!
1061  * Return Max angle of the representation 2D quadratic element as arked polygon
1062  */
1063 double SMESH_DeviceActor::GetQuadraticArcAngle(){
1064   return myGeomFilter->GetQuadraticArcAngle();
1065 }
1066
1067 /*!
1068  * Set point marker enabled
1069  * \param theMarkerEnabled flag to enable/disable point marker
1070  */
1071 void SMESH_DeviceActor::SetMarkerEnabled( bool theMarkerEnabled )
1072 {
1073   myMapper->SetMarkerEnabled( theMarkerEnabled );
1074 }
1075
1076 /*!
1077  * Set point marker enabled
1078  * \param theBallEnabled flag to enable/disable ball drawing
1079  */
1080 void SMESH_DeviceActor::SetBallEnabled( bool theBallEnabled ) {
1081   myMapper->SetBallEnabled( theBallEnabled );
1082 }
1083
1084 /*!
1085  * Set point marker scale factor
1086  * \param theBallScale double value which specifies a scale factor of ball element
1087  */
1088 void SMESH_DeviceActor::SetBallScale( double theBallScale )
1089 {
1090   myMapper->SetBallScale( theBallScale );
1091   myMapper->Modified();
1092 }
1093
1094 /*!
1095  * Set standard point marker
1096  * \param theMarkerType type of the marker
1097  */
1098 void SMESH_DeviceActor::SetMarkerStd( VTK::MarkerType theMarkerType, VTK::MarkerScale theMarkerScale )
1099 {
1100   myMapper->SetMarkerStd( theMarkerType, theMarkerScale );
1101 }
1102
1103 /*!
1104  * Set custom point marker
1105  * \param theMarkerId id of the marker texture
1106  * \param theMarkerTexture marker texture
1107  */
1108 void SMESH_DeviceActor::SetMarkerTexture( int theMarkerId, VTK::MarkerTexture theMarkerTexture )
1109 {
1110   myMapper->SetMarkerTexture( theMarkerId, theMarkerTexture );
1111 }
1112
1113 /*!
1114  * Get type of the point marker
1115  * \return type of the point marker
1116  */
1117 VTK::MarkerType SMESH_DeviceActor::GetMarkerType()
1118 {
1119   return myMapper->GetMarkerType();
1120 }
1121
1122 /*!
1123   Get scale of the point marker
1124   \return scale of the point marker
1125 */
1126 VTK::MarkerScale SMESH_DeviceActor::GetMarkerScale()
1127 {
1128   return myMapper->GetMarkerScale();
1129 }
1130
1131 /*!
1132  * Get texture identifier of the point marker
1133  * \return texture identifier of the point marker
1134  */
1135 int SMESH_DeviceActor::GetMarkerTexture()
1136 {
1137   return myMapper->GetMarkerTexture();
1138 }
1139
1140 /*!
1141  * Get scale factor of ball element
1142  * \return scale factor of ball element
1143  */
1144 double SMESH_DeviceActor::GetBallScale()
1145 {
1146   return myMapper->GetBallScale();
1147 }
1148
1149 void SMESH_DeviceActor::SetCoincident3DAllowed(bool theFlag) {
1150   myGeomFilter->SetAppendCoincident3D(theFlag);
1151 }
1152
1153 bool SMESH_DeviceActor::IsCoincident3DAllowed() const {
1154   return myGeomFilter->GetAppendCoincident3D();
1155 }