source: trunk/source/processes/hadronic/cross_sections/src/G4HadronCrossSections.cc @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 83.1 KB
Line 
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// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29//
30// G4 Hadron Physics class G4HadronCrossSections
31// F.W. Jones, TRIUMF, 03-DEC-96
32//
33// This class encapsulates cross section data and calculations
34// from the Geant3/Gheisha routine GHESIG.
35// The overloaded method MakePhysicsVector can be used to Generate
36// Physics Tables for different processes.
37//
38// Note: this is implemented as a SINGLETON class
39//
40// 27-MAR-97 FWJ: first version for Alpha release
41// 14-APR-97 FWJ: class name changed from G4LCrossSectionData
42//    to G4HadronicCrossSections
43// 23-Apr-97 Johannes Peter Wellisch: debugging, add the residual particles
44// 23-MAY-97 FWJ: corrected problem with neutron cross sections
45// 20-JUN-97 FWJ: added some missing elastic data (e.g. K+) and
46//    fixed a momentum/energy units problem in the physics vectors
47// 14-APR-98 FWJ: rewritten as class G4HadronCrossSections
48//    and adapted to G4CrossSectionDataSet/DataStore class design.
49// 25-JUN-98 FWJ: optimised bin selection.
50// 06-NOV-98 FWJ: added first-order correction for low-energy
51//    inelastic cross sections
52//
53
54
55#include "G4HadronCrossSections.hh"
56#include "G4ios.hh"
57 
58
59// Initialize static pointer for singleton instance
60G4HadronCrossSections* 
61G4HadronCrossSections::theInstance = 0;
62
63
64// Cross section tables from G3.21/GHEISHA routine GHESIG
65
66//---------------------------------------------------------------------
67// Lab Momentum in GeV/c
68//---------------------------------------------------------------------
69G4float G4HadronCrossSections::plab[TSIZE] = {
70       0.00000E+00    , 0.10000    , 0.15000    , 0.20000    , 0.25000    ,
71       0.30000        , 0.35000    , 0.40000    , 0.45000    , 0.50000    ,
72       0.55000        , 0.60000    , 0.65000    , 0.70000    , 0.75000    ,
73       0.80000        , 0.85000    , 0.90000    , 0.95000    ,  1.0000    ,
74        1.1000        ,  1.2000    ,  1.3000    ,  1.4000    ,  1.5000    ,
75        1.6000        ,  1.8000    ,  2.0000    ,  2.2000    ,  2.4000    ,
76        2.6000        ,  2.8000    ,  3.0000    ,  4.0000    ,  5.0000    ,
77        6.0000        ,  8.0000    ,  10.000    ,  20.000    ,  100.00    ,
78        1000.0       
79};
80
81//---------------------------------------------------------------------
82// Elastic scattering on free protons
83//---------------------------------------------------------------------
84
85G4float G4HadronCrossSections::csel[NPARTS][TSIZE] = {
86      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //1
87       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
88       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
89       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
90       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
91       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
92       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
93       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
94       0.00000E+00},
95
96      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //2
97       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
98       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
99       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
100       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
101       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
102       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
103       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
104       0.00000E+00},
105
106      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //3
107       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
108       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
109       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
110       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
111       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
112       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
113       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
114       0.00000E+00},
115
116      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //4
117       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
118       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
119       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
120       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
121       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
122       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
123       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
124       0.00000E+00},
125
126      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //5
127       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
128       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
129       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
130       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
131       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
132       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
133       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
134       0.00000E+00},
135
136      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //6
137       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
138       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
139       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
140       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
141       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
142       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
143       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
144       0.00000E+00},
145
146// Elastic cross section for piplus - p
147      {0.00000E+00,  6.0000    ,  20.000    ,  71.000    ,  155.00    , //7
148        195.00    ,  130.00    ,  78.000    ,  60.000    ,  32.000    ,
149        23.500    ,  18.500    ,  15.000    ,  12.500    ,  10.000    ,
150        9.1000    ,  8.6000    ,  8.8000    ,  9.5000    ,  10.600    ,
151        13.000    ,  15.500    ,  17.100    ,  17.200    ,  16.200    ,
152        15.000    ,  12.300    ,  10.200    ,  9.0000    ,  8.0000    ,
153        7.3000    ,  6.8000    ,  6.5000    ,  5.8000    ,  5.4000    ,
154        5.2000    ,  5.0000    ,  4.9000    ,  3.8000    ,  3.2000    ,
155        3.5000},
156
157      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //8
158       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
159       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
160       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
161       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
162       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
163       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
164       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
165       0.00000E+00},
166
167// Elastic cross section for piminus - p
168      {0.00000E+00,  1.0000    ,  3.0000    ,  8.0000    ,  18.000    , //9
169        25.000    ,  27.500    ,  12.300    ,  10.600    ,  11.000    ,
170        12.500    ,  14.500    ,  17.000    ,  19.400    ,  19.800    ,
171        16.800    ,  14.000    ,  14.800    ,  20.000    ,  26.100    ,
172        19.500    ,  15.000    ,  12.800    ,  11.500    ,  10.500    ,
173        9.8000    ,  8.8000    ,  8.2000    ,  7.8000    ,  7.5000    ,
174        7.2000    ,  7.0000    ,  6.8000    ,  6.1000    ,  5.7000    ,
175        5.4000    ,  4.9000    ,  4.6000    ,  4.0000    ,  3.3000    ,
176        3.5000},
177
178      {10.000    ,  11.200    ,  11.300    ,  11.400    ,  11.500    , //10
179        11.600    ,  11.800    ,  12.000    ,  12.100    ,  12.200    ,
180        12.300    ,  12.400    ,  12.500    ,  12.500    ,  12.500    ,
181        12.400    ,  12.300    ,  12.200    ,  12.000    ,  11.800    ,
182        11.200    ,  11.500    ,  9.9000    ,  9.4000    ,  8.8000    ,
183        8.4000    ,  7.5000    ,  6.9000    ,  6.3000    ,  5.9000    ,
184        5.5000    ,  5.2000    ,  5.0000    ,  4.0000    ,  3.5000    ,
185        3.3000    ,  3.1000    ,  3.1000    ,  3.0000    ,  2.5000    ,
186        3.0000},
187
188      {10.000    ,  11.200    ,  11.300    ,  11.400    ,  11.500    , //11
189        11.600    ,  11.800    ,  12.000    ,  12.100    ,  12.200    ,
190        12.300    ,  12.400    ,  12.500    ,  12.500    ,  12.500    ,
191        12.400    ,  12.300    ,  12.200    ,  12.000    ,  11.800    ,
192        11.200    ,  11.500    ,  9.9000    ,  9.4000    ,  8.8000    ,
193        8.4000    ,  7.5000    ,  6.9000    ,  6.3000    ,  5.9000    ,
194        5.5000    ,  5.2000    ,  5.0000    ,  4.0000    ,  3.5000    ,
195        3.3000    ,  3.1000    ,  3.1000    ,  3.0000    ,  2.5000    ,
196        3.0000},
197
198      {160.83    ,  82.800    ,  58.575    ,  43.683    ,  34.792    , //12
199        28.650    ,  24.367    ,  20.917    ,  18.192    ,  16.300    ,
200        14.608    ,  13.017    ,  12.250    ,  11.700    ,  12.017    ,
201        14.075    ,  15.842    ,  16.433    ,  16.042    ,  15.008    ,
202        12.575    ,  10.708    ,  9.2000    ,  8.0167    ,  7.2833    ,
203        7.0750    ,  6.6333    ,  6.1250    ,  5.6583    ,  5.2750    ,
204        4.9333    ,  4.6250    ,  4.4583    ,  3.7333    ,  3.3833    ,
205        3.1833    ,  2.9833    ,  2.7500    ,  2.3667    ,  2.2000    ,
206        2.6000},
207
208      {300.00    ,  140.00    ,  97.000    ,  70.000    ,  55.000    , //13
209        45.000    ,  37.000    ,  31.000    ,  26.000    ,  23.000    ,
210        20.000    ,  17.000    ,  15.500    ,  14.500    ,  14.700    ,
211        18.500    ,  22.000    ,  23.000    ,  22.500    ,  20.700    ,
212        16.500    ,  14.000    ,  11.500    ,  9.6000    ,  8.6000    ,
213        8.5000    ,  8.3000    ,  7.6000    ,  7.0000    ,  6.4000    ,
214        5.9000    ,  5.5000    ,  5.3000    ,  4.4000    ,  4.1000    ,
215        3.9000    ,  3.7000    ,  3.3000    ,  2.6000    ,  2.5000    ,
216        3.0000},
217
218// Elastic cross section for p-p
219       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //14
220        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
221        20.000    ,  20.000    ,  20.500    ,  21.000    ,  22.000    ,
222        23.000    ,  24.000    ,  24.000    ,  24.400    ,  24.500    ,
223        25.000    ,  25.500    ,  26.000    ,  26.500    ,  27.000    ,
224        27.000    ,  26.000    ,  23.000    ,  21.500    ,  20.000    ,
225        19.000    ,  18.000    ,  17.000    ,  13.000    ,  11.500    ,
226        10.300    ,  9.4000    ,  9.0000    ,  8.8000    ,  7.0000    ,
227        7.5000},
228
229       {200.00    ,  163.00    ,  141.00    ,  120.00    ,  111.00    , //15
230        99.500    ,  92.500    ,  86.500    ,  82.000    ,  78.000    ,
231        74.000    ,  71.000    ,  67.500    ,  65.000    ,  62.500    ,
232        59.700    ,  58.100    ,  56.300    ,  54.700    ,  52.700    ,
233        50.000    ,  48.400    ,  47.000    ,  46.000    ,  45.200    ,
234        42.800    ,  39.200    ,  36.300    ,  32.800    ,  30.400    ,
235        28.100    ,  26.300    ,  24.500    ,  19.250    ,  16.840    ,
236        14.600    ,  12.340    ,  11.210    ,  8.8500    ,  7.5000    ,
237        7.5000}    ,
238
239       {4200.0    ,  440.00    ,  420.00    ,  400.00    ,  230.00    , //16
240        160.00    ,  105.00    ,  80.000    ,  62.000    ,  50.000    ,
241        45.000    ,  41.000    ,  38.000    ,  36.000    ,  35.000    ,
242        34.000    ,  33.000    ,  32.000    ,  31.500    ,  31.000    ,
243        30.500    ,  30.000    ,  29.500    ,  29.000    ,  28.500    ,
244        28.000    ,  26.000    ,  23.000    ,  21.500    ,  20.000    ,
245        19.000    ,  18.000    ,  17.000    ,  13.000    ,  11.500    ,
246        10.300    ,  9.4000    ,  9.0000    ,  8.8000    ,  7.0000    ,
247        7.5000}    ,
248
249       {185.88    ,  133.23    ,  119.37    ,  102.86    ,  93.102    , //17
250        82.752    ,  76.205    ,  71.008    ,  67.366    ,  64.096    ,
251        60.891    ,  58.501    ,  55.735    ,  53.773    ,  51.839    ,
252        49.671    ,  48.485    ,  47.045    ,  45.803    ,  44.306    ,
253        42.623    ,  41.786    ,  41.115    ,  40.630    ,  40.129    ,
254        38.242    ,  35.233    ,  32.662    ,  29.639    ,  27.573    ,
255        25.536    ,  23.948    ,  22.356    ,  17.723    ,  15.614    ,
256        13.653    ,  11.675    ,  10.653    ,  8.6198    ,  7.4464    ,
257        7.4821}    ,
258
259       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //18
260        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
261        20.000    ,  19.067    ,  19.333    ,  19.500    ,  19.833    ,
262        20.567    ,  21.800    ,  22.900    ,  23.869    ,  23.809    ,
263        22.161    ,  21.488    ,  19.732    ,  19.433    ,  19.345    ,
264        19.029    ,  18.121    ,  16.280    ,  15.258    ,  14.280    ,
265        13.644    ,  12.963    ,  12.316    ,  9.5333    ,  8.4333    ,
266        7.5728    ,  6.9696    ,  6.7518    ,  6.6175    ,  5.6000    ,
267        6.1145}    ,
268
269       {157.65    ,  73.701    ,  76.096    ,  68.571    ,  57.305    , //19
270        49.257    ,  43.616    ,  40.024    ,  38.098    ,  36.287    ,
271        34.674    ,  33.105    ,  31.712    ,  30.685    ,  29.613    ,
272        28.602    ,  28.336    ,  28.075    ,  27.786    ,  27.215    ,
273        26.380    ,  26.146    ,  25.108    ,  24.783    ,  24.360    ,
274        23.219    ,  21.431    ,  20.095    ,  18.382    ,  17.267    ,
275        16.100    ,  15.175    ,  14.271    ,  11.573    ,  10.305    ,
276        9.1471    ,  8.0149    ,  7.4349    ,  6.2499    ,  5.8928    ,
277        6.0774}    ,
278
279       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //20
280        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
281        20.000    ,  19.067    ,  19.333    ,  19.500    ,  19.833    ,
282        20.567    ,  21.800    ,  22.900    ,  23.869    ,  23.809    ,
283        22.161    ,  21.488    ,  19.732    ,  19.433    ,  19.345    ,
284        19.029    ,  18.121    ,  16.280    ,  15.258    ,  14.280    ,
285        13.644    ,  12.963    ,  12.316    ,  9.5333    ,  8.4333    ,
286        7.5728    ,  6.9696    ,  6.7518    ,  6.6175    ,  5.6000    ,
287        6.1145}    ,
288
289      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //21
290       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
291       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
292       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
293       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
294       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
295       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
296       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
297       0.00000E+00},
298
299       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //22
300        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
301        20.000    ,  19.067    ,  19.333    ,  19.500    ,  19.833    ,
302        20.567    ,  21.800    ,  22.900    ,  23.869    ,  23.809    ,
303        22.161    ,  21.488    ,  19.732    ,  19.433    ,  19.345    ,
304        19.029    ,  18.121    ,  16.280    ,  15.258    ,  14.280    ,
305        13.644    ,  12.963    ,  12.316    ,  9.5333    ,  8.4333    ,
306        7.5728    ,  6.9696    ,  6.7518    ,  6.6175    ,  5.6000    ,
307        6.1145}    ,
308
309       {185.88    ,  133.23    ,  119.37    ,  102.86    ,  93.102    , //23
310        82.752    ,  76.205    ,  71.008    ,  67.366    ,  64.096    ,
311        60.891    ,  58.104    ,  55.241    ,  53.140    ,  50.934    ,
312        48.660    ,  47.566    ,  46.585    ,  45.581    ,  44.003    ,
313        41.134    ,  39.374    ,  36.878    ,  35.523    ,  34.503    ,
314        32.334    ,  29.365    ,  27.370    ,  24.705    ,  22.921    ,
315        21.229    ,  19.879    ,  18.559    ,  14.625    ,  12.758    ,
316        11.041    ,  9.3440    ,  8.5484    ,  6.7104    ,  6.0000    ,
317        6.1131}    ,
318
319      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //24
320       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
321       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
322       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
323       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
324       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
325       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
326       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
327       0.00000E+00},
328
329       {157.65    ,  73.701    ,  76.096    ,  68.571    ,  57.305    , //25
330        49.257    ,  43.616    ,  40.024    ,  38.098    ,  36.287    ,
331        34.674    ,  33.105    ,  31.712    ,  30.685    ,  29.613    ,
332        28.602    ,  28.336    ,  28.075    ,  27.786    ,  27.215    ,
333        26.380    ,  26.146    ,  25.108    ,  24.783    ,  24.360    ,
334        23.219    ,  21.431    ,  20.095    ,  18.382    ,  17.267    ,
335        16.100    ,  15.175    ,  14.271    ,  11.573    ,  10.305    ,
336        9.1471    ,  8.0149    ,  7.4349    ,  6.2499    ,  5.8928    ,
337        6.0774}    ,
338
339       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //26
340        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
341        20.000    ,  18.133    ,  18.167    ,  18.000    ,  17.667    ,
342        18.133    ,  19.600    ,  21.800    ,  23.338    ,  23.118    ,
343        19.323    ,  17.476    ,  13.464    ,  12.367    ,  11.691    ,
344        11.057    ,  10.242    ,  9.5593    ,  9.0151    ,  8.5591    ,
345        8.2884    ,  7.9253    ,  7.6311    ,  6.0667    ,  5.3667    ,
346        4.8456    ,  4.5392    ,  4.5036    ,  4.4351    ,  4.2000    ,
347        4.7289}    ,
348
349       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //27
350        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
351        20.000    ,  18.133    ,  18.167    ,  18.000    ,  17.667    ,
352        18.133    ,  19.600    ,  21.800    ,  23.338    ,  23.118    ,
353        19.323    ,  17.476    ,  13.464    ,  12.367    ,  11.691    ,
354        11.057    ,  10.242    ,  9.5593    ,  9.0151    ,  8.5591    ,
355        8.2884    ,  7.9253    ,  7.6311    ,  6.0667    ,  5.3667    ,
356        4.8456    ,  4.5392    ,  4.5036    ,  4.4351    ,  4.2000    ,
357        4.7289}    ,
358
359       {157.65    ,  73.701    ,  76.096    ,  68.571    ,  57.305    , //28
360        49.257    ,  43.616    ,  40.024    ,  38.098    ,  36.287    ,
361        34.674    ,  32.708    ,  31.218    ,  30.052    ,  28.707    ,
362        27.591    ,  27.417    ,  27.615    ,  27.564    ,  26.913    ,
363        24.891    ,  23.734    ,  20.871    ,  19.677    ,  18.734    ,
364        17.311    ,  15.563    ,  14.803    ,  13.448    ,  12.615    ,
365        11.794    ,  11.106    ,  10.474    ,  8.4745    ,  7.4498    ,
366        6.5350    ,  5.6835    ,  5.3300    ,  4.3406    ,  4.4464    ,
367        4.7083}    ,
368
369       {143.53    ,  43.935    ,  54.462    ,  51.429    ,  39.407    , //29
370        32.510    ,  27.321    ,  24.532    ,  23.465    ,  22.383    ,
371        21.566    ,  20.209    ,  19.453    ,  18.825    ,  18.046    ,
372        17.562    ,  17.802    ,  18.360    ,  18.667    ,  18.519    ,
373        17.514    ,  17.120    ,  14.985    ,  14.306    ,  13.663    ,
374        12.753    ,  11.596    ,  11.165    ,  10.287    ,  9.7882    ,
375        9.2294    ,  8.7539    ,  8.3300    ,  6.9480    ,  6.2234    ,
376        5.5881    ,  5.0189    ,  4.7733    ,  4.1104    ,  4.3929    ,
377        4.6905}    ,
378
379      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //30
380       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
381       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
382       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
383       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
384       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
385       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
386       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
387       0.00000E+00},
388
389      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //31
390       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
391       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
392       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
393       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
394       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
395       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
396       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
397       0.00000E+00},
398
399      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //32
400       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
401       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
402       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
403       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
404       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
405       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
406       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
407       0.00000E+00},
408
409       {1100.0    ,  115.00    ,  105.00    ,  100.00    ,  56.000    , //33
410        40.000    ,  27.000    ,  22.000    ,  21.000    ,  20.000    ,
411        20.000    ,  18.133    ,  18.167    ,  18.000    ,  17.667    ,
412        18.133    ,  19.600    ,  21.800    ,  23.338    ,  23.118    ,
413        19.323    ,  17.476    ,  13.464    ,  12.367    ,  11.691    ,
414        11.057    ,  10.242    ,  9.5593    ,  9.0151    ,  8.5591    ,
415        8.2884    ,  7.9253    ,  7.6311    ,  6.0667    ,  5.3667    ,
416        4.8456    ,  4.5392    ,  4.5036    ,  4.4351    ,  4.2000    ,
417        4.7289}   ,
418
419       {143.53    ,  43.935    ,  54.462    ,  51.429    ,  39.407    , //34
420        32.510    ,  27.321    ,  24.532    ,  23.465    ,  22.383    ,
421        21.566    ,  20.209    ,  19.453    ,  18.825    ,  18.046    ,
422        17.562    ,  17.802    ,  18.360    ,  18.667    ,  18.519    ,
423        17.514    ,  17.120    ,  14.985    ,  14.306    ,  13.663    ,
424        12.753    ,  11.596    ,  11.165    ,  10.287    ,  9.7882    ,
425        9.2294    ,  8.7539    ,  8.3300    ,  6.9480    ,  6.2234    ,
426        5.5881    ,  5.0189    ,  4.7733    ,  4.1104    ,  4.3929    ,
427        4.6905}    ,
428
429      {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //35
430       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
431       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
432       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
433       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
434       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
435       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
436       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
437       0.00000E+00}
438};
439
440//---------------------------------------------------------------------
441// Inelastic scattering on free protons
442//---------------------------------------------------------------------
443
444G4float G4HadronCrossSections::csin[NPARTS][TSIZE] = {
445
446       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //1
447       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
448       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
449       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
450       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
451       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
452       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
453       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
454       0.00000E+00},
455
456       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //2
457       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
458       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
459       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
460       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
461       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
462       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
463       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
464       0.00000E+00},
465
466       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //3
467       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
468       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
469       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
470       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
471       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
472       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
473       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
474       0.00000E+00},
475
476       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //4
477       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
478       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
479       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
480       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
481       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
482       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
483       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
484       0.00000E+00},
485
486       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //5
487       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
488       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
489       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
490       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
491       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
492       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
493       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
494       0.00000E+00},
495
496       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //6
497       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
498       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
499       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
500       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
501       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
502       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
503       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
504       0.00000E+00},
505
506// Inelastic cross section for piplus - p
507       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //7
508       0.00000E+00, 0.00000E+00, 0.50000    ,  1.2000    ,  1.7000    ,
509        2.2500    ,  3.0000    ,  3.6000    ,  4.5000    ,  5.4000    ,
510        6.3000    ,  8.6000    ,  9.0000    ,  10.000    ,  11.500    ,
511        14.000    ,  17.000    ,  19.500    ,  22.000    ,  24.000    ,
512        21.500    ,  18.500    ,  19.000    ,  20.500    ,  22.200    ,
513        23.000    ,  23.300    ,  23.000    ,  21.000    ,  20.500    ,
514        20.200    ,  20.100    ,  20.000    ,  20.000    ,  20.000    ,
515        21.000},
516
517
518       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //8
519       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
520       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
521       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
522       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
523       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
524       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
525       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
526       0.00000E+00},
527
528// Inelastic cross section for piminus - p
529       {0.00000E+00,  3.0000    ,  9.2000    ,  20.500    ,  36.500    , //9
530        45.000    ,  28.000    ,  19.500    ,  15.500    ,  14.200    ,
531        15.500    ,  17.500    ,  20.000    ,  23.000    ,  26.000    ,
532        20.000    ,  23.000    ,  26.500    ,  32.000    ,  35.000    ,
533        28.500    ,  22.000    ,  22.500    ,  23.500    ,  24.000    ,
534        24.500    ,  26.000    ,  27.500    ,  27.500    ,  27.000    ,
535        26.500    ,  25.500    ,  25.000    ,  23.000    ,  22.500    ,
536        22.200    ,  22.000    ,  22.000    ,  21.200    ,  20.700    ,
537        21.000}    ,
538
539       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //10
540       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
541       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
542       0.50000    ,  1.5000    ,  2.7000    ,  3.8000    ,  4.8000    ,
543        6.5000    ,  7.6000    ,  8.4000    ,  9.0000    ,  9.4000    ,
544        9.8000    ,  10.500    ,  11.000    ,  11.500    ,  11.800    ,
545        12.200    ,  12.400    ,  12.600    ,  13.200    ,  13.500    ,
546        13.700    ,  14.000    ,  14.200    ,  14.500    ,  16.400    ,
547        17.000}    ,
548
549       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //11
550       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
551       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
552       0.50000    ,  1.5000    ,  2.7000    ,  3.8000    ,  4.8000    ,
553        6.5000    ,  7.6000    ,  8.4000    ,  9.0000    ,  9.4000    ,
554        9.8000    ,  10.500    ,  11.000    ,  11.500    ,  11.800    ,
555        12.200    ,  12.400    ,  12.600    ,  13.200    ,  13.500    ,
556        13.700    ,  14.000    ,  14.200    ,  14.500    ,  16.400    ,
557        17.000}    ,
558
559        {266.67    ,  133.33    ,  83.333    ,  57.083    ,  44.500    , //12
560        33.250    ,  24.583    ,  20.833    ,  18.333    ,  16.083    ,
561        15.625    ,  15.083    ,  14.833    ,  15.083    ,  15.833    ,
562        17.042    ,  18.958    ,  20.758    ,  22.533    ,  22.825    ,
563        21.250    ,  18.567    ,  17.767    ,  18.100    ,  19.933    ,
564        20.783    ,  21.225    ,  21.000    ,  20.558    ,  20.258    ,
565        20.017    ,  19.767    ,  19.600    ,  19.183    ,  18.850    ,
566        18.575    ,  18.350    ,  18.175    ,  17.808    ,  17.558    ,
567        19.250}    ,
568
569        {400.00    ,  200.00    ,  120.00    ,  81.000    ,  62.000    , //13
570        47.000    ,  35.000    ,  28.000    ,  24.000    ,  21.000    ,
571        19.500    ,  19.000    ,  18.800    ,  19.000    ,  20.000    ,
572        21.000    ,  23.000    ,  25.000    ,  27.000    ,  27.500    ,
573        25.500    ,  22.000    ,  20.800    ,  21.000    ,  23.000    ,
574        24.000    ,  24.000    ,  23.800    ,  23.000    ,  22.500    ,
575        22.000    ,  21.600    ,  21.400    ,  21.000    ,  20.500    ,
576        20.200    ,  19.800    ,  19.500    ,  18.600    ,  17.500    ,
577        20.000}    ,
578
579// Inelastic cross section for p - p
580       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //14
581       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
582       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
583       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.10000    ,  1.5000    ,
584        7.0000    ,  12.000    ,  17.000    ,  19.500    ,  20.500    ,
585        22.000    ,  23.500    ,  24.800    ,  25.800    ,  26.500    ,
586        27.000    ,  27.500    ,  28.000    ,  30.000    ,  31.000    ,
587        32.000    ,  32.500    ,  32.500    ,  33.000    ,  33.500    ,
588        33.500}    ,
589
590        {1500.0    ,  1160.0    ,  310.00    ,  230.00    ,  178.00    , //15
591        153.00    ,  134.00    ,  124.00    ,  113.00    ,  106.00    ,
592        101.00    ,  96.000    ,  92.000    ,  89.000    ,  87.000    ,
593        84.000    ,  81.000    ,  78.500    ,  76.500    ,  75.000    ,
594        72.000    ,  70.000    ,  68.000    ,  64.500    ,  63.000    ,
595        62.000    ,  61.000    ,  59.500    ,  58.500    ,  56.500    ,
596        56.500    ,  56.000    ,  55.500    ,  52.000    ,  50.000    ,
597        48.000    ,  45.000    ,  44.000    ,  39.200    ,  34.500    ,
598        34.500}    ,
599
600       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //16
601       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
602       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
603       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.10000    ,  1.5000    ,
604        7.0000    ,  12.000    ,  17.000    ,  19.500    ,  20.500    ,
605        22.000    ,  23.500    ,  24.800    ,  25.800    ,  26.500    ,
606        27.000    ,  27.500    ,  28.000    ,  30.000    ,  31.000    ,
607        32.000    ,  32.500    ,  32.500    ,  33.000    ,  33.500    ,
608        34.000}    ,
609
610        {1394.1    ,  948.17    ,  262.43    ,  197.14    ,  149.30    , //17
611        127.25    ,  110.39    ,  101.79    ,  92.834    ,  87.104    ,
612        83.109    ,  79.099    ,  75.965    ,  73.627    ,  72.161    ,
613        69.889    ,  67.595    ,  65.595    ,  64.057    ,  63.054    ,
614        61.377    ,  60.434    ,  59.485    ,  56.970    ,  55.931    ,
615        55.398    ,  54.827    ,  53.538    ,  52.861    ,  51.247    ,
616        51.344    ,  50.992    ,  50.644    ,  47.876    ,  46.358    ,
617        44.887    ,  42.577    ,  41.815    ,  38.180    ,  34.254    ,
618        34.418}    ,
619
620       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //18
621       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
622       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
623       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.97815E-01,  1.4577    ,
624        6.2052    ,  10.112    ,  12.902    ,  14.300    ,  14.688    ,
625        15.505    ,  16.379    ,  17.554    ,  18.309    ,  18.920    ,
626        19.389    ,  19.804    ,  20.284    ,  22.000    ,  22.733    ,
627        23.527    ,  24.097    ,  24.382    ,  24.816    ,  26.800    ,
628        27.719}    ,
629
630        {1182.4    ,  524.50    ,  167.30    ,  131.43    ,  91.895    , //19
631        75.743    ,  63.184    ,  57.376    ,  52.502    ,  49.313    ,
632        47.326    ,  44.762    ,  43.222    ,  42.015    ,  41.221    ,
633        40.244    ,  39.504    ,  39.145    ,  38.860    ,  38.731    ,
634        37.987    ,  37.814    ,  36.326    ,  34.750    ,  33.953    ,
635        33.635    ,  33.349    ,  32.938    ,  32.785    ,  32.092    ,
636        32.373    ,  32.312    ,  32.329    ,  31.261    ,  30.597    ,
637        30.073    ,  29.228    ,  29.182    ,  27.683    ,  27.107    ,
638        27.956}    ,
639
640       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //20
641       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
642       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
643       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.97815E-01,  1.4577    ,
644        6.2052    ,  10.112    ,  12.902    ,  14.300    ,  14.688    ,
645        15.505    ,  16.379    ,  17.554    ,  18.309    ,  18.920    ,
646        19.389    ,  19.804    ,  20.284    ,  22.000    ,  22.733    ,
647        23.527    ,  24.097    ,  24.382    ,  24.816    ,  26.800    ,
648        27.719}    ,
649
650       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //21
651       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
652       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
653       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
654       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
655       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
656       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
657       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
658       0.00000E+00},
659
660       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //22
661       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
662       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
663       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.97815E-01,  1.4577    ,
664        6.2052    ,  10.112    ,  12.902    ,  14.300    ,  14.688    ,
665        15.505    ,  16.379    ,  17.554    ,  18.309    ,  18.920    ,
666        19.389    ,  19.804    ,  20.284    ,  22.000    ,  22.733    ,
667        23.527    ,  24.097    ,  24.382    ,  24.816    ,  26.800    ,
668        27.719}    ,
669
670        {1394.1    ,  948.17    ,  262.43    ,  197.14    ,  149.30    , //23
671        127.25    ,  110.39    ,  101.79    ,  92.834    ,  87.104    ,
672        83.109    ,  78.563    ,  75.292    ,  72.760    ,  70.900    ,
673        68.467    ,  66.314    ,  64.955    ,  63.746    ,  62.623    ,
674        59.233    ,  56.946    ,  53.355    ,  49.810    ,  48.090    ,
675        46.839    ,  45.695    ,  44.863    ,  44.062    ,  42.599    ,
676        42.684    ,  42.328    ,  42.041    ,  39.508    ,  37.880    ,
677        36.299    ,  34.075    ,  33.553    ,  29.723    ,  27.600    ,
678        28.120}    ,
679
680       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //24
681       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
682       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
683       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
684       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
685       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
686       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
687       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
688       0.00000E+00},
689
690        {1182.4    ,  524.50    ,  167.30    ,  131.43    ,  91.895    , //25
691        75.743    ,  63.184    ,  57.376    ,  52.502    ,  49.313    ,
692        47.326    ,  44.762    ,  43.222    ,  42.015    ,  41.221    ,
693        40.244    ,  39.504    ,  39.145    ,  38.860    ,  38.731    ,
694        37.987    ,  37.814    ,  36.326    ,  34.750    ,  33.953    ,
695        33.635    ,  33.349    ,  32.938    ,  32.785    ,  32.092    ,
696        32.373    ,  32.312    ,  32.329    ,  31.261    ,  30.597    ,
697        30.073    ,  29.228    ,  29.182    ,  27.683    ,  27.107    ,
698        27.956}    ,
699
700       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //26
701       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
702       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
703       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.95639E-01,  1.4154    ,
704        5.4104    ,  8.2240    ,  8.8031    ,  9.1000    ,  8.8761    ,
705        9.0095    ,  9.2576    ,  10.307    ,  10.818    ,  11.341    ,
706        11.778    ,  12.108    ,  12.569    ,  14.000    ,  14.467    ,
707        15.054    ,  15.694    ,  16.263    ,  16.632    ,  20.100    ,
708        21.438}    ,
709
710       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //27
711       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
712       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
713       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.95639E-01,  1.4154    ,
714        5.4104    ,  8.2240    ,  8.8031    ,  9.1000    ,  8.8761    ,
715        9.0095    ,  9.2576    ,  10.307    ,  10.818    ,  11.341    ,
716        11.778    ,  12.108    ,  12.569    ,  14.000    ,  14.467    ,
717        15.054    ,  15.694    ,  16.263    ,  16.632    ,  20.100    ,
718        21.438}    ,
719
720        {1182.4    ,  524.50    ,  167.30    ,  131.43    ,  91.895    , //28
721        75.743    ,  63.184    ,  57.376    ,  52.502    ,  49.313    ,
722        47.326    ,  44.225    ,  42.549    ,  41.148    ,  39.960    ,
723        38.822    ,  38.223    ,  38.505    ,  38.549    ,  38.301    ,
724        35.843    ,  34.326    ,  30.196    ,  27.590    ,  26.112    ,
725        25.076    ,  24.217    ,  24.264    ,  23.985    ,  23.445    ,
726        23.713    ,  23.647    ,  23.726    ,  22.892    ,  22.119    ,
727        21.485    ,  20.726    ,  20.921    ,  19.226    ,  20.454    ,
728        21.658}    ,
729
730        {1076.5    ,  312.66    ,  119.74    ,  98.571    ,  63.193    , //29
731        49.990    ,  39.579    ,  35.168    ,  32.335    ,  30.417    ,
732        29.434    ,  27.325    ,  26.514    ,  25.775    ,  25.120    ,
733        24.711    ,  24.818    ,  25.600    ,  26.106    ,  26.355    ,
734        25.220    ,  24.760    ,  21.681    ,  20.060    ,  19.044    ,
735        18.474    ,  18.044    ,  18.301    ,  18.347    ,  18.192    ,
736        18.557    ,  18.639    ,  18.870    ,  18.769    ,  18.478    ,
737        18.372    ,  18.302    ,  18.735    ,  18.206    ,  20.207    ,
738        21.576}    ,
739
740       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //30
741       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
742       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
743       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
744       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
745       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
746       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
747       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
748       0.00000E+00},
749
750       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //31
751       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
752       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
753       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
754       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
755       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
756       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
757       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
758       0.00000E+00},
759
760       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //32
761       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
762       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
763       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
764       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
765       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
766       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
767       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
768       0.00000E+00},
769
770       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //33
771       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
772       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
773       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.95639E-01,  1.4154    ,
774        5.4104    ,  8.2240    ,  8.8031    ,  9.1000    ,  8.8761    ,
775        9.0095    ,  9.2576    ,  10.307    ,  10.818    ,  11.341    ,
776        11.778    ,  12.108    ,  12.569    ,  14.000    ,  14.467    ,
777        15.054    ,  15.694    ,  16.263    ,  16.632    ,  20.100    ,
778        21.438}    ,
779
780        {1076.5    ,  312.66    ,  119.74    ,  98.571    ,  63.193    , //34
781        49.990    ,  39.579    ,  35.168    ,  32.335    ,  30.417    ,
782        29.434    ,  27.325    ,  26.514    ,  25.775    ,  25.120    ,
783        24.711    ,  24.818    ,  25.600    ,  26.106    ,  26.355    ,
784        25.220    ,  24.760    ,  21.681    ,  20.060    ,  19.044    ,
785        18.474    ,  18.044    ,  18.301    ,  18.347    ,  18.192    ,
786        18.557    ,  18.639    ,  18.870    ,  18.769    ,  18.478    ,
787        18.372    ,  18.302    ,  18.735    ,  18.206    ,  20.207    ,
788        21.576}    ,
789
790       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //35
791       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
792       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
793       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
794       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
795       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
796       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
797       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
798       0.00000E+00}
799
800};
801
802//---------------------------------------------------------------------
803// Elastic scattering of pion on nucleus
804//---------------------------------------------------------------------
805
806G4float G4HadronCrossSections::cspiel[3][TSIZE] = {
807
808// Elastic cross section for Pi on Al (mb)
809       {0.00000E+00,  350.00    ,  580.00    ,  600.00    ,  550.00    , //1
810        450.00    ,  410.00    ,  370.00    ,  340.00    ,  230.00    ,
811        220.00    ,  205.00    ,  180.00    ,  155.00    ,  145.00    ,
812        140.00    ,  160.00    ,  195.00    ,  235.00    ,  250.00    ,
813        270.00    ,  280.00    ,  300.00    ,  300.00    ,  290.00    ,
814        285.00    ,  265.00    ,  240.00    ,  230.00    ,  222.00    ,
815        204.00    ,  196.00    ,  190.00    ,  170.00    ,  170.00    ,
816        160.00    ,  150.00    ,  140.00    ,  120.00    ,  80.000    ,
817        80.000},
818
819// Elastic cross section for Pi on Cu (mb)
820       {0.00000E+00,  700.00    ,  1000.0    ,  1200.0    ,  1300.0    , //2
821        1300.0    ,  1250.0    ,  1250.0    ,  1100.0    ,  1000.0    ,
822        940.00    ,  740.00    ,  700.00    ,  670.00    ,  660.00    ,
823        670.00    ,  680.00    ,  700.00    ,  735.00    ,  800.00    ,
824        810.00    ,  820.00    ,  820.00    ,  810.00    ,  800.00    ,
825        800.00    ,  700.00    ,  600.00    ,  500.00    ,  470.00    ,
826        440.00    ,  410.00    ,  380.00    ,  330.00    ,  330.00    ,
827        330.00    ,  330.00    ,  330.00    ,  285.00    ,  240.00    ,
828        240.00},
829
830// Elastic cross section for Pi on Pb (mb)
831       {0.00000E+00,  1700.0    ,  2200.0    ,  2200.0    ,  1800.0    , //3
832        1300.0    ,  1200.0    ,  900.00    ,  900.00    ,  1000.0    ,
833        1100.0    ,  1300.0    ,  1400.0    ,  1420.0    ,  1490.0    ,
834        1560.0    ,  1580.0    ,  1690.0    ,  1795.0    ,  2000.0    ,
835        2070.0    ,  2140.0    ,  2050.0    ,  2010.0    ,  1970.0    ,
836        1880.0    ,  1690.0    ,  1500.0    ,  1420.0    ,  1390.0    ,
837        1350.0    ,  1360.0    ,  1370.0    ,  1280.0    ,  1290.0    ,
838        1295.0    ,  1250.0    ,  1200.0    ,  1050.0    ,  900.00    ,
839        900.00}
840};
841
842//---------------------------------------------------------------------
843// Inelastic scattering of pion on nucleus
844//---------------------------------------------------------------------
845
846G4float G4HadronCrossSections::cspiin[3][TSIZE] = {
847
848// Inelastic cross section for Pi on Al (mb)
849       {0.00000E+00,  200.00    ,  320.00    ,  500.00    ,  600.00    , //1
850        600.00    ,  590.00    ,  530.00    ,  510.00    ,  470.00    ,
851        430.00    ,  425.00    ,  420.00    ,  425.00    ,  425.00    ,
852        430.00    ,  430.00    ,  435.00    ,  435.00    ,  440.00    ,
853        430.00    ,  430.00    ,  420.00    ,  420.00    ,  420.00    ,
854        415.00    ,  415.00    ,  410.00    ,  410.00    ,  408.00    ,
855        406.00    ,  404.00    ,  400.00    ,  380.00    ,  340.00    ,
856        340.00    ,  340.00    ,  340.00    ,  340.00    ,  340.00    ,
857        340.00}    ,
858
859// Inelastic cross section for Pi on Cu (mb)
860       {0.00000E+00,  400.00    ,  800.00    ,  1000.0    ,  1100.0    , //2
861        1200.0    ,  1150.0    ,  1050.0    ,  1000.0    ,  900.00    ,
862        860.00    ,  860.00    ,  850.00    ,  850.00    ,  840.00    ,
863        830.00    ,  820.00    ,  810.00    ,  805.00    ,  800.00    ,
864        800.00    ,  800.00    ,  800.00    ,  800.00    ,  800.00    ,
865        800.00    ,  800.00    ,  800.00    ,  800.00    ,  780.00    ,
866        760.00    ,  740.00    ,  720.00    ,  720.00    ,  700.00    ,
867        690.00    ,  680.00    ,  670.00    ,  665.00    ,  660.00    ,
868        660.00}    ,
869
870// Inelastic cross section for Pi on Pb (mb)
871       {0.00000E+00,  1000.0    ,  1900.0    ,  2600.0    ,  2900.0    , //3
872        3000.0    ,  2800.0    ,  2600.0    ,  2500.0    ,  2300.0    ,
873        2200.0    ,  2000.0    ,  1900.0    ,  1880.0    ,  1860.0    ,
874        1840.0    ,  1820.0    ,  1810.0    ,  1805.0    ,  1800.0    ,
875        1780.0    ,  1760.0    ,  1750.0    ,  1740.0    ,  1730.0    ,
876        1720.0    ,  1710.0    ,  1700.0    ,  1680.0    ,  1660.0    ,
877        1650.0    ,  1640.0    ,  1630.0    ,  1620.0    ,  1610.0    ,
878        1605.0    ,  1600.0    ,  1600.0    ,  1550.0    ,  1500.0    ,
879        1500.0}   
880};
881
882//---------------------------------------------------------------------
883// Elastic scattering of proton on nucleus
884//---------------------------------------------------------------------
885
886G4float G4HadronCrossSections::cspnel[3][TSIZE] = {
887
888// Elastic cross section for P on Al (mb)
889        {2100.0    ,  1800.0    ,  1500.0    ,  1050.0    ,  900.00    , //1
890        950.00    ,  800.00    ,  650.00    ,  570.00    ,  390.00    ,
891        300.00    ,  240.00    ,  230.00    ,  230.00    ,  220.00    ,
892        220.00    ,  225.00    ,  225.00    ,  240.00    ,  240.00    ,
893        290.00    ,  330.00    ,  335.00    ,  350.00    ,  355.00    ,
894        370.00    ,  350.00    ,  330.00    ,  310.00    ,  290.00    ,
895        270.00    ,  265.00    ,  260.00    ,  230.00    ,  210.00    ,
896        210.00    ,  200.00    ,  200.00    ,  190.00    ,  180.00    ,
897        180.00},
898
899// Elastic cross section for P on Cu (mb)
900        {3800.0    ,  2900.0    ,  1850.0    ,  1550.0    ,  1450.0    , //2
901        1520.0    ,  1460.0    ,  1300.0    ,  1140.0    ,  880.00    ,
902        700.00    ,  620.00    ,  540.00    ,  560.00    ,  460.00    ,
903        460.00    ,  470.00    ,  470.00    ,  480.00    ,  480.00    ,
904        580.00    ,  600.00    ,  610.00    ,  620.00    ,  620.00    ,
905        620.00    ,  590.00    ,  580.00    ,  460.00    ,  440.00    ,
906        420.00    ,  400.00    ,  480.00    ,  430.00    ,  380.00    ,
907        380.00    ,  380.00    ,  380.00    ,  380.00    ,  380.00    ,
908        380.00},
909
910// Elastic cross section for P on Pb (mb)
911        {7000.0    ,  6000.0    ,  4500.0    ,  3350.0    ,  2700.0    , //3
912        3000.0    ,  3550.0    ,  3970.0    ,  3280.0    ,  2490.0    ,
913        2100.0    ,  1510.0    ,  1440.0    ,  1370.0    ,  1370.0    ,
914        1370.0    ,  1400.0    ,  1400.0    ,  1420.0    ,  1420.0    ,
915        1440.0    ,  1460.0    ,  1460.0    ,  1450.0    ,  1450.0    ,
916        1470.0    ,  1400.0    ,  1400.0    ,  1380.0    ,  1370.0    ,
917        1360.0    ,  1350.0    ,  1340.0    ,  1330.0    ,  1320.0    ,
918        1310.0    ,  1305.0    ,  1300.0    ,  1300.0    ,  1300.0    ,
919        1300.0}
920};
921
922//---------------------------------------------------------------------
923// Inelastic scattering of proton on nucleus
924//---------------------------------------------------------------------
925
926G4float G4HadronCrossSections::cspnin[3][TSIZE] = {
927
928// Inelastic cross section for P on Al (mb)
929       {0.00000E+00,  200.00    ,  400.00    ,  800.00    ,  800.00    , //1
930        550.00    ,  500.00    ,  450.00    ,  430.00    ,  410.00    ,
931        400.00    ,  390.00    ,  380.00    ,  370.00    ,  370.00    ,
932        370.00    ,  365.00    ,  365.00    ,  360.00    ,  360.00    ,
933        360.00    ,  360.00    ,  365.00    ,  370.00    ,  375.00    ,
934        380.00    ,  400.00    ,  410.00    ,  420.00    ,  430.00    ,
935        440.00    ,  440.00    ,  440.00    ,  440.00    ,  440.00    ,
936        440.00    ,  440.00    ,  440.00    ,  440.00    ,  440.00    ,
937        440.00}    ,
938
939// Inelastic cross section for P on Cu (mb)
940       {0.00000E+00,  400.00    ,  950.00    ,  1050.0    ,  1050.0    , //2
941        980.00    ,  940.00    ,  900.00    ,  860.00    ,  820.00    ,
942        800.00    ,  780.00    ,  760.00    ,  740.00    ,  740.00    ,
943        740.00    ,  730.00    ,  730.00    ,  720.00    ,  720.00    ,
944        720.00    ,  720.00    ,  730.00    ,  740.00    ,  750.00    ,
945        760.00    ,  800.00    ,  820.00    ,  820.00    ,  820.00    ,
946        820.00    ,  820.00    ,  820.00    ,  820.00    ,  820.00    ,
947        820.00    ,  820.00    ,  820.00    ,  820.00    ,  820.00    ,
948        820.00}    ,
949
950// Inelastic cross section for P on Pb (mb)
951       {0.00000E+00, 0.00000E+00,  500.00    ,  1450.0    ,  1700.0    , //3
952        1800.0    ,  1750.0    ,  1730.0    ,  1720.0    ,  1710.0    ,
953        1700.0    ,  1690.0    ,  1660.0    ,  1630.0    ,  1630.0    ,
954        1630.0    ,  1600.0    ,  1600.0    ,  1580.0    ,  1580.0    ,
955        1580.0    ,  1580.0    ,  1600.0    ,  1630.0    ,  1650.0    ,
956        1670.0    ,  1760.0    ,  1800.0    ,  1800.0    ,  1800.0    ,
957        1800.0    ,  1800.0    ,  1800.0    ,  1800.0    ,  1800.0    ,
958        1800.0    ,  1800.0    ,  1800.0    ,  1800.0    ,  1800.0    ,
959        1800.0}   
960};
961
962//---------------------------------------------------------------------
963// Lab kinetic energy in GeV
964//---------------------------------------------------------------------
965G4float G4HadronCrossSections::elab[NELAB] = {
966       0.10000E-03, 0.20000E-03, 0.30000E-03, 0.40000E-03, 0.50000E-03,
967       0.70000E-03, 0.10000E-02, 0.20000E-02, 0.30000E-02, 0.40000E-02,
968       0.50000E-02, 0.70000E-02, 0.10000E-01, 0.15000E-01, 0.20000E-01,
969       0.25000E-01, 0.32700E-01
970};
971
972//---------------------------------------------------------------------
973// Tables for low-energy (< 32.7 MeV) neutrons
974//---------------------------------------------------------------------
975
976// Atomic weight
977G4float G4HadronCrossSections::cnlwat[NCNLW] = {
978        1.0000    ,  16.000    ,  27.000    ,  56.000    ,  59.000    ,
979        64.000    ,  91.000    ,  112.00    ,  119.00    ,  127.00    ,
980        137.00    ,  181.00    ,  207.00    ,  209.00    ,  238.00   
981};
982// Elastic cross section
983G4float G4HadronCrossSections::cnlwel[NCNLW][NELAB] = {
984        {6000.0    ,  5500.0    ,  5200.0    ,  4900.0    ,  4800.0    , //1
985        4400.0    ,  4000.0    ,  2900.0    ,  2200.0    ,  1800.0    ,
986        1400.0    ,  1100.0    ,  900.00    ,  700.00    ,  600.00    ,
987        560.00    ,  520.00}    ,
988        {5400.0    ,  5050.0    ,  4800.0    ,  4600.0    ,  4399.0    , //2
989        4090.0    ,  3700.0    ,  2600.0    ,  1950.0    ,  1600.0    ,
990        1300.0    ,  900.00    ,  700.00    ,  800.00    ,  1050.0    ,
991        1250.0    ,  1320.0}    ,
992        {5500.0    ,  5150.0    ,  4900.0    ,  4699.0    ,  4490.0    , //3
993        4150.0    ,  3750.0    ,  2790.0    ,  2100.0    ,  1650.0    ,
994        1300.0    ,  950.00    ,  800.00    ,  860.00    ,  1000.0    ,
995        1090.0    ,  1080.0}    ,
996        {5499.0    ,  4970.0    ,  4450.0    ,  4080.0    ,  3750.0    , //4
997        3380.0    ,  2900.0    ,  2400.0    ,  2380.0    ,  2350.0    ,
998        2300.0    ,  2100.0    ,  1720.0    ,  1370.0    ,  1200.0    ,
999        1060.0    ,  870.00}    ,
1000        {5399.0    ,  4710.0    ,  4180.0    ,  3760.0    ,  3460.0    , //5
1001        3150.0    ,  2730.0    ,  2270.0    ,  1850.0    ,  1850.0    ,
1002        2130.0    ,  2330.0    ,  2120.0    ,  1640.0    ,  1310.0    ,
1003        1100.0    ,  1050.0}    ,
1004        {5099.0    ,  4405.0    ,  3825.0    ,  3455.0    ,  3125.0    , //6
1005        2695.0    ,  2350.0    ,  1850.0    ,  1580.0    ,  1820.0    ,
1006        2050.0    ,  2210.0    ,  2000.0    ,  1590.0    ,  1310.0    ,
1007        1120.0    ,  1040.0}  ,
1008        {6290.0    ,  5960.0    ,  5640.0    ,  5370.0    ,  5150.0    , //7
1009        4800.0    ,  4250.0    ,  3150.0    ,  2470.0    ,  2100.0    ,
1010        2230.0    ,  2420.0    ,  2450.0    ,  2050.0    ,  1760.0    ,
1011        1550.0    ,  1330.0}  ,
1012        {6885.0    ,  6650.0    ,  6350.0    ,  6150.0    ,  6000.0    , //8
1013        5700.0    ,  5360.0    ,  4250.0    ,  2800.0    ,  1870.0    ,
1014        1810.0    ,  1820.0    ,  2170.0    ,  2450.0    ,  2150.0    ,
1015        1700.0    ,  1390.0}  ,
1016        {6600.0    ,  6500.0    ,  6400.0    ,  6249.0    ,  6190.0    , //9
1017        5950.0    ,  5520.0    ,  4250.0    ,  2750.0    ,  1900.0    ,
1018        1850.0    ,  1950.0    ,  2340.0    ,  2800.0    ,  2540.0    ,
1019        2100.0    ,  1760.0}  ,
1020        {7400.0    ,  7200.0    ,  6999.0    ,  6840.0    ,  6655.0    , //10
1021        6320.0    ,  5820.0    ,  4400.0    ,  2850.0    ,  2000.0    ,
1022        1800.0    ,  1800.0    ,  2150.0    ,  2600.0    ,  2350.0    ,
1023        1950.0    ,  2100.0}  ,
1024        {7900.0    ,  7700.0    ,  7499.0    ,  7390.0    ,  7202.0    , //11
1025        6810.0    ,  6360.0    ,  4920.0    ,  3450.0    ,  2600.0    ,
1026        2200.0    ,  1950.0    ,  2300.0    ,  2800.0    ,  2650.0    ,
1027        2250.0    ,  2050.0}  ,
1028        {7900.0    ,  7750.0    ,  7699.0    ,  7590.0    ,  7450.0    , //12
1029        7200.0    ,  6850.0    ,  5650.0    ,  4400.0    ,  3700.0    ,
1030        3400.0    ,  2800.0    ,  2700.0    ,  3100.0    ,  3250.0    ,
1031        3100.0    ,  2750.0}  ,
1032        {6100.0    ,  5950.0    ,  5750.0    ,  5599.0    ,  5440.0    , //13
1033        5200.0    ,  4800.0    ,  4300.0    ,  5800.0    ,  5750.0    ,
1034        4800.0    ,  3420.0    ,  2650.0    ,  3200.0    ,  3650.0    ,
1035        3500.0    ,  2980.0}  ,
1036        {6100.0    ,  5950.0    ,  5750.0    ,  5599.0    ,  5440.0    , //14
1037        5200.0    ,  4800.0    ,  4300.0    ,  5800.0    ,  5750.0    ,
1038        4800.0    ,  3420.0    ,  2650.0    ,  3200.0    ,  3650.0    ,
1039        3500.0    ,  2980.0}  ,
1040        {6600.0    ,  6350.0    ,  6100.0    ,  5899.0    ,  5690.0    , //15
1041        5300.0    ,  4850.0    ,  4450.0    ,  5650.0    ,  5700.0    ,
1042        4950.0    ,  3850.0    ,  3050.0    ,  3050.0    ,  3460.0    ,
1043        3650.0    ,  3340.0}
1044};
1045// Inelastic cross section
1046G4float G4HadronCrossSections::cnlwin[NCNLW][NELAB] = {
1047       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //1
1048       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
1049       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
1050       0.00000E+00, 0.00000E+00},
1051       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,  1.0000    , //2
1052        10.000    ,  50.000    ,  100.00    ,  200.00    ,  300.00    ,
1053        400.00    ,  600.00    ,  700.00    ,  750.00    ,  700.00    ,
1054        700.00    ,  680.00}  ,
1055       {0.00000E+00, 0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    , //3
1056        50.000    ,  100.00    ,  260.00    ,  450.00    ,  600.00    ,
1057        700.00    ,  800.00    ,  900.00    ,  940.00    ,  900.00    ,
1058        860.00    ,  820.00}  ,
1059        {1.0000    ,  80.000    ,  200.00    ,  320.00    ,  400.00    , //4
1060        520.00    ,  700.00    ,  1000.0    ,  1120.0    ,  1200.0    ,
1061        1200.0    ,  1200.0    ,  1180.0    ,  1130.0    ,  1100.0    ,
1062        1090.0    ,  1080.0}  ,
1063        {1.0000    ,  90.000    ,  220.00    ,  340.00    ,  420.00    , //5
1064        550.00    ,  720.00    ,  1080.0    ,  1300.0    ,  1400.0    ,
1065        1420.0    ,  1420.0    ,  1380.0    ,  1260.0    ,  1190.0    ,
1066        1150.0    ,  1100.0}  ,
1067        {1.0000    ,  95.000    ,  225.00    ,  345.00    ,  425.00    , //6
1068        555.00    ,  750.00    ,  1150.0    ,  1500.0    ,  1680.0    ,
1069        1700.0    ,  1690.0    ,  1550.0    ,  1360.0    ,  1240.0    ,
1070        1180.0    ,  1120.0}  ,
1071        {10.000    ,  140.00    ,  260.00    ,  380.00    ,  450.00    , //7
1072        600.00    ,  750.00    ,  1200.0    ,  1580.0    ,  1800.0    ,
1073        1820.0    ,  1830.0    ,  1800.0    ,  1750.0    ,  1690.0    ,
1074        1650.0    ,  1620.0}  ,
1075        {15.000    ,  150.00    ,  300.00    ,  400.00    ,  500.00    , //8
1076        650.00    ,  840.00    ,  1500.0    ,  2100.0    ,  2130.0    ,
1077        2140.0    ,  2130.0    ,  2080.0    ,  2000.0    ,  1950.0    ,
1078        1900.0    ,  1860.0}  ,
1079       {0.00000E+00, 0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    , //9
1080        150.00    ,  380.00    ,  1000.0    ,  1650.0    ,  2100.0    ,
1081        2100.0    ,  2100.0    ,  2060.0    ,  1950.0    ,  1860.0    ,
1082        1800.0    ,  1740.0}  ,
1083       {0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    ,  45.000    , //10
1084        180.00    ,  380.00    ,  1050.0    ,  1900.0    ,  2300.0    ,
1085        2300.0    ,  2200.0    ,  2150.0    ,  2000.0    ,  1900.0    ,
1086        1800.0    ,  1750.0}  ,
1087       {0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    ,  48.000    , //11
1088        190.00    ,  390.00    ,  1080.0    ,  2000.0    ,  2400.0    ,
1089        2400.0    ,  2300.0    ,  2200.0    ,  2100.0    ,  1950.0    ,
1090        1850.0    ,  1800.0}  ,
1091       {0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    ,  50.000    , //12
1092        200.00    ,  400.00    ,  1100.0    ,  2100.0    ,  2500.0    ,
1093        2500.0    ,  2450.0    ,  2300.0    ,  2100.0    ,  2000.0    ,
1094        1900.0    ,  1850.0}  ,
1095       {0.00000E+00, 0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    , //13
1096        100.00    ,  350.00    ,  900.00    ,  1400.0    ,  2000.0    ,
1097        2300.0    ,  2380.0    ,  2400.0    ,  2300.0    ,  2250.0    ,
1098        2200.0    ,  2120.0}  ,
1099       {0.00000E+00, 0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    , //14
1100        100.00    ,  350.00    ,  900.00    ,  1400.0    ,  2000.0    ,
1101        2300.0    ,  2380.0    ,  2400.0    ,  2300.0    ,  2250.0    ,
1102        2200.0    ,  2120.0}  ,
1103       {0.00000E+00, 0.00000E+00, 0.00000E+00,  1.0000    ,  10.000    , //15
1104        100.00    ,  400.00    ,  950.00    ,  1600.0    ,  2200.0    ,
1105        2550.0    ,  2750.0    ,  2700.0    ,  2600.0    ,  2540.0    ,
1106        2450.0    ,  2360.0}
1107};
1108// Capture cross section indexed by Z
1109G4float G4HadronCrossSections::cscap[100] = {
1110        6.0000    ,  5.7000    ,  5.5000    ,  5.3000    ,  5.2000    ,
1111        5.1000    ,  5.0000    ,  4.9000    ,  4.8000    ,  4.8000    ,
1112        4.8000    ,  4.8000    ,  4.8000    ,  4.8000    ,  4.8000    ,
1113        4.8000    ,  4.9000    ,  5.0000    ,  5.2000    ,  5.5000    ,
1114        6.0000    ,  6.7000    ,  7.5000    ,  8.5000    ,  10.000    ,
1115        12.000    ,  14.500    ,  19.000    ,  26.500    ,  40.000    ,
1116        75.000    ,  120.00    ,  180.00    ,  260.00    ,  360.00    ,
1117        330.00    ,  60.000    ,  7.0000    ,  9.5000    ,  20.000    ,
1118        75.000    ,  140.00    ,  250.00    ,  360.00    ,  480.00    ,
1119        580.00    ,  590.00    ,  500.00    ,  300.00    ,  100.00    ,
1120        200.00    ,  300.00    ,  400.00    ,  470.00    ,  500.00    ,
1121        430.00    ,  100.00    ,  20.000    ,  22.000    ,  40.000    ,
1122        560.00    ,  950.00    ,  1000.0    ,  1000.0    ,  1000.0    ,
1123        990.00    ,  920.00    ,  860.00    ,  790.00    ,  740.00    ,
1124        650.00    ,  600.00    ,  540.00    ,  470.00    ,  440.00    ,
1125        390.00    ,  360.00    ,  340.00    ,  320.00    ,  310.00    ,
1126        280.00    ,  2.0000    ,  2.5000    ,  6.0000    ,  13.000    ,
1127        38.000    ,  65.000    ,  140.00    ,  280.00    ,  300.00    ,
1128        430.00    ,  580.00    ,  650.00    ,  800.00    ,  920.00    ,
1129        1100.0    ,  1250.0    ,  1400.0    ,  1550.0    ,  1700.0   
1130};
1131
1132//---------------------------------------------------------------------
1133// Tables for fission cross sections
1134//---------------------------------------------------------------------
1135
1136// Lab kinetic energy in GeV
1137G4float G4HadronCrossSections::ekfiss[NFISS] = {
1138       0.10000E-03, 0.20000E-03, 0.30000E-03, 0.50000E-03, 0.70000E-03,
1139       0.10000E-02, 0.15000E-02, 0.20000E-02, 0.30000E-02, 0.50000E-02,
1140       0.70000E-02, 0.10000E-01, 0.15000E-01, 0.20000E-01, 0.50000E-01,
1141       0.10000    , 0.20000    , 0.30000    , 0.40000    , 0.50000    ,
1142        1000.0   
1143};
1144// Fission cross sections
1145G4float G4HadronCrossSections::csfiss[4][NFISS] = {
1146        {2600.0    ,  2300.0    ,  2300.0    ,  2100.0    ,  2000.0    , //1
1147        1950.0    ,  1930.0    ,  1900.0    ,  1800.0    ,  1600.0    ,
1148        2100.0    ,  2300.0    , 0.00000E+00, 0.00000E+00, 0.00000E+00,
1149       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
1150       0.00000E+00},
1151        {1850.0    ,  1400.0    ,  1300.0    ,  1150.0    ,  1100.0    , //2
1152        1200.0    ,  1250.0    ,  1300.0    ,  1250.0    ,  1150.0    ,
1153        1600.0    ,  1900.0    , 0.00000E+00, 0.00000E+00, 0.00000E+00,
1154       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
1155       0.00000E+00},
1156        {1700.0    ,  1650.0    ,  1650.0    ,  1700.0    ,  1700.0    , //3
1157        1800.0    ,  1900.0    ,  2000.0    ,  1950.0    ,  1800.0    ,
1158        2150.0    ,  2450.0    , 0.00000E+00, 0.00000E+00, 0.00000E+00,
1159       0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00,
1160       0.00000E+00},
1161       {0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, 0.00000E+00, //4
1162       0.00000E+00,  250.00    ,  550.00    ,  550.00    ,  550.00    ,
1163        550.00    ,  550.00    ,  1000.0    ,  1400.0    ,  1600.0    ,
1164        1500.0    ,  1400.0    ,  1300.0    ,  1200.0    ,  1100.0    ,
1165        1000.0}
1166};
1167
1168G4float G4HadronCrossSections::alpha[NPARTS] = {
1169                     0.7,0.7,0.7,0.7,0.7,0.7,
1170                     0.75,0.75,0.75,
1171                     0.76,0.76,0.76,0.76,
1172                     0.685,0.63,0.685,0.63,0.685,0.63,
1173                     0.685,0.685,0.685,0.63,0.63,0.63,0.685,0.685,0.63,0.63,
1174                     0.7,0.7,0.7,0.685,0.63,0.7
1175};
1176
1177G4float G4HadronCrossSections::alphac[TSIZE] = {
1178                     1.2,1.2,1.2,1.15,0.90,0.91,0.98,1.06,1.10,1.11,
1179                     1.10,1.08,1.05,1.01,0.985,0.962,0.945,0.932,
1180                     0.925,0.920,0.920,0.921,0.922,0.923,0.928,0.931,
1181                     0.940,0.945,0.950,0.955,0.958,0.962,0.965,0.976,
1182                     0.982,0.988,0.992,1.010,1.020,1.030,1.040
1183};
1184
1185G4float G4HadronCrossSections::partel[NPARTS] = {
1186                  0.,0.,0.,0.,0.,0.,
1187                  1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,
1188                  1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,
1189                  1.,1.,1.,1.,1.,1.,1.,1.,1.
1190};
1191G4float G4HadronCrossSections::partin[NPARTS] = {
1192                  0.,0.,0.,0.,0.,0.,
1193                  1.00,0.00,1.05,1.20,1.35,1.30,1.20,1.00,1.30,
1194                  1.00,1.30,1.00,1.30,1.00,1.00,1.00,1.30,1.30,1.30,
1195                  1.00,1.00,1.30,1.30,1.00,1.,1.,1.,1.3,1.
1196};
1197
1198// Enabling flags for corrections for compounds
1199G4int G4HadronCrossSections::icorr[NPARTS] = {
1200                  1,1,1,1,1,1,1,1,1,1,1,1,1,1,
1201                  0,1,0,1,0,1,1,1,0,0,0,1,1,0,0,1,1,1,1,0,0
1202};
1203
1204// Enabling flags for interaction to occur
1205G4int G4HadronCrossSections::intrc[NPARTS] = {
1206                  0,0,0,0,0,0,
1207                  1,0,1,1,1,1,1,1,1,1,1,1,1,1,0,2,2,0,
1208                  1,1,1,1,1,1,1,1,1,1,0
1209};
1210
1211G4float G4HadronCrossSections::csa[4] = {1., 27.00, 63.54, 207.19};
1212
1213G4int G4HadronCrossSections::ipart2[7] = {9, 8, 7, 11, 10, 13, 12};
1214
1215G4bool G4HadronCrossSections::correctInelasticNearZero = 0;
1216
1217
1218G4double G4HadronCrossSections::GetInelasticCrossSection(
1219                          const G4DynamicParticle* particle,
1220                          const G4Element* element)
1221{
1222  if (particle->GetDefinition() != prevParticleDefinition ||
1223      element != prevElement ||
1224      particle->GetKineticEnergy() != prevKineticEnergy) {
1225
1226    G4int nIso = element->GetNumberOfIsotopes();
1227
1228    if (nIso) {
1229      G4double cross_section = 0;
1230      G4IsotopeVector* isoVector = element->GetIsotopeVector();
1231      G4double* abundVector = element->GetRelativeAbundanceVector();
1232      G4double ZZ;
1233      G4double AA;
1234
1235      for (G4int i = 0; i < nIso; i++) {
1236        ZZ = G4double( (*isoVector)[i]->GetZ() );
1237        AA = G4double( (*isoVector)[i]->GetN() );
1238        CalcScatteringCrossSections(particle, ZZ, AA);
1239        cross_section += siginelastic*abundVector[i];
1240      }
1241      siginelastic = cross_section;
1242
1243    } else { 
1244      CalcScatteringCrossSections(particle, element->GetZ(), element->GetN());
1245    }
1246  }
1247  return siginelastic;
1248}
1249
1250
1251G4double G4HadronCrossSections::GetInelasticCrossSection(
1252                          const G4DynamicParticle* particle,
1253                          G4double ZZ, G4double AA)
1254{
1255  prevElement = 0; // force new cross section calculation for next call of
1256                   // GetInelasticCrossSection(particle, element, temp)
1257  CalcScatteringCrossSections(particle, ZZ, AA);
1258  return siginelastic;
1259}
1260
1261
1262G4double G4HadronCrossSections::GetElasticCrossSection(
1263                          const G4DynamicParticle* particle,
1264                          const G4Element* element)
1265{
1266  if (particle->GetDefinition() != prevParticleDefinition ||
1267      element != prevElement ||
1268      particle->GetKineticEnergy() != prevKineticEnergy) {
1269
1270    G4int nIso = element->GetNumberOfIsotopes();
1271
1272    if (nIso) {
1273      G4double cross_section = 0;
1274      G4IsotopeVector* isoVector = element->GetIsotopeVector();
1275      G4double* abundVector = element->GetRelativeAbundanceVector();
1276      G4double ZZ;
1277      G4double AA;
1278
1279      for (G4int i = 0; i < nIso; i++) {
1280        ZZ = G4double( (*isoVector)[i]->GetZ() );
1281        AA = G4double( (*isoVector)[i]->GetN() );
1282        CalcScatteringCrossSections(particle, ZZ, AA);
1283        cross_section += sigelastic*abundVector[i];
1284      }
1285      sigelastic = cross_section;
1286
1287    } else { 
1288      CalcScatteringCrossSections(particle, element->GetZ(), element->GetN());
1289    }
1290  }
1291  return sigelastic;
1292}
1293
1294
1295G4double G4HadronCrossSections::GetElasticCrossSection(
1296                          const G4DynamicParticle* particle,
1297                          G4double ZZ, G4double AA)
1298{
1299  prevElement = 0; // force new cross section calculation for next call of
1300                   // GetElasticCrossSection(particle, element, temp)
1301  CalcScatteringCrossSections(particle, ZZ, AA);
1302  return sigelastic;
1303}
1304
1305
1306// Method to calculate cross sections for all processes.
1307// To facilitate comparison with the original Fortran source, and because
1308// of the interdependence of the elastic and inelastic cross section
1309// calculations, this has not been split into separate processes.
1310
1311void
1312G4HadronCrossSections::CalcScatteringCrossSections(
1313                          const G4DynamicParticle* aParticle,
1314                          G4double ZZ, G4double AA)
1315{
1316   G4double sigel, sigin, sigtot;
1317   G4double xsecel, xsecin=0;
1318   xsecel = 0;
1319   G4double dx, dy, rc, rca, b;
1320   G4double crel, crin;
1321   G4double xspiel, xspiin;
1322
1323   G4int ipart = GetParticleCode(aParticle);
1324   G4double a = AA;
1325   G4double z = ZZ;
1326
1327   if (verboseLevel > 1) {
1328      G4cout << "G4HadronCrossSections: a=" << a << G4endl;
1329      G4cout << "G4HadronCrossSections: z=" << z << G4endl;
1330   }
1331
1332// Ions...
1333
1334   if (ipart >= 30 && ipart <= 32) {
1335
1336      G4double apart=0;
1337      if (ipart == 30) apart = std::pow(2., 1./3.);
1338      else if (ipart == 31) apart = std::pow(3., 1./3.);
1339      else if (ipart == 32) apart = std::pow(4., 1./3.);
1340
1341      G4double term = apart + std::pow(a, 1./3.);
1342      sigin = 49.*term*term;
1343   // Convert cross section from mb to default units
1344      siginelastic = sigin*millibarn;
1345      if(aParticle->GetKineticEnergy() < 6*MeV) siginelastic = 0;
1346      sigelastic = 0.;
1347      return;
1348   }
1349
1350   G4double ek = aParticle->GetKineticEnergy()/GeV;
1351
1352// Low energy neutrons...
1353
1354   if (ipart == 16 && ek <= 0.0327) {
1355
1356     //      G4int je2 = NELAB;
1357     //      for (G4int j = 2; j <= NELAB; j++) {
1358     //         if (ek < elab[j - 1]) {
1359     //            je2 = j;
1360     //            break;
1361     //         }
1362     //      }
1363     //      G4int je1 = je2 -1;
1364     //      je1 = je1 - 1;      // For array indexing
1365     //      je2 = je2 - 1;      // For array indexing
1366
1367      G4int je1 = 0;
1368      G4int je2 = NELAB - 1;
1369      do {
1370         G4int midBin = (je1 + je2)/2;
1371         if (ek < elab[midBin])
1372           je2 = midBin;
1373         else
1374           je1 = midBin;
1375      } while (je2 - je1 > 1); 
1376
1377      G4double delab = elab[je2] - elab[je1];
1378
1379      //      G4int ja2 = NCNLW;
1380      //      for (G4int jj = 2; jj <= NCNLW; jj++) {
1381      //         if (a < cnlwat[jj - 1]) {
1382      //            ja2 = jj;
1383      //            break;
1384      //         }
1385      //      }
1386      //      G4int ja1 = ja2 - 1;
1387      //      ja1 = ja1 - 1;      // For array indexing
1388      //      ja2 = ja2 - 1;      // For array indexing
1389
1390      G4int ja1 = 0;
1391      G4int ja2 = NCNLW - 1;
1392      do {
1393         G4int midBin = (ja1 + ja2)/2;
1394         if (a < cnlwat[midBin])
1395           ja2 = midBin;
1396         else
1397           ja1 = midBin;
1398      } while (ja2 - ja1 > 1); 
1399
1400      G4double dnlwat = cnlwat[ja2] - cnlwat[ja1];
1401
1402// Elastic cross section:
1403      // E interpolation or extrapolation at JA1
1404      dy = cnlwel[ja1][je2] - cnlwel[ja1][je1];
1405      G4double rce = dy/delab;
1406      // A interpolation or extrapolation at JE1
1407      dy = cnlwel[ja2][je1] - cnlwel[ja1][je1];
1408      rca = dy/dnlwat;
1409      b = cnlwel[ja1][je1] - rce*elab[je1] - rca*cnlwat[ja1];
1410      sigelastic = rce*ek + rca*a + b;
1411// Inelastic cross section:
1412      // E interpolation or extrapolation at JA1
1413      dy = cnlwin[ja1][je2] - cnlwin[ja1][je1];
1414      rce = dy/delab;
1415      // A interpolation or extrapolation at JE1
1416      dy = cnlwin[ja2][je1] - cnlwin[ja1][je1];
1417      rca = dy/dnlwat;
1418      b = cnlwin[ja1][je1] - rce*elab[je1] - rca*cnlwat[ja1];
1419      siginelastic = rce*ek + rca*a + b;
1420   // Convert cross sections from mb to default units
1421      sigelastic = sigelastic*millibarn;
1422      siginelastic = siginelastic*millibarn;
1423      return;
1424   }
1425
1426// Remaining particles...
1427
1428// Get momentum bin
1429   G4double p = aParticle->GetTotalMomentum()/GeV;
1430
1431   //   G4int j = TSIZE - 1;
1432   //   for (G4int i = 2; i <= TSIZE; i++) {
1433   //      if (p < plab[i]) {
1434   //         j = i - 1;
1435   //         break;
1436   //      }
1437   //   }
1438   //   j = j - 1;                     // For array indexing
1439
1440   G4int je1 = 0;
1441   G4int je2 = TSIZE - 1;
1442   do {
1443      G4int midBin = (je1 + je2)/2;
1444      if (p < plab[midBin])
1445        je2 = midBin;
1446      else
1447        je1 = midBin;
1448   } while (je2 - je1 > 1); 
1449
1450   G4int ipart1 = ipart - 1;      // For array indexing
1451
1452// Get cross sections for scattering on free protons
1453   dx = plab[je2] - plab[je1];
1454// Elastic cross section
1455   dy = csel[ipart1][je2] - csel[ipart1][je1];
1456   rc = dy/dx;
1457   b = csel[ipart1][je1] - rc*plab[je1];
1458   sigel = rc*p + b;
1459// Inelastic cross section
1460   dy = csin[ipart1][je2] - csin[ipart1][je1];
1461   rc = dy/dx;
1462   b = csin[ipart1][je1] - rc*plab[je1];
1463   sigin = rc*p + b;
1464   if (verboseLevel > 1) {
1465      G4cout << "sigel " << sigel << G4endl;
1466      G4cout << "sigin " << sigin << G4endl;
1467   }
1468   G4double alph = alpha[ipart1];
1469   if (ipart < 14) {
1470      dy = alphac[je2] - alphac[je1];
1471      rc = dy/dx;
1472      b = alphac[je1] - rc*plab[je1];
1473      G4double corfac = rc*p + b;
1474      alph = alph*corfac;
1475      G4int ipart3 = ipart2[ipart - 7];
1476      ipart3 = ipart3 - 1;      // For array indexing
1477      // Elastic cross section
1478      dy = csel[ipart3][je2] - csel[ipart3][je1];
1479      rc = dy/dx;
1480      b = csel[ipart3][je1] - rc*plab[je1];
1481      xsecel = rc*p + b;
1482      // Inelastic cross section
1483      dy = csin[ipart3][je2] - csin[ipart3][je1];
1484      rc = dy/dx;
1485      b = csin[ipart3][je1] - rc*plab[je1];
1486      xsecin = rc*p + b;
1487   }
1488
1489// A-dependence from parameterization...
1490
1491   if (a >= 1.5) {
1492
1493      crel = 1.;
1494      crin = 1.;
1495
1496      G4int i = 3;
1497      if (a < 50.) i = 2;
1498      if (a > 100.) i = 4;
1499      i = i - 1;      // For array indexing
1500
1501// Protons and neutrons
1502      if (ipart == 14 || ipart == 16) {
1503         dy = cspnel[i - 1][je2] - cspnel[i - 1][je1];
1504         rc = dy/dx;
1505         b = cspnel[i - 1][je1] - rc*plab[je1];
1506         xsecel = rc*p + b;
1507         dy = cspnin[i - 1][je2] - cspnin[i - 1][je1];
1508         rc = dy/dx;
1509         b = cspnin[i - 1][je1] - rc*plab[je1];
1510         xsecin = rc*p + b;
1511         // The following is a first-order correction to Gheisha (GHESIG)
1512         // behaviour, where for the lighter elements the inelastic cross
1513         // section is not realistic for particles in the first momentum
1514         // bin.  In the first momentum bin, it better to interpolate
1515         // in K.E.  Subject to further improvements.
1516         if (correctInelasticNearZero && je1 == 0 && i <= 3) {
1517            G4double m0 = aParticle->GetMass()/GeV;
1518            G4double T = std::sqrt(m0*m0 + p*p) - m0;
1519            G4double dx = std::sqrt(m0*m0 + plab[1]*plab[1]) - m0;
1520            rc = dy/dx;
1521            xsecin = rc*T + b;
1522         }
1523
1524         if (sigel >= 0.001) 
1525            crel = xsecel/(0.36*sigel*std::pow(G4double(csa[i]), 1.17));
1526         sigtot = sigel + sigin;
1527         if (sigtot >= 0.001) 
1528            crin = xsecin/(sigtot*std::pow(G4double(csa[i]), alph));
1529      }
1530
1531      else if (ipart < 15) {
1532// Calculate correction factors (crel, crin) from values
1533// on Al, Cu, Pb.  Note that data is only available for pions and protons.
1534         G4double wgch = 0.5;
1535         if (a < 20.) wgch = 0.5 + 0.5*std::exp(-(a - 1.));
1536         sigel = wgch*sigel + (1. - wgch)*xsecel;
1537         sigin = wgch*sigin + (1. - wgch)*xsecin;
1538         
1539// This section not for kaons
1540         if (ipart < 10) {
1541            dy = cspiel[i - 1][je2] - cspiel[i - 1][je1];
1542            rc = dy/dx;
1543            b = cspiel[i - 1][je1] - rc*plab[je1];
1544            xspiel = rc*p + b;
1545            dy = cspiin[i - 1][je2] - cspiin[i - 1][je1];
1546            rc = dy/dx;
1547            b = cspiin[i - 1][je1] - rc*plab[je1];
1548            xspiin = rc*p + b;
1549            if (verboseLevel > 1) {
1550               G4cout << "xspiel " << xspiel << G4endl;
1551               G4cout << "xspiin " << xspiin << G4endl;
1552            }               
1553            if (sigel >= 0.001) 
1554               crel = xspiel/(0.36*sigel*std::pow(G4double(csa[i]),1.17));
1555            sigtot = sigel + sigin;
1556            if (sigtot >= 0.001) 
1557               crin = xspiin/(sigtot*std::pow(G4double(csa[i]), alph));
1558         }
1559      }
1560
1561// Apply correction factors
1562      sigin = crin*(sigin + sigel)*std::pow(a, alph);
1563      sigel = crel*0.36*sigel*std::pow(a, 1.17);
1564      sigel = sigel*partel[ipart1];
1565      sigin = sigin*partin[ipart1];
1566   }
1567
1568// Correction factor for high (p > 100 GeV/c) energies:
1569   G4double corh = 1.;
1570   if (p > 100.) corh = 0.1085736156*std::log(p) + 0.5;
1571
1572   sigel = corh*sigel;
1573   sigin = corh*sigin;
1574   // Convert cross section from mb to default units
1575   sigelastic = sigel*millibarn;
1576   siginelastic = sigin*millibarn;
1577
1578   return;
1579}
1580   
1581
1582G4double
1583G4HadronCrossSections::GetCaptureCrossSection(
1584                          const G4DynamicParticle* aParticle,
1585                          G4double ZZ, G4double /*AA*/)
1586{
1587   if (GetParticleCode(aParticle) != 16)  { return 0.; }
1588   G4double ek = aParticle->GetKineticEnergy()/GeV;
1589   if (ek > 0.0327)  { return 0.; }
1590
1591   G4double ekx = std::max(ek, 1.e-9);
1592   if( ekx != lastEkx )
1593   {
1594     lastEkx = ekx;
1595     lastEkxPower = std::pow(ekx*1.e6, 0.577);
1596   }
1597
1598   G4int izno = static_cast<G4int> (ZZ + 0.01);
1599   if (izno > 100) izno = 100;      // Not in GHESIG
1600   izno = izno - 1;      // For array indexing
1601   G4double sigcap = 11.12*cscap[izno]/lastEkxPower;
1602   // Convert cross section from mb to default units
1603   sigcap = sigcap*millibarn;
1604   return sigcap;
1605}
1606
1607
1608G4double
1609G4HadronCrossSections::GetFissionCrossSection(
1610                          const G4DynamicParticle* aParticle,
1611                          G4double ZZ, G4double AA)
1612{
1613   if (AA < 230.) return 0;
1614
1615   G4double ek = aParticle->GetKineticEnergy()/GeV;     
1616
1617   //   G4int i = NFISS;
1618   //   for (G4int ii = 1; i <= NFISS; i++) {
1619   //      if (ek < ekfiss[ii - 1]) {
1620   //         i = ii;
1621   //         break;
1622   //      }
1623   //   }
1624   //   i = i - 1;      // For array indexing
1625
1626   G4int ie1 = 0;
1627   G4int ie2 = NFISS - 1;
1628   do {
1629      G4int midBin = (ie1 + ie2)/2;
1630      if (ek < ekfiss[midBin])
1631        ie2 = midBin;
1632      else
1633        ie1 = midBin;
1634   } while (ie2 - ie1 > 1); 
1635   G4int i = ie2;
1636   if (ek < ekfiss[0]) i = 0;
1637
1638   G4int j = 4;
1639   if (ek <= 0.01) {
1640      if (ZZ == 92. && std::abs(AA - 233.) < 0.5) j = 1;
1641      else if (ZZ == 92. && std::abs(AA - 235.) < 0.5) j = 2;
1642      else if (ZZ == 94. && std::abs(AA - 239.) < 0.5) j = 3;
1643   }
1644   G4double z43ba;
1645   if (j == 4) {
1646      z43ba = std::pow(ZZ, 4./3.)/AA;
1647      z43ba = std::max(-67. + 38.7*z43ba, 0.);
1648   }
1649   else {
1650      z43ba = 1.;
1651   }
1652   j = j - 1;      // For array indexing
1653   G4double sigfiss = csfiss[j][i]*z43ba;
1654// Convert from mb to cm^2
1655   //         sigfiss = sigfiss*1.e-27;
1656   // Convert cross section from mb to default units
1657   sigfiss = sigfiss*millibarn;
1658   return sigfiss;
1659}
1660
1661
1662G4int
1663G4HadronCrossSections::GetParticleCode(const G4DynamicParticle* aParticle)
1664{
1665  // Returns GHEISHA code for particle
1666  // Case entries ordered by estimated frequency
1667
1668  G4int ipart;
1669
1670  switch( aParticle->GetPDGcode()) {
1671    case 111:
1672      ipart = 8;   // pi0
1673      break;
1674    case 211:
1675      ipart = 7;   // pi+
1676      break;
1677    case -211:
1678      ipart = 9;   // pi-
1679      break;
1680    case 2112:
1681      ipart = 16;  // neutron
1682      break;
1683    case 2212:
1684      ipart = 14;  // proton
1685      break;
1686    case 321:       
1687      ipart = 10;  // K+
1688      break;
1689    case -321:
1690      ipart = 13;  // K-
1691      break;
1692    case 130:
1693      ipart = 12;  // K0L
1694      break;
1695    case 310:     
1696      ipart = 11;  // K0S
1697      break;
1698    case 1000010020:
1699      ipart = 30;  // deuteron
1700      break;
1701    case 1000010030:
1702      ipart = 31;  // triton
1703      break;
1704    case 1000020040:
1705      ipart = 32;  // alpha
1706      break;
1707    case 3122:
1708      ipart = 18;  // lambda
1709      break;
1710    case -2112:
1711      ipart = 17;  // anti-neutron
1712      break;
1713    case -2212:
1714      ipart = 15;  // anti-proton
1715      break;
1716    case -3122:
1717      ipart = 19;  // anti-lambda
1718      break;
1719    case 3222:
1720      ipart = 20;  // sigma+
1721      break;
1722    case 3212:
1723      ipart = 21;  // sigma0
1724      break;
1725    case 3112:
1726      ipart = 22;  // sigma-
1727      break;
1728    case 3322:
1729      ipart = 26;  // xi0
1730      break;
1731    case 3312:
1732      ipart = 27;  // xi-
1733      break;
1734    case 3334:
1735      ipart = 33;  // omega-
1736      break;
1737    case -3222:
1738      ipart = 23;  // anti-sigma+
1739      break;
1740    case -3212:
1741      ipart = 24;  // anti-sigma0
1742      break;
1743    case -3112:
1744      ipart = 25;  // anti-sigma-
1745      break;
1746    case -3322:
1747      ipart = 28;  // anti-xi0
1748      break;
1749    case -3312:
1750      ipart = 29;  // anti-xi-
1751      break;
1752    case -3334:
1753      ipart = 34;  // anti-omega-
1754      break;
1755    default:
1756      G4Exception("G4HadronCrossSections", "007", FatalException,
1757                  "GetParticleCode: unsupported particle");
1758      return 0;
1759  }
1760
1761  return ipart;
1762}
Note: See TracBrowser for help on using the repository browser.