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

Go to the source code of this file.

Functions/Subroutines

subroutine ccomptpath (Eg, p, path)
 
subroutine ccomptea (Eg, Egout, Ee, cosg, cose)
 

Function/Subroutine Documentation

◆ ccomptea()

subroutine ccomptea ( real*8  Eg,
real*8  Egout,
real*8  Ee,
real*8  cosg,
real*8  cose 
)

Definition at line 50 of file ccomptPath.f.

References d, d0, ksamplin(), masele, rndc(), and true.

Referenced by ccompt().

50  implicit none
51 #include "Zmass.h"
52 ! ***********
53  real*8 eg ! input. gamma energy in GeV.
54  real*8 egout ! output. scattered gamma energy in GeV
55  real*8 ee ! output. scattered electron energy in GeV
56  real*8 cosg ! output. cos of scattered gamma angle
57  real*8 cose ! output. cos of scattered electron energy
58 
59  real*8 xmin, a1, a2, x, temp, u, sin2g, cos2e, g
60 ! x = Egout/Eg
61  g =eg/masele
62  xmin = 1.d0/( 1.d0 + 2*g)
63  a1 = - log(xmin)
64  a2 = (1.d0 -xmin*xmin) /2
65  do while (.true.)
66  call rndc(u)
67  if( u .lt. a1/(a1+a2) ) then
68 ! sample form 1/x dx
69  call rndc(u)
70  x = xmin * exp( a1 * u )
71  else
72 ! sample from x dx in (xmin~1)
73  call ksamplin(1.d0, 0.d0, xmin, 1.d0, x)
74  endif
75 ! rejection by (1- xsin^2/(1+x^2)
76  temp = (1 - x)/x/g
77  sin2g = temp*(2.d0-temp)
78  call rndc(u)
79  if(u .lt. (1. - x*sin2g/(1+x*x))) goto 10
80  enddo
81  10 continue
82  egout = eg*x
83  ee = eg- egout + masele
84 ! give cos of phton and electron
85  cosg= 1. - temp
86 ! tan(el)=cot(gm/2)/(1+g);
87 ! cot(t/2)= +-sqrt( (1+cos(t))/(1-cos(t)) ) so that
88  cos2e=(1.d0-cosg) / ( 1.d0-cosg +(1.d0+cosg)/(1.d0+g)**2 )
89 ! electron angle is always 0 to 90 deg.
90  cose=min( max(sqrt(cos2e), 1.d-10), 0.9999999999d0)
91 
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
dE dx *! Nuc Int sampling table d
Definition: cblkMuInt.h:130
const double masele
Definition: Zmass.h:2
subroutine ksamplin(a, b, alpha, beta, x)
Definition: ksampLin.f:15
dE dx *! Nuc Int sampling table g
Definition: cblkMuInt.h:130
! 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:

◆ ccomptpath()

subroutine ccomptpath ( real*8  Eg,
real*8  p,
real*8  path 
)

Definition at line 18 of file ccomptPath.f.

References d0, masele, and rndc().

Referenced by csampgintl().

18  implicit none
19 #include "Zmass.h"
20  real*8 eg ! input. Gamma energy in GeV
21  real*8 p ! output. probability ( number of
22  ! occurence ) of compton scattering per r.l
23  real*8 path ! output. sampled path in r%l
24 
25 
26  real*8 u, g
27 ! tp* 3/8 /cconst= total cross section normalized
28 ! by thomson cross section.
29 ! cconst=3/8 * thomson * n0 * z/a * x0ing
30 !
31 
32  real*8 cconst/2.7429630/
33  save cconst
34 !
35  g =eg/masele
36  if(g .lt. 0.1d0) then
37 ! p=( (5.2*g-2.)*g +1. )*2.66666*media.basearea
38  p=( (5.2*g-2.)*g +1. )*2.66666*cconst
39  else
40  p=( (1. - (g+1.)*2/g**2)*log(g*2+1.) + .5 + 4./g -
41  * 1./(g*2+1.)**2/2 ) /g * cconst
42  endif
43  call rndc(u)
44  path = - log(u)/p
subroutine rndc(u)
Definition: rnd.f:91
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
const double masele
Definition: Zmass.h:2
dE dx *! Nuc Int sampling table g
Definition: cblkMuInt.h:130
Here is the call graph for this function:
Here is the caller graph for this function: