Geant4  10.00.p02
G4GlauberGribovCrossSection.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 // author: V. Grichine
27 //
28 // 17.07.06 V. Grichine - first implementation
29 // 22.01.07 V.Ivanchenko - add interface with Z and A
30 // 05.03.07 V.Ivanchenko - add IfZAApplicable
31 // 11.06.10 V. Grichine - update for antiprotons
32 // 10.11.11 V. Grichine - update for kaons
33 
35 
36 #include "G4PhysicalConstants.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4ParticleTable.hh"
39 #include "G4IonTable.hh"
40 #include "G4ParticleDefinition.hh"
41 #include "G4HadronNucleonXsc.hh"
42 
43 // factory
44 #include "G4CrossSectionFactory.hh"
45 //
47 
49 
51 
52 1.0, 1.0, 1.118517e+00, 1.082002e+00, 1.116171e+00, 1.078747e+00, 1.061315e+00,
53 1.058205e+00, 1.082663e+00, 1.068500e+00, 1.076912e+00, 1.083475e+00, 1.079117e+00,
54 1.071856e+00, 1.071990e+00, 1.073774e+00, 1.079356e+00, 1.081314e+00, 1.082056e+00,
55 1.090772e+00, 1.096776e+00, 1.095828e+00, 1.097678e+00, 1.099157e+00, 1.103677e+00,
56 1.105132e+00, 1.109806e+00, 1.110816e+00, 1.117378e+00, 1.115165e+00, 1.115710e+00,
57 1.111855e+00, 1.110482e+00, 1.110112e+00, 1.106676e+00, 1.108706e+00, 1.105549e+00,
58 1.106318e+00, 1.106242e+00, 1.107672e+00, 1.107342e+00, 1.108119e+00, 1.106655e+00,
59 1.102588e+00, 1.096657e+00, 1.092920e+00, 1.086629e+00, 1.083592e+00, 1.076030e+00,
60 1.083777e+00, 1.089460e+00, 1.086545e+00, 1.079924e+00, 1.082218e+00, 1.077798e+00,
61 1.077062e+00, 1.072825e+00, 1.072241e+00, 1.072104e+00, 1.072490e+00, 1.069829e+00,
62 1.070398e+00, 1.065458e+00, 1.064968e+00, 1.060524e+00, 1.060048e+00, 1.057620e+00,
63 1.056428e+00, 1.055366e+00, 1.055017e+00, 1.052304e+00, 1.051767e+00, 1.049728e+00,
64 1.048745e+00, 1.047399e+00, 1.045876e+00, 1.042972e+00, 1.041824e+00, 1.039993e+00,
65 1.039021e+00, 1.036627e+00, 1.034176e+00, 1.032526e+00, 1.033633e+00, 1.036107e+00,
66 1.037803e+00, 1.031266e+00, 1.032991e+00, 1.033284e+00, 1.035015e+00, 1.033945e+00,
67 1.037075e+00, 1.034721e+00
68 
69 };
70 
72 
73 1.0, 1.0, 1.167421e+00, 1.156250e+00, 1.205364e+00, 1.154225e+00, 1.120391e+00,
74 1.124632e+00, 1.129460e+00, 1.107863e+00, 1.102152e+00, 1.104593e+00, 1.100285e+00,
75 1.098450e+00, 1.092677e+00, 1.101124e+00, 1.106461e+00, 1.115049e+00, 1.123903e+00,
76 1.126661e+00, 1.131259e+00, 1.133949e+00, 1.134185e+00, 1.133767e+00, 1.132813e+00,
77 1.131515e+00, 1.130338e+00, 1.134171e+00, 1.139206e+00, 1.141474e+00, 1.142189e+00,
78 1.140725e+00, 1.140100e+00, 1.139848e+00, 1.137674e+00, 1.138645e+00, 1.136339e+00,
79 1.136439e+00, 1.135946e+00, 1.136431e+00, 1.135702e+00, 1.135703e+00, 1.134113e+00,
80 1.131935e+00, 1.128381e+00, 1.126373e+00, 1.122453e+00, 1.120908e+00, 1.115953e+00,
81 1.115947e+00, 1.114426e+00, 1.111749e+00, 1.106207e+00, 1.107494e+00, 1.103622e+00,
82 1.102576e+00, 1.098816e+00, 1.097889e+00, 1.097306e+00, 1.097130e+00, 1.094578e+00,
83 1.094552e+00, 1.090222e+00, 1.089358e+00, 1.085409e+00, 1.084560e+00, 1.082182e+00,
84 1.080773e+00, 1.079464e+00, 1.078724e+00, 1.076121e+00, 1.075235e+00, 1.073159e+00,
85 1.071920e+00, 1.070395e+00, 1.069503e+00, 1.067525e+00, 1.066919e+00, 1.065779e+00,
86 1.065319e+00, 1.063730e+00, 1.062092e+00, 1.061085e+00, 1.059908e+00, 1.059815e+00,
87 1.059109e+00, 1.051920e+00, 1.051258e+00, 1.049473e+00, 1.048823e+00, 1.045984e+00,
88 1.046435e+00, 1.042614e+00
89 
90 };
91 
93 
94 1.0, 1.0,
95 1.118515e+00, 1.082000e+00, 1.116169e+00, 1.078745e+00, 1.061313e+00, 1.058203e+00,
96 1.082661e+00, 1.068498e+00, 1.076910e+00, 1.083474e+00, 1.079115e+00, 1.071854e+00,
97 1.071988e+00, 1.073772e+00, 1.079355e+00, 1.081312e+00, 1.082054e+00, 1.090770e+00,
98 1.096774e+00, 1.095827e+00, 1.097677e+00, 1.099156e+00, 1.103676e+00, 1.105130e+00,
99 1.109805e+00, 1.110814e+00, 1.117377e+00, 1.115163e+00, 1.115708e+00, 1.111853e+00,
100 1.110480e+00, 1.110111e+00, 1.106674e+00, 1.108705e+00, 1.105548e+00, 1.106317e+00,
101 1.106241e+00, 1.107671e+00, 1.107341e+00, 1.108118e+00, 1.106654e+00, 1.102586e+00,
102 1.096655e+00, 1.092918e+00, 1.086628e+00, 1.083590e+00, 1.076028e+00, 1.083776e+00,
103 1.089458e+00, 1.086543e+00, 1.079923e+00, 1.082216e+00, 1.077797e+00, 1.077061e+00,
104 1.072824e+00, 1.072239e+00, 1.072103e+00, 1.072488e+00, 1.069828e+00, 1.070396e+00,
105 1.065456e+00, 1.064966e+00, 1.060523e+00, 1.060047e+00, 1.057618e+00, 1.056427e+00,
106 1.055365e+00, 1.055016e+00, 1.052303e+00, 1.051766e+00, 1.049727e+00, 1.048743e+00,
107 1.047397e+00, 1.045875e+00, 1.042971e+00, 1.041823e+00, 1.039992e+00, 1.039019e+00,
108 1.036626e+00, 1.034175e+00, 1.032525e+00, 1.033632e+00, 1.036106e+00, 1.037802e+00,
109 1.031265e+00, 1.032990e+00, 1.033283e+00, 1.035014e+00, 1.033944e+00, 1.037074e+00,
110 1.034720e+00
111 
112 };
113 
115 
116 1.0, 1.0,
117 1.167419e+00, 1.156248e+00, 1.205362e+00, 1.154224e+00, 1.120390e+00, 1.124630e+00,
118 1.129459e+00, 1.107861e+00, 1.102151e+00, 1.104591e+00, 1.100284e+00, 1.098449e+00,
119 1.092675e+00, 1.101122e+00, 1.106460e+00, 1.115048e+00, 1.123902e+00, 1.126659e+00,
120 1.131258e+00, 1.133948e+00, 1.134183e+00, 1.133766e+00, 1.132812e+00, 1.131514e+00,
121 1.130337e+00, 1.134170e+00, 1.139205e+00, 1.141472e+00, 1.142188e+00, 1.140724e+00,
122 1.140099e+00, 1.139847e+00, 1.137672e+00, 1.138644e+00, 1.136338e+00, 1.136438e+00,
123 1.135945e+00, 1.136429e+00, 1.135701e+00, 1.135702e+00, 1.134112e+00, 1.131934e+00,
124 1.128380e+00, 1.126371e+00, 1.122452e+00, 1.120907e+00, 1.115952e+00, 1.115946e+00,
125 1.114425e+00, 1.111748e+00, 1.106205e+00, 1.107493e+00, 1.103621e+00, 1.102575e+00,
126 1.098815e+00, 1.097888e+00, 1.097305e+00, 1.097129e+00, 1.094577e+00, 1.094551e+00,
127 1.090221e+00, 1.089357e+00, 1.085408e+00, 1.084559e+00, 1.082181e+00, 1.080772e+00,
128 1.079463e+00, 1.078723e+00, 1.076120e+00, 1.075234e+00, 1.073158e+00, 1.071919e+00,
129 1.070394e+00, 1.069502e+00, 1.067524e+00, 1.066918e+00, 1.065778e+00, 1.065318e+00,
130 1.063729e+00, 1.062091e+00, 1.061084e+00, 1.059907e+00, 1.059814e+00, 1.059108e+00,
131 1.051919e+00, 1.051257e+00, 1.049472e+00, 1.048822e+00, 1.045983e+00, 1.046434e+00,
132 1.042613e+00
133 
134 };
135 
136 
138 
139 1.0, 1.0,
140 1.075927e+00, 1.074407e+00, 1.126098e+00, 1.100127e+00, 1.089742e+00, 1.083536e+00,
141 1.089988e+00, 1.103566e+00, 1.096922e+00, 1.126573e+00, 1.132734e+00, 1.136512e+00,
142 1.136629e+00, 1.133086e+00, 1.132428e+00, 1.129299e+00, 1.125622e+00, 1.126992e+00,
143 1.127840e+00, 1.162670e+00, 1.160392e+00, 1.157864e+00, 1.157227e+00, 1.154627e+00,
144 1.192555e+00, 1.197243e+00, 1.197911e+00, 1.200326e+00, 1.220053e+00, 1.215019e+00,
145 1.211703e+00, 1.209080e+00, 1.204248e+00, 1.203328e+00, 1.198671e+00, 1.196840e+00,
146 1.194392e+00, 1.193037e+00, 1.190408e+00, 1.188583e+00, 1.206127e+00, 1.210028e+00,
147 1.206434e+00, 1.204456e+00, 1.200547e+00, 1.199058e+00, 1.200174e+00, 1.200276e+00,
148 1.198912e+00, 1.213048e+00, 1.207160e+00, 1.208020e+00, 1.203814e+00, 1.202380e+00,
149 1.198306e+00, 1.197002e+00, 1.196027e+00, 1.195449e+00, 1.192563e+00, 1.192135e+00,
150 1.187556e+00, 1.186308e+00, 1.182124e+00, 1.180900e+00, 1.178224e+00, 1.176471e+00,
151 1.174811e+00, 1.173702e+00, 1.170827e+00, 1.169581e+00, 1.167205e+00, 1.165626e+00,
152 1.180244e+00, 1.177626e+00, 1.175121e+00, 1.173903e+00, 1.172192e+00, 1.171128e+00,
153 1.168997e+00, 1.166826e+00, 1.164130e+00, 1.165412e+00, 1.165504e+00, 1.165020e+00,
154 1.158462e+00, 1.158014e+00, 1.156519e+00, 1.156081e+00, 1.153602e+00, 1.154190e+00,
155 1.152974e+00
156 
157 };
158 
160 
161 1.0, 1.0,
162 1.140246e+00, 1.097872e+00, 1.104301e+00, 1.068722e+00, 1.044495e+00, 1.062622e+00,
163 1.047987e+00, 1.037032e+00, 1.035686e+00, 1.042870e+00, 1.052222e+00, 1.065100e+00,
164 1.070480e+00, 1.078286e+00, 1.081488e+00, 1.089713e+00, 1.099105e+00, 1.098003e+00,
165 1.102175e+00, 1.117707e+00, 1.121734e+00, 1.125229e+00, 1.126457e+00, 1.128905e+00,
166 1.137312e+00, 1.126263e+00, 1.126459e+00, 1.115191e+00, 1.116986e+00, 1.117184e+00,
167 1.117037e+00, 1.116777e+00, 1.115858e+00, 1.115745e+00, 1.114489e+00, 1.113993e+00,
168 1.113226e+00, 1.112818e+00, 1.111890e+00, 1.111238e+00, 1.111209e+00, 1.111775e+00,
169 1.110256e+00, 1.109414e+00, 1.107647e+00, 1.106980e+00, 1.106096e+00, 1.107331e+00,
170 1.107849e+00, 1.106407e+00, 1.103426e+00, 1.103896e+00, 1.101756e+00, 1.101031e+00,
171 1.098915e+00, 1.098260e+00, 1.097768e+00, 1.097487e+00, 1.095964e+00, 1.095773e+00,
172 1.093348e+00, 1.092687e+00, 1.090465e+00, 1.089821e+00, 1.088394e+00, 1.087462e+00,
173 1.086571e+00, 1.085997e+00, 1.084451e+00, 1.083798e+00, 1.082513e+00, 1.081670e+00,
174 1.080735e+00, 1.075659e+00, 1.074341e+00, 1.073689e+00, 1.072787e+00, 1.072237e+00,
175 1.071107e+00, 1.069955e+00, 1.064856e+00, 1.065873e+00, 1.065938e+00, 1.065694e+00,
176 1.062192e+00, 1.061967e+00, 1.061180e+00, 1.060960e+00, 1.059646e+00, 1.059975e+00,
177 1.059658e+00
178 
179 };
180 
181 
183 
184 1.0, 1.0,
185 1.075927e+00, 1.077959e+00, 1.129145e+00, 1.102088e+00, 1.089765e+00, 1.083542e+00,
186 1.089995e+00, 1.104895e+00, 1.097154e+00, 1.127663e+00, 1.133063e+00, 1.137425e+00,
187 1.136724e+00, 1.133859e+00, 1.132498e+00, 1.130276e+00, 1.127896e+00, 1.127656e+00,
188 1.127905e+00, 1.164210e+00, 1.162259e+00, 1.160075e+00, 1.158978e+00, 1.156649e+00,
189 1.194157e+00, 1.199177e+00, 1.198983e+00, 1.202325e+00, 1.221967e+00, 1.217548e+00,
190 1.214389e+00, 1.211760e+00, 1.207335e+00, 1.206081e+00, 1.201766e+00, 1.199779e+00,
191 1.197283e+00, 1.195706e+00, 1.193071e+00, 1.191115e+00, 1.208838e+00, 1.212681e+00,
192 1.209235e+00, 1.207163e+00, 1.203451e+00, 1.201807e+00, 1.203283e+00, 1.203388e+00,
193 1.202244e+00, 1.216509e+00, 1.211066e+00, 1.211504e+00, 1.207539e+00, 1.205991e+00,
194 1.202143e+00, 1.200724e+00, 1.199595e+00, 1.198815e+00, 1.196025e+00, 1.195390e+00,
195 1.191137e+00, 1.189791e+00, 1.185888e+00, 1.184575e+00, 1.181996e+00, 1.180229e+00,
196 1.178545e+00, 1.177355e+00, 1.174616e+00, 1.173312e+00, 1.171016e+00, 1.169424e+00,
197 1.184120e+00, 1.181478e+00, 1.179085e+00, 1.177817e+00, 1.176124e+00, 1.175003e+00,
198 1.172947e+00, 1.170858e+00, 1.168170e+00, 1.169397e+00, 1.169304e+00, 1.168706e+00,
199 1.162774e+00, 1.162217e+00, 1.160740e+00, 1.160196e+00, 1.157857e+00, 1.158220e+00,
200 1.157267e+00
201 };
202 
203 
205 
206 1.0, 1.0,
207 1.140246e+00, 1.100898e+00, 1.106773e+00, 1.070289e+00, 1.044514e+00, 1.062628e+00,
208 1.047992e+00, 1.038041e+00, 1.035862e+00, 1.043679e+00, 1.052466e+00, 1.065780e+00,
209 1.070551e+00, 1.078869e+00, 1.081541e+00, 1.090455e+00, 1.100847e+00, 1.098511e+00,
210 1.102226e+00, 1.118865e+00, 1.123143e+00, 1.126904e+00, 1.127785e+00, 1.130444e+00,
211 1.138502e+00, 1.127678e+00, 1.127244e+00, 1.116634e+00, 1.118347e+00, 1.118988e+00,
212 1.118957e+00, 1.118696e+00, 1.118074e+00, 1.117722e+00, 1.116717e+00, 1.116111e+00,
213 1.115311e+00, 1.114745e+00, 1.113814e+00, 1.113069e+00, 1.113141e+00, 1.113660e+00,
214 1.112249e+00, 1.111343e+00, 1.109718e+00, 1.108942e+00, 1.108310e+00, 1.109549e+00,
215 1.110227e+00, 1.108846e+00, 1.106183e+00, 1.106354e+00, 1.104388e+00, 1.103583e+00,
216 1.101632e+00, 1.100896e+00, 1.100296e+00, 1.099873e+00, 1.098420e+00, 1.098082e+00,
217 1.095892e+00, 1.095162e+00, 1.093144e+00, 1.092438e+00, 1.091083e+00, 1.090142e+00,
218 1.089236e+00, 1.088604e+00, 1.087159e+00, 1.086465e+00, 1.085239e+00, 1.084388e+00,
219 1.083473e+00, 1.078373e+00, 1.077136e+00, 1.076450e+00, 1.075561e+00, 1.074973e+00,
220 1.073898e+00, 1.072806e+00, 1.067706e+00, 1.068684e+00, 1.068618e+00, 1.068294e+00,
221 1.065241e+00, 1.064939e+00, 1.064166e+00, 1.063872e+00, 1.062659e+00, 1.062828e+00,
222 1.062699e+00
223 
224 };
225 
226 
228 //
229 
231  : G4VCrossSectionDataSet(Default_Name()),
232 // fUpperLimit(100000*GeV),
233  fLowerLimit(10.*MeV),// fLowerLimit(3*GeV),
234  fRadiusConst(1.08*fermi), // 1.1, 1.3 ?
235  fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
236  fDiffractionXsc(0.0)
237 // , fHadronNucleonXsc(0.0)
238 {
267  theA = G4Alpha::Alpha();
268  theHe3 = G4He3::He3();
269 
270  hnXsc = new G4HadronNucleonXsc();
271 }
272 
274 //
275 //
276 
278 {
279  if (hnXsc) delete hnXsc;
280 }
281 
283 
284 G4bool
286  G4int Z, G4int /*A*/,
287  const G4Element*,
288  const G4Material*)
289 {
290  G4bool applicable = false;
291  // G4int baryonNumber = aDP->GetDefinition()->GetBaryonNumber();
292  G4double kineticEnergy = aDP->GetKineticEnergy();
293 
294  const G4ParticleDefinition* theParticle = aDP->GetDefinition();
295 
296  if ( ( kineticEnergy >= fLowerLimit &&
297  Z > 1 && // >= He
298  ( theParticle == theAProton ||
299  theParticle == theGamma ||
300  theParticle == theKPlus ||
301  theParticle == theKMinus ||
302  theParticle == theK0L ||
303  theParticle == theK0S ||
304  theParticle == theSMinus ||
305  theParticle == theProton ||
306  theParticle == theNeutron ||
307  theParticle == thePiPlus ||
308  theParticle == thePiMinus ) ) ) applicable = true;
309 
310  return applicable;
311 }
312 
314 //
315 // Calculates total and inelastic Xsc, derives elastic as total - inelastic accordong to
316 // Glauber model with Gribov correction calculated in the dipole approximation on
317 // light cone. Gaussian density of point-like nucleons helps to calculate rest integrals of the model.
318 // [1] B.Z. Kopeliovich, nucl-th/0306044 + simplification above
319 
320 G4double
322  G4int Z, G4int A,
323  const G4Isotope*,
324  const G4Element*,
325  const G4Material*)
326 {
327  G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
328  G4double hpInXsc(0.), hnInXsc(0.);
329  G4double R = GetNucleusRadius(A);
330 
331  G4int N = A - Z; // number of neutrons
332  if (N < 0) N = 0;
333 
334  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
335 
336  if( theParticle == theProton ||
337  theParticle == theNeutron ||
338  theParticle == thePiPlus ||
339  theParticle == thePiMinus )
340  {
341  // sigma = GetHadronNucleonXscNS(aParticle, A, Z);
342 
343  sigma = Z*hnXsc->GetHadronNucleonXscNS(aParticle, theProton);
344 
345  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
346 
347  sigma += N*hnXsc->GetHadronNucleonXscNS(aParticle, theNeutron);
348 
349  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
350 
351  cofInelastic = 2.4;
352  cofTotal = 2.0;
353  }
354  else if( theParticle == theKPlus ||
355  theParticle == theKMinus ||
356  theParticle == theK0S ||
357  theParticle == theK0L )
358  {
359  sigma = GetKaonNucleonXscVector(aParticle, A, Z);
360  cofInelastic = 2.2;
361  cofTotal = 2.0;
362  R = 1.3*fermi;
363  R *= std::pow(G4double(A), 0.3333);
364  }
365  else
366  {
367  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
368  cofInelastic = 2.2;
369  cofTotal = 2.0;
370  }
371  // cofInelastic = 2.0;
372 
373  if( A > 1 )
374  {
375  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
376  ratio = sigma/nucleusSquare;
377 
378  xsection = nucleusSquare*std::log( 1. + ratio );
379 
380  xsection *= GetParticleBarCorTot(theParticle, Z);
381 
382  fTotalXsc = xsection;
383 
384 
385 
386  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
387 
388  fInelasticXsc *= GetParticleBarCorIn(theParticle, Z);
389 
391 
392  if(fElasticXsc < 0.) fElasticXsc = 0.;
393 
394  G4double difratio = ratio/(1.+ratio);
395 
396  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
397 
398 
399  // sigma = GetHNinelasticXsc(aParticle, A, Z);
400 
401  sigma = Z*hpInXsc + N*hnInXsc;
402 
403  ratio = sigma/nucleusSquare;
404 
405  fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
406 
407  if (fElasticXsc < 0.) fElasticXsc = 0.;
408  }
409  else // H
410  {
411  fTotalXsc = sigma;
412  xsection = sigma;
413 
414  if ( theParticle != theAProton )
415  {
416  sigma = GetHNinelasticXsc(aParticle, A, Z);
417  fInelasticXsc = sigma;
419  }
420  else
421  {
423  }
424  if (fElasticXsc < 0.) fElasticXsc = 0.;
425 
426  }
427  return xsection;
428 }
429 
431 //
432 // Return single-diffraction/inelastic cross-section ratio
433 
435 GetRatioSD(const G4DynamicParticle* aParticle, G4int A, G4int Z)
436 {
437  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
438  G4double R = GetNucleusRadius(A);
439 
440  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
441 
442  if( theParticle == theProton ||
443  theParticle == theNeutron ||
444  theParticle == thePiPlus ||
445  theParticle == thePiMinus )
446  {
447  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
448  cofInelastic = 2.4;
449  cofTotal = 2.0;
450  }
451  else
452  {
453  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
454  cofInelastic = 2.2;
455  cofTotal = 2.0;
456  }
457  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
458  ratio = sigma/nucleusSquare;
459 
460  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
461 
462  G4double difratio = ratio/(1.+ratio);
463 
464  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
465 
466  if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
467  else ratio = 0.;
468 
469  return ratio;
470 }
471 
473 //
474 // Return suasi-elastic/inelastic cross-section ratio
475 
477 GetRatioQE(const G4DynamicParticle* aParticle, G4int A, G4int Z)
478 {
479  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
480  G4double R = GetNucleusRadius(A);
481 
482  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
483 
484  if( theParticle == theProton ||
485  theParticle == theNeutron ||
486  theParticle == thePiPlus ||
487  theParticle == thePiMinus )
488  {
489  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
490  cofInelastic = 2.4;
491  cofTotal = 2.0;
492  }
493  else
494  {
495  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
496  cofInelastic = 2.2;
497  cofTotal = 2.0;
498  }
499  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
500  ratio = sigma/nucleusSquare;
501 
502  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
503 
504  sigma = GetHNinelasticXsc(aParticle, A, Z);
505  ratio = sigma/nucleusSquare;
506 
507  fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
508 
510  else ratio = 0.;
511  if ( ratio < 0. ) ratio = 0.;
512 
513  return ratio;
514 }
515 
517 //
518 // Returns hadron-nucleon Xsc according to differnt parametrisations:
519 // [2] E. Levin, hep-ph/9710546
520 // [3] U. Dersch, et al, hep-ex/9910052
521 // [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
522 
523 G4double
525  const G4Element* anElement)
526 {
527  G4int At = G4lrint(anElement->GetN()); // number of nucleons
528  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
529 
530  return GetHadronNucleonXsc(aParticle, At, Zt);
531 }
532 
534 //
535 // Returns hadron-nucleon Xsc according to differnt parametrisations:
536 // [2] E. Levin, hep-ph/9710546
537 // [3] U. Dersch, et al, hep-ex/9910052
538 // [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
539 
540 G4double
542  G4int At, G4int /*Zt*/)
543 {
544  G4double xsection;
545 
546  //G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
547 
548  G4double targ_mass = 0.939*GeV; // ~mean neutron and proton ???
549 
550  G4double proj_mass = aParticle->GetMass();
551  G4double proj_momentum = aParticle->GetMomentum().mag();
552  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
553 
554  sMand /= GeV*GeV; // in GeV for parametrisation
555  proj_momentum /= GeV;
556 
557  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
558 
559  G4double aa = At;
560 
561  if(theParticle == theGamma)
562  {
563  xsection = aa*(0.0677*std::pow(sMand,0.0808) + 0.129*std::pow(sMand,-0.4525));
564  }
565  else if(theParticle == theNeutron) // as proton ???
566  {
567  xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
568  }
569  else if(theParticle == theProton)
570  {
571  xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
572  // xsection = At*( 49.51*std::pow(sMand,-0.097) + 0.314*std::log(sMand)*std::log(sMand) );
573  // xsection = At*( 38.4 + 0.85*std::abs(std::pow(log(sMand),1.47)) );
574  }
575  else if(theParticle == theAProton)
576  {
577  xsection = aa*( 21.70*std::pow(sMand,0.0808) + 98.39*std::pow(sMand,-0.4525));
578  }
579  else if(theParticle == thePiPlus)
580  {
581  xsection = aa*(13.63*std::pow(sMand,0.0808) + 27.56*std::pow(sMand,-0.4525));
582  }
583  else if(theParticle == thePiMinus)
584  {
585  // xsection = At*( 55.2*std::pow(sMand,-0.255) + 0.346*std::log(sMand)*std::log(sMand) );
586  xsection = aa*(13.63*std::pow(sMand,0.0808) + 36.02*std::pow(sMand,-0.4525));
587  }
588  else if(theParticle == theKPlus)
589  {
590  xsection = aa*(11.82*std::pow(sMand,0.0808) + 8.15*std::pow(sMand,-0.4525));
591  }
592  else if(theParticle == theKMinus)
593  {
594  xsection = aa*(11.82*std::pow(sMand,0.0808) + 26.36*std::pow(sMand,-0.4525));
595  }
596  else // as proton ???
597  {
598  xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
599  }
600  xsection *= millibarn;
601  return xsection;
602 }
603 
604 
606 //
607 // Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
608 // http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
609 
610 G4double
612  const G4Element* anElement)
613 {
614  G4int At = G4lrint(anElement->GetN()); // number of nucleons
615  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
616 
617  return GetHadronNucleonXscPDG(aParticle, At, Zt);
618 }
619 
620 
621 
622 
624 //
625 // Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
626 // http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
627 // At = number of nucleons, Zt = number of protons
628 
629 G4double
631  G4int At, G4int Zt)
632 {
633  G4double xsection;
634 
635  G4int Nt = At-Zt; // number of neutrons
636  if (Nt < 0) Nt = 0;
637 
638  G4double zz = Zt;
639  G4double aa = At;
640  G4double nn = Nt;
641 
643  GetIonTable()->GetIonMass(Zt, At);
644 
645  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
646 
647  G4double proj_mass = aParticle->GetMass();
648  G4double proj_momentum = aParticle->GetMomentum().mag();
649 
650  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
651 
652  sMand /= GeV*GeV; // in GeV for parametrisation
653 
654  // General PDG fit constants
655 
656  G4double s0 = 5.38*5.38; // in Gev^2
657  G4double eta1 = 0.458;
658  G4double eta2 = 0.458;
659  G4double B = 0.308;
660 
661 
662  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
663 
664 
665  if(theParticle == theNeutron) // proton-neutron fit
666  {
667  xsection = zz*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
668  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
669  xsection += nn*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
670  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
671  }
672  else if(theParticle == theProton)
673  {
674 
675  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
676  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
677 
678  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
679  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
680  }
681  else if(theParticle == theAProton)
682  {
683  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
684  + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
685 
686  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
687  + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
688  }
689  else if(theParticle == thePiPlus)
690  {
691  xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
692  + 19.24*std::pow(sMand,-eta1) - 6.03*std::pow(sMand,-eta2));
693  }
694  else if(theParticle == thePiMinus)
695  {
696  xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
697  + 19.24*std::pow(sMand,-eta1) + 6.03*std::pow(sMand,-eta2));
698  }
699  else if(theParticle == theKPlus || theParticle == theK0L )
700  {
701  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
702  + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
703 
704  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
705  + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
706  }
707  else if(theParticle == theKMinus || theParticle == theK0S )
708  {
709  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
710  + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
711 
712  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
713  + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
714  }
715  else if(theParticle == theSMinus)
716  {
717  xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
718  - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
719  }
720  else if(theParticle == theGamma) // modify later on
721  {
722  xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
723  + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
724 
725  }
726  else // as proton ???
727  {
728  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
729  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
730 
731  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
732  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
733  }
734  xsection *= millibarn; // parametrised in mb
735  return xsection;
736 }
737 
738 
740 //
741 // Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
742 // data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
743 
744 G4double
746  const G4Element* anElement)
747 {
748  G4int At = G4lrint(anElement->GetN()); // number of nucleons
749  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
750 
751  return GetHadronNucleonXscNS(aParticle, At, Zt);
752 }
753 
754 
755 
756 
758 //
759 // Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
760 // data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
761 
762 G4double
764  G4int At, G4int Zt)
765 {
766  G4double xsection(0);
767  // G4double Delta; DHW 19 May 2011: variable set but not used
768  G4double A0, B0;
769  G4double hpXscv(0);
770  G4double hnXscv(0);
771 
772  G4int Nt = At-Zt; // number of neutrons
773  if (Nt < 0) Nt = 0;
774 
775  G4double aa = At;
776  G4double zz = Zt;
777  G4double nn = Nt;
778 
780  GetIonTable()->GetIonMass(Zt, At);
781 
782  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
783 
784  G4double proj_mass = aParticle->GetMass();
785  G4double proj_energy = aParticle->GetTotalEnergy();
786  G4double proj_momentum = aParticle->GetMomentum().mag();
787 
788  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
789 
790  sMand /= GeV*GeV; // in GeV for parametrisation
791  proj_momentum /= GeV;
792  proj_energy /= GeV;
793  proj_mass /= GeV;
794 
795  // General PDG fit constants
796 
797  G4double s0 = 5.38*5.38; // in Gev^2
798  G4double eta1 = 0.458;
799  G4double eta2 = 0.458;
800  G4double B = 0.308;
801 
802 
803  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
804 
805 
806  if(theParticle == theNeutron)
807  {
808  if( proj_momentum >= 373.)
809  {
810  return GetHadronNucleonXscPDG(aParticle,At,Zt);
811  }
812  else if( proj_momentum >= 10.)
813  // if( proj_momentum >= 2.)
814  {
815  // Delta = 1.; // DHW 19 May 2011: variable set but not used
816  // if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
817 
818  if(proj_momentum >= 10.)
819  {
820  B0 = 7.5;
821  A0 = 100. - B0*std::log(3.0e7);
822 
823  xsection = A0 + B0*std::log(proj_energy) - 11
824  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
825  0.93827*0.93827,-0.165); // mb
826  }
827  xsection *= zz + nn;
828  }
829  else
830  {
831  // nn to be pp
832 
833  if( proj_momentum < 0.73 )
834  {
835  hnXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
836  }
837  else if( proj_momentum < 1.05 )
838  {
839  hnXscv = 23 + 40*(std::log(proj_momentum/0.73))*
840  (std::log(proj_momentum/0.73));
841  }
842  else // if( proj_momentum < 10. )
843  {
844  hnXscv = 39.0+
845  75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
846  }
847  // pn to be np
848 
849  if( proj_momentum < 0.8 )
850  {
851  hpXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
852  }
853  else if( proj_momentum < 1.4 )
854  {
855  hpXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
856  }
857  else // if( proj_momentum < 10. )
858  {
859  hpXscv = 33.3+
860  20.8*(std::pow(proj_momentum,2.0)-1.35)/
861  (std::pow(proj_momentum,2.50)+0.95);
862  }
863  xsection = hpXscv*zz + hnXscv*nn;
864  }
865  }
866  else if(theParticle == theProton)
867  {
868  if( proj_momentum >= 373.)
869  {
870  return GetHadronNucleonXscPDG(aParticle,At,Zt);
871  }
872  else if( proj_momentum >= 10.)
873  // if( proj_momentum >= 2.)
874  {
875  // Delta = 1.; DHW 19 May 2011: variable set but not used
876  // if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
877 
878  if(proj_momentum >= 10.)
879  {
880  B0 = 7.5;
881  A0 = 100. - B0*std::log(3.0e7);
882 
883  xsection = A0 + B0*std::log(proj_energy) - 11
884  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
885  0.93827*0.93827,-0.165); // mb
886  }
887  xsection *= zz + nn;
888  }
889  else
890  {
891  // pp
892 
893  if( proj_momentum < 0.73 )
894  {
895  hpXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
896  }
897  else if( proj_momentum < 1.05 )
898  {
899  hpXscv = 23 + 40*(std::log(proj_momentum/0.73))*
900  (std::log(proj_momentum/0.73));
901  }
902  else // if( proj_momentum < 10. )
903  {
904  hpXscv = 39.0+
905  75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
906  }
907  // pn to be np
908 
909  if( proj_momentum < 0.8 )
910  {
911  hnXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
912  }
913  else if( proj_momentum < 1.4 )
914  {
915  hnXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
916  }
917  else // if( proj_momentum < 10. )
918  {
919  hnXscv = 33.3+
920  20.8*(std::pow(proj_momentum,2.0)-1.35)/
921  (std::pow(proj_momentum,2.50)+0.95);
922  }
923  xsection = hpXscv*zz + hnXscv*nn;
924  // xsection = hpXscv*(Zt + Nt);
925  // xsection = hnXscv*(Zt + Nt);
926  }
927  // xsection *= 0.95;
928  }
929  else if( theParticle == theAProton )
930  {
931  // xsection = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
932  // + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
933 
934  // xsection += Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
935  // + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
936 
937  G4double logP = std::log(proj_momentum);
938 
939  if( proj_momentum <= 1.0 )
940  {
941  xsection = zz*(65.55 + 53.84/(proj_momentum+1.e-6) );
942  }
943  else
944  {
945  xsection = zz*( 41.1 + 77.2*std::pow( proj_momentum, -0.68)
946  + 0.293*logP*logP - 1.82*logP );
947  }
948  if ( nn > 0.)
949  {
950  xsection += nn*( 41.9 + 96.2*std::pow( proj_momentum, -0.99) - 0.154*logP);
951  }
952  else // H
953  {
954  fInelasticXsc = 38.0 + 38.0*std::pow( proj_momentum, -0.96)
955  - 0.169*logP*logP;
957  }
958  }
959  else if( theParticle == thePiPlus )
960  {
961  if(proj_momentum < 0.4)
962  {
963  G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
964  hpXscv = Ex3+20.0;
965  }
966  else if( proj_momentum < 1.15 )
967  {
968  G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
969  hpXscv = Ex4+14.0;
970  }
971  else if(proj_momentum < 3.5)
972  {
973  G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
974  G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
975  hpXscv = Ex1+Ex2+27.5;
976  }
977  else // if(proj_momentum > 3.5) // mb
978  {
979  hpXscv = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
980  }
981  // pi+n = pi-p??
982 
983  if(proj_momentum < 0.37)
984  {
985  hnXscv = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
986  }
987  else if(proj_momentum<0.65)
988  {
989  hnXscv = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
990  }
991  else if(proj_momentum<1.3)
992  {
993  hnXscv = 36.1+
994  10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
995  24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
996  }
997  else if(proj_momentum<3.0)
998  {
999  hnXscv = 36.1+0.079-4.313*std::log(proj_momentum)+
1000  3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
1001  1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
1002  }
1003  else // mb
1004  {
1005  hnXscv = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
1006  }
1007  xsection = hpXscv*zz + hnXscv*nn;
1008  }
1009  else if(theParticle == thePiMinus)
1010  {
1011  // pi-n = pi+p??
1012 
1013  if(proj_momentum < 0.4)
1014  {
1015  G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
1016  hnXscv = Ex3+20.0;
1017  }
1018  else if(proj_momentum < 1.15)
1019  {
1020  G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
1021  hnXscv = Ex4+14.0;
1022  }
1023  else if(proj_momentum < 3.5)
1024  {
1025  G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
1026  G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
1027  hnXscv = Ex1+Ex2+27.5;
1028  }
1029  else // if(proj_momentum > 3.5) // mb
1030  {
1031  hnXscv = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
1032  }
1033  // pi-p
1034 
1035  if(proj_momentum < 0.37)
1036  {
1037  hpXscv = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
1038  }
1039  else if(proj_momentum<0.65)
1040  {
1041  hpXscv = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
1042  }
1043  else if(proj_momentum<1.3)
1044  {
1045  hpXscv = 36.1+
1046  10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
1047  24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
1048  }
1049  else if(proj_momentum<3.0)
1050  {
1051  hpXscv = 36.1+0.079-4.313*std::log(proj_momentum)+
1052  3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
1053  1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
1054  }
1055  else // mb
1056  {
1057  hpXscv = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
1058  }
1059  xsection = hpXscv*zz + hnXscv*nn;
1060  }
1061  else if(theParticle == theKPlus)
1062  {
1063  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
1064  + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
1065 
1066  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
1067  + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
1068  }
1069  else if(theParticle == theKMinus)
1070  {
1071  xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
1072  + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
1073 
1074  xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
1075  + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
1076  }
1077  else if(theParticle == theSMinus)
1078  {
1079  xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
1080  - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
1081  }
1082  else if(theParticle == theGamma) // modify later on
1083  {
1084  xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
1085  + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
1086 
1087  }
1088  else // as proton ???
1089  {
1090  xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
1091  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
1092 
1093  xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
1094  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
1095  }
1096  xsection *= millibarn; // parametrised in mb
1097  return xsection;
1098 }
1099 
1100 G4double
1102  G4int At, G4int Zt)
1103 {
1104  G4double Tkin, logTkin, xsc, xscP, xscN;
1105  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
1106 
1107  G4int Nt = At-Zt; // number of neutrons
1108  if (Nt < 0) Nt = 0;
1109 
1110  Tkin = aParticle->GetKineticEnergy(); // Tkin in MeV
1111 
1112  if( Tkin > 70*GeV ) return GetHadronNucleonXscPDG(aParticle,At,Zt);
1113 
1114  logTkin = std::log(Tkin); // Tkin in MeV!!!
1115 
1116  if( theParticle == theKPlus )
1117  {
1118  xscP = hnXsc->GetKpProtonTotXscVector(logTkin);
1119  xscN = hnXsc->GetKpNeutronTotXscVector(logTkin);
1120  }
1121  else if( theParticle == theKMinus )
1122  {
1123  xscP = hnXsc->GetKmProtonTotXscVector(logTkin);
1124  xscN = hnXsc->GetKmNeutronTotXscVector(logTkin);
1125  }
1126  else // K-zero as half of K+ and K-
1127  {
1128  xscP = (hnXsc->GetKpProtonTotXscVector(logTkin)+hnXsc->GetKmProtonTotXscVector(logTkin))*0.5;
1129  xscN = (hnXsc->GetKpNeutronTotXscVector(logTkin)+hnXsc->GetKmNeutronTotXscVector(logTkin))*0.5;
1130  }
1131  xsc = xscP*Zt + xscN*Nt;
1132  return xsc;
1133 }
1135 //
1136 // Returns hadron-nucleon inelastic cross-section based on proper parametrisation
1137 
1138 G4double
1140  const G4Element* anElement)
1141 {
1142  G4int At = G4lrint(anElement->GetN()); // number of nucleons
1143  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
1144 
1145  return GetHNinelasticXsc(aParticle, At, Zt);
1146 }
1147 
1149 //
1150 // Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
1151 
1152 G4double
1154  G4int At, G4int Zt)
1155 {
1156  G4ParticleDefinition* hadron = aParticle->GetDefinition();
1157  G4double sumInelastic;
1158  G4int Nt = At - Zt;
1159  if(Nt < 0) Nt = 0;
1160 
1161  if( hadron == theKPlus )
1162  {
1163  sumInelastic = GetHNinelasticXscVU(aParticle, At, Zt);
1164  }
1165  else
1166  {
1167  //sumInelastic = Zt*GetHadronNucleonXscMK(aParticle, theProton);
1168  // sumInelastic += Nt*GetHadronNucleonXscMK(aParticle, theNeutron);
1169  sumInelastic = G4double(Zt)*GetHadronNucleonXscNS(aParticle, 1, 1);
1170  sumInelastic += G4double(Nt)*GetHadronNucleonXscNS(aParticle, 1, 0);
1171  }
1172  return sumInelastic;
1173 }
1174 
1175 
1177 //
1178 // Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
1179 
1180 G4double
1182  G4int At, G4int Zt)
1183 {
1184  G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
1185  G4int absPDGcode = std::abs(PDGcode);
1186 
1187  G4double Elab = aParticle->GetTotalEnergy();
1188  // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
1189  G4double Plab = aParticle->GetMomentum().mag();
1190  // std::sqrt(Elab * Elab - 0.88);
1191 
1192  Elab /= GeV;
1193  Plab /= GeV;
1194 
1195  G4double LogPlab = std::log( Plab );
1196  G4double sqrLogPlab = LogPlab * LogPlab;
1197 
1198  //G4cout<<"Plab = "<<Plab<<G4endl;
1199 
1200  G4double NumberOfTargetProtons = G4double(Zt);
1201  G4double NumberOfTargetNucleons = G4double(At);
1202  G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
1203 
1204  if(NumberOfTargetNeutrons < 0.0) NumberOfTargetNeutrons = 0.0;
1205 
1206  G4double Xtotal, Xelastic, Xinelastic;
1207 
1208  if( absPDGcode > 1000 ) //------Projectile is baryon --------
1209  {
1210  G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
1211  0.522*sqrLogPlab - 4.51*LogPlab;
1212 
1213  G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
1214  0.513*sqrLogPlab - 4.27*LogPlab;
1215 
1216  G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
1217  0.169*sqrLogPlab - 1.85*LogPlab;
1218 
1219  G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
1220  0.169*sqrLogPlab - 1.85*LogPlab;
1221 
1222  Xtotal = (NumberOfTargetProtons * XtotPP +
1223  NumberOfTargetNeutrons * XtotPN);
1224 
1225  Xelastic = (NumberOfTargetProtons * XelPP +
1226  NumberOfTargetNeutrons * XelPN);
1227  }
1228  else if( PDGcode == 211 ) //------Projectile is PionPlus -------
1229  {
1230  G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) +
1231  0.19 *sqrLogPlab - 0.0 *LogPlab;
1232 
1233  G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) +
1234  0.456*sqrLogPlab - 4.03*LogPlab;
1235 
1236  G4double XelPiP = 0.0 + 11.4*std::pow(Plab,-0.40) +
1237  0.079*sqrLogPlab - 0.0 *LogPlab;
1238 
1239  G4double XelPiN = 1.76 + 11.2*std::pow(Plab,-0.64) +
1240  0.043*sqrLogPlab - 0.0 *LogPlab;
1241 
1242  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1243  NumberOfTargetNeutrons * XtotPiN );
1244 
1245  Xelastic = ( NumberOfTargetProtons * XelPiP +
1246  NumberOfTargetNeutrons * XelPiN );
1247  }
1248  else if( PDGcode == -211 ) //------Projectile is PionMinus -------
1249  {
1250  G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) +
1251  0.456*sqrLogPlab - 4.03*LogPlab;
1252 
1253  G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) +
1254  0.19 *sqrLogPlab - 0.0 *LogPlab;
1255 
1256  G4double XelPiP = 1.76 + 11.2*std::pow(Plab,-0.64) +
1257  0.043*sqrLogPlab - 0.0 *LogPlab;
1258 
1259  G4double XelPiN = 0.0 + 11.4*std::pow(Plab,-0.40) +
1260  0.079*sqrLogPlab - 0.0 *LogPlab;
1261 
1262  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1263  NumberOfTargetNeutrons * XtotPiN );
1264 
1265  Xelastic = ( NumberOfTargetProtons * XelPiP +
1266  NumberOfTargetNeutrons * XelPiN );
1267  }
1268  else if( PDGcode == 111 ) //------Projectile is PionZero -------
1269  {
1270  G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) +
1271  0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
1272  33.0 + 14.0 *std::pow(Plab,-1.36) +
1273  0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
1274 
1275  G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) +
1276  0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
1277  16.4 + 19.3 *std::pow(Plab,-0.42) +
1278  0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1279 
1280  G4double XelPiP =( 0.0 + 11.4*std::pow(Plab,-0.40) +
1281  0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
1282  1.76 + 11.2*std::pow(Plab,-0.64) +
1283  0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1284 
1285  G4double XelPiN =( 1.76 + 11.2*std::pow(Plab,-0.64) +
1286  0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
1287  0.0 + 11.4*std::pow(Plab,-0.40) +
1288  0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1289 
1290  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1291  NumberOfTargetNeutrons * XtotPiN );
1292 
1293  Xelastic = ( NumberOfTargetProtons * XelPiP +
1294  NumberOfTargetNeutrons * XelPiN );
1295  }
1296  else if( PDGcode == 321 ) //------Projectile is KaonPlus -------
1297  {
1298  G4double XtotKP = 18.1 + 0. *std::pow(Plab, 0. ) +
1299  0.26 *sqrLogPlab - 1.0 *LogPlab;
1300  G4double XtotKN = 18.7 + 0. *std::pow(Plab, 0. ) +
1301  0.21 *sqrLogPlab - 0.89*LogPlab;
1302 
1303  G4double XelKP = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1304  0.16 *sqrLogPlab - 1.3 *LogPlab;
1305 
1306  G4double XelKN = 7.3 + 0. *std::pow(Plab,-0. ) +
1307  0.29 *sqrLogPlab - 2.4 *LogPlab;
1308 
1309  Xtotal = ( NumberOfTargetProtons * XtotKP +
1310  NumberOfTargetNeutrons * XtotKN );
1311 
1312  Xelastic = ( NumberOfTargetProtons * XelKP +
1313  NumberOfTargetNeutrons * XelKN );
1314  }
1315  else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------
1316  {
1317  G4double XtotKP = 32.1 + 0. *std::pow(Plab, 0. ) +
1318  0.66 *sqrLogPlab - 5.6 *LogPlab;
1319  G4double XtotKN = 25.2 + 0. *std::pow(Plab, 0. ) +
1320  0.38 *sqrLogPlab - 2.9 *LogPlab;
1321 
1322  G4double XelKP = 7.3 + 0. *std::pow(Plab,-0. ) +
1323  0.29 *sqrLogPlab - 2.4 *LogPlab;
1324 
1325  G4double XelKN = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1326  0.16 *sqrLogPlab - 1.3 *LogPlab;
1327 
1328  Xtotal = ( NumberOfTargetProtons * XtotKP +
1329  NumberOfTargetNeutrons * XtotKN );
1330 
1331  Xelastic = ( NumberOfTargetProtons * XelKP +
1332  NumberOfTargetNeutrons * XelKN );
1333  }
1334  else if( PDGcode == 311 ) //------Projectile is KaonZero ------
1335  {
1336  G4double XtotKP = ( 18.1 + 0. *std::pow(Plab, 0. ) +
1337  0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
1338  32.1 + 0. *std::pow(Plab, 0. ) +
1339  0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
1340 
1341  G4double XtotKN = ( 18.7 + 0. *std::pow(Plab, 0. ) +
1342  0.21 *sqrLogPlab - 0.89*LogPlab + //K+
1343  25.2 + 0. *std::pow(Plab, 0. ) +
1344  0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
1345 
1346  G4double XelKP = ( 5.0 + 8.1*std::pow(Plab,-1.8 )
1347  + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
1348  7.3 + 0. *std::pow(Plab,-0. ) +
1349  0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
1350 
1351  G4double XelKN = ( 7.3 + 0. *std::pow(Plab,-0. ) +
1352  0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
1353  5.0 + 8.1*std::pow(Plab,-1.8 ) +
1354  0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
1355 
1356  Xtotal = ( NumberOfTargetProtons * XtotKP +
1357  NumberOfTargetNeutrons * XtotKN );
1358 
1359  Xelastic = ( NumberOfTargetProtons * XelKP +
1360  NumberOfTargetNeutrons * XelKN );
1361  }
1362  else //------Projectile is undefined, Nucleon assumed
1363  {
1364  G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
1365  0.522*sqrLogPlab - 4.51*LogPlab;
1366 
1367  G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
1368  0.513*sqrLogPlab - 4.27*LogPlab;
1369 
1370  G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
1371  0.169*sqrLogPlab - 1.85*LogPlab;
1372  G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
1373  0.169*sqrLogPlab - 1.85*LogPlab;
1374 
1375  Xtotal = ( NumberOfTargetProtons * XtotPP +
1376  NumberOfTargetNeutrons * XtotPN );
1377 
1378  Xelastic = ( NumberOfTargetProtons * XelPP +
1379  NumberOfTargetNeutrons * XelPN );
1380  }
1381  Xinelastic = Xtotal - Xelastic;
1382 
1383  if( Xinelastic < 0.) Xinelastic = 0.;
1384 
1385  return Xinelastic*= millibarn;
1386 }
1387 
1389 //
1390 //
1391 
1392 G4double
1394  const G4Element* anElement)
1395 {
1396  G4int At = G4lrint(anElement->GetN());
1397  G4double oneThird = 1.0/3.0;
1398  G4double cubicrAt = std::pow(G4double(At), oneThird);
1399 
1400  G4double R; // = fRadiusConst*cubicrAt;
1401  /*
1402  G4double tmp = std::pow( cubicrAt-1., 3.);
1403  tmp += At;
1404  tmp *= 0.5;
1405 
1406  if (At > 20.) // 20.
1407  {
1408  R = fRadiusConst*std::pow (tmp, oneThird);
1409  }
1410  else
1411  {
1412  R = fRadiusConst*cubicrAt;
1413  }
1414  */
1415 
1416  R = fRadiusConst*cubicrAt;
1417 
1418  G4double meanA = 21.;
1419 
1420  G4double tauA1 = 40.;
1421  G4double tauA2 = 10.;
1422  G4double tauA3 = 5.;
1423 
1424  G4double a1 = 0.85;
1425  G4double b1 = 1. - a1;
1426 
1427  G4double b2 = 0.3;
1428  G4double b3 = 4.;
1429 
1430  if (At > 20) // 20.
1431  {
1432  R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
1433  }
1434  else if (At > 3)
1435  {
1436  R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
1437  }
1438  else
1439  {
1440  R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
1441  }
1442  return R;
1443 
1444 }
1446 //
1447 //
1448 
1449 G4double
1451 {
1452  G4double oneThird = 1.0/3.0;
1453  G4double cubicrAt = std::pow(G4double(At), oneThird);
1454 
1455  G4double R; // = fRadiusConst*cubicrAt;
1456 
1457  /*
1458  G4double tmp = std::pow( cubicrAt-1., 3.);
1459  tmp += At;
1460  tmp *= 0.5;
1461 
1462  if (At > 20.)
1463  {
1464  R = fRadiusConst*std::pow (tmp, oneThird);
1465  }
1466  else
1467  {
1468  R = fRadiusConst*cubicrAt;
1469  }
1470  */
1471 
1472  R = fRadiusConst*cubicrAt;
1473 
1474  G4double meanA = 20.;
1475  G4double tauA = 20.;
1476 
1477  if (At > 20) // 20.
1478  {
1479  R *= ( 0.8 + 0.2*std::exp( -(G4double(At) - meanA)/tauA) );
1480  }
1481  else
1482  {
1483  R *= ( 1.0 + 0.1*( 1. - std::exp( (G4double(At) - meanA)/tauA) ) );
1484  }
1485 
1486  return R;
1487 }
1488 
1490 //
1491 //
1492 
1494  const G4double mt ,
1495  const G4double Plab )
1496 {
1497  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1498  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
1499  // G4double Pcm = Plab * mt / Ecm;
1500  // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
1501 
1502  return Ecm ; // KEcm;
1503 }
1504 
1506 //
1507 //
1508 
1510  const G4double mt ,
1511  const G4double Plab )
1512 {
1513  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1514  G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
1515 
1516  return sMand;
1517 }
1518 
1520 //
1521 //
1522 
1524 {
1525  outFile << "G4GlauberGribovCrossSection calculates total, inelastic and\n"
1526  << "elastic cross sections for hadron-nucleus cross sections using\n"
1527  << "the Glauber model with Gribov corrections. It is valid for all\n"
1528  << "targets except hydrogen, and for incident p, pbar, n, sigma-,\n"
1529  << "pi+, pi-, K+, K- and gammas with energies above 3 GeV. This is\n"
1530  << "a cross section component which is to be used to build a cross\n"
1531  << "data set.\n";
1532 }
1533 
1534 //
1535 //
G4double GetHNinelasticXscVU(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetRatioSD(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetKmNeutronTotXscVector(G4double logEnergy)
static const double MeV
Definition: G4SIunits.hh:193
G4_DECLARE_XS_FACTORY(G4GlauberGribovCrossSection)
G4double GetKpProtonTotXscVector(G4double logEnergy)
static G4AntiOmegaMinus * AntiOmegaMinus()
G4double GetKineticEnergy() const
G4double GetParticleBarCorTot(const G4ParticleDefinition *theParticle, G4int Z)
G4double GetTotalEnergy() const
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double GetN() const
Definition: G4Element.hh:134
static const G4double fNeutronBarCorrectionTot[93]
static const G4double a1
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
static G4OmegaMinus * OmegaMinus()
const G4double pi
G4double GetZ() const
Definition: G4Element.hh:131
static G4KaonZeroLong * KaonZeroLong()
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
static const G4double fPionMinusBarCorrectionIn[93]
static G4AntiSigmaPlus * AntiSigmaPlus()
int G4int
Definition: G4Types.hh:78
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)
static const G4double fPionMinusBarCorrectionTot[93]
static const G4double b3
static const G4double fPionPlusBarCorrectionTot[93]
static G4KaonZeroShort * KaonZeroShort()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
G4double GetMass() const
G4double GetKaonNucleonXscVector(const G4DynamicParticle *, G4int At, G4int Zt)
bool G4bool
Definition: G4Types.hh:79
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4AntiXiMinus * AntiXiMinus()
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static const double GeV
Definition: G4SIunits.hh:196
static const G4double b2
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
static const G4double fProtonBarCorrectionTot[93]
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static const G4double A[nN]
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
static G4SigmaMinus * SigmaMinus()
static G4ParticleTable * GetParticleTable()
static G4AntiLambda * AntiLambda()
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4AntiSigmaZero * AntiSigmaZero()
G4double GetKmProtonTotXscVector(G4double logEnergy)
static const double millibarn
Definition: G4SIunits.hh:96
static const G4double b1
static const G4double fNeutronBarCorrectionIn[93]
static G4AntiXiZero * AntiXiZero()
G4double CalculateEcmValue(const G4double, const G4double, const G4double)
virtual void CrossSectionDescription(std::ostream &) const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
G4double GetRatioQE(const G4DynamicParticle *, G4int At, G4int Zt)
static const G4double fPionPlusBarCorrectionIn[93]
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4He3 * He3()
Definition: G4He3.cc:94
const G4double oneThird
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *aDP, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
static const double fermi
Definition: G4SIunits.hh:93
static G4AntiNeutron * AntiNeutron()
G4double GetInelasticHadronNucleonXsc()
G4ThreeVector GetMomentum() const
G4double GetKpNeutronTotXscVector(G4double logEnergy)
static const G4double fProtonBarCorrectionIn[93]