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