Geant4-11
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,
98 G4bool drawn = true):
99 fpPV(pPV),
100 fCopyNo(iCopyNo),
101 fNonCulledDepth(depth),
103 fDrawn(drawn) {}
105 G4int GetCopyNo() const {return fCopyNo;}
107 const G4Transform3D& GetTransform() const {return fTransform;}
108 G4bool GetDrawn() const {return fDrawn;}
109 void SetDrawn(G4bool drawn) {fDrawn = drawn;}
110 G4bool operator< (const G4PhysicalVolumeNodeID& right) const;
111 G4bool operator!= (const G4PhysicalVolumeNodeID& right) const;
113 return !operator!= (right);
114 }
115 private:
121 };
122
123 // Nested class for handling nested parameterisations.
125 public:
127 (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath);
128 const G4ThreeVector& GetTranslation(G4int depth) const;
129 const G4RotationMatrix* GetRotation(G4int depth) const;
130 G4VPhysicalVolume* GetVolume(G4int depth) const;
131 G4VSolid* GetSolid(G4int depth) const;
132 G4int GetReplicaNumber(G4int depth) const;
133 G4int GetHistoryDepth() const {return G4int(fFullPVPath.size());}
134 private:
135 const std::vector<G4PhysicalVolumeNodeID>& fFullPVPath;
136 };
137
138 // Nested struct for encapsulating touchable properties
145 std::vector<G4PhysicalVolumeNodeID> fTouchableBaseFullPVPath;
146 };
147
149 (G4VPhysicalVolume* = 0,
150 G4int requestedDepth = UNLIMITED,
151 const G4Transform3D& modelTransformation = G4Transform3D(),
152 const G4ModelingParameters* = 0,
153 G4bool useFullExtent = false,
154 const std::vector<G4PhysicalVolumeNodeID>& baseFullPVPath =
155 std::vector<G4PhysicalVolumeNodeID>());
156
157 virtual ~G4PhysicalVolumeModel ();
158
160 // The main task of a model is to describe itself to the graphics scene
161 // handler (a object which inherits G4VSceneHandler, which inherits
162 // G4VGraphicsScene).
163
165 // A description which depends on the current state of the model.
166
167 G4String GetCurrentTag () const;
168 // A tag which depends on the current state of the model.
169
171
173
175 {return fpClippingSolid;}
176
178 // Current depth in geom. hierarchy.
179
181 // Initial transformation
182
184 // Current physical volume.
185
187 // Current copy number.
188
190 // Current logical volume.
191
193 // Current material.
194
196 // Current transform.
197
198 const std::vector<G4PhysicalVolumeNodeID>& GetBaseFullPVPath() const
199 {return fBaseFullPVPath;}
200 // Base for this G4PhysicalVolumeModel. It can be empty, which would be the
201 // case for the world volume. But the user may create a G4PhysicalVolumeModel
202 // for a volume anywhere in the tree, in which case the user must provide
203 // the transform (in the constructor) and the base path (SetBaseFullPVPath).
204
205 const std::vector<G4PhysicalVolumeNodeID>& GetFullPVPath() const
206 {return fFullPVPath;}
207 // Vector of physical volume node identifiers for the current
208 // touchable. It is its path in the geometry hierarchy, similar to
209 // the concept of "touchable history" available from the navigator
210 // during tracking.
211
212 const std::vector<G4PhysicalVolumeNodeID>& GetDrawnPVPath() const
213 {return fDrawnPVPath;}
214 // Path of the current drawn (non-culled) touchable in terms of
215 // drawn (non-culled) ancestors. It is a vector of physical volume
216 // node identifiers corresponding to the geometry hierarchy actually
217 // selected, i.e., with "culled" volumes NOT included.
218
220 (const std::vector<G4PhysicalVolumeNodeID>&);
221 // Converts
222
223 const std::map<G4String,G4AttDef>* GetAttDefs() const;
224 // Attribute definitions for current solid.
225
226 std::vector<G4AttValue>* CreateCurrentAttValues() const;
227 // Attribute values for current solid. Each must refer to an
228 // attribute definition in the above map; its name is the key. The
229 // user must test the validity of this pointer (it must be non-zero
230 // and conform to the G4AttDefs, which may be checked with
231 // G4AttCheck) and delete the list after use. See
232 // G4XXXStoredSceneHandler::PreAddSolid for how to access and
233 // G4VTrajectory::ShowTrajectory for an example of the use of
234 // G4Atts.
235
236 void SetRequestedDepth (G4int requestedDepth) {
237 fRequestedDepth = requestedDepth;
238 }
239
240 void SetClippingSolid (G4VSolid* pClippingSolid) {
241 fpClippingSolid = pClippingSolid;
242 }
243
245 fClippingMode = mode;
246 }
247
248 G4bool Validate (G4bool warn);
249 // Validate, but allow internal changes (hence non-const function).
250
251 void Abort () const {fAbort = true;}
252 // Abort all further traversing.
253
254 void CurtailDescent() const {fCurtailDescent = true;}
255 // Curtail descent of current branch.
256
257 void CalculateExtent ();
258
259protected:
260
262 G4int requestedDepth,
263 const G4Transform3D&,
265
267 G4int requestedDepth,
269 G4VSolid*,
270 G4Material*,
271 const G4Transform3D&,
273
274 virtual void DescribeSolid (const G4Transform3D& theAT,
275 G4VSolid* pSol,
276 const G4VisAttributes* pVisAttribs,
277 G4VGraphicsScene& sceneHandler);
278
280 // Data members...
281
282 G4VPhysicalVolume* fpTopPV; // The physical volume.
283 G4String fTopPVName; // ...of the physical volume.
284 G4int fTopPVCopyNo; // ...of the physical volume.
286 // Requested depth of geom. hierarchy search.
287 G4bool fUseFullExtent; // ...if requested.
288 G4Transform3D fTransform; // Initial transform
289 G4int fCurrentDepth; // Current depth of geom. hierarchy.
290 G4VPhysicalVolume* fpCurrentPV; // Current physical volume.
291 G4int fCurrentPVCopyNo; // Current copy number.
292 G4LogicalVolume* fpCurrentLV; // Current logical volume.
293 G4Material* fpCurrentMaterial; // Current material.
294 G4Transform3D fCurrentTransform; // Current transform.
295 std::vector<G4PhysicalVolumeNodeID> fBaseFullPVPath; // Base. May be empty.
296 std::vector<G4PhysicalVolumeNodeID> fFullPVPath; // Starts from base.
297 std::vector<G4PhysicalVolumeNodeID> fDrawnPVPath; // Omits culled volumes.
298 mutable G4bool fAbort; // Abort all further traversing.
299 mutable G4bool fCurtailDescent;// Can be set to curtail descent.
302
303private:
304
305 // Private copy constructor and assigment operator - copying and
306 // assignment not allowed. Keeps CodeWizard happy.
309};
310
311std::ostream& operator<<
312(std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID&);
313
314std::ostream& operator<<
315(std::ostream& os, const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>&);
316
317#endif
HepGeom::Transform3D G4Transform3D
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
std::vector< PVNameCopyNo > PVNameCopyNoPath
const std::vector< G4PhysicalVolumeNodeID > & fFullPVPath
const G4RotationMatrix * GetRotation(G4int depth) const
G4PhysicalVolumeModelTouchable(const std::vector< G4PhysicalVolumeNodeID > &fullPVPath)
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
G4PhysicalVolumeModel(G4VPhysicalVolume *=0, G4int requestedDepth=UNLIMITED, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0, G4bool useFullExtent=false, const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath=std::vector< G4PhysicalVolumeNodeID >())
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
G4bool Validate(G4bool warn)
void SetRequestedDepth(G4int requestedDepth)
G4PhysicalVolumeModel & operator=(const G4PhysicalVolumeModel &)
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 &)
G4PhysicalVolumeModel(const G4PhysicalVolumeModel &)
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 > &)
G4bool transform(G4String &input, const G4String &type)
G4ModelingParameters::PVNameCopyNoPath fTouchablePath
std::vector< G4PhysicalVolumeNodeID > fTouchableBaseFullPVPath