Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/electromagnetic/TestEm1/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 "G4ProductionCutsTable.hh"
42 
43 #include <iomanip>
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 : G4Run(),
49  fDetector(det),
50  fParticle(0), fEkin(0.),
51  fNbOfTraks0(0), fNbOfTraks1(0),
52  fNbOfSteps0(0), fNbOfSteps1(0),
53  fEdep(0),
54  fTrueRange(0.), fTrueRange2(0.),
55  fProjRange(0.), fProjRange2(0.),
56  fTransvDev(0.), fTransvDev2(0.)
57 { }
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 
61 Run::~Run()
62 { }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67 {
68  fParticle = particle;
69  fEkin = energy;
70 }
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  std::map<G4String,G4int>::iterator it = fProcCounter.find(procName);
76  if ( it == fProcCounter.end()) {
77  fProcCounter[procName] = 1;
78  }
79  else {
80  fProcCounter[procName]++;
81  }
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85 
86 void Run::Merge(const G4Run* run)
87 {
88  const Run* localRun = static_cast<const Run*>(run);
89 
90  // pass information about primary particle
91  fParticle = localRun->fParticle;
92  fEkin = localRun->fEkin;
93 
94  // accumulate sums
95  //
96  fNbOfTraks0 += localRun->fNbOfTraks0;
97  fNbOfTraks1 += localRun->fNbOfTraks1;
98  fNbOfSteps0 += localRun->fNbOfSteps0;
99  fNbOfSteps1 += localRun->fNbOfSteps1;
100  fEdep += localRun->fEdep;
101  fTrueRange += localRun->fTrueRange;
102  fTrueRange2 += localRun->fTrueRange2;
103  fProjRange += localRun->fProjRange;
104  fProjRange2 += localRun->fProjRange2;
105  fTransvDev += localRun->fTransvDev;
106  fTransvDev2 += localRun->fTransvDev2;
107 
108  //map: processes count
109  std::map<G4String,G4int>::const_iterator it;
110  for (it = localRun->fProcCounter.begin();
111  it !=localRun->fProcCounter.end(); ++it) {
112 
113  G4String procName = it->first;
114  G4int localCount = it->second;
115  if ( fProcCounter.find(procName) == fProcCounter.end()) {
116  fProcCounter[procName] = localCount;
117  }
118  else {
119  fProcCounter[procName] += localCount;
120  }
121  }
122 
123  G4Run::Merge(run);
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
128 void Run::EndOfRun()
129 {
130  G4int prec = 5, wid = prec + 2;
131  G4int dfprec = G4cout.precision(prec);
132 
133  //run condition
134  //
135  G4String partName = fParticle->GetParticleName();
136  G4Material* material = fDetector->GetMaterial();
137  G4double density = material->GetDensity();
138 
139  G4cout << "\n ======================== run summary ======================\n";
140  G4cout << "\n The run is: " << numberOfEvent << " " << partName << " of "
141  << G4BestUnit(fEkin,"Energy") << " through "
142  << G4BestUnit(fDetector->GetSize(),"Length") << " of "
143  << material->GetName() << " (density: "
144  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
145 
146  if (numberOfEvent == 0) { G4cout.precision(dfprec); return;}
147 
148  G4double dNbOfEvents = double(numberOfEvent);
149  G4cout << "\n total energy deposit: "
150  << G4BestUnit(fEdep/dNbOfEvents, "Energy") << G4endl;
151 
152  //nb of tracks and steps per event
153  //
154  G4cout << "\n nb tracks/event"
155  << " neutral: " << std::setw(wid) << fNbOfTraks0/dNbOfEvents
156  << " charged: " << std::setw(wid) << fNbOfTraks1/dNbOfEvents
157  << "\n nb steps/event"
158  << " neutral: " << std::setw(wid) << fNbOfSteps0/dNbOfEvents
159  << " charged: " << std::setw(wid) << fNbOfSteps1/dNbOfEvents
160  << G4endl;
161 
162  //frequency of processes
163  //
164  G4cout << "\n Process calls frequency :" << G4endl;
165  G4int index = 0;
166  std::map<G4String,G4int>::iterator it;
167  for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
168  G4String procName = it->first;
169  G4int count = it->second;
170  G4String space = " "; if (++index%3 == 0) space = "\n";
171  G4cout << " " << std::setw(20) << procName << "="<< std::setw(7) << count
172  << space;
173  }
174  G4cout << G4endl;
175 
176  //compute true and projected ranges, and transverse dispersion
177  //
180  if (trueRms>0.) trueRms = std::sqrt(trueRms); else trueRms = 0.;
181 
184  if (projRms>0.) projRms = std::sqrt(projRms); else projRms = 0.;
185 
188  if (trvsRms>0.) trvsRms = std::sqrt(trvsRms); else trvsRms = 0.;
189 
190  //compare true range with csda range from PhysicsTables
191  //
192  G4EmCalculator emCalculator;
193  G4double rangeTable = 0.;
194  if (fParticle->GetPDGCharge() != 0.)
195  rangeTable = emCalculator.GetCSDARange(fEkin,fParticle,material);
196 
197  G4cout << "\n---------------------------------------------------------\n";
198  G4cout << " Primary particle : " ;
199  G4cout << "\n true Range = " << G4BestUnit(fTrueRange,"Length")
200  << " rms = " << G4BestUnit(trueRms, "Length");
201 
202  G4cout << "\n proj Range = " << G4BestUnit(fProjRange,"Length")
203  << " rms = " << G4BestUnit(projRms, "Length");
204 
205  G4cout << "\n proj/true = " << fProjRange/fTrueRange;
206 
207  G4cout << "\n transverse dispersion at end = "
208  << G4BestUnit(trvsRms,"Length");
209 
210  G4cout << "\n mass true Range from simulation = "
211  << G4BestUnit(fTrueRange*density, "Mass/Surface")
212  << "\n from PhysicsTable (csda range) = "
213  << G4BestUnit(rangeTable*density, "Mass/Surface");
214  G4cout << "\n---------------------------------------------------------\n";
215  G4cout << G4endl;
216 
217  // remove all contents in fProcCounter
218  fProcCounter.clear();
219 
220  //restore default format
221  G4cout.precision(dfprec);
222 }
223 
224 //....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
G4double GetPDGCharge() const
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
int G4int
Definition: G4Types.hh:78
void CountProcesses(G4String procName)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
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