Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 "G4PhysicsLogVector.hh"
44 #include "G4PhysicsVector.hh"
45 #include "G4MuonMinus.hh"
46 #include "G4MuonPlus.hh"
47 #include "G4NucleiProperties.hh"
48 #include "G4NistManager.hh"
49 #include "G4Log.hh"
50 #include "G4Exp.hh"
51 
52 // factory
53 #include "G4CrossSectionFactory.hh"
54 //
56 
57 G4PhysicsVector* G4KokoulinMuonNuclearXS::theCrossSection[] = {0};
58 
60  :G4VCrossSectionDataSet(Default_Name()),
61  LowestKineticEnergy(1*GeV), HighestKineticEnergy(1*PeV),
62  TotBin(60), CutFixed(0.2*GeV), isInitialized(false), isMaster(false)
63 {}
64 
66 {
67  if (isMaster) {
68  for(G4int i=0; i<MAXZMUN; ++i) {
69  delete theCrossSection[i];
70  theCrossSection[i] = 0;
71  }
72  }
73 }
74 
75 
76 void
78 {
79  outFile << "G4KokoulinMuonNuclearXS provides the total inelastic\n"
80  << "cross section for mu- and mu+ interactions with nuclei.\n"
81  << "R. Kokoulin's approximation of the Borog and Petrukhin double\n"
82  << "differential cross section at high energy and low Q**2 is integrated\n"
83  << "over the muon energy loss to get the total cross section as a\n"
84  << "function of muon kinetic energy\n" ;
85 }
86 
87 
88 G4bool
90  G4int, const G4Material*)
91 {
92  return true;
93 }
94 
95 void
97 {
98  if(!isInitialized) {
99  isInitialized = true;
100  for(G4int i=0; i<MAXZMUN; ++i) {
101  if(theCrossSection[i]) { return; }
102  }
103  isMaster = true;
104  }
105  if(isMaster) { BuildCrossSectionTable(); }
106 }
107 
109 {
110  G4double energy, A, Value;
111  G4int Z;
112 
114  const G4ElementTable* theElementTable = G4Element::GetElementTable();
115  G4NistManager* nistManager = G4NistManager::Instance();
116 
117  for (G4int j = 0; j < nEl; j++) {
118  Z = G4lrint((*theElementTable)[j]->GetZ());
119  A = nistManager->GetAtomicMassAmu(Z);
120  if(Z < MAXZMUN && !theCrossSection[Z]) {
121  theCrossSection[Z] = new G4PhysicsLogVector(LowestKineticEnergy,
122  HighestKineticEnergy,
123  TotBin);
124  for (G4int i = 0; i <= TotBin; ++i) {
125  energy = theCrossSection[Z]->Energy(i);
126  Value = ComputeMicroscopicCrossSection(energy, A);
127  theCrossSection[Z]->PutValue(i,Value);
128  }
129  }
130  }
131 }
132 
133 G4double G4KokoulinMuonNuclearXS::
134 ComputeMicroscopicCrossSection(G4double KineticEnergy, G4double A)
135 {
136  // Calculate cross section (differential in muon incident kinetic energy) by
137  // integrating the double differential cross section over the energy loss
138 
139  static const G4double xgi[] =
140  {0.0199,0.1017,0.2372,0.4083,0.5917,0.7628,0.8983,0.9801};
141  static const G4double wgi[] =
142  {0.0506,0.1112,0.1569,0.1813,0.1813,0.1569,0.1112,0.0506};
143  static const G4double ak1 = 6.9;
144  static const G4double ak2 = 1.0;
145 
147 
148  G4double CrossSection = 0.0;
149  if (KineticEnergy <= CutFixed) return CrossSection;
150 
151  G4double epmin = CutFixed;
152  G4double epmax = KineticEnergy + Mass - 0.5*proton_mass_c2;
153  if (epmax <= epmin) return CrossSection; // NaN bug correction
154 
155  G4double aaa = G4Log(epmin);
156  G4double bbb = G4Log(epmax);
157  G4int kkk = std::max(1,G4int((bbb-aaa)/ak1 +ak2));
158  G4double hhh = (bbb-aaa)/kkk ;
159  G4double epln;
160  G4double ep;
161  G4double x;
162 
163  for (G4int l = 0; l < kkk; ++l) {
164  x = aaa + hhh*l;
165  for (G4int ll = 0; ll < 8; ++ll) {
166  epln=x+xgi[ll]*hhh;
167  ep = G4Exp(epln);
168  CrossSection +=
169  ep*wgi[ll]*ComputeDDMicroscopicCrossSection(KineticEnergy, 0, A, ep);
170  }
171  }
172 
173  CrossSection *= hhh ;
174  if (CrossSection < 0.) { CrossSection = 0.; }
175  return CrossSection;
176 }
177 
181 {
182  // Calculate the double-differential microscopic cross section (in muon
183  // incident kinetic energy and energy loss) using the cross section formula
184  // of R.P. Kokoulin (18/01/98)
185 
186  static const G4double alam2 = 0.400*GeV*GeV;
187  static const G4double alam = 0.632456*GeV;
188  static const G4double coeffn = fine_structure_const/pi;
189 
190  G4double ParticleMass = G4MuonMinus::MuonMinus()->GetPDGMass();
191  G4double TotalEnergy = KineticEnergy + ParticleMass;
192 
193  G4double DCrossSection = 0.;
194 
195  if ((epsilon >= TotalEnergy - 0.5*proton_mass_c2) ||
196  (epsilon <= CutFixed) ) { return DCrossSection; }
197 
198  G4double ep = epsilon/GeV;
199  G4double aeff = 0.22*A+0.78*G4Exp(0.89*G4Log(A)); //shadowing
200  G4double sigph = (49.2+11.1*G4Log(ep)+151.8/std::sqrt(ep))*microbarn;
201 
202  G4double v = epsilon/TotalEnergy;
203  G4double v1 = 1.-v;
204  G4double v2 = v*v;
205  G4double mass2 = ParticleMass*ParticleMass;
206 
207  G4double up = TotalEnergy*TotalEnergy*v1/mass2*(1.+mass2*v2/(alam2*v1));
208  G4double down = 1.+epsilon/alam*(1.+alam/(2.*proton_mass_c2)+epsilon/alam);
209 
210  DCrossSection = coeffn*aeff*sigph/epsilon*
211  (-v1+(v1+0.5*v2*(1.+2.*mass2/alam2))*G4Log(up/down));
212 
213  if (DCrossSection < 0.) { DCrossSection = 0.; }
214  return DCrossSection;
215 }
216 
219  G4int Z, const G4Material*)
220 {
221  return theCrossSection[Z]->Value(aPart->GetKineticEnergy());
222 }
223 
G4double GetKineticEnergy() const
static const G4double ak2
G4double ComputeDDMicroscopicCrossSection(G4double incidentKE, G4double Z, G4double A, G4double epsilon)
virtual void CrossSectionDescription(std::ostream &) const
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
double A(double temperature)
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
G4double GetElementCrossSection(const G4DynamicParticle *particle, G4int Z, const G4Material *)
bool G4bool
Definition: G4Types.hh:79
static const G4double ak1
void PutValue(size_t index, G4double theValue)
float proton_mass_c2
Definition: hepunit.py:275
static constexpr double PeV
Definition: G4SIunits.hh:219
const G4int MAXZMUN
G4double Energy(size_t index) const
G4bool IsElementApplicable(const G4DynamicParticle *particle, G4int Z, const G4Material *)
G4double Value(G4double theEnergy, size_t &lastidx) const
tuple v
Definition: test.py:18
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define G4_DECLARE_XS_FACTORY(cross_section)
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
G4double GetAtomicMassAmu(const G4String &symb) const
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static constexpr double pi
Definition: G4SIunits.hh:75
void BuildPhysicsTable(const G4ParticleDefinition &)
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
G4bool isInitialized()
double epsilon(double density, double temperature)
static constexpr double microbarn
Definition: G4SIunits.hh:107