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

#include <G4VtkPolydataPolyline2DPipeline.hh>

+ Inheritance diagram for G4VtkPolydataPolyline2DPipeline:

Public Member Functions

 G4VtkPolydataPolyline2DPipeline (G4String name, const G4VtkVisContext &vc, const G4VisAttributes *pVA)
 
 ~G4VtkPolydataPolyline2DPipeline () override=default
 
void Modified () override
 
void SetPolydata (const G4Polyline &polyline) override
 
vtkSmartPointer< vtkActor2D > GetActor2D ()
 
void Clear () override
 
virtual void SetPolydata (const G4Polyhedron &polyhedron)
 
virtual void SetPolydata (vtkPolyData *polydata)
 
- 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
 
virtual vtkSmartPointer< vtkActor > GetActor ()
 
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)
 

Protected Attributes

vtkSmartPointer< vtkPolyDataMapper2D > mapper2D
 
vtkSmartPointer< vtkActor2D > actor2D
 
- 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 39 of file G4VtkPolydataPolyline2DPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkPolydataPolyline2DPipeline()

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

Definition at line 48 of file G4VtkPolydataPolyline2DPipeline.cc.

50 : G4VtkPolydataPipeline(nameIn, vcIn)
51{
52 // Set pipeline type
53 SetTypeName(G4String("G4VtkPolydataPolyline2DPipeline"));
54
55 // Get vis attributes
56 G4Color colour = pVisAttributes->GetColour();
57 G4double opacity = colour.GetAlpha();
58 G4double lineWidth = pVisAttributes->GetLineWidth();
59
60 // make 3d actor invisible and remove from renderer
61 GetActor()->SetVisibility(0);
62 vc.fViewer->renderer->RemoveActor(GetActor());
63
64 // Setup 2D mapper
66 mapper2D->SetInputConnection(GetFinalFilter()->GetOutputPort());
68 coordinate->SetCoordinateSystemToNormalizedViewport();
69 mapper2D->SetTransformCoordinate(coordinate);
70
71 // Setup 2D actor
73 GetActor2D()->GetProperty()->SetLineWidth(lineWidth);
74 GetActor2D()->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
75 GetActor2D()->GetProperty()->SetOpacity(opacity);
76 GetActor2D()->SetMapper(mapper2D);
77 GetActor2D()->SetVisibility(1);
78 // GetActor2D()->GetPositionCoordinate()->SetCoordinateSystemToNormalizedDisplay();
79
80 // add to renderer
81 vc.fViewer->renderer->AddActor(GetActor2D());
82}
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< vtkPolyDataAlgorithm > GetFinalFilter()
virtual vtkSmartPointer< vtkActor > GetActor()
G4VtkPolydataPipeline(G4String name, const G4VtkVisContext &vc)
vtkSmartPointer< vtkPolyDataMapper2D > mapper2D
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkPolydataPolyline2DPipeline()

G4VtkPolydataPolyline2DPipeline::~G4VtkPolydataPolyline2DPipeline ( )
overridedefault

Member Function Documentation

◆ Clear()

void G4VtkPolydataPolyline2DPipeline::Clear ( )
overridevirtual

Reimplemented from G4VtkPolydataPipeline.

Definition at line 84 of file G4VtkPolydataPolyline2DPipeline.cc.

85{
86 renderer->RemoveActor(GetActor2D());
88}
vtkSmartPointer< vtkRenderer > renderer

◆ GetActor2D()

vtkSmartPointer< vtkActor2D > G4VtkPolydataPolyline2DPipeline::GetActor2D ( )
inline

Definition at line 50 of file G4VtkPolydataPolyline2DPipeline.hh.

50{ return actor2D; }

Referenced by Clear(), G4VtkPolydataPolyline2DPipeline(), and Modified().

◆ MakeHash()

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

Definition at line 43 of file G4VtkPolydataPolyline2DPipeline.cc.

44{
45 return std::hash<G4VisAttributes>{}(*pVA);
46}

Referenced by G4VtkStore::AddPrimitive().

◆ Modified()

void G4VtkPolydataPolyline2DPipeline::Modified ( )
overridevirtual

Reimplemented from G4VtkPolydataPipeline.

Definition at line 90 of file G4VtkPolydataPolyline2DPipeline.cc.

91{
93 GetActor2D()->Modified();
94}

◆ SetPolydata() [1/3]

void G4VtkPolydataPipeline::SetPolydata ( const G4Polyhedron & polyhedron)
virtual

Reimplemented from G4VtkPolydataPipeline.

Definition at line 68 of file G4VtkPolydataPipeline.cc.

159{
160 G4bool notLastFace;
161 int iVert = 0;
162 do {
163 G4Point3D vertex[4];
164 G4int edgeFlag[4];
165 G4Normal3D normals[4];
166 G4int nEdges;
167 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normals);
168
170 // loop over vertices
171 for (int i = 0; i < nEdges; i++) {
172 polydataPoints->InsertNextPoint(vertex[i].x(), vertex[i].y(), vertex[i].z());
173 poly->InsertNextId(iVert);
174 iVert++;
175 }
176 polydataCells->InsertNextCell(poly);
177
178 } while (notLastFace);
179}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
vtkSmartPointer< vtkCellArray > polydataCells
vtkSmartPointer< vtkPoints > polydataPoints
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const

◆ SetPolydata() [2/3]

void G4VtkPolydataPolyline2DPipeline::SetPolydata ( const G4Polyline & polyline)
overridevirtual

Reimplemented from G4VtkPolydataPipeline.

Definition at line 96 of file G4VtkPolydataPolyline2DPipeline.cc.

97{
98 // Data data
99 const size_t nLines = polyline.size();
100
101 for (size_t i = 0; i < nLines; ++i) {
102 auto id =
103 polydataPoints->InsertNextPoint((polyline[i].x() + 1) / 2.0, (polyline[i].y() + 1) / 2.0, 0);
104
105 if (i < nLines - 1) {
107 line->GetPointIds()->SetId(0, id);
108 line->GetPointIds()->SetId(1, id + 1);
109 polydataCells->InsertNextCell(line);
110 }
111 }
112}

◆ SetPolydata() [3/3]

void G4VtkPolydataPipeline::SetPolydata ( vtkPolyData * polydata)
virtual

Reimplemented from G4VtkPolydataPipeline.

Definition at line 70 of file G4VtkPolydataPipeline.cc.

181 {
182 polydata->DeepCopy(polydataIn);
183}
vtkSmartPointer< vtkPolyData > polydata

Member Data Documentation

◆ actor2D

vtkSmartPointer<vtkActor2D> G4VtkPolydataPolyline2DPipeline::actor2D
protected

◆ mapper2D

vtkSmartPointer<vtkPolyDataMapper2D> G4VtkPolydataPolyline2DPipeline::mapper2D
protected

Definition at line 57 of file G4VtkPolydataPolyline2DPipeline.hh.

Referenced by G4VtkPolydataPolyline2DPipeline().


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