source: trunk/source/visualization/modeling/include/G4PhysicalVolumeModel.hh@ 1155

Last change on this file since 1155 was 1140, checked in by garnier, 16 years ago

update to CVS

File size: 9.1 KB
RevLine 
[834]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//
[1140]27// $Id: G4PhysicalVolumeModel.hh,v 1.35 2009/10/10 14:29:59 allison Exp $
[954]28// GEANT4 tag $Name: $
[834]29//
30//
31// John Allison 31st December 1997.
32//
33// Class Description:
34//
35// Model for physical volumes. It describes a physical volume and its
36// daughters to any desired depth. Note: the "requested depth" is
37// specified in the modeling parameters; enum {UNLIMITED = -1}.
38//
39// For access to base class information, e.g., modeling parameters,
40// use GetModelingParameters() inherited from G4VModel. See Class
41// Description of the base class G4VModel.
42//
43// G4PhysicalVolumeModel assumes the modeling parameters have been set
44// up with meaningful information - default vis attributes and culling
45// policy in particular.
46
47#ifndef G4PHYSICALVOLUMEMODEL_HH
48#define G4PHYSICALVOLUMEMODEL_HH
49
50#include "G4VModel.hh"
[1140]51#include "G4VTouchable.hh"
[834]52
53#include "G4Transform3D.hh"
54#include "G4Plane3D.hh"
55#include <iostream>
56#include <vector>
57#include <map>
58
59class G4VPhysicalVolume;
60class G4LogicalVolume;
61class G4VSolid;
62class G4Material;
63class G4VisAttributes;
64class G4Polyhedron;
65class G4AttDef;
66class G4AttValue;
67
68class G4PhysicalVolumeModel: public G4VModel {
69
70public: // With description
71
72 enum {UNLIMITED = -1};
73
74 enum ClippingMode {subtraction, intersection};
75
76 class G4PhysicalVolumeNodeID {
77 public:
78 G4PhysicalVolumeNodeID
[1140]79 (G4VPhysicalVolume* pPV = 0,
80 G4int iCopyNo = 0,
81 G4int depth = 0,
82 const G4Transform3D& transform = G4Transform3D()):
83 fpPV(pPV),
84 fCopyNo(iCopyNo),
85 fNonCulledDepth(depth),
86 fTransform(transform) {}
[834]87 G4VPhysicalVolume* GetPhysicalVolume() const {return fpPV;}
88 G4int GetCopyNo() const {return fCopyNo;}
89 G4int GetNonCulledDepth() const {return fNonCulledDepth;}
[1140]90 const G4Transform3D& GetTransform() const {return fTransform;}
[834]91 G4bool operator< (const G4PhysicalVolumeNodeID& right) const;
92 private:
93 G4VPhysicalVolume* fpPV;
94 G4int fCopyNo;
95 G4int fNonCulledDepth;
[1140]96 G4Transform3D fTransform;
[834]97 };
98 // Nested class for identifying physical volume nodes.
99
[1140]100 class G4PhysicalVolumeModelTouchable: public G4VTouchable {
101 public:
102 G4PhysicalVolumeModelTouchable
103 (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath);
104 const G4ThreeVector& GetTranslation(G4int depth) const;
105 const G4RotationMatrix* GetRotation(G4int depth) const;
106 G4VPhysicalVolume* GetVolume(G4int depth) const;
107 G4VSolid* GetSolid(G4int depth) const;
108 G4int GetReplicaNumber(G4int depth) const;
109 G4int GetHistoryDepth() const {return fFullPVPath.size();}
110 private:
111 const std::vector<G4PhysicalVolumeNodeID>& fFullPVPath;
112 };
113 // Nested class for handling nested parameterisations.
114
[834]115 G4PhysicalVolumeModel
116 (G4VPhysicalVolume*,
117 G4int requestedDepth = UNLIMITED,
118 const G4Transform3D& modelTransformation = G4Transform3D(),
119 const G4ModelingParameters* = 0,
120 G4bool useFullExtent = false);
121
122 virtual ~G4PhysicalVolumeModel ();
123
124 void DescribeYourselfTo (G4VGraphicsScene&);
125 // The main task of a model is to describe itself to the graphics scene
126 // handler (a object which inherits G4VSceneHandler, which inherits
127 // G4VGraphicsScene).
128
129 G4String GetCurrentDescription () const;
130 // A description which depends on the current state of the model.
131
132 G4String GetCurrentTag () const;
133 // A tag which depends on the current state of the model.
134
135 G4VPhysicalVolume* GetTopPhysicalVolume () const {return fpTopPV;}
136
137 G4int GetRequestedDepth () const {return fRequestedDepth;}
138
139 const G4Polyhedron* GetClippingPolyhedron () const
140 {return fpClippingPolyhedron;}
141
142 G4int GetCurrentDepth() const {return fCurrentDepth;}
143 // Current depth of geom. hierarchy.
144
145 G4VPhysicalVolume* GetCurrentPV() const {return fpCurrentPV;}
146 // Current physical volume.
147
148 G4LogicalVolume* GetCurrentLV() const {return fpCurrentLV;}
149 // Current logical volume.
150
151 G4Material* GetCurrentMaterial() const {return fpCurrentMaterial;}
152 // Current material.
153
154 const std::vector<G4PhysicalVolumeNodeID>& GetDrawnPVPath() const
155 {return fDrawnPVPath;}
156 // Path of the current drawn (non-culled) volume in terms of drawn
157 // (non-culled) ancesters. It is a vector of physical volume node
158 // identifiers corresponding to the geometry hierarchy actually
159 // selected, i.e., not culled. It is guaranteed that volumes are
160 // presented to the scene handlers in top-down hierarchy order,
161 // i.e., ancesters first, mothers before daughters, so the scene
162 // handler can be assured that, if it is building its own scene
163 // graph tree, a mother, if any, will have already been encountered
164 // and there will already be a node in place on which to hang the
165 // current volume.
166
167 const std::map<G4String,G4AttDef>* GetAttDefs() const;
168 // Attribute definitions for current solid.
169
170 std::vector<G4AttValue>* CreateCurrentAttValues() const;
171 // Attribute values for current solid. Each must refer to an
172 // attribute definition in the above map; its name is the key. The
173 // user must test the validity of this pointer (it must be non-zero
174 // and conform to the G4AttDefs, which may be checked with
175 // G4AttCheck) and delete the list after use. See
176 // G4XXXStoredSceneHandler::PreAddSolid for how to access and
177 // G4VTrajectory::ShowTrajectory for an example of the use of
178 // G4Atts.
179
180 void SetRequestedDepth (G4int requestedDepth) {
181 fRequestedDepth = requestedDepth;
182 }
183
184 void SetClippingPolyhedron (const G4Polyhedron* pClippingPolyhedron) {
185 fpClippingPolyhedron = pClippingPolyhedron;
186 }
187
188 void SetClippingMode (ClippingMode mode) {
189 fClippingMode = mode;
190 }
191
192 G4bool Validate (G4bool warn);
193 // Validate, but allow internal changes (hence non-const function).
194
195 void CurtailDescent() {fCurtailDescent = true;}
196
197protected:
198
199 void VisitGeometryAndGetVisReps (G4VPhysicalVolume*,
200 G4int requestedDepth,
201 const G4Transform3D&,
202 G4VGraphicsScene&);
203
204 void DescribeAndDescend (G4VPhysicalVolume*,
205 G4int requestedDepth,
206 G4LogicalVolume*,
207 G4VSolid*,
208 G4Material*,
209 const G4Transform3D&,
210 G4VGraphicsScene&);
211
212 virtual void DescribeSolid (const G4Transform3D& theAT,
213 G4VSolid* pSol,
214 const G4VisAttributes* pVisAttribs,
215 G4VGraphicsScene& sceneHandler);
216
217 void CalculateExtent ();
218
219 /////////////////////////////////////////////////////////
220 // Data members...
221
222 G4VPhysicalVolume* fpTopPV; // The physical volume.
223 G4String fTopPVName; // ...of the physical volume.
224 G4int fTopPVCopyNo; // ...of the physical volume.
225 G4int fRequestedDepth;
226 // Requested depth of geom. hierarchy search.
227 G4bool fUseFullExtent; // ...if requested.
228 G4int fCurrentDepth; // Current depth of geom. hierarchy.
229 G4VPhysicalVolume* fpCurrentPV; // Current physical volume.
230 G4LogicalVolume* fpCurrentLV; // Current logical volume.
231 G4Material* fpCurrentMaterial; // Current material.
232 G4Transform3D* fpCurrentTransform; // Current transform.
233 std::vector<G4PhysicalVolumeNodeID> fFullPVPath;
234 std::vector<G4PhysicalVolumeNodeID> fDrawnPVPath;
235 G4bool fCurtailDescent;// Can be set to curtail descent.
236 const G4Polyhedron*fpClippingPolyhedron;
237 ClippingMode fClippingMode;
238
239private:
240
241 // Private copy constructor and assigment operator - copying and
242 // assignment not allowed. Keeps CodeWizard happy.
243 G4PhysicalVolumeModel (const G4PhysicalVolumeModel&);
244 G4PhysicalVolumeModel& operator = (const G4PhysicalVolumeModel&);
245};
246
247std::ostream& operator<<
248 (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID);
249
250#endif
Note: See TracBrowser for help on using the repository browser.