Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 "HistoManager.hh"
38 #include "PrimaryGeneratorAction.hh"
39 
40 #include "G4Material.hh"
41 #include "G4SystemOfUnits.hh"
42 #include "G4UnitsTable.hh"
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 
47 : G4Run(),
48  fDetector(detector),
49  fParticle(0), fEkin(0.),
50  fEdeposit(0.), fEdeposit2(0.),
51  fTrackLen(0.), fTrackLen2(0.),
52  fProjRange(0.), fProjRange2(0.),
53  fNbOfSteps(0), fNbOfSteps2(0),
54  fStepSize(0.), fStepSize2(0.),
55  fCsdaRange(0.)
56 { }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
60 Run::~Run()
61 { }
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
66 {
67  fParticle = particle;
68  fEkin = energy;
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
73 void Run::AddEdep (G4double e)
74 {
75  fEdeposit += e;
76  fEdeposit2 += e*e;
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82 {
83  fTrackLen += t;
84  fTrackLen2 += t*t;
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90 {
91  fProjRange += x;
92  fProjRange2 += x*x;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
97 void Run::AddStepSize (G4int nb, G4double st)
98 {
99  fNbOfSteps += nb;
100  fNbOfSteps2 += nb*nb;
101  fStepSize += st ;
102  fStepSize2 += st*st;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108 {
109  fCsdaRange = value;
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
115 {
116  return fCsdaRange;
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
121 void Run::Merge(const G4Run* run)
122 {
123  const Run* localRun = static_cast<const Run*>(run);
124 
125  // pass information about primary particle
126  fParticle = localRun->fParticle;
127  fEkin = localRun->fEkin;
128 
129  // accumulate sums
130  fEdeposit += localRun->fEdeposit;
131  fEdeposit2 += localRun->fEdeposit2;
132  fTrackLen += localRun->fTrackLen;
133  fTrackLen2 += localRun->fTrackLen2;
134  fProjRange += localRun->fProjRange;
135  fProjRange2 += localRun->fProjRange2;
136  fNbOfSteps += localRun->fNbOfSteps ;
137  fNbOfSteps2 += localRun->fNbOfSteps2;
138  fStepSize += localRun->fStepSize;
139  fStepSize2 += localRun->fStepSize2;
140 
141  fCsdaRange = localRun->fCsdaRange;
142 
143  G4Run::Merge(run);
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
148 void Run::EndOfRun()
149 {
150  std::ios::fmtflags mode = G4cout.flags();
151  G4cout.setf(std::ios::fixed,std::ios::floatfield);
152  G4int prec = G4cout.precision(2);
153 
154  //run conditions
155  //
156  G4Material* material = fDetector->GetAbsorMaterial();
157  G4double density = material->GetDensity();
158  G4String partName = fParticle->GetParticleName();
159 
160  G4cout << "\n ======================== run summary =====================\n";
161  G4cout
162  << "\n The run is " << numberOfEvent << " "<< partName << " of "
163  << G4BestUnit(fEkin,"Energy") << " through "
164  << G4BestUnit(fDetector->GetAbsorRadius(),"Length") << " of "
165  << material->GetName() << " (density: "
166  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
167 
168  if (numberOfEvent == 0) {
169  G4cout.setf(mode,std::ios::floatfield);
170  G4cout.precision(prec);
171  return;
172  }
173 
174  fEdeposit /= numberOfEvent; fEdeposit2 /= numberOfEvent;
175  G4double rms = fEdeposit2 - fEdeposit*fEdeposit;
176  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
177 
178  G4cout.precision(3);
179  G4cout
180  << "\n Total Energy deposited = " << G4BestUnit(fEdeposit,"Energy")
181  << " +- " << G4BestUnit( rms,"Energy")
182  << G4endl;
183 
184  //compute track length of primary track
185  //
186  fTrackLen /= numberOfEvent; fTrackLen2 /= numberOfEvent;
187  rms = fTrackLen2 - fTrackLen*fTrackLen;
188  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
189 
190  G4cout.precision(3);
191  G4cout
192  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
193  << " +- " << G4BestUnit( rms,"Length");
194 
195  //compare with csda range
196  //
197  G4cout
198  << "\n Range from EmCalculator = " << G4BestUnit(fCsdaRange,"Length")
199  << " (from full dE/dx)" << G4endl;
200 
201 
202  //compute projected range of primary track
203  //
204  fProjRange /= numberOfEvent; fProjRange2 /= numberOfEvent;
205  rms = fProjRange2 - fProjRange*fProjRange;
206  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
207 
208  G4cout
209  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
210  << " +- " << G4BestUnit( rms,"Length")
211  << G4endl;
212 
213  //nb of steps and step size of primary track
214  //
215  G4double dNofEvents = double(numberOfEvent);
216  G4double fNbSteps = fNbOfSteps/dNofEvents,
217  fNbSteps2 = fNbOfSteps2/dNofEvents;
218  rms = fNbSteps2 - fNbSteps*fNbSteps;
219  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
220 
221  G4cout.precision(2);
222  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms;
223 
224  fStepSize /= numberOfEvent; fStepSize2 /= numberOfEvent;
225  rms = fStepSize2 - fStepSize*fStepSize;
226  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
227 
228  G4cout.precision(3);
229  G4cout
230  << "\t Step size= " << G4BestUnit(fStepSize,"Length")
231  << " +- " << G4BestUnit( rms,"Length")
232  << G4endl;
233 
234  // normalize histograms of longitudinal energy profile
235  //
236  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
237  G4int ih = 1;
238  G4double binWidth = analysisManager->GetH1Width(ih)
239  *analysisManager->GetH1Unit(ih);
240  G4double fac = (1./(numberOfEvent*binWidth))*(mm/MeV);
241  analysisManager->ScaleH1(ih,fac);
242 
243  // normalize histogram d(E/E0)/d(r/r0)
244  //
245  ih = 8;
246  binWidth = analysisManager->GetH1Width(ih);
247  fac = 1./(numberOfEvent*binWidth*fEkin);
248  analysisManager->ScaleH1(ih,fac);
249 
250  // reset default formats
251  G4cout.setf(mode,std::ios::floatfield);
252  G4cout.precision(prec);
253 }
254 
255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int numberOfEvent
Definition: G4Run.hh:59
static constexpr double mm
Definition: G4SIunits.hh:115
Run()
Definition: Run.cc:43
const G4String & GetName() const
Definition: G4Material.hh:178
void SetCsdaRange(G4int i, G4double value)
Definition: Run.cc:135
G4double GetDensity() const
Definition: G4Material.hh:180
tuple x
Definition: test.py:50
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4Material * GetAbsorMaterial(G4int i)
string material
Definition: eplot.py:19
void AddEdep(G4double val)
Definition: Run.hh:61
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
Definition: G4Run.hh:46
G4double energy(const ThreeVector &p, const G4double m)
G4double GetCsdaRange()
Definition: Run.cc:114
static const G4double fac
Detector construction class to define materials and geometry.
void EndOfRun()
Definition: Run.cc:147
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void AddTrackLength(G4double t)
Definition: Run.cc:102
void AddStepSize(G4int nb, G4double st)
Definition: Run.cc:118
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
Definition: Run.cc:77
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
virtual void Merge(const G4Run *)
Definition: Run.cc:115
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
~Run()
Definition: Run.cc:72
void AddProjRange(G4double x)
Definition: Run.hh:63