Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4XNNstarTable.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 
28 #include "globals.hh"
29 #include "G4ios.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "G4XNNstarTable.hh"
32 #include "G4PhysicsFreeVector.hh"
33 
34 const G4int G4XNNstarTable::sizeNNstar = 121;
35 
36 // Energies (GeV) corresponding to the cross section table
37 // Units are assigned when filling the PhysicsVector
38 
39 const G4double G4XNNstarTable::energyTable[121] =
40 {
41  0.0,
42  2.014, 2.014, 2.016, 2.018, 2.022, 2.026, 2.031, 2.037, 2.044, 2.052,
43  2.061, 2.071, 2.082, 2.094, 2.107, 2.121, 2.135, 2.151, 2.168, 2.185,
44  2.204, 2.223, 2.244, 2.265, 2.287, 2.311, 2.335, 2.360, 2.386, 2.413,
45  2.441, 2.470, 2.500, 2.531, 2.562, 2.595, 2.629, 2.664, 2.699, 2.736,
46  2.773, 2.812, 2.851, 2.891, 2.933, 2.975, 3.018, 3.062, 3.107, 3.153,
47  3.200, 3.248, 3.297, 3.347, 3.397, 3.449, 3.502, 3.555, 3.610, 3.666,
48  3.722, 3.779, 3.838, 3.897, 3.957, 4.018, 4.081, 4.144, 4.208, 4.273,
49  4.339, 4.406, 4.473, 4.542, 4.612, 4.683, 4.754, 4.827, 4.900, 4.975,
50  5.000, 6.134, 7.269, 8.403, 9.538, 10.672, 11.807, 12.941, 14.076, 15.210,
51  16.345, 17.479, 18.613, 19.748, 20.882, 22.017, 23.151, 24.286, 25.420, 26.555,
52  27.689, 28.824, 29.958, 31.092, 32.227, 33.361, 34.496, 35.630, 36.765, 37.899,
53  39.034, 40.168, 41.303, 42.437, 43.571, 44.706, 45.840, 46.975, 48.109, 49.244
54 };
55 
56 // Cross-sections in mb, from S.A. Bass et al., Prog.Part.Nucl.Phys.41:225-370,1998
57 // Units are assigned when filling the PhysicsVector
58 
59 const G4double G4XNNstarTable::sigmaNN1440[121] =
60 {
61  0.0,
62  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
63  0.000, 0.000, 0.000, 0.001, 0.001, 0.002, 0.003, 0.004, 0.006, 0.009,
64  0.014, 0.020, 0.031, 0.047, 0.072, 0.113, 0.178, 0.266, 0.370, 0.476,
65  0.575, 0.665, 0.744, 0.814, 0.874, 0.926, 0.971, 1.009, 1.040, 1.066,
66  1.087, 1.103, 1.115, 1.124, 1.129, 1.132, 1.131, 1.129, 1.124, 1.117,
67  1.109, 1.099, 1.088, 1.075, 1.062, 1.048, 1.033, 1.017, 1.001, 0.984,
68  0.967, 0.950, 0.932, 0.914, 0.896, 0.879, 0.861, 0.843, 0.825, 0.807,
69  0.790, 0.773, 0.755, 0.738, 0.722, 0.705, 0.689, 0.673, 0.657, 0.641,
70  0.636, 0.453, 0.336, 0.258, 0.204, 0.166, 0.137, 0.115, 0.098, 0.084,
71  0.073, 0.064, 0.057, 0.051, 0.046, 0.041, 0.037, 0.034, 0.031, 0.029,
72  0.026, 0.024, 0.022, 0.021, 0.019, 0.018, 0.017, 0.016, 0.015, 0.014,
73  0.013, 0.013, 0.012, 0.011, 0.011, 0.010, 0.010, 0.009, 0.009, 0.008
74 };
75 
76 const G4double G4XNNstarTable::sigmaNN1520[121] =
77 {
78  0.0,
79  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
80  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.001, 0.001,
81  0.001, 0.002, 0.003, 0.005, 0.008, 0.014, 0.024, 0.043, 0.079, 0.150,
82  0.280, 0.465, 0.662, 0.841, 0.995, 1.125, 1.235, 1.327, 1.410, 1.474,
83  1.525, 1.567, 1.600, 1.625, 1.644, 1.651, 1.657, 1.659, 1.657, 1.651,
84  1.642, 1.629, 1.614, 1.597, 1.578, 1.558, 1.536, 1.512, 1.488, 1.463,
85  1.437, 1.411, 1.384, 1.357, 1.329, 1.302, 1.274, 1.247, 1.219, 1.192,
86  1.165, 1.138, 1.112, 1.086, 1.060, 1.035, 1.010, 0.985, 0.961, 0.938,
87  0.930, 0.652, 0.479, 0.365, 0.287, 0.232, 0.191, 0.160, 0.135, 0.116,
88  0.101, 0.089, 0.078, 0.070, 0.062, 0.056, 0.051, 0.046, 0.042, 0.039,
89  0.036, 0.033, 0.031, 0.028, 0.026, 0.025, 0.023, 0.022, 0.020, 0.019,
90  0.018, 0.017, 0.016, 0.015, 0.015, 0.014, 0.013, 0.013, 0.012, 0.011
91 };
92 
93 const G4double G4XNNstarTable::sigmaNN1535[121] =
94 {
95  0.0,
96  0.000, 0.000, 0.001, 0.001, 0.001, 0.002, 0.002, 0.003,
97  0.004, 0.005, 0.006, 0.008, 0.010, 0.012, 0.015, 0.019,
98  0.024, 0.031, 0.039, 0.052, 0.069, 0.097, 0.145, 0.216,
99  0.298, 0.378, 0.451, 0.513, 0.566, 0.610, 0.646, 0.675,
100  0.699, 0.718, 0.732, 0.742, 0.749, 0.753, 0.754, 0.753,
101  0.751, 0.746, 0.740, 0.733, 0.724, 0.715, 0.705, 0.694,
102  0.683, 0.671, 0.659, 0.647, 0.634, 0.621, 0.608, 0.595,
103  0.582, 0.569, 0.556, 0.543, 0.531, 0.518, 0.506, 0.493,
104  0.481, 0.469, 0.458, 0.446, 0.435, 0.424, 0.413, 0.402,
105  0.399, 0.276, 0.201, 0.153, 0.120, 0.096, 0.079, 0.066,
106  0.056, 0.048, 0.042, 0.037, 0.032, 0.029, 0.026, 0.023,
107  0.021, 0.019, 0.017, 0.016, 0.015, 0.014, 0.013, 0.012,
108  0.011, 0.010, 0.009, 0.009, 0.008, 0.008, 0.007, 0.007,
109  0.007, 0.006, 0.006, 0.006, 0.005, 0.005, 0.005, 0.005
110 };
111 
112 const G4double G4XNNstarTable::sigmaNN1650[121] =
113 {
114  0.0,
115  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
116  0.000, 0.000, 0.000, 0.000, 0.001, 0.001, 0.001, 0.001,
117  0.002, 0.002, 0.002, 0.003, 0.004, 0.005, 0.006, 0.007,
118  0.008, 0.010, 0.012, 0.015, 0.019, 0.023, 0.029, 0.038,
119  0.051, 0.071, 0.104, 0.150, 0.201, 0.249, 0.290, 0.326,
120  0.354, 0.378, 0.397, 0.412, 0.424, 0.434, 0.440, 0.445,
121  0.448, 0.449, 0.449, 0.448, 0.445, 0.442, 0.438, 0.433,
122  0.428, 0.422, 0.416, 0.409, 0.403, 0.395, 0.388, 0.381,
123  0.373, 0.366, 0.358, 0.350, 0.343, 0.335, 0.327, 0.320,
124  0.312, 0.305, 0.298, 0.291, 0.284, 0.277, 0.270, 0.263,
125  0.261, 0.182, 0.133, 0.102, 0.080, 0.064, 0.053, 0.044,
126  0.037, 0.032, 0.028, 0.024, 0.022, 0.019, 0.017, 0.015,
127  0.014, 0.013, 0.012, 0.011, 0.010, 0.009, 0.008, 0.008,
128  0.007, 0.007, 0.006, 0.006, 0.006, 0.005, 0.005, 0.005,
129  0.004, 0.004, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003
130 };
131 
132 const G4double G4XNNstarTable::sigmaNN1675[121] =
133 {
134  0.0,
135  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
136  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
137  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.001,
138  0.001, 0.001, 0.002, 0.002, 0.003, 0.005, 0.009, 0.015,
139  0.026, 0.048, 0.095, 0.189, 0.324, 0.463, 0.586, 0.691,
140  0.780, 0.855, 0.919, 0.972, 1.016, 1.052, 1.081, 1.105,
141  1.123, 1.136, 1.145, 1.151, 1.153, 1.153, 1.149, 1.144,
142  1.136, 1.127, 1.116, 1.104, 1.090, 1.076, 1.061, 1.045,
143  1.028, 1.011, 0.993, 0.975, 0.957, 0.939, 0.921, 0.902,
144  0.884, 0.865, 0.847, 0.828, 0.810, 0.792, 0.775, 0.757,
145  0.751, 0.538, 0.399, 0.307, 0.242, 0.196, 0.162, 0.136,
146  0.115, 0.099, 0.086, 0.076, 0.067, 0.060, 0.053, 0.048,
147  0.044, 0.040, 0.036, 0.033, 0.031, 0.028, 0.026, 0.024,
148  0.023, 0.021, 0.020, 0.019, 0.018, 0.016, 0.016, 0.015,
149  0.014, 0.013, 0.013, 0.012, 0.011, 0.011, 0.010, 0.010
150 };
151 
152 const G4double G4XNNstarTable::sigmaNN1680[121] =
153 {
154  0.0,
155  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
156  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
157  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
158  0.000, 0.001, 0.001, 0.001, 0.002, 0.003, 0.006, 0.010,
159  0.018, 0.035, 0.073, 0.156, 0.294, 0.446, 0.580, 0.693,
160  0.788, 0.867, 0.933, 0.988, 1.033, 1.070, 1.100, 1.124,
161  1.142, 1.155, 1.163, 1.168, 1.170, 1.168, 1.164, 1.158,
162  1.149, 1.139, 1.127, 1.114, 1.100, 1.085, 1.068, 1.052,
163  1.034, 1.016, 0.998, 0.979, 0.960, 0.941, 0.922, 0.903,
164  0.884, 0.865, 0.846, 0.827, 0.809, 0.790, 0.772, 0.754,
165  0.748, 0.533, 0.394, 0.301, 0.238, 0.192, 0.158, 0.133,
166  0.113, 0.097, 0.084, 0.074, 0.065, 0.058, 0.052, 0.047,
167  0.042, 0.039, 0.035, 0.032, 0.030, 0.028, 0.026, 0.024,
168  0.022, 0.021, 0.019, 0.018, 0.017, 0.016, 0.015, 0.014,
169  0.013, 0.013, 0.012, 0.012, 0.011, 0.010, 0.010, 0.010
170 };
171 
172 const G4double G4XNNstarTable::sigmaNN1700[121] =
173 {
174  0.0,
175  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
176  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
177  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
178  0.000, 0.001, 0.001, 0.001, 0.002, 0.003, 0.005, 0.008,
179  0.013, 0.022, 0.038, 0.070, 0.140, 0.245, 0.347, 0.431,
180  0.500, 0.556, 0.601, 0.637, 0.666, 0.689, 0.708, 0.721,
181  0.731, 0.738, 0.742, 0.743, 0.743, 0.740, 0.736, 0.730,
182  0.723, 0.716, 0.707, 0.697, 0.687, 0.676, 0.665, 0.654,
183  0.642, 0.630, 0.617, 0.605, 0.593, 0.580, 0.568, 0.555,
184  0.543, 0.531, 0.519, 0.507, 0.495, 0.483, 0.471, 0.460,
185  0.456, 0.322, 0.236, 0.180, 0.142, 0.114, 0.094, 0.079,
186  0.067, 0.057, 0.050, 0.044, 0.039, 0.034, 0.031, 0.028,
187  0.025, 0.023, 0.021, 0.019, 0.018, 0.016, 0.015, 0.014,
188  0.013, 0.012, 0.011, 0.011, 0.010, 0.009, 0.009, 0.008,
189  0.008, 0.008, 0.007, 0.007, 0.007, 0.006, 0.006, 0.006
190 };
191 
192 const G4double G4XNNstarTable::sigmaNN1710[121] =
193 {
194  0.0,
195  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
196  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
197  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
198  0.000, 0.000, 0.000, 0.001, 0.001, 0.001, 0.002, 0.003,
199  0.005, 0.008, 0.013, 0.025, 0.052, 0.096, 0.144, 0.185,
200  0.219, 0.246, 0.269, 0.288, 0.304, 0.316, 0.326, 0.334,
201  0.340, 0.344, 0.347, 0.349, 0.349, 0.349, 0.348, 0.346,
202  0.343, 0.340, 0.336, 0.332, 0.328, 0.323, 0.319, 0.313,
203  0.308, 0.303, 0.297, 0.292, 0.286, 0.280, 0.275, 0.269,
204  0.263, 0.257, 0.252, 0.246, 0.241, 0.235, 0.230, 0.224,
205  0.223, 0.158, 0.117, 0.090, 0.071, 0.057, 0.047, 0.040,
206  0.034, 0.029, 0.025, 0.022, 0.019, 0.017, 0.016, 0.014,
207  0.013, 0.012, 0.011, 0.010, 0.009, 0.008, 0.008, 0.007,
208  0.007, 0.006, 0.006, 0.005, 0.005, 0.005, 0.005, 0.004,
209  0.004, 0.004, 0.004, 0.003, 0.003, 0.003, 0.003, 0.003
210 };
211 
212 const G4double G4XNNstarTable::sigmaNN1720[121] =
213 {
214  0.0,
215  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
216  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
217  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
218  0.001, 0.001, 0.002, 0.002, 0.003, 0.005, 0.007, 0.010,
219  0.015, 0.023, 0.036, 0.061, 0.106, 0.174, 0.250, 0.321,
220  0.382, 0.434, 0.477, 0.513, 0.543, 0.568, 0.589, 0.605,
221  0.617, 0.627, 0.634, 0.639, 0.641, 0.642, 0.641, 0.639,
222  0.635, 0.630, 0.625, 0.618, 0.611, 0.603, 0.595, 0.586,
223  0.576, 0.567, 0.557, 0.547, 0.537, 0.527, 0.516, 0.506,
224  0.496, 0.485, 0.475, 0.465, 0.455, 0.444, 0.435, 0.425,
225  0.421, 0.302, 0.224, 0.172, 0.136, 0.117, 0.091, 0.076,
226  0.065, 0.056, 0.049, 0.043, 0.038, 0.034, 0.030, 0.027,
227  0.025, 0.022, 0.020, 0.019, 0.017, 0.016, 0.015, 0.014,
228  0.013, 0.012, 0.011, 0.010, 0.010, 0.009, 0.009, 0.008,
229  0.008, 0.007, 0.007, 0.007, 0.006, 0.006, 0.006, 0.006
230 };
231 
232 const G4double G4XNNstarTable::sigmaNN1900[121] =
233 {
234  0.0,
235  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
236  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
237  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.001,
238  0.001, 0.001, 0.001, 0.002, 0.002, 0.003, 0.003, 0.005,
239  0.006, 0.008, 0.010, 0.014, 0.019, 0.026, 0.037, 0.054,
240  0.074, 0.094, 0.114, 0.131, 0.147, 0.161, 0.173, 0.184,
241  0.193, 0.201, 0.208, 0.213, 0.218, 0.221, 0.224, 0.226,
242  0.228, 0.229, 0.229, 0.229, 0.229, 0.228, 0.227, 0.225,
243  0.223, 0.221, 0.219, 0.217, 0.214, 0.212, 0.209, 0.206,
244  0.203, 0.200, 0.197, 0.194, 0.190, 0.187, 0.184, 0.181,
245  0.180, 0.137, 0.106, 0.083, 0.067, 0.056, 0.047, 0.039,
246  0.034, 0.029, 0.026, 0.023, 0.020, 0.018, 0.016, 0.015,
247  0.013, 0.012, 0.011, 0.010, 0.010, 0.009, 0.008, 0.008,
248  0.007, 0.007, 0.006, 0.006, 0.006, 0.005, 0.005, 0.005,
249  0.004, 0.004, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003
250 };
251 
252 const G4double G4XNNstarTable::sigmaNN1990[121] =
253 {
254  0.0,
255  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
256  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
257  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
258  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
259  0.000, 0.000, 0.001, 0.001, 0.002, 0.003, 0.005, 0.009,
260  0.017, 0.030, 0.051, 0.076, 0.102, 0.127, 0.149, 0.168,
261  0.185, 0.199, 0.211, 0.221, 0.230, 0.237, 0.243, 0.247,
262  0.251, 0.254, 0.256, 0.258, 0.258, 0.259, 0.258, 0.258,
263  0.257, 0.255, 0.254, 0.252, 0.249, 0.247, 0.245, 0.242,
264  0.239, 0.236, 0.233, 0.230, 0.226, 0.223, 0.220, 0.216,
265  0.215, 0.167, 0.131, 0.104, 0.085, 0.070, 0.059, 0.050,
266  0.043, 0.038, 0.033, 0.029, 0.026, 0.023, 0.021, 0.019,
267  0.017, 0.016, 0.014, 0.013, 0.012, 0.011, 0.010, 0.010,
268  0.009, 0.009, 0.008, 0.007, 0.007, 0.007, 0.006, 0.006,
269  0.006, 0.005, 0.005, 0.005, 0.005, 0.004, 0.004, 0.004
270 };
271 
272 const G4double G4XNNstarTable::sigmaNN2090[121] =
273 {
274  0.0,
275  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
276  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
277  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
278  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
279  0.000, 0.001, 0.001, 0.001, 0.002, 0.002, 0.003, 0.005,
280  0.007, 0.011, 0.016, 0.024, 0.036, 0.053, 0.071, 0.089,
281  0.106, 0.120, 0.133, 0.142, 0.151, 0.158, 0.164, 0.169,
282  0.172, 0.175, 0.178, 0.179, 0.180, 0.180, 0.180, 0.180,
283  0.179, 0.178, 0.176, 0.175, 0.173, 0.171, 0.169, 0.166,
284  0.164, 0.162, 0.159, 0.156, 0.154, 0.151, 0.148, 0.145,
285  0.144, 0.107, 0.081, 0.063, 0.050, 0.041, 0.034, 0.028,
286  0.024, 0.021, 0.018, 0.016, 0.014, 0.013, 0.011, 0.010,
287  0.009, 0.008, 0.008, 0.007, 0.006, 0.006, 0.006, 0.005,
288  0.005, 0.004, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003,
289  0.003, 0.003, 0.003, 0.003, 0.002, 0.002, 0.002, 0.002
290 };
291 
292 const G4double G4XNNstarTable::sigmaNN2190[121] =
293 {
294  0.0,
295  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
296  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
297  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
298  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.001,
299  0.001, 0.002, 0.003, 0.005, 0.009, 0.015, 0.024, 0.036,
300  0.050, 0.064, 0.078, 0.090, 0.100, 0.110, 0.118, 0.125,
301  0.131, 0.136, 0.140, 0.143, 0.146, 0.149, 0.150, 0.152,
302  0.153, 0.153, 0.153, 0.153, 0.153, 0.152, 0.152, 0.151,
303  0.150, 0.148, 0.147, 0.145, 0.144, 0.142, 0.140, 0.139,
304  0.138, 0.110, 0.087, 0.069, 0.056, 0.047, 0.039, 0.033,
305  0.029, 0.025, 0.022, 0.019, 0.017, 0.015, 0.014, 0.013,
306  0.012, 0.010, 0.009, 0.009, 0.008, 0.007, 0.007, 0.006,
307  0.006, 0.006, 0.005, 0.005, 0.005, 0.004, 0.004, 0.004,
308  0.004, 0.003, 0.003, 0.003, 0.003, 0.003, 0.003, 0.003
309 };
310 
311 const G4double G4XNNstarTable::sigmaNN2220[121] =
312 {
313  0.0,
314  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
315  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
316  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
317  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
318  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
319  0.001, 0.001, 0.001, 0.003, 0.005, 0.008, 0.015, 0.026,
320  0.041, 0.058, 0.074, 0.089, 0.101, 0.113, 0.122, 0.131,
321  0.138, 0.144, 0.149, 0.153, 0.156, 0.159, 0.161, 0.163,
322  0.164, 0.165, 0.165, 0.165, 0.165, 0.165, 0.164, 0.163,
323  0.162, 0.161, 0.160, 0.158, 0.157, 0.155, 0.153, 0.151,
324  0.150, 0.121, 0.096, 0.077, 0.062, 0.052, 0.043, 0.037,
325  0.032, 0.028, 0.024, 0.022, 0.019, 0.017, 0.015, 0.014,
326  0.013, 0.011, 0.011, 0.010, 0.009, 0.008, 0.008, 0.007,
327  0.007, 0.006, 0.006, 0.005, 0.005, 0.005, 0.005, 0.004,
328  0.004, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003, 0.003
329 };
330 
331 const G4double G4XNNstarTable::sigmaNN2250[121] =
332 {
333  0.0,
334  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
335  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
336  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
337  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
338  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
339  0.000, 0.001, 0.001, 0.002, 0.003, 0.005, 0.009, 0.016,
340  0.027, 0.043, 0.060, 0.076, 0.090, 0.103, 0.114, 0.123,
341  0.131, 0.138, 0.144, 0.149, 0.153, 0.156, 0.159, 0.161,
342  0.163, 0.164, 0.165, 0.165, 0.165, 0.165, 0.165, 0.164,
343  0.163, 0.162, 0.161, 0.159, 0.158, 0.156, 0.154, 0.152,
344  0.152, 0.122, 0.096, 0.077, 0.062, 0.051, 0.043, 0.037,
345  0.031, 0.027, 0.024, 0.022, 0.019, 0.017, 0.015, 0.014,
346  0.012, 0.011, 0.010, 0.009, 0.009, 0.008, 0.008, 0.007,
347  0.007, 0.006, 0.006, 0.005, 0.005, 0.005, 0.004, 0.004,
348  0.004, 0.004, 0.004, 0.003, 0.003, 0.003, 0.003, 0.003
349 };
350 
351 
353 {
354  xMap["N(1440)0"] = (G4double*) sigmaNN1440;
355  xMap["N(1440)+"] = (G4double*) sigmaNN1440;
356 
357  xMap["N(1520)0"] = (G4double*) sigmaNN1520;
358  xMap["N(1520)+"] = (G4double*) sigmaNN1520;
359 
360  xMap["N(1535)0"] = (G4double*) sigmaNN1535;
361  xMap["N(1535)+"] = (G4double*) sigmaNN1535;
362 
363  xMap["N(1650)0"] = (G4double*) sigmaNN1650;
364  xMap["N(1650)+"] = (G4double*) sigmaNN1650;
365 
366  xMap["N(1675)0"] = (G4double*) sigmaNN1675;
367  xMap["N(1675)+"] = (G4double*) sigmaNN1675;
368 
369  xMap["N(1680)0"] = (G4double*) sigmaNN1680;
370  xMap["N(1680)+"] = (G4double*) sigmaNN1680;
371 
372  xMap["N(1700)0"] = (G4double*) sigmaNN1700;
373  xMap["N(1700)+"] = (G4double*) sigmaNN1700;
374 
375  xMap["N(1710)0"] = (G4double*) sigmaNN1710;
376  xMap["N(1710)+"] = (G4double*) sigmaNN1710;
377 
378  xMap["N(1720)0"] = (G4double*) sigmaNN1720;
379  xMap["N(1720)+"] = (G4double*) sigmaNN1720;
380 
381  xMap["N(1900)0"] = (G4double*) sigmaNN1900;
382  xMap["N(1900)+"] = (G4double*) sigmaNN1900;
383 
384  xMap["N(1990)0"] = (G4double*) sigmaNN1990;
385  xMap["N(1990)+"] = (G4double*) sigmaNN1990;
386 
387  xMap["N(2090)0"] = (G4double*) sigmaNN2090;
388  xMap["N(2090)+"] = (G4double*) sigmaNN2090;
389 
390  xMap["N(2190)0"] = (G4double*) sigmaNN2190;
391  xMap["N(2190)+"] = (G4double*) sigmaNN2190;
392 
393  xMap["N(2220)0"] = (G4double*) sigmaNN2220;
394  xMap["N(2220)+"] = (G4double*) sigmaNN2220;
395 
396  xMap["N(2250)0"] = (G4double*) sigmaNN2250;
397  xMap["N(2250)+"] = (G4double*) sigmaNN2250;
398 }
399 
400 
402 { }
403 
404 
406 {
407  // NOTE: the returned pointer is owned by the client
408 
409  if (xMap.find(particleName) != xMap.end())
410  {
411  // Cross section table for the requested particle available in the Map
412  G4PhysicsFreeVector* sigmaVector = new G4PhysicsFreeVector(sizeNNstar);
413  std::map <G4String, G4double*, std::less<G4String> >::const_iterator iter;
414  G4double* sigmaPointer = 0;
415  for (iter = xMap.begin(); iter != xMap.end(); ++iter)
416  {
417  G4String str = (*iter).first;
418  if (str == particleName)
419  {
420  sigmaPointer = (*iter).second;
421  }
422  }
423  G4int i;
424  for (i=0; i<sizeNNstar; i++)
425  {
426  G4double value = *(sigmaPointer + i) * millibarn;
427  G4double energy = energyTable[i] * GeV;
428  sigmaVector->PutValue(i,energy,value);
429  }
430  return sigmaVector;
431  }
432  else
433  // No cross section table for the requested particle is available in the Map
434  return 0;
435 }
436 
437 
438 
G4int first(char) const
void PutValue(size_t index, G4double energy, G4double dataValue)
int G4int
Definition: G4Types.hh:78
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual const G4PhysicsVector * CrossSectionTable(const G4String &particleName) const
virtual ~G4XNNstarTable()
G4double energy(const ThreeVector &p, const G4double m)
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
static constexpr double millibarn
Definition: G4SIunits.hh:106