COSMOS v7.655  COSMOSv7655
(AirShowerMC)
cinteElec.f
Go to the documentation of this file.
1  subroutine cinteelec
2  implicit none
3 !---- include 'Ztrack.h'
4 #include "Ztrack.h"
5 !---- include 'Ztrackv.h'
6 #include "Ztrackv.h"
7  character*70 msg
8 !
9  if(intinfarray(processno)%process .eq. 'brems') then
10  call cbrems
11  elseif(intinfarray(processno)%process .eq. 'mscat') then
12  call cknockon(-1)
13  elseif(intinfarray(processno)%process .eq. 'bscat') then
14  call cknockon(1)
15  elseif(intinfarray(processno)%process .eq. 'anihi') then
16  call canihi
17  elseif(intinfarray(processno)%process .eq. 'mbrem') then
18  call cmbrem
19  else
20 
21  write(msg, *) ' process#=', processno,
22  * 'and process for electron=',
23  * intinfarray(processno)%process, ' undef'
24  call cerrormsg(msg, 0)
25  endif
26  end
27 !
28 ! Pwork assumes that higher energy particles are normally
29 ! in the last part, and stacked later from the last to
30 ! save the stack area. Due to this fact, putting ptlcs
31 ! in Pwork must be reversed (realized from v6.10)
32 ! *****************
33  subroutine cbrems
34 ! *****************
35  implicit none
36 
37 #include "Zcode.h"
38 #include "Ztrackp.h"
39 #include "Ztrack.h"
40 #include "Ztrackv.h"
41 
42 
43  real*8 beg, den, cvh2den
44  real*8 theta, sint, cost, sn, cs
45  type(coord)::dc1
46  type(coord)::w
47  type(track)::aTrack
48 
49 ! get brem gamma energy
50  if(lpmeffect .and.
51  * trackbefmove%p%fm%p(4) .gt. lpmbrememin) then ! use BefMove%p
52 ! den = cthick2den(TrackBefMove.pos.depth)
53  den = cvh2den(trackbefmove%pos%height) ! this should be better
54  call cbremerglpm(trackbefmove%p%fm%p(4), den, beg)
55  else
56  call cbremsenergy(movedtrack%p%fm%p(4), beg)
57  endif
58 
59  atrack = movedtrack
60 
61 !
62  atrack%p%fm%p(4) = movedtrack%p%fm%p(4) -beg
63  call ce2p(atrack)
64 !
65 ! see !! note above
66  pwork(nproduced+2) = atrack%p
67 
68 
69  atrack%p%fm%p(4) = beg
70  call cmkptc(kphoton, kcasg, 0, atrack%p)
71 ! brems angle at low energy
72  if( trackbefmove%p%fm%p(4) .lt. 10. ) then
73  call cbremang( trackbefmove%p%fm%p(4), trackbefmove%p%mass,
74  * beg, targetatomicn, theta)
75  ! v7.646; v8.03 revived
76 ! call cBremAng(1, TrackBefMove%p%fm%p(4), beg, theta) removed
77  if(theta .lt. 0.03d0) then
78  sint = theta
79  cost = 1.- theta**2 / 2
80  else
81  sint = sin(theta)
82  cost = cos(theta)
83  endif
84  call kcossn(cs,sn)
85  w%r(1) = cs*sint
86  w%r(2) = sn*sint
87  w%r(3) = cost
88  call ctransvectz(movedtrack%vec%w, w, dc1)
89  call csetdircos(dc1, atrack)
90  endif
91  call ce2p(atrack)
92  pwork(nproduced+1) = atrack%p
93  nproduced = nproduced + 2
94  end
95 ! *****************
96  subroutine cknockon(cg)
97 ! ***************** moller and bhabha scattering
98  implicit none
99 !---- include '../Particle/Zcode.h'
100 #include "Zcode.h"
101 !---- include 'Ztrack.h'
102 #include "Ztrack.h"
103 !---- include 'Ztrackv.h'
104 #include "Ztrackv.h"
105 #include "Zelemagp.h"
106 !
107  integer cg
108  real*8 e1, er, tmp, cos1, cosr, sine, cs, sn, sinr
109  type(coord)::dc
110  type(coord)::dc1
111  type(coord)::dcr
112 !
113  type(track)::aTrack
114 !
115  if(cg .eq. -1) then
116  call cmollerea(movedtrack%p%fm%p(4), recoilkinemine,
117  * e1, er, cos1, cosr)
118  else
119  call cbhabhaea(movedtrack%p%fm%p(4), recoilkinemine,
120  * e1, er, cos1, cosr)
121  endif
122 
123  tmp=1.d0-cos1*cos1
124  if(tmp .lt. 0.d0) then
125  tmp=0.d0
126  cos1=1.d0
127  endif
128  sine=sqrt(tmp)
129  call kcossn(cs,sn)
130  dc%r(1) = cs*sine
131  dc%r(2) = sn*sine
132  dc%r(3) = cos1
133 ! convert angle
134  call ctransvectz(movedtrack%vec%w, dc, dc1)
135  atrack = movedtrack
136  atrack%p%fm%p(4) = e1
137  call csetdircos(dc1, atrack)
138  call ce2p(atrack)
139 
140  pwork(nproduced + 2) = atrack%p
141 ! knock on electron
142  tmp=1.d0-cosr*cosr
143  if(tmp .lt. 0.d0) then
144  tmp=0.d0
145  cosr=1.d0
146  endif
147  sinr=sqrt(tmp)
148  dc%r(1) = -cs*sinr
149  dc%r(2) = -sn*sinr
150  dc%r(3) = cosr
151  call ctransvectz(movedtrack%vec%w, dc, dcr)
152  atrack%p%fm%p(4) = er
153  call csetdircos(dcr, atrack)
154 !
155  if(cg .eq. 1) then
156  call cmkptc(kelec, 0, -1, atrack%p)
157  endif
158 
159  call ce2p(atrack)
160 
161  pwork(nproduced + 1) = atrack%p
162 
163  nproduced = nproduced + 2
164  end
165 ! *****************
166  subroutine canihi
167 ! *****************
168  implicit none
169 !---- include '../Particle/Zcode.h'
170 #include "Zcode.h"
171 !---- include 'Ztrack.h'
172 #include "Ztrack.h"
173 !---- include 'Ztrackv.h'
174 #include "Ztrackv.h"
175 
176  real*8 e1, er, cos1, cosr, tmp, sine
177  real*8 cs, sn, sinr
178  type(track)::aTrack
179  type(coord)::dc
180  type(coord)::dc1
181  type(coord)::dcr
182 !
183  call canihiea(movedtrack%p%fm%p(4), e1, er, cos1, cosr)
184 !
185  tmp=1.d0-cos1*cos1
186  if(tmp .lt. 0.d0) then
187  tmp=0.d0
188  cos1=-1.d0
189  endif
190  sine=sqrt(tmp)
191  call kcossn(cs,sn)
192  dc%r(1) = cs*sine
193  dc%r(2) = sn*sine
194  dc%r(3) = cos1
195 !
196  call ctransvectz(movedtrack%vec%w, dc, dc1)
197 ! higher gamma
198  atrack = movedtrack
199  call csetdircos(dc1, atrack)
200  call cmkptc(kphoton, kcasg, 0, atrack%p)
201  atrack%p%fm%p(4) = e1
202  call ce2p(atrack)
203 
204 
205  nproduced = nproduced + 1
206  pwork(nproduced) = atrack%p
207 ! call cpush(aTrack)
208 ! low gamma
209  tmp=1.d0-cosr*cosr
210  if(tmp .lt. 0.d0) then
211  tmp=0.d0
212  cosr=-1.d0
213  endif
214  sinr=sqrt(tmp)
215  dc%r(1) = -cs*sinr
216  dc%r(2) = -sn*sinr
217  dc%r(3) = cosr
218  call ctransvectz(movedtrack%vec%w, dc, dcr)
219  call csetdircos(dcr, atrack)
220  atrack%p%fm%p(4) = er
221  call ce2p(atrack)
222 
223 
224  nproduced = nproduced + 1
225  pwork(nproduced) = atrack%p
226 ! call cpush(aTrack)
227  end
228 
229 ! *****************
230  subroutine cmbrem
231 ! magnetic brems
232 ! *****************
233  implicit none
234 #include "Zcode.h"
235 #include "Ztrack.h"
236 #include "Ztrackv.h"
237 
238  real*8 x, beg
239 
240  type(track)::aTrack
241 
242 ! get brem gamma energy
243  call cmbreme(upsilon, x)
244  beg = x * movedtrack%p%fm%p(4) ! gamma energy
245 
246  atrack = movedtrack
247  atrack%p%fm%p(4) = movedtrack%p%fm%p(4) -beg
248  call ce2p(atrack)
249  pwork( nproduced + 2 ) = atrack%p
250 ! stack gamma (probably lower energy than e)
251  atrack%p%fm%p(4) = beg
252  call cmkptc(kphoton, kcasg, 0, atrack%p)
253  call ce2p(atrack)
254  pwork(nproduced + 1) = atrack%p
255 
256  nproduced = nproduced + 2
257  end
258 
259 
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
subroutine cbremsenergy(Ee, Eg)
Definition: cbremsPath.f:42
! for length to thickness conversion or v v ! integer maxnodes real Hinf ! This is used when making table for dim simulation ! The slant length for vertical height to km is ! divided by LenStep m steps ! It can cover the slant length of about km for cos
Definition: Zatmos.h:8
subroutine ce2p(aTrack)
Definition: ce2p.f:5
const int kphoton
Definition: Zcode.h:6
Definition: Ztrack.h:44
subroutine cbremang(e, m, eg, z, teta)
Definition: cPBA.f:7
max ptcl codes in the kelec
Definition: Zcode.h:2
! Parameters used for hadronic cascade shower is generated newline ! For you may give as as or em quick generation of AS for heavy primaries is tried See chookASbyH f character *Generate2 don t touch this for skeleton flesh use integer MagBrem no magnetic bremsstrahlung is considered newline ! if and Ee energy loss due to magnetic brems is considered newline ! if and Ee real sampling of gamma is performed WaitRatio ! must be made small so that WaitRatio *E0 sim MagBremEmin integer MagPair no magnetic pair creation is considered newline ! if and Eg real sampling is the LPM effect is considered when Ee LpmBremEmin for electrons and ! Eg LpmPairEmin for gamma rays real *MagBremEmin E magnetic bremsstrahlung by electrons may be considered if not considered at all newline total energy loss due to brems is considered newline gamma energy is sampled actually newline ! If upsilon(Ee/m *B/Bcr) is small
max ptcl codes in the kseethru ! subcode integer kcasg
Definition: Zcode.h:2
subroutine cmollerea(ein, w, es, er, coss, cosr)
Definition: cmollerPath.f:49
subroutine canihi
Definition: cinteElec.f:167
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
subroutine ctransvectz(zax, dir1, dir2)
Definition: ctransVectZ.f:21
subroutine cbrems
Definition: cinteElec.f:34
subroutine csetdircos(dc, aTrack)
Definition: cgetZenith.f:4
subroutine cknockon(cg)
Definition: cinteElec.f:97
subroutine cmbreme(up, x)
Definition: cmBremE.f:2
subroutine kcossn(cs, sn)
Definition: kcossn.f:13
subroutine cbremerglpm(ee, rho, eg)
Definition: cbremErgLPM.f:20
subroutine cinteelec
Definition: cinteElec.f:2
subroutine cmkptc(code, subcode, charge, p)
Definition: cmkptc.f:15
subroutine canihiea(ein, eg1, eg2, cos1, cos2)
Definition: canihiPath.f:60
subroutine cbhabhaea(ein, w, es, er, coss, cosr)
Definition: cbhabhaPath.f:46
Definition: Zcoord.h:43
subroutine cmbrem
Definition: cinteElec.f:231