Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HistoManager.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 //---------------------------------------------------------------------------
32 //
33 // ClassName: HistoManager
34 //
35 //
36 // Author: V.Ivanchenko 30/01/01
37 //
38 // Modified:
39 // 04.06.2006 Adoptation of Hadr01 (V.Ivanchenko)
40 // 16.11.2006 Add beamFlag (V.Ivanchenko)
41 //
42 //----------------------------------------------------------------------------
43 //
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 #include "HistoManager.hh"
49 #include "G4Track.hh"
50 #include "G4Step.hh"
51 #include "G4UnitsTable.hh"
52 #include "G4Neutron.hh"
53 #include "G4Proton.hh"
54 #include "G4AntiProton.hh"
55 #include "G4Gamma.hh"
56 #include "G4Electron.hh"
57 #include "G4Positron.hh"
58 #include "G4MuonPlus.hh"
59 #include "G4MuonMinus.hh"
60 #include "G4PionPlus.hh"
61 #include "G4PionMinus.hh"
62 #include "G4PionZero.hh"
63 #include "G4KaonPlus.hh"
64 #include "G4KaonMinus.hh"
65 #include "G4KaonZeroShort.hh"
66 #include "G4KaonZeroLong.hh"
67 #include "G4Deuteron.hh"
68 #include "G4Triton.hh"
69 #include "G4He3.hh"
70 #include "G4Alpha.hh"
71 #include "Histo.hh"
72 #include "globals.hh"
73 #include "G4SystemOfUnits.hh"
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 
77 HistoManager* HistoManager::fManager = 0;
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82 {
83  if(!fManager) {
84  static HistoManager manager;
85  fManager = &manager;
86  }
87  return fManager;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {
94  fVerbose= 0;
95  fNSlices = 300;
96  fNBinsE = 100;
97  fNHisto = 25;
98  fLength = 300.*mm;
99  fEdepMax = 1.0*GeV;
100  fBeamFlag = true;
101  fHistoBooked = false;
102  fHisto = new Histo();
103  fHisto->SetVerbose(fVerbose);
104  fNeutron = G4Neutron::Neutron();
105  fPrimaryKineticEnergy = 0.0;
106  fPrimaryDef = 0;
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110 
112 {
113  delete fHisto;
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117 
119 {
120  fHistoBooked = true;
121  fHisto->Add1D("1","Energy deposition (MeV/mm/event) in the target",
122  fNSlices,0.0,fLength/mm,MeV/mm);
123  fHisto->Add1D("2","Log10 Energy (MeV) of gammas",fNBinsE,-5.,5.,1.0);
124  fHisto->Add1D("3","Log10 Energy (MeV) of electrons",fNBinsE,-5.,5.,1.0);
125  fHisto->Add1D("4","Log10 Energy (MeV) of positrons",fNBinsE,-5.,5.,1.0);
126  fHisto->Add1D("5","Log10 Energy (MeV) of protons",fNBinsE,-5.,5.,1.0);
127  fHisto->Add1D("6","Log10 Energy (MeV) of neutrons",fNBinsE,-5.,5.,1.0);
128  fHisto->Add1D("7","Log10 Energy (MeV) of charged pions",fNBinsE,-4.,6.,1.0);
129  fHisto->Add1D("8","Log10 Energy (MeV) of pi0",fNBinsE,-4.,6.,1.0);
130  fHisto->Add1D("9","Log10 Energy (MeV) of charged kaons",fNBinsE,-4.,6.,1.0);
131  fHisto->Add1D("10","Log10 Energy (MeV) of neutral kaons",fNBinsE,-4.,6.,1.0);
132  fHisto->Add1D("11","Log10 Energy (MeV) of deuterons and tritons",fNBinsE,-5.,5.,1.0);
133  fHisto->Add1D("12","Log10 Energy (MeV) of He3 and alpha",fNBinsE,-5.,5.,1.0);
134  fHisto->Add1D("13","Log10 Energy (MeV) of Generic Ions",fNBinsE,-5.,5.,1.0);
135  fHisto->Add1D("14","Log10 Energy (MeV) of muons",fNBinsE,-4.,6.,1.0);
136  fHisto->Add1D("15","log10 Energy (MeV) of side-leaked neutrons",fNBinsE,-5.,5.,1.0);
137  fHisto->Add1D("16","log10 Energy (MeV) of forward-leaked neutrons",fNBinsE,-5.,5.,1.0);
138  fHisto->Add1D("17","log10 Energy (MeV) of backward-leaked neutrons",fNBinsE,-5.,5.,1.0);
139  fHisto->Add1D("18","log10 Energy (MeV) of leaking protons",fNBinsE,-4.,6.,1.0);
140  fHisto->Add1D("19","log10 Energy (MeV) of leaking charged pions",fNBinsE,-4.,6.,1.0);
141  fHisto->Add1D("20","Log10 Energy (MeV) of pi+",fNBinsE,-4.,6.,1.0);
142  fHisto->Add1D("21","Log10 Energy (MeV) of pi-",fNBinsE,-4.,6.,1.0);
143  fHisto->Add1D("22","Energy deposition in the target normalized to beam energy",
144  110,0.0,1.1,1.0);
145  fHisto->Add1D("23","EM energy deposition in the target normalized to beam energy",
146  110,0.0,1.1,1.0);
147  fHisto->Add1D("24","Pion energy deposition in the target normalized to beam energy",
148  110,0.0,1.1,1.0);
149  fHisto->Add1D("25","Proton energy deposition in the target normalized to beam energy",
150  110,0.0,1.1,1.0);
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154 
156 {
157  fAbsZ0 = -0.5*fLength;
158  fNevt = 0;
159  fNelec = 0;
160  fNposit = 0;
161  fNgam = 0;
162  fNstep = 0;
163  fNprot_leak = 0;
164  fNpiofNleak = 0;
165  fNions = 0;
166  fNdeut = 0;
167  fNalpha = 0;
168  fNkaons = 0;
169  fNmuons = 0;
170  fNcpions = 0;
171  fNpi0 = 0;
172  fNneutron = 0;
173  fNproton = 0;
174  fNaproton = 0;
175  fNneu_forw = 0;
176  fNneu_leak = 0;
177  fNneu_back = 0;
178 
179  fEdepSum = 0.0;
180  fEdepSum2 = 0.0;
181 
182  if(!fHistoBooked) { bookHisto(); }
183  fHisto->Book();
184 
185  if(fVerbose > 0) {
186  G4cout << "HistoManager: Histograms are booked and run has been started"
187  <<G4endl;
188  }
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192 
194 {
195 
196  G4cout << "HistoManager: End of run actions are started" << G4endl;
197 
198  // Average values
199  G4cout<<"========================================================"<<G4endl;
200 
201  G4double x = (G4double)fNevt;
202  if(fNevt > 0) { x = 1.0/x; }
203 
204  G4double xe = x*(G4double)fNelec;
205  G4double xg = x*(G4double)fNgam;
206  G4double xp = x*(G4double)fNposit;
207  G4double xs = x*(G4double)fNstep;
208  G4double xn = x*(G4double)fNneutron;
209  G4double xpn = x*(G4double)fNproton;
210  G4double xap = x*(G4double)fNaproton;
211  G4double xnf = x*(G4double)fNneu_forw;
212  G4double xnb = x*(G4double)fNneu_leak;
213  G4double xnbw= x*(G4double)fNneu_back;
214  G4double xpl = x*(G4double)fNprot_leak;
215  G4double xal = x*(G4double)fNpiofNleak;
216  G4double xpc = x*(G4double)fNcpions;
217  G4double xp0 = x*(G4double)fNpi0;
218  G4double xpk = x*(G4double)fNkaons;
219  G4double xpm = x*(G4double)fNmuons;
220  G4double xid = x*(G4double)fNdeut;
221  G4double xia = x*(G4double)fNalpha;
222  G4double xio = x*(G4double)fNions;
223 
224  fEdepSum *= x;
225  fEdepSum2 *= x;
226  fEdepSum2 -= fEdepSum*fEdepSum;
227  if(fEdepSum2 > 0.0) { fEdepSum2 = std::sqrt(fEdepSum2); }
228  else { fEdepSum2 = 0.0; }
229 
230  G4cout << "Beam particle "
231  << fPrimaryDef->GetParticleName() <<G4endl;
232  G4cout << "Beam Energy(MeV) "
233  << fPrimaryKineticEnergy/MeV <<G4endl;
234  G4cout << "Number of events " << fNevt <<G4endl;
235  G4cout << std::setprecision(4) << "Average energy deposit (MeV) " << fEdepSum/MeV
236  << " RMS(MeV) " << fEdepSum2/MeV << G4endl;
237  G4cout << std::setprecision(4) << "Average number of steps " << xs << G4endl;
238  G4cout << std::setprecision(4) << "Average number of gamma " << xg << G4endl;
239  G4cout << std::setprecision(4) << "Average number of e- " << xe << G4endl;
240  G4cout << std::setprecision(4) << "Average number of e+ " << xp << G4endl;
241  G4cout << std::setprecision(4) << "Average number of neutrons " << xn << G4endl;
242  G4cout << std::setprecision(4) << "Average number of protons " << xpn << G4endl;
243  G4cout << std::setprecision(4) << "Average number of antiprotons " << xap << G4endl;
244  G4cout << std::setprecision(4) << "Average number of pi+ & pi- " << xpc << G4endl;
245  G4cout << std::setprecision(4) << "Average number of pi0 " << xp0 << G4endl;
246  G4cout << std::setprecision(4) << "Average number of kaons " << xpk << G4endl;
247  G4cout << std::setprecision(4) << "Average number of muons " << xpm << G4endl;
248  G4cout << std::setprecision(4) << "Average number of deuterons+tritons " << xid << G4endl;
249  G4cout << std::setprecision(4) << "Average number of He3+alpha " << xia << G4endl;
250  G4cout << std::setprecision(4) << "Average number of ions " << xio << G4endl;
251  G4cout << std::setprecision(4) << "Average number of forward neutrons " << xnf << G4endl;
252  G4cout << std::setprecision(4) << "Average number of reflected neutrons " << xnb << G4endl;
253  G4cout << std::setprecision(4) << "Average number of leaked neutrons " << xnbw << G4endl;
254  G4cout << std::setprecision(4) << "Average number of proton leak " << xpl << G4endl;
255  G4cout << std::setprecision(4) << "Average number of pion leak " << xal << G4endl;
256  G4cout<<"========================================================"<<G4endl;
257  G4cout<<G4endl;
258 
259  // normalise histograms
260  for(G4int i=0; i<fNHisto; i++) {
261  fHisto->ScaleH1(i,x);
262  }
263 
264  fHisto->Save();
265 }
266 
267 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
268 
270 {
271  fEdepEvt = 0.0;
272  fEdepEM = 0.0;
273  fEdepPI = 0.0;
274  fEdepP = 0.0;
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278 
280 {
281  fEdepSum += fEdepEvt;
282  fEdepSum2 += fEdepEvt*fEdepEvt;
283  fHisto->Fill(21,fEdepEvt/fPrimaryKineticEnergy,1.0);
284  fHisto->Fill(22,fEdepEM/fPrimaryKineticEnergy,1.0);
285  fHisto->Fill(23,fEdepPI/fPrimaryKineticEnergy,1.0);
286  fHisto->Fill(24,fEdepP/fPrimaryKineticEnergy,1.0);
287 }
288 
289 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
290 
291 void HistoManager::ScoreNewTrack(const G4Track* track)
292 {
293  const G4ParticleDefinition* pd = track->GetDefinition();
294  G4String name = pd->GetParticleName();
295  G4double e = track->GetKineticEnergy();
296 
297  // Primary track
298  if(0 == track->GetParentID()) {
299 
300  fNevt++;
301  fPrimaryKineticEnergy = e;
302  fPrimaryDef = pd;
304  if(1 < fVerbose)
305  G4cout << "### Primary " << name
306  << " kinE(MeV)= " << e/MeV
307  << "; m(MeV)= " << pd->GetPDGMass()/MeV
308  << "; pos(mm)= " << track->GetPosition()/mm
309  << "; dir= " << track->GetMomentumDirection()
310  << G4endl;
311 
312  // Secondary track
313  } else {
314  if(1 < fVerbose)
315  G4cout << "=== Secondary " << name
316  << " kinE(MeV)= " << e/MeV
317  << "; m(MeV)= " << pd->GetPDGMass()/MeV
318  << "; pos(mm)= " << track->GetPosition()/mm
319  << "; dir= " << track->GetMomentumDirection()
320  << G4endl;
321  e = std::log10(e/MeV);
322  if(pd == G4Gamma::Gamma()) {
323  fNgam++;
324  fHisto->Fill(1,e,1.0);
325  } else if ( pd == G4Electron::Electron()) {
326  fNelec++;
327  fHisto->Fill(2,e,1.0);
328  } else if ( pd == G4Positron::Positron()) {
329  fNposit++;
330  fHisto->Fill(3,e,1.0);
331  } else if ( pd == G4Proton::Proton()) {
332  fNproton++;
333  fHisto->Fill(4,e,1.0);
334  } else if ( pd == fNeutron) {
335  fNneutron++;
336  fHisto->Fill(5,e,1.0);
337  } else if ( pd == G4AntiProton::AntiProton()) {
338  fNaproton++;
339  } else if ( pd == G4PionPlus::PionPlus() ) {
340  fNcpions++;
341  fHisto->Fill(6,e,1.0);
342  fHisto->Fill(19,e,1.0);
343 
344  } else if ( pd == G4PionMinus::PionMinus()) {
345  fNcpions++;
346  fHisto->Fill(6,e,1.0);
347  fHisto->Fill(20,e,1.0);
348 
349  } else if ( pd == G4PionZero::PionZero()) {
350  fNpi0++;
351  fHisto->Fill(7,e,1.0);
352  } else if ( pd == G4KaonPlus::KaonPlus() || pd == G4KaonMinus::KaonMinus()) {
353  fNkaons++;
354  fHisto->Fill(8,e,1.0);
355  } else if ( pd == G4KaonZeroShort::KaonZeroShort() || pd == G4KaonZeroLong::KaonZeroLong()) {
356  fNkaons++;
357  fHisto->Fill(9,e,1.0);
358  } else if ( pd == G4Deuteron::Deuteron() || pd == G4Triton::Triton()) {
359  fNdeut++;
360  fHisto->Fill(10,e,1.0);
361  } else if ( pd == G4He3::He3() || pd == G4Alpha::Alpha()) {
362  fNalpha++;
363  fHisto->Fill(11,e,1.0);
364  } else if ( pd->GetParticleType() == "nucleus") {
365  fNions++;
366  fHisto->Fill(12,e,1.0);
367  } else if ( pd == G4MuonPlus::MuonPlus() || pd == G4MuonMinus::MuonMinus()) {
368  fNmuons++;
369  fHisto->Fill(13,e,1.0);
370  }
371  }
372 }
373 
374 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
375 
377 {
378  fNstep++;
379  G4double fEdep = step->GetTotalEnergyDeposit();
380  if(fEdep >= DBL_MIN) {
381  const G4Track* track = step->GetTrack();
382 
383  G4ThreeVector pos =
384  (step->GetPreStepPoint()->GetPosition() +
385  step->GetPostStepPoint()->GetPosition())*0.5;
386 
387  G4double z = pos.z() - fAbsZ0;
388 
389  // scoring
390  fEdepEvt += fEdep;
391  fHisto->Fill(0,z,fEdep);
392  const G4ParticleDefinition* pd = track->GetDefinition();
393 
394  if(pd == G4Gamma::Gamma() || pd == G4Electron::Electron()
395  || pd == G4Positron::Positron()) {
396  fEdepEM += fEdep;
397  } else if ( pd == G4PionPlus::PionPlus() || pd == G4PionMinus::PionMinus()) {
398  fEdepPI += fEdep;
399  } else if ( pd == G4Proton::Proton() || pd == G4AntiProton::AntiProton()) {
400  fEdepP += fEdep;
401  }
402 
403  if(1 < fVerbose) {
404  G4cout << "HistoManager::AddEnergy: e(keV)= " << fEdep/keV
405  << "; z(mm)= " << z/mm
406  << "; step(mm)= " << step->GetStepLength()/mm
407  << " by " << pd->GetParticleName()
408  << " E(MeV)= " << track->GetKineticEnergy()/MeV
409  << G4endl;
410  }
411  }
412 }
413 
414 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
415 
417 {
418  const G4ParticleDefinition* pd = track->GetDefinition();
419  G4double e = std::log10(track->GetKineticEnergy()/MeV);
420 
421  G4ThreeVector pos = track->GetPosition();
422  G4ThreeVector dir = track->GetMomentumDirection();
423  G4double x = pos.x();
424  G4double y = pos.y();
425  G4double z = pos.z();
426 
427  G4bool isLeaking = false;
428 
429  // Forward
430  if(z > -fAbsZ0 && dir.z() > 0.0) {
431  isLeaking = true;
432  if(pd == fNeutron) {
433  ++fNneu_forw;
434  fHisto->Fill(15,e,1.0);
435  } else isLeaking = true;
436 
437  // Backward
438  } else if (z < fAbsZ0 && dir.z() < 0.0) {
439  isLeaking = true;
440  if(pd == fNeutron) {
441  ++fNneu_back;
442  fHisto->Fill(16,e,1.0);
443  } else isLeaking = true;
444 
445  // Side
446  } else if (std::abs(z) <= -fAbsZ0 && x*dir.x() + y*dir.y() > 0.0) {
447  isLeaking = true;
448  if(pd == fNeutron) {
449  ++fNneu_leak;
450  fHisto->Fill(14,e,1.0);
451  } else isLeaking = true;
452  }
453 
454  // protons and pions
455  if(isLeaking) {
456  if(pd == G4Proton::Proton()) {
457  fHisto->Fill(17,e,1.0);
458  ++fNprot_leak;
459  } else if (pd == G4PionPlus::PionPlus() || pd == G4PionMinus::PionMinus()) {
460  fHisto->Fill(18,e,1.0);
461  ++fNpiofNleak;
462  }
463  }
464 }
465 
466 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
467 
469 {
470  fVerbose = val;
471  fHisto->SetVerbose(val);
472 }
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
475 
477 {
478  fHisto->Fill(id, x, w);
479 }
480 
481 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
482