source: trunk/source/processes/hadronic/cross_sections/src/G4NucleonNuclearCrossSection.cc @ 1316

Last change on this file since 1316 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

File size: 29.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// author: Vladimir.Grichine@cern.ch
27//
28// Implements data from: Barashenkov V.S., Nucleon-Nucleus Cross Section,
29// Preprint JINR P2-89-770, p. 12, Dubna 1989 (scanned version from KEK)
30// Based on G. Folger version of G4PiNuclearCrossSection class
31//
32// Modified:
33// 05.03.07 V.Ivanchenko - add IfZAApplicable, remove "debug"
34// 11.09.09 V.Ivanchenko - fixed bug in interpolation
35//
36
37#include "G4NucleonNuclearCrossSection.hh"
38#include "G4Neutron.hh"
39#include "G4Proton.hh"
40
41// Group 1: He, Be, C for 44 energies 
42
43const G4double G4NucleonNuclearCrossSection::e1[44] =     
44{
45  0.014, 0.015, 0.017, .02, 0.022, 0.025, 0.027, 0.03, 0.035, .04, 0.045, 0.05, .06, 0.07, 
46  .08, 0.09,  .1, .12, .14, .15, .16, .18, .20, .25, .30, .35, .4 , 0.5, 0.6, 0.7,  0.8, 
47  0.9,   1, 1.5,   2,   3,   5,  7, 10,   
48  20,   50,  100,  500, 1000
49};
50
51const G4double G4NucleonNuclearCrossSection::he_m_t[44] =   
52{ 
53  1090, 1020, 915, 800, 710, 640, 600, 560, 500, 440, 390, 360, 295, 256, 220, 192, 
54  168, 136, 120, 116, 114, 110, 107, 104, 106, 108, 110, 120, 126, 135, 140, 144, 146, 
55  148, 152, 150, 146, 142, 138, 132, 129, 126, 127, 128 
56};
57const G4double G4NucleonNuclearCrossSection::he_m_in[44] = 
58{ 
59  0, 5, 10, 20, 35, 55, 70, 80, 90, 105, 115, 115, 100, 90, 86, 84, 84, 82, 80, 80, 80, 80, 
60  79, 78, 80, 84, 88, 94, 100, 105, 108, 108, 108, 112, 114, 114, 112, 110, 108, 106, 104, 
61  101, 102, 102
62};
63const G4double G4NucleonNuclearCrossSection::he_p_in[44] = 
64{ 
65  0, 2, 3, 13, 30, 50, 65, 77, 90, 105, 115, 115, 100, 90, 86, 84, 84, 82, 80, 80, 80, 80, 
66  79, 78, 80, 84, 88, 94, 100, 105, 108, 108, 108, 112, 114, 114, 112, 110, 108, 106, 104, 
67  101, 102, 102
68};
69
70const G4double G4NucleonNuclearCrossSection::be_m_t[44] = 
71{
72  1490, 1460, 1400, 1350, 1270, 1200, 1160, 1100, 1000, 910, 810, 740, 625, 575, 455, 406, 
73  365, 310, 275, 262, 255, 240, 235, 225, 225, 230, 238, 252, 270, 282, 288, 290, 294, 303, 
74  303, 300, 292, 284, 277, 267, 263, 264, 268, 268 
75};
76const G4double G4NucleonNuclearCrossSection::be_m_in[44] =
77{ 
78  650, 640, 617, 595, 555, 520, 495, 470, 430, 385, 350, 320, 270, 250, 210, 190, 185, 178, 
79  175, 175, 175, 175, 175, 170, 170, 172, 176, 184, 194, 200, 209, 213, 214, 216, 216, 212, 
80  210, 210, 210, 210, 210, 210, 210, 210
81};
82const G4double G4NucleonNuclearCrossSection::be_p_in[44] =
83{ 
84  490, 540, 580, 545, 525, 495, 470, 450, 420, 370, 340, 310, 262, 242, 205, 185, 180, 175, 
85  172, 175, 175, 175, 175, 170, 170, 172, 176, 184, 194, 200, 209, 213, 214, 216, 216, 212, 
86  210, 210, 210, 210, 210, 210, 210, 210
87};
88
89const G4double G4NucleonNuclearCrossSection::c_m_t[44] = 
90{
91  1240, 1370, 1450, 1455, 1445, 1385, 1345, 1290, 1210, 1110, 1020, 940, 800, 700, 604, 530, 
92  475, 396, 350, 336, 320, 303, 294, 280, 280, 286, 296, 314, 330, 344, 356, 360, 364, 384, 
93  388, 384, 364, 352, 344, 330, 324, 324, 332, 332
94};
95const G4double G4NucleonNuclearCrossSection::c_m_in[44] =
96{
97  590, 570, 542, 510, 500, 460, 445, 430, 395, 380, 350, 330, 295, 270, 255, 240, 228, 222, 
98  216, 216, 210, 210, 210, 208, 210, 214, 216, 228, 240, 248, 254, 257, 260, 262, 260, 256, 
99  252, 252, 250, 250, 248, 248, 248, 248
100};
101const G4double G4NucleonNuclearCrossSection::c_p_in[44] =
102{ 
103  310, 330, 400, 440, 450, 435, 430, 420, 385, 370, 340, 320, 288, 263, 249, 234, 222, 216, 
104  210, 211, 205, 208, 210, 208, 210, 214, 216, 228, 240, 248, 254, 257, 260, 262, 260, 256, 
105  252, 252, 250, 250, 248, 248, 248, 248
106};
107
108// Group 2: N, O, Na for 44 energies (e1=e2)
109
110const G4double G4NucleonNuclearCrossSection::e2[44] =   
111{
112 0.014, 0.015, 0.017, .02, 0.022, 0.025, 0.027, 0.03, 0.035, .04, 0.045, 0.05, .06, 0.07, 
113  .08, 0.09,  .1, .12, .14, .15, .16, .18, .20, .25, .30, .35, .4 , 0.5, 0.6, 0.7,  0.8, 
114  0.9,   1, 1.5,   2,   3,   5,  7, 10,   
115  20,   50,  100,  500, 1000 
116};
117
118const G4double G4NucleonNuclearCrossSection::n_m_t[44] = 
119{
120  1420,1480, 1537, 1550, 1525, 1500, 1480, 1425, 1340, 1260, 1175, 1090, 930, 805, 690, 612, 
121  552, 462, 402, 384, 372, 350, 345, 326, 324, 328, 336, 356, 372, 388, 400, 408, 415, 430, 
122  435, 432, 415, 402, 390, 375, 367, 370, 382, 385
123};
124const G4double G4NucleonNuclearCrossSection::n_m_in[44] =
125{
126  680, 665, 625, 580, 562, 525, 510, 485, 450, 435, 410, 387, 340, 310, 290, 280, 276, 274, 
127  260, 258, 254, 247, 245, 240, 240, 244, 250, 260, 268, 275, 280, 285, 290, 295, 300, 294, 
128  292, 290, 285, 285, 282, 282, 282, 282
129};
130const G4double G4NucleonNuclearCrossSection::n_p_in[44] =
131{ 
132  420, 440, 470, 490, 497, 500, 480, 462, 440, 425, 400, 377, 333, 303, 284, 274, 270, 268, 
133  254, 252, 247, 245, 245, 240, 240, 244, 250, 260, 268, 275, 280, 285, 290, 295, 300, 294, 
134  292, 290, 285, 285, 282, 282, 282, 282
135};
136
137const G4double G4NucleonNuclearCrossSection::o_m_t[44] = 
138{
139  1520, 1570, 1630, 1660, 1647, 1623, 1595, 1555, 1475, 1395, 1290, 1207, 1035, 925, 816, 
140  720, 645, 540, 462, 438, 415, 392, 378, 362, 361, 381, 390, 403, 417, 440, 460, 470, 
141  479, 498, 504, 498, 477, 457, 443, 427, 420, 425, 429, 430
142};
143const G4double G4NucleonNuclearCrossSection::o_m_in[44] = 
144{
145  750, 740, 700, 650, 620, 575, 555, 530, 505, 462, 435, 420, 375, 345, 320, 310, 300, 293, 
146  288, 282, 282, 280, 276, 270, 271, 275, 280, 290, 295, 304, 310, 315, 318, 332, 335, 330, 
147  323, 320, 317, 315, 315, 315, 315, 315
148};
149const G4double G4NucleonNuclearCrossSection::o_p_in[44] = 
150{
151  460, 485, 510, 535, 537, 532, 520, 500, 460, 432, 405, 390, 350, 320, 310, 304, 293, 287, 
152  283, 279, 279, 278, 276, 270, 271, 275, 280, 290, 295, 304, 310, 315, 318, 332, 335, 330, 
153  323, 320, 317, 315, 315, 315, 315, 315
154};
155
156const G4double G4NucleonNuclearCrossSection::na_m_t[44] = 
157{
158  1570, 1620, 1695, 1730, 1750, 1760, 1755, 1740, 1710, 1643, 1560, 1480, 1343, 1220, 1073, 
159  953, 860, 720, 618, 582, 546, 522, 504, 484, 492, 500, 512, 538, 560, 586, 608, 622, 632, 
160  660, 668, 664, 640, 616, 596, 568, 568, 568, 568, 568 
161};
162const G4double G4NucleonNuclearCrossSection::na_m_in[44] =
163{
164  960, 930, 890, 822, 790, 750, 725, 686, 620, 600, 575, 540, 497, 450, 414, 390, 380, 372, 
165  354, 360, 355, 354, 350, 350, 350, 356, 364, 384, 392, 400, 408, 410, 420, 408, 412, 420, 
166  411, 409, 407, 403, 400, 400, 400, 400
167};
168const G4double G4NucleonNuclearCrossSection::na_p_in[44] =
169{
170  600, 617, 660, 675, 680, 680, 670, 650, 575, 550, 525, 490, 450, 420, 385, 367, 360, 350, 
171  350, 350, 345, 347, 350, 350, 350, 356, 364, 384, 392, 400, 408, 410, 420, 408, 412, 420, 
172  411, 409, 407, 403, 400, 400, 400, 400
173};
174
175// Al, Si, Ca for 45 energies
176
177const G4double G4NucleonNuclearCrossSection::e3[45] =     
178{
179  0.014, 0.015, 0.016, 0.017, .02, 0.022, 0.025, 0.027, 0.03, 0.035, .04, 0.045, 0.05, 
180  .06, 0.07, .08, 0.09,  .1, .12, .14, .15, .16, .18, .20, .25, .30, .35, .4 , 0.5, 0.6, 
181  0.7,  0.8,  0.9,   1, 1.5,   2,   3,   5,  7, 10,  20,   50,  100,  500, 1000 
182};
183
184const G4double G4NucleonNuclearCrossSection::al_m_t[45] = 
185{ 
186  1735, 1750, 1760, 1795, 1830, 1855, 1885, 1895, 1900, 1870, 1835, 1785, 1710, 1522, 
187  1350, 1212, 1080, 972, 816, 720, 678, 642, 600, 567, 558, 560, 578, 592, 616, 644, 
188  672, 688, 708, 720, 736, 754, 736, 706, 680, 672, 646, 632, 632, 632, 632
189};
190const G4double G4NucleonNuclearCrossSection::al_m_in[45] = 
191{
192  1000, 990, 975, 950, 905, 875, 825, 800, 762, 690, 652, 610, 570, 495, 480, 456, 444, 
193  432, 420, 420, 420, 420, 410, 410, 400, 402, 404, 408, 424, 438, 448, 450, 454, 456, 
194  472, 480, 466, 456, 452, 448, 444, 440, 440, 440, 440
195};
196const G4double G4NucleonNuclearCrossSection::al_p_in[45] = 
197{
198  650, 682, 690, 715, 750, 762, 750, 740, 720, 655, 617, 575, 540, 470, 455, 532, 420, 
199  408, 400, 403, 403, 408, 406, 404, 400, 402, 404, 408, 424, 438, 448, 450, 454, 456, 
200  472, 480, 466, 456, 452, 448, 444, 440, 440, 440, 440 
201};
202
203const G4double G4NucleonNuclearCrossSection::si_m_t[45] = 
204{ 
205  1810, 1833, 1850, 1872, 1920, 1950, 1995, 2020, 2035, 2000, 1930, 1850, 1760, 1570, 
206  1400, 1255, 1110, 1008, 846, 742, 696, 671, 623, 588, 584, 584, 602, 618, 645, 679, 
207  708, 727, 746, 757, 769, 782, 771, 734, 710, 698, 672, 654, 650, 650, 650 
208};
209const G4double G4NucleonNuclearCrossSection::si_m_in[45] = 
210{
211  1060, 1035, 1015, 990, 935, 900, 860, 830, 790, 725, 665, 630, 600, 520, 504, 486, 
212  470, 456, 444, 432, 432, 432, 418, 418, 415, 412, 416, 422, 440, 460, 472, 476, 479, 
213  480, 492, 496, 488, 472, 472, 464, 460, 452, 448, 448, 448 
214};
215const G4double G4NucleonNuclearCrossSection::si_p_in[45] = 
216{
217  670, 700, 725, 750, 780, 780, 770, 757, 735, 690, 635, 585, 570, 490, 475, 460, 446, 
218  431, 423, 425, 425, 425, 425, 422, 422, 412, 416, 422, 440, 460, 472, 476, 479, 
219  480, 492, 496, 488, 472, 472, 464, 460, 452, 448, 448, 448
220};
221
222const G4double G4NucleonNuclearCrossSection::ca_m_t[45] = 
223{ 
224  2180, 2130, 2095, 2075, 2115, 2150, 2220, 2250, 2300, 2365, 2360, 2280, 2180, 2000, 
225  1805, 1650, 1500, 1340, 1140, 990, 940, 890, 825, 790, 770, 773, 787, 800, 830, 870, 
226  905, 930, 950, 965, 990, 1002, 990, 965, 945, 925, 892, 860, 860, 860, 860 
227};
228const G4double G4NucleonNuclearCrossSection::ca_m_in[45] = 
229{
230  1240, 1225, 1200, 1180, 1125, 1090, 1045, 1020, 980, 925, 880, 825, 770, 680, 640, 
231  620, 615, 600, 580, 565, 560, 560, 560, 550, 535, 530, 540, 550, 570, 595, 610, 615, 
232  620, 622, 629, 630, 620, 612, 607, 592, 587, 580, 580, 580, 580 
233};
234const G4double G4NucleonNuclearCrossSection::ca_p_in[45] = 
235{
236  770, 800, 823, 850, 900, 925, 935, 920, 895, 835, 800, 750, 715, 640, 605, 590, 588, 
237  573, 555, 543, 540, 540, 540, 535, 530, 530, 540, 550, 570, 595, 610, 615, 
238  620, 622, 629, 630, 620, 612, 607, 592, 587, 580, 580, 580, 580 
239};
240
241// Fe, Cu, Mo for 47 energies
242
243const G4double G4NucleonNuclearCrossSection::e4[47] =     
244{
245  0.014, 0.015, 0.017, .02, 0.022, 0.025, 0.027, 0.03, 0.033, 0.035, 0.037, .04, 0.045, 
246  0.05, 0.055, .06, 0.07, .08, 0.09,  .1, .12, .14, .15, .16, .18, .20, .25, .30, .35, 
247  .4 , 0.5, 0.6, 0.7,  0.8,  0.9,   1, 1.5,   2,   3,   5,  7, 10,   
248  20,   50,  100,  500, 1000
249};
250
251const G4double G4NucleonNuclearCrossSection::fe_m_t[47] = 
252{
253  2580, 2490, 2370, 2282, 2275, 2285, 2320, 2370, 2432, 2445, 2460, 2485, 2530, 2540, 
254  2517, 2480, 2290, 2110, 1940, 1790, 1510, 1290, 1220, 1150, 1070, 1030, 1013, 1020, 
255  1030, 1043, 1075, 1110, 1133, 1163, 1185, 1225, 1252, 1260, 1260, 1233, 1207, 1185, 
256  1140, 1110, 1110, 1110, 1110
257};
258const G4double G4NucleonNuclearCrossSection::fe_m_in[47] = 
259{
260  1440, 1433, 1390, 1325, 1280, 1260, 1215, 1180, 1140, 1110, 1080, 1040, 990, 955, 920, 
261  885, 835, 800, 780, 765, 750, 725, 720, 720, 710, 700, 700, 700, 712, 705, 735, 750, 
262  765, 775, 780, 795, 810, 813, 810, 784, 757, 743, 735, 720, 720, 720, 720 
263};
264const G4double G4NucleonNuclearCrossSection::fe_p_in[47] = 
265{
266  900, 960, 1070, 1090, 1115, 1120, 1115, 1080, 1045, 1025, 1000, 960, 900, 885, 865, 790, 
267  765, 740, 720, 700, 697, 697, 697, 697, 695, 690, 688, 690, 712, 705, 735, 750, 
268  765, 775, 780, 795, 810, 813, 810, 784, 757, 743, 735, 720, 720, 720, 720 
269};
270
271const G4double G4NucleonNuclearCrossSection::cu_m_t[47] = 
272{
273  2920, 2800, 2615, 2480, 2455, 2430, 2440, 2460, 2500, 2530, 2560, 2615, 2690, 2720, 
274  2700, 2645, 2500, 2320, 2140, 1970, 1670, 1460, 1380, 1285, 1200, 1160, 1140, 1147, 
275  1163, 1170, 1200, 1237, 1265, 1285, 1305, 1328, 1375, 1390, 1395, 1370, 1335, 1315, 
276  1270, 1230, 1230, 1230, 1230 
277};
278const G4double G4NucleonNuclearCrossSection::cu_m_in[47] = 
279{
280  1540, 1535, 1500, 1445, 1407, 1380, 1330, 1300, 1285, 1270, 1240, 1190, 1090, 1010, 
281  940, 920, 860, 835, 820, 810, 800, 780, 775, 770, 760, 760, 758, 765, 765, 770, 795, 
282  810, 825, 830, 840, 848, 870, 870, 868, 840, 825, 810, 803, 795, 795, 795, 795 
283};
284const G4double G4NucleonNuclearCrossSection::cu_p_in[47] = 
285{
286  935, 1000, 1060, 1190, 1220, 1250, 1240, 1210, 1150, 1130, 1115, 1050, 985, 950, 890, 
287  870, 820, 800, 785, 780, 770, 750, 745, 740, 735, 735, 745, 760, 762, 770, 795, 
288  810, 825, 830, 840, 848, 870, 870, 868, 840, 825, 810, 803, 795, 795, 795, 795 
289};
290
291const G4double G4NucleonNuclearCrossSection::mo_m_t[47] = 
292{
293  4150, 4040, 3800, 3490, 3300, 3060, 2960, 2845, 2785, 2820, 2850, 2980, 3170, 3230, 
294  3270, 3280, 3225, 3075, 2895, 2710, 2355, 2060, 1925, 1800, 1630, 1560, 1540, 1550, 
295  1570, 1590, 1650, 1685, 1715, 1740, 1760, 1780, 1850, 1880, 1858, 1815, 1790, 1782, 
296  1720, 1690, 1690, 1690, 1690 
297};
298const G4double G4NucleonNuclearCrossSection::mo_m_in[47] = 
299{
300  1790, 1775, 1740, 1680, 1640, 1580, 1550, 1510, 1460, 1440, 1418, 1380, 1330, 1280, 
301  1240, 1200, 1155, 1140, 1110, 1110, 1080, 1065, 1050, 1050, 1025, 1020, 1015, 1020, 
302  1022, 1026, 1060, 1085, 1100, 1110, 1120, 1127, 1150, 1160, 1140, 1100, 1085, 1080, 
303  1070, 1070, 1070, 1070, 1070 
304};
305const G4double G4NucleonNuclearCrossSection::mo_p_in[47] = 
306{
307  1025, 1080, 1190, 1380, 1440, 1495, 1475, 1420, 1350, 1310, 1300, 1290, 1250, 1200, 
308  1170, 1130, 1095, 1060, 1040, 1022, 1020, 1016, 1016, 1016, 1016, 1012, 1005, 1005, 
309  1005, 1010, 1060, 1085, 1100, 1110, 1120, 1127, 1150, 1160, 1140, 1100, 1085, 1080, 
310  1070, 1070, 1070, 1070, 1070
311};
312
313// Cd, Sn, W for 48 energies
314
315const G4double G4NucleonNuclearCrossSection::e5[48] =     
316{
317   0.014, 0.015, 0.017, 0.018, .02, 0.022, 0.025, 0.027, 0.03, 0.033, 0.035, .04, 
318   0.045, 0.05, 0.055, .06, .065, 0.07, .08, 0.09,  .1, .12, .14, .15, .16, .18, 
319   .20, .25, .30, .35, .4 , 0.5, 0.6, 0.7,  0.8, 0.9, 1, 1.5, 2, 3, 5, 7, 10,  20, 
320   50, 100, 500, 1000
321};
322
323const G4double G4NucleonNuclearCrossSection::cd_m_t[48] = 
324{
325  4420, 4280, 4170, 4070, 3860, 3680, 3420, 3280, 3125, 3060, 3080, 3190, 3350, 3445, 
326  3510, 3540, 3560, 3550, 3460, 3300, 3030, 2640, 2340, 2190, 2070, 1950, 1770, 1732, 
327  1740, 1760, 1780, 1832, 1885, 1925, 1945, 1960, 1980, 2070, 2080, 2065, 2040, 2022, 
328  1980, 1940, 1870, 1870, 1870, 1870 
329};
330const G4double G4NucleonNuclearCrossSection::cd_m_in[48]= 
331{
332  1920, 1910, 1880, 1860, 1840, 1800, 1760, 1720, 1675, 1630, 1600, 1520, 1465, 1420, 
333  1390, 1340, 1310, 1280, 1275, 1235, 1225, 1200, 1170, 1170, 1170, 1165, 1145, 1140, 
334  1140, 1135, 1160, 1180, 1220, 1240, 1250, 1260, 1265, 1270, 1275, 1250, 1222, 1222, 
335  1220, 1215, 1190, 1190, 1190, 1190 
336};
337const G4double G4NucleonNuclearCrossSection::cd_p_in[48] = 
338{
339  1020, 1100, 1225, 1290, 1440, 1520, 1575, 1560, 1518, 1460, 1420, 1400, 1365, 1340, 
340  1300, 1280, 1260, 1200, 1190, 1160, 1125, 1125, 1125, 1125, 1125, 1125, 1120, 1120, 
341  1120, 1118, 1146, 1180, 1220, 1240, 1250, 1260, 1265, 1270, 1275, 1250, 1222, 1222, 
342  1220, 1215, 1190, 1190, 1190, 1190 
343};
344
345const G4double G4NucleonNuclearCrossSection::sn_m_t[48] = 
346{
347  4420, 4400, 4260, 4150, 3980, 3770, 3530, 3370, 3245, 3180, 3170, 3260, 3400, 3500, 
348  3560, 3610, 3650, 3680, 3580, 3390, 3190, 2760, 2430, 2295, 2175, 1990, 1880, 1810, 
349  1820, 1840, 1865, 1940, 1985, 2020, 2040, 2060, 2080, 2160, 2185, 2180, 2110, 2105, 
350  2080, 2050, 1980, 1980, 1980, 1980 
351};
352const G4double G4NucleonNuclearCrossSection::sn_m_in[48] = 
353{
354  1945, 1940, 1905, 1890, 1860, 1830, 1780, 1755, 1717, 1680, 1645, 1570, 1500, 1455, 
355  1410, 1370, 1340, 1320, 1290, 1285, 1260, 1240, 1235, 1212, 1200, 1200, 1200, 1190, 
356  1190, 1200, 1210, 1240, 1270, 1285, 1300, 1300, 1310, 1320, 1320, 1290, 1240, 1240, 
357  1240, 1240, 1240, 1240, 1240, 1240 
358};
359const G4double G4NucleonNuclearCrossSection::sn_p_in[48] = 
360{ 
361  1020, 1080, 1270, 1335, 1465, 1505, 1610, 1610, 1550, 1535, 1500, 1440, 1407, 1370, 
362  1340, 1300, 1285, 1260, 1230, 1215, 1200, 1180, 1170, 1170, 1165, 1165, 1170, 1165, 
363  1165, 1183, 1195, 1240, 1270, 1285, 1300, 1300, 1310, 1320, 1320, 1290, 1240, 1240, 
364  1240, 1240, 1240, 1240, 1240, 1240
365};
366
367const G4double G4NucleonNuclearCrossSection::w_m_t[48] =   
368{
369  5320, 5430, 5480, 5450, 5330, 5190, 4960, 4790, 4550, 4340, 4200, 4070, 4000, 4030, 
370  4125, 4220, 4270, 4390, 4440, 4360, 4200, 3800, 3380, 3200, 3040, 2790, 2660, 2575, 
371  2575, 2600, 2640, 2690, 2755, 2790, 2812, 2837, 2850, 2950, 3000, 2970, 2940, 2910, 
372  2880, 2820, 2730, 2730, 2730, 2730 
373};
374const G4double G4NucleonNuclearCrossSection::w_m_in[48] = 
375{
376  2440, 2400, 2370, 2350, 2310, 2270, 2220, 2195, 2150, 2100, 2070, 2010, 1945, 1900, 
377  1850, 1820, 1780, 1760, 1730, 1720, 1680, 1680, 1660, 1660, 1650, 1650, 1640, 1640, 
378  1612, 1615, 1625, 1640, 1700, 1720, 1730, 1740, 1750, 1780, 1780, 1750, 1740, 1735, 
379  1710, 1695, 1680, 1680, 1680, 1680 
380};
381const G4double G4NucleonNuclearCrossSection::w_p_in[48] = 
382{ 
383  950,  1020, 1240, 1400, 1560, 1670, 1760, 1830, 1850, 1855, 1870, 1840, 1800, 1770, 
384  1740, 1715, 1680, 1670, 1650, 1620, 1610, 1600, 1600, 1600, 1600, 1600, 1600, 1595, 
385  1585, 1595, 1615, 1640, 1700, 1720, 1730, 1740, 1750, 1780, 1780, 1750, 1740, 1735, 
386  1710, 1695, 1680, 1680, 1680, 1680
387};
388
389// Pb, U for 46 energies
390
391const G4double G4NucleonNuclearCrossSection::e6[46] =     
392{
393  0.014, 0.015, 0.017, 0.019, .02, 0.022, 0.025, 0.027, 0.03, 0.035, .04, 0.045, 0.05, 
394  0.055, .06, 0.07, .08, 0.09, .1, .12, .14, .15, .16, .18, .20, .25, .30, .35, .4 , 
395  0.5, 0.6, 0.7,  0.8, 0.9, 1, 1.5, 2, 3, 5, 7, 10, 20, 50, 100, 500, 1000
396};
397
398const G4double G4NucleonNuclearCrossSection::pb_m_t[46] = 
399{
400  5300, 5440, 5720, 5880, 5765, 5745, 5480, 5280, 4970, 4550, 4390, 4300, 4265, 4325, 
401  4450, 4540, 4740, 4710, 4600, 4100, 3660, 3480, 3300, 3000, 2890, 2865, 2855, 2850, 
402  2865, 2920, 2955, 3000, 3030, 3060, 3105, 3240, 3290, 3270, 3240, 3180, 3090, 3060, 
403  2970, 2970, 2970, 2970 
404
405};
406const G4double G4NucleonNuclearCrossSection::pb_m_in[46] = 
407{
408  2580, 2550, 2505, 2462, 2460, 2435, 2380, 2355, 2280, 2180, 2170, 2130, 2080, 2035, 
409  1980, 1940, 1900, 1870, 1840, 1800, 1800, 1800, 1780, 1760, 1760, 1740, 1730, 1725, 
410  1740, 1785, 1815, 1835, 1860, 1890, 1895, 1920, 1920, 1890, 1850, 1835, 1830, 1830, 
411  1830, 1830, 1830, 1830 
412};
413const G4double G4NucleonNuclearCrossSection::pb_p_in[46] = 
414{ 
415  900,  1060, 1200, 1420, 1515, 1620, 1750, 1800, 1915, 2030, 1960, 1940, 1910, 1860, 
416  1840, 1780, 1770, 1760, 1740, 1720, 1725, 1740, 1740, 1730, 1720, 1700, 1710, 1720, 
417  1730, 1740, 1815, 1835, 1860, 1890, 1895, 1920, 1920, 1890, 1850, 1835, 1830, 1830, 
418  1830, 1830, 1830, 1830
419};
420
421const G4double G4NucleonNuclearCrossSection::u_m_t[46] =   
422{
423  5800, 5940, 6160, 6345, 6360, 6350, 6170, 6020, 5760, 5350, 4990, 4800, 4710, 4690, 
424  4760, 5040, 5190, 5200, 5080, 4600, 4120, 3920, 3720, 3420, 3240, 3150, 3160, 3180, 
425  3210, 3240, 3280, 3350, 3390, 3435, 3480, 3560, 3585, 3580, 3540, 3500, 3470, 3410, 
426  3335, 3335, 3335, 3335   
427};
428const G4double G4NucleonNuclearCrossSection::u_m_in[46] = 
429{
430  2820, 2770, 2700, 2660, 2645, 2620, 2580, 2550, 2515, 2450, 2390, 2320, 2260, 2225, 
431  2200, 2140, 2080, 2060, 2040, 2000, 1980, 1965, 1960, 1930, 1920, 1890, 1905, 1920, 
432  1945, 1970, 1985, 2010, 2040, 2070, 2080, 2090, 2095, 2080, 2063, 2060, 2050, 2040, 
433  2005, 2005, 2005, 2005 
434};
435const G4double G4NucleonNuclearCrossSection::u_p_in[46] = 
436{ 
437  800,  900,  1100, 1300, 1410, 1510, 1680, 1800, 2000, 2200, 2080, 2060, 2035, 2100, 
438  2030, 2030, 2000, 1960, 1960, 1960, 1940, 1925, 1920, 1905, 1890, 1860, 1880, 1910, 
439  1930, 1945, 1985, 2010, 2040, 2070, 2080, 2090, 2095, 2080, 2063, 2060, 2050, 2040, 
440  2005, 2005, 2005, 2005
441};
442
443using namespace std;
444
445//////////////////////////////////////////////////////////////////////////////////
446//
447//
448
449G4NucleonNuclearCrossSection::G4NucleonNuclearCrossSection()
450{
451  theNeutron = G4Neutron::Neutron();
452  theProton  = G4Proton::Proton();
453 
454  // He, Be, C
455   
456   thePimData.push_back(new G4PiData(he_m_t, he_m_in, e1, 44));
457   thePipData.push_back(new G4PiData(he_m_t, he_p_in, e1, 44));
458
459   thePimData.push_back(new G4PiData(be_m_t, be_m_in, e1, 44));
460   thePipData.push_back(new G4PiData(be_m_t, be_p_in, e1, 44));
461
462   thePimData.push_back(new G4PiData(c_m_t,  c_m_in,  e1, 44));
463   thePipData.push_back(new G4PiData(c_m_t,  c_p_in,  e1, 44));
464
465   // N, O, Na
466
467   thePimData.push_back(new G4PiData(n_m_t,  n_m_in,  e2, 44));
468   thePipData.push_back(new G4PiData(n_m_t,  n_p_in,  e2, 44));
469
470   thePimData.push_back(new G4PiData(o_m_t,  o_m_in,  e2, 44));
471   thePipData.push_back(new G4PiData(o_m_t,  o_p_in,  e2, 44));
472
473   thePimData.push_back(new G4PiData(na_m_t, na_m_in, e2, 44));
474   thePipData.push_back(new G4PiData(na_m_t, na_p_in, e2, 44));
475
476   // Al, Si, Ca
477
478   thePimData.push_back(new G4PiData(al_m_t, al_m_in, e3, 45));
479   thePipData.push_back(new G4PiData(al_m_t, al_p_in, e3, 45));
480
481   thePimData.push_back(new G4PiData(si_m_t, si_m_in, e3, 45));
482   thePipData.push_back(new G4PiData(si_m_t, si_p_in, e3, 45));
483
484   thePimData.push_back(new G4PiData(ca_m_t, ca_m_in, e3, 45));
485   thePipData.push_back(new G4PiData(ca_m_t, ca_p_in, e3, 45));
486
487   // Fe, Cu, Mo
488
489   thePimData.push_back(new G4PiData(fe_m_t, fe_m_in, e4, 47));
490   thePipData.push_back(new G4PiData(fe_m_t, fe_p_in, e4, 47));
491
492   thePimData.push_back(new G4PiData(cu_m_t, cu_m_in, e4, 47));
493   thePipData.push_back(new G4PiData(cu_m_t, cu_p_in, e4, 47));
494
495   thePimData.push_back(new G4PiData(mo_m_t, mo_m_in, e4, 47));
496   thePipData.push_back(new G4PiData(mo_m_t, mo_p_in, e4, 47));
497
498   // Cd, Sn, W
499
500   thePimData.push_back(new G4PiData(cd_m_t, cd_m_in, e5, 48));
501   thePipData.push_back(new G4PiData(cd_m_t, cd_p_in, e5, 48));
502
503   thePimData.push_back(new G4PiData(sn_m_t, sn_m_in, e5, 48));
504   thePipData.push_back(new G4PiData(sn_m_t, sn_p_in, e5, 48));
505
506   thePimData.push_back(new G4PiData(w_m_t,  w_m_in,  e5, 48));
507   thePipData.push_back(new G4PiData(w_m_t,  w_p_in,  e5, 48));
508
509   // Pb, U
510
511   thePimData.push_back(new G4PiData(pb_m_t, pb_m_in, e6, 46));
512   thePipData.push_back(new G4PiData(pb_m_t, pb_p_in, e6, 46));
513
514   thePimData.push_back(new G4PiData(u_m_t,  u_m_in,  e6, 46));
515   thePipData.push_back(new G4PiData(u_m_t,  u_p_in,  e6, 46));
516
517   theZ.push_back(2); // He
518   theZ.push_back(4); // Be
519   theZ.push_back(6); // C
520   theZ.push_back(7); // N
521   theZ.push_back(8); // O
522   theZ.push_back(11); // Na
523   theZ.push_back(13); // Al
524   theZ.push_back(14); // Si
525   theZ.push_back(20); // Ca
526   theZ.push_back(26); // Fe
527   theZ.push_back(29); // Cu
528   theZ.push_back(42); // Mo
529   theZ.push_back(48); // Cd
530   theZ.push_back(50); // Sn
531   theZ.push_back(74); // W
532   theZ.push_back(82); // Pb
533   theZ.push_back(92); // U
534
535}
536
537///////////////////////////////////////////////////////////////////////////////////
538//
539//
540
541G4NucleonNuclearCrossSection::~G4NucleonNuclearCrossSection()
542{
543   std::for_each(thePimData.begin(), thePimData.end(), G4PiData::Delete());
544   std::for_each(thePipData.begin(), thePipData.end(), G4PiData::Delete());
545}
546
547////////////////////////////////////////////////////////////////////////////
548//
549//
550
551G4double G4NucleonNuclearCrossSection::
552GetCrossSection( const G4DynamicParticle* aParticle, 
553                 const G4Element* anElement,
554                       G4double                )
555
556{
557  return GetIsoZACrossSection(aParticle, anElement->GetZ(), anElement->GetN(), 0.);
558}
559
560////////////////////////////////////////////////////////////////////////////
561//
562//
563
564G4bool G4NucleonNuclearCrossSection::IsApplicable(const G4DynamicParticle* aParticle, 
565                                                  const G4Element* anElement)
566{
567  return IsZAApplicable(aParticle, anElement->GetZ(), anElement->GetN()); 
568}
569
570////////////////////////////////////////////////////////////////////////////
571//
572//
573
574G4bool G4NucleonNuclearCrossSection::IsZAApplicable(const G4DynamicParticle* aParticle, 
575                                                    G4double Z, G4double)
576{
577  G4bool result = false;
578  if(aParticle->GetDefinition() == theNeutron ) result = true;
579  if(aParticle->GetDefinition() == theProton)   result = true;
580  if(Z < 1.5)                                   result = false;
581  if(aParticle->GetKineticEnergy() > 999.9*GeV) result = false;
582  return result;
583}
584
585////////////////////////////////////////////////////////////////////////////
586//
587//
588
589G4double G4NucleonNuclearCrossSection::
590GetIsoZACrossSection( const G4DynamicParticle* aParticle, 
591                      G4double  zElement, G4double, G4double  )
592{
593   G4double kineticEnergy = aParticle->GetKineticEnergy();
594 
595   G4double result = 0;
596   G4int Z = G4int(zElement + 0.5);
597
598   // G4cout<<"Z = "<<Z<<G4endl;
599
600   size_t it = 0;
601
602   while( it < theZ.size() && Z > theZ[it] ) {++it;}
603   if(it >= theZ.size()) it = theZ.size() - 1; 
604
605   if( Z > theZ[it] ) Z = theZ[it]; 
606   G4int Z1, Z2;
607   G4double x1, x2, xt1, xt2;
608
609   if(aParticle->GetDefinition() == theNeutron )
610   {
611     if( theZ[it] == Z )
612     {
613       result = thePimData[it]->ReactionXSection(kineticEnergy);
614       fTotalXsc = thePimData[it]->TotalXSection(kineticEnergy);
615     }
616     else
617     {
618       x1 = thePimData[it-1]->ReactionXSection(kineticEnergy);
619       xt1 = thePimData[it-1]->TotalXSection(kineticEnergy);
620       Z1 = theZ[it-1];
621       x2 = thePimData[it]->ReactionXSection(kineticEnergy);
622       xt2 = thePimData[it]->TotalXSection(kineticEnergy);
623       Z2 = theZ[it];
624
625       result = Interpolate(Z1, Z2, Z, x1, x2);
626       fTotalXsc = Interpolate(Z1, Z2, Z, xt1, xt2);
627     }
628   }
629   else if(aParticle->GetDefinition() == theProton)
630   {
631     if( theZ[it] == Z )
632     {
633       // at high energies, when no data for proton, use neutron
634
635       std::vector<G4PiData *> * theData = &thePimData;
636
637       if( thePipData[it]->AppliesTo(kineticEnergy) )
638       {
639         theData = &thePipData;
640       }
641       result = theData->operator[](it)->ReactionXSection(kineticEnergy);
642       fTotalXsc = theData->operator[](it)->TotalXSection(kineticEnergy);
643
644     }
645     else
646     {
647       std::vector<G4PiData*>* theLData = &thePimData;
648
649       if(thePipData[it-1]->AppliesTo(kineticEnergy))
650       {
651         theLData = &thePipData;
652       }
653       std::vector<G4PiData *> * theHData = &thePimData;
654
655       if( thePipData[it]->AppliesTo(kineticEnergy) )
656       {
657         theHData = &thePipData;
658       }
659       x1 = theLData->operator[](it-1)->ReactionXSection(kineticEnergy);
660       xt1 = theLData->operator[](it-1)->TotalXSection(kineticEnergy);
661       Z1 = theZ[it-1];
662       x2 = theHData->operator[](it)->ReactionXSection(kineticEnergy);
663       xt2 = theHData->operator[](it)->TotalXSection(kineticEnergy);
664       Z2 = theZ[it];
665
666       result = Interpolate(Z1, Z2, Z, x1, x2);
667       fTotalXsc = Interpolate(Z1, Z2, Z, xt1, xt2);
668     }
669   }
670   fElasticXsc = fTotalXsc - result;
671   if( fElasticXsc < 0.) fElasticXsc = 0.;
672
673   return result;
674}
675
676/////////////////////////////////////////////////////////////////////////////
677//
678//
679
680G4double G4NucleonNuclearCrossSection::
681Interpolate(G4int Z1, G4int Z2, G4int Z, G4double x1, G4double x2)
682{ 
683//   Nucleon numbers obtained from G4NistManager G4 8.0
684
685  static G4double alpha = 2./3.;
686
687  static const G4double A[92] = 
688  {
689    1.0001, 4.0000, 6.9241, 9.0000, 10.801, 12.011, 14.004, 16.004, 19.000, 20.188,
690    23.000, 24.320, 27.000, 28.109, 31.000, 32.094, 35.484, 39.985, 39.135, 40.116,
691    45.000, 47.918, 50.998, 52.055, 55.000, 55.910, 59.000, 58.760, 63.617, 65.468,
692    69.798, 72.691, 75.000, 79.042, 79.986, 83.887, 85.557, 87.710, 89.000, 91.318,
693    93.000, 96.025, 98.000, 101.16, 103.00, 106.51, 107.96, 112.51, 114.91, 118.81,
694    121.86, 127.70, 127.00, 131.39, 133.00, 137.42, 139.00, 140.21, 141.00, 144.32,
695    145.00, 150.45, 152.04, 157.33, 159.00, 162.57, 165.00, 167.32, 169.00, 173.10,
696    175.03, 178.54, 181.00, 183.89, 186.25, 190.27, 192.25, 195.11, 197.00, 200.63,
697    204.41, 207.24, 209.00, 209.00, 210.00, 222.00, 223.00, 226.00, 227.00, 232.00,
698    231.00, 237.98
699  };                     
700  static G4bool NeedInit = true;
701                     
702  static G4double A75[92];
703
704  if ( NeedInit )
705  {
706    for (G4int i=0; i<92; ++i)
707    {
708      A75[i] = std::pow(A[i], alpha); // interpolate by square ~ A^(2/3)
709    }
710    NeedInit=false;
711  }
712
713  // for tabulated data, cross section scales with A^(2/3)
714  G4double r1 = x1 / A75[Z1-1] * A75[Z-1];
715  G4double r2 = x2 / A75[Z2-1] * A75[Z-1];
716  G4double result = 0.5*(r1+r2);
717
718  // More precise average
719  if(Z1 != Z2) {
720    G4double alp1 = (A[Z-1] - A[Z1-1]);
721    G4double alp2 = (A[Z2-1] - A[Z-1]);
722    result = (r1*alp2 + r2*alp1)/(alp1 + alp2);
723  }
724  //       G4cout << "x1/2, z1/2 z" <<x1<<" "<<x2<<" "<<Z1<<" "<<Z2<<" "<<Z<<G4endl;
725  //       G4cout << "res1/2 " << r1 <<" " << r2 <<" " << result<< G4endl;
726  return result;
727}
Note: See TracBrowser for help on using the repository browser.