Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
extended/electromagnetic/TestEm8/src/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 //
26 /// \file electromagnetic/TestEm8/src/HistoManager.cc
27 /// \brief Implementation of the HistoManager class
28 //
29 // $Id: HistoManager.cc 79241 2014-02-20 16:03:35Z gcosmo $
30 //
31 //---------------------------------------------------------------------------
32 //
33 // ClassName: HistoManager
34 //
35 // Author: V.Ivanchenko 01.09.2010
36 //
37 //----------------------------------------------------------------------------
38 //
39 
40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42 
43 #include "HistoManager.hh"
44 #include "G4UnitsTable.hh"
45 #include "Histo.hh"
46 #include "G4Step.hh"
47 #include "G4LossTableManager.hh"
48 #include "G4ElectronIonPair.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4PhysicalConstants.hh"
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53 
54 HistoManager* HistoManager::fManager = 0;
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57 
59 {
60  if(!fManager) {
61  fManager = new HistoManager();
62  }
63  return fManager;
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
69  : fHisto(0),
70  fElIonPair(0)
71 {
72  fNHisto = 2;
73  fVerbose = 1;
74  fMaxEnergy = 100.*keV;
75  fBinsE = 100;
76  fBinsCluster = 200;
77 
78  // normalisation to PAI
79  fFactorALICE = 325;
80 
81  // normalisation to Opt0
82  //fFactorALICE = 275;
83 
84  fEvt = 0.0;
85  fTotStepGas = 0.0;
86  fTotCluster = 0.0;
87  fMeanCluster= 0.0;
88  fOverflow = 0.0;
89 
90  fTotEdep = 0.0;
91  fStepGas = 0;
92  fCluster = 0;
93 
94  fHistoBooked = false;
95 
96  fHisto = new Histo();
98 }
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101 
103 {
104  delete fHisto;
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108 
110 {
111  // initilise scoring
112  fEvt = 0.0;
113  fTotStepGas = 0.0;
114  fTotCluster = 0.0;
115  fMeanCluster= 0.0;
116  fOverflow = 0.0;
117 
118  if(fHisto->IsActive() && !fHistoBooked) {
119 
120  fHisto->Add1D("10","Energy deposition in detector (keV)",
121  fBinsE,0.0,fMaxEnergy/keV,1.0);
122  fHisto->Add1D("11","Number of primary clusters",
123  fBinsCluster,-0.5,fBinsCluster-0.5,1.0);
124  fHisto->Add1D("12","Energy deposition in detector (ADC)",
125  200,0.0,2000,1.0);
126 
127  fHisto->Activate(0, true);
128  fHisto->Activate(1, true);
129  fHisto->Activate(2, true);
130  }
131  fHisto->Book();
132 
133  fEgas.resize(fBinsE,0.0);
134  fEdep.reset();
135 
136  if(fVerbose > 0) {
137  G4cout << "HistoManager: Histograms are booked and run has been started"
138  << G4endl;
139  G4cout << " BinsCluster= " << fBinsCluster << " BinsE= " << fBinsE
140  << " Emax(keV)= " << fMaxEnergy/keV << G4endl;
141  }
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145 
147 {
148  G4double norm = fEvt;
149  if(fEvt > 0.0) { norm = 1.0/norm; }
150 
151  fTotStepGas *= norm;
152  fTotCluster *= norm;
153  fMeanCluster *= norm;
154  fOverflow *= norm;
155 
156  G4double y1 = fEdep.mean();
157  G4double y2 = sqrt(fEdep.rms());
158 
159  G4double de = fMaxEnergy/G4double(fBinsE);
160  G4double x1 = -de*0.5;
161 
162  G4cout << " ================== run summary =====================" << G4endl;
163  G4int prec = G4cout.precision(5);
164  G4cout << " End of Run TotNbofEvents = "
165  << G4int(fEvt) << G4endl;
166  G4cout << " Energy(keV) per ADC channel = "
167  << 1.0/(keV*fFactorALICE) << G4endl;
168  /*
169  G4double p1 = 1*GeV;
170  G4double p2 = 3*GeV;
171  G4double mass = proton_mass_c2;
172  G4cout << sqrt(p1*p1 + mass*mass) - mass << " "
173  << sqrt(p2*p2 + mass*mass) - mass << G4endl;
174  */
175  G4cout << G4endl;
176  G4cout << " Mean energy deposit in absorber = " <<
177  y1/keV << " +- " << y2*std::sqrt(norm)/keV << " keV; ";
178  if(y1 > 0.0) { G4cout << " RMS/Emean = " << y2/y1; }
179  G4cout << G4endl;
180  G4cout << " Mean number of steps in absorber= "
181  << fTotStepGas << "; mean number of ion-clusters = "
182  << fTotCluster << " MeanCluster= " << fMeanCluster
183  << G4endl;
184  G4cout << G4endl;
185 
186  G4cout << " ====== Energy deposit distribution Noverflows= " << fOverflow
187  << " ====== " << G4endl ;
188  G4cout << " bin nb Elow entries normalized " << G4endl;
189 
190  std::ofstream fileOut("distribution.out", std::ios::out );
191  fileOut.setf( std::ios::scientific, std::ios::floatfield );
192 
193  x1 = 0.0;
194 
195  fileOut << fBinsE << G4endl;
196 
197  for(G4int j=0; j<fBinsE; ++j)
198  {
199  G4cout << std::setw(5) << j << std::setw(10) << x1/keV
200  << std::setw(12) << fEgas[j] << std::setw(12) << fEgas[j]*norm
201  << G4endl ;
202  fileOut << x1/keV << "\t" << fEgas[j] << G4endl;
203  x1 += de;
204  }
205  G4cout.precision(prec);
206 
207  // normalise histograms
208  if(fHisto->IsActive() && !fHistoBooked) {
209  fHisto->ScaleH1(0,norm);
210  fHisto->ScaleH1(1,norm);
211  // fHisto->ScaleH1(2,0.05);
212  fHisto->ScaleH1(2,0.1);
213  fHisto->Save();
214  }
215  G4cout << " ================== run end ==========================" << G4endl;
216 }
217 
218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
219 
221 {
222  fEvt += 1.0;
223  fTotEdep = 0.0;
224  fStepGas = 0;
225  fCluster = 0;
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229 
231 {
232  fTotStepGas += fStepGas;
233  fTotCluster += fCluster;
234 
235  G4int idx = G4int(fTotEdep*fBinsE/fMaxEnergy);
236  if(idx < 0) { idx = 0; }
237  if(idx >= fBinsE) { fOverflow += 1.0; }
238  else { fEgas[idx] += 1.0; }
239 
240  // fill histo
241  fHisto->Fill(0,fTotEdep/keV,1.0);
242  fHisto->Fill(1,fCluster,1.0);
243  fHisto->Fill(2,fTotEdep*fFactorALICE/keV,1.0);
244 
245  fEdep.fill(fTotEdep, 1.0);
246 }
247 
248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
249 
250 void HistoManager::AddEnergy(G4double edep, const G4Step* step)
251 {
252  if(1 < fVerbose) {
253  G4cout << "HistoManager::AddEnergy: e(keV)= " << edep/keV
254  << G4endl;
255  }
256  fTotEdep += edep;
257  if(step) {
258  if(1 == step->GetTrack()->GetTrackID()) { fStepGas += 1.0; }
259 
260  fMeanCluster += fElIonPair->MeanNumberOfIonsAlongStep(step);
261  fCluster += fElIonPair->SampleNumberOfIonsAlongStep(step);
262  }
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266 
G4double MeanNumberOfIonsAlongStep(const G4ParticleDefinition *, const G4Material *, G4double edepTotal, G4double edepNIEL=0.0)
static G4LossTableManager * Instance()
void fill(G4double x, G4double weight=1.)
Definition: G4StatDouble.cc:57
void AddEnergy(G4double edep, const G4Step *)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4int SampleNumberOfIonsAlongStep(const G4Step *)
G4ElectronIonPair * ElectronIonPair()
Definition: G4Step.hh:76
G4int GetTrackID() const
G4double mean() const
Definition: G4StatDouble.cc:77
#define G4endl
Definition: G4ios.hh:61
G4double rms()
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const