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

#include <G4VtkTextPipeline.hh>

+ Inheritance diagram for G4VtkTextPipeline:

Public Member Functions

 G4VtkTextPipeline (const G4Text &text, const G4VtkVisContext &vc, const G4VisAttributes *pVisAttributes)
 
 ~G4VtkTextPipeline () override=default
 
void Enable () override
 
void Disable () override
 
void Print () override
 
void Modified () override
 
void Clear () override
 
virtual vtkSmartPointer< vtkBillboardTextActor3D > GetActor ()
 
virtual void SetText (const G4String &text)
 
- 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 G4Text &text, const G4VtkVisContext &vc, const G4VisAttributes *pVA)
 

Protected Attributes

vtkSmartPointer< vtkBillboardTextActor3D > 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 38 of file G4VtkTextPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkTextPipeline()

G4VtkTextPipeline::G4VtkTextPipeline ( const G4Text & text,
const G4VtkVisContext & vc,
const G4VisAttributes * pVisAttributes )

Definition at line 55 of file G4VtkTextPipeline.cc.

57 : G4VVtkPipeline(text.GetText().c_str(), "G4VtkTextPipeline", vcIn, false, vcIn.fViewer->renderer)
58{
59 G4double x = text.GetPosition().x();
60 G4double y = text.GetPosition().y();
61 G4double z = text.GetPosition().z();
62
63 G4Color colour = pVA->GetColour();
64 G4double opacity = colour.GetAlpha();
65
67 actor->SetInput(text.GetText().c_str());
68 actor->SetPosition(x, y, z);
69 actor->GetTextProperty()->SetFontSize(text.GetScreenSize());
70 actor->GetTextProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
71 actor->GetTextProperty()->SetOpacity(opacity);
72
73 vc.fViewer->renderer->AddActor(actor);
74}
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
G4String GetText() const
G4double GetScreenSize() const
G4Point3D GetPosition() const
G4VtkVisContext vc
vtkSmartPointer< vtkBillboardTextActor3D > actor
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkTextPipeline()

G4VtkTextPipeline::~G4VtkTextPipeline ( )
overridedefault

Member Function Documentation

◆ Clear()

void G4VtkTextPipeline::Clear ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 96 of file G4VtkTextPipeline.cc.

97{
98 renderer->RemoveActor(actor);
100}
virtual void Clear()
vtkSmartPointer< vtkRenderer > renderer

◆ Disable()

void G4VtkTextPipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 81 of file G4VtkTextPipeline.cc.

82{
83 actor->SetVisibility(0);
84}

◆ Enable()

void G4VtkTextPipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 76 of file G4VtkTextPipeline.cc.

77{
78 actor->SetVisibility(1);
79}

◆ GetActor()

virtual vtkSmartPointer< vtkBillboardTextActor3D > G4VtkTextPipeline::GetActor ( )
inlinevirtual

Definition at line 52 of file G4VtkTextPipeline.hh.

52{ return actor; }

◆ MakeHash()

std::size_t G4VtkTextPipeline::MakeHash ( const G4Text & text,
const G4VtkVisContext & vc,
const G4VisAttributes * pVA )
static

Definition at line 38 of file G4VtkTextPipeline.cc.

40{
41 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
42
43 G4double x = text.GetPosition().x();
44 G4double y = text.GetPosition().y();
45 G4double z = text.GetPosition().z();
46
47 std::hash_combine(hash, std::hash<G4String>{}(text.GetText()));
48 std::hash_combine(hash, std::hash<G4double>{}(x));
49 std::hash_combine(hash, std::hash<G4double>{}(y));
50 std::hash_combine(hash, std::hash<G4double>{}(z));
51
52 return hash;
53}
void hash_combine(std::size_t)

Referenced by G4VtkStore::AddPrimitive().

◆ Modified()

void G4VtkTextPipeline::Modified ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 91 of file G4VtkTextPipeline.cc.

92{
94}
virtual void Modified()

◆ Print()

void G4VtkTextPipeline::Print ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 86 of file G4VtkTextPipeline.cc.

87{
89}
virtual void Print()

◆ SetText()

void G4VtkTextPipeline::SetText ( const G4String & text)
virtual

Definition at line 102 of file G4VtkTextPipeline.cc.

103{
104 actor->SetInput(text.c_str());
105}

Member Data Documentation

◆ actor

vtkSmartPointer<vtkBillboardTextActor3D> G4VtkTextPipeline::actor
protected

Definition at line 60 of file G4VtkTextPipeline.hh.

Referenced by Clear(), Disable(), Enable(), G4VtkTextPipeline(), GetActor(), and SetText().


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