Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkPolydataPolygonPipeline Class Reference

#include <G4VtkPolydataPolygonPipeline.hh>

+ Inheritance diagram for G4VtkPolydataPolygonPipeline:

Public Member Functions

 G4VtkPolydataPolygonPipeline (G4String name, const G4VtkVisContext &vc, const G4VisAttributes *pVA)
 
 ~G4VtkPolydataPolygonPipeline () override=default
 
- Public Member Functions inherited from G4VtkPolydataPipeline
 G4VtkPolydataPipeline (G4String name, const G4VtkVisContext &vc)
 
 ~G4VtkPolydataPipeline () override=default
 
void AddFilter (vtkSmartPointer< vtkPolyDataAlgorithm > f)
 
vtkSmartPointer< vtkPolyDataAlgorithm > GetFilter (G4int iFilter)
 
G4int GetNumberOfFilters ()
 
vtkSmartPointer< vtkPolyDataAlgorithm > GetFinalFilter ()
 
void Enable () override
 
void Disable () override
 
void Print () override
 
void Modified () override
 
void Clear () override
 
virtual vtkSmartPointer< vtkActor > GetActor ()
 
virtual void SetPolydata (const G4Polyhedron &polyhedron)
 
virtual void SetPolydata (const G4Polyline &polyline)
 
virtual void SetPolydata (vtkPolyData *polydata)
 
virtual void SetPolydataData (const G4Point3D &p)
 
virtual void SetPolydataData (double x, double y, double z)
 
virtual void SetActorTransform (G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22)
 
virtual void SetActorColour (G4double r, G4double g, G4double b, G4double a)
 
virtual vtkSmartPointer< vtkPolyData > GetPolydata ()
 
virtual G4doubleGetBounds ()
 
- Public Member Functions inherited from G4VVtkPipeline
 G4VVtkPipeline ()
 
 G4VVtkPipeline (G4String nameIn, G4String typeIn, const G4VtkVisContext &vcIn, G4bool disableParentIn=false, vtkSmartPointer< vtkRenderer > rendererIn=nullptr)
 
virtual ~G4VVtkPipeline ()
 
void SetDisableParent (G4bool disableParentIn)
 
bool GetDisableParent ()
 
void SetName (G4String nameIn)
 
const G4String GetName ()
 
void SetTypeName (G4String typeNameIn)
 
G4String GetTypeName ()
 
G4VtkVisContextGetVtkVisContext ()
 
void AddChildPipeline (G4VVtkPipeline *child)
 
G4VVtkPipelineGetChildPipeline (G4int iPipeline)
 
G4int GetNumberOfChildPipelines ()
 
std::vector< G4VVtkPipeline * > GetChildPipelines ()
 
void ClearChildPipeline ()
 

Static Public Member Functions

static std::size_t MakeHash (const G4VisAttributes *va)
 
- Static Public Member Functions inherited from G4VtkPolydataPipeline
static std::size_t MakeHash (const G4Polyhedron &p, const G4VtkVisContext &vc)
 

Additional Inherited Members

- Protected Attributes inherited from G4VtkPolydataPipeline
vtkSmartPointer< vtkPoints > polydataPoints
 
vtkSmartPointer< vtkCellArray > polydataCells
 
vtkSmartPointer< vtkPolyData > polydata
 
std::vector< vtkSmartPointer< vtkPolyDataAlgorithm > > filters
 
vtkSmartPointer< vtkPolyDataMapper > mapper
 
vtkSmartPointer< vtkActor > actor
 
- Protected Attributes inherited from G4VVtkPipeline
G4String name
 
G4String type
 
G4bool disableParent
 
std::vector< G4VVtkPipeline * > childPipelines
 
vtkSmartPointer< vtkRenderer > renderer
 
G4VtkVisContext vc
 

Detailed Description

Definition at line 35 of file G4VtkPolydataPolygonPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkPolydataPolygonPipeline()

G4VtkPolydataPolygonPipeline::G4VtkPolydataPolygonPipeline ( G4String name,
const G4VtkVisContext & vc,
const G4VisAttributes * pVA )

Definition at line 43 of file G4VtkPolydataPolygonPipeline.cc.

46 : G4VtkPolydataPipeline(nameIn, vcIn)
47{
48 // Set pipeline type
49 SetTypeName(G4String("G4VtkPolydataSpherePipeline"));
50
51 // Get vis attributes
52 G4Color colour = pVisAttributes->GetColour();
53 G4double opacity = colour.GetAlpha();
54 G4double lineWidth = pVisAttributes->GetLineWidth();
55
56 // Draw in world coordinates.
59 polygonSource->SetNumberOfSides(4);
60 polygonSource->SetRadius(vc.fSize);
61
62 // vertex glyph filter
64 glyphFilter->SetInputData(polydata);
65 AddFilter(glyphFilter);
66
67 // Setup actor and mapper
68 mapper->SetInputConnection(glyphFilter->GetOutputPort());
69
70 GetActor()->GetProperty()->SetLineWidth(lineWidth);
71 GetActor()->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
72 GetActor()->GetProperty()->SetOpacity(opacity);
73 GetActor()->SetVisibility(1);
74
75 GetActor()->GetProperty()->SetRenderPointsAsSpheres(true);
76 GetActor()->GetProperty()->SetPointSize(vc.fSize * 5);
77
78 vc.fViewer->renderer->AddActor(GetActor());
79}
double G4double
Definition G4Types.hh:83
static G4bool GetColour(const G4String &key, G4Colour &result)
Definition G4Colour.cc:155
G4double GetBlue() const
Definition G4Colour.hh:154
G4double GetAlpha() const
Definition G4Colour.hh:155
G4double GetRed() const
Definition G4Colour.hh:152
G4double GetGreen() const
Definition G4Colour.hh:153
void SetTypeName(G4String typeNameIn)
G4VtkVisContext vc
vtkSmartPointer< vtkPolyData > polydata
vtkSmartPointer< vtkPolyDataMapper > mapper
virtual vtkSmartPointer< vtkActor > GetActor()
void AddFilter(vtkSmartPointer< vtkPolyDataAlgorithm > f)
G4VtkPolydataPipeline(G4String name, const G4VtkVisContext &vc)
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkPolydataPolygonPipeline()

G4VtkPolydataPolygonPipeline::~G4VtkPolydataPolygonPipeline ( )
overridedefault

Member Function Documentation

◆ MakeHash()

std::size_t G4VtkPolydataPolygonPipeline::MakeHash ( const G4VisAttributes * va)
static

Definition at line 38 of file G4VtkPolydataPolygonPipeline.cc.

39{
40 return std::hash<G4VisAttributes>{}(*pVA);
41}

The documentation for this class was generated from the following files: