VTKFIG  0.20.20181114
Easy VTK based in situ visualization
vtkfigGridView.cxx
Go to the documentation of this file.
1 #include <vtkSliderRepresentation2D.h>
2 #include <vtkProperty2D.h>
3 #include <vtkAlgorithmOutput.h>
4 #include <vtkTextProperty.h>
5 #include <vtkRectilinearGrid.h>
6 #include <vtkGeometryFilter.h>
7 #include <vtkRectilinearGridGeometryFilter.h>
8 #include <vtkUnstructuredGrid.h>
9 #include <vtkUnstructuredGridGeometryFilter.h>
10 #include <vtkPlane.h>
11 #include <vtkCutter.h>
12 #include <vtkImplicitBoolean.h>
13 #include <vtkOutlineFilter.h>
14 #include <vtkCubeAxesActor2D.h>
15 #include <vtkAppendPolyData.h>
16 #include <vtkAssignAttribute.h>
17 #include <vtkCamera.h>
18 #include <vtkTextActor.h>
19 #include <vtkCornerAnnotation.h>
20 #include <vtkCoordinate.h>
21 #include <vtkExtractEdges.h>
22 #include <vtkTransformPolyDataFilter.h>
23 #include <vtkTransformFilter.h>
24 
25 
26 
27 #include "vtkfigGridView.h"
28 #include "config.h"
29 
30 
31 namespace vtkfig
32 {
33 
34 
38  {
39  cutgeometry=vtkSmartPointer<vtkExtractGeometry>::New();
40  cutgeometry->SetImplicitFunction(planeZ);
41  bcutgeometry=vtkSmartPointer<vtkExtractGeometry>::New();
42  bcutgeometry->SetImplicitFunction(planeZ);
43  }
44 
45  int GridView::RTProcessPlaneMove(const std::string plane,int idim, int dx, int dy, bool & edit,
46  vtkSmartPointer<vtkPlane> planeXYZ )
47  {
48  if (edit)
49  {
50  double trans_origin[3];
51  planeXYZ->GetOrigin(trans_origin);
52  double planepos=trans_origin[idim];
53  double delta=trans_bounds[2*idim+1]-trans_bounds[2*idim+0];
54  planepos+=(0.01+1.0e-8)*((double)dx)*delta;
55  planepos=std::min(planepos,trans_bounds[2*idim+1]+1.0e-10*delta);
56  planepos=std::max(planepos,trans_bounds[2*idim+0]-1.0e-10*delta);
57  trans_origin[idim]=planepos;
58  planeXYZ->SetOrigin(trans_origin);
59  RTShowPlanePos(planeXYZ,plane,idim);
60  planeXYZ->Modified();
61  cutgeometry->Modified();
62  bcutgeometry->Modified();
63  return 1;
64  }
65  return 0;
66  }
67 
68 
69  void GridView::RTProcessMove(int dx, int dy)
70  {
71  if (RTProcessPlaneMove("x",0,dx,dy,edit.x_plane, planeX)) return;
72  if (RTProcessPlaneMove("y",1,dx,dy,edit.y_plane, planeY)) return;
73  if (RTProcessPlaneMove("z",2,dx,dy,edit.z_plane, planeZ)) return;
74  }
75 
76 
77  void GridView::RTShowPlanePos(vtkSmartPointer<vtkPlane> planeXYZ, const std::string plane, int idim)
78  {
79  double trans_origin[3];
80  planeXYZ->GetOrigin(trans_origin);
81  double planepos=data_bounds[2*idim]+ (data_bounds[2*idim+1]-data_bounds[2*idim])*(trans_origin[idim]-trans_bounds[2*idim])/(trans_bounds[2*idim+1]-trans_bounds[2*idim]);
82  RTMessage("plane_" + plane+"=" + std::to_string(planepos));
83  }
84 
85 
87  const std::string plane,
88  int idim,
89  const std::string key,
90  bool & edit,
91  vtkSmartPointer<vtkPlane> planeXYZ)
92  {
93 
94  if (!edit && key==plane)
95  {
96  edit=true;
97  RTShowPlanePos(planeXYZ,plane,idim);
98  cutgeometry->SetImplicitFunction(planeXYZ);
99  cutgeometry->Modified();
100  bcutgeometry->SetImplicitFunction(planeXYZ);
101  bcutgeometry->Modified();
102  return 1;
103  }
104 
105  if (edit&& key=="Escape")
106  {
107  edit=false;
108  return 1;
109  }
110 
111  return 0;
112  }
113 
114 
115  void GridView::RTProcessKey(const std::string key)
116  {
117 
118  if (key=="slash")
119  {
120  cutgeometry->SetExtractInside(!cutgeometry->GetExtractInside());
121  cutgeometry->Modified();
122  bcutgeometry->SetExtractInside(!bcutgeometry->GetExtractInside());
123  bcutgeometry->Modified();
124  return;
125  }
126 
127  if (key=="A")
128  {
129  if (axes)
130  {
131  int vis=axes->GetVisibility();
132  vis=!vis;
133  axes->SetVisibility(vis);
134  if (outline)
135  outline->SetVisibility(vis);
136  }
137  return;
138  }
139 
140  if (key=="O")
141  {
142  if (splot)
143  {
144  int vis=splot->GetVisibility();
145  vis=!vis;
146  splot->SetVisibility(vis);
147  }
148  return;
149  }
150 
151 
152 
153  if (key=="C" && cellplot)
154  {
155  int vis=cellplot->GetVisibility();
156  vis=!vis;
157  cellplot->SetVisibility(vis);
158  if (edgeplot)
159  edgeplot->SetVisibility(vis);
160 
161  if (cbar)
162  cbar->SetVisibility(vis);
163  return;
164  }
165 
166  if (key=="B" && bcellplot)
167  {
168  int vis=bcellplot->GetVisibility();
169  vis=!vis;
170  bcellplot->SetVisibility(vis);
171  if (bedgeplot)
172  bedgeplot->SetVisibility(vis);
173  if (bcbar)
174  bcbar->SetVisibility(vis);
175  return;
176  }
177 
178 
179  if (
180  (edit.x_plane||edit.y_plane|| edit.z_plane)
181  &&(key=="x"|| key=="y"|| key=="z")
182  )
183  {
184  edit.x_plane=false;
185  edit.y_plane=false;
186  edit.z_plane=false;
187  }
188 
189  if (RTProcessPlaneKey("x",0,key,edit.x_plane, planeX)) return;
190  if (RTProcessPlaneKey("y",1,key,edit.y_plane, planeY)) return;
191  if (RTProcessPlaneKey("z",2,key,edit.z_plane, planeZ)) return;
192 
193 
194  }
195 
196 
199 
200  template <class DATA, class FILTER>
202  {
203  RTCalcTransform();
204  {
205  auto cr=vtkDoubleArray::SafeDownCast(DATA::SafeDownCast(data_producer->GetOutputDataObject(0))->GetCellData()->GetAbstractArray("cellregions"));
206 
207  auto scalar = vtkSmartPointer<vtkAssignAttribute>::New();
208  if (cr)
209  {
210  scalar->Assign("cellregions",vtkDataSetAttributes::SCALARS,vtkAssignAttribute::CELL_DATA);
211  scalar->SetInputConnection(data_producer->GetOutputPort());
212  }
213 
214  auto geometry=vtkSmartPointer<FILTER>::New();
215  if (cr)
216  geometry->SetInputConnection(scalar->GetOutputPort());
217  else
218  geometry->SetInputConnection(data_producer->GetOutputPort());
219 
220 
221 
222  auto transgeometry=vtkSmartPointer<vtkTransformPolyDataFilter>::New();
223  transgeometry->SetInputConnection(geometry->GetOutputPort());
224  transgeometry->SetTransform(transform);
225 
226  auto cells = vtkSmartPointer<vtkPolyDataMapper>::New();
227  cells->SetInputConnection(transgeometry->GetOutputPort());
228 
229  if (cr)
230  {
231  cells->UseLookupTableScalarRangeOn();
232  cells->SetLookupTable(cell_lut);
233  }
234  else
235  cells->ScalarVisibilityOff();
236 
237 #ifdef VTK_HAS_MAPPER_IMMEDIATE_RENDERING_ON
238  cells->ImmediateModeRenderingOn();
239 #endif
240  cellplot = vtkSmartPointer<vtkActor>::New();
241  cellplot->SetMapper(cells);
242  if (!cr)
243  cellplot->GetProperty()->SetColor(0.9,0.9,0.9);
245 
246 
247 
249  {
250  // Extract edges is slow for large datasets so it is better to
251  // plot cells twice: once with, once without wireframe.
252 
253  // auto edges= vtkSmartPointer<vtkExtractEdges>::New();
254  // edges->SetInputConnection(transgeometry->GetOutputPort());
255  // auto emapper = vtkSmartPointer<vtkPolyDataMapper>::New();
256  // emapper->SetInputConnection(edges->GetOutputPort());
257  // emapper->ScalarVisibilityOff();
258 
259  // edgeplot = vtkSmartPointer<vtkActor>::New();
260  // edgeplot->GetProperty()->SetColor(0,0,0);
261  // edgeplot->SetMapper(emapper);
262  // Figure::RTAddActor(edgeplot);
263 
264 
265  auto celledges = vtkSmartPointer<vtkPolyDataMapper>::New();
266  celledges->SetInputConnection(transgeometry->GetOutputPort());
267  celledges->ScalarVisibilityOff();
268 #ifdef VTK_HAS_MAPPER_IMMEDIATE_RENDERING_ON
269  celledges->ImmediateModeRenderingOn();
270 #endif
271  auto celledgeplot = vtkSmartPointer<vtkActor>::New();
272  celledgeplot->SetMapper(celledges);
273  celledgeplot->GetProperty()->SetColor(0,0,0);
274  celledgeplot->GetProperty()->SetRepresentationToWireframe();
275  Figure::RTAddActor(celledgeplot);
276  }
277 
278  if (cr && state.show_grid_colorbar)
279  {
280  cbar=BuildColorBar(cells);
281  cbar->SetTitle("C ");
282  cbar->SetLabelFormat(" %-2.0f ");
283  double range[2];
284  cell_lut->GetTableRange(range);
285  cbar->SetNumberOfLabels((int)(range[1]-range[0]+1));
287  }
288  }
289 
290 
291  auto boundary_data=vtkDataSet::SafeDownCast(boundary_data_producer->GetOutputDataObject(0));
292  if (boundary_data)
293  {
294 
295  auto bcr=vtkDoubleArray::SafeDownCast(boundary_data->GetCellData()->GetAbstractArray("boundarycellregions"));
296  if (bcr)
297  {
298  auto bscalar = vtkSmartPointer<vtkAssignAttribute>::New();
299  bscalar->Assign("boundarycellregions",vtkDataSetAttributes::SCALARS,vtkAssignAttribute::CELL_DATA);
300  bscalar->SetInputConnection(boundary_data_producer->GetOutputPort());
301 
302  auto bgeometry=vtkSmartPointer<FILTER>::New();
303  bgeometry->SetInputConnection(bscalar->GetOutputPort());
304 
305  auto transbgeometry=vtkSmartPointer<vtkTransformPolyDataFilter>::New();
306  transbgeometry->SetInputConnection(bgeometry->GetOutputPort());
307  transbgeometry->SetTransform(transform);
308 
309  auto bcells = vtkSmartPointer<vtkPolyDataMapper>::New();
310  bcells->SetInputConnection(transbgeometry->GetOutputPort());
311 
312  bcells->UseLookupTableScalarRangeOn();
313  bcells->SetLookupTable(bface_lut);
314  bcells->ScalarVisibilityOn();
315 #ifdef VTK_HAS_MAPPER_IMMEDIATE_RENDERING_ON
316  bcells->ImmediateModeRenderingOn();
317 #endif
318  bcellplot = vtkSmartPointer<vtkActor>::New();
319  bcellplot->SetMapper(bcells);
320  bcellplot->GetProperty()->SetLineWidth(5);
322 
323  if( state.show_grid_colorbar)
324  {
325  bcbar=BuildColorBar(bcells,1);
326  bcbar->SetTitle("B ");
327  bcbar->SetLabelFormat(" %-2.0f ");
328  double brange[2];
329  bface_lut->GetTableRange(brange);
330  bcbar->SetNumberOfLabels((int)(brange[1]-brange[0]+1));
332  }
333  }
334  }
335  }
336 
337 
340 
341  template <class DATA,class FILTER>
343  {
344  RTCalcTransform();
345  {
346  double range[2];
347  auto cr=vtkDoubleArray::SafeDownCast(DATA::SafeDownCast(data_producer->GetOutputDataObject(0))->GetCellData()->GetAbstractArray("cellregions"));
348 
349  auto scalar = vtkSmartPointer<vtkAssignAttribute>::New();
350  if (cr)
351  {
352  scalar->Assign("cellregions",vtkDataSetAttributes::SCALARS,vtkAssignAttribute::CELL_DATA);
353  scalar->SetInputConnection(data_producer->GetOutputPort());
354  cr->GetRange(range);
355 
356  cell_lut->SetTableRange(range[0],range[1]);
357  cell_lut->Modified();
358  }
359 
360 
361  planeX->SetOrigin(trans_center);
362  planeY->SetOrigin(trans_center);
363  planeZ->SetOrigin(trans_center);
364 
365  auto transgeometry=vtkSmartPointer<vtkTransformFilter>::New();
366  transgeometry->SetTransform(transform);
367 
368 
369  if (cr)
370  transgeometry->SetInputConnection(scalar->GetOutputPort());
371  else
372  transgeometry->SetInputConnection(data_producer->GetOutputPort());
373 
374 
375 
376 
377  cutgeometry->SetInputConnection(transgeometry->GetOutputPort());
378 
379  auto cutpolydata=vtkSmartPointer<vtkGeometryFilter>::New();
380  cutpolydata->SetInputConnection(cutgeometry->GetOutputPort());
381 
382  auto cells = vtkSmartPointer<vtkPolyDataMapper>::New();
383  cells->SetInputConnection(cutpolydata->GetOutputPort());
384 
385  if (cr)
386  {
387  cells->UseLookupTableScalarRangeOn();
388  cells->SetLookupTable(cell_lut);
389  }
390  else
391  cells->ScalarVisibilityOff();
392 #ifdef VTK_HAS_MAPPER_IMMEDIATE_RENDERING_ON
393  cells->ImmediateModeRenderingOn();
394 #endif
395  cellplot = vtkSmartPointer<vtkActor>::New();
396  cellplot->SetMapper(cells);
397  if (!cr)
398  cellplot->GetProperty()->SetColor(0.9,0.9,0.9);
400 
401 
402  // Extract edges is slow, plot cell polydata with wireframe instead
403  // auto edges= vtkSmartPointer<vtkExtractEdges>::New();
404  // edges->SetInputConnection(cutgeometry->GetOutputPort());
405  // auto emapper = vtkSmartPointer<vtkPolyDataMapper>::New();
406  // emapper->SetInputConnection(edges->GetOutputPort());
407  // emapper->ScalarVisibilityOff();
408  // edgeplot = vtkSmartPointer<vtkActor>::New();
409  // edgeplot->GetProperty()->SetColor(0,0,0);
410  // edgeplot->SetMapper(emapper);
411  // Figure::RTAddActor(edgeplot);
412 
413 
414  auto celledges = vtkSmartPointer<vtkPolyDataMapper>::New();
415  celledges->SetInputConnection(cutpolydata->GetOutputPort());
416  celledges->ScalarVisibilityOff();
417 #ifdef VTK_HAS_MAPPER_IMMEDIATE_RENDERING_ON
418  celledges->ImmediateModeRenderingOn();
419 #endif
420  auto celledgeplot = vtkSmartPointer<vtkActor>::New();
421  celledgeplot->SetMapper(celledges);
422  celledgeplot->GetProperty()->SetColor(0,0,0);
423  celledgeplot->GetProperty()->SetRepresentationToWireframe();
424  Figure::RTAddActor(celledgeplot);
425 
426 
427  if ( cr)
428  {
429  cbar=BuildColorBar(cells);
430  cbar->SetLabelFormat(" %-2.0f ");
431  cbar->SetNumberOfLabels((int)(range[1]-range[0]+1));
432  cbar->SetTitle("c ");
434  }
435  }
436 
437  auto boundary_data=vtkDataSet::SafeDownCast(boundary_data_producer->GetOutputDataObject(0));
438  if (boundary_data)
439  {
440  double brange[2];
441  auto bcr=vtkDoubleArray::SafeDownCast(boundary_data->GetCellData()->GetAbstractArray("boundarycellregions"));
442 
443  if (bcr)
444  {
445  auto bscalar = vtkSmartPointer<vtkAssignAttribute>::New();
446  bscalar->Assign("boundarycellregions",vtkDataSetAttributes::SCALARS,vtkAssignAttribute::CELL_DATA);
447  bscalar->SetInputDataObject(boundary_data);
448  bcr->GetRange(brange);
449 
450  bface_lut->SetTableRange(brange[0],brange[1]);
451  bface_lut->Modified();
452  auto bgeometry=vtkSmartPointer<FILTER>::New();
453  bgeometry->SetInputConnection(bscalar->GetOutputPort());
454 
455  auto transbgeometry=vtkSmartPointer<vtkTransformPolyDataFilter>::New();
456  transbgeometry->SetInputConnection(bgeometry->GetOutputPort());
457  transbgeometry->SetTransform(transform);
458 
459 
460  bcutgeometry->SetInputConnection(transbgeometry->GetOutputPort());
461 
462  auto bcutpolydata=vtkSmartPointer<vtkGeometryFilter>::New();
463  bcutpolydata->SetInputConnection(bcutgeometry->GetOutputPort());
464 
465 
466 
467  auto bcells = vtkSmartPointer<vtkPolyDataMapper>::New();
468  bcells->SetInputConnection(bcutpolydata->GetOutputPort());
469 
470  bcells->UseLookupTableScalarRangeOn();
471  bcells->SetLookupTable(bface_lut);
472  bcells->ScalarVisibilityOn();
473 #ifdef VTK_HAS_MAPPER_IMMEDIATE_RENDERING_ON
474  bcells->ImmediateModeRenderingOn();
475 #endif
476  bcellplot = vtkSmartPointer<vtkActor>::New();
477  bcellplot->SetMapper(bcells);
479 
480 
481  auto bedges= vtkSmartPointer<vtkExtractEdges>::New();
482  bedges->SetInputConnection(bcutgeometry->GetOutputPort());
483  auto bemapper = vtkSmartPointer<vtkPolyDataMapper>::New();
484  bemapper->SetInputConnection(bedges->GetOutputPort());
485  bemapper->ScalarVisibilityOff();
486  bedgeplot = vtkSmartPointer<vtkActor>::New();
487  bedgeplot->GetProperty()->SetColor(0,0,0);
488  bedgeplot->SetMapper(bemapper);
490 
491 
492  bcbar=BuildColorBar(bcells,1);
493  bcbar->SetTitle("b ");
494  bcbar->SetLabelFormat(" %-2.0f ");
495  bcbar->SetNumberOfLabels((int)(brange[1]-brange[0]+1));
497  }
498  }
499  }
500 
501 
505  {
507  {
508  if (state.spacedim==2)
509  this->RTBuildVTKPipeline2D<vtkUnstructuredGrid,vtkGeometryFilter>();
510  else
511  this->RTBuildVTKPipeline3D<vtkUnstructuredGrid,vtkGeometryFilter>();
512  return;
513  }
514  else if (state.datatype==DataSet::DataType::RectilinearGrid)
515  {
516  if (state.spacedim==2)
517  this->RTBuildVTKPipeline2D<vtkRectilinearGrid,vtkRectilinearGridGeometryFilter>();
518  else
519  this->RTBuildVTKPipeline3D<vtkRectilinearGrid,vtkRectilinearGridGeometryFilter>();
520  return;
521  }
522  }
523 }
vtkSmartPointer< vtkActor > outline
Definition: vtkfigFigure.h:237
double trans_bounds[6]
Definition: vtkfigFigure.h:532
void RTMessage(std::string msg)
Show string in message field.
void RTBuildVTKPipeline2D()
2D Filter
void RTAddActor2D(vtkSmartPointer< vtkActor2D > prop)
Add vtk Actor to renderer showing figure.
vtkSmartPointer< vtkScalarBarActor > bcbar
Definition: vtkfigFigure.h:294
vtkSmartPointer< vtkActor > bcellplot
int RTProcessPlaneKey(const std::string plane, int idim, const std::string key, bool &edit, vtkSmartPointer< vtkPlane > planeXYZ)
Base class for all figures.
Definition: vtkfigFigure.h:52
void RTShowPlanePos(vtkSmartPointer< vtkPlane > planeXYZ, const std::string plane, int idim)
vtkSmartPointer< vtkActor > splot
Definition: vtkfigFigure.h:238
vtkSmartPointer< vtkCubeAxesActor2D > axes
Definition: vtkfigFigure.h:236
vtkSmartPointer< vtkExtractGeometry > cutgeometry
struct vtkfig::Figure::@1 edit
edit state
vtkSmartPointer< vtkLookupTable > bface_lut
Definition: vtkfigFigure.h:276
vtkSmartPointer< vtkActor > bedgeplot
vtkSmartPointer< vtkTrivialProducer > data_producer
Data producer for grid dataset.
Definition: vtkfigFigure.h:247
vtkSmartPointer< vtkActor > cellplot
vtkSmartPointer< vtkTransform > transform
Definition: vtkfigFigure.h:200
void RTCalcTransform()
Calculate transformation to unit cube This shall be applied to all data.
vtkSmartPointer< vtkActor > edgeplot
void RTProcessKey(const std::string key) override final
Process keyboard and mouse move events.
struct vtkfig::Figure::@0 state
figure state
vtkSmartPointer< vtkExtractGeometry > bcutgeometry
void RTBuildVTKPipeline3D()
3D Filter
vtkSmartPointer< vtkPlane > planeY
Definition: vtkfigFigure.h:209
int RTProcessPlaneMove(const std::string plane, int idim, int dx, int dy, bool &edit, vtkSmartPointer< vtkPlane > planeXYZ)
vtkSmartPointer< vtkTrivialProducer > boundary_data_producer
Data producer for boundary grid dataset.
Definition: vtkfigFigure.h:250
vtkSmartPointer< vtkLookupTable > cell_lut
Definition: vtkfigFigure.h:275
void RTAddActor(vtkSmartPointer< vtkActor > prop)
Add vtk Actor to renderer showing figure.
vtkSmartPointer< vtkPlane > planeX
Plane equations for plane sections.
Definition: vtkfigFigure.h:208
vtkSmartPointer< vtkScalarBarActor > cbar
Definition: vtkfigFigure.h:293
vtkSmartPointer< vtkScalarBarActor > BuildColorBar(vtkSmartPointer< vtkPolyDataMapper > mapper, int irank=0)
Definition: vtkfigTools.cxx:88
vtkSmartPointer< vtkPlane > planeZ
Definition: vtkfigFigure.h:210
double trans_center[3]
Definition: vtkfigFigure.h:533
GridView()
Constructor.
void RTProcessMove(int dx, int dy) override final
double data_bounds[6]
Definition: vtkfigFigure.h:529
void RTBuildVTKPipeline() override final
Generic access to filter.