Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RToEConvForGamma.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 //
26 //
27 // $Id$
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file/ History:
32 // 5 Oct. 2002, H.Kuirashige : Structure created based on object model
33 // --------------------------------------------------------------
34 
35 #include "G4RToEConvForGamma.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4Material.hh"
39 #include "G4PhysicsLogVector.hh"
40 
41 #include "G4ios.hh"
42 #include "G4SystemOfUnits.hh"
43 
45 {
47  if (theParticle ==0) {
48 #ifdef G4VERBOSE
49  if (GetVerboseLevel()>0) {
50  G4cout << " G4RToEConvForGamma::G4RToEConvForGamma() ";
51  G4cout << " Gamma is not defined !!" << G4endl;
52  }
53 #endif
54  }
55 }
56 
58 {
59 }
60 
61 
62 // ***********************************************************************
63 // ******************* BuildAbsorptionLengthVector ***********************
64 // ***********************************************************************
66  const G4Material* aMaterial,
67  G4RangeVector* absorptionLengthVector )
68 {
69  // fill the absorption length vector for this material
70  // absorption length is defined here as
71  //
72  // absorption length = 5./ macroscopic absorption cross section
73  //
74  const G4CrossSectionTable* aCrossSectionTable = (G4CrossSectionTable*)(theLossTable);
75  const G4ElementVector* elementVector = aMaterial->GetElementVector();
76  const G4double* atomicNumDensityVector = aMaterial->GetAtomicNumDensityVector();
77 
78  // fill absorption length vector
79  G4int NumEl = aMaterial->GetNumberOfElements();
80  G4double absorptionLengthMax = 0.0;
81  for (size_t ibin=0; ibin<size_t(TotBin); ibin++) {
82  G4double SIGMA = 0. ;
83  for (size_t iel=0; iel<size_t(NumEl); iel++) {
84  G4int IndEl = (*elementVector)[iel]->GetIndex();
85  SIGMA += atomicNumDensityVector[iel]*
86  (*((*aCrossSectionTable)[IndEl]))[ibin];
87  }
88  // absorption length=5./SIGMA
89  absorptionLengthVector->PutValue(ibin, 5./SIGMA);
90  if (absorptionLengthMax < 5./SIGMA ) absorptionLengthMax = 5./SIGMA;
91  }
92 }
93 
94 
95 
96 // ***********************************************************************
97 // ********************** ComputeCrossSection ****************************
98 // ***********************************************************************
100  G4double KineticEnergy) const
101 {
102  // Compute the "absorption" cross section of the photon "absorption"
103  // cross section means here the sum of the cross sections of the
104  // pair production, Compton scattering and photoelectric processes
105  static G4double Z;
106  const G4double t1keV = 1.*keV;
107  const G4double t200keV = 200.*keV;
108  const G4double t100MeV = 100.*MeV;
109 
110  static G4double s200keV, s1keV;
111  static G4double tmin, tlow;
112  static G4double smin, slow;
113  static G4double cmin, clow, chigh;
114  // compute Z dependent quantities in the case of a new AtomicNumber
115  if(std::abs(AtomicNumber-Z)>0.1) {
116  Z = AtomicNumber;
117  G4double Zsquare = Z*Z;
118  G4double Zlog = std::log(Z);
119  G4double Zlogsquare = Zlog*Zlog;
120 
121  s200keV = (0.2651-0.1501*Zlog+0.02283*Zlogsquare)*Zsquare;
122  tmin = (0.552+218.5/Z+557.17/Zsquare)*MeV;
123  smin = (0.01239+0.005585*Zlog-0.000923*Zlogsquare)*std::exp(1.5*Zlog);
124  cmin=std::log(s200keV/smin)/(std::log(tmin/t200keV)*std::log(tmin/t200keV));
125  tlow = 0.2*std::exp(-7.355/std::sqrt(Z))*MeV;
126  slow = s200keV*std::exp(0.042*Z*std::log(t200keV/tlow)*std::log(t200keV/tlow));
127  s1keV = 300.*Zsquare;
128  clow =std::log(s1keV/slow)/std::log(tlow/t1keV);
129 
130  chigh=(7.55e-5-0.0542e-5*Z)*Zsquare*Z/std::log(t100MeV/tmin);
131  }
132 
133  // calculate the cross section (using an approximate empirical formula)
134  G4double xs;
135  if ( KineticEnergy<tlow ) {
136  if(KineticEnergy<t1keV) xs = slow*std::exp(clow*std::log(tlow/t1keV));
137  else xs = slow*std::exp(clow*std::log(tlow/KineticEnergy));
138 
139  } else if ( KineticEnergy<t200keV ) {
140  xs = s200keV
141  * std::exp(0.042*Z*std::log(t200keV/KineticEnergy)*std::log(t200keV/KineticEnergy));
142 
143  } else if( KineticEnergy<tmin ){
144  xs = smin
145  * std::exp(cmin*std::log(tmin/KineticEnergy)*std::log(tmin/KineticEnergy));
146 
147  } else {
148  xs = smin + chigh*std::log(KineticEnergy/tmin);
149 
150  }
151  return xs * barn;
152 }
153