Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkTextPipeline.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Created by Stewart Boogert on 05/03/2023.
27//
28
29#include "G4VtkTextPipeline.hh"
30
31#include "G4Text.hh"
32#include "G4VisAttributes.hh"
33#include "G4VtkViewer.hh"
34
35#include <vtkBillboardTextActor3D.h>
36#include <vtkTextProperty.h>
37
38std::size_t G4VtkTextPipeline::MakeHash(const G4Text& text, const G4VtkVisContext& /*vc*/,
39 const G4VisAttributes* pVA)
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}
54
56 const G4VisAttributes* pVA)
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}
75
77{
78 actor->SetVisibility(1);
79}
80
82{
83 actor->SetVisibility(0);
84}
85
90
95
97{
98 renderer->RemoveActor(actor);
100}
101
103{
104 actor->SetInput(text.c_str());
105}
double G4double
Definition G4Types.hh:83
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
virtual void Clear()
virtual void Print()
virtual void Modified()
vtkSmartPointer< vtkRenderer > renderer
G4VtkVisContext vc
const G4Colour & GetColour() const
virtual void SetText(const G4String &text)
static std::size_t MakeHash(const G4Text &text, const G4VtkVisContext &vc, const G4VisAttributes *pVA)
vtkSmartPointer< vtkBillboardTextActor3D > actor
void Clear() override
void Enable() override
void Modified() override
void Print() override
void Disable() override
G4VtkTextPipeline(const G4Text &text, const G4VtkVisContext &vc, const G4VisAttributes *pVisAttributes)
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer
void hash_combine(std::size_t)