COSMOS v7.655  COSMOSv7655
(AirShowerMC)
clen2thick.f File Reference

Go to the source code of this file.

Functions/Subroutines

real *8 function clen2thick (z, cosz, leng)
 
real *8 function clen2thickap (z, cosz, leng)
 
real *8 function clen2thickex (z, cosz, leng, n)
 

Function/Subroutine Documentation

◆ clen2thick()

real*8 function clen2thick ( real*8  z,
real*8  cosz,
real*8  leng 
)

Definition at line 27 of file clen2thick.f.

27 !
28 ! z: real*8. input. vertical height in m.
29 ! cosz: real*8. input. cos of zenith angle at z.
30 ! leng: real*8. input. length along cosz direction.
31 ! function value. thickness of air in kg/m2.
32 !
33  implicit none
34 #include "Ztrack.h"
35 #include "Ztrackp.h"
36 #include "Ztrackv.h"
37 #include "Zatmos.h"
38 
39  real*8 z, cosz, leng, clen2thickex, clen2thickap
40  real*8 clen2thickta
41 
42  if(usetbl .and. z .lt. htop) then
44  elseif(exactthick) then
45  clen2thick = clen2thickex(z, cosz, leng, 10)
46  elseif(abs(cosz) .gt. 0.6) then
47  clen2thick = clen2thickap(z, cosz, leng)
48  else
49  clen2thick = clen2thickex(z, cosz, leng, 10)
50  endif
nodes z
real *8 function clen2thickap(z, cosz, leng)
Definition: clen2thick.f:80
real *8 function clen2thick(z, cosz, leng)
Definition: clen2thick.f:27
integer leng
Definition: interface2.h:1
real *8 function clen2thickta(z, leng)
Definition: cl2tTblA.f:5
real *8 function clen2thickex(z, cosz, leng, n)
Definition: clen2thick.f:131

◆ clen2thickap()

real*8 function clen2thickap ( real*8  z,
real*8  cosz,
real*8  leng 
)

Definition at line 80 of file clen2thick.f.

References cnewh().

80 !
81 ! z: real*8. input. vertical height in m.
82 ! cosz: real*8. input. cos of zenith angle at z.
83 ! leng: real*8. input. length along cosz direction.
84 ! function value. thickness of air in kg/m2.
85 !
86  implicit none
87 !---- include 'Zearth.h'
88 #include "Zearth.h"
89 
90  real*8 z, cosz, leng
91 !
92 
93  real*8 cosm, rs, re, cvh2thick, cnewcos, ze, t1, t2,
94  * cnewh
95 
96 
97 
98  rs = z + eradius
99  re = cnewh(rs, cosz, leng)
100  ze = re - eradius
101  cosm = cnewcos(rs, cosz, leng/2)
102  t1 = cvh2thick(z)
103  t2 = cvh2thick(ze)
104  clen2thickap = (t2 - t1)/cosm
nodes z
real *8 function clen2thickap(z, cosz, leng)
Definition: clen2thick.f:80
integer leng
Definition: interface2.h:1
real *8 function cnewcos(H, cost, L)
Definition: catmosutil.f:70
real *8 function cvh2thick(z)
Definition: ciniSegAtoms.f:93
real *8 function cnewh(H, cost, L)
Definition: catmosutil.f:98
Here is the call graph for this function:

◆ clen2thickex()

real*8 function clen2thickex ( real*8  z,
real*8  cosz,
real*8  leng,
integer  n 
)

Definition at line 131 of file clen2thick.f.

References catmosrho(), k16pgaussleg(), and true.

131 !
132 ! z: real*8. input. vertical height in m.
133 ! cosz: real*8. input. cos of zenith angle at z.
134 ! leng: real*8. input. length along cosz direction.
135 ! n: integer.input. how many points for Gauss-Legendre integ.
136 ! function value. thickness of air in kg/m2.
137 !
138  implicit none
139 !---- include 'Zearth.h'
140 #include "Zearth.h"
141 
142  real*8 z, cosz, leng
143  integer n
144 !
145  real*8 seg/1000./, ans, a, b, inte
146 !
147  external catmosrho
148  real*8 catmosrho
149 
150  common /ccatmosrho/ coss, rs
151  real*8 coss, rs
152 !/////
153  real(8),parameter::eps=1.d-5
154  real(8):: error
155  integer icon
156 !///////
157 
158  inte = 0.
159  b =0.
160  coss = cosz
161  rs = z + eradius
162  do while (.true.)
163  a = b
164  b = min(leng, a + seg)
165 ! call kdexpIntF(catmosrho, a, b, eps, ans, error, icon)
166  call k16pgaussleg(catmosrho, a, b, n, ans)
167  inte = inte + ans
168 ! write(*, *) 'a =', a, ' b=',b, ' ans=',ans
169  if(b .eq. leng) goto 10
170  enddo
171  10 continue
172  clen2thickex = inte
nodes z
real *8 function catmosrho(x)
Definition: catmosrho.f:5
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 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
integer leng
Definition: interface2.h:1
subroutine k16pgaussleg(func, a, b, n, ans)
Definition: k16pGaussLeg.f:20
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
nodes a
dE dx *! Nuc Int sampling table b
Definition: cblkMuInt.h:130
real *8 function clen2thickex(z, cosz, leng, n)
Definition: clen2thick.f:131
integer n
Definition: Zcinippxc.h:1
Here is the call graph for this function: