COSMOS v7.655  COSMOSv7655
(AirShowerMC)
modpdgxs Module Reference

Functions/Subroutines

real(8) function, public csofstu (p, Mp, Mt)
 
real(8) function, public csigmat (a, b, p)
 

Variables

real(8), parameter pi = 3.14159265358979323846d0
 
real(8), parameter hbarc =197.327d-3
 
real(8), parameter, public m = 2.15
 
real(8), parameter eta1 = 0.462
 
real(8), parameter eta2 = 0.55
 
real(8), parameter lb = pi*(hbarc/M)**2*10.0d0
 
real(8), parameter delta = 0.003056
 
real(8), parameter lambda = 1.630
 

Function/Subroutine Documentation

◆ csigmat()

real(8) function, public modpdgxs::csigmat ( character(*), intent(in)  a,
character(*), intent(in)  b,
real(8), intent(in)  p 
)

Definition at line 31 of file cpdgXs.f.

References csofstu(), d, d0, delta, eta1, eta2, lambda, lb, m, masd, maskc, masn, masp, maspic, masrho, and massigmam.

Referenced by cgpxs1(), ckmptotxs(), ckpntotxs(), ckpptotxs(), cnptotxs(), cpbarptotxs(), cpimptotxs(), cpipptotxs(), and cpptotxs().

31 ! total cross-section by pdg (ed. 2012)
32  implicit none
33 #include "Zmass.h"
34  character(*),intent(in):: a ! one of "p", "pb", "pi+","pi-"
35  ! "pi0"
36  ! "K+","K-" "K0" "Sig", "g"
37  character(*),intent(in):: b ! one of "p", "n", "d"
38  real(8),intent(in)::p ! momentum of a in GeV
39  real(8):: ans ! total cross section in mb.
40 
41  real(8):: lbt, sm, s
42 
43  real(8):: z, y1, y2
44 
45  lbt = lb
46  select case(a)
47  case('p', 'pb')
48  select case(b)
49  case('p')
50  z = 34.71
51  y1 = 12.72
52  y2 = 7.35
53  sm = (masp *2 + m)**2
54  s = csofstu(p,masp, masp)
55  case('n')
56  z = 35.00
57  y1 = 12.19
58  y2 = 6.62
59  sm = (masp + masn + m)**2
60  s = csofstu(p, masp, masn)
61  case('d')
62  z = 65.02
63  y1= 29.04
64  y2= 14.9
65  sm = (masp + masd + m)**2
66  s = csofstu(p, masp, masd)
67  lbt =lbt*lambda
68  case default
69  write(0,*)
70  * 'b ', b, ' for csigma(a,b) is invalid'
71  stop
72  end select
73  case('pi+', 'pi-', 'pi0')
74  select case(b)
75  case('p','n')
76  z = 19.02
77  y1 = 9.22
78  y2= 1.75
79  sm = (maspic + masn + m)**2
80  s = csofstu(p, maspic, masp)
81  case('d')
82  z = 37.06
83  y1 = 18.28
84  y2 = 0.34
85  sm = (maspic + masd + m)**2
86  s = csofstu(p, maspic, masd)
87  lbt =lbt*lambda
88  case default
89  write(0,*)
90  * 'b ', b, ' for csigma(a,b) is invalid'
91  write(0,*) ' a is ',a
92  stop
93  end select
94  case('K+', 'K-', 'K0')
95  select case(b)
96  case('p')
97  z = 16.55
98  y1 = 4.02
99  y2 = 3.39
100  sm = (maskc + masp + m)**2
101  s = csofstu(p, maskc, masp)
102  case('n')
103  z = 16.49
104  y1 = 3.44
105  y2 = 1.82
106  sm = (maskc + masn + m)**2
107  s = csofstu(p, maskc, masn)
108  case('d')
109  z = 32.34
110  y1 = 7.33
111  y2 = 5.59
112  sm = (maskc + masd + m)**2
113  s = csofstu(p, maskc, masd)
114  lbt =lbt*lambda
115  case default
116  write(0,*)
117  * 'b ', b, ' for csigma(a,b) is invalid'
118  write(0,*) ' a is ',a
119  stop
120  end select
121  case('Sig')
122  select case(b)
123  case('p')
124  z = 34.9
125  y1 = -55.0
126  y2 = -57.0
127  sm = (massigmam + masp + m)**2
128  s = csofstu(p, massigmam, masp)
129  case default
130  write(0,*)
131  * 'b ', b, ' for csigma(a,b) is invalid'
132  write(0,*) ' a is ',a
133  stop
134  end select
135  case('g')
136  select case(b)
137  case('p','n')
138  z = 34.71
139  y1 = 0.0128
140  y2= 0.
141  sm = (masp+ masp + m)**2
142  s = csofstu(p, masp, masp)
143  ans = (z + lbt* (log(s/sm))**2)* delta
144  sm = (masp + m)**2
145  s = csofstu(p, 0.d0, masp)
146  ans = ans + y1*(sm/s)**eta1
147  return
148  case('d')
149  z = 65.02
150  y1 = 0.0128
151  y2= 0.
152  sm = (masp + masd + m)**2
153  s = csofstu(p, masp, masd)
154  ans = (z + lambda* lbt* (log(s/sm))**2 )* delta
155  sm = ( masd + m)**2
156  s = csofstu(p, 0.d0, masd)
157  ans = ans + y1*(sm/s)**eta1
158  return
159  case('g')
160  z = 34.71
161  y1 = -0.034d-4
162  y2 = 0.
163  sm = (masp + masp + m)**2
164  s = csofstu(p, masp, masp)
165  ans = (z + lbt* (log(s/sm))**2 )* delta**2
166  sm = ( masrho + m)**2
167  s = csofstu(p, 0.d0, masrho)
168  ans = ans + y1*(sm/s)**eta1
169  return
170  case default
171  write(0,*)
172  * 'b ', b, ' for csigma(a,b) is invalid'
173  write(0,*) ' a is ',a
174  stop
175  end select
176  case default
177  write(0,*)
178  * 'a ', a, ' for csigma(a,b) is invalid'
179  write(0,*) ' b is ',b
180  stop
181  end select
182  ans = z + lbt* (log(s/sm))**2 + y1*(sm/s)**eta1
183  if( a == "pb" .or. a =="pi-" .or. a == "K-") then
184  ans = ans + y2*(sm/s)**eta2
185  else
186  ans = ans - y2*(sm/s)**eta2
187  endif
nodes z
atmos%rho(atmos%nodes) **exp(-(z-atmos%z(atmos%nodes))/Hinf) elseif(z .lt. atmos%z(1)) then ans=atmos%rho(1) **exp((atmos%z(1) -z)/atmos%H(1)) else call kdwhereis(z, atmos%nodes, atmos%z, 1, i) a=atmos%a(i) if(a .ne. 0.d0) then ans=atmos%rho(i) **(1+a *(z-atmos%z(i))/atmos%H(i)) **(-1.0d0-1.d0/a) else ans=*atmos%rho(i) *exp(-(z-atmos%z(i))/atmos%H(i)) endif endif ! zsave=z ! endif cvh2den=ans end ! ---------------------------------- real *8 function cvh2temp(z) implicit none ! vettical height to temperatur(Kelvin) real *8 z ! input. vertical height in m ! output is temperature of the atmospher in Kelvin real *8 ans integer i if(z .gt. atmos%z(atmos%nodes)) then ans=atmos%T(atmos%nodes) elseif(z .lt. atmos%z(1)) then ans=atmos%T(1)+atmos%b(1) *(z - atmos%z(1)) else call kdwhereis(z, atmos%nodes, atmos%z, 1, i) ans=atmos%T(i)+atmos%b(i) *(z-atmos%z(i)) endif cvh2temp=ans end !--------------------------------------------- real *8 function cthick2h(t) implicit none real *8 t ! input. air thickness in kg/m^2 real *8 logt, ans integer i real *8 dod0, fd, a logt=log(t) if(t .ge. atmos%cumd(1)) then ans=atmos%z(1) - *(logt - atmos%logcumd(1)) *atmos%H(1) elseif(t .le. atmos%cumd(atmos%nodes)) then ans=atmos%z(atmos%nodes) - *Hinf *log(t/atmos%cumd(atmos%nodes)) else call kdwhereis(t, atmos%nodes, atmos%cumd, 1, i) ! i is such that X(i) > x >=x(i+1) ans
masn
Definition: Zmass.h:5
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
masd
Definition: Zmass.h:5
maskc
Definition: Zmass.h:5
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
masrho
Definition: Zmass.h:5
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
real(8) function, public csofstu(p, Mp, Mt)
Definition: cpdgXs.f:21
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data m
Definition: ZlatfitBD.h:35
real(4), save a
Definition: cNRLAtmos.f:20
masp
Definition: Zmass.h:5
real(4), save b
Definition: cNRLAtmos.f:21
massigmam
Definition: Zmass.h:5
maspic
Definition: Zmass.h:5
Here is the call graph for this function:
Here is the caller graph for this function:

◆ csofstu()

real(8) function, public modpdgxs::csofstu ( real(8), intent(in)  p,
real(8), intent(in)  Mp,
real(8), intent(in)  Mt 
)

Definition at line 21 of file cpdgXs.f.

Referenced by cnpelaxs(), cnptotxs(), cpippelaxs(), cppelaxs(), cpptotxs(), and csigmat().

21 ! s = (E+Mt)^2 - p**2 = 2*E*Mt + Mt**2 + Mp**2
22  implicit none
23  real(8),intent(in):: p ! momentum in GeV
24  real(8),intent(in):: mp ! proj. mass
25  real(8),intent(in):: mt ! target mass
26  real(8):: ans
27  ans = (sqrt(p**2+mp**2)*2 + mt)*mt + mp**2
atmos%rho(atmos%nodes) **exp(-(z-atmos%z(atmos%nodes))/Hinf) elseif(z .lt. atmos%z(1)) then ans=atmos%rho(1) **exp((atmos%z(1) -z)/atmos%H(1)) else call kdwhereis(z, atmos%nodes, atmos%z, 1, i) a=atmos%a(i) if(a .ne. 0.d0) then ans=atmos%rho(i) **(1+a *(z-atmos%z(i))/atmos%H(i)) **(-1.0d0-1.d0/a) else ans=*atmos%rho(i) *exp(-(z-atmos%z(i))/atmos%H(i)) endif endif ! zsave=z ! endif cvh2den=ans end ! ---------------------------------- real *8 function cvh2temp(z) implicit none ! vettical height to temperatur(Kelvin) real *8 z ! input. vertical height in m ! output is temperature of the atmospher in Kelvin real *8 ans integer i if(z .gt. atmos%z(atmos%nodes)) then ans=atmos%T(atmos%nodes) elseif(z .lt. atmos%z(1)) then ans=atmos%T(1)+atmos%b(1) *(z - atmos%z(1)) else call kdwhereis(z, atmos%nodes, atmos%z, 1, i) ans=atmos%T(i)+atmos%b(i) *(z-atmos%z(i)) endif cvh2temp=ans end !--------------------------------------------- real *8 function cthick2h(t) implicit none real *8 t ! input. air thickness in kg/m^2 real *8 logt, ans integer i real *8 dod0, fd, a logt=log(t) if(t .ge. atmos%cumd(1)) then ans=atmos%z(1) - *(logt - atmos%logcumd(1)) *atmos%H(1) elseif(t .le. atmos%cumd(atmos%nodes)) then ans=atmos%z(atmos%nodes) - *Hinf *log(t/atmos%cumd(atmos%nodes)) else call kdwhereis(t, atmos%nodes, atmos%cumd, 1, i) ! i is such that X(i) > x >=x(i+1) ans
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
Here is the caller graph for this function:

Variable Documentation

◆ delta

real(8), parameter delta = 0.003056
private

Definition at line 16 of file cpdgXs.f.

Referenced by csigmat().

16  real(8),parameter:: delta = 0.003056, lambda = 1.630

◆ eta1

real(8), parameter eta1 = 0.462
private

Definition at line 14 of file cpdgXs.f.

Referenced by csigmat().

14  real(8),parameter:: eta1 = 0.462, eta2 = 0.55

◆ eta2

real(8), parameter eta2 = 0.55
private

Definition at line 14 of file cpdgXs.f.

Referenced by csigmat().

◆ hbarc

real(8), parameter hbarc =197.327d-3
private

Definition at line 12 of file cpdgXs.f.

12  real(8),parameter:: hbarc=197.327d-3 ! GeV%Fermi
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130

◆ lambda

real(8), parameter lambda = 1.630
private

Definition at line 16 of file cpdgXs.f.

Referenced by csigmat().

◆ lb

real(8), parameter lb = pi*(hbarc/M)**2*10.0d0
private

Definition at line 15 of file cpdgXs.f.

Referenced by csigmat().

15  real(8),parameter:: lb = pi*(hbarc/m)**2*10.0d0 ! mb
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
! constants thru Cosmos real * pi
Definition: Zglobalc.h:2
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data m
Definition: ZlatfitBD.h:35

◆ m

real(8), parameter, public m = 2.15

Definition at line 13 of file cpdgXs.f.

Referenced by ckmpelaxs(), ckpnelaxs(), ckppelaxs(), cpimpelaxs(), cpippelaxs(), and csigmat().

13  real(8),parameter:: m = 2.15 ! GeV
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data m
Definition: ZlatfitBD.h:35

◆ pi

real(8), parameter pi = 3.14159265358979323846d0
private

Definition at line 11 of file cpdgXs.f.

11  real(8),parameter:: pi= 3.14159265358979323846d0
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
! constants thru Cosmos real * pi
Definition: Zglobalc.h:2