VTKFIG  0.20.20181114
Easy VTK based in situ visualization
template<class DATA >
void vtkfig::Figure::RTBuildDomainPipeline0 ( vtkSmartPointer< vtkRenderer >  renderer)
protected

Duck typing interface allowing to handle different VTK datatypes with the same code.

Todo:
add filter to extract boundary edges

if boundary cell color set, use this one!

!!!

Definition at line 829 of file vtkfigFigure.cxx.

830  {
831  RTCalcTransform();
832 
833  auto geometry=vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
834  geometry->SetInputConnection(data_producer->GetOutputPort());
835 
836  auto transgeometry=vtkSmartPointer<vtkTransformFilter>::New();
837  transgeometry->SetInputConnection(geometry->GetOutputPort());
838  transgeometry->SetTransform(transform);
839 
840 
842  if (state.show_domain_boundary && state.spacedim==3)
843  {
844  vtkSmartPointer<vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
845  mapper->SetInputConnection(transgeometry->GetOutputPort());
846  splot=vtkSmartPointer<vtkActor>::New();
847  if (state.spacedim==3)
848  {
849  splot->GetProperty()->SetOpacity(state.domain_opacity);
850  splot->GetProperty()->SetColor(state.domain_surface_color);
851  }
852  else // TODO add filter to extract boundary edges
853  {
854  splot->GetProperty()->SetOpacity(1.0);
855  splot->GetProperty()->SetColor(0,0,0);
856  }
857 
858  splot->SetMapper(mapper);
860  }
861 
862 
863  if (state.show_domain_box&& state.spacedim==3)
864  {
865  // create outline
866  vtkSmartPointer<vtkOutlineFilter>outlinefilter = vtkSmartPointer<vtkOutlineFilter>::New();
867  outlinefilter->SetInputConnection(transgeometry->GetOutputPort());
868  vtkSmartPointer<vtkPolyDataMapper> outlineMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
869  outlineMapper->SetInputConnection(outlinefilter->GetOutputPort());
870  outline = vtkSmartPointer<vtkActor>::New();
871  outline->SetMapper(outlineMapper);
872  outline->GetProperty()->SetColor(0, 0, 0);
874  }
875 
876 
877  if (state.show_domain_axes)
878  {
879  axes=vtkSmartPointer<vtkCubeAxesActor2D>::New();
880  double unscaled_data_bounds[6];
881  unscaled_data_bounds[0]=data_bounds[0]/this->state.coordinate_scale_factor_xyz[0];
882  unscaled_data_bounds[1]=data_bounds[1]/this->state.coordinate_scale_factor_xyz[0];
883  unscaled_data_bounds[2]=data_bounds[2]/this->state.coordinate_scale_factor_xyz[1];
884  unscaled_data_bounds[3]=data_bounds[3]/this->state.coordinate_scale_factor_xyz[1];
885  unscaled_data_bounds[4]=data_bounds[4]/this->state.coordinate_scale_factor_xyz[2];
886  unscaled_data_bounds[5]=data_bounds[5]/this->state.coordinate_scale_factor_xyz[2];
887  axes->SetRanges(unscaled_data_bounds);
888  axes->SetUseRanges(1);
889  axes->SetInputConnection(transgeometry->GetOutputPort());
890  axes->GetProperty()->SetColor(0, 0, 0);
891  axes->SetFontFactor(1.5);
892  axes->SetCornerOffset(0);
893  axes->SetNumberOfLabels(3);
894  axes->SetInertia(100);
895  axes->SetLabelFormat("%6.2g");
896 
898  axes->SetCamera(renderer->GetActiveCamera());
899 
900  if (state.spacedim==2)
901  {
902  axes->SetXLabel("");
903  axes->SetYLabel("");
904  axes->ZAxisVisibilityOff();
905  }
906  else
907  {
908  axes->SetXLabel("x");
909  axes->SetYLabel("y");
910  axes->SetZLabel("z");
911  }
912  auto textprop=axes->GetAxisLabelTextProperty();
913  textprop->ItalicOff();
914  textprop->BoldOn();
915  textprop->SetFontFamilyToCourier();
916  textprop->SetColor(0,0,0);
917 
918  textprop=axes->GetAxisTitleTextProperty();
919  textprop->ItalicOff();
920  textprop->BoldOn();
921  textprop->SetFontFamilyToCourier();
922  textprop->SetColor(0,0,0);
923 
925  }
926  }
vtkSmartPointer< vtkActor > outline
Definition: vtkfigFigure.h:237
void RTAddActor2D(vtkSmartPointer< vtkActor2D > prop)
Add vtk Actor to renderer showing figure.
vtkSmartPointer< vtkActor > splot
Definition: vtkfigFigure.h:238
vtkSmartPointer< vtkCubeAxesActor2D > axes
Definition: vtkfigFigure.h:236
vtkSmartPointer< vtkTrivialProducer > data_producer
Data producer for grid dataset.
Definition: vtkfigFigure.h:247
vtkSmartPointer< vtkTransform > transform
Definition: vtkfigFigure.h:200
void RTCalcTransform()
Calculate transformation to unit cube This shall be applied to all data.
struct vtkfig::Figure::@0 state
figure state
void RTAddActor(vtkSmartPointer< vtkActor > prop)
Add vtk Actor to renderer showing figure.
double data_bounds[6]
Definition: vtkfigFigure.h:529

+ Here is the call graph for this function: