Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VVtkPipeline.hh
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#ifndef G4VVTKPIPELINE_HH
27#define G4VVTKPIPELINE_HH
28
29#include "G4Polyhedron.hh"
30#include "G4String.hh"
31#include "G4Types.hh"
32#include "G4VisAttributes.hh"
33#include "G4VtkVisContext.hh"
34
35#include <vtkSmartPointer.h>
36#include <vtkRenderer.h>
37
38namespace std
39{
40inline void hash_combine(std::size_t) {}
41
42template<typename T, typename... Rest>
43inline void hash_combine(std::size_t& seed, const T& v, Rest... rest)
44{
45 std::hash<T> hasher;
46 seed ^= hasher(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
47 std::hash_combine(seed, rest...);
48}
49
50template<>
51struct hash<G4String>
52{
53 std::size_t operator()(const G4String& strng) const
54 {
55 using std::hash;
56 using std::size_t;
57
58 std::size_t h = 0;
59
60 for (char const& c : strng) {
62 }
63
64 return h;
65 }
66};
67
68template<>
69struct hash<G4VisAttributes>
70{
71 std::size_t operator()(const G4VisAttributes& va) const
72 {
73 using std::hash;
74 using std::size_t;
75
76 std::size_t h = 0;
77
84 std::hash_combine(h, static_cast<int>(va.GetLineStyle()));
85
86 return h;
87 }
88};
89
90template<>
91struct hash<G4Polyhedron>
92{
93 std::size_t operator()(const G4Polyhedron& ph) const
94 {
95 using std::hash;
96 using std::size_t;
97
98 G4bool notLastFace;
99 G4Point3D vertex[4];
100 G4int edgeFlag[4];
101 G4Normal3D normals[4];
102 G4int nEdges;
103
104 std::size_t h = 0;
105
106 do {
107 notLastFace = ph.GetNextFacet(nEdges, vertex, edgeFlag, normals);
108
109 for (int i = 0; i < nEdges; i++) {
110 std::size_t hx = std::hash<double>()(vertex[i].x());
111 std::size_t hy = std::hash<double>()(vertex[i].y());
112 std::size_t hz = std::hash<double>()(vertex[i].z());
113 std::hash_combine(h, hx);
114 std::hash_combine(h, hy);
115 std::hash_combine(h, hz);
116 }
117 } while (notLastFace);
118
119 return h;
120 }
121};
122} // namespace std
123
125{
126 public:
127 G4VVtkPipeline() : name("none"), type("G4VVtkPipeline"), disableParent(false), renderer(nullptr)
128 {}
129 G4VVtkPipeline(G4String nameIn, G4String typeIn, const G4VtkVisContext& vcIn,
130 G4bool disableParentIn = false,
131 vtkSmartPointer<vtkRenderer> rendererIn = nullptr)
132 : name(nameIn), type(typeIn), disableParent(disableParentIn), renderer(rendererIn), vc(vcIn)
133 {}
135 {
136 for (auto i : childPipelines)
137 delete i;
138 }
139
140 virtual void Enable() = 0;
141 virtual void Disable() = 0;
142
143 virtual void Print()
144 {
145 G4cout << "children (";
146 for (auto i : childPipelines)
147 G4cout << i->GetName() << "-" << i->GetTypeName() << ",";
148 G4cout << ")" << G4endl;
149 };
150 virtual void Modified()
151 {
152 for (auto i : childPipelines)
153 i->Modified();
154 }
155 virtual void Clear()
156 {
157 for (auto i : childPipelines)
158 i->Clear();
159 }
160
161 void SetDisableParent(G4bool disableParentIn) { disableParent = disableParentIn; }
163
164 void SetName(G4String nameIn) { name = nameIn; }
165 const G4String GetName() { return name; }
166
167 void SetTypeName(G4String typeNameIn) { type = typeNameIn; }
169
171
173 {
174 childPipelines.push_back(child);
175 if (child->GetDisableParent()) {
176 Disable();
177 }
178 }
179 G4VVtkPipeline* GetChildPipeline(G4int iPipeline) { return childPipelines[iPipeline]; }
181 std::vector<G4VVtkPipeline*> GetChildPipelines() { return childPipelines; }
183
184 protected:
188 std::vector<G4VVtkPipeline*> childPipelines;
191};
192
193#endif // G4VVTKPIPELINE_HH
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
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
virtual void Enable()=0
virtual void Clear()
virtual void Print()
std::vector< G4VVtkPipeline * > childPipelines
G4String GetTypeName()
void SetTypeName(G4String typeNameIn)
virtual void Modified()
G4VVtkPipeline * GetChildPipeline(G4int iPipeline)
G4int GetNumberOfChildPipelines()
virtual void Disable()=0
vtkSmartPointer< vtkRenderer > renderer
void AddChildPipeline(G4VVtkPipeline *child)
G4VtkVisContext vc
const G4String GetName()
void SetDisableParent(G4bool disableParentIn)
void SetName(G4String nameIn)
std::vector< G4VVtkPipeline * > GetChildPipelines()
G4VtkVisContext & GetVtkVisContext()
G4VVtkPipeline(G4String nameIn, G4String typeIn, const G4VtkVisContext &vcIn, G4bool disableParentIn=false, vtkSmartPointer< vtkRenderer > rendererIn=nullptr)
virtual ~G4VVtkPipeline()
G4bool IsDaughtersInvisible() const
LineStyle GetLineStyle() const
const G4Colour & GetColour() const
G4bool IsVisible() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const
void hash_combine(std::size_t)
std::size_t operator()(const G4Polyhedron &ph) const
std::size_t operator()(const G4String &strng) const
std::size_t operator()(const G4VisAttributes &va) const