COSMOS v7.655  COSMOSv7655
(AirShowerMC)
ctotx.f
Go to the documentation of this file.
1 ! ctotx: x A collsion total cross section
2 !
3  subroutine ctotx(pj, A, xs)
4  implicit none
5 #include "Zptcl.h"
6 #include "Zcode.h"
7 #include "Zmass.h"
8 #include "Zevhnp.h"
9  type(ptcl)::pj !input projectile particle
10  real(8),intent(in):: A ! input. effective target mass number
11  real(8),intent(out):: xs ! output. in mb. total xsection.
12 ! ctotx0 gives small Xs so we renormalize by
13 ! using cPDGsigmaTotpA(A) at @200 GeV
14  real(8):: cPDGsigmaTotpA
15  type(ptcl):: proton
16  save proton
17  real(8),save:: xst
18  logical,save::first=.true.
19  real(8),save::Asave=-1.
20  real(8),save::ratio
21 
22 
23  call ctotx0(pj, a, xs)
24  if( a /= 1.0d0 ) then
25  if(first) then
26  call cmkptc(knuc, -1, 1, proton)
27  proton.fm.p(4) = 200. + masp
28  first = .false.
29  endif
30  if( a /= asave ) then
31  call ctotx0(proton, a, xst)
32  ratio = cpdgsigmatotpa(a)/xst
33  asave = a
34  endif
35  xs = xs *ratio
36  endif
37  end
38 
39  subroutine ctotx0(pj, A, xs)
40  implicit none
41 #include "Zptcl.h"
42 #include "Zcode.h"
43 #include "Zevhnp.h"
44  type(ptcl)::pj !input projectile particle
45  real*8 A ! input. effective target mass number
46  real*8 xs ! output. in mb. total xsection.
47 
48  real*8 p
49 
50  real(8)::shp, shn
51  if( pj.fm.p(4) .le. pj.mass) then
52  if(pj.code .eq. knuc .and. pj.subcode .eq. antip) then
53  xs = largexs
54  elseif(pj.code .eq. kelec .and. pj.charge .eq. 1) then
55  xs = largexs
56  else
57  xs = smallxs
58  endif
59  else
60  p = sqrt(pj.fm.p(4)**2 - pj.mass**2)
61 ! p = max(p, 0.1d0)
62 ! if( p .ge. 20.) then
63 ! call cerrorMsg('Momentum is >20 GeV/c for ctotx',0)
64 ! endif
65 
66  if(pj.code .eq. knuc) then
67  if(pj.charge .eq. 1) then
68 ! proton
69 ! shp = ctotpp1(p)
70  call cpptotxs(p, shp)
71 
72  if(a .gt. 1.) then
73  call cnptotxs(p, shn)
74 ! use average of pp,pn
75 ! shp =( ctotpn1(p) + shp)/2.0
76  shp = (shp + shn)/2.0
77  call cxp2xaxsec(a, shp, xs)
78  else
79  xs = shp
80  endif
81  elseif(pj.charge .eq. -1) then
82 ! pbar
83 ! shp = ctotpbarp1(p)
84  call cpbarptotxs(p, shp)
85  if( a .gt. 1.) then
86 ! shp = (shp + ctotpbarn1(p))/2.0
87  call cnbarptotxs(p, shn)
88  shp = (shp + shn)/2.
89  call cxp2xaxsec(a, shp, xs)
90  else
91  xs = shp
92  endif
93  elseif(pj.subcode .eq. antip) then
94 ! anti-neutron; assume the same one as pbar
95 ! shp = ctotpbarp1(p)
96  call cnbarptotxs(p, shp)
97  if( a .gt. 1.) then
98  call cxp2xaxsec(a, shp, xs)
99  else
100  xs = shp
101  endif
102  else
103 ! neutron
104 ! shp = ctotnp1(p)
105  call cnptotxs(p, shp)
106  if(a .gt. 1.) then
107  call cxp2xaxsec(a, shp, xs)
108  else
109  xs = shp
110  endif
111  endif
112  elseif(pj.code .eq. kpion) then
113  if(pj.charge .eq. 1) then
114 ! shp = ctotpiPp1(p)
115  call cpipptotxs(p, shp)
116  elseif(pj.charge .eq. -1) then
117 ! shp = ctotpiMp1(p)
118  call cpimptotxs(p, shp)
119  else
120 ! at low energy, pi0 would not interact. any would be o.k
121 ! shp = ctotpiMp1(p)
122  call cpimptotxs(p, shp)
123  endif
124  if(a .ne. 1.0) then
125  call cxp2xaxsec(a, shp, xs)
126  else
127  xs = shp
128  endif
129  elseif(pj.code .eq. kkaon) then
130  if(pj.charge .eq. 1) then
131 ! shp = ctotkPp1(p)
132  call ckpptotxs(p, shp)
133  if(a .gt. 1.) then
134  call ckpntotxs(p, shn)
135  shp = (shp + shn)/2.0
136  call cxp2xaxsec(a, shp, xs)
137  else
138  xs = shp
139  endif
140  elseif(pj.charge .eq. -1) then
141 ! shp = ctotkMp1(p)
142  call ckmptotxs(p, shp)
143  if(a .gt. 1.) then
144  call ckmntotxs(p, shn)
145  shp = (shp + shn)/2.
146  call cxp2xaxsec(a, shp, xs)
147  else
148  xs = shp
149  endif
150  else
151 ! k0; don't worry so much
152 ! xs =( ctotkMp1(p) +ctotkPp1(p))/2.0
153  call ckmptotxs(p, shp)
154  call ckmntotxs(p, shn)
155  shp = (shp+ shn)/2.
156  if(a .gt. 1.) then
157  call cxp2xaxsec(a, shp, xs)
158  else
159  xs = shp
160  endif
161  endif
162  elseif(pj.code .eq. kgzai .or. pj.code .eq. ksigma .or.
163  * pj.code .eq. kbomega .or. pj.code .eq. klambda .or.
164  * pj.code .eq. klambdac) then
165 ! don't worry, not used almost at all; use proton
166 ! shp = ctotpp1(p)
167  call cpptotxs(p, shp)
168  call cxp2xaxsec(a, shp, xs)
169  elseif( pj.code == kgnuc ) then
170  write(0,*) 'Sorry: ctotx is not usable for heavy ions'
171  write(0,*) ' only inelastic xs is used for heavy ions'
172  write(0,*) ' so that you may use cinela or cAAXsec2'
173  stop
174  else
175 ! use pion
176 ! shp =ctotpiPp1(p)
177  call cpipptotxs(p, shp)
178  call cxp2xaxsec(a, shp, xs)
179  endif
180  endif
181  end
max ptcl codes in the kgzai
Definition: Zcode.h:2
subroutine cpbarptotxs(p, xs)
Definition: cpbarpTotXs.f:9
subroutine cpptotxs(p, xs)
Definition: cppTotXs.f:3
max ptcl codes in the kgnuc
Definition: Zcode.h:2
subroutine cnptotxs(p, xs)
Definition: cnpTotXs.f:3
max ptcl codes in the klambdac
Definition: Zcode.h:2
subroutine ckpptotxs(p, xs)
Definition: ckppTotXs.f:3
max ptcl codes in the kkaon
Definition: Zcode.h:2
max ptcl codes in the kelec
Definition: Zcode.h:2
block data cblkElemag data *AnihiE ! Eposi< 1 TeV, anihilation considered *X0/365.667/, ! radiation length of air in kg/m2 *Ecrit/81.e-3/, ! critical energy of air in GeV *MaxComptonE/1./, ! compton is considered below 1 GeV *MaxPhotoE/1.e-3/, ! above this, PhotoElectric effect neg. *MinPhotoProdE/153.e-3/, ! below 153 MeV, no gp --> hadrons ! scattering const not MeV *Knockon true
Definition: cblkElemag.h:7
subroutine cxp2xaxsec(a, xsxp, xsxa)
Definition: cxp2xAXsec.f:30
subroutine ckpntotxs(p, xs)
Definition: ckpnTotXs.f:3
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
*Zfirst p fm *Zfirst p Zfirst p code
Definition: ZavoidUnionMap.h:1
subroutine ctotx(pjin, Ain, xs)
Definition: ctotx.f:7
subroutine ckmptotxs(p, xs)
Definition: ckmpTotXs.f:3
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
max ptcl codes in the klambda
Definition: Zcode.h:2
********************block data cblkHeavy ********************integer j data *HeavyG2symbol *data *HeavyG2code knuc
Definition: cblkHeavy.h:7
subroutine ckmntotxs(p, xs)
Definition: ckmnTotXs.f:2
subroutine cpipptotxs(p, xs)
Definition: cpippTotXs.f:3
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p charge
Definition: ZavoidUnionMap.h:1
block data cblkElemag data *AnihiE ! Eposi< 1 TeV, anihilation considered *X0/365.667/, ! radiation length of air in kg/m2 *Ecrit/81.e-3/, ! critical energy of air in GeV *MaxComptonE/1./, ! compton is considered below 1 GeV *MaxPhotoE/1.e-3/, ! above this, PhotoElectric effect neg. *MinPhotoProdE/153.e-3/, ! below 153 MeV, no gp --> hadrons ! scattering const not MeV *Knockon ! knockon is considered Obsolete *PhotoProd false
Definition: cblkElemag.h:7
subroutine cmkptc(code, subcode, charge, p)
Definition: cmkptc.f:15
masp
Definition: Zmass.h:5
Definition: Zptcl.h:75
max ptcl codes in the kseethru ! subcode integer antip
Definition: Zcode.h:2
*Zfirst p fm *Zfirst p mass
Definition: ZavoidUnionMap.h:1
max ptcl codes in the kpion
Definition: Zcode.h:2
max ptcl codes in the ksigma
Definition: Zcode.h:2
max ptcl codes in the kbomega
Definition: Zcode.h:2
subroutine cpimptotxs(p, xs)
Definition: cpimpTotXs.f:3
*Zfirst p fm *Zfirst p Zfirst p Zfirst p subcode
Definition: ZavoidUnionMap.h:1
subroutine ctotx0(pj, A, Z, xs)
Definition: ctotx.f:66
subroutine cnbarptotxs(p, xs)
Definition: cnbarpTotXs.f:3