Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
examples/extended/electromagnetic/TestEm6/src/RunAction.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: RunAction.cc 83428 2014-08-21 15:46:01Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "G4SystemOfUnits.hh"
35 #include "G4PhysicalConstants.hh"
36 #include "G4EmCalculator.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4ParticleDefinition.hh"
39 #include "G4Positron.hh"
40 #include "G4AnnihiToMuPair.hh"
41 #include "G4eeToHadrons.hh"
42 #include "G4eeToHadronsModel.hh"
43 #include "G4eBremsstrahlung.hh"
44 
45 #include "RunAction.hh"
46 #include "G4Run.hh"
47 #include "G4RunManager.hh"
48 #include "G4MuonMinus.hh"
49 
50 #include "Randomize.hh"
51 
52 #include <sstream>
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
57  : G4UserRunAction(),fDetector(det),fProcCounter(0),fAnalysis(0),fMat(0)
58 {
59  fMinE = 40*GeV;
60  fMaxE = 10000*GeV;
61  fnBin = 10000;
62  fNtColId[0] = fNtColId[1] = fNtColId[2] = 0;
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
68 {}
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
72 void RunAction::BeginOfRunAction(const G4Run* aRun)
73 {
74  G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
75 
76  //get material
77  //
79  G4cout << "###RunAction::BeginOfRunAction Material:"
80  << fMat->GetName() << G4endl;
81 
83 
84  fAnalysis = G4AnalysisManager::Instance();
85 
86  // Open an output file
87  //
88  std::stringstream tmp;
89  tmp << "testem6_" << aRun->GetRunID();
90  G4String fileName = tmp.str();
91  fAnalysis->OpenFile(fileName);
92  fAnalysis->SetVerboseLevel(2);
93  G4String extension = fAnalysis->GetFileType();
94  fileName = fileName + "." + extension;
95 
96  // Creating histograms 1,2,3,4,5,6
97  //
98  fAnalysis->SetFirstHistoId(1);
99  fAnalysis->CreateH1("h1","1/(1+(theta+[g]+)**2)",100, 0 ,1.);
100  fAnalysis->CreateH1("h2","log10(theta+ [g]+)", 100,-3.,1.);
101  fAnalysis->CreateH1("h3","log10(theta- [g]-)", 100,-3.,1.);
102  fAnalysis->CreateH1("h4","log10(theta+ [g]+ -theta- [g]-)", 100,-3.,1.);
103  fAnalysis->CreateH1("h5","xPlus" ,100,0.,1.);
104  fAnalysis->CreateH1("h6","xMinus",100,0.,1.);
105 
106  //creating histogram 7,8,9,10,11 (CrossSectionPerAtom)
107  //
108  G4double minBin = std::log10(fMinE/GeV);
109  G4double maxBin = std::log10(fMaxE/GeV);
110  fAnalysis->CreateH1("h7","CrossSectionPerAtom of AnnihiToMuMu (microbarn)",
111  fnBin,minBin,maxBin);
112  fAnalysis->CreateH1("h8",
113  "CrossSectionPerAtom of AnnihiToTwoGamma (microbarn)",fnBin,minBin,maxBin);
114  fAnalysis->CreateH1("h9","CrossSectionPerAtom of AnnihiToHadrons (microbarn)",
115  fnBin,minBin,maxBin);
116  fAnalysis->CreateH1("h10",
117  "Theoretical CrossSectionPerAtom of AnnihiToTwoGamma (microbarn)",
118  fnBin,minBin,maxBin);
119  fAnalysis->CreateH1("h11",
120  "Theoretical CrossSectionPerAtom of AnnihiToMuMu (microbarn)",
121  fnBin,minBin,maxBin);
122 
123  //creating histogram 12,13,14,15,16(CrossSectionPerVolume)
124  //
125  fAnalysis->CreateH1("h12","CrossSectionPerVol of Bremsstraulung (1/mm) ",
126  fnBin,minBin,maxBin);
127  fAnalysis->CreateH1("h13","CrossSectionPerVol of Ionization (1/mm)",
128  fnBin,minBin,maxBin);
129  fAnalysis->CreateH1("h14","CrossSectionPerVol of AnnihiToMuMu (1/mm)",
130  fnBin,minBin,maxBin);
131  fAnalysis->CreateH1("h15","CrossSectionPerVol of AnnihiToTwoGamma (1/mm)",
132  fnBin,minBin,maxBin);
133  fAnalysis->CreateH1("h16","CrossSectionPerVol of AnnihiToHadrons (1/mm)",
134  fnBin,minBin,maxBin);
135 
136  //creating histogram 17 (R ratio)
137  fAnalysis->CreateH1("h17","R : eeToHadr/eeToMu",fnBin,minBin,maxBin);
138 
139  G4cout << "\n----> Histogram file is opened in " << fileName << G4endl;
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
144 void RunAction::CountProcesses(G4String procName)
145 {
146  //does the process already encounted ?
147  //
148  size_t nbProc = fProcCounter->size();
149  size_t i = 0;
150  while ((i<nbProc)&&((*fProcCounter)[i]->GetName()!=procName)) i++;
151  if (i == nbProc) fProcCounter->push_back( new OneProcessCount(procName));
152 
153  (*fProcCounter)[i]->Count();
154 }
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
157 
158 void RunAction::EndOfRunAction(const G4Run*)
159 {
160  // show Rndm status : not needed
161  //
162  //CLHEP::HepRandom::showEngineStatus();
163 
164  //total number of process calls
165  //
166  G4double countTot = 0.;
167  G4cout << "\n Number of process calls --->";
168  for (size_t i=0; i< fProcCounter->size();i++) {
169  G4String procName = (*fProcCounter)[i]->GetName();
170  if (procName != "Transportation") {
171  G4int count = (*fProcCounter)[i]->GetCounter();
172  G4cout << "\t" << procName << " : " << count;
173  countTot += count;
174  }
175  }
176 
177  //instance EmCalculator
178  //
179  G4EmCalculator emCal;
180  emCal.SetVerbose(0);
181 
182  //get positron
183  //
184  G4String positronName = "e+";
187 
188  //process name
189  //
190  G4String annihilName = "annihil";
191  G4String annihiToMuName = "AnnihiToMuPair";
192  G4String annihiToHadrName = "ee2hadr";
193  G4String BremName = "eBrem";
194  G4String IoniName = "eIoni";
195 
196 
197  //get AnnihiToMuPair
198  //
199  G4AnnihiToMuPair* annihiToMu =
200  reinterpret_cast<G4AnnihiToMuPair*>(emCal.FindProcess(positron,
201  annihiToMuName));
202 
203  G4double de, x, energy;
204  G4double crs_annihil, crs_annihiToMu, crs_annihiToHadr;
205  G4double crsVol_annihil, crsVol_annihiToMu, crsVol_annihiToHadr,
206  crsVol_Brem, crsVol_Ioni;
207  G4double crs_annihil_theory, crs_annihiToMu_theory, RR;
208  G4double X1, X2;
209 
210  G4double atomicZ, atomicA;
211 
212  G4double Mu; //for muon mass
213  G4double Me; //for electron mass
214  G4double Re; //for classical electron radius
215  G4double Ru; //for classical muon radius
216  G4double Eth; //for energy threshould to annihiToMu
217 
218  //parameters for ComputeCrossSection
219  //
220  G4double minBin = std::log10(fMinE/GeV);
221  G4double maxBin = std::log10(fMaxE/GeV);
222  de = (maxBin - minBin)/G4double(fnBin);
223  x = minBin - de*0.5;
224  atomicZ = 1.;
225  atomicA = 2.;
226 
227  //set parameters for theory
228  //
230  Mu = muon->GetPDGMass();
231  Me = electron_mass_c2;
233  Ru = Re*Me/Mu;
234  Eth = 2*Mu*Mu/Me-Me;
235 
236  //Compute CrossSections and Fill histgrams
237  //
238  for(int i=0;i<fnBin;i++){
239  x += de;
240  energy=std::pow(10,x)*GeV;
241 
242  //CrossSectionPerAtom
243  //
244  crs_annihiToMu = annihiToMu->ComputeCrossSectionPerAtom(energy,atomicZ);
245  fAnalysis->FillH1(7,x,crs_annihiToMu/microbarn);
246 
247  crs_annihil =
248  emCal.ComputeCrossSectionPerAtom(energy,positron,annihilName,
249  atomicZ,atomicA);
250  fAnalysis->FillH1(8,x,crs_annihil/microbarn);
251 
252  crs_annihiToHadr =
253  emCal.ComputeCrossSectionPerAtom(energy,positron,annihiToHadrName,
254  atomicZ,atomicA);
255  fAnalysis->FillH1(9,x,crs_annihiToHadr/microbarn);
256 
257  //CrossSectionPerVolume
258  //
259  crsVol_Brem =
260  emCal.ComputeCrossSectionPerVolume(energy,positron,BremName,fMat,100*keV);
261  fAnalysis->FillH1(12,x,crsVol_Brem*mm);
262 
263  crsVol_Ioni =
264  emCal.ComputeCrossSectionPerVolume(energy,positron,IoniName,fMat,100*keV);
265  fAnalysis->FillH1(13,x,crsVol_Ioni*mm);
266 
267  crsVol_annihiToMu = annihiToMu->CrossSectionPerVolume(energy,fMat);
268  fAnalysis->FillH1(14,x,crsVol_annihiToMu*mm);
269 
270  crsVol_annihil =
271  emCal.ComputeCrossSectionPerVolume(energy,positron,annihilName,fMat);
272  fAnalysis->FillH1(15,x,crsVol_annihil*mm);
273 
274  crsVol_annihiToHadr =
275  emCal.ComputeCrossSectionPerVolume(energy,positron,annihiToHadrName,fMat);
276  fAnalysis->FillH1(16,x,crsVol_annihiToHadr*mm);
277 
278  //R ratio
279  //
280  RR = 0.0;
281  if(crsVol_annihiToMu != 0) RR = crsVol_annihiToHadr/crsVol_annihiToMu;
282  fAnalysis->FillH1(17,x,RR);
283 
284  //Theoretical calculation
285  //
286  X1 = energy/Me;
287  if(X1>1 && i%1000==0){
288  crs_annihil_theory = atomicZ*pi*Re*Re*
289  ( (X1*X1+4*X1+1)*G4Log(X1+std::sqrt(X1*X1-1))/(X1*X1-1)
290  -(X1+3)/std::sqrt(X1*X1-1) )/(X1+1);
291  fAnalysis->FillH1(10,x,crs_annihil_theory/microbarn);
292  }
293 
294  X2 = Eth/energy;
295  if(X2<1 && i%1000==0){
296  crs_annihiToMu_theory = atomicZ*pi*Ru*Ru/3*X2*(1+X2/2)*std::sqrt(1-X2);
297  fAnalysis->FillH1(11,x,crs_annihiToMu_theory/microbarn);
298  }
299 
300  //if(i%1000==0)G4cout <<"###energy:" << energy << "/crs_ToMuMu:"
301  // << crs_annihiToMu << "/crs_ToTwoGamma:"<< crs_annihil
302  // <<"/crs_ToToHadr:"<<crs_annihiToHadr<< G4endl;
303  }
304 
305  fAnalysis->Write();
306  fAnalysis->CloseFile();
307 
308  delete fAnalysis;
309  G4cout << G4endl;
310 }
311 
312 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
Float_t x
Definition: compare.C:6
static G4ParticleTable * GetParticleTable()
static constexpr double keV
Definition: G4SIunits.hh:216
static constexpr double mm
Definition: G4SIunits.hh:115
std::vector< OneProcessCount * > ProcessesCount
#define G4endl
Definition: G4ios.hh:61
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
Float_t tmp
G4double GetPDGMass() const
static constexpr double classic_electr_radius
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Double_t X2
const G4String & GetName() const
Definition: G4Material.hh:179
double G4double
Definition: G4Types.hh:76
Double_t X1
static constexpr double electron_mass_c2
G4double CrossSectionPerVolume(G4double PositronEnergy, const G4Material *)
G4int GetRunID() const
Definition: G4Run.hh:76
double energy
Definition: plottest35.C:25
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
Definition: G4Run.hh:46
G4double ComputeCrossSectionPerAtom(G4double PositronEnergy, G4double AtomicZ)
static constexpr double microbarn
Definition: G4SIunits.hh:107
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double pi
Definition: G4SIunits.hh:75
void SetVerbose(G4int val)
static constexpr double GeV
Definition: G4SIunits.hh:217
Simple detector construction with a box volume placed in a world.