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 "PrimaryGeneratorAction.hh"
36 #include "HistoManager.hh"
37 
38 #include "G4SystemOfUnits.hh"
39 #include "G4UnitsTable.hh"
40 #include "G4PhysicalConstants.hh"
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43 
44 Run::Run()
45 : G4Run(),
46  fParticle(0), fEkin(0.),
47  fDecayCount(0), fTimeCount(0), fPrimaryTime(0.),
48  fTimeWindow1(0.), fTimeWindow2(0.)
49 {
50  fEkinTot[0] = fPbalance[0] = fEventTime[0] = fEvisEvent[0] = 0. ;
51  fEkinTot[1] = fPbalance[1] = fEventTime[1] = fEvisEvent[1] = DBL_MAX;
52  fEkinTot[2] = fPbalance[2] = fEventTime[2] = fEvisEvent[2] = 0. ;
53 }
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
56 Run::~Run()
57 { }
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 
62 {
63  fParticle = particle;
64  fEkin = energy;
65 }
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70 {
71  std::map<G4String, ParticleData>::iterator it = fParticleDataMap.find(name);
72  if ( it == fParticleDataMap.end()) {
73  fParticleDataMap[name] = ParticleData(1, Ekin, Ekin, Ekin);
74  }
75  else {
76  ParticleData& data = it->second;
77  data.fCount++;
78  data.fEmean += Ekin;
79  //update min max
80  G4double emin = data.fEmin;
81  if (Ekin < emin) data.fEmin = Ekin;
82  G4double emax = data.fEmax;
83  if (Ekin > emax) data.fEmax = Ekin;
84  }
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90 {
91  fTimeWindow1 = t1;
92  fTimeWindow2 = t2;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98  G4bool life2, G4bool decay)
99 {
100  std::map<G4String, ActivityData>::iterator it = fActivityMap.find(name);
101  if ( it == fActivityMap.end()) {
102  G4int n1(0), n2(0), nd(0);
103  if(life1) n1 = 1;
104  if(life2) n2 = 1;
105  if(decay) nd = 1;
106  fActivityMap[name] = ActivityData(n1, n2, nd);
107  }
108  else {
109  ActivityData& data = it->second;
110  if(life1) data.fNlife_t1++;
111  if(life2) data.fNlife_t2++;
112  if(decay) data.fNdecay_t1t2++;
113  }
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117 
119 {
120  fDecayCount++;
121  fEkinTot[0] += Ekin;
122  //update min max
123  if (fDecayCount == 1) fEkinTot[1] = fEkinTot[2] = Ekin;
124  if (Ekin < fEkinTot[1]) fEkinTot[1] = Ekin;
125  if (Ekin > fEkinTot[2]) fEkinTot[2] = Ekin;
126 
127  fPbalance[0] += Pbal;
128  //update min max
129  if (fDecayCount == 1) fPbalance[1] = fPbalance[2] = Pbal;
130  if (Pbal < fPbalance[1]) fPbalance[1] = Pbal;
131  if (Pbal > fPbalance[2]) fPbalance[2] = Pbal;
132 }
133 
134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135 
137 {
138  fTimeCount++;
139  fEventTime[0] += time;
140  if (fTimeCount == 1) fEventTime[1] = fEventTime[2] = time;
141  if (time < fEventTime[1]) fEventTime[1] = time;
142  if (time > fEventTime[2]) fEventTime[2] = time;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146 
148 {
149  fPrimaryTime += ptime;
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
155 {
156  fEvisEvent[0] += Evis;
157  if (fTimeCount == 1) fEvisEvent[1] = fEvisEvent[2] = Evis;
158  if (Evis < fEvisEvent[1]) fEvisEvent[1] = Evis;
159  if (Evis > fEvisEvent[2]) fEvisEvent[2] = Evis;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163 
164 void Run::Merge(const G4Run* run)
165 {
166  const Run* localRun = static_cast<const Run*>(run);
167 
168  //primary particle info
169  //
170  fParticle = localRun->fParticle;
171  fEkin = localRun->fEkin;
172 
173  // accumulate sums
174  //
175  fDecayCount += localRun->fDecayCount;
176  fTimeCount += localRun->fTimeCount;
177  fPrimaryTime += localRun->fPrimaryTime;
178 
179  fEkinTot[0] += localRun->fEkinTot[0];
180  fPbalance[0] += localRun->fPbalance[0];
181  fEventTime[0] += localRun->fEventTime[0];
182  fEvisEvent[0] += localRun->fEvisEvent[0];
183 
184  G4double min,max;
185  min = localRun->fEkinTot[1]; max = localRun->fEkinTot[2];
186  if (fEkinTot[1] > min) fEkinTot[1] = min;
187  if (fEkinTot[2] < max) fEkinTot[2] = max;
188  //
189  min = localRun->fPbalance[1]; max = localRun->fPbalance[2];
190  if (fPbalance[1] > min) fPbalance[1] = min;
191  if (fPbalance[2] < max) fPbalance[2] = max;
192  //
193  min = localRun->fEventTime[1]; max = localRun->fEventTime[2];
194  if (fEventTime[1] > min) fEventTime[1] = min;
195  if (fEventTime[2] < max) fEventTime[2] = max;
196  //
197  min = localRun->fEvisEvent[1]; max = localRun->fEvisEvent[2];
198  if (fEvisEvent[1] > min) fEvisEvent[1] = min;
199  if (fEvisEvent[2] < max) fEvisEvent[2] = max;
200 
201  //maps
202  std::map<G4String,ParticleData>::const_iterator itn;
203  for (itn = localRun->fParticleDataMap.begin();
204  itn != localRun->fParticleDataMap.end(); ++itn) {
205 
206  G4String name = itn->first;
207  const ParticleData& localData = itn->second;
208  if ( fParticleDataMap.find(name) == fParticleDataMap.end()) {
209  fParticleDataMap[name]
210  = ParticleData(localData.fCount,
211  localData.fEmean,
212  localData.fEmin,
213  localData.fEmax);
214  }
215  else {
216  ParticleData& data = fParticleDataMap[name];
217  data.fCount += localData.fCount;
218  data.fEmean += localData.fEmean;
219  G4double emin = localData.fEmin;
220  if (emin < data.fEmin) data.fEmin = emin;
221  G4double emax = localData.fEmax;
222  if (emax > data.fEmax) data.fEmax = emax;
223  }
224  }
225 
226  //activity
227  fTimeWindow1 = localRun->fTimeWindow1;
228  fTimeWindow2 = localRun->fTimeWindow2;
229 
230  std::map<G4String,ActivityData>::const_iterator ita;
231  for (ita = localRun->fActivityMap.begin();
232  ita != localRun->fActivityMap.end(); ++ita) {
233 
234  G4String name = ita->first;
235  const ActivityData& localData = ita->second;
236  if ( fActivityMap.find(name) == fActivityMap.end()) {
237  fActivityMap[name]
238  = ActivityData(localData.fNlife_t1,
239  localData.fNlife_t2,
240  localData.fNdecay_t1t2);
241  } else {
242  ActivityData& data = fActivityMap[name];
243  data.fNlife_t1 += localData.fNlife_t1;
244  data.fNlife_t2 += localData.fNlife_t2;
245  data.fNdecay_t1t2 += localData.fNdecay_t1t2;
246  }
247  }
248 
249  G4Run::Merge(run);
250 }
251 
252 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
253 
254 void Run::EndOfRun()
255 {
256  G4int nbEvents = numberOfEvent;
257  G4String partName = fParticle->GetParticleName();
258 
259  G4cout << "\n ======================== run summary ======================";
260  G4cout << "\n The run was " << nbEvents << " " << partName << " of "
261  << G4BestUnit(fEkin,"Energy");
262  G4cout << "\n ===========================================================\n";
263  G4cout << G4endl;
264  if (nbEvents == 0) { return; }
265 
266  G4int prec = 4, wid = prec + 2;
267  G4int dfprec = G4cout.precision(prec);
268 
269  //particle count
270  //
271  G4cout << " Nb of generated particles: \n" << G4endl;
272 
273  std::map<G4String,ParticleData>::iterator it;
274  for (it = fParticleDataMap.begin(); it != fParticleDataMap.end(); it++) {
275  G4String name = it->first;
276  ParticleData data = it->second;
277  G4int count = data.fCount;
278  G4double eMean = data.fEmean/count;
279  G4double eMin = data.fEmin;
280  G4double eMax = data.fEmax;
281 
282  G4cout << " " << std::setw(15) << name << ": " << std::setw(7) << count
283  << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy")
284  << "\t( " << G4BestUnit(eMin, "Energy")
285  << " --> " << G4BestUnit(eMax, "Energy")
286  << ")" << G4endl;
287  }
288 
289  //energy momentum balance
290  //
291 
292  if (fDecayCount > 0) {
293  G4double Ebmean = fEkinTot[0]/fDecayCount;
294  G4double Pbmean = fPbalance[0]/fDecayCount;
295 
296  G4cout << "\n Ekin Total (Q single decay): mean = "
297  << std::setw(wid) << G4BestUnit(Ebmean, "Energy")
298  << "\t( " << G4BestUnit(fEkinTot[1], "Energy")
299  << " --> " << G4BestUnit(fEkinTot[2], "Energy")
300  << ")" << G4endl;
301 
302  G4cout << "\n Momentum balance (excluding gamma desexcitation): mean = "
303  << std::setw(wid) << G4BestUnit(Pbmean, "Energy")
304  << "\t( " << G4BestUnit(fPbalance[1], "Energy")
305  << " --> " << G4BestUnit(fPbalance[2], "Energy")
306  << ")" << G4endl;
307  }
308 
309  //total time of life
310  //
311  if (fTimeCount > 0) {
312  G4double Tmean = fEventTime[0]/fTimeCount;
313  G4double halfLife = Tmean*std::log(2.);
314 
315  G4cout << "\n Total time of life (full chain): mean = "
316  << std::setw(wid) << G4BestUnit(Tmean, "Time")
317  << " half-life = "
318  << std::setw(wid) << G4BestUnit(halfLife, "Time")
319  << " ( " << G4BestUnit(fEventTime[1], "Time")
320  << " --> " << G4BestUnit(fEventTime[2], "Time")
321  << ")" << G4endl;
322  }
323 
324  //total visible Energy
325  //
326  if (fTimeCount > 0) {
327  G4double Evmean = fEvisEvent[0]/fTimeCount;
328 
329  G4cout << "\n Total visible energy (full chain) : mean = "
330  << std::setw(wid) << G4BestUnit(Evmean, "Energy")
331  << " ( " << G4BestUnit(fEvisEvent[1], "Energy")
332  << " --> " << G4BestUnit(fEvisEvent[2], "Energy")
333  << ")" << G4endl;
334  }
335 
336  //activity of primary ion
337  //
338  G4double pTimeMean = fPrimaryTime/nbEvents;
339  G4double molMass = fParticle->GetAtomicMass()*g/mole;
340  G4double nAtoms_perUnitOfMass = Avogadro/molMass;
341  G4double Activity_perUnitOfMass = 0.0;
342  if (pTimeMean > 0.0)
343  { Activity_perUnitOfMass = nAtoms_perUnitOfMass/pTimeMean; }
344 
345  G4cout << "\n Activity of " << partName << " = "
346  << std::setw(wid) << Activity_perUnitOfMass*g/becquerel
347  << " Bq/g (" << Activity_perUnitOfMass*g/curie
348  << " Ci/g) \n"
349  << G4endl;
350 
351 
352  //activities in time window
353  //
354  if (fTimeWindow2 > 0.) {
355  G4double dt = fTimeWindow2 - fTimeWindow1;
356  G4cout << " Activities in time window [t1, t2] = ["
357  << G4BestUnit(fTimeWindow1, "Time") << ", "
358  << G4BestUnit(fTimeWindow2, "Time") << "] (delta time = "
359  << G4BestUnit(dt, "Time") << ") : \n" << G4endl;
360 
361  std::map<G4String,ActivityData>::iterator ita;
362  for (ita = fActivityMap.begin(); ita != fActivityMap.end(); ita++) {
363  G4String name = ita->first;
364  ActivityData data = ita->second;
365  G4int n1 = data.fNlife_t1;
366  G4int n2 = data.fNlife_t2;
367  G4int ndecay = data.fNdecay_t1t2;
368  G4double actv = ndecay/dt;
369 
370  G4cout << " " << std::setw(15) << name << ": "
371  << " n(t1) = " << std::setw(7) << n1
372  << "\tn(t2) = " << std::setw(7) << n2
373  << "\t decays = " << std::setw(7) << ndecay
374  << " ---> <actv> = " << G4BestUnit(actv, "Activity") << "\n";
375  }
376  }
377  G4cout << G4endl;
378 
379  //normalize histograms
380  //
381  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
382  G4double factor = 100./nbEvents;
383  analysisManager->ScaleH1(1,factor);
384  analysisManager->ScaleH1(2,factor);
385  analysisManager->ScaleH1(3,factor);
386  analysisManager->ScaleH1(4,factor);
387  analysisManager->ScaleH1(5,factor);
388 
389  // remove all contents in fParticleDataMap
390  //
391  fParticleDataMap.clear();
392  fActivityMap.clear();
393 
394  // restore default precision
395  //
396  G4cout.precision(dfprec);
397 }
398 
399 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
const XML_Char * name
Definition: expat.h:151
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int numberOfEvent
Definition: G4Run.hh:59
void PrimaryTiming(G4double)
Definition: Run.cc:147
void CountInTimeWindow(G4String, G4bool, G4bool, G4bool)
Definition: Run.cc:97
Run()
Definition: Run.cc:43
G4int first(char) const
static constexpr double becquerel
Definition: G4SIunits.hh:291
void EventTiming(G4double)
Definition: Run.cc:136
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void Balance(G4double)
Definition: Run.cc:133
int G4int
Definition: G4Types.hh:78
float Avogadro
Definition: hepunit.py:253
const G4String & GetParticleName() const
void EvisEvent(G4double)
Definition: Run.cc:154
void SetTimeWindow(G4double, G4double)
Definition: Run.cc:89
const XML_Char const XML_Char * data
Definition: expat.h:268
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
bool G4bool
Definition: G4Types.hh:79
Definition: G4Run.hh:46
G4int GetAtomicMass() const
static const G4double emax
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
tuple t1
Definition: plottest35.py:33
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void EndOfRun()
Definition: Run.cc:147
#define G4endl
Definition: G4ios.hh:61
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
Definition: Run.cc:77
static constexpr double curie
Definition: G4SIunits.hh:292
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
#define DBL_MAX
Definition: templates.hh:83
virtual void Merge(const G4Run *)
Definition: Run.cc:115
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
static constexpr double mole
Definition: G4SIunits.hh:286
~Run()
Definition: Run.cc:72
void ParticleCount(G4String, G4double)
Definition: Run.cc:114