Geant4  10.01
Run.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: Run.cc 71376 2013-06-14 07:44:50Z maire $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "Run.hh"
35 #include "DetectorConstruction.hh"
36 #include "PrimaryGeneratorAction.hh"
37 #include "HistoManager.hh"
38 
39 #include "G4EmCalculator.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4UnitsTable.hh"
42 
43 #include <iomanip>
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 : G4Run(),
49  fDetector(det),
50  fParticle(0), fEkin(0.)
51 {
58  fMscThetaCentral = 0.;
59 
60  fNbGamma = fNbElect = fNbPosit = 0;
61 
62  fTransmit[0] = fTransmit[1] = fReflect[0] = fReflect[1] = 0;
63 
64  fMscEntryCentral = 0;
65 
66  fEnergyLeak[0] = fEnergyLeak[1] = fEnergyLeak2[0] = fEnergyLeak2[1] = 0.;
67  }
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
71 Run::~Run()
72 { }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75 
77 {
78  fParticle = particle;
79  fEkin = energy;
80 
81  //compute theta0
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
87 void Run::Merge(const G4Run* run)
88 {
89  const Run* localRun = static_cast<const Run*>(run);
90 
91  // pass information about primary particle
92  fParticle = localRun->fParticle;
93  fEkin = localRun->fEkin;
94 
96 
97  // accumulate sums
98  //
99  fEnergyDeposit += localRun->fEnergyDeposit;
100  fEnergyDeposit2 += localRun->fEnergyDeposit2;
101  fTrakLenCharged += localRun->fTrakLenCharged;
102  fTrakLenCharged2 += localRun->fTrakLenCharged2;
103  fTrakLenNeutral += localRun->fTrakLenNeutral;
104  fTrakLenNeutral2 += localRun->fTrakLenNeutral2;
105  fNbStepsCharged += localRun->fNbStepsCharged;
106  fNbStepsCharged2 += localRun->fNbStepsCharged2;
107  fNbStepsNeutral += localRun->fNbStepsNeutral;
108  fNbStepsNeutral2 += localRun->fNbStepsNeutral2;
109 
110  fNbGamma += localRun->fNbGamma;
111  fNbElect += localRun->fNbElect;
112  fNbPosit += localRun->fNbPosit;
113 
114  fTransmit[0] += localRun->fTransmit[0];
115  fTransmit[1] += localRun->fTransmit[1];
116  fReflect[0] += localRun->fReflect[0];
117  fReflect[1] += localRun->fReflect[1];
118 
119  fMscEntryCentral += localRun->fMscEntryCentral;
120 
121  fEnergyLeak[0] += localRun->fEnergyLeak[0];
122  fEnergyLeak[1] += localRun->fEnergyLeak[1];
123  fEnergyLeak2[0] += localRun->fEnergyLeak2[0];
124  fEnergyLeak2[1] += localRun->fEnergyLeak2[1];
125 
126  G4Run::Merge(run);
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
131 void Run::EndOfRun()
132 {
133  // compute mean and rms
134  //
135  G4int TotNbofEvents = numberOfEvent;
136  if (TotNbofEvents == 0) return;
137 
138  G4double EnergyBalance = fEnergyDeposit + fEnergyLeak[0] + fEnergyLeak[1];
139  EnergyBalance /= TotNbofEvents;
140 
141  fEnergyDeposit /= TotNbofEvents; fEnergyDeposit2 /= TotNbofEvents;
143  if (rmsEdep>0.) rmsEdep = std::sqrt(rmsEdep/TotNbofEvents);
144  else rmsEdep = 0.;
145 
146  fTrakLenCharged /= TotNbofEvents; fTrakLenCharged2 /= TotNbofEvents;
148  if (rmsTLCh>0.) rmsTLCh = std::sqrt(rmsTLCh/TotNbofEvents);
149  else rmsTLCh = 0.;
150 
151  fTrakLenNeutral /= TotNbofEvents; fTrakLenNeutral2 /= TotNbofEvents;
153  if (rmsTLNe>0.) rmsTLNe = std::sqrt(rmsTLNe/TotNbofEvents);
154  else rmsTLNe = 0.;
155 
156  fNbStepsCharged /= TotNbofEvents; fNbStepsCharged2 /= TotNbofEvents;
158  if (rmsStCh>0.) rmsStCh = std::sqrt(rmsTLCh/TotNbofEvents);
159  else rmsStCh = 0.;
160 
161  fNbStepsNeutral /= TotNbofEvents; fNbStepsNeutral2 /= TotNbofEvents;
163  if (rmsStNe>0.) rmsStNe = std::sqrt(rmsTLCh/TotNbofEvents);
164  else rmsStNe = 0.;
165 
166  G4double Gamma = (G4double)fNbGamma/TotNbofEvents;
167  G4double Elect = (G4double)fNbElect/TotNbofEvents;
168  G4double Posit = (G4double)fNbPosit/TotNbofEvents;
169 
170  G4double transmit[2];
171  transmit[0] = 100.*fTransmit[0]/TotNbofEvents;
172  transmit[1] = 100.*fTransmit[1]/TotNbofEvents;
173 
174  G4double reflect[2];
175  reflect[0] = 100.*fReflect[0]/TotNbofEvents;
176  reflect[1] = 100.*fReflect[1]/TotNbofEvents;
177 
178  G4double rmsMsc = 0., tailMsc = 0.;
179  if (fMscEntryCentral > 0) {
182  if (rmsMsc > 0.) { rmsMsc = std::sqrt(rmsMsc); }
183  if(fTransmit[1] > 0.0) {
184  tailMsc = 100.- (100.*fMscEntryCentral)/(2*fTransmit[1]);
185  }
186  }
187 
188  fEnergyLeak[0] /= TotNbofEvents; fEnergyLeak2[0] /= TotNbofEvents;
189  G4double rmsEl0 = fEnergyLeak2[0] - fEnergyLeak[0]*fEnergyLeak[0];
190  if (rmsEl0>0.) rmsEl0 = std::sqrt(rmsEl0/TotNbofEvents);
191  else rmsEl0 = 0.;
192 
193  fEnergyLeak[1] /= TotNbofEvents; fEnergyLeak2[1] /= TotNbofEvents;
194  G4double rmsEl1 = fEnergyLeak2[1] - fEnergyLeak[1]*fEnergyLeak[1];
195  if (rmsEl1>0.) rmsEl1 = std::sqrt(rmsEl1/TotNbofEvents);
196  else rmsEl1 = 0.;
197 
198 
199  //Stopping Power from input Table.
200  //
203  G4double density = material->GetDensity();
204  G4String partName = fParticle->GetParticleName();
205 
206  G4EmCalculator emCalculator;
207  G4double dEdxTable = 0., dEdxFull = 0.;
208  if (fParticle->GetPDGCharge()!= 0.) {
209  dEdxTable = emCalculator.GetDEDX(fEkin,fParticle,material);
210  dEdxFull = emCalculator.ComputeTotalDEDX(fEkin,fParticle,material);
211  }
212  G4double stopTable = dEdxTable/density;
213  G4double stopFull = dEdxFull /density;
214 
215  //Stopping Power from simulation.
216  //
217  G4double meandEdx = fEnergyDeposit/length;
218  G4double stopPower = meandEdx/density;
219 
220  G4cout << "\n ======================== run summary ======================\n";
221 
222  G4int prec = G4cout.precision(3);
223 
224  G4cout << "\n The run was " << TotNbofEvents << " " << partName << " of "
225  << G4BestUnit(fEkin,"Energy") << " through "
226  << G4BestUnit(length,"Length") << " of "
227  << material->GetName() << " (density: "
228  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
229 
230  G4cout.precision(4);
231 
232  G4cout << "\n Total energy deposit in absorber per event = "
233  << G4BestUnit(fEnergyDeposit,"Energy") << " +- "
234  << G4BestUnit(rmsEdep, "Energy")
235  << G4endl;
236 
237  G4cout << "\n -----> Mean dE/dx = " << meandEdx/(MeV/cm) << " MeV/cm"
238  << "\t(" << stopPower/(MeV*cm2/g) << " MeV*cm2/g)"
239  << G4endl;
240 
241  G4cout << "\n From formulas :" << G4endl;
242  G4cout << " restricted dEdx = " << dEdxTable/(MeV/cm) << " MeV/cm"
243  << "\t(" << stopTable/(MeV*cm2/g) << " MeV*cm2/g)"
244  << G4endl;
245 
246  G4cout << " full dEdx = " << dEdxFull/(MeV/cm) << " MeV/cm"
247  << "\t(" << stopFull/(MeV*cm2/g) << " MeV*cm2/g)"
248  << G4endl;
249 
250  G4cout << "\n Leakage : primary = "
251  << G4BestUnit(fEnergyLeak[0],"Energy") << " +- "
252  << G4BestUnit(rmsEl0, "Energy")
253  << " secondaries = "
254  << G4BestUnit(fEnergyLeak[1],"Energy") << " +- "
255  << G4BestUnit(rmsEl1, "Energy")
256  << G4endl;
257 
258  G4cout << " Energy balance : edep + eleak = "
259  << G4BestUnit(EnergyBalance,"Energy")
260  << G4endl;
261 
262  G4cout << "\n Total track length (charged) in absorber per event = "
263  << G4BestUnit(fTrakLenCharged,"Length") << " +- "
264  << G4BestUnit(rmsTLCh, "Length") << G4endl;
265 
266  G4cout << " Total track length (neutral) in absorber per event = "
267  << G4BestUnit(fTrakLenNeutral,"Length") << " +- "
268  << G4BestUnit(rmsTLNe, "Length") << G4endl;
269 
270  G4cout << "\n Number of steps (charged) in absorber per event = "
271  << fNbStepsCharged << " +- " << rmsStCh << G4endl;
272 
273  G4cout << " Number of steps (neutral) in absorber per event = "
274  << fNbStepsNeutral << " +- " << rmsStNe << G4endl;
275 
276  G4cout << "\n Number of secondaries per event : Gammas = " << Gamma
277  << "; electrons = " << Elect
278  << "; positrons = " << Posit << G4endl;
279 
280  G4cout << "\n Number of events with the primary particle transmitted = "
281  << transmit[1] << " %" << G4endl;
282 
283  G4cout << " Number of events with at least 1 particle transmitted "
284  << "(same charge as primary) = " << transmit[0] << " %" << G4endl;
285 
286  G4cout << "\n Number of events with the primary particle reflected = "
287  << reflect[1] << " %" << G4endl;
288 
289  G4cout << " Number of events with at least 1 particle reflected "
290  << "(same charge as primary) = " << reflect[0] << " %" << G4endl;
291 
292  // compute width of the Gaussian central part of the MultipleScattering
293  //
294  G4cout << "\n MultipleScattering:"
295  << "\n rms proj angle of transmit primary particle = "
296  << rmsMsc/mrad << " mrad (central part only)" << G4endl;
297 
298  G4cout << " computed theta0 (Highland formula) = "
299  << ComputeMscHighland()/mrad << " mrad" << G4endl;
300 
301  G4cout << " central part defined as +- "
302  << fMscThetaCentral/mrad << " mrad; "
303  << " Tail ratio = " << tailMsc << " %" << G4endl;
304 
305  // normalize histograms
306  //
307  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
308 
309  G4int ih = 1;
310  G4double binWidth = analysisManager->GetH1Width(ih);
311  G4double unit = analysisManager->GetH1Unit(ih);
312  G4double fac = unit/(TotNbofEvents*binWidth);
313  analysisManager->ScaleH1(ih,fac);
314 
315  ih = 10;
316  binWidth = analysisManager->GetH1Width(ih);
317  unit = analysisManager->GetH1Unit(ih);
318  fac = unit/(TotNbofEvents*binWidth);
319  analysisManager->ScaleH1(ih,fac);
320 
321  ih = 12;
322  analysisManager->ScaleH1(ih,1./TotNbofEvents);
323 
324  // reset default precision
325  G4cout.precision(prec);
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
329 
331 {
332  //compute the width of the Gaussian central part of the MultipleScattering
333  //projected angular distribution.
334  //Eur. Phys. Jour. C15 (2000) page 166, formule 23.9
335 
338  if (t < DBL_MIN) return 0.;
339 
340  G4double T = fEkin;
342  G4double z = std::abs(fParticle->GetPDGCharge()/eplus);
343 
344  G4double bpc = T*(T+2*M)/(T+M);
345  G4double teta0 = 13.6*MeV*z*std::sqrt(t)*(1.+0.038*std::log(t))/bpc;
346  return teta0;
347 }
348 
349 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static const double cm
Definition: G4SIunits.hh:106
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int numberOfEvent
Definition: G4Run.hh:59
static const double MeV
Definition: G4SIunits.hh:193
G4double fEkin
Definition: Run.hh:65
static const double cm2
Definition: G4SIunits.hh:107
G4double fMscThetaCentral
Definition: Run.hh:113
G4double fTrakLenCharged
Definition: Run.hh:108
G4double fNbStepsNeutral
Definition: Run.hh:111
G4double fTrakLenNeutral2
Definition: Run.hh:109
Run()
Definition: Run.cc:43
G4int fTransmit[2]
Definition: Run.hh:116
G4double z
Definition: TRTMaterials.hh:39
const G4String & GetName() const
Definition: G4Material.hh:176
G4double ComputeMscHighland()
Definition: Run.cc:330
G4double GetDensity() const
Definition: G4Material.hh:178
G4double fNbStepsCharged2
Definition: Run.hh:110
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double fTrakLenCharged2
Definition: Run.hh:108
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double fEnergyDeposit2
Definition: Run.hh:107
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4int fReflect[2]
Definition: Run.hh:116
G4double density
Definition: TRTMaterials.hh:39
static const double prec
Definition: RanecuEngine.cc:58
G4double fNbStepsCharged
Definition: Run.hh:110
G4GLOB_DLL std::ostream G4cout
G4Material * GetAbsorberMaterial()
G4double fTrakLenNeutral
Definition: Run.hh:109
G4int fNbElect
Definition: Run.hh:115
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:61
G4double fMscProjecTheta2
Definition: Run.hh:112
G4int fNbPosit
Definition: Run.hh:115
static const double mrad
Definition: G4SIunits.hh:131
G4double GetRadlen() const
Definition: G4Material.hh:218
G4double fEnergyLeak[2]
Definition: Run.hh:119
G4int fMscEntryCentral
Definition: Run.hh:117
G4double fEnergyLeak2[2]
Definition: Run.hh:119
G4double GetPDGMass() const
G4int fNbGamma
Definition: Run.hh:115
G4double energy(const ThreeVector &p, const G4double m)
static const double g
Definition: G4SIunits.hh:162
#define DBL_MIN
Definition: templates.hh:75
DetectorConstruction * fDetector
Definition: Run.hh:63
static const G4double fac
Detector construction class to demonstrate various ways of placement.
void EndOfRun()
Definition: Run.cc:147
#define G4endl
Definition: G4ios.hh:61
G4ParticleDefinition * fParticle
Definition: Run.hh:64
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
Definition: Run.cc:77
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
static const double eplus
Definition: G4SIunits.hh:178
std::vector< G4double > fEnergyDeposit[MaxAbsor]
Definition: Run.hh:82
G4double GetPDGCharge() const
G4double fNbStepsNeutral2
Definition: Run.hh:111
virtual void Merge(const G4Run *)
Definition: Run.cc:115
~Run()
Definition: Run.cc:72
G4double fMscProjecTheta
Definition: Run.hh:112