Geant4  10.01.p03
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 
37 #include "EventAction.hh"
38 #include "HistoManager.hh"
39 #include "PrimaryGeneratorAction.hh"
40 
41 #include "G4Material.hh"
42 #include "G4Event.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4UnitsTable.hh"
45 #include <iomanip>
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
50 : G4Run(),
51  fDetector(detector),
52  fParticle(0), fEkin(0.),
53  fEdeposit(0.), fEdeposit2(0.),
54  fTrackLen(0.), fTrackLen2(0.),
55  fProjRange(0.), fProjRange2(0.),
56  fNbOfSteps(0), fNbOfSteps2(0),
57  fStepSize(0.), fStepSize2(0.)
58 {
59  for (G4int i=0; i<3; ++i) fStatus[i] = 0;
60  for (G4int i=0; i<MaxAbsor; ++i) {
61  fCsdaRange[i] = 0.; fXfrontNorm[i] = 0.;
62  }
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
67 Run::~Run()
68 { }
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
73 {
74  fParticle = particle;
75  fEkin = energy;
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
80 void Run::AddEdep (G4double e)
81 {
82  fEdeposit += e;
83  fEdeposit2 += e*e;
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 {
90  fTrackLen += t;
91  fTrackLen2 += t*t;
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95 
97 {
98  fProjRange += x;
99  fProjRange2 += x*x;
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103 
105 {
106  fNbOfSteps += nb;
107  fNbOfSteps2 += nb*nb;
108  fStepSize += st ;
109  fStepSize2 += st*st;
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
115 {
116  fStatus[i]++ ;
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
122 {
123  fCsdaRange[i] = value;
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
129 {
130  fXfrontNorm[i] = value;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136 {
137  return fCsdaRange[i];
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
143 {
144  return fXfrontNorm[i];
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
149 void Run::Merge(const G4Run* run)
150 {
151  const Run* localRun = static_cast<const Run*>(run);
152 
153  // pass information about primary particle
154  fParticle = localRun->fParticle;
155  fEkin = localRun->fEkin;
156 
157  // accumulate sums
158  fEdeposit += localRun->fEdeposit;
159  fEdeposit2 += localRun->fEdeposit2;
160  fTrackLen += localRun->fTrackLen;
161  fTrackLen2 += localRun->fTrackLen2;
162  fProjRange += localRun->fProjRange;
163  fProjRange2 += localRun->fProjRange2;
164  fNbOfSteps += localRun->fNbOfSteps ;
165  fNbOfSteps2 += localRun->fNbOfSteps2;
166  fStepSize += localRun->fStepSize;
167  fStepSize2 += localRun->fStepSize2;
168 
169  G4int nbOfAbsor = fDetector->GetNbOfAbsor();
170  for (G4int i=1; i<=nbOfAbsor; ++i) {
171  fCsdaRange[i] = localRun->fCsdaRange[i];
172  fXfrontNorm[i] = localRun->fXfrontNorm[i];
173  }
174 
175  for (G4int i=0; i<3; ++i) fStatus[i] += localRun->fStatus[i];
176 
177  G4Run::Merge(run);
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181 
182 void Run::EndOfRun()
183 {
184  std::ios::fmtflags mode = G4cout.flags();
185  G4cout.setf(std::ios::fixed,std::ios::floatfield);
186  G4int prec = G4cout.precision(2);
187 
188  //run conditions
189  //
190  G4String partName = fParticle->GetParticleName();
191  G4int nbOfAbsor = fDetector->GetNbOfAbsor();
192 
193  G4cout << "\n ======================== run summary =====================\n";
194  G4cout
195  << "\n The run is " << numberOfEvent << " "<< partName << " of "
196  << G4BestUnit(fEkin,"Energy")
197  << " through " << nbOfAbsor << " absorbers: \n";
198  for (G4int i=1; i<= nbOfAbsor; i++) {
199  G4Material* material = fDetector->GetAbsorMaterial(i);
200  G4double thickness = fDetector->GetAbsorThickness(i);
201  G4double density = material->GetDensity();
202  G4cout << std::setw(20) << G4BestUnit(thickness,"Length") << " of "
203  << material->GetName() << " (density: "
204  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
205  }
206 
207  if (numberOfEvent == 0) {
208  G4cout.setf(mode,std::ios::floatfield);
209  G4cout.precision(prec);
210  return;
211  }
212 
215  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
216 
217  G4cout.precision(3);
218  G4cout
219  << "\n Total Energy deposited = " << G4BestUnit(fEdeposit,"Energy")
220  << " +- " << G4BestUnit( rms,"Energy")
221  << G4endl;
222 
223  //compute track length of primary track
224  //
226  rms = fTrackLen2 - fTrackLen*fTrackLen;
227  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
228 
229  G4cout.precision(3);
230  G4cout
231  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
232  << " +- " << G4BestUnit( rms,"Length");
233 
234  //compare with csda range
235  //
236  G4int NbOfAbsor = fDetector->GetNbOfAbsor();
237  if (NbOfAbsor == 1) {
238  G4cout
239  << "\n Range from EmCalculator = " << G4BestUnit(fCsdaRange[1],"Length")
240  << " (from full dE/dx)" << G4endl;
241  }
242 
243  //compute projected range of primary track
244  //
247  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
248 
249  G4cout
250  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
251  << " +- " << G4BestUnit( rms,"Length")
252  << G4endl;
253 
254  //nb of steps and step size of primary track
255  //
256  G4double dNofEvents = double(numberOfEvent);
257  G4double fNbSteps = fNbOfSteps/dNofEvents,
258  fNbSteps2 = fNbOfSteps2/dNofEvents;
259  rms = fNbSteps2 - fNbSteps*fNbSteps;
260  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
261 
262  G4cout.precision(2);
263  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms;
264 
266  rms = fStepSize2 - fStepSize*fStepSize;
267  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
268 
269  G4cout.precision(3);
270  G4cout
271  << "\t Step size= " << G4BestUnit(fStepSize,"Length")
272  << " +- " << G4BestUnit( rms,"Length")
273  << G4endl;
274 
275  //transmission coefficients
276  //
277  G4double absorbed = 100.*fStatus[0]/dNofEvents;
278  G4double transmit = 100.*fStatus[1]/dNofEvents;
279  G4double reflected = 100.*fStatus[2]/dNofEvents;
280 
281  G4cout.precision(2);
282  G4cout
283  << "\n absorbed = " << absorbed << " %"
284  << " transmit = " << transmit << " %"
285  << " reflected = " << reflected << " %" << G4endl;
286 
287  // normalize histograms of longitudinal energy profile
288  //
289  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
290  G4int ih = 1;
291  G4double binWidth = analysisManager->GetH1Width(ih);
292  G4double fac = (1./(numberOfEvent*binWidth))*(mm/MeV);
293  analysisManager->ScaleH1(ih,fac);
294 
295  ih = 8;
296  binWidth = analysisManager->GetH1Width(ih);
297  fac = (1./(numberOfEvent*binWidth))*(g/(MeV*cm2));
298  analysisManager->ScaleH1(ih,fac);
299 
300  // reset default formats
301  G4cout.setf(mode,std::ios::floatfield);
302  G4cout.precision(prec);
303 }
304 
305 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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
Run()
Definition: Run.cc:43
G4double fTrackLen
Definition: Run.hh:76
G4int fNbOfSteps
Definition: Run.hh:78
const G4String & GetName() const
Definition: G4Material.hh:178
void SetCsdaRange(G4int i, G4double value)
Definition: Run.cc:121
G4double GetDensity() const
Definition: G4Material.hh:180
const G4int MaxAbsor
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double fProjRange2
Definition: Run.hh:78
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4Material * GetAbsorMaterial(G4int i)
void AddEdep(G4double val)
Definition: Run.hh:61
G4double fStepSize2
Definition: Run.hh:79
G4double density
Definition: TRTMaterials.hh:39
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
G4double fEdeposit
Definition: Run.hh:75
void SetXfrontNorm(G4int i, G4double value)
Definition: Run.cc:128
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:61
G4double fProjRange
Definition: Run.hh:78
G4double fEdeposit2
Definition: Run.hh:75
G4double GetAbsorThickness(G4int i)
G4int fNbOfSteps2
Definition: Run.hh:78
G4double GetXfrontNorm(G4int i)
Definition: Run.cc:142
G4double energy(const ThreeVector &p, const G4double m)
G4double GetCsdaRange()
Definition: Run.cc:114
static const double g
Definition: G4SIunits.hh:162
void AddTrackStatus(G4int i)
Definition: Run.cc:114
G4double fStepSize
Definition: Run.hh:79
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
G4double fCsdaRange[MaxAbsor]
Definition: Run.hh:82
void AddTrackLength(G4double t)
Definition: Run.cc:88
void AddStepSize(G4int nb, G4double st)
Definition: Run.cc:104
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
Definition: Run.cc:77
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
G4double fXfrontNorm[MaxAbsor]
Definition: Run.hh:83
G4int fStatus[3]
Definition: Run.hh:80
virtual void Merge(const G4Run *)
Definition: Run.cc:115
static const double mm
Definition: G4SIunits.hh:102
G4double fTrackLen2
Definition: Run.hh:76
~Run()
Definition: Run.cc:72
void AddProjRange(G4double x)
Definition: Run.hh:63