Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
examples/extended/electromagnetic/TestEm15/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/TestEm15/src/RunAction.cc
27 /// \brief Implementation of the RunAction class
28 //
29 // $Id: RunAction.cc 67268 2013-02-13 11:38:40Z ihrivnac $
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 
45 #include "Randomize.hh"
46 #include "G4SystemOfUnits.hh"
47 #include <iomanip>
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
52  : G4UserRunAction(),fDetector(det), fPrimary(prim), fProcCounter(0),
53  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
72  ////G4RunManager::GetRunManager()->SetRandomNumberStore(true);
74 
75  fProcCounter = new ProcessesCount;
76  fTotalCount = 0;
77 
78  fTruePL = fTruePL2 = fGeomPL = fGeomPL2 = 0.;
79  fLDispl = fLDispl2 = fPsiSpa = fPsiSpa2 = 0.;
80  fTetPrj = fTetPrj2 = 0.;
81  fPhiCor = fPhiCor2 = 0.;
82 
83  //histograms
84  //
85  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
86  if ( analysisManager->IsActive() ) {
87  analysisManager->OpenFile();
88  }
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 
94 {
95  //does the process already encounted ?
96  size_t nbProc = fProcCounter->size();
97  size_t i = 0;
98  while ((i<nbProc)&&((*fProcCounter)[i]->GetName()!=procName)) i++;
99  if (i == nbProc) fProcCounter->push_back( new OneProcessCount(procName));
100 
101  (*fProcCounter)[i]->Count();
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
106 void RunAction::EndOfRunAction(const G4Run* aRun)
107 {
108  G4int NbOfEvents = aRun->GetNumberOfEvent();
109  if (NbOfEvents == 0) return;
110 
111  G4int prec = G4cout.precision(5);
112 
113  G4Material* material = fDetector->GetMaterial();
114  G4double density = material->GetDensity();
115 
116  G4ParticleDefinition* particle =
117  fPrimary->GetParticleGun()->GetParticleDefinition();
118  G4String Particle = particle->GetParticleName();
120  G4cout << "\n The run consists of " << NbOfEvents << " "<< Particle << " of "
121  << G4BestUnit(energy,"Energy") << " through "
122  << G4BestUnit(fDetector->GetBoxSize(),"Length") << " of "
123  << material->GetName() << " (density: "
124  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
125 
126  //frequency of processes
127  G4cout << "\n Process calls frequency --->";
128  for (size_t i=0; i< fProcCounter->size();i++) {
129  G4String procName = (*fProcCounter)[i]->GetName();
130  G4int count = (*fProcCounter)[i]->GetCounter();
131  G4cout << "\t" << procName << " = " << count;
132  }
133 
134  if (fTotalCount == 0) return;
135 
136  //compute path length and related quantities
137  //
138  G4double MeanTPL = fTruePL /fTotalCount;
139  G4double MeanTPL2 = fTruePL2/fTotalCount;
140  G4double rmsTPL = std::sqrt(std::fabs(MeanTPL2 - MeanTPL*MeanTPL));
141 
142  G4double MeanGPL = fGeomPL /fTotalCount;
143  G4double MeanGPL2 = fGeomPL2/fTotalCount;
144  G4double rmsGPL = std::sqrt(std::fabs(MeanGPL2 - MeanGPL*MeanGPL));
145 
146  G4double MeanLaD = fLDispl /fTotalCount;
147  G4double MeanLaD2 = fLDispl2/fTotalCount;
148  G4double rmsLaD = std::sqrt(std::fabs(MeanLaD2 - MeanLaD*MeanLaD));
149 
150  G4double MeanPsi = fPsiSpa /(fTotalCount);
151  G4double MeanPsi2 = fPsiSpa2/(fTotalCount);
152  G4double rmsPsi = std::sqrt(std::fabs(MeanPsi2 - MeanPsi*MeanPsi));
153 
154  G4double MeanTeta = fTetPrj /(2*fTotalCount);
155  G4double MeanTeta2 = fTetPrj2/(2*fTotalCount);
156  G4double rmsTeta = std::sqrt(std::fabs(MeanTeta2 - MeanTeta*MeanTeta));
157 
158  G4double MeanCorrel = fPhiCor /(fTotalCount);
159  G4double MeanCorrel2 = fPhiCor2/(fTotalCount);
160  G4double rmsCorrel = std::sqrt(std::fabs(MeanCorrel2-MeanCorrel*MeanCorrel));
161 
162  G4cout << "\n\n truePathLength :\t" << G4BestUnit(MeanTPL,"Length")
163  << " +- " << G4BestUnit( rmsTPL,"Length")
164  << "\n geomPathLength :\t" << G4BestUnit(MeanGPL,"Length")
165  << " +- " << G4BestUnit( rmsGPL,"Length")
166  << "\n lateralDisplac :\t" << G4BestUnit(MeanLaD,"Length")
167  << " +- " << G4BestUnit( rmsLaD,"Length")
168  << "\n Psi :\t" << MeanPsi/mrad << " mrad"
169  << " +- " << rmsPsi /mrad << " mrad"
170  << " (" << MeanPsi/deg << " deg"
171  << " +- " << rmsPsi /deg << " deg)"
172  << G4endl;
173 
174  G4cout << "\n Theta_plane :\t" << rmsTeta/mrad << " mrad"
175  << " (" << rmsTeta/deg << " deg)"
176  << "\n phi correlation:\t" << MeanCorrel
177  << " +- " << rmsCorrel
178  << " (std::cos(phi_pos - phi_dir))"
179  << G4endl;
180 
181 
182  //cross check from G4EmCalculator
183  //
184  G4cout << "\n Verification from G4EmCalculator. \n";
185 
186  G4EmCalculator emCal;
187 
188  //get transport mean free path (for multiple scattering)
189  G4double MSmfp = emCal.GetMeanFreePath(energy,particle,"msc",material);
190 
191  //get range from restricted dedx
192  G4double range = emCal.GetRangeFromRestricteDEDX(energy,particle,material);
193 
194  //effective facRange
195  G4double efFacrange = MeanTPL/std::max(MSmfp, range);
196  if (MeanTPL/range >= 0.99) efFacrange = 1.;
197 
198  G4cout << "\n transport mean free path :\t" << G4BestUnit(MSmfp,"Length")
199  << "\n range from restrict dE/dx:\t" << G4BestUnit(range,"Length")
200  << "\n ---> effective facRange :\t" << efFacrange
201  << G4endl;
202 
203  G4cout << "\n compute theta0 from Highland :\t"
204  << ComputeMscHighland(MeanTPL)/mrad << " mrad"
205  << " (" << ComputeMscHighland(MeanTPL)/deg << " deg)"
206  << G4endl;
207 
208  //restore default format
209  G4cout.precision(prec);
210 
211  // delete and remove all contents in fProcCounter
212  while (fProcCounter->size()>0){
213  OneProcessCount* aProcCount=fProcCounter->back();
214  fProcCounter->pop_back();
215  delete aProcCount;
216  }
217  delete fProcCounter;
218 
219  //save histograms
220  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
221  if ( analysisManager->IsActive() ) {
222  analysisManager->Write();
223  analysisManager->CloseFile();
224  }
225 
226  // show Rndm status
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231 
233 {
234  //compute the width of the Gaussian central part of the MultipleScattering
235  //projected angular distribution.
236  //Eur. Phys. Jour. C15 (2000) page 166, formule 23.9
237 
238  G4double t = pathLength/(fDetector->GetMaterial()->GetRadlen());
239  if (t < DBL_MIN) return 0.;
240 
241  G4ParticleGun* particle = fPrimary->GetParticleGun();
242  G4double T = particle->GetParticleEnergy();
243  G4double M = particle->GetParticleDefinition()->GetPDGMass();
244  G4double z = std::abs(particle->GetParticleDefinition()->GetPDGCharge()/eplus);
245 
246  G4double bpc = T*(T+2*M)/(T+M);
247  G4double teta0 = 13.6*MeV*z*std::sqrt(t)*(1.+0.038*std::log(t))/bpc;
248  return teta0;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< OneProcessCount * > ProcessesCount
G4double z
Definition: TRTMaterials.hh:39
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
G4double GetMeanFreePath(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
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4int GetRunID() const
Definition: G4Run.hh:76
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
G4double GetRadlen() const
Definition: G4Material.hh:218
static void showEngineStatus()
Definition: Random.cc:203
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define DBL_MIN
Definition: templates.hh:75
G4ParticleDefinition * GetParticleDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double GetParticleEnergy() const