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 
46 Run::Run(const DetectorConstruction* detector)
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 { }
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
59 Run::~Run()
60 { }
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
65 {
66  fParticle = particle;
67  fEkin = energy;
68 }
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
72 void Run::AddEdep (G4double e)
73 {
74  fEdeposit += e;
75  fEdeposit2 += e*e;
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
81 {
82  fTrackLen += t;
83  fTrackLen2 += t*t;
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 {
90  fProjRange += x;
91  fProjRange2 += x*x;
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95 
96 void Run::AddStepSize (G4int nb, G4double st)
97 {
98  fNbOfSteps += nb;
99  fNbOfSteps2 += nb*nb;
100  fStepSize += st ;
101  fStepSize2 += st*st;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
106 void Run::Merge(const G4Run* run)
107 {
108  const Run* localRun = static_cast<const Run*>(run);
109 
110  // pass information about primary particle
111  fParticle = localRun->fParticle;
112  fEkin = localRun->fEkin;
113 
114  // accumulate sums
115  fEdeposit += localRun->fEdeposit;
116  fEdeposit2 += localRun->fEdeposit2;
117  fTrackLen += localRun->fTrackLen;
118  fTrackLen2 += localRun->fTrackLen2;
119  fProjRange += localRun->fProjRange;
120  fProjRange2 += localRun->fProjRange2;
121  fNbOfSteps += localRun->fNbOfSteps ;
122  fNbOfSteps2 += localRun->fNbOfSteps2;
123  fStepSize += localRun->fStepSize;
124  fStepSize2 += localRun->fStepSize2;
125 
126  G4Run::Merge(run);
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
131 void Run::EndOfRun()
132 {
133  std::ios::fmtflags mode = G4cout.flags();
134  G4cout.setf(std::ios::fixed,std::ios::floatfield);
135  G4int prec = G4cout.precision(2);
136 
137  //run conditions
138  //
139  G4Material* material = fDetector->GetAbsorMaterial();
140  G4double density = material->GetDensity();
141  G4String partName = fParticle->GetParticleName();
142 
143  G4cout << "\n ======================== run summary =====================\n";
144  G4cout
145  << "\n The run is " << numberOfEvent << " "<< partName << " of "
146  << G4BestUnit(fEkin,"Energy") << " through a sphere of radius "
147  << G4BestUnit(fDetector->GetAbsorRadius(),"Length") << "of "
148  << material->GetName() << " (density: "
149  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
150 
151  if (numberOfEvent == 0) {
152  G4cout.setf(mode,std::ios::floatfield);
153  G4cout.precision(prec);
154  return;
155  }
156 
157  fEdeposit /= numberOfEvent; fEdeposit2 /= numberOfEvent;
158  G4double rms = fEdeposit2 - fEdeposit*fEdeposit;
159  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
160 
161  G4cout.precision(3);
162  G4cout
163  << "\n Total Energy deposited = " << G4BestUnit(fEdeposit,"Energy")
164  << " +- " << G4BestUnit( rms,"Energy")
165  << G4endl;
166 
167  G4double sValue=fEdeposit/fDetector->GetAbsorMass();
168  G4double rmsSValue=rms/fDetector->GetAbsorMass();
169 
170  G4cout.precision(3);
171  G4cout
172  << "\n S value = " << sValue/gray << " Gy/Bq.s "
173  << " +- " << rmsSValue/gray
174  << " Gy/Bq.s "
175  << G4endl;
176 
177  //compute track length of primary track
178  //
179  fTrackLen /= numberOfEvent; fTrackLen2 /= numberOfEvent;
180  rms = fTrackLen2 - fTrackLen*fTrackLen;
181  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
182 
183  G4cout.precision(3);
184  G4cout
185  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
186  << " +- " << G4BestUnit( rms,"Length");
187 
188  //compute projected range of primary track
189  //
190  fProjRange /= numberOfEvent; fProjRange2 /= numberOfEvent;
191  rms = fProjRange2 - fProjRange*fProjRange;
192  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
193 
194  G4cout
195  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
196  << " +- " << G4BestUnit( rms,"Length")
197  << G4endl;
198 
199  //nb of steps and step size of primary track
200  //
201  G4double dNofEvents = double(numberOfEvent);
202  G4double fNbSteps = fNbOfSteps/dNofEvents,
203  fNbSteps2 = fNbOfSteps2/dNofEvents;
204  rms = fNbSteps2 - fNbSteps*fNbSteps;
205  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
206 
207  G4cout.precision(2);
208  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms
209  << G4endl;
210 
211  fStepSize /= numberOfEvent; fStepSize2 /= numberOfEvent;
212  rms = fStepSize2 - fStepSize*fStepSize;
213  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
214 
215  G4cout.precision(3);
216  G4cout
217  << "\n Step size = " << G4BestUnit(fStepSize,"Length")
218  << " +- " << G4BestUnit( rms,"Length")
219  << G4endl;
220 
221  // normalize histograms of longitudinal energy profile
222  //
223  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
224  G4int ih = 1;
225  G4double binWidth = analysisManager->GetH1Width(ih);
226  G4double fac = (1./(numberOfEvent*binWidth))*(mm/MeV);
227  analysisManager->ScaleH1(ih,fac);
228 
229  // reset default formats
230  G4cout.setf(mode,std::ios::floatfield);
231  G4cout.precision(prec);
232 
233  //output file
234  FILE *myFile;
235  myFile = fopen ("s.txt","a");
236  fprintf (myFile, "%e %e %e %e \n", fDetector->GetAbsorRadius()/nm,fEkin/eV,
237  sValue/gray, rmsSValue/gray );
238  fclose (myFile);
239 
240 }
241 
242 //....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
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
static constexpr double gray
Definition: G4SIunits.hh:309
Definition: G4Run.hh:46
static constexpr double eV
Definition: G4SIunits.hh:215
G4double GetAbsorMass() const
static constexpr double nm
Definition: G4SIunits.hh:112
G4double energy(const ThreeVector &p, const G4double m)
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