Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkSceneHandler.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//
27//
28//
29// John Allison 5th April 2001
30// A template for a simplest possible graphics driver.
31//?? Lines or sections marked like this require specialisation for your driver.
32
33#ifndef G4VTKSCENEHANDLER_HH
34#define G4VTKSCENEHANDLER_HH
35
36
37#include "G4VSceneHandler.hh"
38
39#include <map>
40#include <vector>
41
42#include "G4VtkViewer.hh"
43#pragma GCC diagnostic push
44#pragma GCC diagnostic ignored "-Wextra-semi"
45#include "vtkSmartPointer.h"
46#include "vtkPolyData.h"
47#include "vtkLine.h"
48#include "vtkNamedColors.h"
49#include "vtkProperty.h"
50#include "vtkTextProperty.h"
51#include "vtkPolyDataMapper.h"
52#include "vtkActor.h"
53#include "vtkFollower.h"
54#include "vtkTextActor.h"
55#include "vtkBillboardTextActor3D.h"
56#include "vtkRegularPolygonSource.h"
57#include "vtkSphereSource.h"
58#include "vtkFeatureEdges.h"
59#include "vtkCleanPolyData.h"
60#include "vtkPolyDataNormals.h"
61#include "vtkGlyph2D.h"
62#include "vtkGlyph3D.h"
63#include "vtkVertexGlyphFilter.h"
64#include "vtkMatrix4x4.h"
65#include "vtkDoubleArray.h"
66#include "vtkPointData.h"
67#include "vtkTriangleFilter.h"
68#include "vtkTensorGlyph.h"
69#include "vtkTensorGlyphColor.h"
70#include "vtkScalarsToColors.h"
71#include "vtkImageData.h"
72#include "vtkVolume.h"
73#include "vtkOpenGLGPUVolumeRayCastMapper.h"
74#pragma GCC diagnostic pop
75
77 friend class G4VtkViewer;
78
79public:
80 G4VtkSceneHandler(G4VGraphicsSystem& system, const G4String& name);
81 virtual ~G4VtkSceneHandler() = default;
82
83 ////////////////////////////////////////////////////////////////
84 // Required implementation of pure virtual functions...
85
86 void AddPrimitive(const G4Polyline&);
87 void AddPrimitive(const G4Text&);
88 void AddPrimitive(const G4Circle&);
89 void AddPrimitive(const G4Square&);
90 void AddPrimitive(const G4Polyhedron&);
93
94 // Further optional AddPrimitive methods. Explicitly invoke base
95 // class methods if not otherwise defined to avoid warnings about
96 // hiding of base class methods.
97 void AddPrimitive(const G4Polymarker& polymarker)
98 {
99#ifdef G4VTKDEBUG
100 G4cout << "=================================" << G4endl;
101 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polymarker& polymarker) called." << G4endl;
102#endif
104 }
105 /*
106 void AddPrimitive(const G4Scale& scale)
107 {
108#ifdef G4VTKDEBUG
109 G4cout << "=================================" << G4endl;
110 G4cout << "G4VtkSceneHandler::AddScale(const G4Scale& scale) called." <<
111G4endl; #endif G4VSceneHandler::AddPrimitive(scale);
112 }
113 */
114
115 void AddSolid(const G4Box& box);
116 void AddCompound(const G4Mesh& mesh);
117
118 void Modified();
119 void Clear();
120
121 protected:
122 static G4int fSceneIdCount; // Counter for Vtk scene handlers.
123
124 // maps for polylines
125 std::map<std::size_t, const G4VisAttributes*> polylineVisAttributesMap;
126 std::map<std::size_t, vtkSmartPointer<vtkPoints>> polylineDataMap;
127 std::map<std::size_t, vtkSmartPointer<vtkCellArray>> polylineLineMap;
128 std::map<std::size_t, vtkSmartPointer<vtkPolyData>> polylinePolyDataMap;
129 std::map<std::size_t, vtkSmartPointer<vtkPolyDataMapper>> polylinePolyDataMapperMap;
130 std::map<std::size_t, vtkSmartPointer<vtkActor>> polylinePolyDataActorMap;
131
132 // maps for circles
133 std::map<std::size_t, const G4VisAttributes*> circleVisAttributesMap;
134 std::map<std::size_t, vtkSmartPointer<vtkPoints>> circleDataMap;
135 std::map<std::size_t, vtkSmartPointer<vtkPolyData>> circlePolyDataMap;
136 std::map<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter>> circleFilterMap;
137 std::map<std::size_t, vtkSmartPointer<vtkPolyDataMapper>> circlePolyDataMapperMap;
138 std::map<std::size_t, vtkSmartPointer<vtkActor>> circlePolyDataActorMap;
139
140 // maps for squares
141 std::map<std::size_t, const G4VisAttributes*> squareVisAttributesMap;
142 std::map<std::size_t, vtkSmartPointer<vtkPoints>> squareDataMap;
143 std::map<std::size_t, vtkSmartPointer<vtkPolyData>> squarePolyDataMap;
144 std::map<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter>> squareFilterMap;
145 std::map<std::size_t, vtkSmartPointer<vtkPolyDataMapper>> squarePolyDataMapperMap;
146 std::map<std::size_t, vtkSmartPointer<vtkActor>> squarePolyDataActorMap;
147
148 // map for polyhedra
149 std::map<std::size_t, const G4VisAttributes*> polyhedronVisAttributesMap;
150 std::map<std::size_t, vtkSmartPointer<vtkPoints>> polyhedronDataMap;
151 std::map<std::size_t, vtkSmartPointer<vtkCellArray>> polyhedronPolyMap;
152 std::map<std::size_t, vtkSmartPointer<vtkPolyData>> polyhedronPolyDataMap;
153 std::map<std::size_t, std::size_t> polyhedronPolyDataCountMap;
154
155 // map for polyhedra instances
156 std::map<std::size_t, vtkSmartPointer<vtkPoints>> instancePositionMap;
157 std::map<std::size_t, vtkSmartPointer<vtkDoubleArray>> instanceRotationMap;
158 std::map<std::size_t, vtkSmartPointer<vtkDoubleArray>> instanceColoursMap;
159 std::map<std::size_t, vtkSmartPointer<vtkPolyData>> instancePolyDataMap;
160 std::map<std::size_t, vtkSmartPointer<vtkTensorGlyphColor>> instanceTensorGlyphMap;
161 std::map<std::size_t, vtkSmartPointer<vtkActor>> instanceActorMap;
162
163
164private:
165#ifdef G4VTKDEBUG
166 void PrintThings();
167#endif
168
169};
170
171#endif
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Box.hh:56
Definition: G4Mesh.hh:48
Definition: G4Text.hh:72
virtual void AddPrimitive(const G4Polyline &)=0
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > polylinePolyDataMapperMap
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polylineLineMap
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > circleFilterMap
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > squareFilterMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > instancePolyDataMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > circlePolyDataMap
std::map< std::size_t, const G4VisAttributes * > polylineVisAttributesMap
void AddPrimitive(const G4Polymarker &polymarker)
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceColoursMap
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceRotationMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > instancePositionMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > circleDataMap
void AddPrimitiveBakedTransform(const G4Polyhedron &)
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polylinePolyDataMap
static G4int fSceneIdCount
std::map< std::size_t, vtkSmartPointer< vtkPoints > > squareDataMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polyhedronPolyDataMap
void AddPrimitive(const G4Polyline &)
std::map< std::size_t, std::size_t > polyhedronPolyDataCountMap
std::map< std::size_t, const G4VisAttributes * > squareVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > circlePolyDataMapperMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > instanceActorMap
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polyhedronPolyMap
std::map< std::size_t, const G4VisAttributes * > polyhedronVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > squarePolyDataActorMap
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > squarePolyDataMapperMap
virtual ~G4VtkSceneHandler()=default
void AddSolid(const G4Box &box)
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > squarePolyDataMap
void AddPrimitiveTensorGlyph(const G4Polyhedron &)
std::map< std::size_t, vtkSmartPointer< vtkTensorGlyphColor > > instanceTensorGlyphMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > polylinePolyDataActorMap
std::map< std::size_t, const G4VisAttributes * > circleVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > circlePolyDataActorMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polylineDataMap
void AddCompound(const G4Mesh &mesh)
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polyhedronDataMap