source: trunk/examples/extended/electromagnetic/TestEm15/src/RunAction.cc @ 1337

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

tag geant4.9.4 beta 1 + modifs locales

File size: 8.8 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 2006/06/29 16:47:13 gunter 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
43#include "Randomize.hh"
44#include <iomanip>
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47
48RunAction::RunAction(DetectorConstruction* det, PrimaryGeneratorAction* prim,
49                     HistoManager* histo)
50  : detector(det), primary(prim), ProcCounter(0), histoManager(histo)
51{ }
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54
55RunAction::~RunAction()
56{ }
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60void RunAction::BeginOfRunAction(const G4Run* aRun)
61{ 
62  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
63 
64  // save Rndm status
65  G4RunManager::GetRunManager()->SetRandomNumberStore(false);
66  CLHEP::HepRandom::showEngineStatus();
67
68  ProcCounter = new ProcessesCount;
69  totalCount = 0;
70 
71  truePL = truePL2 = geomPL = geomPL2 = 0.;
72  lDispl = lDispl2 = psiSpa = psiSpa2 = 0.;
73  tetPrj = tetPrj2 = 0.;
74  phiCor = phiCor2 = 0.;
75 
76  histoManager->book();
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
81void RunAction::CountProcesses(G4String procName)
82{
83   //does the process  already encounted ?
84   size_t nbProc = ProcCounter->size();
85   size_t i = 0;
86   while ((i<nbProc)&&((*ProcCounter)[i]->GetName()!=procName)) i++;
87   if (i == nbProc) ProcCounter->push_back( new OneProcessCount(procName));
88
89   (*ProcCounter)[i]->Count();
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
94void RunAction::EndOfRunAction(const G4Run* aRun)
95{
96  G4int NbOfEvents = aRun->GetNumberOfEvent();
97  if (NbOfEvents == 0) return;
98 
99  G4int  prec = G4cout.precision(5);
100   
101  G4Material* material = detector->GetMaterial();
102  G4double density = material->GetDensity();
103   
104  G4ParticleDefinition* particle = 
105                            primary->GetParticleGun()->GetParticleDefinition();
106  G4String Particle = particle->GetParticleName();   
107  G4double energy = primary->GetParticleGun()->GetParticleEnergy();
108  G4cout << "\n The run consists of " << NbOfEvents << " "<< Particle << " of "
109         << G4BestUnit(energy,"Energy") << " through " 
110         << G4BestUnit(detector->GetBoxSize(),"Length") << " of "
111         << material->GetName() << " (density: " 
112         << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
113 
114  //frequency of processes
115  G4cout << "\n Process calls frequency --->";
116  for (size_t i=0; i< ProcCounter->size();i++) {
117     G4String procName = (*ProcCounter)[i]->GetName();
118     G4int    count    = (*ProcCounter)[i]->GetCounter(); 
119     G4cout << "\t" << procName << " = " << count;
120  }
121 
122  if (totalCount == 0) return;
123 
124  //compute path length and related quantities
125  //
126  G4double MeanTPL  = truePL /totalCount;     
127  G4double MeanTPL2 = truePL2/totalCount;     
128  G4double rmsTPL   = std::sqrt(std::fabs(MeanTPL2 - MeanTPL*MeanTPL));
129 
130  G4double MeanGPL  = geomPL /totalCount;     
131  G4double MeanGPL2 = geomPL2/totalCount;     
132  G4double rmsGPL   = std::sqrt(std::fabs(MeanGPL2 - MeanGPL*MeanGPL));
133 
134  G4double MeanLaD  = lDispl /totalCount;     
135  G4double MeanLaD2 = lDispl2/totalCount;     
136  G4double rmsLaD   = std::sqrt(std::fabs(MeanLaD2 - MeanLaD*MeanLaD));
137 
138  G4double MeanPsi  = psiSpa /(totalCount);     
139  G4double MeanPsi2 = psiSpa2/(totalCount);     
140  G4double rmsPsi   = std::sqrt(std::fabs(MeanPsi2 - MeanPsi*MeanPsi));
141 
142  G4double MeanTeta  = tetPrj /(2*totalCount);     
143  G4double MeanTeta2 = tetPrj2/(2*totalCount);     
144  G4double rmsTeta   = std::sqrt(std::fabs(MeanTeta2 - MeanTeta*MeanTeta));
145 
146  G4double MeanCorrel  = phiCor /(totalCount);     
147  G4double MeanCorrel2 = phiCor2/(totalCount);     
148  G4double rmsCorrel = std::sqrt(std::fabs(MeanCorrel2-MeanCorrel*MeanCorrel));
149           
150  G4cout << "\n\n truePathLength :\t" << G4BestUnit(MeanTPL,"Length")
151         << " +- "                    << G4BestUnit( rmsTPL,"Length")
152         <<   "\n geomPathLength :\t" << G4BestUnit(MeanGPL,"Length")
153         << " +- "                    << G4BestUnit( rmsGPL,"Length")
154         <<   "\n lateralDisplac :\t" << G4BestUnit(MeanLaD,"Length")
155         << " +- "                    << G4BestUnit( rmsLaD,"Length")
156         <<   "\n Psi            :\t" << MeanPsi/mrad << " mrad"
157         << " +- "                    << rmsPsi /mrad << " mrad"
158         <<   "  ("                   << MeanPsi/deg  << " deg"
159         << " +- "                    << rmsPsi /deg  << " deg)"
160         << G4endl;
161                         
162  G4cout <<   "\n Theta_plane    :\t" << rmsTeta/mrad << " mrad"
163         <<   "  ("                   << rmsTeta/deg  << " deg)"
164         <<   "\n phi correlation:\t" << MeanCorrel
165         << " +- "                    << rmsCorrel
166         << "  (std::cos(phi_pos - phi_dir))"           
167         << G4endl;
168         
169
170  //cross check from G4EmCalculator
171  //
172  G4cout << "\n Verification from G4EmCalculator. \n";
173 
174  G4EmCalculator emCal;
175 
176  //get transport mean free path (for multiple scattering)
177  G4double MSmfp = emCal.GetMeanFreePath(energy,particle,"msc",material);
178   
179  //get range from restricted dedx
180  G4double range = emCal.GetRangeFromRestricteDEDX(energy,particle,material);
181 
182  //effective facRange
183  G4double efFacrange = MeanTPL/std::max(MSmfp, range);
184  if (MeanTPL/range >= 0.99) efFacrange = 1.;
185
186  G4cout << "\n transport mean free path :\t" << G4BestUnit(MSmfp,"Length")
187         << "\n range from restrict dE/dx:\t" << G4BestUnit(range,"Length")
188         << "\n ---> effective facRange  :\t" << efFacrange
189         << G4endl;
190
191  G4cout << "\n compute theta0 from Highland :\t"
192         << ComputeMscHighland(MeanTPL)/mrad << " mrad" 
193         << "  (" << ComputeMscHighland(MeanTPL)/deg << " deg)" 
194         << G4endl;
195                         
196  //restore default format       
197  G4cout.precision(prec);         
198
199  // delete and remove all contents in ProcCounter
200  while (ProcCounter->size()>0){
201    OneProcessCount* aProcCount=ProcCounter->back();
202    ProcCounter->pop_back();
203    delete aProcCount;
204  }
205  delete ProcCounter;
206 
207  histoManager->save();
208
209  // show Rndm status
210  CLHEP::HepRandom::showEngineStatus();
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214
215G4double RunAction::ComputeMscHighland(G4double pathLength)
216{
217 //compute the width of the Gaussian central part of the MultipleScattering
218 //projected angular distribution.
219 //Eur. Phys. Jour. C15 (2000) page 166, formule 23.9
220
221 G4double t = pathLength/(detector->GetMaterial()->GetRadlen());
222 if (t < DBL_MIN) return 0.;
223
224 G4ParticleGun* particle = primary->GetParticleGun();
225 G4double T = particle->GetParticleEnergy();
226 G4double M = particle->GetParticleDefinition()->GetPDGMass();
227 G4double z = std::abs(particle->GetParticleDefinition()->GetPDGCharge()/eplus);
228
229 G4double bpc = T*(T+2*M)/(T+M);
230 G4double teta0 = 13.6*MeV*z*std::sqrt(t)*(1.+0.038*std::log(t))/bpc;
231 return teta0;
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.