source: trunk/examples/extended/electromagnetic/TestEm14/src/RunAction.cc@ 1353

Last change on this file since 1353 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

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.5 2010/04/05 18:02:39 maire Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
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#include "HistoManager.hh"
37
38#include "G4Run.hh"
39#include "G4RunManager.hh"
40#include "G4UnitsTable.hh"
41#include "G4EmCalculator.hh"
42#include "G4Gamma.hh"
43
44#include "Randomize.hh"
45#include <iomanip>
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
49RunAction::RunAction(DetectorConstruction* det, PrimaryGeneratorAction* prim,
50 HistoManager* histo)
51 : detector(det), primary(prim), histoManager(histo)
52{ }
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55
56RunAction::~RunAction()
57{ }
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61void RunAction::BeginOfRunAction(const G4Run* aRun)
62{
63 G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
64
65 // save Rndm status
66 G4RunManager::GetRunManager()->SetRandomNumberStore(false);
67 CLHEP::HepRandom::showEngineStatus();
68
69 totalCount = 0;
70 sumTrack = sumTrack2 = 0.;
71 eTransfer = 0.;
72
73 histoManager->book();
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77
78void RunAction::EndOfRunAction(const G4Run* aRun)
79{
80 G4int NbOfEvents = aRun->GetNumberOfEvent();
81 if (NbOfEvents == 0) return;
82
83 G4int prec = G4cout.precision(5);
84
85 G4Material* material = detector->GetMaterial();
86 G4double density = material->GetDensity();
87 G4int survive = 0;
88
89 G4ParticleDefinition* particle =
90 primary->GetParticleGun()->GetParticleDefinition();
91 G4String Particle = particle->GetParticleName();
92 G4double energy = primary->GetParticleGun()->GetParticleEnergy();
93 G4cout << "\n The run consists of " << NbOfEvents << " "<< Particle << " of "
94 << G4BestUnit(energy,"Energy") << " through "
95 << G4BestUnit(detector->GetSize(),"Length") << " of "
96 << material->GetName() << " (density: "
97 << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
98
99 //frequency of processes
100 G4cout << "\n Process calls frequency --->";
101 std::map<G4String,G4int>::iterator it;
102 for (it = procCounter.begin(); it != procCounter.end(); it++) {
103 G4String procName = it->first;
104 G4int count = it->second;
105 G4cout << "\t" << procName << " = " << count;
106 if (procName == "Transportation") survive = count;
107 }
108
109 if (survive > 0) {
110 G4cout << "\n\n Nb of incident particles surviving after "
111 << G4BestUnit(detector->GetSize(),"Length") << " of "
112 << material->GetName() << " : " << survive << G4endl;
113 }
114
115 if (totalCount == 0) totalCount = 1; //force printing anyway
116
117 //compute mean free path and related quantities
118 //
119 G4double MeanFreePath = sumTrack /totalCount;
120 G4double MeanTrack2 = sumTrack2/totalCount;
121 G4double rms = std::sqrt(std::fabs(MeanTrack2 - MeanFreePath*MeanFreePath));
122 G4double CrossSection = 1./MeanFreePath;
123 G4double massicMFP = MeanFreePath*density;
124 G4double massicCS = 1./massicMFP;
125
126 G4cout << "\n\n MeanFreePath:\t" << G4BestUnit(MeanFreePath,"Length")
127 << " +- " << G4BestUnit( rms,"Length")
128 << "\tmassic: " << G4BestUnit(massicMFP, "Mass/Surface")
129 << "\n CrossSection:\t" << CrossSection*cm << " cm^-1 "
130 << "\t\t\tmassic: " << G4BestUnit(massicCS, "Surface/Mass")
131 << G4endl;
132
133 //compute energy transfer coefficient
134 //
135 G4double MeanTransfer = eTransfer/totalCount;
136 G4double massTransfCoef = massicCS*MeanTransfer/energy;
137
138 G4cout << "\n mean energy of charged secondaries: " << G4BestUnit(MeanTransfer, "Energy")
139 << "\tmass_energy_transfer coef: " << G4BestUnit(massTransfCoef, "Surface/Mass")
140 << G4endl;
141
142 //check cross section from G4EmCalculator
143 //
144 G4cout << "\n Verification : "
145 << "crossSections from G4EmCalculator \n";
146
147 G4EmCalculator emCalculator;
148 G4double sumc = 0.0;
149 for (it = procCounter.begin(); it != procCounter.end(); it++) {
150 G4String procName = it->first;
151 G4double massSigma =
152 emCalculator.GetCrossSectionPerVolume(energy,particle,
153 procName,material)/density;
154 if (particle == G4Gamma::Gamma())
155 massSigma =
156 emCalculator.ComputeCrossSectionPerVolume(energy,particle,
157 procName,material)/density;
158 sumc += massSigma;
159 G4cout << "\t" << procName << "= "
160 << G4BestUnit(massSigma, "Surface/Mass");
161 }
162 G4cout << "\ttotal= "
163 << G4BestUnit(sumc, "Surface/Mass") << G4endl;
164
165 //restore default format
166 G4cout.precision(prec);
167
168 // remove all contents in procCounter
169 procCounter.clear();
170
171 histoManager->save();
172
173 // show Rndm status
174 CLHEP::HepRandom::showEngineStatus();
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.