Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/hadronic/Hadr04/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 "G4UnitsTable.hh"
40 #include "G4SystemOfUnits.hh"
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43 
45 : G4Run(),
46  fDetector(det), fParticle(0), fEkin(0.),
47  fNbStep1(0), fNbStep2(0),
48  fTrackLen1(0.), fTrackLen2(0.),
49  fTime1(0.),fTime2(0.)
50 { }
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 
53 Run::~Run()
54 { }
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
59 {
60  fParticle = particle;
61  fEkin = energy;
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
66 void Run::CountProcesses(const G4VProcess* process)
67 {
68  G4String procName = process->GetProcessName();
69  std::map<G4String,G4int>::iterator it = fProcCounter.find(procName);
70  if ( it == fProcCounter.end()) {
71  fProcCounter[procName] = 1;
72  }
73  else {
74  fProcCounter[procName]++;
75  }
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
81 {
82  std::map<G4String, ParticleData>::iterator it = fParticleDataMap.find(name);
83  if ( it == fParticleDataMap.end()) {
84  fParticleDataMap[name] = ParticleData(1, Ekin, Ekin, Ekin);
85  }
86  else {
87  ParticleData& data = it->second;
88  data.fCount++;
89  data.fEmean += Ekin;
90  //update min max
91  G4double emin = data.fEmin;
92  if (Ekin < emin) data.fEmin = Ekin;
93  G4double emax = data.fEmax;
94  if (Ekin > emax) data.fEmax = Ekin;
95  }
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99 
100 void Run::SumTrackLength(G4int nstep1, G4int nstep2,
101  G4double trackl1, G4double trackl2,
102  G4double time1, G4double time2)
103 {
104  fNbStep1 += nstep1; fNbStep2 += nstep2;
105  fTrackLen1 += trackl1; fTrackLen2 += trackl2;
106  fTime1 += time1; fTime2 += time2;
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110 
111 void Run::Merge(const G4Run* run)
112 {
113  const Run* localRun = static_cast<const Run*>(run);
114 
115  //primary particle info
116  //
117  fParticle = localRun->fParticle;
118  fEkin = localRun->fEkin;
119 
120  // accumulate sums
121  //
122  fNbStep1 += localRun->fNbStep1;
123  fNbStep2 += localRun->fNbStep2;
124  fTrackLen1 += localRun->fTrackLen1;
125  fTrackLen2 += localRun->fTrackLen2;
126  fTime1 += localRun->fTime1;
127  fTime2 += localRun->fTime2;
128 
129  //map: processes count
130  std::map<G4String,G4int>::const_iterator itp;
131  for ( itp = localRun->fProcCounter.begin();
132  itp != localRun->fProcCounter.end(); ++itp ) {
133 
134  G4String procName = itp->first;
135  G4int localCount = itp->second;
136  if ( fProcCounter.find(procName) == fProcCounter.end()) {
137  fProcCounter[procName] = localCount;
138  }
139  else {
140  fProcCounter[procName] += localCount;
141  }
142  }
143 
144  //map: created particles count
145  std::map<G4String,ParticleData>::const_iterator itn;
146  for (itn = localRun->fParticleDataMap.begin();
147  itn != localRun->fParticleDataMap.end(); ++itn) {
148 
149  G4String name = itn->first;
150  const ParticleData& localData = itn->second;
151  if ( fParticleDataMap.find(name) == fParticleDataMap.end()) {
153  = ParticleData(localData.fCount,
154  localData.fEmean,
155  localData.fEmin,
156  localData.fEmax);
157  }
158  else {
159  ParticleData& data = fParticleDataMap[name];
160  data.fCount += localData.fCount;
161  data.fEmean += localData.fEmean;
162  G4double emin = localData.fEmin;
163  if (emin < data.fEmin) data.fEmin = emin;
164  G4double emax = localData.fEmax;
165  if (emax > data.fEmax) data.fEmax = emax;
166  }
167  }
168 
169  G4Run::Merge(run);
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173 
174 void Run::EndOfRun()
175 {
176  G4int prec = 5, wid = prec + 2;
177  G4int dfprec = G4cout.precision(prec);
178 
179  //run condition
180  //
181  G4Material* material = fDetector->GetMaterial();
182  G4double density = material->GetDensity();
183 
184  G4String Particle = fParticle->GetParticleName();
185  G4cout << "\n The run is " << numberOfEvent << " "<< Particle << " of "
186  << G4BestUnit(fEkin,"Energy") << " through "
187  << G4BestUnit(0.5*(fDetector->GetSize()),"Length") << " of "
188  << material->GetName() << " (density: "
189  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
190 
191  if (numberOfEvent == 0) { G4cout.precision(dfprec); return;}
192 
193  //frequency of processes
194  //
195  G4cout << "\n Process calls frequency :" << G4endl;
196  G4int survive = 0;
197  std::map<G4String,G4int>::iterator it;
198  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
199  G4String procName = it->first;
200  G4int count = it->second;
201  G4cout << "\t" << procName << "= " << count;
202  if (procName == "Transportation") survive = count;
203  }
204  G4cout << G4endl;
205 
206  if (survive > 0) {
207  G4cout << "\n Nb of incident particles surviving after "
208  << G4BestUnit(0.5*(fDetector->GetSize()),"Length") << " of "
209  << fDetector->GetMaterial()->GetName() << " : " << survive << G4endl;
210  }
211 
212  // total track length of incident neutron
213  //
214  G4cout << "\n Parcours of incident neutron:";
215 
216  G4double meanCollision1 = (G4double)fNbStep1/numberOfEvent;
217  G4double meanCollision2 = (G4double)fNbStep2/numberOfEvent;
218  G4double meanCollisTota = meanCollision1 + meanCollision2;
219 
220  G4cout << "\n nb of collisions E>1*eV= " << meanCollision1
221  << " E<1*eV= " << meanCollision2
222  << " total= " << meanCollisTota;
223 
224  G4double meanTrackLen1 = fTrackLen1/numberOfEvent;
225  G4double meanTrackLen2 = fTrackLen2/numberOfEvent;
226  G4double meanTrackLtot = meanTrackLen1 + meanTrackLen2;
227 
228  G4cout
229  << "\n track length E>1*eV= " << G4BestUnit(meanTrackLen1,"Length")
230  << " E<1*eV= " << G4BestUnit(meanTrackLen2, "Length")
231  << " total= " << G4BestUnit(meanTrackLtot, "Length");
232 
233  G4double meanTime1 = fTime1/numberOfEvent;
234  G4double meanTime2 = fTime2/numberOfEvent;
235  G4double meanTimeTo = meanTime1 + meanTime2;
236 
237  G4cout
238  << "\n time of flight E>1*eV= " << G4BestUnit(meanTime1,"Time")
239  << " E<1*eV= " << G4BestUnit(meanTime2, "Time")
240  << " total= " << G4BestUnit(meanTimeTo, "Time") << G4endl;
241 
242  //particles count
243  //
244  G4cout << "\n List of generated particles:" << G4endl;
245 
246  std::map<G4String,ParticleData>::iterator itn;
247  for (itn = fParticleDataMap.begin(); itn != fParticleDataMap.end(); itn++) {
248  G4String name = itn->first;
249  ParticleData data = itn->second;
250  G4int count = data.fCount;
251  G4double eMean = data.fEmean/count;
252  G4double eMin = data.fEmin;
253  G4double eMax = data.fEmax;
254 
255  G4cout << " " << std::setw(13) << name << ": " << std::setw(7) << count
256  << " Emean = " << std::setw(wid) << G4BestUnit(eMean, "Energy")
257  << "\t( " << G4BestUnit(eMin, "Energy")
258  << " --> " << G4BestUnit(eMax, "Energy")
259  << ")" << G4endl;
260  }
261 
262  //normalize histograms
266 
267  //remove all contents in fProcCounter, fCount
268  fProcCounter.clear();
269  fParticleDataMap.clear();
270 
271  //restore default format
272  G4cout.precision(dfprec);
273 }
274 
275 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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)
DetectorConstruction * fDetector
G4int numberOfEvent
Definition: G4Run.hh:59
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
static const G4double emax
void SumTrackLength(G4int, G4int, G4double, G4double, G4double, G4double)
void ParticleCount(G4String, G4double)
const XML_Char const XML_Char * data
Definition: expat.h:268
const G4String & GetName() const
Definition: G4Material.hh:179
double G4double
Definition: G4Types.hh:76
double energy
Definition: plottest35.C:25
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
virtual void Merge(const G4Run *)
std::map< G4String, G4int > fProcCounter
Definition: G4Run.hh:46
G4ParticleDefinition * fParticle
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
void CountProcesses(G4String procName)
G4GLOB_DLL std::ostream G4cout
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
G4int first(char) const
Simple detector construction with a box volume placed in a world.
G4double GetDensity() const
Definition: G4Material.hh:181