Geant4  10.03.p01
 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 // $Id: G4PhysicalVolumeModel.hh 99668 2016-09-30 08:02:13Z gcosmo $
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 
74 class G4VPhysicalVolume;
75 class G4LogicalVolume;
76 class G4VSolid;
77 class G4Material;
78 class G4VisAttributes;
79 class G4AttDef;
80 class G4AttValue;
81 
83 
84 public: // 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 
158  const G4VSolid* GetClippingSolid () const
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 
174  // Current transform.
175 
176  const std::vector<G4PhysicalVolumeNodeID>& GetFullPVPath() const
177  {return fFullPVPath;}
178  // Vector of physical volume node identifiers for the current
179  // touchable. It is its path in the geometry hierarchy, similar to
180  // the concept of "touchable history" available from the navigator
181  // during tracking.
182 
183  const std::vector<G4PhysicalVolumeNodeID>& GetDrawnPVPath() const
184  {return fDrawnPVPath;}
185  // Path of the current drawn (non-culled) touchable in terms of
186  // drawn (non-culled) ancesters. It is a vector of physical volume
187  // node identifiers corresponding to the geometry hierarchy actually
188  // selected, i.e., not culled.
189 
190  const std::map<G4String,G4AttDef>* GetAttDefs() const;
191  // Attribute definitions for current solid.
192 
193  std::vector<G4AttValue>* CreateCurrentAttValues() const;
194  // Attribute values for current solid. Each must refer to an
195  // attribute definition in the above map; its name is the key. The
196  // user must test the validity of this pointer (it must be non-zero
197  // and conform to the G4AttDefs, which may be checked with
198  // G4AttCheck) and delete the list after use. See
199  // G4XXXStoredSceneHandler::PreAddSolid for how to access and
200  // G4VTrajectory::ShowTrajectory for an example of the use of
201  // G4Atts.
202 
203  void SetBaseFullPVPath
204  (const std::vector<G4PhysicalVolumeNodeID>&
205  baseFullPVPath) {
206  fBaseFullPVPath = baseFullPVPath;
207  fFullPVPath = baseFullPVPath;
208  }
209 
210  void SetRequestedDepth (G4int requestedDepth) {
211  fRequestedDepth = requestedDepth;
212  }
213 
214  void SetClippingSolid (G4VSolid* pClippingSolid) {
215  fpClippingSolid = pClippingSolid;
216  }
217 
219  fClippingMode = mode;
220  }
221 
222  G4bool Validate (G4bool warn);
223  // Validate, but allow internal changes (hence non-const function).
224 
225  void Abort () const {fAbort = true;}
226  // Abort all further traversing.
227 
228  void CurtailDescent() const {fCurtailDescent = true;}
229  // Curtail descent of current branch.
230 
231 protected:
232 
234  G4int requestedDepth,
235  const G4Transform3D&,
237 
239  G4int requestedDepth,
241  G4VSolid*,
242  G4Material*,
243  const G4Transform3D&,
245 
246  virtual void DescribeSolid (const G4Transform3D& theAT,
247  G4VSolid* pSol,
248  const G4VisAttributes* pVisAttribs,
249  G4VGraphicsScene& sceneHandler);
250 
251  void CalculateExtent ();
252 
254  // Data members...
255 
256  G4VPhysicalVolume* fpTopPV; // The physical volume.
257  G4String fTopPVName; // ...of the physical volume.
258  G4int fTopPVCopyNo; // ...of the physical volume.
260  // Requested depth of geom. hierarchy search.
261  G4bool fUseFullExtent; // ...if requested.
262  G4int fCurrentDepth; // Current depth of geom. hierarchy.
263  G4VPhysicalVolume* fpCurrentPV; // Current physical volume.
264  G4LogicalVolume* fpCurrentLV; // Current logical volume.
265  G4Material* fpCurrentMaterial; // Current material.
266  G4Transform3D* fpCurrentTransform; // Current transform.
267  std::vector<G4PhysicalVolumeNodeID> fBaseFullPVPath;
268  std::vector<G4PhysicalVolumeNodeID> fFullPVPath;
269  std::vector<G4PhysicalVolumeNodeID> fDrawnPVPath;
270  mutable G4bool fAbort; // Abort all further traversing.
271  mutable G4bool fCurtailDescent;// Can be set to curtail descent.
274 
275 private:
276 
277  // Private copy constructor and assigment operator - copying and
278  // assignment not allowed. Keeps CodeWizard happy.
280  G4PhysicalVolumeModel& operator = (const G4PhysicalVolumeModel&);
281 };
282 
283 std::ostream& operator<<
284 (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID&);
285 
286 std::ostream& operator<<
287 (std::ostream& os, const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>&);
288 
289 #endif
G4bool operator<(const G4PhysicalVolumeNodeID &right) const
void SetClippingSolid(G4VSolid *pClippingSolid)
const std::map< G4String, G4AttDef > * GetAttDefs() const
G4Material * GetCurrentMaterial() const
const G4ThreeVector & GetTranslation(G4int depth) const
G4bool Validate(G4bool warn)
G4Transform3D * fpCurrentTransform
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
int G4int
Definition: G4Types.hh:78
G4VPhysicalVolume * fpCurrentPV
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
G4Transform3D * GetCurrentTransform() const
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath
bool G4bool
Definition: G4Types.hh:79
const G4VSolid * GetClippingSolid() const
G4String GetCurrentDescription() const
HepGeom::Transform3D G4Transform3D
G4PhysicalVolumeNodeID(G4VPhysicalVolume *pPV=0, G4int iCopyNo=0, G4int depth=0, const G4Transform3D &transform=G4Transform3D(), G4bool drawn=true)
void SetBaseFullPVPath(const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath)
virtual void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
G4VPhysicalVolume * GetTopPhysicalVolume() const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
G4VPhysicalVolume * fpTopPV
G4PhysicalVolumeModelTouchable(const std::vector< G4PhysicalVolumeNodeID > &fullPVPath)
std::vector< G4AttValue > * CreateCurrentAttValues() const
void SetClippingMode(ClippingMode mode)
G4String GetCurrentTag() const
G4PhysicalVolumeModel(G4VPhysicalVolume *=0, G4int requestedDepth=UNLIMITED, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0, G4bool useFullExtent=false)
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
void VisitGeometryAndGetVisReps(G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
G4LogicalVolume * GetCurrentLV() const
G4VPhysicalVolume * GetCurrentPV() const
void DescribeAndDescend(G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
const G4RotationMatrix * GetRotation(G4int depth) const
void SetRequestedDepth(G4int requestedDepth)
void DescribeYourselfTo(G4VGraphicsScene &)