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