Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
examples/extended/electromagnetic/TestEm0/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 107324 2017-11-08 16:35:06Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "RunAction.hh"
35 #include "DetectorConstruction.hh"
36 #include "PrimaryGeneratorAction.hh"
37 
38 #include "G4Run.hh"
39 #include "G4ProcessManager.hh"
40 #include "G4LossTableManager.hh"
41 #include "G4UnitsTable.hh"
42 #include "G4EmCalculator.hh"
43 #include "G4PhysicalConstants.hh"
44 #include "G4SystemOfUnits.hh"
45 #include "G4Electron.hh"
46 #include "G4Positron.hh"
47 
48 #include <vector>
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
53 :G4UserRunAction(),fDetector(det), fPrimary(kin)
54 { }
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
59 { }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
64 {
65  //set precision for printing
66  G4int prec = G4cout.precision(6);
67 
68  //instanciate EmCalculator
69  G4EmCalculator emCal;
70  // emCal.SetVerbose(2);
71 
72  // get particle
75  G4String partName = particle->GetParticleName();
76  G4double charge = particle->GetPDGCharge();
78 
79  // get material
80  const G4Material* material = fDetector->GetMaterial();
81  G4String matName = material->GetName();
82  G4double density = material->GetDensity();
83  G4double radl = material->GetRadlen();
84 
85  G4cout << "\n " << partName << " ("
86  << G4BestUnit(energy,"Energy") << ") in "
87  << material->GetName() << " (density: "
88  << G4BestUnit(density,"Volumic Mass") << "; radiation length: "
89  << G4BestUnit(radl, "Length") << ")" << G4endl;
90 
91  // get cuts
92  GetCuts();
93  if (charge != 0.) {
94  G4cout << "\n Range cuts : \t gamma "
95  << std::setw(8) << G4BestUnit(fRangeCut[0],"Length")
96  << "\t e- " << std::setw(8) << G4BestUnit(fRangeCut[1],"Length");
97  G4cout << "\n Energy cuts : \t gamma "
98  << std::setw(8) << G4BestUnit(fEnergyCut[0],"Energy")
99  << "\t e- " << std::setw(8) << G4BestUnit(fEnergyCut[1],"Energy")
100  << G4endl;
101  }
102 
103  // max energy transfert
104  if (charge != 0.) {
105  G4double Mass_c2 = particle->GetPDGMass();
106  G4double moverM = electron_mass_c2/Mass_c2;
107  G4double gamM1 = energy/Mass_c2, gam = gamM1 + 1., gamP1 = gam + 1.;
108  G4double Tmax = energy;
109  if(particle == G4Electron::Electron()) {
110  Tmax *= 0.5;
111  } else if(particle != G4Positron::Positron()) {
112  Tmax = (2*electron_mass_c2*gamM1*gamP1)/(1.+2*gam*moverM+moverM*moverM);
113  }
114  G4double range = emCal.GetCSDARange(Tmax,G4Electron::Electron(),material);
115 
116  G4cout << "\n Max_energy _transferable : " << G4BestUnit(Tmax,"Energy")
117  << " (" << G4BestUnit(range,"Length") << ")" << G4endl;
118  }
119 
120  // get processList and extract EM processes (but not MultipleScattering)
121  G4ProcessVector* plist = particle->GetProcessManager()->GetProcessList();
122  G4String procName;
123  G4double cut;
124  std::vector<G4String> emName;
125  std::vector<G4double> enerCut;
126  size_t length = plist->size();
127  for (size_t j=0; j<length; j++) {
128  procName = (*plist)[j]->GetProcessName();
129  cut = fEnergyCut[1];
130  if ((procName == "eBrem")||(procName == "muBrems")) cut = fEnergyCut[0];
131  if (((*plist)[j]->GetProcessType() == fElectromagnetic) &&
132  (procName != "msc")) {
133  emName.push_back(procName);
134  enerCut.push_back(cut);
135  }
136  }
137 
138  // write html documentation, if requested
139  char* htmlDocName = getenv("G4PhysListName"); // file name
140  char* htmlDocDir = getenv("G4PhysListDocDir"); // directory
141  if (htmlDocName && htmlDocDir) {
143  }
144 
145  // print list of processes
146  G4cout << "\n processes : ";
147  for (size_t j=0; j<emName.size();j++)
148  G4cout << "\t" << std::setw(13) << emName[j] << "\t";
149  G4cout << "\t" << std::setw(13) <<"total";
150 
151  //compute cross section per atom (only for single material)
152  if (material->GetNumberOfElements() == 1) {
153  G4double Z = material->GetZ();
154  G4double A = material->GetA();
155 
156  std::vector<G4double> sigma0;
157  G4double sig, sigtot = 0.;
158 
159  for (size_t j=0; j<emName.size();j++) {
160  sig = emCal.ComputeCrossSectionPerAtom
161  (energy,particle,emName[j],Z,A,enerCut[j]);
162  sigtot += sig;
163  sigma0.push_back(sig);
164  }
165  sigma0.push_back(sigtot);
166 
167  G4cout << "\n \n cross section per atom : ";
168  for (size_t j=0; j<sigma0.size();j++) {
169  G4cout << "\t" << std::setw(13) << G4BestUnit(sigma0[j], "Surface");
170  }
171  G4cout << G4endl;
172  }
173 
174  //get cross section per volume
175  std::vector<G4double> sigma0;
176  std::vector<G4double> sigma1;
177  std::vector<G4double> sigma2;
178  G4double Sig, SigtotComp = 0., Sigtot = 0.;
179 
180  for (size_t j=0; j<emName.size();j++) {
181  Sig = emCal.ComputeCrossSectionPerVolume
182  (energy,particle,emName[j],material,enerCut[j]);
183  SigtotComp += Sig;
184  sigma0.push_back(Sig);
185  Sig = emCal.GetCrossSectionPerVolume(energy,particle,emName[j],material);
186  Sigtot += Sig;
187  sigma1.push_back(Sig);
188  sigma2.push_back(Sig/density);
189  }
190  sigma0.push_back(SigtotComp);
191  sigma1.push_back(Sigtot);
192  sigma2.push_back(Sigtot/density);
193 
194  //print cross sections
195  G4cout << "\n \n compCrossSectionPerVolume : ";
196  for (size_t j=0; j<sigma0.size();j++) {
197  G4cout << "\t" << std::setw(13) << sigma0[j]*cm << " cm^-1";
198  }
199  G4cout << "\n cross section per volume : ";
200  for (size_t j=0; j<sigma1.size();j++) {
201  G4cout << "\t" << std::setw(13) << sigma1[j]*cm << " cm^-1";
202  }
203 
204  G4cout << "\n cross section per mass : ";
205  for (size_t j=0; j<sigma2.size();j++) {
206  G4cout << "\t" << std::setw(13)
207  << G4BestUnit(sigma2[j], "Surface/Mass");
208  }
209 
210  //print mean free path
211 
213 
214  G4cout << "\n \n mean free path : ";
215  for (size_t j=0; j<sigma1.size();j++) {
216  lambda = DBL_MAX;
217  if (sigma1[j] > 0.) lambda = 1/sigma1[j];
218  G4cout << "\t" << std::setw(13) << G4BestUnit( lambda, "Length");
219  }
220 
221  //mean free path (g/cm2)
222  G4cout << "\n (g/cm2) : ";
223  for (size_t j=0; j<sigma2.size();j++) {
224  lambda = DBL_MAX;
225  if (sigma2[j] > 0.) lambda = 1/sigma2[j];
226  G4cout << "\t" << std::setw(13) << G4BestUnit( lambda, "Mass/Surface");
227  }
228  G4cout << G4endl;
229 
230  if (charge == 0.) {
231  G4cout.precision(prec);
232  G4cout << "\n-----------------------------------------------------------\n"
233  << G4endl;
234  return;
235  }
236 
237  //get stopping power
238  std::vector<G4double> dedx1;
239  std::vector<G4double> dedx2;
240  G4double dedx, dedxtot = 0.;
241  size_t nproc = emName.size();
242 
243  for (size_t j=0; j<nproc; j++) {
244  dedx = emCal.ComputeDEDX(energy,particle,emName[j],material,enerCut[j]);
245  dedxtot += dedx;
246  dedx1.push_back(dedx);
247  dedx2.push_back(dedx/density);
248  }
249  dedx1.push_back(dedxtot);
250  dedx2.push_back(dedxtot/density);
251 
252  //print stopping power
253  G4cout << "\n \n restricted dE/dx : ";
254  for (size_t j=0; j<=nproc; j++) {
255  G4cout << "\t" << std::setw(13)
256  << G4BestUnit(dedx1[j],"Energy/Length");
257  }
258 
259  G4cout << "\n (MeV/g/cm2) : ";
260  for (size_t j=0; j<=nproc; j++) {
261  G4cout << "\t" << std::setw(13)
262  << G4BestUnit(dedx2[j],"Energy*Surface/Mass");
263  }
264  dedxtot = 0.;
265 
266  for (size_t j=0; j<nproc; j++) {
267  dedx = emCal.ComputeDEDX(energy,particle,emName[j],material,energy);
268  dedxtot += dedx;
269  dedx1[j] = dedx;
270  dedx2[j] = dedx/density;
271  }
272  dedx1[nproc] = dedxtot;
273  dedx2[nproc] = dedxtot/density;
274 
275  //print stopping power
276  G4cout << "\n \n unrestricted dE/dx : ";
277  for (size_t j=0; j<=nproc; j++) {
278  G4cout << "\t" << std::setw(13) << G4BestUnit(dedx1[j],"Energy/Length");
279  }
280 
281  G4cout << "\n (MeV/g/cm2) : ";
282  for (size_t j=0; j<=nproc; j++) {
283  G4cout << "\t" << std::setw(13)
284  << G4BestUnit(dedx2[j],"Energy*Surface/Mass");
285  }
286 
287  //get range from restricted dedx
288  G4double range1 = emCal.GetRangeFromRestricteDEDX(energy,particle,material);
289  G4double range2 = range1*density;
290 
291  //print range
292  G4cout << "\n \n range from restrict dE/dx: "
293  << "\t" << std::setw(8) << G4BestUnit(range1,"Length")
294  << " (" << std::setw(8) << G4BestUnit(range2,"Mass/Surface") << ")";
295 
296  //get range from full dedx
298  if(energy < EmaxTable) {
299  G4double Range1 = emCal.GetCSDARange(energy,particle,material);
300  G4double Range2 = Range1*density;
301 
302  G4cout << "\n range from full dE/dx : "
303  << "\t" << std::setw(8) << G4BestUnit(Range1,"Length")
304  << " (" << std::setw(8) << G4BestUnit(Range2,"Mass/Surface") << ")";
305  }
306 
307  //get transport mean free path (for multiple scattering)
308  G4double MSmfp1 = emCal.GetMeanFreePath(energy,particle,"msc",material);
309  G4double MSmfp2 = MSmfp1*density;
310 
311  //print transport mean free path
312  G4cout << "\n \n transport mean free path : "
313  << "\t" << std::setw(8) << G4BestUnit(MSmfp1,"Length")
314  << " (" << std::setw(8) << G4BestUnit(MSmfp2,"Mass/Surface") << ")";
315 
316  if (particle == G4Electron::Electron()) CriticalEnergy();
317 
318  G4cout << "\n-------------------------------------------------------------\n";
319  G4cout << G4endl;
320 
321  // reset default precision
322  G4cout.precision(prec);
323 }
324 
325 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
326 
327 void RunAction::EndOfRunAction(const G4Run* )
328 { }
329 
330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
331 
332 #include "G4ProductionCutsTable.hh"
333 
334 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
335 
336 void RunAction::GetCuts()
337 {
338  G4ProductionCutsTable* theCoupleTable =
340 
341  size_t numOfCouples = theCoupleTable->GetTableSize();
342  const G4MaterialCutsCouple* couple = 0;
343  G4int index = 0;
344  for (size_t i=0; i<numOfCouples; i++) {
345  couple = theCoupleTable->GetMaterialCutsCouple(i);
346  if (couple->GetMaterial() == fDetector->GetMaterial()) {index = i; break;}
347  }
348 
349  fRangeCut[0] =
350  (*(theCoupleTable->GetRangeCutsVector(idxG4GammaCut)))[index];
351  fRangeCut[1] =
352  (*(theCoupleTable->GetRangeCutsVector(idxG4ElectronCut)))[index];
353  fRangeCut[2] =
354  (*(theCoupleTable->GetRangeCutsVector(idxG4PositronCut)))[index];
355 
356  fEnergyCut[0] =
357  (*(theCoupleTable->GetEnergyCutsVector(idxG4GammaCut)))[index];
358  fEnergyCut[1] =
359  (*(theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut)))[index];
360  fEnergyCut[2] =
361  (*(theCoupleTable->GetEnergyCutsVector(idxG4PositronCut)))[index];
362 
363 }
364 
365 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
366 
368 {
369  // compute e- critical energy (Rossi definition) and Moliere radius.
370  // Review of Particle Physics - Eur. Phys. J. C3 (1998) page 147
371  //
372  G4EmCalculator emCal;
373 
374  const G4Material* material = fDetector->GetMaterial();
375  const G4double radl = material->GetRadlen();
376  G4double ekin = 5*MeV;
377  G4double deioni;
378  G4double err = 1., errmax = 0.001;
379  G4int iter = 0 , itermax = 10;
380  while (err > errmax && iter < itermax) {
381  iter++;
382  deioni = radl*
383  emCal.ComputeDEDX(ekin,G4Electron::Electron(),"eIoni",material);
384  err = std::abs(deioni - ekin)/ekin;
385  ekin = deioni;
386  }
387  G4cout << "\n \n critical energy (Rossi) : "
388  << "\t" << std::setw(8) << G4BestUnit(ekin,"Energy");
389 
390  //Pdg formula (only for single material)
391  G4double pdga[2] = { 610*MeV, 710*MeV };
392  G4double pdgb[2] = { 1.24, 0.92 };
393  G4double EcPdg = 0.;
394 
395  if (material->GetNumberOfElements() == 1) {
396  G4int istat = 0;
397  if (material->GetState() == kStateGas) istat = 1;
398  G4double Zeff = material->GetZ() + pdgb[istat];
399  EcPdg = pdga[istat]/Zeff;
400  G4cout << "\t\t\t (from Pdg formula : "
401  << std::setw(8) << G4BestUnit(EcPdg,"Energy") << ")";
402  }
403 
404  const G4double Es = 21.2052*MeV;
405  G4double rMolier1 = Es/ekin, rMolier2 = rMolier1*radl;
406  G4cout << "\n Moliere radius : "
407  << "\t" << std::setw(8) << rMolier1 << " X0 "
408  << "= " << std::setw(8) << G4BestUnit(rMolier2,"Length");
409 
410  if (material->GetNumberOfElements() == 1) {
411  G4double rMPdg = radl*Es/EcPdg;
412  G4cout << "\t (from Pdg formula : "
413  << std::setw(8) << G4BestUnit(rMPdg,"Length") << ")";
414  }
415 }
416 
417 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
static const double prec
Definition: RanecuEngine.cc:58
static constexpr double MeV
Definition: G4SIunits.hh:214
#define G4endl
Definition: G4ios.hh:61
G4double GetRadlen() const
Definition: G4Material.hh:221
const G4String & GetParticleName() const
G4double GetPDGCharge() const
G4ProcessVector * GetProcessList() const
G4double GetPDGMass() const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
const G4String & GetName() const
Definition: G4Material.hh:179
Float_t Z
G4double GetA() const
Definition: G4Material.cc:642
double G4double
Definition: G4Types.hh:76
G4double MaxEnergyForCSDARange() const
G4State GetState() const
Definition: G4Material.hh:182
const G4ParticleDefinition const G4Material *G4double range
static constexpr double electron_mass_c2
double energy
Definition: plottest35.C:25
double A(double temperature)
static G4Positron * Positron()
Definition: G4Positron.cc:94
The primary generator action class with particle gun.
G4double GetZ() const
Definition: G4Material.cc:629
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
static G4Electron * Electron()
Definition: G4Electron.cc:94
Definition: G4Run.hh:46
#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)
static G4ProductionCutsTable * GetProductionCutsTable()
int G4int
Definition: G4Types.hh:78
G4ProcessManager * GetProcessManager() const
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double GetParticleEnergy() const
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
const std::vector< G4double > * GetRangeCutsVector(size_t pcIdx) const
G4int size() const
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const G4Material * GetMaterial() const
G4ParticleDefinition * GetParticleDefinition() const
#define DBL_MAX
Definition: templates.hh:83
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
static G4EmParameters * Instance()
Simple detector construction with a box volume placed in a world.
size_t GetNumberOfElements() const
Definition: G4Material.hh:187
G4double GetDensity() const
Definition: G4Material.hh:181