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

Go to the source code of this file.

Functions/Subroutines

subroutine cbremang (e, m, eg, z, teta)
 
subroutine cpairang (e, m, teta)
 
subroutine cpbang (a, d, u)
 

Function/Subroutine Documentation

◆ cbremang()

subroutine cbremang ( real*8  e,
real*8  m,
real*8  eg,
real*8  z,
real*8  teta 
)

Definition at line 7 of file cPBA.f.

References cpbang(), parameter(), and true.

Referenced by cbrems().

7  implicit none
8  real*8 e ! input. energy of electron/positron GeV or...
9  real*8 m ! input. electron mass (in the same unit of e)
10  real*8 eg ! input. brems gamma energy. GeV.
11  real*8 z ! input. matter z.
12  real*8 teta ! output. sampled angle of photons relative to
13  ! the inciedent electron
14  real*8 d, pi, hpi, maxu, u
15  real*8 a/0.625d0/
16 !
17  parameter( pi=3.14159265d0, hpi= pi/2)
18 
19  d = 0.13d0*(0.8d0 + 1.3d0/z)*
20  * (100.d0+ 1.d0/e) * (1.+ eg/e)
21 
22 !
23  maxu = e*hpi/m
24  do while (.true.)
25  call cpbang(a, d, u)
26  if(u .lt. maxu) goto 10
27  enddo
28  10 continue
29 !
30  teta = u*m/e
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
nodes z
dE dx *! Nuc Int sampling table e
Definition: cblkMuInt.h:130
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
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
subroutine cpbang(a, d, u)
Definition: cPBA.f:61
! constants thru Cosmos real * pi
Definition: Zglobalc.h:2
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cpairang()

subroutine cpairang ( real*8  e,
real*8  m,
real*8  teta 
)

Definition at line 34 of file cPBA.f.

References cpbang(), parameter(), and true.

Referenced by cpair().

34  implicit none
35 ! samples polar angle of an electron or positron
36 ! from pair creation. It must be a smaller eneregy
37 ! one.
38 ! The angular distribution is approximated as
39 ! { uexp(-au) + duexp(-3au) } du
40 ! where a = 0.625 and d=27.0. This is employed in
41 ! GEANT.
42  real*8 e ! input. energy of e+/e_ (smaller energy)
43  real*8 m ! input. mass of the electron. must be in the same unit as
44  ! e
45 
46  real*8 teta ! output. sampled angle in radian. < pi/2
47 !
48  real*8 a/0.625d0/, d/27.0d0/, pi, hpi, maxu, u
49  parameter( pi=3.14159265d0, hpi= pi/2)
50 !
51  maxu = e*hpi/m
52  do while (.true.)
53  call cpbang(a, d, u)
54  if(u .lt. maxu) goto 10
55  enddo
56  10 continue
57 !
58  teta = u*m/e
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
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
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
subroutine cpbang(a, d, u)
Definition: cPBA.f:61
! constants thru Cosmos real * pi
Definition: Zglobalc.h:2
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cpbang()

subroutine cpbang ( real*8  a,
real*8  d,
real*8  u 
)

Definition at line 61 of file cPBA.f.

References d0, and rndc().

Referenced by cbremang(), and cpairang().

61  implicit none
62 ! sample u from
63 ! { uexp(-au) + duexp(-3au) } du
64 !
65  real*8 a ! input. > 0. see above.
66  real*8 d ! iput. >= 0 see above.
67 !
68  real*8 u ! output. sampled variable u.
69 ! the integration gives
70 !
71 ! 1/a^2 + d/9/a^2 so that the weight is
72 ! 1 : d/9
73 ! 9: d = 9/(9+d): d/(9+d)
74 !
75  real*8 x, x1, x2
76  call rndc(x1)
77  call rndc(x2)
78  u = -log(x1*x2) ! this obeys uexp(-u)du
79  call rndc(x)
80  if(x .lt. 9./(9.+d)) then
81 ! use u exp(-au)du
82  u = u/a
83  else
84 ! use u exp(-3au) du
85  u = u/3.d0/a
86  endif
block data include Zlatfit h c fitting region data x1(1)/0.03/
subroutine rndc(u)
Definition: rnd.f:91
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
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: