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

Go to the source code of this file.

Functions/Subroutines

subroutine chncol (pj, tg, a, ntp, icon)
 
subroutine chncol2 (pj, tg, a, ntp, icon)
 
subroutine cdcycp (a, nin, n)
 
subroutine cpikcd (a, ntp)
 
subroutine cnnbdc (a, b, n)
 
subroutine cddbdc (a, b, n)
 
subroutine cspipm (pj, x, a)
 
subroutine cskchg (pj, x, a)
 

Function/Subroutine Documentation

◆ cdcycp()

subroutine cdcycp ( type(ptcl), dimension(*)  a,
integer  nin,
integer  n 
)

Definition at line 158 of file chncol.f.

References cddbdc(), cnnbdc(), kddb, and knnb.

Referenced by chncol2().

158 ! ****************************************************
159  implicit none
160 !---- include '../../Zptcl.h'
161 #include "Zptcl.h"
162 !---- include '../../Zcode.h'
163 #include "Zcode.h"
164  integer nin, n
165  type(ptcl):: a(*)
166  integer i, k, nx
167  type(ptcl):: p(2)
168  n=nin
169  do i=1, nin
170  k=a(i)%code
171  if(k .eq. knnb .or. k .eq. kddb) then
172  if(k .eq. knnb) then
173  call cnnbdc(a(i), p, nx)
174  elseif(k .eq. kddb) then
175  call cddbdc(a(i), p, nx)
176  endif
177 ! put n or D in the i-th pos. and append c antiptcl
178 ! at the bottome and increase n
179  a(i) = p(1)
180  a(n +1 ) = p(2)
181  n=n+1
182  endif
183  enddo
subroutine cddbdc(a, b, n)
Definition: chncol.f:259
nodes i
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
subroutine cnnbdc(a, b, n)
Definition: chncol.f:231
max ptcl codes in the knnb
Definition: Zcode.h:2
real(4), save a
Definition: cNRLAtmos.f:20
Definition: Zptcl.h:75
max ptcl codes in the kddb
Definition: Zcode.h:2
integer n
Definition: Zcinippxc.h:1
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cddbdc()

subroutine cddbdc ( type(ptcl a,
type(ptcl), dimension(2)  b,
integer  n 
)

Definition at line 259 of file chncol.f.

References c2bdcy(), cmkptc(), kd0, kd0b, kdmes, and rndc().

Referenced by cdcycp(), and cinteddb().

259  implicit none
260 !---- include '../../Zptcl.h'
261 #include "Zptcl.h"
262 !---- include '../../Zcode.h'
263 #include "Zcode.h"
264  integer n
265  type(ptcl):: a, b(2)
266 !
267  integer ic
268  real*8 u
269 !
270  call rndc(u)
271  ic=int(u*3)-1 ! charge
272  if(ic .eq. 0) then
273  call cmkptc(kdmes, kd0, ic, b(1))
274  call cmkptc(kdmes, kd0b, -ic, b(2))
275  else
276  call cmkptc(kdmes, 0, ic, b(1))
277  call cmkptc(kdmes, 0, -ic, b(2))
278  endif
279  call c2bdcy(a, b(1), b(2))
280  n=2
max ptcl codes in the kdmes
Definition: Zcode.h:2
subroutine rndc(u)
Definition: rnd.f:91
subroutine c2bdcy(p, p1, p2)
Definition: c2bdcy.f:44
max ptcl codes in the kseethru ! subcode integer kd0
Definition: Zcode.h:2
real(4), save a
Definition: cNRLAtmos.f:20
subroutine cmkptc(code, subcode, charge, p)
Definition: cmkptc.f:15
real(4), save b
Definition: cNRLAtmos.f:21
max ptcl codes in the kseethru ! subcode integer kd0b
Definition: Zcode.h:2
Definition: Zptcl.h:75
integer n
Definition: Zcinippxc.h:1
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chncol()

subroutine chncol ( type(ptcl pj,
type(ptcl tg,
type(ptcl), dimension(*)  a,
integer  ntp,
integer  icon 
)

Definition at line 11 of file chncol.f.

References chncol2(), rndc(), and true.

Referenced by chacoladhoc().

11 ! pj: type ptcl. Input. to give the incident
12 ! projectile hadron in the Lab. system.
13 ! tg: type ptcl. Input. to give the target
14 ! nucleon in the Lab. system.
15 ! a: array of type ptcl. Output. to get
16 ! produced particles.
17 ! ntp: integer. Output. to get the total number of produced
18 ! particles in a.
19 ! icon: Integer. Output. if 0, particle generation is ok
20 ! if non 0, you may give up generation.
21 !
22 ! Note: target recoil is put in a(ntp-1),
23 ! projectile recoil is put in a(ntp)
24 !
25  implicit none
26 
27 #include "Zptcl.h"
28 #include "Zevhnv.h"
29 
30 !
31  type(ptcl):: pj, tg, a(*)
32  integer ntp, icon
33  integer i
34  real*8 u, xmax
35 
36 ! real*8 BigXRejCnst/1./, BigXRejPw/2./
37 ! default
38  real*8 bigxrejcnst/.4/, bigxrejpw/2.2/
39 ! enhance
40 ! real*8 BigXRejCnst/.8/, BigXRejPw/2.2/
41 ! denhance
42 ! real*8 BigXRejCnst/.01/, BigXRejPw/1./
43 ! real*8 BigXRejCnst/0.0003/, BigXRejPw/3.5/
44 
45  do while (.true.)
46  call chncol2(pj, tg, a, ntp, icon)
47  if(icon .ne. 0) goto 100
48 ! find max secondary energy.
49  xmax = 0.
50  do i = 1, ntp-2
51  if(xmax .lt. a(i)%fm%p(4)) xmax = a(i)%fm%p(4)
52  enddo
53  xmax = xmax/pj%fm%p(4)
54 ! if very large x appears, discard it
55 ! with some probabilty to adjust
56 ! the x-distribution.
57 !
58 !c if(xmax .gt. 0.4) then
59 !c call rndc(u)
60 !c if(u .lt. (1. - (xmax-0.40)/0.6 )**2 ) then
61 !c goto 100
62 !c endif
63 !c else
64 !c goto 100
65 !c endif
66  call rndc(u)
67  if(u .lt.
68  * ( bigxrejcnst/(bigxrejcnst +
69  * (xmax/(1.0-xmax))**bigxrejpw) ) ) goto 100
70 
71  enddo
72  100 continue
subroutine chncol2(pj, tg, a, ntp, icon)
Definition: chncol.f:84
nodes i
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), save a
Definition: cNRLAtmos.f:20
Definition: Zptcl.h:75
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chncol2()

subroutine chncol2 ( type(ptcl pj,
type(ptcl tg,
type(ptcl), dimension(*)  a,
integer  ntp,
integer  icon 
)

Definition at line 84 of file chncol.f.

References cconsvchg(), cdcycp(), cgnlp(), cibst1(), cpikcd(), crot3mom(), cs2lp(), and true.

Referenced by chncol().

84 ! pj: type ptcl. Input. to give the incident
85 ! projectile hadron in the Lab. system.
86 ! tg: type ptcl. Input. to give the target
87 ! nucleon in the Lab. system.
88 ! a: array of type ptcl. Output. to get
89 ! produced particles.
90 ! ntp: integer. Output. to get the total number of produced
91 ! particles in a.
92 ! icon: Integer. Output. if 0, particle generation is ok
93 ! if non 0, you may give up generation.
94 !
95 ! Note: target recoil is put in a(ntp-1),
96 ! projectile recoil is put in a(ntp)
97 !
98  implicit none
99 
100 #include "Zptcl.h"
101 #include "Zevhnv.h"
102 !
103 !
104  type(ptcl):: pj, tg, a(*)
105  integer ntp, icon, jcon, nfin
106  integer i, outc
107  type(ptcl):: pjin
108 ! use pt=0 incident
109  pjin = pj
110  pjin%fm%p(1) = 0.
111  pjin%fm%p(2) = 0.
112  pjin%fm%p(3) = sqrt(pjin%fm%p(4)**2 - pjin%mass**2)
113 !
114 
115 ! *** until loop**until succeed
116  do while (.true.)
117 ! generate 2 leading ptcls Zevhnv become ready for use.
118  call cs2lp(pjin, tg, icon)
119  if(icon .ne. 0) then
120 ! neglect this event; because of very low energy
121 ! no particle can be produced
122  ntp = 0
123  jcon = 0
124  else
125 ! generation of particles other than the leadings.
126  call cgnlp(a, ntp, jcon)
127  endif
128 
129  if(jcon .eq. 0 ) goto 10
130  enddo
131  10 continue
132 ! ptcls in 'a' should have cms energy, here.
133 ! give final ptcl charge/ subcode for
134 ! pi+-, k0,k0~,k+,k-
135  call cpikcd(a, ntp)
136  outc =rtglab%charge+rpjlab%charge - pj%charge- tg%charge
137  call cconsvchg(outc, a, ntp, icon)
138  if(icon .eq. 0) then
139 ! boost to lab system
140  do i=1, ntp
141  call cibst1(i, cmsp, a(i), a(i))
142  enddo
143 ! decay of "composit ptcls" (nn~, dd~)
144  call cdcycp(a, ntp, nfin)
145  ntp=nfin
146 ! store recoils in a
147  a(ntp+1) = rtglab
148  a(ntp+2) = rpjlab
149  ntp = ntp +2
150 ! -------------- rotate
151  call crot3mom(pj, a, ntp)
152  endif
153 
subroutine cgnlp(a, ntp, icon)
Definition: cgnlp.f:4
nodes i
subroutine cpikcd(a, ntp)
Definition: chncol.f:188
subroutine cibst1(init, p1, p2, po)
Definition: cibst1.f:29
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 cdcycp(a, nin, n)
Definition: chncol.f:158
subroutine crot3mom(p, a, n)
Definition: crot3mom.f:2
subroutine cconsvchg(outc, a, ntp, icon)
Definition: cconsvChg.f:2
subroutine cs2lp(proj, trgt, icon)
Definition: cs2lp.f:5
real(4), save a
Definition: cNRLAtmos.f:20
Definition: Zptcl.h:75
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cnnbdc()

subroutine cnnbdc ( type(ptcl a,
type(ptcl), dimension(2)  b,
integer  n 
)

Definition at line 231 of file chncol.f.

References c2bdcy(), cmkptc(), kneutron, kneutronb, knuc, and rndc().

Referenced by cdcycp(), and cintennb().

231  implicit none
232 !---- include '../../Zptcl.h'
233 #include "Zptcl.h"
234 !---- include '../../Zcode.h'
235 #include "Zcode.h"
236  integer n
237  type(ptcl):: a, b(2)
238 !
239  integer ic
240  real*8 u
241 !
242  call rndc(u)
243  ic=u*2 ! charge
244  if(ic .eq. 0) then
245  call cmkptc(knuc, kneutron, ic, b(1))
246  call cmkptc(knuc, kneutronb, -ic, b(2))
247  else
248  call cmkptc(knuc, 0, ic, b(1))
249  call cmkptc(knuc, 0, -ic, b(2))
250  endif
251  call c2bdcy(a, b(1), b(2))
252  n=2
subroutine rndc(u)
Definition: rnd.f:91
max ptcl codes in the kseethru ! subcode integer kneutronb
Definition: Zcode.h:2
max ptcl codes in the kseethru ! subcode integer kneutron
Definition: Zcode.h:2
********************block data cblkHeavy ********************integer j data *HeavyG2symbol *data *HeavyG2code knuc
Definition: cblkHeavy.h:7
subroutine c2bdcy(p, p1, p2)
Definition: c2bdcy.f:44
real(4), save a
Definition: cNRLAtmos.f:20
subroutine cmkptc(code, subcode, charge, p)
Definition: cmkptc.f:15
real(4), save b
Definition: cNRLAtmos.f:21
Definition: Zptcl.h:75
integer n
Definition: Zcinippxc.h:1
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cpikcd()

subroutine cpikcd ( type(ptcl), dimension(ntp)  a,
integer  ntp 
)

Definition at line 188 of file chncol.f.

References cskchg(), cspipm(), kkaon, and kpion.

Referenced by chncol2().

188 !
189 !*******************************************************
190  implicit none
191 !---- include '../../Zptcl.h'
192 #include "Zptcl.h"
193 !---- include '../../Zcode.h'
194 #include "Zcode.h"
195 !---- include '../Zevhnv.h'
196 #include "Zevhnv.h"
197  integer ntp
198  type(ptcl):: a(ntp)
199 !
200  integer i, k
201  real*8 x
202  if(ntp .eq. 1) then
203 ! nothing to do. charge assignment should have
204 ! been done in c1pion
205  else
206  do i=1, ntp
207  x = a(i)%fm%p(3)/pjcms%fm%p(3)
208  k = a(i)%code
209  if(k .eq. kpion .and. a(i)%charge .ne. 0) then
210  if(x .gt. 0.) then
211  call cspipm(pjcms, x, a(i))
212  else
213  call cspipm(tgcms, -x, a(i))
214  endif
215  elseif(k .eq. kkaon) then
216 ! set kaon charge
217  if(x .gt. 0.) then
218  call cskchg(pjcms, x, a(i))
219  else
220  call cskchg(tgcms, -x, a(i))
221  endif
222  endif
223  enddo
224  endif
nodes i
subroutine cspipm(pj, x, a)
Definition: chncol.f:292
max ptcl codes in the kkaon
Definition: Zcode.h:2
subroutine cskchg(pj, x, a)
Definition: chncol.f:333
real(4), save a
Definition: cNRLAtmos.f:20
Definition: Zptcl.h:75
max ptcl codes in the kpion
Definition: Zcode.h:2
! 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:

◆ cskchg()

subroutine cskchg ( type(ptcl pj,
real*8  x,
type(ptcl a 
)

Definition at line 333 of file chncol.f.

References cmkptc(), d0, k0l, k0s, kkaon, regptcl, and rndc().

Referenced by cpikcd().

333 ! ************************************************
334  implicit none
335 
336 #include "Zptcl.h"
337 #include "Zcode.h"
338 
339  type(ptcl):: pj, a
340  real*8 x
341 !
342  real*8 f, u, up
343  if(a%charge .ne. 0) then
344 ! k+/k- for p incident
345  if(x .lt. .3d0) then
346  f=exp(4.36d0*x)
347  elseif(x .lt. .6d0) then
348  f=3.7d0*exp(6.5d0*(x-.3d0))
349  else
350  f=27.d0*exp(11.3d0*(x-.6d0))
351  endif
352  if(pj%charge .eq. 1) then
353  up=f/(1.d0+f)
354  elseif(pj%charge .eq. -1) then
355  up=1./(1.d0+f)
356  else
357  if(pj%subcode .eq. 0) then
358  up=0.5d0
359  elseif(pj%subcode .eq. regptcl) then
360  up=1.d0/(1.d0+f)
361  else
362  up=f/(1.d0+f)
363  endif
364  endif
365  call rndc(u)
366  if(u .lt. up) then
367  call cmkptc(kkaon, 0, 1, a)
368  else
369  call cmkptc(kkaon, 0, -1, a)
370  endif
371  else
372 ! k0
373  call rndc(u)
374  if(u .lt. .50d0) then
375  call rndc(u)
376  if(u .lt. 0.5) then
377 ! k0 short
378  call cmkptc(kkaon, k0s, 0, a)
379  else
380  call cmkptc(kkaon, -k0s, 0, a)
381  endif
382  else
383  call rndc(u)
384  if(u .lt. 0.5) then
385  call cmkptc(kkaon, k0l, 0, a)
386  else
387  call cmkptc(kkaon, -k0l, 0, a)
388  endif
389  endif
390  endif
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data h g *is for param c g data up(2, 1)/7.0d0/
max ptcl codes in the kseethru ! subcode integer k0l
Definition: Zcode.h:2
max ptcl codes in the kseethru ! subcode integer k0s
Definition: Zcode.h:2
max ptcl codes in the kkaon
Definition: Zcode.h:2
subroutine rndc(u)
Definition: rnd.f:91
max ptcl codes in the kseethru ! subcode integer regptcl
Definition: Zcode.h:2
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
real(4), save a
Definition: cNRLAtmos.f:20
subroutine cmkptc(code, subcode, charge, p)
Definition: cmkptc.f:15
Definition: Zptcl.h:75
! 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
dE dx *! Nuc Int sampling table f
Definition: cblkMuInt.h:130
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cspipm()

subroutine cspipm ( type(ptcl pj,
real*8  x,
type(ptcl a 
)

Definition at line 292 of file chncol.f.

References antip, d0, regptcl, and rndc().

Referenced by cpikcd().

292 
293  implicit none
294 
295 #include "Zptcl.h"
296 #include "Zcode.h"
297 
298  type(ptcl):: pj, a
299  real*8 x
300 !
301  real*8 f, up, u
302 !
303  if(x .lt. .6d0) then
304  f = 1.d0+3.33d0*x
305  else
306  f = 3.0d0* exp( 4.4d0*(x-.6d0))
307  endif
308  if(pj%charge .eq. 1) then
309  up=f/(1.d0+f)
310  elseif(pj%charge .eq. -1) then
311  up=1.d0/(1.d0+f)
312  else
313  if(pj%subcode .eq. 0) then
314  up=0.5
315  elseif(pj%subcode .eq. regptcl) then
316  up=1.d0/(1.d0+f)
317  else
318  up=f/(1.d0+f)
319  endif
320  endif
321  call rndc(u)
322  if(u .lt. up) then
323  a%charge = 1
324  a%subcode = regptcl
325  else
326  a%charge = -1
327  a%subcode = antip
328  endif
block data include Zlatfit h c fitting region data data data data data d0 data data d0 data data h g *is for param c g data up(2, 1)/7.0d0/
subroutine rndc(u)
Definition: rnd.f:91
max ptcl codes in the kseethru ! subcode integer regptcl
Definition: Zcode.h:2
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
real(4), save a
Definition: cNRLAtmos.f:20
Definition: Zptcl.h:75
max ptcl codes in the kseethru ! subcode integer antip
Definition: Zcode.h:2
! 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
dE dx *! Nuc Int sampling table f
Definition: cblkMuInt.h:130
Here is the call graph for this function:
Here is the caller graph for this function: