Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/radioactivedecay/Activation/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 "DetectorConstruction.hh"
36 #include "PrimaryGeneratorAction.hh"
37 #include "HistoManager.hh"
38 
39 #include "G4Threading.hh"
40 #include "G4AutoLock.hh"
41 #include "G4UnitsTable.hh"
42 #include "G4SystemOfUnits.hh"
43 
44 // mutex in a file scope
45 
46 namespace {
47  //Mutex to lock updating the global ion map
48  G4Mutex ionIdMapMutex = G4MUTEX_INITIALIZER;
49 }
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 
53 std::map<G4String,G4int> Run::fgIonMap;
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
59 : G4Run(),
60  fDetector(det), fParticle(nullptr), fEkin(0.)
61 {
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
68 Run::~Run()
69 { }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
73 void Run::Merge(std::map<G4String, ParticleData>& destinationMap,
74  const std::map<G4String, ParticleData>& sourceMap) const
75 {
76  for ( const auto& particleData : sourceMap ) {
77  G4String name = particleData.first;
78  const ParticleData& localData = particleData.second;
79  if ( destinationMap.find(name) == destinationMap.end()) {
80  destinationMap[name]
81  = ParticleData(localData.fCount,
82  localData.fEmean,
83  localData.fEmin,
84  localData.fEmax,
85  localData.fTmean);
86  }
87  else {
88  ParticleData& data = destinationMap[name];
89  data.fCount += localData.fCount;
90  data.fEmean += localData.fEmean;
91  G4double emin = localData.fEmin;
92  if (emin < data.fEmin) data.fEmin = emin;
93  G4double emax = localData.fEmax;
94  if (emax > data.fEmax) data.fEmax = emax;
95  data.fTmean = localData.fTmean;
96  }
97  }
98 }
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101 
103 {
104  fParticle = particle;
105  fEkin = energy;
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109 
110 void Run::CountProcesses(const G4VProcess* process)
111 {
112  G4String procName = process->GetProcessName();
113  std::map<G4String,G4int>::iterator it = fProcCounter.find(procName);
114  if ( it == fProcCounter.end()) {
115  fProcCounter[procName] = 1;
116  }
117  else {
118  fProcCounter[procName]++;
119  }
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123 
125 {
126  std::map<G4String, ParticleData>::iterator it = fParticleDataMap1.find(name);
127  if ( it == fParticleDataMap1.end()) {
128  fParticleDataMap1[name] = ParticleData(1, Ekin, Ekin, Ekin, meanLife);
129  }
130  else {
131  ParticleData& data = it->second;
132  data.fCount++;
133  data.fEmean += Ekin;
134  //update min max
135  G4double emin = data.fEmin;
136  if (Ekin < emin) data.fEmin = Ekin;
137  G4double emax = data.fEmax;
138  if (Ekin > emax) data.fEmax = Ekin;
139  data.fTmean = meanLife;
140  }
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144 
146 {
147  fEnergyDeposit += edep;
148  fEnergyDeposit2 += edep*edep;
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152 
153 void Run::AddEflow(G4double eflow)
154 {
155  fEnergyFlow += eflow;
156  fEnergyFlow2 += eflow*eflow;
157 }
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
161 {
162  std::map<G4String, ParticleData>::iterator it = fParticleDataMap2.find(name);
163  if ( it == fParticleDataMap2.end()) {
164  fParticleDataMap2[name] = ParticleData(1, Ekin, Ekin, Ekin, -1*ns);
165  }
166  else {
167  ParticleData& data = it->second;
168  data.fCount++;
169  data.fEmean += Ekin;
170  //update min max
171  G4double emin = data.fEmin;
172  if (Ekin < emin) data.fEmin = Ekin;
173  G4double emax = data.fEmax;
174  if (Ekin > emax) data.fEmax = Ekin;
175  data.fTmean = -1*ns;
176  }
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
182 {
183  G4AutoLock lock(&ionIdMapMutex);
184  // updating the global ion map needs to be locked
185 
186  std::map<G4String,G4int>::const_iterator it = fgIonMap.find(ionName);
187  if ( it == fgIonMap.end()) {
188  fgIonMap[ionName] = fgIonId;
189  if (fgIonId < kMaxHisto2) fgIonId++;
190  }
191  return fgIonMap[ionName];
192 }
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195 
196 void Run::Merge(const G4Run* run)
197 {
198  const Run* localRun = static_cast<const Run*>(run);
199 
200  //primary particle info
201  //
202  fParticle = localRun->fParticle;
203  fEkin = localRun->fEkin;
204 
205  // accumulate sums
206  //
207  fEnergyDeposit += localRun->fEnergyDeposit;
208  fEnergyDeposit2 += localRun->fEnergyDeposit2;
209  fEnergyFlow += localRun->fEnergyFlow;
210  fEnergyFlow2 += localRun->fEnergyFlow2;
211 
212  //map: processes count
213  for ( const auto& procCounter : localRun->fProcCounter ) {
214  G4String procName = procCounter.first;
215  G4int localCount = procCounter.second;
216  if ( fProcCounter.find(procName) == fProcCounter.end()) {
217  fProcCounter[procName] = localCount;
218  }
219  else {
220  fProcCounter[procName] += localCount;
221  }
222  }
223 
224  //map: created particles count
226 
227  //map: particles flux count
229 
230  G4Run::Merge(run);
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
234 
235 void Run::EndOfRun()
236 {
237  G4int prec = 5, wid = prec + 2;
238  G4int dfprec = G4cout.precision(prec);
239 
240  //run condition
241  //
242  G4Material* material = fDetector->GetAbsorMaterial();
243  G4double density = material->GetDensity();
244 
245  G4String Particle = fParticle->GetParticleName();
246  G4cout << "\n The run is " << numberOfEvent << " "<< Particle << " of "
247  << G4BestUnit(fEkin,"Energy") << " through "
248  << G4BestUnit(fDetector->GetAbsorThickness(),"Length") << " of "
249  << material->GetName() << " (density: "
250  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
251 
252  if (numberOfEvent == 0) { G4cout.precision(dfprec); return;}
253 
254  //frequency of processes
255  //
256  G4cout << "\n Process calls frequency :" << G4endl;
257  G4int index = 0;
258  for ( const auto& procCounter : fProcCounter ) {
259  G4String procName = procCounter.first;
260  G4int count = procCounter.second;
261  G4String space = " "; if (++index%3 == 0) space = "\n";
262  G4cout << " " << std::setw(20) << procName << "="<< std::setw(7) << count
263  << space;
264  }
265  G4cout << G4endl;
266 
267  //particles count
268  //
269  G4cout << "\n List of generated particles:" << G4endl;
270 
271  for ( const auto& particleData : fParticleDataMap1 ) {
272  G4String name = particleData.first;
273  ParticleData data = particleData.second;
274  G4int count = data.fCount;
275  G4double eMean = data.fEmean/count;
276  G4double eMin = data.fEmin;
277  G4double eMax = data.fEmax;
278  G4double meanLife = data.fTmean;
279 
280  G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
281  << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy")
282  << "\t( " << G4BestUnit(eMin, "Energy")
283  << " --> " << G4BestUnit(eMax, "Energy") << ")";
284  if (meanLife >= 0.)
285  G4cout << "\tmean life = " << G4BestUnit(meanLife, "Time") << G4endl;
286  else G4cout << "\tstable" << G4endl;
287  }
288 
289  // compute mean Energy deposited and rms
290  //
291  G4int TotNbofEvents = numberOfEvent;
292  fEnergyDeposit /= TotNbofEvents; fEnergyDeposit2 /= TotNbofEvents;
294  if (rmsEdep>0.) rmsEdep = std::sqrt(rmsEdep);
295  else rmsEdep = 0.;
296 
297  G4cout << "\n Mean energy deposit per event = "
298  << G4BestUnit(fEnergyDeposit,"Energy") << "; rms = "
299  << G4BestUnit(rmsEdep, "Energy")
300  << G4endl;
301 
302  // compute mean Energy flow and rms
303  //
304  fEnergyFlow /= TotNbofEvents; fEnergyFlow2 /= TotNbofEvents;
306  if (rmsEflow>0.) rmsEflow = std::sqrt(rmsEflow);
307  else rmsEflow = 0.;
308 
309  G4cout << " Mean energy flow per event = "
310  << G4BestUnit(fEnergyFlow,"Energy") << "; rms = "
311  << G4BestUnit(rmsEflow, "Energy")
312  << G4endl;
313 
314  //particles flux
315  //
316  G4cout << "\n List of particles emerging from the target :" << G4endl;
317 
318  for ( const auto& particleData : fParticleDataMap2 ) {
319  G4String name = particleData.first;
320  ParticleData data = particleData.second;
321  G4int count = data.fCount;
322  G4double eMean = data.fEmean/count;
323  G4double eMin = data.fEmin;
324  G4double eMax = data.fEmax;
325  G4double Eflow = data.fEmean/TotNbofEvents;
326 
327  G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
328  << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy")
329  << "\t( " << G4BestUnit(eMin, "Energy")
330  << " --> " << G4BestUnit(eMax, "Energy")
331  << ") \tEflow/event = " << G4BestUnit(Eflow, "Energy") << G4endl;
332  }
333 
334  //histogram Id for populations
335  //
336  G4cout << "\n histo Id for populations :" << G4endl;
337 
338  // Update the histogram titles according to the ion map
339  // and print new titles
340  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
341  for ( const auto& ionMapElement : fgIonMap ) {
342  G4String ionName = ionMapElement.first;
343  G4int h1Id = ionMapElement.second;
344  // print new titles
345  G4cout << " " << std::setw(20) << ionName << " id = "<< std::setw(3) << h1Id
346  << G4endl;
347 
348  // update histogram ids
349  if ( ! analysisManager->GetH1(h1Id) ) continue;
350  // Skip inactive histograms, this is not necessary
351  // but it makes the code safe wrt modifications in future
352  G4String title = analysisManager->GetH1Title(h1Id);
353  title = ionName + title;
354  analysisManager->SetH1Title(h1Id, title);
355  }
356  G4cout << G4endl;
357 
358  //normalize histograms
359  G4int ih = 2;
360  G4double binWidth = analysisManager->GetH1Width(ih);
361  G4double fac = (1./(numberOfEvent*binWidth))*(mm/MeV);
362  analysisManager->ScaleH1(ih,fac);
363 
364  for (ih=14; ih<24; ih++) {
365  binWidth = analysisManager->GetH1Width(ih);
366  G4double unit = analysisManager->GetH1Unit(ih);
367  fac = (second/(binWidth*unit));
368  analysisManager->ScaleH1(ih,fac);
369  }
370 
371  //remove all contents in fProcCounter, fCount
372  fProcCounter.clear();
373  fParticleDataMap1.clear();
374  fParticleDataMap2.clear();
375  fgIonMap.clear();
376 
377  //restore default format
378  G4cout.precision(dfprec);
379 }
380 
381 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
const XML_Char * name
Definition: expat.h:151
static const double prec
Definition: RanecuEngine.cc:58
std::map< G4String, ParticleData > fParticleDataMap1
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
DetectorConstruction * fDetector
static constexpr double MeV
Definition: G4SIunits.hh:214
std::vector< G4double > fEnergyDeposit[kMaxAbsor]
G4int numberOfEvent
Definition: G4Run.hh:59
static constexpr double mm
Definition: G4SIunits.hh:115
void ParticleFlux(G4String, G4double)
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
void AddEflow(G4double eflow)
static const G4double emax
std::vector< G4double > fEnergyFlow
static constexpr double second
Definition: G4SIunits.hh:157
void ParticleCount(G4String, G4double)
const XML_Char const XML_Char * data
Definition: expat.h:268
const G4String & GetName() const
Definition: G4Material.hh:179
static std::map< G4String, G4int > fgIonMap
double G4double
Definition: G4Types.hh:76
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:88
double energy
Definition: plottest35.C:25
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
virtual void Merge(const G4Run *)
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
std::map< G4String, G4int > fProcCounter
Double_t edep
Definition: G4Run.hh:46
G4ParticleDefinition * fParticle
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
std::map< G4String, ParticleData > fParticleDataMap2
int G4int
Definition: G4Types.hh:78
void CountProcesses(G4String procName)
static const G4double fac
G4GLOB_DLL std::ostream G4cout
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int first(char) const
#define ns
Definition: xmlparse.cc:614
Simple detector construction with a box volume placed in a world.
G4double GetDensity() const
Definition: G4Material.hh:181
std::mutex G4Mutex
Definition: G4Threading.hh:84