Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
28//
29// John Allison 31st December 1997.
30//
31// Class Description:
32//
33// Model for physical volumes. It describes a physical volume and its
34// daughters to any desired depth. Note: the "requested depth" is
35// specified in the modeling parameters; enum {UNLIMITED = -1}.
36//
37// For access to base class information, e.g., modeling parameters,
38// use GetModelingParameters() inherited from G4VModel. See Class
39// Description of the base class G4VModel.
40//
41// G4PhysicalVolumeModel assumes the modeling parameters have been set
42// up with meaningful information - default vis attributes and culling
43// policy in particular.
44//
45// The volumes are unpacked and sent to the scene handler. They are,
46// in effect, "touchables" as defined by G4TouchableHistory. A
47// touchable is defined not only by its physical volume name and copy
48// number but also by its position in the geometry hierarchy.
49//
50// It is guaranteed that touchables are presented to the scene handler
51// in top-down hierarchy order, i.e., ancestors first, mothers before
52// daughters, so the scene handler can be assured that, if it is
53// building its own scene graph tree, a mother, if any, will have
54// already been encountered and there will already be a node in place
55// on which to hang the current volume. But be aware that the
56// visibility and culling policy might mean that some touchables are
57// not passed to the scene handler so the drawn tree might have
58// missing layers. GetFullPVPath allows you to know the full
59// hierarchy.
60
61#ifndef G4PHYSICALVOLUMEMODEL_HH
62#define G4PHYSICALVOLUMEMODEL_HH
63
64#include "G4VModel.hh"
66
67#include "G4VTouchable.hh"
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
90 // Nested class for identifying physical volume nodes.
92 public:
94 (G4VPhysicalVolume* pPV = 0,
95 G4int iCopyNo = 0,
96 G4int depth = 0,
97 const G4Transform3D& transform = G4Transform3D(),
98 G4bool drawn = true):
99 fpPV(pPV),
100 fCopyNo(iCopyNo),
101 fNonCulledDepth(depth),
102 fTransform(transform),
103 fDrawn(drawn) {}
104 G4VPhysicalVolume* GetPhysicalVolume() const {return fpPV;}
105 G4int GetCopyNo() const {return fCopyNo;}
106 G4int GetNonCulledDepth() const {return fNonCulledDepth;}
107 const G4Transform3D& GetTransform() const {return fTransform;}
108 G4bool GetDrawn() const {return fDrawn;}
110 void SetCopyNo (G4int n) {fCopyNo = n;}
111 void SetNonCulledDepth(G4int d) {fNonCulledDepth = d;}
112 void SetTransform (const G4Transform3D& t) {fTransform = t;}
113 void SetDrawn (G4bool b) {fDrawn = b;}
114 G4bool operator< (const G4PhysicalVolumeNodeID& right) const;
115 G4bool operator!= (const G4PhysicalVolumeNodeID& right) const;
117 return !operator!= (right);
118 }
119 private:
120 G4VPhysicalVolume* fpPV;
121 G4int fCopyNo;
122 G4int fNonCulledDepth;
123 G4Transform3D fTransform;
124 G4bool fDrawn;
125 };
126
127 // Nested class for handling nested parameterisations.
129 public:
131 (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath);
132 const G4ThreeVector& GetTranslation(G4int depth) const;
133 const G4RotationMatrix* GetRotation(G4int depth) const;
134 G4VPhysicalVolume* GetVolume(G4int depth) const;
135 G4VSolid* GetSolid(G4int depth) const;
136 G4int GetReplicaNumber(G4int depth) const;
137 G4int GetHistoryDepth() const {return G4int(fFullPVPath.size());}
138 private:
139 const std::vector<G4PhysicalVolumeNodeID>& fFullPVPath;
140 };
141
142 // Nested struct for encapsulating touchable properties
149 std::vector<G4PhysicalVolumeNodeID> fTouchableBaseFullPVPath;
150 std::vector<G4PhysicalVolumeNodeID> fTouchableFullPVPath;
151 };
152
154 (G4VPhysicalVolume* = 0,
155 G4int requestedDepth = UNLIMITED,
156 const G4Transform3D& modelTransformation = G4Transform3D(),
157 const G4ModelingParameters* = 0,
158 G4bool useFullExtent = false,
159 const std::vector<G4PhysicalVolumeNodeID>& baseFullPVPath =
160 std::vector<G4PhysicalVolumeNodeID>());
161
162 virtual ~G4PhysicalVolumeModel ();
163
165 // The main task of a model is to describe itself to the graphics scene
166 // handler (a object which inherits G4VSceneHandler, which inherits
167 // G4VGraphicsScene).
168
170 // A description which depends on the current state of the model.
171
172 G4String GetCurrentTag () const;
173 // A tag which depends on the current state of the model.
174
176
178
180 {return fpClippingSolid;}
181
183 // Current depth in geom. hierarchy.
184
186 // Initial transformation
187
189 // Current physical volume.
190
192 // Current copy number.
193
195 // Current logical volume.
196
198 // Current material.
199
201 // Current transform.
202
203 const std::vector<G4PhysicalVolumeNodeID>& GetBaseFullPVPath() const
204 {return fBaseFullPVPath;}
205 // Base for this G4PhysicalVolumeModel. It can be empty, which would be the
206 // case for the world volume. But the user may create a G4PhysicalVolumeModel
207 // for a volume anywhere in the tree, in which case the user must provide
208 // the transform (in the constructor) and the base path (SetBaseFullPVPath).
209
210 const std::vector<G4PhysicalVolumeNodeID>& GetFullPVPath() const
211 {return fFullPVPath;}
212 // Vector of physical volume node identifiers for the current
213 // touchable. It is its path in the geometry hierarchy, similar to
214 // the concept of "touchable history" available from the navigator
215 // during tracking.
216
217 const std::vector<G4PhysicalVolumeNodeID>& GetDrawnPVPath() const
218 {return fDrawnPVPath;}
219 // Path of the current drawn (non-culled) touchable in terms of
220 // drawn (non-culled) ancestors. It is a vector of physical volume
221 // node identifiers corresponding to the geometry hierarchy actually
222 // selected, i.e., with "culled" volumes NOT included.
223
225 (const std::vector<G4PhysicalVolumeNodeID>&);
226 // Converts
227
228 const std::map<G4String,G4AttDef>* GetAttDefs() const;
229 // Attribute definitions for current solid.
230
231 std::vector<G4AttValue>* CreateCurrentAttValues() const;
232 // Attribute values for current solid. Each must refer to an
233 // attribute definition in the above map; its name is the key. The
234 // user must test the validity of this pointer (it must be non-zero
235 // and conform to the G4AttDefs, which may be checked with
236 // G4AttCheck) and delete the list after use. See
237 // G4XXXStoredSceneHandler::PreAddSolid for how to access and
238 // G4VTrajectory::ShowTrajectory for an example of the use of
239 // G4Atts.
240
241 void SetRequestedDepth (G4int requestedDepth) {
242 fRequestedDepth = requestedDepth;
243 }
244
245 void SetClippingSolid (G4VSolid* pClippingSolid) {
246 fpClippingSolid = pClippingSolid;
247 }
248
250 fClippingMode = mode;
251 }
252
253 G4bool Validate (G4bool warn);
254 // Validate, but allow internal changes (hence non-const function).
255
256 void Abort () const {fAbort = true;}
257 // Abort all further traversing.
258
259 void CurtailDescent() const {fCurtailDescent = true;}
260 // Curtail descent of current branch.
261
262 void CalculateExtent ();
263
264protected:
265
267 G4int requestedDepth,
268 const G4Transform3D&,
270
272 G4int requestedDepth,
274 G4VSolid*,
275 G4Material*,
276 const G4Transform3D&,
278
279 virtual void DescribeSolid (const G4Transform3D& theAT,
280 G4VSolid* pSol,
281 const G4VisAttributes* pVisAttribs,
282 G4VGraphicsScene& sceneHandler);
283
284 /////////////////////////////////////////////////////////
285 // Data members...
286
287 G4VPhysicalVolume* fpTopPV; // The physical volume.
288 G4String fTopPVName; // ...of the physical volume.
289 G4int fTopPVCopyNo; // ...of the physical volume.
291 // Requested depth of geom. hierarchy search.
292 G4bool fUseFullExtent; // ...if requested.
293 G4Transform3D fTransform; // Initial transform
294 G4int fCurrentDepth; // Current depth of geom. hierarchy.
295 G4VPhysicalVolume* fpCurrentPV; // Current physical volume.
296 G4int fCurrentPVCopyNo; // Current copy number.
297 G4LogicalVolume* fpCurrentLV; // Current logical volume.
298 G4Material* fpCurrentMaterial; // Current material.
299 G4Transform3D fCurrentTransform; // Current transform.
300 std::vector<G4PhysicalVolumeNodeID> fBaseFullPVPath; // Base. May be empty.
301 std::vector<G4PhysicalVolumeNodeID> fFullPVPath; // Starts from base.
302 std::vector<G4PhysicalVolumeNodeID> fDrawnPVPath; // Omits culled volumes.
303 mutable G4bool fAbort; // Abort all further traversing.
304 mutable G4bool fCurtailDescent;// Can be set to curtail descent.
307
308private:
309
310 // Private copy constructor and assigment operator - copying and
311 // assignment not allowed. Keeps CodeWizard happy.
313 G4PhysicalVolumeModel& operator = (const G4PhysicalVolumeModel&);
314};
315
316std::ostream& operator<<
317(std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID&);
318
319std::ostream& operator<<
320(std::ostream& os, const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>&);
321
322#endif
HepGeom::Transform3D G4Transform3D
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
std::vector< PVNameCopyNo > PVNameCopyNoPath
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
G4bool operator<(const G4PhysicalVolumeNodeID &right) const
G4bool operator!=(const G4PhysicalVolumeNodeID &right) const
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
const G4VSolid * GetClippingSolid() const
const G4Transform3D & GetCurrentTransform() 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
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)
const G4Transform3D & GetTransformation() const
void DescribeYourselfTo(G4VGraphicsScene &)
G4LogicalVolume * GetCurrentLV() const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath
G4VPhysicalVolume * GetTopPhysicalVolume() const
const std::vector< G4PhysicalVolumeNodeID > & GetBaseFullPVPath() const
G4VPhysicalVolume * fpCurrentPV
G4Material * GetCurrentMaterial() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath(const std::vector< G4PhysicalVolumeNodeID > &)
G4ModelingParameters::PVNameCopyNoPath fTouchablePath
std::vector< G4PhysicalVolumeNodeID > fTouchableFullPVPath
std::vector< G4PhysicalVolumeNodeID > fTouchableBaseFullPVPath