COSMOS v7.655  COSMOSv7655
(AirShowerMC)
Gencol3.f
Go to the documentation of this file.
1 #include "Zintmodel.h"
2 #include "ZcosmosBD.h"
3  implicit none
4 #include "Zcode.h"
5 #include "Zptcl.h"
6 #include "Ztrackp.h"
7  include "Zprivate.h"
8 
9  integer i, nev, j, ntp, eof, ntpout
10  type(ptcl):: w(maxn)
11 ! write(*,*) ' this is dummy'
12  integer,save:: Jdecay(klast)
13 !///////////
14 ! real(8),save:: tau=0.3e-10
15  real(8),save:: tau=0.95e-10 ! k0s decay
16  jdecay(:) = 0
17 ! make Jdecay non zero to see if it is needed to decay,if
18 ! the life is < tau.
19 ! data jcode/kbomega, kdmes, keta, kgzai, klambda, klambdac,
20 ! * kpion, ksigma, kkaon/
21 ! kpion--> pi0 (pi+/- is never decayed)
22 ! kkaon--> k0s (only k0s is examined)
23  jdecay(kbomega) = 1
24  jdecay(kdmes) = 1
25  jdecay(kgzai) = 1
26  jdecay(klambda) = 1
27  jdecay(klambdac) = 1
28  jdecay(ksigma) = 1
29  jdecay(kkaon) = 1
30  jdecay(keta) = 1
31  jdecay(krho) = 1
32  jdecay(kphi) = 1
33  jdecay(komega) = 1
34  jdecay(kpion) = 1
35 !////////////
36  call init
37  do j = 1, nevent
38  if(inpfileno .gt. 0) then
39  call readinpfile(eof)
40  if(eof .eq. 1) then
41  write(0,*)
42  * ' number of events generated is ',j-1
43  goto 100
44  endif
45  call formpjtg(0)
46  endif
47  call gencol(w, ntp)
48 !//////////
49  call cmydecay(jdecay, tau, w, ntp, ntpout)
50  ntp = ntpout
51 !////////////
52  call cutbyangle(w, ntp, ntp)
53  call sortbyke(w, ntp) ! sort by kinetic energy
54  if(trace .gt. 0) then
55  call outtrace(j, w, ntp)
56  endif
57  call outresul(w, ntp)
58  enddo
59  write(0,*)
60  * ' number of events generated is ',nevent
61  100 continue
62  end
63 
64  subroutine init
65  implicit none
66 #include "Zptcl.h"
67 #include "Zcode.h"
68 #include "Zevhnp.h"
69 #include "Zevhnv.h"
70 #include "Zmass.h"
71 #include "Zmanager.h"
72 #include "Zmanagerp.h"
73 #include "Ztrackp.h"
74 
75  include "Zprivate.h"
76  character*200 input, file
77  character*20 uid
78  integer klena, icon, eof
79 
80  external cblkmanager
81  external cblkevhnp
82 
83  call creadparam(5)
84  if(tracedir .eq. ' ') then
85  call cgetloginn(uid)
86  tracedir = '/tmp/'//uid(1:klena(uid))
87  endif
88 
89  if(desteventno(2) .eq. 0) then
90  nevent =abs( desteventno(1) )
91  else
92  nevent = abs( desteventno(2) )
93  endif
94  call cmkseed(initrn(1), initrn)
95  call rnd1i(initrn) ! random number init.
96  call cquhookr(1, wzmin)
97  call cquhookr(2, wzmax)
98  call cquhookr(3, trackl)
99  call cquhooki(1, outzero)
100 !
101 ! make projectile going +z
102  call cquhookc(1, input)
103 
104  if(input(1:4) .eq. "file") then
105  read(input(5:10), *) inpfileno
106  xyz=index(input, "xyz")
107  call cquhookc(2, input)
108 
109  file = ' '
110  file=input(1:klena(input))
111 
112  call copenfw2(inpfileno, file, 1, icon)
113  if(icon .ne. 1) then
114  write(0,*)
115  * ' file=', file, ' cannot be opened in Gencol'
116  stop 9999
117  endif
118  call readinpfile(eof)
119 ! once rewind to read successively for event generation
120  rewind inpfileno
121  else
122  inpfileno=0
123  read(input, *)
124  * pjcode, pjsub, pjchg, pjpx, pjpy, pjpz
125  call cquhookc(2, input)
126  read(input, *)
127  * tgcode, tgsub, tgchg, tgpx, tgpy, tgpz
128  call cquhookc(3, input)
129  if(input .ne. ' ') then
130  read(input, *) xpos, ypos, zpos
131  xyz = 1
132  else
133  xyz = 0
134  endif
135  endif
136 
137  call formpjtg(1)
138 
139  call cfixprefix('configDummy')
140  call csetcosorepi('epics')
141  if( index( intmodel,'qgsjet1') .ne. 0 ) then
142 #ifdef QGSJET1
143  call qgs01init
144  activemdl = 'qgsjet1'
145 #else
146  write(0,*) 'to use qgsjet1, define it in Zintmodel%h'
147 #endif
148  elseif(index(intmodel, 'sibyll') .ne. 0 ) then
149 #ifdef SIBYLL
150  call sibyllinit
151  activemdl = 'sibyll'
152 #else
153  write(0,*) 'to use sibyll, define it in Zintmodel%h'
154 #endif
155  else
156  call cintmodels('epics')
157  call cfixmodel( plab )
158  endif
159 
160  write(0, *) 'Active int. model=',activemdl
161  write(0, *) ' equiv. lab E=', plab%fm%p(4)
162  if(xyz .eq. 0) then
163 ! write(*, '(a)') '# mulsubKEdir user '
164  else
165 ! write(*, '(a)') '# mulsubKExyzdir user '
166  endif
167 ! write(*, '(a)') '#--------------------------------'
168  end
169  subroutine readinpfile(eof)
170  implicit none
171 #include "Zptcl.h"
172  include "Zprivate.h"
173 
174  integer eof ! output . data read--> 0
175  ! eof reached --> 1
176  read(inpfileno,*, end=100)
177  * pjcode, pjsub, pjchg, pjpx, pjpy, pjpz
178  read(inpfileno,*, end=100)
179  * tgcode, tgsub, tgchg, tgpx, tgpy, tgpz
180  if(xyz .gt. 0 ) then
181  read(inpfileno,*, end=100) xpos, ypos, zpos
182  endif
183  eof = 0
184  return
185  100 continue
186  eof = 1
187  end
188 ! *******************
189  subroutine formpjtg(confirm)
190 ! ******************
191  implicit none
192 #include "Zptcl.h"
193 #include "Zcode.h"
194 #include "Zevhnp.h"
195 #include "Zevhnv.h"
196 #include "Zmass.h"
197 #include "Zmanager.h"
198 #include "Zmanagerp.h"
199 #include "Ztrackp.h"
200 
201  include "Zprivate.h"
202 
203  integer confirm ! input. if 0, root s is not printed.
204  ! else printed
205  real*8 roots, s
206 ! form projectile and target
207 
208  call cmkptc(pjcode, pjsub, pjchg, pj)
209 
210  if(pjcode .ne. kgnuc) then
211  pjmnum = 1
212  else
213  pjmnum = pjsub
214  endif
215  pj%fm%p(1) = pjpx*pjmnum
216  pj%fm%p(2) = pjpy*pjmnum
217  pj%fm%p(3) = pjpz*pjmnum
218  pj%fm%p(4) =
219  * sqrt(pj%fm%p(1)**2 + pj%fm%p(2)**2 + pj%fm%p(3)**2
220  * + pj%mass**2)
221 
222 ! make taget (rest of moving -z)
223  call cmkptc(tgcode, tgsub, tgchg, tg)
224  if(tgcode .ne. kgnuc) then
225  tgmnum = 1
226  else
227  tgmnum = tgsub
228  endif
229  tg%fm%p(1) = tgpx*tgmnum
230  tg%fm%p(2) = tgpy*tgmnum
231  tg%fm%p(3) = tgpz*tgmnum
232  tg%fm%p(4) =
233  * sqrt(tg%fm%p(1)**2 + tg%fm%p(2)**2 + tg%fm%p(3)**2
234  * + tg%mass**2)
235 !
236  if(tgpx .eq. 0. .and. tgpy .eq. 0. .and.
237  * tgpz .eq. 0.) then
238 ! target is at rest;
239  s= 2*pj%fm%p(4)*tg%mass +tg%mass**2 + pj%mass**2
240  else
241 ! by general formula
242  s = (pj%fm%p(4)+tg%fm%p(4))**2 -
243  * (pjpx+tgpx)**2 + (pjpy+tgpy)**2+(pjpz+tgpz)**2
244  endif
245  roots = sqrt(s)
246  if(confirm .ne. 0) then
247  write(0, *) ' roots/2=', roots/2
248  endif
249 !c boost to target rest system
250  call cbst1(1, tg, pj, plab)
251  end
252 ! ************************
253  subroutine outresul(a, ntp)
254  implicit none
255 #include "Zptcl.h"
256 #include "Zcode.h"
257 #include "Zevhnp.h"
258 #include "Zevhnv.h"
259 #include "Zmass.h"
260 #include "Zmanagerp.h"
261 #include "Zmanager.h"
262 #include "Ztrackp.h"
263  include "Zprivate.h"
264 ! general process information; only for dpmjet3
265  INTEGER IPROCE,IDNODF,IDIFR1,IDIFR2,IDDPOM,IPRON
266  COMMON /poprcs/ iproce,idnodf,idifr1,idifr2,iddpom,ipron(15,4)
267 ! IPROCE
268 ! 1 non-diffractive inelastic
269 ! 2 elestic
270 ! 3 quasi elestic vector meson prod. (photon)
271 ! 4 central diffraction
272 ! 5 single diff. ptcl 1
273 ! 6 // ptcl 2
274 ! 7 double diff.
275 ! 8 direct photo-hadron
276 ! For moore detail, see manual in Documents/CosmicRays/phojetShort.pdf
277 ! say, IDIFR1 classifies IPROCE=5
278 
279  integer ntp
280  type(ptcl):: a(ntp)
281  integer i, j
282 
283  real*8 p, wx, wy, wz, pt, y, eta
284  integer npzm, npzp, Ncht, Nchpzm, Nchpzp
285  integer nw, diffcode(20)
286 !//////////////
287  integer icon
288  type(ptcl):: pcms, pl
289  type(fmom):: gc
290 !/////////////
291  npzm = 0
292  npzp = 0
293  nchpzm = 0
294  nchpzp = 0
295 
296  do j = 1, ntp
297  i = indx(j)
298 ! pt = sqrt( a(i).fm.p(1)**2 + a(i).fm.p(2)**2)
299 ! p= sqrt( pt**2 + a(i).fm.p(3)**2 )
300  if( a(i)%fm%p(3) .gt. 0.) then
301  npzp = npzp + 1
302  if(a(i)%charge .ne. 0) then
303  nchpzp = nchpzp +1
304  endif
305  else
306  npzm = npzm + 1
307  if(a(i)%charge .ne. 0) then
308  nchpzm = nchpzm + 1
309  endif
310  endif
311  enddo
312  ncht = nchpzm + nchpzp
313  call getdiffcode(nw, diffcode)
314  write(*,'("h ",i3, 6i6)' )
315  * diffcode(1), ntp, npzm, npzp, ncht, nchpzm, nchpzp
316 !////////
317 ! call cgetcm(tg, plab, gc, icon)
318 ! write(0,*) ' gcbeta=', gc
319 !//////////////////
320 
321  do j = 1, ntp
322  i = indx(j)
323  pt = sqrt( a(i)%fm%p(1)**2 + a(i)%fm%p(2)**2)
324 ! p= sqrt( pt**2 + a(i).fm.p(3)**2 )
325 !
326 ! wx = a(i).fm.p(1)/p
327 ! wy = a(i).fm.p(2)/p
328 ! wz = a(i).fm.p(3)/p
329 ! if(xyz .eq. 0) then
330 ! write(*,'(3i3,g14.5,3f17.13, i5)')
331 ! * a(i).code, a(i).subcode, a(i).charge,
332 ! * a(i).fm.p(4)-a(i).mass, wx, wy, wz, j
333 ! else
334 ! write(*,'(3i3,g14.5,1p3E11.3,0p3f17.13, i5)')
335 ! * a(i).code, a(i).subcode, a(i).charge,
336 ! * a(i).fm.p(4)-a(i).mass, xpos, ypos, zpos,
337 ! * wx, wy, wz, j
338 !
339 ! endif
340  call cyeta(a(i), y, eta)
341 !/////////////
342 ! call clorez(gc, a(i), pl)
343 ! write(*,'(2i4, 1p, g13.4, 0p, f10.4, 1p,6g13.4)')
344 ! * a(i).code, a(i).charge, pt, eta, pl.fm.p(4)-pl.mass,
345 ! * a(i).fm.p(1:4)
346 !//////////////////
347  write(*,'(2i4, 1p, g13.4, 0p, f10.4,1p, g13.4)')
348  * a(i)%code, a(i)%charge, pt, eta, a(i)%fm%p(4)
349  enddo
350  if(ntp .gt. 0 .or. outzero .eq. 0) then
351  write(*, *)
352  endif
353  end
354  subroutine gencol(a, ntp)
355  implicit none
356 #include "Zptcl.h"
357 #include "Zcode.h"
358 #include "Zevhnv.h"
359 #include "Zevhnp.h"
360 #include "Zmanagerp.h"
361  include "Zprivate.h"
362  type(ptcl):: a(*)
363 ! projectile and target information (both befor
364 ! and after collision ) in different system.
365 !
366  integer ntp
367  integer j
368  integer tZ, tA
369  real*8 xs
370 !
371  if( tg%code .eq. knuc ) then
372  ta = 1
373  elseif( tg%code .eq. kgnuc ) then
374  ta = tg%subcode
375  else
376  write(0,*) ' target code=', tg%code, 'invalid'
377  stop 9999
378  endif
379  tz = tg%charge
380 
381 #ifdef QGSJET1
382  if(activemdl .eq. 'qgsjet1') then
383  call qgs01event(plab, ta, tz, a, ntp)
384  endif
385 #elif defined SIBYLL
386  if(activemdl .eq. 'sibyll') then
387  call sibyllevent(plab, ta, tz, a, ntp)
388  endif
389 #else
390  if(activemdl .eq. 'qgsjet2' ) then
391  call cxsecqgs(plab, ta, xs)
392  endif
393  call chacol(plab, ta, tz, a, ntp)
394 #endif
395  do j = 1, ntp
396 ! boost to target mooving system
397  call cibst1(j, tg, a(j), a(j))
398  enddo
399  end
400 
401 
402  subroutine cutbyangle(a, ntp0, ntp)
403  implicit none
404 #include "Zptcl.h"
405 #include "Zcode.h"
406 #include "Zevhnv.h"
407 #include "Zevhnp.h"
408 #include "Zmanagerp.h"
409  include "Zprivate.h"
410  type(ptcl):: a(*)
411  integer ntp0 ! input. number of ptcls. in a
412  integer ntp ! output. could be the same as ntp0
413  integer j
414  integer i
415  real*8 p, wz
416  j = 0
417  do i = 1, ntp0
418  p = a(i)%fm%p(1)**2 + a(i)%fm%p(2)**2 +
419  * a(i)%fm%p(3)**2
420  p = sqrt(p)
421  wz = a(i)%fm%p(3)/p
422  if( wz .ge. wzmin .and. wz .le. wzmax ) then
423  j = j + 1
424  a(j)=a(i)
425  endif
426  enddo
427  ntp = j
428  end
429  subroutine sortbyke(a, ntp)
430  implicit none
431 #include "Zptcl.h"
432 #include "Zcode.h"
433 
434  include "Zprivate.h"
435  integer ntp
436  type(ptcl):: a(*)
437 ! projectile and target information (both befor
438 ! and after collision ) in different system.
439 !
440 
441  integer i
442  do i = 1, ntp
443  ke(i) = a(i)%fm%p(4) - a(i)%mass
444  enddo
445  call kqsortd(ke, indx, ntp)
446  call ksortinv(indx, ntp)
447 ! ke( indx(1) ) is the highest energy
448  end
449  subroutine outtrace(nev, a, ntp)
450  implicit none
451 #include "Zptcl.h"
452 #include "Zcode.h"
453 #include "Zevhnp.h"
454 #include "Zevhnv.h"
455 #include "Zmass.h"
456 #include "Zmanagerp.h"
457 #include "Zmanager.h"
458 #include "Ztrackp.h"
459  include "Zprivate.h"
460 
461  integer ntp, nev
462  type(ptcl):: a(ntp)
463  integer i, j, leng, icon, klena
464  real*8 p, wx, wy, wz
465  real x1, y1, z1, x2, y2, z2
466  character*100 tracefile
467 
468  write(tracefile, *) tracedir(1:klena(tracedir))//'/trace', nev
469  call kseblk(tracefile, ' ', leng)
470  call copenfw(tracedev, tracefile(1:leng), icon)
471  if(icon .ne. 0) then
472  call cerrormsg('tracefile could not be opened',0)
473  endif
474  do j = 1, ntp
475  i = indx(j)
476  p= sqrt( a(i)%fm%p(1)**2 + a(i)%fm%p(2)**2
477  * + a(i)%fm%p(3)**2 )
478  wx = a(i)%fm%p(1)/p
479  wy = a(i)%fm%p(2)/p
480  wz = a(i)%fm%p(3)/p
481  if(xyz .eq. 0) then
482  x1 = 0.
483  y1 = 0.
484  z1 = 0.
485  else
486  x1 = xpos
487  y1 = ypos
488  z1 = zpos
489  endif
490  x2 = x1 + wx*trackl
491  y2 = y1 + wy*trackl
492  z2 = z1 + wz*trackl
493  write(tracedev,'(3g14.5, i3, g14.4, i3, i2)')
494  * x1, y1, z1,
495  * a(i)%code, a(i)%fm%p(4) - a(i)%mass, a(i)%charge,
496  * 0
497  write(tracedev, '(3g14.5, i3, g14.4, i3, g14.4)' )
498  * x2, y2, z2,
499  * a(i)%code, a(i)%fm%p(4) - a(i)%mass, a(i)%charge,
500  * trackl
501  write(tracedev, *)
502  write(tracedev, *)
503  enddo
504  close(tracedev)
505  end
506 
subroutine cmkseed(dummy, seed)
Definition: cmkSeed.f:3
subroutine cerrormsg(msg, needrtn)
Definition: cerrorMsg.f:4
subroutine csetcosorepi(from)
Definition: cintModels.f:252
subroutine sortbyke(a, ntp)
Definition: Gencol.f:754
integer maxn LabEquivE real outresul integer pjchg integer tgchg integer xyz real pjpz real tgpy
Definition: Zprivate.h:8
max ptcl codes in the kgzai
Definition: Zcode.h:2
subroutine cbst1(init, p1, p2, po)
Definition: cbst1.f:54
dE dx *! Nuc Int sampling table e
Definition: cblkMuInt.h:130
integer maxn LabEquivE real outresul integer pjchg integer tgchg integer inpfileno
Definition: Zprivate.h:8
max ptcl codes in the kdmes
Definition: Zcode.h:2
max ptcl codes in the kgnuc
Definition: Zcode.h:2
subroutine cintmodels(from)
Definition: cintModels.f:3
subroutine cfixmodel(aPtcl)
Definition: cfixModel.f:2
max ptcl codes in the kphi
Definition: Zcode.h:2
latitude latitude this system is used *****************************************************************! type coord sequence union map real z z in m endmap xyz map real ! latitude in deg is to the north ! longitude in deg is to the east *h ! height in m endmap llh map real ! polar angle ! azimuthal angle *radius ! radial distance endmap sph endunion character *sys ! which system xyz
Definition: Zcoord.h:25
max ptcl codes in the klambdac
Definition: Zcode.h:2
max ptcl codes in the kkaon
Definition: Zcode.h:2
subroutine cibst1(init, p1, p2, po)
Definition: cibst1.f:29
subroutine readinpfile(eof)
Definition: Gencol.f:422
integer maxn LabEquivE real outresul integer pjcode
Definition: Zprivate.h:8
subroutine outtrace(nev, a, ntp)
Definition: Gencol.f:772
subroutine cmydecay(Jdecay, tau, a, nin, nout)
Definition: cmydecay.f:4
subroutine qgs01init
Definition: qgs01init.f:2
max ptcl codes in the komega
Definition: Zcode.h:2
Definition: Zptcl.h:72
integer maxn wzmax
Definition: Zprivate.h:3
********************block data cblkHeavy ********************integer j data *HeavyG2symbol *data *HeavyG2code kiron data j
Definition: cblkHeavy.h:36
subroutine copenfw2(io, fnin, form, icon)
Definition: copenf.f:205
subroutine getdiffcode(nw, difcode)
Definition: getDiffCode.f:3
subroutine kqsortd(A, ORD, N)
Definition: kqsortd.f:23
subroutine cquhookr(i, rv)
Definition: cqUHookr.f:2
integer maxn LabEquivE real outresul integer pjchg integer tgcode
Definition: Zprivate.h:8
subroutine init
Definition: Gencol.f:62
subroutine ksortinv(idx, n)
Definition: ksortinv.f:2
subroutine creadparam(io)
Definition: creadParam.f:5
max ptcl codes in the klambda
Definition: Zcode.h:2
integer maxn LabEquivE real outresul integer pjchg integer tgsub
Definition: Zprivate.h:8
subroutine gencol(a, ntp)
Definition: Gencol.f:585
integer maxn LabEquivE real outresul integer pjchg integer tgchg integer xyz real pjpz real * tgpx
Definition: Zprivate.h:8
integer maxn LabEquivE real * ke(maxn) integer indx(maxn) integer nevent integer outzero
integer maxn trackl
Definition: Zprivate.h:3
********************block data cblkHeavy ********************integer j data *HeavyG2symbol *data *HeavyG2code knuc
Definition: cblkHeavy.h:7
max ptcl codes in the krho
Definition: Zcode.h:2
subroutine cfixprefix(dsn)
Definition: cintModels.f:209
integer maxn LabEquivE real outresul integer pjchg integer tgchg integer xyz real pjpy
Definition: Zprivate.h:8
nodes a
subroutine formpjtg(confirm)
Definition: Gencol.f:443
subroutine copenfw(io, fnin, icon)
Definition: copenf.f:122
subroutine cutbyangle(a, ntp0, ntp)
Definition: Gencol.f:727
subroutine cquhooki(i, iv)
Definition: cqUHookr.f:15
subroutine cmkptc(code, subcode, charge, p)
Definition: cmkptc.f:15
max ptcl codes in the keta
Definition: Zcode.h:2
subroutine outresul(a, ntp)
Definition: Gencol0.f:229
integer maxn LabEquivE real outresul integer pjsub
Definition: Zprivate.h:8
subroutine cquhookc(i, cv)
Definition: cqUHookr.f:28
max ptcl codes in the klast
Definition: Zcode.h:2
Definition: Zptcl.h:75
max ptcl codes in the kpion
Definition: Zcode.h:2
integer maxn LabEquivE real outresul integer pjchg integer tgchg integer xyz real * pjpx
Definition: Zprivate.h:8
max ptcl codes in the ksigma
Definition: Zcode.h:2
subroutine chacol(pj, ia, iz, a, ntp)
Definition: chAcol.f:3
subroutine cyeta(p, y, eta)
Definition: cyeta.f:5
subroutine kseblk(text, c, lc)
Definition: kseblk.f:18
max ptcl codes in the kbomega
Definition: Zcode.h:2
subroutine cgetloginn(userid)
Definition: cgetLoginN.f:29