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