Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicalVolumeMassScene.cc
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//
28//
29// John Allison 10th August 1998.
30// An artificial scene to find physical volumes.
31
33
34#include "G4VSolid.hh"
35#include "G4Vector3D.hh"
37#include "G4LogicalVolume.hh"
38#include "G4Polyhedron.hh"
39#include "G4Material.hh"
41#include "G4UnitsTable.hh"
42
44(G4PhysicalVolumeModel* pPVModel):
45 fpPVModel (pPVModel),
46 fVolume (0.),
47 fMass (0.),
48 fpLastPV (0),
49 fPVPCount (0),
50 fLastDepth (0),
51 fLastDensity (0.)
52{}
53
55
57{
58 fVolume = 0.;
59 fMass = 0.;
60 fpLastPV = 0;
61 fPVPCount = 0;
62 fLastDepth = 0;
63 fLastDensity = 0.;
64 fDensityStack.clear();
65}
66
67void G4PhysicalVolumeMassScene::ProcessVolume (const G4VSolid& solid)
68{
69 G4int currentDepth = fpPVModel->GetCurrentDepth();
70 G4VPhysicalVolume* pCurrentPV = fpPVModel->GetCurrentPV();
71 //G4LogicalVolume* pCurrentLV = fpPVModel->GetCurrentLV();
72 G4Material* pCurrentMaterial = fpPVModel->GetCurrentMaterial();
73
74 if (pCurrentPV != fpLastPV) {
75 fpLastPV = pCurrentPV;
76 fPVPCount = 0;
77 }
78
79 G4double currentVolume = ((G4VSolid&)solid).GetCubicVolume();
80 G4double currentDensity = pCurrentMaterial? pCurrentMaterial->GetDensity() : 0.;
81 /* Using G4Polyhedron... (gives slightly different answers on Tubs, e.g.).
82 G4Polyhedron* pPolyhedron = solid.GetPolyhedron();
83 if (!pPolyhedron) {
84 G4cout <<
85 "G4PhysicalVolumeMassScene::AccrueMass: WARNING:"
86 "\n No G4Polyhedron for" << solid.GetEntityType() <<
87 ". \"" << solid.GetName() << "\" will not be accounted."
88 "\n It will be as though not there, i.e., the density as its mother."
89 "\n Its daughters will still be found and accounted."
90 << G4endl;
91 currentVolume = 0.;
92 currentDensity = 0.;
93 }
94 */
95
96 if (currentDepth == 0) fVolume = currentVolume;
97
98 if (currentDepth > fLastDepth) {
99 fDensityStack.push_back (fLastDensity);
100 } else if (currentDepth < fLastDepth) {
101 fDensityStack.pop_back();
102 }
103 fLastDepth = currentDepth;
104 fLastDensity = currentDensity;
105 G4double motherDensity = 0.;
106 if (currentDepth > 0) motherDensity = fDensityStack.back();
107
108 G4double subtractedMass = currentVolume * motherDensity;
109 G4double addedMass = currentVolume * currentDensity;
110 fMass -= subtractedMass;
111 fMass += addedMass;
112 /* Debug
113 G4cout << "current vol = "
114 << G4BestUnit (currentVolume,"Volume")
115 << ", current density = "
116 << G4BestUnit (currentDensity, "Volumic Mass")
117 << ", mother density = "
118 << G4BestUnit (motherDensity, "Volumic Mass")
119 << G4endl;
120 G4cout << "Subtracted mass = " << G4BestUnit (subtractedMass, "Mass")
121 << ", added mass = " << G4BestUnit (addedMass, "Mass")
122 << ", new mass = " << G4BestUnit (fMass, "Mass")
123 << G4endl;
124 */
125 if (fMass < 0.) {
126 G4cout <<
127 "G4PhysicalVolumeMassScene::AccrueMass: WARNING:"
128 "\n Mass going negative for \""
129 << pCurrentPV->GetName() <<
130 "\", copy "
131 << pCurrentPV->GetCopyNo() <<
132 ". Larger than mother?"
133 << G4endl;
134 }
135}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetDensity() const
Definition: G4Material.hh:178
G4PhysicalVolumeMassScene(G4PhysicalVolumeModel *)
G4VPhysicalVolume * GetCurrentPV() const
G4Material * GetCurrentMaterial() const
virtual G4int GetCopyNo() const =0
const G4String & GetName() const