Geant4  10.01.p01
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: RunAction.cc 86443 2014-11-12 09:40:56Z gcosmo $
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 "G4ParticleDefinition.hh"
39 
40 #include "G4Run.hh"
41 #include "G4RunManager.hh"
42 #include "G4UnitsTable.hh"
43 #include "G4EmCalculator.hh"
44 
45 #include "G4Gamma.hh"
46 #include "G4Electron.hh"
47 #include "G4Positron.hh"
48 
49 #include "G4SystemOfUnits.hh"
50 #include "G4PhysicalConstants.hh"
51 #include <iomanip>
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56  : detector(det), primary(prim), ProcCounter(0), fAnalysisManager(0)
57 {
62 
63  BookHisto();
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
69 {}
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
73 void RunAction::BeginOfRunAction(const G4Run* aRun)
74 {
75  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
76 
77  // save Rndm status
78  // G4RunManager::GetRunManager()->SetRandomNumberStore(false);
79  // CLHEP::HepRandom::showEngineStatus();
80 
81  if (ProcCounter) delete ProcCounter;
83  totalEventCount = 0;
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90 
92  G4double kinEnergy, G4double costheta,
93  G4double phi,
94  G4double longitudinalPolarization)
95 {
96  G4int id = -1;
97  if (particle == fGamma) {
98  photonStats.FillData(kinEnergy, costheta, longitudinalPolarization);
99  if(fAnalysisManager) { id = 1; }
100  } else if (particle == fElectron) {
101  electronStats.FillData(kinEnergy, costheta, longitudinalPolarization);
102  if(fAnalysisManager) { id = 5; }
103  } else if (particle == fPositron) {
104  positronStats.FillData(kinEnergy, costheta, longitudinalPolarization);
105  if(fAnalysisManager) { id = 9; }
106  }
107  if(id > 0) {
108  fAnalysisManager->FillH1(id,kinEnergy,1.0);
109  fAnalysisManager->FillH1(id+1,costheta,1.0);
110  fAnalysisManager->FillH1(id+2,phi,1.0);
111  fAnalysisManager->FillH1(id+3,longitudinalPolarization,1.0);
112  }
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 
118 {
119  // Always creating analysis manager
120  fAnalysisManager = G4AnalysisManager::Instance();
121  fAnalysisManager->SetActivation(true);
122  fAnalysisManager->SetVerboseLevel(1);
123 
124  // Open file histogram file
125  fAnalysisManager->OpenFile("pol01");
126  fAnalysisManager->SetFirstHistoId(1);
127 
128  // Creating an 1-dimensional histograms in the root directory of the tree
129  const G4String id[] = { "h1", "h2", "h3", "h4", "h5",
130  "h6", "h7", "h8", "h9", "h10", "h11", "h12"};
131  const G4String title[] =
132  { "Gamma Energy distribution", //1
133  "Gamma Cos(Theta) distribution", //2
134  "Gamma Phi angular distribution", //3
135  "Gamma longitudinal Polarization", //4
136  "Electron Energy distribution", //5
137  "Electron Cos(Theta) distribution", //6
138  "Electron Phi angular distribution", //7
139  "Electron longitudinal Polarization", //8
140  "Positron Energy distribution", //9
141  "Positron Cos(Theta) distribution", //10
142  "Positron Phi angular distribution", //11
143  "Positron longitudinal Polarization" //12
144  };
145  G4double vmin, vmax;
146  G4int nbins = 120;
147  for(int i=0; i<12; ++i) {
148  G4int j = i - i/4*4;
149  if(0==j) { vmin = 0.; vmax = 12.*MeV; }
150  else if(1==j) { vmin = -1.; vmax = 1.; }
151  else if(2==j) { vmin = 0.; vmax = pi; }
152  else { vmin = -1.5; vmax = 1.5; }
153  G4int ih = fAnalysisManager->CreateH1(id[i],title[i],nbins,vmin,vmax);
154  fAnalysisManager->SetH1Activation(ih, false);
155  }
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
161 {
162  if(fAnalysisManager) {
163  G4double norm = 1.0/G4double(nevents);
164  for(int i=0; i<12; ++i) {
165  fAnalysisManager->ScaleH1(i, norm);
166  }
167  fAnalysisManager->Write();
168  fAnalysisManager->CloseFile();
169  delete fAnalysisManager;
170  fAnalysisManager = 0;
171  }
172 }
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175 
176 void RunAction::CountProcesses(G4String procName)
177 {
178  // is the process already counted ?
179  // *AS* change to std::map?!
180  size_t nbProc = ProcCounter->size();
181  size_t i = 0;
182  while ((i<nbProc)&&((*ProcCounter)[i]->GetName()!=procName)) i++;
183  if (i == nbProc) ProcCounter->push_back( new OneProcessCount(procName));
184 
185  (*ProcCounter)[i]->Count();
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 
190 void RunAction::EndOfRunAction(const G4Run* aRun)
191 {
192  G4int NbOfEvents = aRun->GetNumberOfEvent();
193  if (NbOfEvents == 0) return;
194 
195  G4int prec = G4cout.precision(5);
196 
197  G4Material* material = detector->GetMaterial();
198  G4double density = material->GetDensity();
199 
200  G4ParticleDefinition* particle =
202  G4String Particle = particle->GetParticleName();
204  G4cout << "\n The run consists of " << NbOfEvents << " "<< Particle << " of "
205  << G4BestUnit(energy,"Energy") << " through "
206  << G4BestUnit(detector->GetBoxSizeZ(),"Length") << " of "
207  << material->GetName() << " (density: "
208  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
209 
210  //frequency of processes
211  G4cout << "\n Process calls frequency --->\n";
212  for (size_t i=0; i< ProcCounter->size();i++) {
213  G4String procName = (*ProcCounter)[i]->GetName();
214  G4int count = (*ProcCounter)[i]->GetCounter();
215  G4cout << "\t" << procName << " = " << count<<"\n";
216  }
217 
218  if (totalEventCount == 0) return;
219 
220  G4cout<<" Gamma: \n";
222  G4cout<<" Electron: \n";
224  G4cout<<" Positron: \n";
226 
227  //cross check from G4EmCalculator
228  // G4cout << "\n Verification from G4EmCalculator. \n";
229  // G4EmCalculator emCal;
230 
231  //restore default format
232  G4cout.precision(prec);
233 
234  // write out histograms
235  SaveHisto(NbOfEvents);
236 
237  // show Rndm status
238  CLHEP::HepRandom::showEngineStatus();
239 }
240 
241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242 
244 {
245  ++totalEventCount;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252 
254  : currentNumber(0),
255  totalNumber(0), totalNumber2(0),
256  sumEnergy(0), sumEnergy2(0),
257  sumPolarization(0), sumPolarization2(0),
258  sumCosTheta(0), sumCosTheta2(0)
259 {}
260 
261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
262 
264 {}
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267 
269 {
270  totalNumber+=currentNumber;
271  totalNumber2+=currentNumber*currentNumber;
272  currentNumber=0;
273 }
274 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
275 
277  G4double costheta,
278  G4double longitudinalPolarization)
279 {
280  ++currentNumber;
281  sumEnergy+=kinEnergy;
282  sumEnergy2+=kinEnergy*kinEnergy;
283  sumPolarization+=longitudinalPolarization;
284  sumPolarization2+=longitudinalPolarization*longitudinalPolarization;
285  sumCosTheta+=costheta;
286  sumCosTheta2+=costheta*costheta;
287 }
288 
289 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
290 
292 {
293  G4cout<<"Mean Number per Event :"
294  <<G4double(totalNumber)/G4double(totalNumberOfEvents)<<"\n";
295  if (totalNumber==0) totalNumber=1;
296  G4double energyMean=sumEnergy/totalNumber;
297  G4double energyRms=std::sqrt(sumEnergy2/totalNumber-energyMean*energyMean);
298  G4cout<<"Mean Energy :"<< G4BestUnit(energyMean,"Energy")
299  <<" +- "<<G4BestUnit(energyRms,"Energy")<<"\n";
300  G4double polarizationMean=sumPolarization/totalNumber;
301  G4double polarizationRms=
302  std::sqrt(sumPolarization2/totalNumber-polarizationMean*polarizationMean);
303  G4cout<<"Mean Polarization :"<< polarizationMean
304  <<" +- "<<polarizationRms<<"\n";
305 }
306 
307 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
308 
310 {
311  currentNumber=0;
312  totalNumber=totalNumber2=0;
313  sumEnergy=sumEnergy2=0;
314  sumPolarization=sumPolarization2=0;
315  sumCosTheta=sumCosTheta2=0;
316 }
317 
318 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< OneProcessCount * > ProcessesCount
static const double MeV
Definition: G4SIunits.hh:193
void BeginOfRunAction(const G4Run *)
Definition: RunAction.cc:57
const G4ParticleDefinition * fPositron
Definition: RunAction.hh:93
const G4String & GetName() const
Definition: G4Material.hh:178
const G4double pi
G4double GetDensity() const
Definition: G4Material.hh:180
G4int totalEventCount
Definition: RunAction.hh:101
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
ParticleStatistics electronStats
Definition: RunAction.hh:104
ParticleStatistics photonStats
Definition: RunAction.hh:103
G4double density
Definition: TRTMaterials.hh:39
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
void FillData(G4double kinEnergy, G4double costheta, G4double longitudinalPolarization)
Definition: RunAction.cc:276
ParticleStatistics positronStats
Definition: RunAction.hh:105
G4int GetNumberOfEvent() const
Definition: G4Run.hh:79
void EndOfRunAction(const G4Run *)
Definition: RunAction.cc:260
void PrintResults(G4int totalNumberOfEvents)
Definition: RunAction.cc:291
G4int GetRunID() const
Definition: G4Run.hh:76
Definition: G4Run.hh:46
void EventFinished()
Definition: RunAction.cc:243
void CountProcesses(G4String)
Definition: RunAction.cc:91
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void BookHisto()
Definition: RunAction.cc:68
const G4ParticleDefinition * fElectron
Definition: RunAction.hh:92
~RunAction()
Definition: RunAction.cc:52
ProcessesCount * ProcCounter
Definition: RunAction.hh:97
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4AnalysisManager * fAnalysisManager
Definition: RunAction.hh:77
G4ParticleGun * GetParticleGun()
G4double energy(const ThreeVector &p, const G4double m)
PrimaryGeneratorAction * primary
Definition: RunAction.hh:60
G4ParticleDefinition * GetParticleDefinition() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
DetectorConstruction * detector
Definition: RunAction.hh:59
Detector construction class to demonstrate various ways of placement.
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void FillData(const G4ParticleDefinition *particle, G4double kinEnergy, G4double costheta, G4double phi, G4double longitudinalPolarization)
Definition: RunAction.cc:91
const G4ParticleDefinition * fGamma
Definition: RunAction.hh:91
G4double GetParticleEnergy() const
void SaveHisto(G4int nevents)
Definition: RunAction.cc:160