Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicalVolumeModel.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// $Id$
28//
29//
30// John Allison 31st December 1997.
31//
32// Class Description:
33//
34// Model for physical volumes. It describes a physical volume and its
35// daughters to any desired depth. Note: the "requested depth" is
36// specified in the modeling parameters; enum {UNLIMITED = -1}.
37//
38// For access to base class information, e.g., modeling parameters,
39// use GetModelingParameters() inherited from G4VModel. See Class
40// Description of the base class G4VModel.
41//
42// G4PhysicalVolumeModel assumes the modeling parameters have been set
43// up with meaningful information - default vis attributes and culling
44// policy in particular.
45//
46// The volumes are unpacked and sent to the scene handler. They are,
47// in effect, "touchables" as defined by G4TouchableHistory. A
48// touchable is defined not only by its physical volume name and copy
49// number but also by its position in the geometry hierarchy.
50//
51// It is guaranteed that touchables are presented to the scene handler
52// in top-down hierarchy order, i.e., ancesters first, mothers before
53// daughters, so the scene handler can be assured that, if it is
54// building its own scene graph tree, a mother, if any, will have
55// already been encountered and there will already be a node in place
56// on which to hang the current volume. But be aware that the
57// visibility and culling policy might mean that some touchables are
58// not passed to the scene handler so the drawn tree might have
59// missing layers. GetFullPVPath allows you to know the full
60// hierarchy.
61
62#ifndef G4PHYSICALVOLUMEMODEL_HH
63#define G4PHYSICALVOLUMEMODEL_HH
64
65#include "G4VModel.hh"
66#include "G4VTouchable.hh"
67
68#include "G4Transform3D.hh"
69#include "G4Plane3D.hh"
70#include <iostream>
71#include <vector>
72#include <map>
73
75class G4LogicalVolume;
76class G4VSolid;
77class G4Material;
78class G4VisAttributes;
79class G4AttDef;
80class G4AttValue;
81
83
84public: // With description
85
86 enum {UNLIMITED = -1};
87
89
91 public:
93 (G4VPhysicalVolume* pPV = 0,
94 G4int iCopyNo = 0,
95 G4int depth = 0,
96 const G4Transform3D& transform = G4Transform3D(),
97 G4bool drawn = true):
98 fpPV(pPV),
99 fCopyNo(iCopyNo),
100 fNonCulledDepth(depth),
101 fTransform(transform),
102 fDrawn(drawn) {}
103 G4VPhysicalVolume* GetPhysicalVolume() const {return fpPV;}
104 G4int GetCopyNo() const {return fCopyNo;}
105 G4int GetNonCulledDepth() const {return fNonCulledDepth;}
106 const G4Transform3D& GetTransform() const {return fTransform;}
107 G4bool GetDrawn() const {return fDrawn;}
108 void SetDrawn(G4bool drawn) {fDrawn = drawn;}
109 G4bool operator< (const G4PhysicalVolumeNodeID& right) const;
110 private:
111 G4VPhysicalVolume* fpPV;
112 G4int fCopyNo;
113 G4int fNonCulledDepth;
114 G4Transform3D fTransform;
115 G4bool fDrawn;
116 };
117 // Nested class for identifying physical volume nodes.
118
120 public:
122 (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath);
123 const G4ThreeVector& GetTranslation(G4int depth) const;
124 const G4RotationMatrix* GetRotation(G4int depth) const;
125 G4VPhysicalVolume* GetVolume(G4int depth) const;
126 G4VSolid* GetSolid(G4int depth) const;
127 G4int GetReplicaNumber(G4int depth) const;
128 G4int GetHistoryDepth() const {return fFullPVPath.size();}
129 private:
130 const std::vector<G4PhysicalVolumeNodeID>& fFullPVPath;
131 };
132 // Nested class for handling nested parameterisations.
133
135 (G4VPhysicalVolume* = 0,
136 G4int requestedDepth = UNLIMITED,
137 const G4Transform3D& modelTransformation = G4Transform3D(),
138 const G4ModelingParameters* = 0,
139 G4bool useFullExtent = false);
140
141 virtual ~G4PhysicalVolumeModel ();
142
144 // The main task of a model is to describe itself to the graphics scene
145 // handler (a object which inherits G4VSceneHandler, which inherits
146 // G4VGraphicsScene).
147
149 // A description which depends on the current state of the model.
150
151 G4String GetCurrentTag () const;
152 // A tag which depends on the current state of the model.
153
155
157
159 {return fpClippingSolid;}
160
162 // Current depth of geom. hierarchy.
163
165 // Current physical volume.
166
168 // Current logical volume.
169
171 // Current material.
172
173 const std::vector<G4PhysicalVolumeNodeID>& GetFullPVPath() const
174 {return fFullPVPath;}
175 // Vector of physical volume node identifiers for the current
176 // touchable. It is its path in the geometry hierarchy, similar to
177 // the concept of "touchable history" available from the navigator
178 // during tracking.
179
180 const std::vector<G4PhysicalVolumeNodeID>& GetDrawnPVPath() const
181 {return fDrawnPVPath;}
182 // Path of the current drawn (non-culled) touchable in terms of
183 // drawn (non-culled) ancesters. It is a vector of physical volume
184 // node identifiers corresponding to the geometry hierarchy actually
185 // selected, i.e., not culled.
186
187 const std::map<G4String,G4AttDef>* GetAttDefs() const;
188 // Attribute definitions for current solid.
189
190 std::vector<G4AttValue>* CreateCurrentAttValues() const;
191 // Attribute values for current solid. Each must refer to an
192 // attribute definition in the above map; its name is the key. The
193 // user must test the validity of this pointer (it must be non-zero
194 // and conform to the G4AttDefs, which may be checked with
195 // G4AttCheck) and delete the list after use. See
196 // G4XXXStoredSceneHandler::PreAddSolid for how to access and
197 // G4VTrajectory::ShowTrajectory for an example of the use of
198 // G4Atts.
199
201 (const std::vector<G4PhysicalVolumeNodeID>&
202 baseFullPVPath) {
203 fBaseFullPVPath = baseFullPVPath;
204 fFullPVPath = baseFullPVPath;
205 }
206
207 void SetRequestedDepth (G4int requestedDepth) {
208 fRequestedDepth = requestedDepth;
209 }
210
211 void SetClippingSolid (G4VSolid* pClippingSolid) {
212 fpClippingSolid = pClippingSolid;
213 }
214
216 fClippingMode = mode;
217 }
218
219 G4bool Validate (G4bool warn);
220 // Validate, but allow internal changes (hence non-const function).
221
223
224protected:
225
227 G4int requestedDepth,
228 const G4Transform3D&,
230
232 G4int requestedDepth,
234 G4VSolid*,
235 G4Material*,
236 const G4Transform3D&,
238
239 virtual void DescribeSolid (const G4Transform3D& theAT,
240 G4VSolid* pSol,
241 const G4VisAttributes* pVisAttribs,
242 G4VGraphicsScene& sceneHandler);
243
244 void CalculateExtent ();
245
246 /////////////////////////////////////////////////////////
247 // Data members...
248
249 G4VPhysicalVolume* fpTopPV; // The physical volume.
250 G4String fTopPVName; // ...of the physical volume.
251 G4int fTopPVCopyNo; // ...of the physical volume.
253 // Requested depth of geom. hierarchy search.
254 G4bool fUseFullExtent; // ...if requested.
255 G4int fCurrentDepth; // Current depth of geom. hierarchy.
256 G4VPhysicalVolume* fpCurrentPV; // Current physical volume.
257 G4LogicalVolume* fpCurrentLV; // Current logical volume.
258 G4Material* fpCurrentMaterial; // Current material.
259 G4Transform3D* fpCurrentTransform; // Current transform.
260 std::vector<G4PhysicalVolumeNodeID> fBaseFullPVPath;
261 std::vector<G4PhysicalVolumeNodeID> fFullPVPath;
262 std::vector<G4PhysicalVolumeNodeID> fDrawnPVPath;
263 G4bool fCurtailDescent;// Can be set to curtail descent.
266
267private:
268
269 // Private copy constructor and assigment operator - copying and
270 // assignment not allowed. Keeps CodeWizard happy.
272 G4PhysicalVolumeModel& operator = (const G4PhysicalVolumeModel&);
273};
274
275std::ostream& operator<<
276 (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID);
277
278#endif
HepGeom::Transform3D G4Transform3D
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
const G4RotationMatrix * GetRotation(G4int depth) const
const G4ThreeVector & GetTranslation(G4int depth) const
G4PhysicalVolumeNodeID(G4VPhysicalVolume *pPV=0, G4int iCopyNo=0, G4int depth=0, const G4Transform3D &transform=G4Transform3D(), G4bool drawn=true)
G4bool operator<(const G4PhysicalVolumeNodeID &right) const
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
const G4VSolid * GetClippingSolid() const
void DescribeAndDescend(G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
G4VPhysicalVolume * GetCurrentPV() const
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
void SetClippingSolid(G4VSolid *pClippingSolid)
G4VPhysicalVolume * fpTopPV
void VisitGeometryAndGetVisReps(G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
std::vector< G4AttValue > * CreateCurrentAttValues() const
G4Transform3D * fpCurrentTransform
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
G4bool Validate(G4bool warn)
void SetRequestedDepth(G4int requestedDepth)
virtual void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
G4String GetCurrentDescription() const
void SetClippingMode(ClippingMode mode)
void DescribeYourselfTo(G4VGraphicsScene &)
G4LogicalVolume * GetCurrentLV() const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath
G4VPhysicalVolume * GetTopPhysicalVolume() const
void SetBaseFullPVPath(const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath)
G4VPhysicalVolume * fpCurrentPV
G4Material * GetCurrentMaterial() const
const std::map< G4String, G4AttDef > * GetAttDefs() const