Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
examples/extended/electromagnetic/TestEm14/src/RunAction.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 /// \file electromagnetic/TestEm14/src/RunAction.cc
27 /// \brief Implementation of the RunAction class
28 //
29 // $Id: RunAction.cc 68313 2013-03-21 18:15:21Z maire $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "RunAction.hh"
35 
36 #include "DetectorConstruction.hh"
37 #include "PrimaryGeneratorAction.hh"
38 #include "HistoManager.hh"
39 
40 #include "G4Run.hh"
41 #include "G4RunManager.hh"
42 #include "G4UnitsTable.hh"
43 #include "G4EmCalculator.hh"
44 #include "G4Gamma.hh"
45 
46 #include "G4SystemOfUnits.hh"
47 #include "Randomize.hh"
48 #include <iomanip>
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
53  : G4UserRunAction(),fDetector(det), fPrimary(prim), fHistoManager(0)
54 {
55  fHistoManager = new HistoManager();
56 }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61 {
62  delete fHistoManager;
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
67 void RunAction::BeginOfRunAction(const G4Run* aRun)
68 {
69  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
70 
71  // save Rndm status
74 
75  fTotalCount = 0;
76  fSumTrack = fSumTrack2 = 0.;
77  fEnTransfer = 0.;
78 
79  //histograms
80  //
81  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
82  if ( analysisManager->IsActive() ) {
83  analysisManager->OpenFile();
84  }
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
89 void RunAction::EndOfRunAction(const G4Run* aRun)
90 {
91  G4int NbOfEvents = aRun->GetNumberOfEvent();
92  if (NbOfEvents == 0) return;
93 
94  G4int prec = G4cout.precision(5);
95 
96  G4Material* material = fDetector->GetMaterial();
97  G4double density = material->GetDensity();
98  G4int survive = 0;
99 
100  G4ParticleDefinition* particle =
101  fPrimary->GetParticleGun()->GetParticleDefinition();
102  G4String Particle = particle->GetParticleName();
104  G4cout << "\n The run consists of " << NbOfEvents << " "<< Particle << " of "
105  << G4BestUnit(energy,"Energy") << " through "
106  << G4BestUnit(fDetector->GetSize(),"Length") << " of "
107  << material->GetName() << " (density: "
108  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
109 
110  //frequency of processes
111  G4cout << "\n Process calls frequency --->";
112  std::map<G4String,G4int>::iterator it;
113  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
114  G4String procName = it->first;
115  G4int count = it->second;
116  G4cout << "\t" << procName << " = " << count;
117  if (procName == "Transportation") survive = count;
118  }
119 
120  if (survive > 0) {
121  G4cout << "\n\n Nb of incident particles surviving after "
122  << G4BestUnit(fDetector->GetSize(),"Length") << " of "
123  << material->GetName() << " : " << survive << G4endl;
124  }
125 
126  if (fTotalCount == 0) fTotalCount = 1; //force printing anyway
127 
128  //compute mean free path and related quantities
129  //
130  G4double MeanFreePath = fSumTrack /fTotalCount;
131  G4double MeanTrack2 = fSumTrack2/fTotalCount;
132  G4double rms = std::sqrt(std::fabs(MeanTrack2 - MeanFreePath*MeanFreePath));
133  G4double CrossSection = 1./MeanFreePath;
134  G4double massicMFP = MeanFreePath*density;
135  G4double massicCS = 1./massicMFP;
136 
137  G4cout << "\n\n MeanFreePath:\t" << G4BestUnit(MeanFreePath,"Length")
138  << " +- " << G4BestUnit( rms,"Length")
139  << "\tmassic: " << G4BestUnit(massicMFP, "Mass/Surface")
140  << "\n CrossSection:\t" << CrossSection*cm << " cm^-1 "
141  << "\t\t\tmassic: " << G4BestUnit(massicCS, "Surface/Mass")
142  << G4endl;
143 
144  //compute energy transfer coefficient
145  //
146  G4double MeanTransfer = fEnTransfer/fTotalCount;
147  G4double massTransfCoef = massicCS*MeanTransfer/energy;
148 
149  G4cout << "\n mean energy of charged secondaries: "
150  << G4BestUnit(MeanTransfer, "Energy")
151  << "\tmass_energy_transfer coef: "
152  << G4BestUnit(massTransfCoef, "Surface/Mass")
153  << G4endl;
154 
155  //check cross section from G4EmCalculator
156  //
157  G4cout << "\n Verification : "
158  << "crossSections from G4EmCalculator \n";
159 
160  G4EmCalculator emCalculator;
161  G4double sumc = 0.0;
162  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
163  G4String procName = it->first;
164  G4double massSigma =
165  emCalculator.GetCrossSectionPerVolume(energy,particle,
166  procName,material)/density;
167  if (particle == G4Gamma::Gamma())
168  massSigma =
169  emCalculator.ComputeCrossSectionPerVolume(energy,particle,
170  procName,material)/density;
171  sumc += massSigma;
172  G4cout << "\t" << procName << "= "
173  << G4BestUnit(massSigma, "Surface/Mass");
174  }
175  G4cout << "\ttotal= "
176  << G4BestUnit(sumc, "Surface/Mass") << G4endl;
177 
178  //restore default format
179  G4cout.precision(prec);
180 
181  // remove all contents in fProcCounter
182  fProcCounter.clear();
183 
184  //save histograms
185  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
186  if ( analysisManager->IsActive() ) {
187  analysisManager->Write();
188  analysisManager->CloseFile();
189  }
190 
191  // show Rndm status
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4int first(char) const
const G4String & GetName() const
Definition: G4Material.hh:176
G4double GetDensity() const
Definition: G4Material.hh:178
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void SetRandomNumberStore(G4bool flag)
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
string material
Definition: eplot.py:19
G4double density
Definition: TRTMaterials.hh:39
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4int GetNumberOfEvent() const
Definition: G4Run.hh:79
G4int GetRunID() const
Definition: G4Run.hh:76
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static void showEngineStatus()
Definition: Random.cc:203
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4ParticleDefinition * GetParticleDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetParticleEnergy() const