Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Run.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 //
28 //
29 // $Id: Run.cc 71376 2013-06-14 07:44:50Z maire $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "Run.hh"
35 #include "DetectorConstruction.hh"
36 #include "PrimaryGeneratorAction.hh"
37 
38 #include "G4UnitsTable.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "G4EmCalculator.hh"
41 #include "G4Gamma.hh"
42 
43 #include <iomanip>
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 : G4Run(),
49  fDetector(det),
50  fParticle(0), fEkin(0.),
51  fTotalCount(0), fSumTrack(0.), fSumTrack2(0.), fEnTransfer(0.)
52 { }
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
56 Run::~Run()
57 { }
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 
62 {
63  fParticle = particle;
64  fEkin = energy;
65 }
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
68 void Run::CountProcesses(G4String procName)
69 {
70  std::map<G4String,G4int>::iterator it = fProcCounter.find(procName);
71  if ( it == fProcCounter.end()) {
72  fProcCounter[procName] = 1;
73  }
74  else {
75  fProcCounter[procName]++;
76  }
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
81 void Run::SumTrack (G4double track)
82 {
83  fTotalCount++;
84  fSumTrack += track;
85  fSumTrack2 += track*track;
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89 
91 {
92  fEnTransfer += energy;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
97 void Run::Merge(const G4Run* run)
98 {
99  const Run* localRun = static_cast<const Run*>(run);
100 
101  // pass information about primary particle
102  fParticle = localRun->fParticle;
103  fEkin = localRun->fEkin;
104 
105  //map: processes count
106  std::map<G4String,G4int>::const_iterator it;
107  for (it = localRun->fProcCounter.begin();
108  it !=localRun->fProcCounter.end(); ++it) {
109 
110  G4String procName = it->first;
111  G4int localCount = it->second;
112  if ( fProcCounter.find(procName) == fProcCounter.end()) {
113  fProcCounter[procName] = localCount;
114  }
115  else {
116  fProcCounter[procName] += localCount;
117  }
118  }
119 
120  fTotalCount += localRun->fTotalCount;
121  fSumTrack += localRun->fSumTrack;
122  fSumTrack2 += localRun->fSumTrack2;
123  fEnTransfer += localRun->fEnTransfer;
124 
125  G4Run::Merge(run);
126 }
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129 
130 void Run::EndOfRun()
131 {
132  G4int prec = 5;
133  G4int dfprec = G4cout.precision(prec);
134 
135  //run condition
136  //
137  G4String partName = fParticle->GetParticleName();
138  G4Material* material = fDetector->GetMaterial();
139  G4double density = material->GetDensity();
140  G4double tickness = fDetector->GetSize();
141 
142  G4cout << "\n ======================== run summary ======================\n";
143  G4cout << "\n The run is: " << numberOfEvent << " " << partName << " of "
144  << G4BestUnit(fEkin,"Energy") << " through "
145  << G4BestUnit(tickness,"Length") << " of "
146  << material->GetName() << " (density: "
147  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
148 
149  //frequency of processes
150  G4int survive = 0;
151  G4cout << "\n Process calls frequency --->";
152  std::map<G4String,G4int>::iterator it;
153  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
154  G4String procName = it->first;
155  G4int count = it->second;
156  G4cout << "\t" << procName << " = " << count;
157  if (procName == "Transportation") survive = count;
158  }
159 
160  if (survive > 0) {
161  G4cout << "\n\n Nb of incident particles surviving after "
162  << G4BestUnit(fDetector->GetSize(),"Length") << " of "
163  << material->GetName() << " : " << survive << G4endl;
164  }
165 
166  if (fTotalCount == 0) fTotalCount = 1; //force printing anyway
167 
168  //compute mean free path and related quantities
169  //
170  G4double MeanFreePath = fSumTrack /fTotalCount;
171  G4double MeanTrack2 = fSumTrack2/fTotalCount;
172  G4double rms = std::sqrt(std::fabs(MeanTrack2 - MeanFreePath*MeanFreePath));
173  G4double CrossSection = 1./MeanFreePath;
174  G4double massicMFP = MeanFreePath*density;
175  G4double massicCS = 1./massicMFP;
176 
177  G4cout << "\n\n MeanFreePath:\t" << G4BestUnit(MeanFreePath,"Length")
178  << " +- " << G4BestUnit( rms,"Length")
179  << "\tmassic: " << G4BestUnit(massicMFP, "Mass/Surface")
180  << "\n CrossSection:\t" << CrossSection*cm << " cm^-1 "
181  << "\t\t\tmassic: " << G4BestUnit(massicCS, "Surface/Mass")
182  << G4endl;
183 
184  //compute energy transfer coefficient
185  //
186  G4double MeanTransfer = fEnTransfer/fTotalCount;
187  G4double massTransfCoef = massicCS*MeanTransfer/fEkin;
188 
189  G4cout << "\n mean energy of charged secondaries: "
190  << G4BestUnit(MeanTransfer, "Energy")
191  << "\tmass_energy_transfer coef: "
192  << G4BestUnit(massTransfCoef, "Surface/Mass")
193  << G4endl;
194 
195  //check cross section from G4EmCalculator
196  //
197  G4cout << "\n Verification : "
198  << "crossSections from G4EmCalculator \n";
199 
200  G4EmCalculator emCalculator;
201  G4double sumc = 0.0;
202  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
203  G4String procName = it->first;
204  G4double massSigma =
205  emCalculator.GetCrossSectionPerVolume(fEkin,fParticle,
206  procName,material)/density;
207  if (fParticle == G4Gamma::Gamma())
208  massSigma =
209  emCalculator.ComputeCrossSectionPerVolume(fEkin,fParticle,
210  procName,material)/density;
211  sumc += massSigma;
212  G4cout << "\t" << procName << "= "
213  << G4BestUnit(massSigma, "Surface/Mass");
214  }
215  G4cout << "\ttotal= "
216  << G4BestUnit(sumc, "Surface/Mass") << G4endl;
217 
218  // remove all contents in fProcCounter
219  fProcCounter.clear();
220 
221  //restore default format
222  G4cout.precision(dfprec);
223 }
224 
225 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int numberOfEvent
Definition: G4Run.hh:59
void CountProcesses(G4String procName)
Definition: Run.cc:72
Run()
Definition: Run.cc:43
G4int first(char) const
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetDensity() const
Definition: G4Material.hh:180
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
string material
Definition: eplot.py:19
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
void SumTrack(G4double track)
Definition: Run.cc:81
static constexpr double cm
Definition: G4SIunits.hh:119
Definition: G4Run.hh:46
void SumeTransf(G4double energy)
Definition: Run.cc:90
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double energy(const ThreeVector &p, const G4double m)
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
Detector construction class to define materials and geometry.
void EndOfRun()
Definition: Run.cc:147
#define G4endl
Definition: G4ios.hh:61
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
Definition: Run.cc:77
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
virtual void Merge(const G4Run *)
Definition: Run.cc:115
~Run()
Definition: Run.cc:72