Geant4  10.03.p03
 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  fPenetration(0.), fPenetration2(0.),
54  fNbOfSteps(0), fNbOfSteps2(0),
55  fStepSize(0.), fStepSize2(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 
98 {
99  fPenetration += x;
100  fPenetration2 += x*x;
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104 
105 void Run::AddStepSize (G4int nb, G4double st)
106 {
107  fNbOfSteps += nb;
108  fNbOfSteps2 += nb*nb;
109  fStepSize += st ;
110  fStepSize2 += st*st;
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
115 void Run::Merge(const G4Run* run)
116 {
117  const Run* localRun = static_cast<const Run*>(run);
118 
119  // pass information about primary particle
120  fParticle = localRun->fParticle;
121  fEkin = localRun->fEkin;
122 
123  // accumulate sums
124  fEdeposit += localRun->fEdeposit;
125  fEdeposit2 += localRun->fEdeposit2;
126  fTrackLen += localRun->fTrackLen;
127  fTrackLen2 += localRun->fTrackLen2;
128  fProjRange += localRun->fProjRange;
129  fProjRange2 += localRun->fProjRange2;
130  fPenetration += localRun->fPenetration;
131  fPenetration2 += localRun->fPenetration2;
132  fNbOfSteps += localRun->fNbOfSteps ;
133  fNbOfSteps2 += localRun->fNbOfSteps2;
134  fStepSize += localRun->fStepSize;
135  fStepSize2 += localRun->fStepSize2;
136 
137  G4Run::Merge(run);
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
142 void Run::EndOfRun()
143 {
144  std::ios::fmtflags mode = G4cout.flags();
145  G4cout.setf(std::ios::fixed,std::ios::floatfield);
146  G4int prec = G4cout.precision(2);
147 
148  //run conditions
149  //
150  G4Material* material = fDetector->GetAbsorMaterial();
151  G4double density = material->GetDensity();
152  G4String partName = fParticle->GetParticleName();
153 
154  G4cout << "\n ======================= run summary ====================\n";
155  G4cout
156  << "\n The run is " << numberOfEvent << " "<< partName << " of "
157  << G4BestUnit(fEkin,"Energy") << " through a sphere of radius "
158  << G4BestUnit(fDetector->GetAbsorRadius(),"Length") << "of "
159  << material->GetName() << " (density: "
160  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
161 
162  if (numberOfEvent == 0) {
163  G4cout.setf(mode,std::ios::floatfield);
164  G4cout.precision(prec);
165  return;
166  }
167 
168  //compute track length of primary track
169  //
170  fTrackLen /= numberOfEvent; fTrackLen2 /= numberOfEvent;
171  G4double rmsTrack = fTrackLen2 - fTrackLen*fTrackLen;
172 
173  if (rmsTrack>0.) rmsTrack = std::sqrt(rmsTrack); else rmsTrack = 0.;
174 
175  G4cout.precision(3);
176  G4cout
177  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
178  << " +- "
179  << G4BestUnit( rmsTrack,"Length");
180 
181  //compute projected range of primary track
182  //
183  fProjRange /= numberOfEvent; fProjRange2 /= numberOfEvent;
184  G4double rmsProj = fProjRange2 - fProjRange*fProjRange;
185  if (rmsProj>0.) rmsProj = std::sqrt(rmsProj); else rmsProj = 0.;
186 
187  G4cout
188  << "\n Projected range = "
189  << G4BestUnit(fProjRange,"Length")
190  << " +- " << G4BestUnit( rmsProj,"Length");
191 
192  //compute penetration of primary track
193  //
194  fPenetration /= numberOfEvent; fPenetration2 /= numberOfEvent;
195  G4double rmsPene = fPenetration2 - fPenetration*fPenetration;
196  if (rmsPene>0.) rmsPene = std::sqrt(rmsPene); else rmsPene = 0.;
197 
198  G4cout
199  << "\n Penetration = "
200  << G4BestUnit(fPenetration,"Length")
201  << " +- " << G4BestUnit( rmsPene,"Length")
202  << G4endl;
203 
204  //
205 
206  //output file
207  FILE *myFile;
208  myFile = fopen ("range.txt","a");
209  fprintf (myFile, "%e %e %e %e %e %e %e\n",
210  fEkin/eV,
211  fTrackLen/nm,
212  rmsTrack/nm,
213  fProjRange/nm,
214  rmsProj/nm,
215  fPenetration/nm,
216  rmsPene/nm
217  );
218  fclose (myFile);
219 
220  // reset default formats
221  G4cout.setf(mode,std::ios::floatfield);
222  G4cout.precision(prec);
223 
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int numberOfEvent
Definition: G4Run.hh:59
Run()
Definition: Run.cc:43
const G4String & GetName() const
Definition: G4Material.hh:178
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
Definition: G4Run.hh:46
static constexpr double eV
Definition: G4SIunits.hh:215
static constexpr double nm
Definition: G4SIunits.hh:112
G4double energy(const ThreeVector &p, const G4double m)
void AddPenetration(G4double x)
Definition: Run.cc:97
Detector construction class to define materials and geometry.
void EndOfRun()
Definition: Run.cc:147
#define G4endl
Definition: G4ios.hh:61
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
~Run()
Definition: Run.cc:72
void AddProjRange(G4double x)
Definition: Run.hh:63