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

Go to the source code of this file.

Functions/Subroutines

subroutine cpairerglpm (eg, rho, ee)
 

Function/Subroutine Documentation

◆ cpairerglpm()

subroutine cpairerglpm ( real*8  eg,
real*8  rho,
real*8  ee 
)

Definition at line 20 of file cpairErgLPM.f.

References d0, masele, parameter(), rndc(), and true.

Referenced by cpair().

20  implicit none
21 #include "Zmass.h"
22  real*8 eg !input. gamma energy in GeV
23  real*8 rho !input. air density in kg/m3
24  real*8 ee !output. sampled electron or positron energy in GeV
25 !
26  real*8 u, v, temp, s, w, u1, v1, w1, x
28 
29  real*8 a1, a2, suma, a1prob
30  parameter( a1 = 2./3., a2 = 1./9.,
31  * suma = a1 + a2, a1prob = a1/suma)
32 !
33  do while (.true.)
34  call rndc(u)
35  if(u .lt. a1prob) then
36 ! sampling by dv
37  call rndc(v)
38 ! rejection
39 ! get s
40  s = cmigdsforpair(eg, rho, v)
41  temp = cmigdgzai(s) * (cmigdg(s) +cmigdpsi(s))/2.
42  call rndc(u)
43  if(u .lt. temp) then
44 ! accepted
45  ee =max(min(v, 1.d0-masele/eg),masele/eg) * eg
46  goto 100
47  endif
48  else
49 ! sampling by 12(v-0.5)dv
50  call rndc(u)
51  call rndc(v)
52  call rndc(w)
53 ! get one which is most far from 1/2
54  u1 = abs(u - 0.5)
55  v1 = abs(v - 0.5)
56  w1 = abs(w - 0.5)
57  x = max(u1, v1, w1)
58  if(x .eq. u1) then
59  v = u
60  elseif(x .eq. w1) then
61  v = w
62  endif
63 ! rejection, get s
64  s = cmigdsforpair(eg, rho, v)
65  temp = cmigdgzai(s) * cmigdpsi(s)
66  call rndc(u)
67  if(u .lt. temp) then
68 ! accepted
69  ee =max(min(v, 1.d0-masele/eg), masele/eg) * eg
70  goto 100
71  endif
72  endif
73  enddo
74  100 continue
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
real *8 function cmigdsforpair(eg, rho, v)
Definition: cmigdFunc.f:214
real *8 function cmigdg(s)
Definition: cmigdFunc.f:52
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 rndc(u)
Definition: rnd.f:91
real(4), dimension(:), allocatable, save temp
Definition: cNRLAtmos.f:29
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
real *8 function cmigdgzai(s)
Definition: cmigdFunc.f:146
const double masele
Definition: Zmass.h:2
real *8 function cmigdpsi(s)
Definition: cmigdFunc.f:103
! 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: