Geant4_10
G4KokoulinMuonNuclearXS.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 // $Id: $
27 //
28 // Author: D.H. Wright
29 // Date: 1 February 2011
30 //
31 // Modified:
32 //
33 // 19 Aug 2011, V.Ivanchenko move to new design and make x-section per element
34 
35 // Description: use Kokoulin's parameterized calculation of virtual
36 // photon production cross section and conversion to
37 // real photons.
38 
40 
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4PhysicsTable.hh"
44 #include "G4PhysicsLogVector.hh"
45 
46 
48  :G4VCrossSectionDataSet("KokoulinMuonNuclearXS"), theCrossSectionTable(0),
49  LowestKineticEnergy(1*GeV), HighestKineticEnergy(1*PeV),
50  TotBin(60), CutFixed(0.2*GeV)
51 {
53 }
54 
56 {
57  if (theCrossSectionTable) {
58  theCrossSectionTable->clearAndDestroy();
59  delete theCrossSectionTable;
60  }
61 }
62 
63 G4bool
65  G4int, const G4Material*)
66 {
67  return true;
68 }
69 
70 
72 {
73  G4double lowEdgeEnergy;
74  G4PhysicsLogVector* ptrVector;
75 
77  const G4ElementTable* theElementTable = G4Element::GetElementTable();
78 
79  theCrossSectionTable = new G4PhysicsTable(nEl);
80 
81  G4double AtomicNumber;
82  G4double AtomicWeight;
83  G4double Value;
84 
85  for (G4int j = 0; j < nEl; j++) {
86  ptrVector = new G4PhysicsLogVector(LowestKineticEnergy,
87  HighestKineticEnergy, TotBin);
88  AtomicNumber = (*theElementTable)[j]->GetZ();
89  AtomicWeight = (*theElementTable)[j]->GetA();
90 
91  for (G4int i = 0; i <= TotBin; ++i) {
92  lowEdgeEnergy = ptrVector->Energy(i);
93  Value = ComputeMicroscopicCrossSection(lowEdgeEnergy,
94  AtomicNumber, AtomicWeight);
95  ptrVector->PutValue(i,Value);
96  }
97 
98  theCrossSectionTable->insertAt(j, ptrVector);
99  }
100 
101  // Build (Z,element) look-up table (for use in GetZandACrossSection)
102  for (G4int j = 0; j < nEl; j++) {
103  G4int ZZ = G4int((*theElementTable)[j]->GetZ());
104  zelMap[ZZ] = j;
105  }
106 }
107 
108 G4double G4KokoulinMuonNuclearXS::
109 ComputeMicroscopicCrossSection(G4double KineticEnergy,
110  G4double AtomicNumber, G4double AtomicWeight)
111 {
112  // Calculate cross section (differential in muon incident kinetic energy) by
113  // integrating the double differential cross section over the energy loss
114 
115  static const G4double xgi[] = {0.0199,0.1017,0.2372,0.4083,0.5917,0.7628,0.8983,0.9801};
116  static const G4double wgi[] = {0.0506,0.1112,0.1569,0.1813,0.1813,0.1569,0.1112,0.0506};
117  static const G4double ak1 = 6.9;
118  static const G4double ak2 = 1.0;
119 
121 
122  G4double CrossSection = 0.0;
123  if (AtomicNumber < 1.) return CrossSection;
124  if (KineticEnergy <= CutFixed) return CrossSection;
125 
126  G4double epmin = CutFixed;
127  G4double epmax = KineticEnergy + Mass - 0.5*proton_mass_c2;
128  if (epmax <= epmin) return CrossSection; // NaN bug correction
129 
130  G4double aaa = std::log(epmin) ;
131  G4double bbb = std::log(epmax) ;
132  G4int kkk = G4int((bbb-aaa)/ak1 +ak2) ;
133  G4double hhh = (bbb-aaa)/kkk ;
134  G4double epln;
135  G4double ep;
136  G4double x;
137 
138  for (G4int l = 0; l < kkk; l++) {
139  x = aaa + hhh*l;
140  for (G4int ll = 0; ll < 8; ll++) {
141  epln=x+xgi[ll]*hhh;
142  ep = std::exp(epln);
143  CrossSection += ep*wgi[ll]*ComputeDDMicroscopicCrossSection(KineticEnergy,
144  AtomicNumber, AtomicWeight, ep);
145  }
146  }
147 
148  CrossSection *= hhh ;
149  if (CrossSection < 0.) CrossSection = 0.;
150  return CrossSection;
151 }
152 
155  G4double, G4double AtomicWeight,
156  G4double epsilon)
157 {
158  // Calculate the double-differential microscopic cross section (in muon
159  // incident kinetic energy and energy loss) using the cross section formula
160  // of R.P. Kokoulin (18/01/98)
161 
162  static const G4double alam2 = 0.400*GeV*GeV;
163  static const G4double alam = 0.632456*GeV;
164  static const G4double coeffn = fine_structure_const/pi;
165 
166  G4double ParticleMass = G4MuonMinus::MuonMinus()->GetPDGMass();
167  G4double TotalEnergy = KineticEnergy + ParticleMass;
168 
169  G4double DCrossSection = 0.;
170 
171  if ((epsilon >= TotalEnergy - 0.5*proton_mass_c2) ||
172  (epsilon <= CutFixed) ) return DCrossSection;
173 
174  G4double ep = epsilon/GeV;
175  G4double a = AtomicWeight/(g/mole);
176  G4double aeff = 0.22*a+0.78*std::exp(0.89*std::log(a)); //shadowing
177  G4double sigph = (49.2+11.1*std::log(ep)+151.8/std::sqrt(ep))*microbarn;
178 
179  G4double v = epsilon/TotalEnergy;
180  G4double v1 = 1.-v;
181  G4double v2 = v*v;
182  G4double mass2 = ParticleMass*ParticleMass;
183 
184  G4double up = TotalEnergy*TotalEnergy*v1/mass2*(1.+mass2*v2/(alam2*v1));
185  G4double down = 1.+epsilon/alam*(1.+alam/(2.*proton_mass_c2)+epsilon/alam);
186 
187  DCrossSection = coeffn*aeff*sigph/epsilon*
188  (-v1+(v1+0.5*v2*(1.+2.*mass2/alam2))*std::log(up/down));
189 
190  if (DCrossSection < 0.) DCrossSection = 0.;
191  return DCrossSection;
192 }
193 
196  G4int ZZ, const G4Material*)
197 {
198  G4int j = zelMap[ZZ];
199  return (*theCrossSectionTable)[j]->Value(aPart->GetKineticEnergy());
200 }
201 
tuple a
Definition: test.py:11
G4double GetKineticEnergy() const
G4double ComputeDDMicroscopicCrossSection(G4double incidentKE, G4double, G4double AtomicWeight, G4double epsilon)
tuple x
Definition: test.py:50
int microbarn
Definition: hepunit.py:41
int G4int
Definition: G4Types.hh:78
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
G4double GetElementCrossSection(const G4DynamicParticle *particle, G4int Z, const G4Material *)
bool G4bool
Definition: G4Types.hh:79
void PutValue(size_t index, G4double theValue)
float proton_mass_c2
Definition: hepunit.py:275
G4double Energy(size_t index) const
G4bool IsElementApplicable(const G4DynamicParticle *particle, G4int Z, const G4Material *)
tuple v
Definition: test.py:18
G4double GetPDGMass() const
void insertAt(size_t, G4PhysicsVector *)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
TH1F * hhh
Definition: AddMC.C:5
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
void clearAndDestroy()