Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/electromagnetic/TestEm13/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 
38 #include "G4UnitsTable.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "G4EmCalculator.hh"
41 #include "G4Gamma.hh"
42 
43 #include <iomanip>
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 : G4Run(),
49  fDetector(det),
50  fParticle(0), fEkin(0.)
51 { }
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
55 Run::~Run()
56 { }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61 {
62  fParticle = particle;
63  fEkin = energy;
64 }
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
67 void Run::CountProcesses(G4String procName)
68 {
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 
80 void Run::Merge(const G4Run* run)
81 {
82  const Run* localRun = static_cast<const Run*>(run);
83 
84  // pass information about primary particle
85  fParticle = localRun->fParticle;
86  fEkin = localRun->fEkin;
87 
88  //map: processes count
89  std::map<G4String,G4int>::const_iterator it;
90  for (it = localRun->fProcCounter.begin();
91  it !=localRun->fProcCounter.end(); ++it) {
92 
93  G4String procName = it->first;
94  G4int localCount = it->second;
95  if ( fProcCounter.find(procName) == fProcCounter.end()) {
96  fProcCounter[procName] = localCount;
97  }
98  else {
99  fProcCounter[procName] += localCount;
100  }
101  }
102 
103  G4Run::Merge(run);
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
108 void Run::EndOfRun()
109 {
110  G4int prec = 5;
111  G4int dfprec = G4cout.precision(prec);
112 
113  //run condition
114  //
115  G4String partName = fParticle->GetParticleName();
116  G4Material* material = fDetector->GetMaterial();
117  G4double density = material->GetDensity();
118  G4double tickness = fDetector->GetSize();
119 
120  G4cout << "\n ======================== run summary ======================\n";
121  G4cout << "\n The run is: " << numberOfEvent << " " << partName << " of "
122  << G4BestUnit(fEkin,"Energy") << " through "
123  << G4BestUnit(tickness,"Length") << " of "
124  << material->GetName() << " (density: "
125  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
126 
127  //frequency of processes
128  G4int totalCount = 0;
129  G4int survive = 0;
130  G4cout << "\n Process calls frequency --->";
131  std::map<G4String,G4int>::iterator it;
132  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
133  G4String procName = it->first;
134  G4int count = it->second;
135  totalCount += count;
136  G4cout << "\t" << procName << " = " << count;
137  if (procName == "Transportation") survive = count;
138  }
139  G4cout << G4endl;
140 
141  if (totalCount == 0) { G4cout.precision(dfprec); return;};
142  G4double ratio = double(survive)/totalCount;
143 
144  G4cout << "\n Nb of incident particles unaltered after "
145  << G4BestUnit(tickness,"Length") << " of "
146  << material->GetName() << " : " << survive
147  << " over " << totalCount << " incident particles."
148  << " Ratio = " << 100*ratio << " %" << G4endl;
149 
150  if (ratio == 0.) return;
151 
152  //compute cross section and related quantities
153  //
154  G4double CrossSection = - std::log(ratio)/tickness;
155  G4double massicCS = CrossSection/density;
156 
157  G4cout << " ---> CrossSection per volume:\t" << CrossSection*cm << " cm^-1 "
158  << "\tCrossSection per mass: " << G4BestUnit(massicCS, "Surface/Mass")
159  << G4endl;
160 
161  //check cross section from G4EmCalculator
162  //
163  G4cout << "\n Verification from G4EmCalculator: \n";
164  G4EmCalculator emCalculator;
165  G4double sumc = 0.0;
166  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
167  G4String procName = it->first;
168  G4double massSigma =
170  procName,material)/density;
171  if (fParticle == G4Gamma::Gamma())
172  massSigma =
174  procName,material)/density;
175  sumc += massSigma;
176  if (procName != "Transportation")
177  G4cout << "\t" << procName << "= "
178  << G4BestUnit(massSigma, "Surface/Mass");
179  }
180  G4cout << "\ttotal= "
181  << G4BestUnit(sumc, "Surface/Mass") << G4endl;
182 
183  //expected ratio of transmitted particles
184  G4double Ratio = std::exp(-sumc*density*tickness);
185  G4cout << "\tExpected ratio of transmitted particles= "
186  << 100*Ratio << " %" << G4endl;
187 
188  // remove all contents in fProcCounter
189  fProcCounter.clear();
190 
191  //restore default format
192  G4cout.precision(dfprec);
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static const double prec
Definition: RanecuEngine.cc:58
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 G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4String & GetName() const
Definition: G4Material.hh:179
double G4double
Definition: G4Types.hh:76
double energy
Definition: plottest35.C:25
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
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
int G4int
Definition: G4Types.hh:78
void CountProcesses(G4String procName)
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
static constexpr double cm
Definition: G4SIunits.hh:119
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