Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkSceneHandler.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//
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#include "G4VtkSceneHandler.hh"
34
35#include "G4Box.hh"
36#include "G4Circle.hh"
37#include "G4LogicalVolume.hh"
39#include "G4Material.hh"
40#include "G4Mesh.hh"
42#include "G4Polyhedron.hh"
43#include "G4Polyline.hh"
44#include "G4PseudoScene.hh"
45#include "G4Square.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4Text.hh"
48#include "G4UnitsTable.hh"
50#include "G4VPhysicalVolume.hh"
51#include "G4VtkStore.hh"
52#include "G4VtkVisContext.hh"
53
54#include <vtkColorTransferFunction.h>
55#include <vtkContourValues.h>
56#include <vtkPiecewiseFunction.h>
57#include <vtkVolumeProperty.h>
58
59#include <cstdlib>
60
62// Counter for XXX scene handlers.
63
65 : G4VSceneHandler(system, fSceneIdCount++, name), polyhedronPipelineType(G4String("tensor"))
66{}
67
69{
70#ifdef G4VTKDEBUG
71 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline&)" << G4endl;
72#endif
73 auto vc = MakeDefaultVisContext();
74
76 transientStore.AddPrimitive(polyline, vc);
77 else
78 store.AddPrimitive(polyline, vc);
79}
80
82{
83#ifdef G4VTKDEBUG
84 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Text& text)" << G4endl;
85#endif
86
87 auto vc = MakeDefaultVisContext();
88
91 else
92 store.AddPrimitive(text, vc);
93}
94
96{
97#ifdef G4VTKDEBUG
98 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle)" << G4endl;
99#endif
100
101 auto vc = MakeDefaultVisContext();
103 vc.fSize = GetMarkerSize(circle, sizeType);
104
106 transientStore.AddPrimitive(circle, vc);
107 else
108 store.AddPrimitive(circle, vc);
109}
110
112{
113#ifdef G4VTKDEBUG
114 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Square& square)" << G4endl;
115#endif
116
117 auto vc = MakeDefaultVisContext();
119 vc.fSize = GetMarkerSize(square, sizeType);
120
122 transientStore.AddPrimitive(square, vc);
123 else
124 store.AddPrimitive(square, vc);
125}
126
128{
129#ifdef G4VTKDEBUG
130 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron)" << G4endl;
131#endif
132
133 auto vc = MakeDefaultVisContext();
134 auto visAtt = vc.fViewer->GetApplicableVisAttributes(polyhedron.GetVisAttributes());
135 auto colour = visAtt->GetColour();
136
137 vc.fDrawingStyle = GetDrawingStyle(visAtt);
138 vc.alpha = colour.GetAlpha();
139 vc.red = colour.GetRed();
140 vc.green = colour.GetGreen();
141 vc.blue = colour.GetBlue();
142
143 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
144 if (pPVModel != nullptr) {
145 vc.fDepth = pPVModel->GetCurrentDepth();
146 vc.fDescription = pPVModel->GetCurrentDescription();
147 }
148
150 if (polyhedronPipelineType == "tensor")
152 else if (polyhedronPipelineType == "append")
153 transientStore.AddPrimitiveAppend(polyhedron, vc);
154 else if (polyhedronPipelineType == "bake")
156 else if (polyhedronPipelineType == "separate")
157 transientStore.AddPrimitiveSeparate(polyhedron, vc);
158 }
159 else {
160 if (polyhedronPipelineType == "tensor")
161 store.AddPrimitiveTensorGlyph(polyhedron, vc);
162 else if (polyhedronPipelineType == "append")
163 store.AddPrimitiveAppend(polyhedron, vc);
164 else if (polyhedronPipelineType == "bake")
165 store.AddPrimitiveTransformBake(polyhedron, vc);
166 else if (polyhedronPipelineType == "separate")
167 store.AddPrimitiveSeparate(polyhedron, vc);
168 }
169}
170
172{
173#ifdef G4VTKDEBUG
174 G4cout << "G4VtkSceneHandler::Modified()" << G4endl;
175#endif
176
177 store.Modified();
179}
180
182{
183#ifdef G4VTKDEBUG
184 G4cout << "G4VtkSceneHandler::ClearStore()" << G4endl;
185#endif
186 store.Clear();
187}
188
190{
191#ifdef G4VTKDEBUG
192 G4cout << "G4VtkSceneHandler::ClearTransientStore()" << G4endl;
193#endif
195}
196
198{
201
202 if (fpVisAttribs != nullptr) {
204 vc.red = c.GetRed();
205 vc.green = c.GetGreen();
206 vc.blue = c.GetBlue();
207 vc.alpha = c.GetAlpha();
208 vc.fDrawingStyle = fpViewer->GetViewParameters().GetDrawingStyle();
209 }
210
211 return vc;
212}
213
215{
217
218 return;
219
220 const G4VModel* pv_model = GetModel();
221 if (pv_model == nullptr) {
222 return;
223 }
224
225 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
226 if (pPVModel == nullptr) {
227 return;
228 }
229
230 //-- debug information
231#ifdef G4VTKDEBUG
232 G4VPhysicalVolume* pv = pPVModel->GetCurrentPV();
234 G4cout << "name=" << box.GetName() << " volumeType=" << pv->VolumeType()
235 << " pvName=" << pv->GetName() << " lvName=" << lv->GetName()
236 << " multiplicity=" << pv->GetMultiplicity() << " isparametrised=" << pv->IsParameterised()
237 << " isreplicated=" << pv->IsReplicated()
238 << " parametrisation=" << pv->GetParameterisation()
239 << G4endl
240
241 G4Material* mat = pPVModel->GetCurrentMaterial();
242 G4String name = mat->GetName();
243 G4double dens = mat->GetDensity() / (g / cm3);
244 G4int copyNo = pPVModel->GetCurrentPV()->GetCopyNo();
245 G4int depth = pPVModel->GetCurrentDepth();
246 G4cout << " name : " << box.GetName() << G4endl;
247 G4cout << " copy no.: " << copyNo << G4endl;
248 G4cout << " depth : " << depth << G4endl;
249 G4cout << " density : " << dens << " [g/cm3]" << G4endl;
250 G4cout << " location: " << pPVModel->GetCurrentPV()->GetObjectTranslation() << G4endl;
251 G4cout << " Multiplicity : " << pPVModel->GetCurrentPV()->GetMultiplicity() << G4endl;
252 G4cout << " Is replicated? : " << pPVModel->GetCurrentPV()->IsReplicated() << G4endl;
253 G4cout << " Is parameterised? : " << pPVModel->GetCurrentPV()->IsParameterised() << G4endl;
254 G4cout << " top phys. vol. name : " << pPVModel->GetTopPhysicalVolume()->GetName() << G4endl;
255#endif
256}
257
259{
260#ifdef G4VTKDEBUG
261 G4cout << "G4VtkSceneHandler::AddCompound> mesh type " << mesh.GetMeshType() << " "
263#endif
264
266 {
267 auto vc = MakeDefaultVisContext();
268
270 transientStore.AddCompound(mesh, vc);
271 else
272 store.AddCompound(mesh, vc);
273 }
274 else {
276 }
277}
278
280
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Definition G4Box.hh:56
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
const G4String & GetName() const
MeshType GetMeshType() const
Definition G4Mesh.hh:75
virtual G4bool IsReplicated() const =0
virtual EVolume VolumeType() const =0
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4bool IsParameterised() const =0
G4VModel * GetModel() const
G4Transform3D fObjectTransformation
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
const G4VisAttributes * fpVisAttribs
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
virtual void AddSolid(const G4Box &)
void StandardSpecialMeshRendering(const G4Mesh &)
G4String GetName() const
const G4ViewParameters & GetViewParameters() const
SMROption GetSpecialMeshRenderingOption() const
DrawingStyle GetDrawingStyle() const
const G4Colour & GetColour() const
const G4VisAttributes * GetVisAttributes() const
void ClearTransientStore() override
void AddCompound(const G4Mesh &mesh) override
void AddPrimitive(const G4Polyline &) override
void SetPolyhedronPipeline(const G4String &str)
void ClearStore() override
void AddSolid(const G4Box &box) override
G4VtkSceneHandler(G4VGraphicsSystem &system, const G4String &name)
G4VtkVisContext MakeDefaultVisContext()
void AddPrimitiveSeparate(const G4Polyhedron &polyhedron, const G4VtkVisContext &vc)
void Modified()
Definition G4VtkStore.cc:90
void AddPrimitiveTensorGlyph(const G4Polyhedron &polyhedron, const G4VtkVisContext &vc)
void Clear()
void AddPrimitiveTransformBake(const G4Polyhedron &polyhedron, const G4VtkVisContext &vc)
void AddPrimitiveAppend(const G4Polyhedron &polyhedron, const G4VtkVisContext &vc)
void AddCompound(const G4Mesh &mesh, const G4VtkVisContext &vc)
void AddPrimitive(const G4Polyline &polyline, const G4VtkVisContext &vc)