Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RunAction.cc
Go to the documentation of this file.
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 99373 2016-09-20 07:13:41Z 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 "G4UnitsTable.hh"
41 #include "G4EmCalculator.hh"
42 #include "G4Electron.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
73  G4ParticleDefinition* particle = fPrimary->GetParticleGun()
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  // print list of processes
139  G4cout << "\n processes : ";
140  for (size_t j=0; j<emName.size();j++)
141  G4cout << "\t" << std::setw(13) << emName[j] << "\t";
142  G4cout << "\t" << std::setw(13) <<"total";
143 
144  //compute cross section per atom (only for single material)
145  if (material->GetNumberOfElements() == 1) {
146  G4double Z = material->GetZ();
147  G4double A = material->GetA();
148 
149  std::vector<G4double> sigma0;
150  G4double sig, sigtot = 0.;
151 
152  for (size_t j=0; j<emName.size();j++) {
153  sig = emCal.ComputeCrossSectionPerAtom
154  (energy,particle,emName[j],Z,A,enerCut[j]);
155  sigtot += sig;
156  sigma0.push_back(sig);
157  }
158  sigma0.push_back(sigtot);
159 
160  G4cout << "\n \n cross section per atom : ";
161  for (size_t j=0; j<sigma0.size();j++) {
162  G4cout << "\t" << std::setw(13) << G4BestUnit(sigma0[j], "Surface");
163  }
164  G4cout << G4endl;
165  }
166 
167  //get cross section per volume
168  std::vector<G4double> sigma0;
169  std::vector<G4double> sigma1;
170  std::vector<G4double> sigma2;
171  G4double Sig, SigtotComp = 0., Sigtot = 0.;
172 
173  for (size_t j=0; j<emName.size();j++) {
174  Sig = emCal.ComputeCrossSectionPerVolume
175  (energy,particle,emName[j],material,enerCut[j]);
176  SigtotComp += Sig;
177  sigma0.push_back(Sig);
178  Sig = emCal.GetCrossSectionPerVolume(energy,particle,emName[j],material);
179  Sigtot += Sig;
180  sigma1.push_back(Sig);
181  sigma2.push_back(Sig/density);
182  }
183  sigma0.push_back(SigtotComp);
184  sigma1.push_back(Sigtot);
185  sigma2.push_back(Sigtot/density);
186 
187  //print cross sections
188  G4cout << "\n \n compCrossSectionPerVolume : ";
189  for (size_t j=0; j<sigma0.size();j++) {
190  G4cout << "\t" << std::setw(13) << sigma0[j]*cm << " cm^-1";
191  }
192  G4cout << "\n cross section per volume : ";
193  for (size_t j=0; j<sigma1.size();j++) {
194  G4cout << "\t" << std::setw(13) << sigma1[j]*cm << " cm^-1";
195  }
196 
197  G4cout << "\n cross section per mass : ";
198  for (size_t j=0; j<sigma2.size();j++) {
199  G4cout << "\t" << std::setw(13)
200  << G4BestUnit(sigma2[j], "Surface/Mass");
201  }
202 
203  //print mean free path
204 
206 
207  G4cout << "\n \n mean free path : ";
208  for (size_t j=0; j<sigma1.size();j++) {
209  lambda = DBL_MAX;
210  if (sigma1[j] > 0.) lambda = 1/sigma1[j];
211  G4cout << "\t" << std::setw(13) << G4BestUnit( lambda, "Length");
212  }
213 
214  //mean free path (g/cm2)
215  G4cout << "\n (g/cm2) : ";
216  for (size_t j=0; j<sigma2.size();j++) {
217  lambda = DBL_MAX;
218  if (sigma2[j] > 0.) lambda = 1/sigma2[j];
219  G4cout << "\t" << std::setw(13) << G4BestUnit( lambda, "Mass/Surface");
220  }
221  G4cout << G4endl;
222 
223  if (charge == 0.) {
224  G4cout.precision(prec);
225  G4cout << "\n-----------------------------------------------------------\n"
226  << G4endl;
227  return;
228  }
229 
230  //get stopping power
231  std::vector<G4double> dedx1;
232  std::vector<G4double> dedx2;
233  G4double dedx, dedxtot = 0.;
234  size_t nproc = emName.size();
235 
236  for (size_t j=0; j<nproc; j++) {
237  dedx = emCal.ComputeDEDX(energy,particle,emName[j],material,enerCut[j]);
238  dedxtot += dedx;
239  dedx1.push_back(dedx);
240  dedx2.push_back(dedx/density);
241  }
242  dedx1.push_back(dedxtot);
243  dedx2.push_back(dedxtot/density);
244 
245  //print stopping power
246  G4cout << "\n \n restricted dE/dx : ";
247  for (size_t j=0; j<=nproc; j++) {
248  G4cout << "\t" << std::setw(13)
249  << G4BestUnit(dedx1[j],"Energy/Length");
250  }
251 
252  G4cout << "\n (MeV/g/cm2) : ";
253  for (size_t j=0; j<=nproc; j++) {
254  G4cout << "\t" << std::setw(13)
255  << G4BestUnit(dedx2[j],"Energy*Surface/Mass");
256  }
257  dedxtot = 0.;
258 
259  for (size_t j=0; j<nproc; j++) {
260  dedx = emCal.ComputeDEDX(energy,particle,emName[j],material,energy);
261  dedxtot += dedx;
262  dedx1[j] = dedx;
263  dedx2[j] = dedx/density;
264  }
265  dedx1[nproc] = dedxtot;
266  dedx2[nproc] = dedxtot/density;
267 
268  //print stopping power
269  G4cout << "\n \n unrestricted dE/dx : ";
270  for (size_t j=0; j<=nproc; j++) {
271  G4cout << "\t" << std::setw(13) << G4BestUnit(dedx1[j],"Energy/Length");
272  }
273 
274  G4cout << "\n (MeV/g/cm2) : ";
275  for (size_t j=0; j<=nproc; j++) {
276  G4cout << "\t" << std::setw(13)
277  << G4BestUnit(dedx2[j],"Energy*Surface/Mass");
278  }
279 
280  //get range from restricted dedx
281  G4double range1 = emCal.GetRangeFromRestricteDEDX(energy,particle,material);
282  G4double range2 = range1*density;
283 
284  //print range
285  G4cout << "\n \n range from restrict dE/dx: "
286  << "\t" << std::setw(8) << G4BestUnit(range1,"Length")
287  << " (" << std::setw(8) << G4BestUnit(range2,"Mass/Surface") << ")";
288 
289  //get range from full dedx
290  if(energy < GeV) {
291  G4double Range1 = emCal.GetCSDARange(energy,particle,material);
292  G4double Range2 = Range1*density;
293 
294  G4cout << "\n range from full dE/dx : "
295  << "\t" << std::setw(8) << G4BestUnit(Range1,"Length")
296  << " (" << std::setw(8) << G4BestUnit(Range2,"Mass/Surface") << ")";
297  }
298 
299  //get transport mean free path (for multiple scattering)
300  G4double MSmfp1 = emCal.GetMeanFreePath(energy,particle,"msc",material);
301  G4double MSmfp2 = MSmfp1*density;
302 
303  //print transport mean free path
304  G4cout << "\n \n transport mean free path : "
305  << "\t" << std::setw(8) << G4BestUnit(MSmfp1,"Length")
306  << " (" << std::setw(8) << G4BestUnit(MSmfp2,"Mass/Surface") << ")";
307 
308  if (particle == G4Electron::Electron()) CriticalEnergy();
309 
310  G4cout << "\n-------------------------------------------------------------\n";
311  G4cout << G4endl;
312 
313  // reset default precision
314  G4cout.precision(prec);
315 }
316 
317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
318 
319 void RunAction::EndOfRunAction(const G4Run* )
320 { }
321 
322 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
323 
324 #include "G4ProductionCutsTable.hh"
325 
326 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
327 
328 void RunAction::GetCuts()
329 {
330  G4ProductionCutsTable* theCoupleTable =
332 
333  size_t numOfCouples = theCoupleTable->GetTableSize();
334  const G4MaterialCutsCouple* couple = 0;
335  G4int index = 0;
336  for (size_t i=0; i<numOfCouples; i++) {
337  couple = theCoupleTable->GetMaterialCutsCouple(i);
338  if (couple->GetMaterial() == fDetector->GetMaterial()) {index = i; break;}
339  }
340 
341  fRangeCut[0] =
342  (*(theCoupleTable->GetRangeCutsVector(idxG4GammaCut)))[index];
343  fRangeCut[1] =
344  (*(theCoupleTable->GetRangeCutsVector(idxG4ElectronCut)))[index];
345  fRangeCut[2] =
346  (*(theCoupleTable->GetRangeCutsVector(idxG4PositronCut)))[index];
347 
348  fEnergyCut[0] =
349  (*(theCoupleTable->GetEnergyCutsVector(idxG4GammaCut)))[index];
350  fEnergyCut[1] =
351  (*(theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut)))[index];
352  fEnergyCut[2] =
353  (*(theCoupleTable->GetEnergyCutsVector(idxG4PositronCut)))[index];
354 
355 }
356 
357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
358 
360 {
361  // compute e- critical energy (Rossi definition) and Moliere radius.
362  // Review of Particle Physics - Eur. Phys. J. C3 (1998) page 147
363  //
364  G4EmCalculator emCal;
365 
366  const G4Material* material = fDetector->GetMaterial();
367  const G4double radl = material->GetRadlen();
368  G4double ekin = 5*MeV;
369  G4double deioni;
370  G4double err = 1., errmax = 0.001;
371  G4int iter = 0 , itermax = 10;
372  while (err > errmax && iter < itermax) {
373  iter++;
374  deioni = radl*
375  emCal.ComputeDEDX(ekin,G4Electron::Electron(),"eIoni",material);
376  err = std::abs(deioni - ekin)/ekin;
377  ekin = deioni;
378  }
379  G4cout << "\n \n critical energy (Rossi) : "
380  << "\t" << std::setw(8) << G4BestUnit(ekin,"Energy");
381 
382  //Pdg formula (only for single material)
383  G4double pdga[2] = { 610*MeV, 710*MeV };
384  G4double pdgb[2] = { 1.24, 0.92 };
385  G4double EcPdg = 0.;
386 
387  if (material->GetNumberOfElements() == 1) {
388  G4int istat = 0;
389  if (material->GetState() == kStateGas) istat = 1;
390  G4double Zeff = material->GetZ() + pdgb[istat];
391  EcPdg = pdga[istat]/Zeff;
392  G4cout << "\t\t\t (from Pdg formula : "
393  << std::setw(8) << G4BestUnit(EcPdg,"Energy") << ")";
394  }
395 
396  const G4double Es = 21.2052*MeV;
397  G4double rMolier1 = Es/ekin, rMolier2 = rMolier1*radl;
398  G4cout << "\n Moliere radius : "
399  << "\t" << std::setw(8) << rMolier1 << " X0 "
400  << "= " << std::setw(8) << G4BestUnit(rMolier2,"Length");
401 
402  if (material->GetNumberOfElements() == 1) {
403  G4double rMPdg = radl*Es/EcPdg;
404  G4cout << "\t (from Pdg formula : "
405  << std::setw(8) << G4BestUnit(rMPdg,"Length") << ")";
406  }
407 }
408 
409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
The primary generator action class with particle gun.
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4double GetZ() const
Definition: G4Material.cc:623
void BeginOfRunAction(const G4Run *)
Definition: RunAction.cc:57
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetDensity() const
Definition: G4Material.hh:180
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
string material
Definition: eplot.py:19
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
const G4ParticleDefinition const G4Material *G4double range
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
void EndOfRunAction(const G4Run *)
Definition: RunAction.cc:260
static constexpr double cm
Definition: G4SIunits.hh:119
void CriticalEnergy()
Definition: RunAction.cc:300
Definition: G4Run.hh:46
float electron_mass_c2
Definition: hepunit.py:274
const std::vector< G4double > * GetRangeCutsVector(size_t pcIdx) const
G4double GetRadlen() const
Definition: G4Material.hh:220
~RunAction()
Definition: RunAction.cc:52
G4double GetA() const
Definition: G4Material.cc:636
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4int size() const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4ProcessManager * GetProcessManager() const
G4ParticleGun * GetParticleGun()
G4double energy(const ThreeVector &p, const G4double m)
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4ParticleDefinition * GetParticleDefinition() const
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
Detector construction class to define materials and geometry.
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4State GetState() const
Definition: G4Material.hh:181
double G4double
Definition: G4Types.hh:76
void GetCuts()
Definition: RunAction.cc:269
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double GetParticleEnergy() const
G4ProcessVector * GetProcessList() const
const G4Material * GetMaterial() const