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

#include <G4VtkText2DPipeline.hh>

+ Inheritance diagram for G4VtkText2DPipeline:

Public Member Functions

 G4VtkText2DPipeline (const G4Text &text, const G4VtkVisContext &vc, const G4VisAttributes *pVisAttributes)
 
 ~G4VtkText2DPipeline () override=default
 
void Enable () override
 
void Disable () override
 
void Print () override
 
void Modified () override
 
void Clear () override
 
virtual vtkSmartPointer< vtkTextActor > 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< vtkTextActor > 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 G4VtkText2DPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkText2DPipeline()

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

Definition at line 53 of file G4VtkText2DPipeline.cc.

55 : G4VVtkPipeline(text.GetText().c_str(), "G4VtkText2DPipeline", vcIn, false, vcIn.fViewer->renderer)
56{
57 G4double x = text.GetPosition().x();
58 G4double y = text.GetPosition().y();
59
60 G4Color colour = pVA->GetColour();
61 G4double opacity = colour.GetAlpha();
62
64 actor->SetInput(text.GetText().c_str());
65 actor->GetPositionCoordinate()->SetCoordinateSystemToNormalizedViewport();
66 // actor->SetTextScaleModeToViewport();
67 actor->SetPosition((x + 1.) / 2.0, (y + 1.) / 2.);
68 actor->GetTextProperty()->SetFontSize(text.GetScreenSize());
69 actor->GetTextProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
70 actor->GetTextProperty()->SetOpacity(opacity);
71
72 switch (text.GetLayout()) {
74 actor->GetTextProperty()->SetJustificationToLeft();
75 break;
77 actor->GetTextProperty()->SetJustificationToCentered();
78 break;
80 actor->GetTextProperty()->SetJustificationToRight();
81 break;
82 }
83
84 vc.fViewer->renderer->AddActor(actor);
85}
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
Layout GetLayout() const
G4String GetText() const
@ centre
Definition G4Text.hh:76
@ right
Definition G4Text.hh:76
@ left
Definition G4Text.hh:76
G4double GetScreenSize() const
G4Point3D GetPosition() const
G4VtkVisContext vc
vtkSmartPointer< vtkTextActor > actor
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkText2DPipeline()

G4VtkText2DPipeline::~G4VtkText2DPipeline ( )
overridedefault

Member Function Documentation

◆ Clear()

void G4VtkText2DPipeline::Clear ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 107 of file G4VtkText2DPipeline.cc.

108{
109 renderer->RemoveViewProp(actor);
110
112}
virtual void Clear()
vtkSmartPointer< vtkRenderer > renderer

◆ Disable()

void G4VtkText2DPipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 92 of file G4VtkText2DPipeline.cc.

93{
94 actor->SetVisibility(0);
95}

◆ Enable()

void G4VtkText2DPipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 87 of file G4VtkText2DPipeline.cc.

88{
89 actor->SetVisibility(1);
90}

◆ GetActor()

virtual vtkSmartPointer< vtkTextActor > G4VtkText2DPipeline::GetActor ( )
inlinevirtual

Definition at line 52 of file G4VtkText2DPipeline.hh.

52{ return actor; }

◆ MakeHash()

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

Definition at line 38 of file G4VtkText2DPipeline.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
46 std::hash_combine(hash, std::hash<G4String>{}(text.GetText()));
47 std::hash_combine(hash, std::hash<G4double>{}(x));
48 std::hash_combine(hash, std::hash<G4double>{}(y));
49
50 return hash;
51}
void hash_combine(std::size_t)

◆ Modified()

void G4VtkText2DPipeline::Modified ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 102 of file G4VtkText2DPipeline.cc.

103{
105}
virtual void Modified()

◆ Print()

void G4VtkText2DPipeline::Print ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 97 of file G4VtkText2DPipeline.cc.

98{
100}
virtual void Print()

◆ SetText()

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

Definition at line 114 of file G4VtkText2DPipeline.cc.

115{
116 actor->SetInput(text.c_str());
117}

Member Data Documentation

◆ actor

vtkSmartPointer<vtkTextActor> G4VtkText2DPipeline::actor
protected

Definition at line 60 of file G4VtkText2DPipeline.hh.

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


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