Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VTreeSceneHandler.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 base class for a scene handler to dump geometry hierarchy.
31// Based on a provisional G4VTreeGraphicsScene (was in modeling).
32
34
35#include "G4VSolid.hh"
37#include "G4VPhysicalVolume.hh"
38#include "G4LogicalVolume.hh"
39
41// Counter for Tree scene handlers.
42
44 const G4String& name):
45 G4VSceneHandler(system, fSceneIdCount++, name),
46 fpCurrentObjectTransformation (0)
47{}
48
50
52 G4VSceneHandler::BeginModeling(); // Required: see G4VSceneHandler.hh.
53}
54
56 fDrawnLVStore.clear();
57 G4VSceneHandler::EndModeling(); // Required: see G4VSceneHandler.hh.
58}
59
61(const G4Transform3D& objectTransformation, const G4VisAttributes& visAttribs)
62{
63 G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
64
65 G4PhysicalVolumeModel* pPVModel =
66 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
67 if (!pPVModel) return; // Not from a G4PhysicalVolumeModel.
68
69 // This call comes from a G4PhysicalVolumeModel, drawnPVPath is
70 // the path of the current drawn (non-culled) volume in terms of
71 // drawn (non-culled) ancestors. Each node is identified by a
72 // PVNodeID object, which is a physical volume and copy number. It
73 // is a vector of PVNodeIDs corresponding to the geometry hierarchy
74 // actually selected, i.e., not culled.
76 typedef std::vector<PVNodeID> PVPath;
77 const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
78 //G4int currentDepth = pPVModel->GetCurrentDepth();
79 //G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
80 //G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
81 //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
82
83 // Actually, it is enough to store the logical volume of current
84 // physical volume...
85 fDrawnLVStore.insert
86 (drawnPVPath.back().GetPhysicalVolume()->GetLogicalVolume());
87
88 // Find mother. ri points to drawn mother, if any.
89 PVPath::const_reverse_iterator ri = ++drawnPVPath.rbegin();
90 if (ri != drawnPVPath.rend()) {
91 // This volume has a mother.
92 G4LogicalVolume* drawnMotherLV =
93 ri->GetPhysicalVolume()->GetLogicalVolume();
94 if (fDrawnLVStore.find(drawnMotherLV) != fDrawnLVStore.end()) {
95 // Mother previously encountered. Add this volume to
96 // appropriate node in scene graph tree.
97 // ...
98 } else {
99 // Mother not previously encountered. Shouldn't happen, since
100 // G4PhysicalVolumeModel sends volumes as it encounters them,
101 // i.e., mothers before daughters, in its descent of the
102 // geometry tree. Error!
103 G4cerr << "ERROR: G4VTreeSceneHandler::PreAddSolid: Mother "
104 << ri->GetPhysicalVolume()->GetName()
105 << ':' << ri->GetCopyNo()
106 << " not previously encountered."
107 "\nShouldn't happen! Please report to visualization coordinator."
108 << G4endl;
109 // Continue anyway. Add to root of scene graph tree.
110 // ...
111 }
112 } else {
113 // This volume has no mother. Must be a top level un-culled
114 // volume. Add to root of scene graph tree.
115 // ...
116 }
117}
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
virtual void BeginModeling()
virtual void EndModeling()
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
G4VTreeSceneHandler(G4VGraphicsSystem &system, const G4String &name)
std::set< G4LogicalVolume * > fDrawnLVStore
virtual void EndModeling()
void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
virtual void BeginModeling()