Geant4  10.02.p01
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  fNbInelastic(0), fNbInelastic2(0),
51  fEdeposit(0.), fEdeposit2(0.),
52  fTrackLen(0.), fTrackLen2(0.),
53  fProjRange(0.), fProjRange2(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 
74 {
75  fNbInelastic += nb;
76  fNbInelastic2 += nb*nb;
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
81 void Run::AddEdep (G4double e)
82 {
83  fEdeposit += e;
84  fEdeposit2 += e*e;
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90 {
91  fTrackLen += t;
92  fTrackLen2 += t*t;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98 {
99  fProjRange += x;
100  fProjRange2 += 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  fNbInelastic += localRun->fNbInelastic;
125  fNbInelastic2 += localRun->fNbInelastic2;
126  fEdeposit += localRun->fEdeposit;
127  fEdeposit2 += localRun->fEdeposit2;
128  fTrackLen += localRun->fTrackLen;
129  fTrackLen2 += localRun->fTrackLen2;
130  fProjRange += localRun->fProjRange;
131  fProjRange2 += localRun->fProjRange2;
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 
169 
171  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
172 
173  G4cout.precision(3);
174  G4cout
175  << "\n Nb of ionisations = " << fNbInelastic
176  << " +- " << rms
177  << G4endl;
178 
179  G4cout.precision(3);
180  G4cout
181  << "\n w = " << G4BestUnit((fEkin)/fNbInelastic,"Energy")
182  << " +- " << G4BestUnit((fEkin)*rms/(fNbInelastic*fNbInelastic),"Energy")
183  << G4endl;
184 
185  //output file
186  if(fNbInelastic>0.)
187  {
188  FILE *myFile;
189  myFile = fopen ("wvalue.txt","a");
190  fprintf (myFile, "%e %e %e %e %e \n", fEkin/eV, fNbInelastic, rms, fEkin/eV/fNbInelastic,
191  (fEkin/eV)*rms/(fNbInelastic*fNbInelastic) );
192  fclose (myFile);
193  }
194  //
195 
197  rms = fEdeposit2 - fEdeposit*fEdeposit;
198  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
199 
200  G4cout.precision(3);
201  G4cout
202  << "\n Total Energy deposited = " << G4BestUnit(fEdeposit,"Energy")
203  << " +- " << G4BestUnit( rms,"Energy")
204  << G4endl;
205 
206  //compute track length of primary track
207  //
209  rms = fTrackLen2 - fTrackLen*fTrackLen;
210  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
211 
212  G4cout.precision(3);
213  G4cout
214  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
215  << " +- " << G4BestUnit( rms,"Length");
216 
217  //compute projected range of primary track
218  //
221  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
222 
223  G4cout
224  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
225  << " +- " << G4BestUnit( rms,"Length")
226  << G4endl;
227 
228  //nb of steps and step size of primary track
229  //
230  G4double dNofEvents = double(numberOfEvent);
231  G4double fNbSteps = fNbOfSteps/dNofEvents,
232  fNbSteps2 = fNbOfSteps2/dNofEvents;
233  rms = fNbSteps2 - fNbSteps*fNbSteps;
234  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
235 
236  G4cout.precision(2);
237  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms
238  << G4endl;
239 
241  rms = fStepSize2 - fStepSize*fStepSize;
242  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
243 
244  G4cout.precision(3);
245  G4cout
246  << "\n Step size = " << G4BestUnit(fStepSize,"Length")
247  << " +- " << G4BestUnit( rms,"Length")
248  << G4endl;
249 
250  // normalize histograms of longitudinal energy profile
251  //
252  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
253  G4int ih = 1;
254  G4double binWidth = analysisManager->GetH1Width(ih);
255  G4double fac = (1./(numberOfEvent*binWidth))*(mm/MeV);
256  analysisManager->ScaleH1(ih,fac);
257 
258  // reset default formats
259  G4cout.setf(mode,std::ios::floatfield);
260  G4cout.precision(prec);
261 
262 }
263 
264 //....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:211
G4double fNbInelastic2
Definition: Run.hh:72
G4double fEkin
Definition: Run.hh:65
G4double fNbInelastic
Definition: Run.hh:72
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
G4double GetDensity() const
Definition: G4Material.hh:180
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double fProjRange2
Definition: Run.hh:78
int G4int
Definition: G4Types.hh:78
void AddInelastic(G4int nb)
Definition: Run.cc:73
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
Definition: G4Run.hh:46
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:66
G4double fProjRange
Definition: Run.hh:78
static const double eV
Definition: G4SIunits.hh:212
G4double fEdeposit2
Definition: Run.hh:75
G4int fNbOfSteps2
Definition: Run.hh:78
G4double energy(const ThreeVector &p, const G4double m)
G4double fStepSize
Definition: Run.hh:79
const G4double x[NPOINTSGL]
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 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
virtual void Merge(const G4Run *)
Definition: Run.cc:115
static const double mm
Definition: G4SIunits.hh:114
G4double fTrackLen2
Definition: Run.hh:76
~Run()
Definition: Run.cc:72
void AddProjRange(G4double x)
Definition: Run.hh:63