Geant4 Cross Reference |
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 38 namespace std 39 { 40 inline void hash_combine(std::size_t) {} 41 42 template<typename T, typename... Rest> 43 inline 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 50 template<> 51 struct 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) { 61 std::hash_combine(h, c); 62 } 63 64 return h; 65 } 66 }; 67 68 template<> 69 struct 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 78 std::hash_combine(h, va.IsVisible()); 79 std::hash_combine(h, va.IsDaughtersInvisible()); 80 std::hash_combine(h, va.GetColour().GetRed()); 81 std::hash_combine(h, va.GetColour().GetGreen()); 82 std::hash_combine(h, va.GetColour().GetBlue()); 83 std::hash_combine(h, va.GetColour().GetAlpha()); 84 std::hash_combine(h, static_cast<int>(va.GetLineStyle())); 85 86 return h; 87 } 88 }; 89 90 template<> 91 struct 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 124 class G4VVtkPipeline 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 {} 134 virtual ~G4VVtkPipeline() 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; } 162 bool GetDisableParent() { return disableParent; } 163 164 void SetName(G4String nameIn) { name = nameIn; } 165 const G4String GetName() { return name; } 166 167 void SetTypeName(G4String typeNameIn) { type = typeNameIn; } 168 G4String GetTypeName() { return type; } 169 170 G4VtkVisContext& GetVtkVisContext() { return vc; } 171 172 void AddChildPipeline(G4VVtkPipeline* child) 173 { 174 childPipelines.push_back(child); 175 if (child->GetDisableParent()) { 176 Disable(); 177 } 178 } 179 G4VVtkPipeline* GetChildPipeline(G4int iPipeline) { return childPipelines[iPipeline]; } 180 G4int GetNumberOfChildPipelines() { return (G4int)childPipelines.size(); } 181 std::vector<G4VVtkPipeline*> GetChildPipelines() { return childPipelines; } 182 void ClearChildPipeline() { childPipelines.clear(); } 183 184 protected: 185 G4String name; 186 G4String type; 187 G4bool disableParent; 188 std::vector<G4VVtkPipeline*> childPipelines; 189 vtkSmartPointer<vtkRenderer> renderer; 190 G4VtkVisContext vc; 191 }; 192 193 #endif // G4VVTKPIPELINE_HH 194