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 // Created by Stewart Boogert on 13/11/2021. 27 // 28 29 #ifndef G4VTKSTORE_HH 30 #define G4VTKSTORE_HH 31 32 #include "G4Transform3D.hh" 33 #include "G4ViewParameters.hh" 34 35 #include <vtkActor.h> 36 #include <vtkActor2D.h> 37 #include <vtkCellArray.h> 38 #include <vtkCellData.h> 39 #include <vtkCleanPolyData.h> 40 #include <vtkClipClosedSurface.h> 41 #include <vtkClipPolyData.h> 42 #include <vtkDataSetMapper.h> 43 #include <vtkDoubleArray.h> 44 #include <vtkFeatureEdges.h> 45 #include <vtkNamedColors.h> 46 #include <vtkPlane.h> 47 #include <vtkPlaneCollection.h> 48 #include <vtkPointData.h> 49 #include <vtkPoints.h> 50 #include <vtkPolyData.h> 51 #include <vtkPolyDataAlgorithm.h> 52 #include <vtkPolyDataMapper.h> 53 #include <vtkPolyDataMapper2D.h> 54 #include <vtkPolyDataNormals.h> 55 #include <vtkProperty.h> 56 #include <vtkRenderer.h> 57 #include <vtkStripper.h> 58 #include <vtkStructuredGrid.h> 59 #include <vtkTensorGlyphColor.h> 60 #include <vtkTriangleFilter.h> 61 #include <vtkVertexGlyphFilter.h> 62 #include <vtkSmartPointer.h> 63 64 #include <map> 65 #include <vector> 66 67 class G4VViewer; 68 class G4VtkViewer; 69 class G4Polyline; 70 class G4Text; 71 class G4Circle; 72 class G4Square; 73 class G4Polyhedron; 74 class G4VisAttributes; 75 class G4VtkVisContext; 76 class G4Mesh; 77 78 class G4VtkPolydataPipeline; 79 class G4VtkPolydataPolylinePipeline; 80 class G4VtkPolydataPolyline2DPipeline; 81 class G4VtkPolydataSpherePipeline; 82 class G4VtkPolydataPolygonPipeline; 83 class G4VtkPolydataInstanceTensorPipeline; 84 class G4VtkPolydataInstanceAppendPipeline; 85 class G4VtkPolydataInstanceBakePipeline; 86 class G4VtkClipCloseSurfacePipeline; 87 class G4VtkImagePipeline; 88 class G4VtkTextPipeline; 89 class G4VtkText2DPipeline; 90 class G4VtkUnstructuredGridPipeline; 91 92 template<typename> 93 class vtkSmartPointer; 94 95 class G4VtkTensorGlyphPolydataPipeline; 96 97 class G4VtkStore 98 { 99 public: 100 G4VtkStore(G4String name); 101 ~G4VtkStore(); 102 103 void Modified(); 104 void Clear(); 105 void ClearNonG4(); 106 void Print(); 107 108 void AddPrimitive(const G4Polyline& polyline, const G4VtkVisContext& vc); 109 void AddPrimitive(const G4Text& text, const G4VtkVisContext& vc); 110 void AddPrimitive(const G4Circle& circle, const G4VtkVisContext& vc); 111 void AddPrimitive(const G4Square& square, const G4VtkVisContext& vc); 112 void AddPrimitiveSeparate(const G4Polyhedron& polyhedron, const G4VtkVisContext& vc); 113 void AddPrimitiveTensorGlyph(const G4Polyhedron& polyhedron, const G4VtkVisContext& vc); 114 void AddPrimitiveAppend(const G4Polyhedron& polyhedron, const G4VtkVisContext& vc); 115 void AddPrimitiveTransformBake(const G4Polyhedron& polyhedron, const G4VtkVisContext& vc); 116 void AddCompound(const G4Mesh& mesh, const G4VtkVisContext& vc); 117 118 void UpdatePlanePipelines(G4String name, G4String type, const G4Plane3D); 119 120 void AddClipper(G4String name, const G4Plane3D& plane); 121 void UpdateClipper(G4String name, const G4Plane3D& plane); 122 void RemoveClipper(G4String name); 123 124 void AddCutter(G4String name, const G4Plane3D& plane); 125 void UpdateCutter(G4String name, const G4Plane3D& plane); 126 void RemoveCutter(G4String name); 127 128 void AddNonG4ObjectImage(const G4String& fileName, const G4VtkVisContext& vc); 129 void AddNonG4ObjectPolydata(const G4String fileName, const G4VtkVisContext& vc); 130 131 void GetBounds(G4double maxBound[6]); 132 133 void AddToRenderer(vtkRenderer* renderer); 134 135 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPolylinePipeline>>& GetPolylinePipeMap() 136 { 137 return polylinePipeMap; 138 } 139 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPolyline2DPipeline>>& GetPolyline2DPipeMap() 140 { 141 return polyline2DPipeMap; 142 } 143 std::map<std::size_t, std::shared_ptr<G4VtkPolydataSpherePipeline>>& GetPolydataSpherePipeMap() 144 { 145 return circlePipeMap; 146 } 147 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPolygonPipeline>>& 148 GetPolydataPolygonPipeMap() 149 { 150 return squarePipeMap; 151 } 152 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPipeline>>& GetSeparatePipeMap() 153 { 154 return separatePipeMap; 155 } 156 std::map<std::size_t, std::shared_ptr<G4VtkPolydataInstanceTensorPipeline>>& GetTensorPipeMap() 157 { 158 return tensorGlyphPipeMap; 159 } 160 std::map<std::size_t, std::shared_ptr<G4VtkPolydataInstanceAppendPipeline>>& GetAppendPipeMap() 161 { 162 return appendPipeMap; 163 } 164 std::map<std::size_t, std::shared_ptr<G4VtkPolydataInstanceBakePipeline>>& GetBakePipeMap() 165 { 166 return bakePipeMap; 167 } 168 169 private: 170 G4String name; 171 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPolylinePipeline>> polylinePipeMap; 172 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPolyline2DPipeline>> polyline2DPipeMap; 173 std::map<std::size_t, std::shared_ptr<G4VtkPolydataSpherePipeline>> circlePipeMap; 174 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPolygonPipeline>> squarePipeMap; 175 std::map<std::size_t, std::shared_ptr<G4VtkTextPipeline>> textPipeMap; 176 std::map<std::size_t, std::shared_ptr<G4VtkText2DPipeline>> text2DPipeMap; 177 178 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPipeline>> separatePipeMap; 179 std::map<std::size_t, std::shared_ptr<G4VtkPolydataInstanceTensorPipeline>> tensorGlyphPipeMap; 180 std::map<std::size_t, std::shared_ptr<G4VtkPolydataInstanceAppendPipeline>> appendPipeMap; 181 std::map<std::size_t, std::shared_ptr<G4VtkPolydataInstanceBakePipeline>> bakePipeMap; 182 183 std::map<G4String, std::shared_ptr<G4VtkUnstructuredGridPipeline>> ugridPipeMap; 184 185 std::map<G4String, std::shared_ptr<G4VtkImagePipeline>> imagePipeMap; 186 std::map<std::size_t, std::shared_ptr<G4VtkPolydataPipeline>> sideloadPipeMap; 187 }; 188 189 #endif // G4VTKSTORE_HH 190