source: trunk/source/visualization/modeling/include/G4PhysicalVolumeMassScene.hh@ 1177

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

update to CVS

File size: 5.3 KB
Line 
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,v 1.10 2009/10/21 14:17:33 allison Exp $
28// GEANT4 tag $Name: $
29//
30//
31// John Allison 12th September 2004
32
33// Class Description:
34//
35// Calculates the mass of a geometry tree taking into account daughters
36// up to the depth specified in the G4PhysicalVolumeModel. Culling is
37// ignored so that all volumes are seen.
38//
39// Do not use this for a "parallel world" for which materials are not
40// defined. Use only for the material world.
41//
42// The calculation is quite tricky, since it involves subtracting the
43// mass of that part of the mother that is occupied by each daughter and
44// then adding the mass of the daughter, and so on down the heirarchy.
45//
46// Usage for a given G4PhysicalVolumeModel* pvModel:
47// G4PhysicalVolumeMassScene massScene(pvModel);
48// pvModel->DescribeYourselfTo (massScene);
49// G4double volume = massScene.GetVolume();
50// G4double mass = massScene.GetMass();
51// massScene.Reset();
52// See, for example, G4ASCIITreeSceneHandler::EndModeling().
53
54#ifndef G4PHYSICALVOLUMEMASSSCENE_HH
55#define G4PHYSICALVOLUMEMASSSCENE_HH
56
57#include "G4VGraphicsScene.hh"
58
59#include "G4Box.hh"
60#include "G4Cons.hh"
61#include "G4Tubs.hh"
62#include "G4Trd.hh"
63#include "G4Trap.hh"
64#include "G4Sphere.hh"
65#include "G4Para.hh"
66#include "G4Torus.hh"
67#include "G4Polycone.hh"
68#include "G4Polyhedra.hh"
69#include <deque>
70
71class G4VPhysicalVolume;
72class G4LogicalVolume;
73class G4PhysicalVolumeModel;
74class G4Material;
75
76class G4PhysicalVolumeMassScene: public G4VGraphicsScene {
77
78public:
79 G4PhysicalVolumeMassScene (G4PhysicalVolumeModel*);
80 virtual ~G4PhysicalVolumeMassScene ();
81
82public: // With description
83
84 G4double GetVolume () const {return fVolume;}
85 // Overall volume.
86
87 G4double GetMass () const {return fMass;}
88 // Mass of whole tree, i.e., accounting for all daughters.
89
90 void Reset ();
91 // Reset for subsequent re-use.
92
93public:
94
95 // Force execution of AccrueMass for all solids...
96 void PreAddSolid (const G4Transform3D&, const G4VisAttributes&) {}
97 void PostAddSolid () {}
98 void AddSolid (const G4Box& s) {AccrueMass (s);}
99 void AddSolid (const G4Cons & s) {AccrueMass (s);}
100 void AddSolid (const G4Tubs& s) {AccrueMass (s);}
101 void AddSolid (const G4Trd& s) {AccrueMass (s);}
102 void AddSolid (const G4Trap& s) {AccrueMass (s);}
103 void AddSolid (const G4Sphere& s) {AccrueMass (s);}
104 void AddSolid (const G4Para& s) {AccrueMass (s);}
105 void AddSolid (const G4Torus& s) {AccrueMass (s);}
106 void AddSolid (const G4Polycone& s) {AccrueMass (s);}
107 void AddSolid (const G4Polyhedra& s) {AccrueMass (s);}
108 void AddSolid (const G4VSolid& s) {AccrueMass (s);}
109 void AddCompound (const G4VTrajectory&) {}
110 void AddCompound (const G4VHit&) {}
111 void AddCompound (const G4THitsMap<G4double>&) {}
112
113 ////////////////////////////////////////////////////////////////
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 virtual void AddPrimitive (const G4NURBS&) {}
128
129private:
130 void AccrueMass (const G4VSolid&);
131 G4PhysicalVolumeModel* fpPVModel;
132 G4double fVolume;
133 G4double fMass;
134 G4VPhysicalVolume* fpLastPV;
135 G4int fPVPCount;
136 G4int fLastDepth;
137 G4double fLastDensity;
138 std::deque<G4double> fDensityStack;
139};
140
141#endif
Note: See TracBrowser for help on using the repository browser.