Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/electromagnetic/TestEm11/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 
37 #include "EventAction.hh"
38 #include "HistoManager.hh"
39 #include "PrimaryGeneratorAction.hh"
40 
41 #include "G4Material.hh"
42 #include "G4Event.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4UnitsTable.hh"
45 #include <iomanip>
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
50 : G4Run(),
51  fDetector(detector),
52  fParticle(0), fEkin(0.),
53  fTrackLen(0.), fTrackLen2(0.),
54  fProjRange(0.), fProjRange2(0.),
55  fNbOfSteps(0), fNbOfSteps2(0),
56  fStepSize(0.), fStepSize2(0.)
57 {
58  for (G4int i=0; i<3; ++i) { fStatus[i] = 0; fTotEdep[i] = 0.; }
59  fTotEdep[1] = joule;
60  for (G4int i=0; i<kMaxAbsor; ++i) {
61  fEdeposit[i] = 0.; fEmin[i] = joule; fEmax[i] = 0.;
62  fCsdaRange[i] = 0.; fXfrontNorm[i] = 0.;
63  }
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
68 Run::~Run()
69 { }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  fParticle = particle;
76  fEkin = energy;
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82 {
83  if (e > 0.) {
84  fEdeposit[i] += e;
85  if (e < fEmin[i]) fEmin[i] = e;
86  if (e > fEmax[i]) fEmax[i] = e;
87  }
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {
94  if (e > 0.) {
95  fTotEdep[0] += e;
96  if (e < fTotEdep[1]) fTotEdep[1] = e;
97  if (e > fTotEdep[2]) fTotEdep[2] = e;
98  }
99 }
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101 
103 {
104  fTrackLen += t;
105  fTrackLen2 += t*t;
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109 
111 {
112  fProjRange += x;
113  fProjRange2 += x*x;
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117 
119 {
120  fNbOfSteps += nb;
121  fNbOfSteps2 += nb*nb;
122  fStepSize += st ;
123  fStepSize2 += st*st;
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
129 {
130  fStatus[i]++ ;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136 {
137  fCsdaRange[i] = value;
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
143 {
144  fXfrontNorm[i] = value;
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150 {
151  return fCsdaRange[i];
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
157 {
158  return fXfrontNorm[i];
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
163 void Run::Merge(const G4Run* run)
164 {
165  const Run* localRun = static_cast<const Run*>(run);
166 
167  // pass information about primary particle
168  fParticle = localRun->fParticle;
169  fEkin = localRun->fEkin;
170 
171  // accumulate sums
172  fTrackLen += localRun->fTrackLen;
173  fTrackLen2 += localRun->fTrackLen2;
174  fProjRange += localRun->fProjRange;
175  fProjRange2 += localRun->fProjRange2;
176  fNbOfSteps += localRun->fNbOfSteps ;
177  fNbOfSteps2 += localRun->fNbOfSteps2;
178  fStepSize += localRun->fStepSize;
179  fStepSize2 += localRun->fStepSize2;
180 
181  G4int nbOfAbsor = fDetector->GetNbOfAbsor();
182  for (G4int i=1; i<=nbOfAbsor; ++i) {
183  fEdeposit[i] += localRun->fEdeposit[i];
184  fCsdaRange[i] = localRun->fCsdaRange[i];
185  fXfrontNorm[i] = localRun->fXfrontNorm[i];
186  // min, max
187  G4double min,max;
188  min = localRun->fEmin[i]; max = localRun->fEmax[i];
189  if (fEmin[i] > min) fEmin[i] = min;
190  if (fEmax[i] < max) fEmax[i] = max;
191  }
192 
193  for (G4int i=0; i<3; ++i) fStatus[i] += localRun->fStatus[i];
194 
195  // total Edep
196  fTotEdep[0] += localRun->fTotEdep[0];
197  G4double min,max;
198  min = localRun->fTotEdep[1]; max = localRun->fTotEdep[2];
199 if (fTotEdep[1] > min) fTotEdep[1] = min;
200 if (fTotEdep[2] < max) fTotEdep[2] = max;
201 
202  G4Run::Merge(run);
203 }
204 
205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206 
207 void Run::EndOfRun()
208 {
209  std::ios::fmtflags mode = G4cout.flags();
210  G4cout.setf(std::ios::fixed,std::ios::floatfield);
211  G4int prec = G4cout.precision(2);
212 
213  //run conditions
214  //
215  G4String partName = fParticle->GetParticleName();
216  G4int nbOfAbsor = fDetector->GetNbOfAbsor();
217 
218  G4cout << "\n ======================== run summary =====================\n";
219  G4cout
220  << "\n The run is " << numberOfEvent << " "<< partName << " of "
221  << G4BestUnit(fEkin,"Energy")
222  << " through " << nbOfAbsor << " absorbers: \n";
223  for (G4int i=1; i<= nbOfAbsor; i++) {
224  G4Material* material = fDetector->GetAbsorMaterial(i);
225  G4double thickness = fDetector->GetAbsorThickness(i);
226  G4double density = material->GetDensity();
227  G4cout << std::setw(5) << i
228  << std::setw(10) << G4BestUnit(thickness,"Length") << " of "
229  << material->GetName() << " (density: "
230  << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
231  }
232 
233  if (numberOfEvent == 0) {
234  G4cout.setf(mode,std::ios::floatfield);
235  G4cout.precision(prec);
236  return;
237  }
238 
239  G4cout.precision(3);
240  G4double rms (0);
241 
242  for (G4int i=1; i<= nbOfAbsor; i++) {
243  fEdeposit[i] /= numberOfEvent;
244 
245  G4cout
246  << "\n Edep in absorber " << i << " = "
247  << G4BestUnit(fEdeposit[i],"Energy")
248  << "\t(" << G4BestUnit(fEmin[i], "Energy")
249  << "-->" << G4BestUnit(fEmax[i], "Energy")
250  << ")";
251  }
252  G4cout << G4endl;
253 
254  if (nbOfAbsor > 1) {
255  fTotEdep[0] /= numberOfEvent;
256  G4cout
257  << "\n Edep in all absorbers = " << G4BestUnit(fTotEdep[0],"Energy")
258  << "\t(" << G4BestUnit(fTotEdep[1], "Energy")
259  << "-->" << G4BestUnit(fTotEdep[2], "Energy")
260  << ")" << G4endl;
261  }
262 
263  //compute track length of primary track
264  //
266  rms = fTrackLen2 - fTrackLen*fTrackLen;
267  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
268 
269  G4cout.precision(3);
270  G4cout
271  << "\n Track length of primary track = " << G4BestUnit(fTrackLen,"Length")
272  << " +- " << G4BestUnit( rms,"Length");
273 
274  //compare with csda range
275  //
276  G4int NbOfAbsor = fDetector->GetNbOfAbsor();
277  if (NbOfAbsor == 1) {
278  G4cout
279  << "\n Range from EmCalculator = " << G4BestUnit(fCsdaRange[1],"Length")
280  << " (from full dE/dx)" << G4endl;
281  }
282 
283  //compute projected range of primary track
284  //
287  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
288 
289  G4cout
290  << "\n Projected range = " << G4BestUnit(fProjRange,"Length")
291  << " +- " << G4BestUnit( rms,"Length")
292  << G4endl;
293 
294  //nb of steps and step size of primary track
295  //
296  G4double dNofEvents = double(numberOfEvent);
297  G4double fNbSteps = fNbOfSteps/dNofEvents,
298  fNbSteps2 = fNbOfSteps2/dNofEvents;
299  rms = fNbSteps2 - fNbSteps*fNbSteps;
300  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
301 
302  G4cout.precision(2);
303  G4cout << "\n Nb of steps of primary track = " << fNbSteps << " +- " << rms;
304 
306  rms = fStepSize2 - fStepSize*fStepSize;
307  if (rms>0.) rms = std::sqrt(rms); else rms = 0.;
308 
309  G4cout.precision(3);
310  G4cout
311  << "\t Step size= " << G4BestUnit(fStepSize,"Length")
312  << " +- " << G4BestUnit( rms,"Length")
313  << G4endl;
314 
315  //transmission coefficients
316  //
317  G4double absorbed = 100.*fStatus[0]/dNofEvents;
318  G4double transmit = 100.*fStatus[1]/dNofEvents;
319  G4double reflected = 100.*fStatus[2]/dNofEvents;
320 
321  G4cout.precision(2);
322  G4cout
323  << "\n absorbed = " << absorbed << " %"
324  << " transmit = " << transmit << " %"
325  << " reflected = " << reflected << " %" << G4endl;
326 
327  // normalize histograms of longitudinal energy profile
328  //
329  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
330  G4int ih = 1;
331  G4double binWidth = analysisManager->GetH1Width(ih)
332  *analysisManager->GetH1Unit(ih);
333  G4double fac = (1./(numberOfEvent*binWidth))*(mm/MeV);
334  analysisManager->ScaleH1(ih,fac);
335 
336  ih = 8;
337  binWidth = analysisManager->GetH1Width(ih);
338  fac = (1./(numberOfEvent*binWidth))*(g/(MeV*cm2));
339  analysisManager->ScaleH1(ih,fac);
340 
341  // reset default formats
342  G4cout.setf(mode,std::ios::floatfield);
343  G4cout.precision(prec);
344 }
345 
346 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static const double prec
Definition: RanecuEngine.cc:58
void SetPrimary(G4ParticleDefinition *particle, G4double energy)
DetectorConstruction * fDetector
static constexpr double MeV
Definition: G4SIunits.hh:214
G4int numberOfEvent
Definition: G4Run.hh:59
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
const G4String & GetName() const
Definition: G4Material.hh:179
static constexpr double g
Definition: G4SIunits.hh:183
double G4double
Definition: G4Types.hh:76
static constexpr double cm2
Definition: G4SIunits.hh:120
if(nlines<=0)
const XML_Char int const XML_Char * value
Definition: expat.h:331
void SetXfrontNorm(G4int i, G4double value)
double energy
Definition: plottest35.C:25
virtual void Merge(const G4Run *)
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
Definition: G4Run.hh:46
G4ParticleDefinition * fParticle
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
static const G4double fac
static constexpr double joule
Definition: G4SIunits.hh:204
G4GLOB_DLL std::ostream G4cout
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
void SetCsdaRange(G4int i, G4double value)
void AddStepSize(G4int nb, G4double st)
Simple detector construction with a box volume placed in a world.
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetDensity() const
Definition: G4Material.hh:181