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

Go to the source code of this file.

Functions/Subroutines

subroutine ksamppeang (ain, cost)
 
real *8 function ksamppeangf (x)
 

Function/Subroutine Documentation

◆ ksamppeang()

subroutine ksamppeang ( real*8  ain,
real*8  cost 
)

Definition at line 16 of file ksampPEang.f.

References d, kbinchop(), and rndc().

Referenced by cphotoee().

16  implicit none
17 ! This program samples cost=x from the following distribution
18 ! (1-x**2)/(a-x)**4 dx
19 ! where a>1. This is related to the cosine angle of the
20 ! ejected photo-electron at the photo-electric effect.
21 ! (angle is relative to the incident photon).
22 ! This distribution is valid for few keV to *10 keV
23 ! photons.
24 ! a= ((alfa*Z)**2 +2Ee + Eg**2)/2Egsqrt(2Ee)
25 ! where alfa: fin structure const.
26 ! Z: effective atomic number
27 ! Ee: energy of electron /Me
28 ! Eg: // photon /Me
29 !
30  real*8 ain !input "a" value > 1
31  real*8 cost ! output. sampled cos angle
32 
33  real*8 x1, x2, x, eps, ans
34  integer icon
35  real*8 ksamppeangf
36  external ksamppeangf
37  real*8 norm, u, peakpos, am1
38 
39  real*8 a, a1, a2, a3, c1, un
40  common /cksamppeang/a, a1, a2, a3, c1, un
41 
42  a = ain
43  c1 = (1.0-a**2)/3.0
44  a1 = (a+1.0)
45  a2 = a1*a1
46  a3 = a2*a1
47  am1=(a-1.0)
48  norm = c1*(1.0/am1**3 - 1.0/a3) + a*(1./am1**2 -1.0/a2)
49  * - (1.0/am1 - 1.0/a1)
50  call rndc(u)
51  un = norm*u
52 ! peak pos. of the distribution
53  peakpos =(sqrt(a**2 + 8)-a)/2.0
54 
55  x = peakpos
56  x1=-1.
57  x2 =1.0
58  eps = 1.d-4
59  call kbinchop(ksamppeangf,
60  * x1, x2, x, eps, cost, icon)
61 ! write(0,*) ' icon=',icon, ' x=',ans
common ZdedxAir norm
Definition: ZdedxAir.h:2
block data include Zlatfit h c fitting region data x1(1)/0.03/
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
subroutine kbinchop(f, x1, x2, x, eps, ans, icon)
Definition: kbinChop.f:20
subroutine rndc(u)
Definition: rnd.f:91
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
real *8 function ksamppeangf(x)
Definition: ksampPEang.f:66
real(4), save a
Definition: cNRLAtmos.f:20
block data include Zlatfit h c fitting region data x2(1)/0.5/data x1(2)/0.3/
! structure defining a particle at production ! Basic idea of what is to be contained in ! the particle structue is that dynamical ones should be included those derivable from the particle code ! is not included ******************************************************type fmom momentum sequence union map real e endmap map real * x
Definition: Zptcl.h:21
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ksamppeangf()

real*8 function ksamppeangf ( real*8  x)

Definition at line 66 of file ksampPEang.f.

66  implicit none
67  real*8 x
68 
69  real*8 a, a1, a2, a3, c1, un
70  common /cksamppeang/a, a1, a2, a3, c1, un
71 
72  real*8 ax
73 
74  ax = a-x
75  ksamppeangf = c1*(a3-ax**3) +a*ax*a1*(a2-ax**2)
76  * - ax**2*a2*(a1-ax) - un*ax**3*a3
77 
real *8 function ksamppeangf(x)
Definition: ksampPEang.f:66
real(4), save a
Definition: cNRLAtmos.f:20
! structure defining a particle at production ! Basic idea of what is to be contained in ! the particle structue is that dynamical ones should be included those derivable from the particle code ! is not included ******************************************************type fmom momentum sequence union map real e endmap map real * x
Definition: Zptcl.h:21