COSMOS v7.655  COSMOSv7655
(AirShowerMC)
cgpxs1.f
Go to the documentation of this file.
1 ! gp xsection gg is gp/330
2  subroutine cgpxs1(Eg, xs)
3  use modpdgxs, only : csigmat
4  implicit none
5 #include "Zmass.h"
6  real*8 Eg ! input. momentum /Energy of gamma in GeV
7  real*8 xs ! output. total np cross section in mb
8  integer np, i, m
9  real*8 error, s, rts
10 
11  parameter(np=66, m=5)
12  real*8 roots(np), mb(np)
13  data ( roots(i), i= 1 , np )/
14  1 1.0705 , 1.1016 , 1.1256 ,
15  2 1.1501 , 1.1625 , 1.1878 , 1.2049 ,
16  3 1.2179 , 1.2355 , 1.2534 , 1.2715 ,
17  4 1.2898 , 1.3179 , 1.3274 , 1.3369 ,
18  5 1.3660 , 1.3758 , 1.4158 , 1.4570 ,
19  6 1.4887 , 1.5431 , 1.5823 , 1.6167 ,
20  7 1.6518 , 1.6877 , 1.7121 , 1.7306 ,
21  8 1.7556 , 1.7873 , 1.8460 , 1.9133 ,
22  9 1.9761 , 2.0046 , 2.0778 , 2.1769 ,
23  a 2.2563 , 2.3981 , 2.5672 , 2.7187 ,
24  b 2.9208 , 3.1718 , 3.4692 , 3.7809 ,
25  c 4.1801 , 4.6049 , 5.1094 , 5.8132 ,
26  d 6.6377 , 7.6611 , 8.7165 , 9.7759 ,
27  e 11.243 , 13.070 , 15.468 , 17.917 ,
28  f 20.239 , 23.360 , 26.865 , 32.956 ,
29  g 39.995 , 52.520 , 68.968 , 93.871 ,
30  h 116.40 , 151.76 , 199.28
31  * /
32  data ( mb(i), i= 1 , np )/
33  1 0.34130e-01, 0.79199e-01, 0.11944 ,
34  2 0.25327 , 0.33197 , 0.45290 , 0.52636 ,
35  3 0.56460 , 0.51589 , 0.43072 , 0.32208 ,
36  4 0.25321 , 0.20721 , 0.18008 , 0.16455 ,
37  5 0.18372 , 0.20513 , 0.22674 , 0.25316 ,
38  6 0.27704 , 0.25313 , 0.22670 , 0.21133 ,
39  7 0.22442 , 0.23831 , 0.21130 , 0.19308 ,
40  8 0.17642 , 0.16282 , 0.15958 , 0.16280 ,
41  9 0.16117 , 0.15483 , 0.14725 , 0.14288 ,
42  a 0.14003 , 0.13451 , 0.13316 , 0.13050 ,
43  b 0.12789 , 0.12533 , 0.12657 , 0.12529 ,
44  c 0.12278 , 0.12154 , 0.12030 , 0.11907 ,
45  d 0.11904 , 0.11901 , 0.11779 , 0.11777 ,
46  e 0.11774 , 0.11537 , 0.11418 , 0.11530 ,
47  f 0.11879 , 0.12116 , 0.12235 , 0.12478 ,
48  g 0.12726 , 0.13108 , 0.13501 , 0.14186 ,
49  h 0.14468 , 0.15357 , 0.15977
50  * /
51 
52  save
53 ! s = (E+Mp)^2-E^2 = 2EMp+Mp^2
54  s = (2*eg+masp)*masp
55  rts = sqrt(s)
56 !//////
57 ! write(*,*) Eg, s, rts
58 !/////////
59  if( rts .lt. 1.08) then
60  xs = 0.
61  elseif(rts .lt. 15.) then
62  call
63  * kpolintplogxyfe(roots, 1, mb, 1, np, m, 3, rts, xs, error)
64  else
65  xs = csigmat('g', 'p', eg)
66  endif
67  end
68  subroutine cggxsec(Eg,xs)
69  implicit none
70  real(8),intent(in):: Eg ! GeV
71  real(8),intent(out):: xs ! xs in mb
72  call cgpxs1(eg, xs)
73  xs = xs /330.
74  end subroutine cggxsec
75 
76 !c These are x-sections near threshold.
77 ! implicit none
78 ! real*8 e, eg, xs
79 !
80 !c to test cgppi0, cgppip, cgppi2, cgppi3
81 ! do e=2.2, 3.7, .02
82 ! eg = 10.**e/1000.
83 ! call cgppi0(e, xs)
84 ! write(*, *) eg, xs
85 ! enddo
86 ! write(*, *)
87 ! do e=2.2, 3.7, .02
88 ! eg = 10.**e/1000.
89 ! call cgppip(e, xs)
90 ! write(*, *) eg, xs
91 ! enddo
92 ! write(*, *)
93 ! do e=2.2, 3.7, .02
94 ! eg = 10.**e/1000.
95 ! call cgppi2(e, xs)
96 ! write(*, *) eg, xs
97 ! enddo
98 ! write(*, *)
99 ! do e=2.2, 3.7, .02
100 ! eg = 10.**e/1000.
101 ! call cgppi3(e, xs)
102 ! write(*, *) eg, xs
103 ! enddo
104 ! end
105  subroutine cgppi0(egl10, xs)
106 ! egl10: input. log10(Eg/MeV)
107 ! xs: output. xsection in micro-barn for gp-->p+pi0
108  implicit none
109  real*8 xs, egl10
110  integer i
111  real*8 xs1(107), xs2(100), xs3(109), xs4(92)
112 ! log10(Eg/MeV) range step 0.01
113  real*8 e11/2.2328224/, e12/3.2981243/, eps/1./
114  real*8 e21/2.3053560/, e22/3.2981243/
115  real*8 e31/2.6146288/, e32/3.7040262/
116  real*8 e41/2.8246450/
117 
118 !
119 ! gp-->ppi0
120  data ( xs1(i),i= 1, 72)/
121  1 0.0, 2.4, 5.8, 8.9, 11.6, 12.8, 13.6, 14.4,
122  2 17.2, 20.9, 25.3, 30.3, 35.9, 42.4, 48.5, 53.4,
123  3 60.9, 71.8, 85.8, 120.9, 150.4, 204.4, 224.5, 241.6,
124  4 251.3, 252.4, 248.9, 233.7, 210.9, 196.1, 186.7, 170.5,
125  5 155.6, 145.1, 140.3, 134.9, 124.5, 111.0, 100.6, 91.3,
126  6 86.8, 77.9, 66.8, 55.4, 46.5, 41.0, 40.0, 39.5,
127  7 40.0, 40.4, 40.8, 41.2, 41.6, 42.0, 42.8, 41.7,
128  8 40.3, 38.8, 37.3, 35.5, 34.2, 33.5, 32.7, 32.8,
129  9 32.4, 32.0, 31.6, 31.2, 30.6, 30.0, 29.7, 28.9/
130  data ( xs1(i),i= 73, 107)/
131  1 27.9, 26.7, 25.3, 23.4, 21.8, 20.3, 18.8, 16.8,
132  2 15.4, 14.2, 13.2, 12.3, 11.9, 11.2, 10.5, 9.8,
133  3 9.0, 8.0, 7.3, 6.8, 6.2, 5.7, 5.3, 4.8,
134  4 4.8, 4.7, 4.7, 4.7, 4.5, 4.1, 4.0, 2.5,
135  5 1.2, 0.7, 0.7/
136 !
137 
138 ! gp-->npi+
139  data ( xs2(i),i= 1, 72)/
140  1 0.0, 59.9, 82.4, 105.5, 127.1, 147.2, 160.8, 166.5,
141  2 175.5, 185.1, 191.6, 194.8, 197.7, 200.1, 201.8, 203.0,
142  3 203.2, 203.4, 202.3, 198.6, 195.0, 193.8, 190.1, 159.7,
143  4 152.0, 143.0, 130.8, 119.4, 109.8, 101.7, 92.2, 81.5,
144  5 73.8, 67.7, 64.9, 61.5, 58.4, 56.5, 56.5, 59.2,
145  6 62.0, 64.1, 68.6, 76.9, 84.7, 90.1, 93.4, 96.1,
146  7 97.6, 98.0, 97.3, 97.0, 96.0, 94.7, 93.0, 91.0,
147  8 88.8, 86.4, 83.6, 80.6, 78.0, 74.3, 69.7, 63.7,
148  9 58.1, 52.8, 47.4, 43.3, 41.2, 38.2, 35.0, 31.2/
149  data ( xs2(i),i= 73, 100)/
150  1 28.3, 25.7, 23.6, 22.1, 20.5, 19.1, 18.5, 17.3,
151  2 16.0, 14.6, 12.8, 11.6, 10.6, 9.8, 9.2, 8.7,
152  3 8.7, 8.7, 8.6, 8.6, 8.2, 7.5, 6.9, 6.2,
153  4 5.7, 4.6, 3.2, 1.3/
154 !
155 ! gp-->ppi+pi-
156  data ( xs3(i),i= 1, 72)/
157  1 0.0, 0.7, 3.0, 5.7, 8.9, 12.9, 16.6, 21.1,
158  2 23.0, 30.9, 44.3, 50.9, 57.5, 63.4, 67.4, 69.9,
159  3 70.7, 71.3, 72.6, 72.4, 72.4, 71.6, 70.3, 68.5,
160  4 67.4, 66.4, 65.3, 64.3, 63.3, 62.3, 61.2, 60.2,
161  5 59.1, 58.0, 56.9, 55.8, 54.8, 53.8, 52.9, 51.9,
162  6 51.0, 50.0, 49.4, 48.4, 47.4, 46.4, 45.4, 44.3,
163  7 43.2, 42.2, 40.7, 39.6, 38.5, 37.5, 36.4, 35.4,
164  8 34.4, 33.4, 32.4, 31.4, 30.4, 29.4, 28.5, 27.6,
165  9 26.7, 25.8, 24.9, 24.0, 22.3, 21.5, 20.7, 20.1/
166  data ( xs3(i),i= 73, 109)/
167  1 19.5, 19.1, 19.2, 19.3, 19.4, 19.9, 19.6, 19.7,
168  2 19.8, 20.0, 19.8, 19.7, 19.1, 18.8, 18.5, 18.3,
169  3 18.0, 17.8, 17.6, 17.1, 16.9, 16.8, 16.7, 16.9,
170  4 17.1, 17.4, 17.8, 18.2, 18.8, 19.4, 20.4, 21.1,
171  5 21.8, 22.4, 23.0, 23.5, 23.9/
172 !
173 
174 ! e42/3.7369471/
175 ! gp-->ppi+pi- x (x=pi+,pi-,,,)
176  data ( xs4(i),i= 1, 72)/
177  1 0.0, 2.6, 5.6, 8.7, 11.8, 14.9, 18.1, 21.2,
178  2 24.4, 28.2, 31.1, 34.0, 36.7, 39.2, 41.6, 43.9,
179  3 45.9, 48.0, 50.1, 52.0, 53.9, 55.6, 57.3, 58.4,
180  4 60.0, 61.6, 63.2, 64.7, 66.2, 67.6, 68.8, 70.2,
181  5 71.6, 73.0, 74.4, 75.8, 77.2, 77.9, 79.3, 80.7,
182  6 82.3, 84.0, 85.9, 87.8, 89.9, 93.2, 95.4, 97.4,
183  7 99.3, 101.2, 103.0, 104.6, 106.7, 108.0, 109.2, 109.9,
184  8 110.3, 111.7, 112.0, 112.3, 111.9, 112.3, 112.3, 112.0,
185  9 111.6, 110.9, 110.1, 109.2, 108.2, 107.0, 105.7, 104.0/
186  data ( xs4(i),i= 73, 92)/
187  1 102.5, 100.9, 99.3, 97.6, 95.7, 93.9, 92.1, 90.3,
188  2 88.5, 86.1, 84.4, 82.8, 81.4, 81.1, 80.3, 79.1,
189  3 77.7, 75.9, 73.9, 71.6/
190 !
191  goto 10
192 ! **********************
193  entry cgppip(egl10, xs)
194 ! **********************
195  goto 20
196 ! **********************
197  entry cgppi2(egl10, xs)
198 ! **********************
199  goto 30
200 ! **********************
201  entry cgppi3(egl10, xs)
202 ! **********************
203  goto 40
204 !
205 
206  10 continue
207  if(egl10 .lt. e11) then
208  xs=0.
209  elseif(egl10 .lt. e12) then
210  call kintp3(xs1, 1, 107, e11, 0.01d0, egl10, xs)
211  else
212  xs=0.
213  endif
214  return
215 !
216  20 continue
217  if(egl10 .lt. e21) then
218  xs=0.
219  elseif(egl10 .lt. e22) then
220  call kintp3(xs2, 1, 100, e21, 0.01d0, egl10, xs)
221  else
222  xs=0.
223  endif
224  return
225 !
226  30 continue
227  if(egl10 .lt. e31) then
228  xs=0.
229  elseif(egl10 .lt. e32) then
230  call kintp3(xs3, 1, 109, e31, 0.01d0, egl10, xs)
231  else
232  xs=0.
233  endif
234  return
235 !
236  40 continue
237  if(egl10 .lt. e41) then
238  xs=0.
239  else
240  call kintp3(xs4, 1, 92, e41, 0.01d0, egl10, xs)
241  endif
242  end
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
dE dx *! Nuc Int sampling table e
Definition: cblkMuInt.h:130
subroutine cgppi0(egl10, xs)
Definition: cgpxs1.f:106
real(8) function, public csigmat(a, b, p)
Definition: cpdgXs.f:31
subroutine cggxsec(Eg, xs)
Definition: cgpxs1.f:69
subroutine kpolintplogxyfe(xa, xstep, ya, ystep, nt, m, logxy, x, y, error)
Definition: kpolintp.f:22
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
nodes a
dE dx *! Nuc Int sampling table b
Definition: cblkMuInt.h:130
dE dx *! Nuc Int sampling table h
Definition: cblkMuInt.h:130
dE dx *! Nuc Int sampling table g
Definition: cblkMuInt.h:130
subroutine kintp3(f, intv, n, x1, h, x, ans)
Definition: kintp3.f:19
masp
Definition: Zmass.h:5
subroutine cgpxs1(Eg, xs)
Definition: cgpxs1.f:3
Definition: cpdgXs.f:1
dE dx *! Nuc Int sampling table f
Definition: cblkMuInt.h:130
dE dx *! Nuc Int sampling table c
Definition: cblkMuInt.h:130