Geant4  10.02.p01
G4PhysicalVolumeMassScene.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: G4PhysicalVolumeMassScene.hh 66373 2012-12-18 09:41:34Z gcosmo $
28 //
29 //
30 // John Allison 12th September 2004
31 
32 // Class Description:
33 //
34 // Calculates the mass of a geometry tree taking into account daughters
35 // up to the depth specified in the G4PhysicalVolumeModel. Culling is
36 // ignored so that all volumes are seen.
37 //
38 // Do not use this for a "parallel world" for which materials are not
39 // defined. Use only for the material world.
40 //
41 // The calculation is quite tricky, since it involves subtracting the
42 // mass of that part of the mother that is occupied by each daughter and
43 // then adding the mass of the daughter, and so on down the heirarchy.
44 //
45 // Usage for a given G4PhysicalVolumeModel* pvModel:
46 // G4PhysicalVolumeMassScene massScene(pvModel);
47 // pvModel->DescribeYourselfTo (massScene);
48 // G4double volume = massScene.GetVolume();
49 // G4double mass = massScene.GetMass();
50 // massScene.Reset();
51 // See, for example, G4ASCIITreeSceneHandler::EndModeling().
52 
53 #ifndef G4PHYSICALVOLUMEMASSSCENE_HH
54 #define G4PHYSICALVOLUMEMASSSCENE_HH
55 
56 #include "G4VGraphicsScene.hh"
57 
58 #include "G4Box.hh"
59 #include "G4Cons.hh"
60 #include "G4Tubs.hh"
61 #include "G4Trd.hh"
62 #include "G4Trap.hh"
63 #include "G4Sphere.hh"
64 #include "G4Para.hh"
65 #include "G4Torus.hh"
66 #include "G4Polycone.hh"
67 #include "G4Polyhedra.hh"
68 #include <deque>
69 
70 class G4VPhysicalVolume;
71 class G4LogicalVolume;
73 class G4Material;
74 
76 
77 public:
79  virtual ~G4PhysicalVolumeMassScene ();
80 
81 public: // With description
82 
83  G4double GetVolume () const {return fVolume;}
84  // Overall volume.
85 
86  G4double GetMass () const {return fMass;}
87  // Mass of whole tree, i.e., accounting for all daughters.
88 
89  void Reset ();
90  // Reset for subsequent re-use.
91 
92 public:
93 
94  // Force execution of AccrueMass for all solids...
95  void PreAddSolid (const G4Transform3D&, const G4VisAttributes&) {}
96  void PostAddSolid () {}
97  void AddSolid (const G4Box& solid) {AccrueMass (solid);}
98  void AddSolid (const G4Cons & solid) {AccrueMass (solid);}
99  void AddSolid (const G4Tubs& solid) {AccrueMass (solid);}
100  void AddSolid (const G4Trd& solid) {AccrueMass (solid);}
101  void AddSolid (const G4Trap& solid) {AccrueMass (solid);}
102  void AddSolid (const G4Sphere& solid) {AccrueMass (solid);}
103  void AddSolid (const G4Para& solid) {AccrueMass (solid);}
104  void AddSolid (const G4Torus& solid) {AccrueMass (solid);}
105  void AddSolid (const G4Polycone& solid) {AccrueMass (solid);}
106  void AddSolid (const G4Polyhedra& solid) {AccrueMass (solid);}
107  void AddSolid (const G4VSolid& solid) {AccrueMass (solid);}
108  void AddCompound (const G4VTrajectory&) {}
109  void AddCompound (const G4VHit&) {}
110  void AddCompound (const G4VDigi&) {}
112 
114  // Functions not used but required by the abstract interface.
115 
116  virtual void BeginPrimitives (const G4Transform3D&) {}
117  virtual void EndPrimitives () {}
118  virtual void BeginPrimitives2D (const G4Transform3D&) {}
119  virtual void EndPrimitives2D () {}
120  virtual void AddPrimitive (const G4Polyline&) {}
121  virtual void AddPrimitive (const G4Scale&) {}
122  virtual void AddPrimitive (const G4Text&) {}
123  virtual void AddPrimitive (const G4Circle&) {}
124  virtual void AddPrimitive (const G4Square&) {}
125  virtual void AddPrimitive (const G4Polymarker&) {}
126  virtual void AddPrimitive (const G4Polyhedron&) {}
127 
128 private:
129  void AccrueMass (const G4VSolid&);
137  std::deque<G4double> fDensityStack;
138 };
139 
140 #endif
virtual void AddPrimitive(const G4Square &)
Definition: G4Para.hh:77
Definition: G4Text.hh:73
G4PhysicalVolumeMassScene(G4PhysicalVolumeModel *)
void AddSolid(const G4Box &solid)
Definition: G4Box.hh:64
virtual void AddPrimitive(const G4Text &)
void AddSolid(const G4Polyhedra &solid)
Definition: G4Tubs.hh:85
void AddSolid(const G4Para &solid)
Definition: G4VHit.hh:48
Definition: G4Trd.hh:72
virtual void AddPrimitive(const G4Polyline &)
int G4int
Definition: G4Types.hh:78
void PreAddSolid(const G4Transform3D &, const G4VisAttributes &)
void AddCompound(const G4VTrajectory &)
void AddCompound(const G4THitsMap< G4double > &)
void AddSolid(const G4Trd &solid)
void AddSolid(const G4Cons &solid)
void AddSolid(const G4Sphere &solid)
virtual void BeginPrimitives2D(const G4Transform3D &)
Definition: G4Cons.hh:83
HepGeom::Transform3D G4Transform3D
virtual void BeginPrimitives(const G4Transform3D &)
void AddSolid(const G4Trap &solid)
virtual void AddPrimitive(const G4Polymarker &)
virtual void AddPrimitive(const G4Scale &)
void AddSolid(const G4Tubs &solid)
virtual void AddPrimitive(const G4Circle &)
void AddSolid(const G4VSolid &solid)
virtual void AddPrimitive(const G4Polyhedron &)
double G4double
Definition: G4Types.hh:76
std::deque< G4double > fDensityStack
void AddSolid(const G4Torus &solid)
void AddSolid(const G4Polycone &solid)