Geant4 Cross Reference |
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" 65 #include "G4ModelingParameters.hh" 66 67 #include "G4VTouchable.hh" 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 82 class G4PhysicalVolumeModel: public G4VModel { 83 84 public: // With description 85 86 enum {UNLIMITED = -1}; 87 88 enum ClippingMode {subtraction, intersection}; 89 90 // Nested class for identifying physical volume nodes. 91 class G4PhysicalVolumeNodeID { 92 public: 93 G4PhysicalVolumeNodeID 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;} 109 void SetPhysicalVolume(G4VPhysicalVolume* v) {fpPV = v;} 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; 116 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. 128 class G4PhysicalVolumeModelTouchable: public G4VTouchable { 129 public: 130 G4PhysicalVolumeModelTouchable 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 143 struct TouchableProperties { 144 TouchableProperties(): fpTouchablePV(nullptr), fCopyNo(0) {} 145 G4ModelingParameters::PVNameCopyNoPath fTouchablePath; 146 G4VPhysicalVolume* fpTouchablePV; 147 G4int fCopyNo; 148 G4Transform3D fTouchableGlobalTransform; 149 std::vector<G4PhysicalVolumeNodeID> fTouchableBaseFullPVPath; 150 std::vector<G4PhysicalVolumeNodeID> fTouchableFullPVPath; 151 }; 152 153 G4PhysicalVolumeModel 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 164 void DescribeYourselfTo (G4VGraphicsScene&); 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 169 G4String GetCurrentDescription () const; 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 175 G4VPhysicalVolume* GetTopPhysicalVolume () const {return fpTopPV;} 176 177 G4int GetRequestedDepth () const {return fRequestedDepth;} 178 179 const G4VSolid* GetClippingSolid () const 180 {return fpClippingSolid;} 181 182 G4int GetCurrentDepth() const {return fCurrentDepth;} 183 // Current depth in geom. hierarchy. 184 185 const G4Transform3D& GetTransformation() const {return fTransform;} 186 // Initial transformation 187 188 G4VPhysicalVolume* GetCurrentPV() const {return fpCurrentPV;} 189 // Current physical volume. 190 191 G4int GetCurrentPVCopyNo() const {return fCurrentPVCopyNo;} 192 // Current copy number. 193 194 G4LogicalVolume* GetCurrentLV() const {return fpCurrentLV;} 195 // Current logical volume. 196 197 G4Material* GetCurrentMaterial() const {return fpCurrentMaterial;} 198 // Current material. 199 200 const G4Transform3D& GetCurrentTransform() const {return fCurrentTransform;} 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 224 static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath 225 (const std::vector<G4PhysicalVolumeNodeID>&); 226 // Converts to PVNameCopyNoPath 227 228 static G4String GetPVNamePathString(const std::vector<G4PhysicalVolumeNodeID>&); 229 // Converts to path string, e.g., " World 0 Envelope 0 Shape1 0". 230 // Note leading space character. 231 232 const std::map<G4String,G4AttDef>* GetAttDefs() const; 233 // Attribute definitions for current solid. 234 235 std::vector<G4AttValue>* CreateCurrentAttValues() const; 236 // Attribute values for current solid. Each must refer to an 237 // attribute definition in the above map; its name is the key. The 238 // user must test the validity of this pointer (it must be non-zero 239 // and conform to the G4AttDefs, which may be checked with 240 // G4AttCheck) and delete the list after use. See 241 // G4XXXStoredSceneHandler::PreAddSolid for how to access and 242 // G4VTrajectory::ShowTrajectory for an example of the use of 243 // G4Atts. 244 245 const std::map<G4int,G4int>& GetNumberOfTouchables() const {return fNTouchables;} 246 // Number of touchables drawn at each depth. 247 248 G4int GetTotalTouchables () {return fTotalTouchables;} 249 // Total numbere of touchables. 250 251 void SetRequestedDepth (G4int requestedDepth) { 252 fRequestedDepth = requestedDepth; 253 } 254 255 void SetClippingSolid (G4VSolid* pClippingSolid) { 256 fpClippingSolid = pClippingSolid; 257 } 258 259 void SetClippingMode (ClippingMode mode) { 260 fClippingMode = mode; 261 } 262 263 G4bool Validate (G4bool warn); 264 // Validate, but allow internal changes (hence non-const function). 265 266 void Abort () const {fAbort = true;} 267 // Abort all further traversing. 268 269 void CurtailDescent() const {fCurtailDescent = true;} 270 // Curtail descent of current branch. 271 272 void CalculateExtent (); 273 274 protected: 275 276 void VisitGeometryAndGetVisReps (G4VPhysicalVolume*, 277 G4int requestedDepth, 278 const G4Transform3D&, 279 G4VGraphicsScene&); 280 281 void DescribeAndDescend (G4VPhysicalVolume*, 282 G4int requestedDepth, 283 G4LogicalVolume*, 284 G4VSolid*, 285 G4Material*, 286 const G4Transform3D&, 287 G4VGraphicsScene&); 288 289 virtual void DescribeSolid (const G4Transform3D& theAT, 290 G4VSolid* pSol, 291 const G4VisAttributes* pVisAttribs, 292 G4VGraphicsScene& sceneHandler); 293 294 ///////////////////////////////////////////////////////// 295 // Data members... 296 297 G4VPhysicalVolume* fpTopPV; // The physical volume. 298 G4String fTopPVName; // ...of the physical volume. 299 G4int fTopPVCopyNo; // ...of the physical volume. 300 G4int fRequestedDepth; 301 // Requested depth of geom. hierarchy search. 302 G4bool fUseFullExtent; // ...if requested. 303 G4Transform3D fTransform; // Initial transform 304 G4int fCurrentDepth; // Current depth of geom. hierarchy. 305 G4VPhysicalVolume* fpCurrentPV; // Current physical volume. 306 G4int fCurrentPVCopyNo; // Current copy number. 307 G4LogicalVolume* fpCurrentLV; // Current logical volume. 308 G4Material* fpCurrentMaterial; // Current material. 309 G4Transform3D fCurrentTransform; // Current transform. 310 std::vector<G4PhysicalVolumeNodeID> fBaseFullPVPath; // Base. May be empty. 311 std::vector<G4PhysicalVolumeNodeID> fFullPVPath; // Starts from base. 312 std::vector<G4PhysicalVolumeNodeID> fDrawnPVPath; // Omits culled volumes. 313 mutable G4bool fAbort; // Abort all further traversing. 314 mutable G4bool fCurtailDescent;// Can be set to curtail descent. 315 G4VSolid* fpClippingSolid; 316 ClippingMode fClippingMode; 317 G4int fNClippers; // No of clipping/cutting solids - only 0 or 1 allowed 318 std::map<G4int,G4int> fNTouchables; // No of touchables at each depth 319 G4int fTotalTouchables; 320 321 private: 322 323 // Private copy constructor and assigment operator - copying and 324 // assignment not allowed. Keeps CodeWizard happy. 325 G4PhysicalVolumeModel (const G4PhysicalVolumeModel&); 326 G4PhysicalVolumeModel& operator = (const G4PhysicalVolumeModel&); 327 }; 328 329 std::ostream& operator<< 330 (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID&); 331 332 std::ostream& operator<< 333 (std::ostream& os, const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>&); 334 335 #endif 336