source: trunk/examples/extended/electromagnetic/TestEm13/src/RunAction.cc@ 811

Last change on this file since 811 was 807, checked in by garnier, 17 years ago

update

File size: 6.7 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// $Id: RunAction.cc,v 1.6 2006/06/29 16:44:47 gunter Exp $
27// GEANT4 tag $Name: $
28//
29//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
31
32#include "RunAction.hh"
33
34#include "DetectorConstruction.hh"
35#include "PrimaryGeneratorAction.hh"
36
37#include "G4Run.hh"
38#include "G4RunManager.hh"
39#include "G4UnitsTable.hh"
40#include "G4EmCalculator.hh"
41#include "G4Gamma.hh"
42
43#include "Randomize.hh"
44#include <iomanip>
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47
48RunAction::RunAction(DetectorConstruction* det, PrimaryGeneratorAction* prim)
49 : detector(det), primary(prim), ProcCounter(0)
50{ }
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
54RunAction::~RunAction()
55{ }
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
59void RunAction::BeginOfRunAction(const G4Run* aRun)
60{
61 G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
62
63 // save Rndm status
64 G4RunManager::GetRunManager()->SetRandomNumberStore(false);
65 CLHEP::HepRandom::showEngineStatus();
66
67 ProcCounter = new ProcessesCount;
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
72void RunAction::CountProcesses(G4String procName)
73{
74 //does the process already encounted ?
75 size_t nbProc = ProcCounter->size();
76 size_t i = 0;
77 while ((i<nbProc)&&((*ProcCounter)[i]->GetName()!=procName)) i++;
78 if (i == nbProc) ProcCounter->push_back( new OneProcessCount(procName));
79
80 (*ProcCounter)[i]->Count();
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
85void RunAction::EndOfRunAction(const G4Run* aRun)
86{
87 G4int NbOfEvents = aRun->GetNumberOfEvent();
88 if (NbOfEvents == 0) return;
89
90 G4int prec = G4cout.precision(5);
91
92 G4Material* material = detector->GetMaterial();
93 G4double density = material->GetDensity();
94 G4double tickness = detector->GetSize();
95
96 G4ParticleDefinition* particle =
97 primary->GetParticleGun()->GetParticleDefinition();
98 G4String Particle = particle->GetParticleName();
99 G4double energy = primary->GetParticleGun()->GetParticleEnergy();
100 G4cout << "\n The run consists of " << NbOfEvents << " "<< Particle << " of "
101 << G4BestUnit(energy,"Energy") << " through "
102 << G4BestUnit(tickness,"Length") << " of "
103 << material->GetName() << " (density: "
104 << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
105
106 //frequency of processes
107 G4int totalCount = 0;
108 G4int survive = 0;
109 G4cout << "\n Process calls frequency --->";
110 for (size_t i=0; i< ProcCounter->size();i++) {
111 G4String procName = (*ProcCounter)[i]->GetName();
112 G4int count = (*ProcCounter)[i]->GetCounter();
113 totalCount += count;
114 G4cout << "\t" << procName << " = " << count;
115 if (procName == "Transportation") survive = count;
116 }
117 G4cout << G4endl;
118 if (totalCount == 0) return;
119
120 G4double ratio = double(survive)/totalCount;
121
122 G4cout << "\n Nb of incident particles unaltered after "
123 << G4BestUnit(tickness,"Length") << " of "
124 << material->GetName() << " : " << survive
125 << " over " << totalCount << " incident particles."
126 << " Ratio = " << 100*ratio << " %" << G4endl;
127
128 if (ratio == 0.) return;
129
130 //compute cross section and related quantities
131 //
132 G4double CrossSection = std::log(1./ratio)/tickness;
133 G4double massicCS = CrossSection/density;
134
135 G4cout << " ---> CrossSection per volume:\t" << CrossSection*cm << " cm^-1 "
136 << "\tCrossSection per mass: " << G4BestUnit(massicCS, "Surface/Mass")
137 << G4endl;
138
139 //check cross section from G4EmCalculator
140 //
141 G4cout << "\n Verification from G4EmCalculator: \n";
142 G4EmCalculator emCalculator;
143 G4double sumc = 0.0;
144 for (size_t i=0; i< ProcCounter->size();i++) {
145 G4String procName = (*ProcCounter)[i]->GetName();
146 G4double massSigma =
147 emCalculator.GetCrossSectionPerVolume(energy,particle,
148 procName,material)/density;
149 if (particle == G4Gamma::Gamma())
150 massSigma =
151 emCalculator.ComputeCrossSectionPerVolume(energy,particle,
152 procName,material)/density;
153 sumc += massSigma;
154 if (procName != "Transportation")
155 G4cout << "\t" << procName << "= "
156 << G4BestUnit(massSigma, "Surface/Mass");
157 }
158 G4cout << "\ttotal= "
159 << G4BestUnit(sumc, "Surface/Mass") << G4endl;
160
161 //expected ratio of transmitted particles
162 G4double Ratio = std::exp(-sumc*density*tickness);
163 G4cout << "\tExpected ratio of transmitted particles= "
164 << 100*Ratio << " %" << G4endl;
165
166 //restore default format
167 G4cout.precision(prec);
168
169 // delete and remove all contents in ProcCounter
170 while (ProcCounter->size()>0){
171 OneProcessCount* aProcCount=ProcCounter->back();
172 ProcCounter->pop_back();
173 delete aProcCount;
174 }
175 delete ProcCounter;
176
177 // show Rndm status
178 CLHEP::HepRandom::showEngineStatus();
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.