COSMOS v7.655  COSMOSv7655
(AirShowerMC)
csoftenPiK.f
Go to the documentation of this file.
1 
2  module softenpik
3 ! This part is common to Cosmos(modifyX) and Cosmos/Util/Gencol
4 ! Used to soften the x-distribution.
5  implicit none
6 
7 ! The parameters are listed above the --------------- line
8 ! mode=0) don't do any cut/modification of X
9 ! mode=1) discard all pi/K with x>xth.
10 ! However, if the incident is pi/K the highest
11 ! energy pi/K, resp, is untouched.
12 ! mode=2) discard only pi0 with x> xth
13 ! mode=3) soften pi/K/eta spectrum;
14 ! If their x is > x>xth, this is applied.
15 ! (except for the leading pi/K when
16 ! incident is pi/K).
17 ! mode=-3 probability of softening depends on the user's algorithm
18 ! The user must modify cJudegeSplit.
19 ! trial version method
20 ! If we find a pi/K/eta with x>xth,
21 ! apply softening with the prob. of (x-xth)**a (a=0.1 default).
22 ! Softening is done by splitting the particle into two by
23 ! calling csplitMeson. The sum of two particle is made to
24 ! be the same as the original particle. However, momentum
25 ! cannot be conserved: the two particle will have the same
26 ! Pt as the original one and 3 memeta are adjuested to
27 ! satisfy P^2 + m^2 =E^2. Let the two particles' energy be
28 ! E1 and E2. If E2>E1, we split E2 into two again.
29 ! (At present, if E1> E2, we don't split E1).
30 ! (so, x>xth pi/K/eta may split into two or three).
31 ! This method results in some unnatural small bump at x<xth.
32 !
33 
34  integer,save:: mode = 3
35  real(8),save:: xth = 0.05 ! over this x, we apply cut or modification
36  real(8),save:: pw = 0.5 ! softening is determined by
37 ! if( (x-xth)**pw < u ) cycle
38 ! call csplitMeson(pstore(i), E1, E2, icon)
39 ! so if pw in (x-xth)**pw is smaller, stronger softening
40  real(8),save:: repeat = 2.5 ! apply the simple softenning repeat times
41  ! for one event if the event has pi/K/eta with X> Xth.
42  ! the odd number is probabilistic. If negative,
43  ! positive number is taken and understood as poisson
44  ! average.
45 
46  real(8),save:: e0th = 500. ! 500 GeV over this we apply cut/mod.
47  ! if make this 10 TeV, effect after shower max
48  ! decreases even for E1ry=10^19 eV.
49  ! This is in lab. frame
50  integer,save:: fwbw = 3 ! Used in Gencol; modification is applied to
51  ! 1--> forward only 2--> backword only
52  ! 3--> both ; However, 1 is used when
53  ! Cosmos output (in cms) is the target, irespectively of
54  ! fwbw.
55  logical,save:: special=.false. ! This is used olny in modifyX of
56  ! Cosmos. but not used in Gencol.
57  ! In Cosomos, it is used to see the x-dist. of
58  ! pi, K, eta at the first interaction;
59  ! For that, make this .ture. and
60  ! ***set next in parm*** (don't forget to restore
61  ! them to the default if special=.false. is reset.)
62  ! =====================================
63  ! Generate="em" (don't use "em/as")
64  ! EminObs = 8*xxxx,
65  ! =====================================
66  ! where xxxx is a value (GeV) little bit smaller then
67  ! the 1ry energy. E.g 0.99999e8 if primary is 10^17 eV.
68  !
69  ! output will look like
70  ! xd 4 1 12.3e-4
71  ! xd 4 0 2.89e-1
72  ! where xd is id, next two are code and charge
73  ! last one is X. in lab. frame.
74  ! For special=.false., i.e, to generate air showers,
75  ! you should have
76  ! =====================================
77  ! Generate="em/as"
78  ! EminObs = 500e-6,
79  ! =====================================
80 
81  logical,save:: leadingpik=.true. ! for Pi/K incident case,
82  ! treat highest energy Pi/K of the same
83  ! charge as leading Pi/K and don't
84  ! apply the method here
85  ! If used in Cosmos and
86  ! if make this f, effect increase a bit esp. at
87  ! after max.
88  ! leadingPiK=T --> Xmax ==> 35 g/cm^2
89  ! =F Xmax ==> 25 g/cm^2
90  logical,save:: usexincms=.false. !This is used only in Cosmos.
91  ! In Cosmos, generated particle's E is in Lob. so
92  ! x= E/E0 is also in Lob. If this is .false., softening by
93  ! csoftenPiK is applied to the Lab.X.
94  ! if this is .true., softening is performed after
95  ! converting the energy into cms system, and then
96  ! re-boosted to Lab.
97  !
98  real(8),save:: k1u=0.25 !
99  real(8),save:: xthl=1.d0
100  real(8),save:: xthu=10.d0
101  real(8),save:: rejpw=0.
102  logical,save:: modified ! not input. used in Cosmos
103  real(8),save:: e0lab=1.e12 ! not input. current E0 in lab.
104  ! frame. Used in Cosmos
105 
106 ! next one is about 0.6 of maxn in Epics/Util/Gencol/Zprivate.h
107 ! but this cannot inherit maxn
108  integer,parameter::half=12000
109 !--------------------------------------------------------------------
110 
111  contains
112 
113  subroutine cjudgesplit(pj, x, split)
114  implicit none
115 ! judge if this x should be split
116 #include "Zptcl.h"
117  type(ptcl):: pj ! input. projectile. you may use this
118  real(8),intent(in):: x ! KE ratio
119  logical,intent(out):: split ! .true. or .false.. if true, split meson
120 
121  real(8):: u
122 
123  if( mode == 3) then
124  call rndc(u) ! uniform rn in (0,1)
125 ! if pw in (x-xth)**pw is smaller, stronger softening
126  split = (x-xth)**pw > u
127  elseif(mode == -3 ) then
128 ! comment out next lines and supply your judgement
129 ! and fix split (.true. or .false.)
130 ! you may use xth, pw etc which are defined before "contains"
131  write(0,*) ' you have to supply your owon code here'
132  write(0,*) ' to judge if this particle should be '
133  write(0,*) ' split or not'
134  write(0,*)
135  * ' Place is cJudgeSplit; first subroutine in csoftenPiK'
136  stop
137  endif
138 
139  end subroutine cjudgesplit
140 
141  subroutine csoftenpik(inciptcl, pstore, nin, nout)
142 ! inciptcl made a collsion and generated nin particles in pstore
143 ! particle information may be at CMS or Lab.
144 ! If we want to make the x-dist. softer,
145 ! A) For Cosmos (air target), we may use the Lab system
146 ! (for large X, the distribution would be the same
147 ! as CMS)
148  implicit none
149 #include "Zcode.h"
150 #include "Zptcl.h"
151 #include "Zmanagerp.h"
152  type(ptcl):: inciptcl ! input. incident ptcl at collision
153 !
154  integer,intent(in)::nin ! # of ptcls in pstore
155  type(ptcl):: pstore(*) ! input/output. ptcls to be softened.
156  ! Softened ones are put here
157  integer,intent(out)::nout ! # of ptcls re-stored in pstore
158  logical,save::first=.true.
159  integer:: jcon, i
160  integer:: nrepeat
161  real(8):: u
162  integer:: exist
163  integer:: nnow
164 
165  if( first ) then
166  call copenf(tempdev,
167  * "$COSMOSTOP/UserHook/modifyX/softenparam.dat", jcon)
168  if( jcon /= 0 ) then
169  write(0,*) 'Data file for csoftenPiK could not be found'
170  write(0,*) 'in $COSMOSTOP/Util/Data/softenpiK.dat'
171  write(0,*) 'Forgot to set COSMOSTOP ?'
172  stop
173  else
174  call creadsoftenpara(tempdev)
175  endif
176  call cwritesoftenpara(0)
177  first = .false.
178  endif
179  if( repeat < 0. ) then
180  call kpoisn(-repeat, nrepeat)
181  else
182  call rndc(u)
183  nrepeat = repeat
184  if( u < repeat-nrepeat) then
185  nrepeat = nrepeat + 1
186  endif
187  endif
188 
189  nnow = nin
191  k1u = ck1u(e0lab)
192  xthu =cxthu(e0lab)
193  do i = 1, nrepeat
194  call csoftenpik0(inciptcl, pstore, nnow, nout, exist)
195  if( exist == 0 ) exit
196  nnow = nout
197  enddo
198  do i = 1, 1
199  call csoftenpik1(inciptcl, pstore, nnow, nout, exist)
200  if( exist == 0 ) exit
201  nnow = nout
202  enddo
203  end subroutine csoftenpik
204 
205 
206  subroutine csoftenpik0(inciptcl, pstore, nin, nout, exist)
207 ! inciptcl made a collision and generated nin particles in pstore
208 ! particle. information may be at CMS or Lab.
209 ! If we want to make the x-dist. softer,
210 ! A) For Cosmos (air target), we may use the Lab system
211 ! (for large X, the distribution would be the same
212 ! as CMS)
213  implicit none
214 #include "Zcode.h"
215 #include "Zptcl.h"
216 
217  type(ptcl):: inciptcl ! input. incident ptcl at collision
218 !
219  integer,intent(in)::nin ! # of ptcls in pstore
220  type(ptcl):: pstore(*) ! input/output. ptcls to be softened.
221  ! Softened ones are put here
222  integer,intent(out)::nout ! # of ptcls re-stored in pstore
223  integer,intent(out)::exist ! # of particcls with X> Xth
224 ! ----------------------------------------------
225  type(ptcl):: work(2)
226  logical ok
227  integer::maxi
228  real(8)::maxE
229  integer::i, j
230  real(8):: x, E0, Ex, u, u1, u2
231  real(8):: E1, E2
232  real(8)::temp
233  integer:: nc ! counting # of ptcls in pstore
234  integer:: nw, icon
235  logical:: split
236 
237  e0 = inciptcl.fm.p(4) - inciptcl.mass
238  nout = nin
239  exist = 0
240  if( e0lab < e0th) return !!!!!
241  nc = nin
242 ! if pi/K is incident, regards the highest one as leading
243 ! and avoid to modify
244  if(leadingpik .and. (inciptcl.code == kpion .or.
245  * inciptcl.code == kkaon ) ) then
246 ! find max energy ptcl index and E; eta cannot be incident
247  call cgetmaxptcl(pstore, nc, inciptcl.code, inciptcl.charge,
248  * maxi, maxe)
249  else
250  maxi = 0
251  endif
252 
253  do i = 1, nin
254  x = (pstore(i).fm.p(4)-pstore(i).mass)/e0
255  if( abs(mode) == 3 ) then
256  if( i /= maxi ) then
257 ! softening
258  if(x > xth ) then
259 ! modify high X of pi/K/eta
260  if(pstore(i).code == kpion .or.
261  * pstore(i).code == kkaon
262  * .or. pstore(i).code == keta ) then
263  nw = 0
264  exist = exist + 1
265 ! judge to split this
266  call cjudgesplit(inciptcl, x, split)
267  if( .not. split ) cycle
268 
269  call csplitmeson(pstore(i), e1, e2, icon)
270  if(icon == 0) then
271  modified = .true.
272  ! split; if impossilbe, do nothing
273  nw = nw + 1
274  work(nw) = pstore(i)
275  work(nw).fm.p(4) = e1 + work(nw).mass
276  call cae2p(work(nw)) ! adjust momentum
277  nw = nw + 1
278  work(nw) = pstore(i)
279  work(nw).fm.p(4) = e2 + work(nw).mass
280  call cae2p(work(nw))
281 ! original one is replaced by E1
282  pstore(i) = work(1)
283 ! others are appended to pstore.
284  do j = 2, nw ! altough always nw=2
285  nc = nc + 1
286  pstore(nc) = work(j)
287  enddo
288  endif
289  endif
290  endif
291  endif
292  elseif( mode == 1 ) then
293  if( i /= maxi ) then
294 ! discard all pi/K/eta with x>xth (except leanding)
295  if(pstore(i).code == kpion .or. pstore(i).code == kkaon
296  * .or. pstore(i).code == keta ) then
297  if( x > xth ) then
298 ! to neglect, put mass as E
299  pstore(i).fm.p(4) = pstore(i).mass
300  pstore(i).fm.p(1:3) = 0.
301  modified = .true.
302  endif
303  endif
304  endif
305  elseif( mode == 2 ) then
306  if( i /= maxi) then
307 ! discard only pi0/eta
308  if( ( pstore(i).code == kpion .and.
309  * pstore(i).charge == 0 ) .or.
310  * pstore(i).code == keta ) then
311  if(x> xth) then
312  pstore(i).fm.p(4) = pstore(i).mass
313  pstore(i).fm.p(1:3) = 0. ! put zero energy
314  modified = .true.
315  endif
316  endif
317  endif
318  elseif( mode == 0 ) then
319 ! do nothing
320  else
321  write(0,*) ' mode err=',mode, ' in csoftenPiK0'
322  stop
323  endif
324  enddo
325  if(abs(mode) == 3) then
326  nout = nc
327  endif
328  end subroutine csoftenpik0
329 
330 
331  subroutine csoftenpik1(inciptcl, pstore, nin, nout, exist)
332 ! inciptcl made a collision and generated nin particles in pstore
333 ! particle. information may be at CMS or Lab.
334 ! If we want to make the x-dist. softer,
335 ! A) For Cosmos (air target), we may use the Lab system
336 ! (for large X, the distribution would be the same
337 ! as CMS)
338  implicit none
339 #include "Zcode.h"
340 #include "Zptcl.h"
341 
342  type(ptcl):: inciptcl ! input. incident ptcl at collision
343 !
344  integer,intent(in)::nin ! # of ptcls in pstore
345  type(ptcl):: pstore(*) ! input/output. ptcls to be softened.
346  ! Softened ones are put here
347  integer,intent(out)::nout ! # of ptcls re-stored in pstore
348  integer,intent(out)::exist ! # of particcls with X> Xth
349 ! ----------------------------------------------
350  type(ptcl):: work(2)
351  logical ok
352  integer::maxi
353  real(8)::maxE
354  integer::i, j
355  real(8):: x, E0, Ex, u, u1, u2
356  real(8):: E1, E2
357  real(8)::temp
358  integer:: nc ! counting # of ptcls in pstore
359  integer:: nw, icon
360  logical:: split
361 
362  real(8):: xcent, centval
363  xcent = sqrt(xth/xthu*xth*xthl)
364  centval = crejk1(xcent)
365 
366 
367  e0 = inciptcl.fm.p(4) - inciptcl.mass
368  nout = nin
369  exist = 0
370  if( e0lab < e0th) return !!!!!
371  nc = nin
372 ! if pi/K is incident, regards the highest one as leading
373 ! and avoid to modify
374  if(leadingpik .and. (inciptcl.code == kpion .or.
375  * inciptcl.code == kkaon ) ) then
376 ! find max energy ptcl index and E; eta cannot be incident
377  call cgetmaxptcl(pstore, nc, inciptcl.code, inciptcl.charge,
378  * maxi, maxe)
379  else
380  maxi = 0
381  endif
382 
383  do i = 1, nin
384  x = (pstore(i).fm.p(4)-pstore(i).mass)/e0
385  if( x > xth*xthl ) cycle
386  if( abs(mode) == 3 ) then
387  if( i /= maxi ) then
388 ! softening
389  if(x > xth/xthu) then
390 ! modify high X of pi/K/eta
391  if(pstore(i).code == kpion .or.
392  * pstore(i).code == kkaon
393  * .or. pstore(i).code == keta ) then
394  nw = 0
395  exist = exist + 1
396 ! judge to split this
397 !! call cJudgeSplit(inciptcl, x, split)
398  call rndc(u)
399  temp = crejk1(x)
400 ! split = u < 0.25*(temp/centval)
401 ! split = u < 0.125*(temp/centval)
402  split = u < k1u*(temp/centval)
403  if( .not. split ) cycle
404 
405  call csplitmeson(pstore(i), e1, e2, icon)
406  if(icon == 0) then
407  modified = .true.
408  ! split; if impossilbe, do nothing
409  nw = nw + 1
410  work(nw) = pstore(i)
411  work(nw).fm.p(4) = e1 + work(nw).mass
412  call cae2p(work(nw)) ! adjust momentum
413  nw = nw + 1
414  work(nw) = pstore(i)
415  work(nw).fm.p(4) = e2 + work(nw).mass
416  call cae2p(work(nw))
417 ! original one is replaced by E1
418  pstore(i) = work(1)
419 ! others are appended to pstore.
420  do j = 2, nw ! altough always nw=2
421  nc = nc + 1
422  pstore(nc) = work(j)
423  enddo
424  endif
425  endif
426  endif
427  endif
428  elseif( mode == 1 ) then
429  if( i /= maxi ) then
430 ! discard all pi/K/eta with x>xth (except leanding)
431  if(pstore(i).code == kpion .or. pstore(i).code == kkaon
432  * .or. pstore(i).code == keta ) then
433  if( x > xth ) then
434 ! to neglect, put mass as E
435  pstore(i).fm.p(4) = pstore(i).mass
436  pstore(i).fm.p(1:3) = 0.
437  modified = .true.
438  endif
439  endif
440  endif
441  elseif( mode == 2 ) then
442  if( i /= maxi) then
443 ! discard only pi0/eta
444  if( ( pstore(i).code == kpion .and.
445  * pstore(i).charge == 0 ) .or.
446  * pstore(i).code == keta ) then
447  if(x> xth) then
448  pstore(i).fm.p(4) = pstore(i).mass
449  pstore(i).fm.p(1:3) = 0. ! put zero energy
450  modified = .true.
451  endif
452  endif
453  endif
454  elseif( mode == 0 ) then
455 ! do nothing
456  else
457  write(0,*) ' mode err=',mode, ' in csoftenPiK0'
458  stop
459  endif
460  enddo
461  if(abs(mode) == 3) then
462  nout = nc
463  endif
464  end subroutine csoftenpik1
465 
466  subroutine csplitmeson(p, E1, E2, icon)
467  implicit none
468 #include "Zcode.h"
469 #include "Zptcl.h"
470  type(ptcl):: p ! input, a high energy ptcl
471  real(8),intent(out):: E1 ! split ptcl energy KE
472  real(8),intent(out):: E2 ! split ptcl
473  integer,intent(out):: icon ! 0--> split ok, 1--> no split
474  real(8):: u, Emin, Em
475  logical ok
476  integer::count
477 
478  ok = .false.
479 
480  count = 0
481  do while(.not. ok)
482  call rndc(u)
483  u = u*(1.-xth) + xth
484  e1 = u*( p.fm.p(4) - p.mass)
485  e2 = (p.fm.p(4) -p.mass) - e1
486 ! if(E1 > p.mass .and.
487 ! * E2 > p.mass ) then
488  ok = .true.
489 ! else
490 ! count = count + 1
491 ! if(count > 20) then
492 ! icon = 1
493 ! return
494 ! endif
495 ! endif
496  enddo
497  icon = 0
498  end subroutine csplitmeson
499 
500  subroutine cae2p( pc )
501 ! adjust mpmentum by refering to changed E
502 ! keep the Pt same, if possible
503  implicit none
504 #include "Zcode.h"
505 #include "Zptcl.h"
506 
507 
508  type(ptcl):: pc
509 
510  real(8):: E, P2, cf
511 
512  e = pc.fm.p(4)
513 !
514 ! Pt^2 + Pz^2 +m^2= E^2
515 ! so new Pz =sqrt( E^2-Pt^2-m^2)
516 ! the sign is the same as original one
517  p2 =e**2- ( pc.fm.p(1)**2 + pc.fm.p(2)**2 + pc.mass**2)
518  if( p2 > pc.mass**2 ) then
519  pc.fm.p(3) = sign(sqrt(p2), pc.fm.p(3))
520  else
521 ! keep the direction and shirnk the magnitude of p
522  p2 = pc.fm.p(1)**2 + pc.fm.p(2)**2 + pc.fm.p(3)**2
523  if( e <= pc.mass .or. p2 == 0. ) then
524  pc.fm.p(1:3) = 0.
525  else
526  cf = sqrt( (e**2 - pc.mass**2) / p2 )
527  pc.fm.p(1:3) = pc.fm.p(1:3)*cf
528  endif
529  endif
530  end subroutine cae2p
531  subroutine cgetmaxptcl(pstore, nin, pcode, pcharge, maxi, maxE)
532 ! get max energy ptcl with the same code / charge as incident
533 ! if meson is incident, most probably, it is leading.
534  implicit none
535 #include "Zcode.h"
536 #include "Zptcl.h"
537 
538  integer,intent(in):: nin ! # of ptcls in pstore
539  type(ptcl):: pstore(nin)
540 
541  integer(2),intent(in):: pcode ! incident code
542  integer(2),intent(in):: pcharge ! // charge
543  integer,intent(out):: maxi ! index of maxE in pstore //
544  real(8),intent(out):: maxE ! max Energy with the same code/charg as
545  ! incident. if there is no such, 0
546 
547 
548  integer i
549 
550  maxi = 0
551  maxe = 0.
552  do i = 1, nin
553  if( pstore(i).code /= pcode ) cycle
554  if( pstore(i).charge /= pcharge ) cycle
555  if( maxe < pstore(i).fm.p(4) ) then
556  maxe = pstore(i).fm.p(4)
557  maxi = i
558  endif
559  enddo
560  end subroutine cgetmaxptcl
561 
562  subroutine csoftenfe(inci, fwbwin, a, nin, nout)
563 ! front end of softening when it is to be done at CMS
564 ! x is defined simply by E/E0,
565  implicit none
566 #include "Zptcl.h"
567  type(ptcl):: inci ! incident ptcl (in cms); symmetric case
568  integer,intent(in):: fwbwin ! modification is applied to
569  ! 1--> forward only 2--> backword only
570  ! 3--> both ; However, 1 is used when
571  ! Cosmos output (in cms) is the target, independently of
572  ! fwbw.
573  type(ptcl):: a(*) ! array containing ptcl info.
574  integer,intent(in)::nin ! # of ptcls in w
575  integer,intent(out)::nout ! # of ptcls put in w after modification
576 
577  type(ptcl):: w1(half) ! working array
578  type(ptcl):: w2(half) ! working array
579 
580  integer::i, nc1, nc2, ncout1, ncout2
581 
582 ! do modification extract Pz>0 and Pz<0
583  nc1 = 0
584  nc2 = 0
585  do i = 1, nin
586  if(a(i).fm.p(3) > 0.) then
587  nc1 = nc1 + 1
588  w1(nc1) = a(i)
589  else
590  nc2 = nc2 + 1
591  w2(nc2) = a(i)
592  endif
593  enddo
594  if( nc1 > 0 .and. ibits(fwbwin,0,1)>0 ) then ! LSB pos=0
595  call csoftenpik(inci, w1, nc1, ncout1)
596  else
597  ncout1 = nc1
598  endif
599  if( nc2 > 0 .and. ibits(fwbwin,1,1)>0 ) then ! 2nd bit pos=1
600  call csoftenpik(inci, w2, nc2, ncout2)
601  else
602  ncout2 = nc2
603  endif
604 
605  nout = 0
606  do i = 1, ncout1
607  nout = nout + 1
608  a(nout) = w1(i)
609  enddo
610  do i = 1, ncout2
611  nout = nout + 1
612  a(nout) = w2(i)
613  enddo
614  end subroutine csoftenfe
615 
616  subroutine creadsoftenpara(io)
617  implicit none
618  integer,intent(in):: io ! logical dev. #
619  character*24 vname
620  character*100 vvalue
621 
622 
623  call cskipsep(io)
624  do while( cgetparmn(io, vname, vvalue ) )
625  select case(vname)
626  case('mode')
627  call creadparai(vvalue, mode)
628  case('xth')
629  call creadparar(vvalue, xth)
630  case('E0th')
631  call creadparar(vvalue, e0th)
632  case('fwbw')
633  call creadparai(vvalue, fwbw)
634  case('pw')
635  call creadparar(vvalue, pw)
636  case('repeat')
637  call creadparar(vvalue, repeat)
638  case('special')
639  call creadparal(vvalue, special)
640  case('leadingPiK')
641  call creadparal(vvalue, leadingpik)
642  case('useXinCMS')
643  call creadparal(vvalue, usexincms)
644  case('k1u')
645  call creadparar(vvalue, k1u)
646  case('xthl')
647  call creadparar(vvalue, xthl)
648  case('xthu')
649  call creadparar(vvalue, xthu)
650  case('rejpw')
651  call creadparar(vvalue, rejpw)
652  end select
653  enddo
654  end subroutine creadsoftenpara
655 ! *************
656  subroutine cwritesoftenpara(io)
657  implicit none
658  integer,intent(in):: io
659 
660  write(io,*)'----------------------'
661  call cwriteparai(io,'mode', mode)
662  call cwriteparar(io,'xth', xth)
663  call cwriteparar(io,'E0th', e0th)
664  call cwriteparai(io,'fwbw', fwbw)
665  call cwriteparar(io,'pw', pw)
666  call cwriteparar(io,'repeat',repeat)
667  call cwriteparal(io,'special',special)
668  call cwriteparal(io,'leadingPiK',leadingpik)
669  call cwriteparal(io,'useXinCMS', usexincms)
670  call cwriteparar(io,'k1u', k1u)
671  call cwriteparar(io,'xthl', xthl)
672  call cwriteparar(io,'xthu', xthu)
673  call cwriteparar(io,'rejpw', rejpw)
674  end subroutine cwritesoftenpara
675 
676 
677  subroutine cskipsep(io)
678  implicit none
679  integer io
680  character(10) sep
681  do while (.true.)
682  read(io, '(a)') sep
683  if(sep(2:10) == '---------') exit
684  enddo
685  end subroutine cskipsep
686 ! ************************* real*8 data
687  subroutine creadparar(vvalue, x)
688  implicit none
689  integer io
690  character*(*) vvalue
691  real*8 x
692 ! read(vvalue, *) x, x
693  read(vvalue, *) x
694  end subroutine creadparar
695  subroutine creadparar2(vvalue, x, n)
696  implicit none
697  integer io
698  character*(*) vvalue
699  integer n
700  real*8 x(n)
701  read(vvalue, *) x
702  end subroutine creadparar2
703 
704 ! ************************* complex data
705  subroutine creadparacx(vvalue, c)
706  implicit none
707  character*(*) vvalue
708  complex*8 c
709  read( vvalue, *) c
710  end subroutine creadparacx
711 ! ************************ integer data
712  subroutine creadparai(vvalue, i)
713  implicit none
714  character*(*) vvalue
715  integer i
716  read(vvalue, *) i
717  end subroutine creadparai
718 ! ************************* character data
719  subroutine creadparac(vvalue, cha)
720  implicit none
721  character*(*) vvalue
722  character*(*) cha
723  read(vvalue, *) cha
724  end subroutine creadparac
725 ! ***************************** logical data
726  subroutine creadparal(vvalue, logi)
727  implicit none
728  character*(*) vvalue
729  logical logi
730  read(vvalue, *) logi
731  end subroutine creadparal
732 ! ---------------------------------------------
733  subroutine cwriteparar(io, vname, x)
734  implicit none
735  integer io
736  character*(*) vname
737  real*8 x
738 
739  write(io, *) ' ', vname,' ', x,' /'
740  end subroutine cwriteparar
741  subroutine cwriteparar2(io, vname, x, n)
742  implicit none
743  integer io
744  integer n ! arra size of x
745  character*(*) vname
746  real*8 x(n)
747 
748  write(io,*) ' ', vname,' ', x,' /'
749  end subroutine cwriteparar2
750 
751  subroutine cwriteparacx(io, vname, c)
752  implicit none
753  integer io
754  character*(*) vname
755  complex*8 c
756  write(io, *) ' ', vname,' ', c,' /'
757  end subroutine cwriteparacx
758 
759  subroutine cwriteparai(io, vname, i)
760  implicit none
761  integer io
762  character*(*) vname
763  integer i
764 
765  write(io, *) ' ', vname,' ', i,' /'
766  end subroutine cwriteparai
767 
768  subroutine cwriteparac(io, vname, cha)
769  implicit none
770  integer io
771  character*(*) vname
772  character*(*) cha
773  integer klena
774  character*2 qmk/" '"/ ! '
775  if(klena(cha) .gt. 0) then
776  write(io, *) ' ', vname, qmk, cha(1:klena(cha)),
777  * qmk,' /'
778  else
779  write(io, *) ' ', vname, qmk, ' ', qmk, ' /'
780  endif
781  end subroutine cwriteparac
782  subroutine cwriteparal(io, vname, logi)
783  implicit none
784  integer io
785  character*(*) vname
786  logical logi
787 
788  write(io, *) ' ', vname,' ', logi,' /'
789  end subroutine cwriteparal
790 
791  function crejk1( x ) result(ans)
792  implicit none
793  real(8),intent(in):: x
794  real(8):: ans
795  ans = ( log(x/(xth/xthu)) * log((xth*xthl)/x) )**rejpw
796  end function crejk1
797 
798  function cxthu( E0 ) result(ans)
799  implicit none
800  real(8),intent(in):: E0 ! lab E0 in Gev
801  real(8):: ans ! xthu
802  ans = 12.0*(e0/1.e8)**0.1
803  end function cxthu
804 
805  function ck1u( E0 ) result(ans)
806  implicit none
807  real(8),intent(in):: E0 ! lab E0 in Gev
808  real(8):: ans ! k1u
809  ans = 0.3*(e0/1.e8)**0.09
810  end function ck1u
811 
812  function cgetparmn( io, vname, vv ) result(ans)
813 ! get parameter variable name and given value(s) from io
814  implicit none
815  integer io
816  character*(*) vname, vv ! output
817  logical ans
818 
819  integer linel
820  parameter( linel = 100)
821  character*(linel) line
822  integer loc, loc2
823  vname = " "
824  do while(.true.)
825  read(io, '(a)', end=100 ) line
826  if( line(1:1) .eq. " " .and. line(2:2) .ne. " ") then
827  loc = index( line(3:linel), " ") + 2
828  vname = line(2:loc-1)
829  loc2 = index( line, "/")
830  if(loc2 .eq. 0 ) then
831  write(0,* ) ' "/" is missing in the input data file '
832  write(0,*) ' The line is: ', line
833  stop 1234
834  endif
835  vv = line(loc+1:linel) ! some data containes '/' such as '../../Media' so put all
836  ! data.
837  goto 50
838  endif
839  enddo
840  50 continue
841  ans = .true.
842  return
843  100 continue
844  ans =.false.
845  end function cgetparmn
846 
847  function csoftenfixxth(E0) result(xth)
848  implicit none
849  real(8),intent(in)::E0 ! proton/pi/K incident E. in GeV
850  real(8):: xth !
851 ! fix the xth above which softening is performed
852 ! at 10^12 eV: 0.1
853 ! 10^13 0.063
854 ! 10^14 0.04
855 ! 10^16 0.01585
856 ! 10^17 0.01
857 ! 10^19 0.004
858  xth = 0.1d0/(e0/1000.d0)**0.2
859  end function csoftenfixxth
860 
861  end module softenpik
integer npitbl real *nx parameter(n=101, npitbl=46, nx=n-1) real *8 uconst
subroutine csoftenpik(inciptcl, pstore, nin, nout)
Definition: csoftenPiK.f:142
subroutine cwriteparai(io, vname, i)
real(8), save rejpw
Definition: csoftenPiK.f:101
subroutine creadparac(vvalue, cha)
logical, save special
Definition: csoftenPiK.f:55
integer, save fwbw
Definition: csoftenPiK.f:50
logical, save modified
Definition: csoftenPiK.f:102
logical, save usexincms
Definition: csoftenPiK.f:90
subroutine creadparal(vvalue, logi)
subroutine creadparar2(vvalue, x, n)
max ptcl codes in the kkaon
Definition: Zcode.h:2
real(8) function cxthu(E0)
Definition: csoftenPiK.f:799
logical, save leadingpik
Definition: csoftenPiK.f:81
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
logical function cgetparmn(io, vname, vv)
integer, save mode
Definition: csoftenPiK.f:34
subroutine cjudgesplit(pj, x, split)
Definition: csoftenPiK.f:114
subroutine rndc(u)
Definition: rnd.f:91
subroutine creadsoftenpara(io)
********************block data cblkHeavy ********************integer j data *HeavyG2symbol p
Definition: cblkHeavy.h:7
real(8), save xthl
Definition: csoftenPiK.f:99
subroutine kpoisn(am, n)
Definition: kpoisn.f:25
*Zfirst p fm *Zfirst p Zfirst p code
Definition: ZavoidUnionMap.h:1
real(8), save pw
Definition: csoftenPiK.f:36
subroutine creadparar(vvalue, x)
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
real(8), save e0th
Definition: csoftenPiK.f:46
subroutine copenf(io, fnin, icon)
Definition: copenf.f:8
real(8) function csoftenfixxth(E0)
Definition: csoftenPiK.f:848
real(8), save repeat
Definition: csoftenPiK.f:40
subroutine cwriteparacx(io, vname, c)
real(8) function ck1u(E0)
Definition: csoftenPiK.f:806
subroutine cwriteparal(io, vname, logi)
real(8) function crejk1(x)
Definition: csoftenPiK.f:792
subroutine creadparacx(vvalue, c)
subroutine cwriteparar2(io, vname, x, n)
subroutine cae2p(pc)
Definition: csoftenPiK.f:501
real(8), save k1u
Definition: csoftenPiK.f:98
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p charge
Definition: ZavoidUnionMap.h:1
subroutine csoftenfe(inci, fwbwin, a, nin, nout)
Definition: csoftenPiK.f:563
subroutine csplitmeson(p, E1, E2, icon)
Definition: csoftenPiK.f:467
subroutine csoftenpik1(inciptcl, pstore, nin, nout, exist)
Definition: csoftenPiK.f:332
subroutine cgetmaxptcl(pstore, nin, pcode, pcharge, maxi, maxE)
Definition: csoftenPiK.f:532
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 ! knockon is considered Obsolete *PhotoProd false
Definition: cblkElemag.h:7
real(8), save xthu
Definition: csoftenPiK.f:100
max ptcl codes in the keta
Definition: Zcode.h:2
subroutine cwriteparac(io, vname, cha)
subroutine csoftenpik0(inciptcl, pstore, nin, nout, exist)
Definition: csoftenPiK.f:207
subroutine creadparai(vvalue, i)
Definition: Zptcl.h:75
subroutine cwriteparar(io, vname, x)
*Zfirst p fm *Zfirst p mass
Definition: ZavoidUnionMap.h:1
max ptcl codes in the kpion
Definition: Zcode.h:2
real(8), save e0lab
Definition: csoftenPiK.f:103
integer, parameter half
Definition: csoftenPiK.f:108
subroutine cwritesoftenpara(io)
subroutine cskipsep(io)
real(8), save xth
Definition: csoftenPiK.f:35