Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
c----------------------------------------------------------------------
c The two elementary fuctions of our approach, the profile fuction G
c and the Phi exponent G~, are here referred to as Gpro and Gtilde.
c Both functions can be written as
c
c G = \sum_type alp * xp**bet * xm**betp
c
c The parameters alp, bet, betp are calculated in GfunParK (k-mode,
c b is takento the one of pair k) or GfunPar (b-mode: arbitrary b) as
c
c Gpro: bet = betD' + epsGp + gamD*b**2 + epsG -alppar
c bet' = betD'' + epsGt + gamD*b**2 + epsG -alppar
c
c alp = alpD*f * s**(betD+gamD*b**2+epsG) * exp(-b**2/delD)
c
c Gtilde: bet~ = bet + 1
c bet~' = bet' + 1
c
c alp~ = alp * gam(bet~) * gam(bet~')
c * gam(1+alppro) * gam(1+alptar)
c * gam(1+alppro+bet~) * gam(1+alptar+bet~')
c * (1+epsGt') * (1+epsGt')
c
c The parameters depend implicitely on s.
c
c In the program we use om1 = Gpro
c (they differ by a constant which is actually one)
c and om5 = om1 * 0.5
c All functions related to om1 are called om1... .
c
c The inclusive Pomeron distributions are
c
c PomInc(xp,xm) = Gpro(xp,xm) * (1-xp)**alppro * (1-xm)**alptar
c
c----------------------------------------------------------------------
c----------------------------------------------------------------------
subroutine GfunParK(irea) !---MC---
c----------------------------------------------------------------------
c calculates parameters alp,bet,betp of the G functions (k-mode)
c and the screening exponents epsilongp(k,i), epsilongt(k,i), epsilongs(k,i)
c----------------------------------------------------------------------
c Gpro parameters written to /comtilde/atilde(,)btildep(,),btildepp(,)
c Gtilde parameters written to /cgtilde/atildg(,)btildgp(,),btildgpp(,)
c two subscripts: first=type, second=collision k
c Certain pieces of this routine are only done if irea is <= or = zero.
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incems'
include 'epos.incpar'
double precision atildg,btildgp,btildgpp
common/cgtilde/atildg(idxD0:idxD1,kollmx)
*,btildgp(idxD0:idxD1,kollmx),btildgpp(idxD0:idxD1,kollmx)
double precision utgam2,coefgdp,coefgdt
parameter(nbkbin=40)
common /kfitd/ xkappafit(nclegy,nclha,nclha,nbkbin),xkappa,bkbin
common /cgtilnu/ cfbetpnp,cfbetppnp,cfbetpnm,cfbetppnm,cfalpro
&,cfaltar,cfbpap,cfbpam,cfbppap,cfbppam
double precision cfbetpnp,cfbetppnp,cfbetpnm,cfbetppnm,cfalpro
&,cfaltar,cfbpap,cfbpam,cfbppap,cfbppam,gamv,eps
parameter (eps=1.d-25)
parameter(nxeps=20,nyeps=32)
common/cxeps1/w(0:nxeps,nyeps),y1(nyeps),y2(nyeps)
common/cxeps2/db,b1,b2
common/geom/rmproj,rmtarg,bmax,bkmx
common/nucl3/phi,bimp
common /cncl/xproj(mamx),yproj(mamx),zproj(mamx)
* ,xtarg(mamx),ytarg(mamx),ztarg(mamx)
b1=0.03
b2=bkmx*1.2
db=(b2-b1)/nyeps
call utprj('GfunParK ',ish,ishini,10)
cfbetpnp=0.d0
cfbetppnp=0.d0
cfbetpnm=0.d0
cfbetppnm=0.d0
cfalpro=dble(ucfpro)
cfaltar=dble(ucftar)
cfbpap=1.d0
cfbppap=1.d0
cfbpam=1.d0
cfbppam=1.d0
alpfom=0.
sy=engy*engy
do k=1,koll
do i=ntymi,ntymx
atilde(i,k)=0.d0
btildep(i,k)=0.d0
btildepp(i,k)=0.d0
enddo
do i=idxD0,idxD1
atildg(i,k)=0.d0
btildgp(i,k)=0.d0
btildgpp(i,k)=0.d0
enddo
enddo
! calculate collision number according to Glauber -----------------------
bglaub=sqrt(sigine/10./pi)
nglevt=0
do ko=1,koll
r=bk(ko)
if(r.le.bglaub)nglevt=nglevt+1
enddo
! Z_parton_projectile (zparpro) and Z_parton_target (zpartar)-----------
ztav=0.
zpav=0.
if(iscreen.eq.1.or.isplit.eq.1)then
b2x=b2xscr
alpfom=alpfomi*fegypp
rho0p=conrho(1,0.)
rho0t=conrho(2,0.)
bcut=0.
do k=1,koll
ip=iproj(k)
it=itarg(k)
!....... targ partons seen by proj
if(lproj(ip).gt.0)then
absb=max(1.e-9,bk(k))
b2=absb*absb
zkp=fegypp*exp(-b2/b2x)
zpartar(k)=min(zkp,epscrx)
if(lproj3(ip).gt.1)then
do lt=1,lproj3(ip)
kp=kproj3(ip,lt)
if(kp.ne.k)then
ikt=itarg(kp)
rtarg=sqrt(xtarg(ikt)**2+ytarg(ikt)**2+ztarg(ikt)**2)
rho=conrho(2,rtarg)/rho0t
fegyAA=epscrw*fscro(engy/egyscr,rho)
absb=max(1.e-9,abs(bk(kp))-bcut)
b2=absb*absb
zkp=fegyAA*exp(-b2/b2x)
zpartar(k)=zpartar(k)+min(zkp,epscrx)
endif
c alpfom=max(alpfom,dble(zpartar(k)))
enddo
endif
else
zpartar(k)=0.
endif
ztav=ztav+zpartar(k)
!...........proj partons seen by targ
if(ltarg(it).gt.0)then
absb=max(1.e-9,bk(k))
b2=absb*absb
zkt=fegypp*exp(-b2/b2x)
zparpro(k)=min(zkt,epscrx)
if(ltarg3(it).gt.1)then
do lp=1,ltarg3(it)
kt=ktarg3(it,lp)
if(kt.ne.k)then
ikp=iproj(kt)
rproj=sqrt(xproj(ikp)**2+yproj(ikp)**2+zproj(ikp)**2)
rho=conrho(1,rproj)/rho0p
fegyAA=epscrw*fscro(engy/egyscr,rho)
absb=max(1.e-9,abs(bk(kt))-bcut)
b2=absb*absb
zkt=fegyAA*exp(-b2/b2x)
zparpro(k)=zparpro(k)+min(zkt,epscrx)
endif
c alpfom=max(alpfom,dble(zparpro(k)))
enddo
endif
else
zparpro(k)=0.
endif
zpav=zpav+zparpro(k)
xzcutpar(k)=dble(exp(-min(50.,xzcut*(zparpro(k)+zpartar(k)))))
enddo
else ! no screening
do k=1,koll
zparpro(k)=0.
zpartar(k)=0.
xzcutpar(k)=1d0
enddo
endif
c calculation of epsilongp epsilongt
if(iscreen.eq.1)then
!ip=0
do k=1,koll
!ipp=ip
epsG=0.
epsilongs(k,0)=0.
epsilongs(k,1)=0.
ip=iproj(k) !...........projectile
if(lproj(ip).gt.0)then
x=zparpro(k)
epsilongs(k,0)=sign(abs(epscrd)*x
& ,epscrd)
epsilongp(k,0)=max(-betDp(0,iclpro,icltar)-0.95+alppar,
& epscrs*x)
c & min(epscrx,epscrs*x))
if(sy.ge.sfshlim)then
epsilongp(k,1)=max(-betDp(1,iclpro,icltar)-0.95+alppar,
& epscrh*x)
c & min(epscrx,epscrh*x))
else
epsilongp(k,1)=epsilongp(k,0)
c epsilongs(k,1)=epsilongs(k,0)
endif
gammaV(k)=1.d0
else
epsilongp(k,0)=0.
epsilongp(k,1)=0.
gammaV(k)=1.d0
endif
it=itarg(k) !...........target
if(ltarg(it).gt.0)then
x=zpartar(k)
epsilongs(k,1)=sign(abs(epscrd)*x
& ,epscrd)
epsilongt(k,0)=max(-betDpp(0,iclpro,icltar)-0.95+alppar,
& epscrs*x)
c & min(epscrx,epscrs*x))
if(sy.ge.sfshlim)then
epsilongt(k,1)=max(-betDpp(1,iclpro,icltar)-0.95+alppar,
& epscrh*x)
c & min(epscrx,epscrh*x))
else
epsilongt(k,1)=epsilongt(k,0)
c epsilongs(k,1)=epsilongs(k,0)
endif
cc gammaV(k)=gammaV(k)
else
epsilongt(k,0)=0.
epsilongt(k,1)=0.
gammaV(k)=gammaV(k)
endif
enddo
else ! no screening
do k=1,koll
epsilongs(k,0)=0.
epsilongs(k,1)=0.
epsilongp(k,0)=0.
epsilongp(k,1)=0.
epsilongt(k,0)=0.
epsilongt(k,1)=0.
gammaV(k)=1.d0
enddo
endif
!..............alpha beta betap for Gtilde (->PhiExpo).......................
imax=idxD1
if(iomega.eq.2)imax=1
do k=1,koll
b=bk(k)
b2=bk(k)*bk(k)
ip=iproj(k)
it=itarg(k)
if(b.lt.(nbkbin-1)*bkbin)then
ibk=int(bk(k)/bkbin)+1
if(isetcs.gt.1.and.iclegy.lt.iclegy2)then
egy0=egylow*egyfac**float(iclegy-1)
xkappa1=xkappafit(iclegy,iclpro,icltar,ibk)
* +(bk(k)-bkbin*float(ibk-1))/bkbin
* *(xkappafit(iclegy,iclpro,icltar,ibk+1)
* -xkappafit(iclegy,iclpro,icltar,ibk))
xkappa2=xkappafit(iclegy+1,iclpro,icltar,ibk)
* +(bk(k)-bkbin*float(ibk-1))/bkbin
* *(xkappafit(iclegy+1,iclpro,icltar,ibk+1)
* -xkappafit(iclegy+1,iclpro,icltar,ibk))
xkappa=xkappa1+(xkappa2-xkappa1)/log(egyfac)
* *(log(engy)-log(egy0))
xkappa=facmc*xkappa
else
xkappa=xkappafit(iclegy,iclpro,icltar,ibk)
* +(bk(k)-bkbin*float(ibk-1))/bkbin
* *(xkappafit(iclegy,iclpro,icltar,ibk+1)
* -xkappafit(iclegy,iclpro,icltar,ibk))
xkappa=facmc*xkappa
endif
else
xkappa=1.
endif
gfactorp=1.!+(gfactor-1)*exp(-5*b/gwidth/bglaub)
gfactort=1.!+(gfactor-1)*exp(-5*b/gwidth/bglaub)
do i=idxDmin,imax
gamV=gammaV(k)
c if(i.lt.2)then
c if(i.eq.0)then
c epsG=epsilongs(k,i)
if(i.eq.2)then
if(epscrd.lt.0.)then
epsG=epsilongs(k,0)+epsilongs(k,1)
epsGp=0.
epsGt=0.
else
epsG=0.
epsGp=epsilongs(k,0)
epsGt=epsilongs(k,1)
endif
else
epsG=0.
epsGp=0.
epsGt=0.
endif
gamb=gamD(i,iclpro,icltar)*b2
atildg(i,k)=dble(alpD(i,iclpro,icltar))
* *cfalpro*cfaltar
* *gamv
c if(i.eq.0) atildg(i,k)=atildg(i,k)
atildg(i,k)=atildg(i,k)
* *dble(xkappa*xkappa)
if(i.lt.2)then
atildg(i,k)=atildg(i,k)*dble(
* chad(iclpro)*chad(icltar)
* *exp(-b2/delD(i,iclpro,icltar))
* *sy**(betD(i,iclpro,icltar)+gamb+epsG)
* *gfactorp *gfactort)
epsGp=epsilongp(k,i)
epsGt=epsilongt(k,i)
btildgp(i,k)=dble(betDp(i,iclpro,icltar)
* +epsGp
* +gamb-alppar)+1.d0
btildgpp(i,k)=dble(betDpp(i,iclpro,icltar)
* +epsGt
* +gamb-alppar)+1.d0
else
absb=abs(bk(k))-bmxdif(iclpro,icltar)
b2a=absb*absb
atildg(i,k)=atildg(i,k)*dble(
* sy**(betD(i,iclpro,icltar)+epsG)
* *exp(-b2a/delD(i,iclpro,icltar)))
btildgp(i,k)=dble(betDp(i,iclpro,icltar)-alppar+epsGp)+1.d0
btildgpp(i,k)=dble(betDpp(i,iclpro,icltar)-alppar+epsGt)+1.d0
endif
coefgdp=utgam2(1.d0+dble(alplea(iclpro))+btildgp(i,k))
coefgdt=utgam2(1.d0+dble(alplea(icltar))+btildgpp(i,k))
atildg(i,k)=atildg(i,k)
* *utgam2(btildgp(i,k))*utgam2(btildgpp(i,k))
* /coefgdp/coefgdt
!...........prepare plot in xepsilon
if(irea.eq.0)then
kk=max(1,int((bk(k)-b1)/db)+1)
if(i.lt.2)then
if(i.eq.0)w(0,kk)=w(0,kk)+1
if(i.eq.0)w(1,kk)=w(1,kk)+epsGp
if(i.eq.0)w(2,kk)=w(2,kk)+epsGt
c w(3+i,kk)=w(3+i,kk)+abs(epsG)
!...5-8 soft ... 9-12 semi
w(5+4*i,kk)=w(5+4*i,kk)
* +betDp(i,iclpro,icltar) !prj eff
* +epsGp+gamb
w(6+4*i,kk)=w(6+4*i,kk)
* +betDpp(i,iclpro,icltar) !tgt eff
* +epsGt+gamb
w(7+4*i,kk)=w(7+4*i,kk)
* +betDp(i,iclpro,icltar) !prj unscr
* +gamb
w(8+4*i,kk)=w(8+4*i,kk)
* +betDpp(i,iclpro,icltar) !tgt unscr
* +gamb
if(i.eq.0)w(13,kk)=w(13,kk)+zparpro(k)
if(i.eq.0)w(14,kk)=w(14,kk)+zpartar(k)
else
if(epscrd.lt.0.)then
w(3,kk)=w(3,kk)+abs(epsG)
else
w(3,kk)=w(3,kk)+epsGp
w(4,kk)=w(4,kk)+epsGt
endif
endif
endif
!................
enddo
enddo
!...........................................................................
zppevt=zpav/koll
zptevt=ztav/koll
if(irea.eq.0)then
ktot=int(bimp)+1
if(ktot.le.nyeps)then
w(15,ktot)=w(15,ktot)+zppevt
w(16,ktot)=w(16,ktot)+zptevt
w(17,ktot)=w(17,ktot)+1
endif
n=1+int(float(nglevt)/(0.1*maproj*matarg))*(nyeps-1)
if(nglevt.ge.1.and.n.ge.1.and.n.le.nyeps)then
w(18,n)=w(18,n)+zppevt
w(19,n)=w(19,n)+zptevt
w(20,n)=w(20,n)+1
endif
endif
!........alpha beta betap for Gpro...........................................
if(irea.le.0)then
do k=1,koll
ip=iproj(k)
it=itarg(k)
b2=bk(k)*bk(k)
imax=ntymx
if(iomega.eq.2)imax=1
do i=ntymin,imax
c if(i.lt.2)then
c if(i.eq.0)then
c epsG=epsilongs(k,i)
if(i.eq.2)then
if(epscrd.lt.0.)then
epsG=epsilongs(k,0)+epsilongs(k,1)
epsGp=0.
epsGt=0.
else
epsG=0.
epsGp=epsilongs(k,0)
epsGt=epsilongs(k,1)
endif
else
epsG=0.
epsGp=0.
epsGt=0.
endif
gamb=gamD(i,iclpro,icltar)*b2
atilde(i,k)=dble(alpD(i,iclpro,icltar))
if(i.lt.2)then
atilde(i,k)=atilde(i,k)*dble(
* exp(-b2/delD(i,iclpro,icltar))
* *sy**(betD(i,iclpro,icltar)
* +gamb+epsG)
* *chad(iclpro)*chad(icltar))
epsGp=epsilongp(k,i)
epsGt=epsilongt(k,i)
btildep(i,k)=dble(betDp(i,iclpro,icltar)
* +epsGp
* +gamb-alppar)
btildepp(i,k)=dble(betDpp(i,iclpro,icltar)
* +epsGt
* +gamb-alppar)
else
absb=abs(bk(k))-bmxdif(iclpro,icltar)
b2a=absb*absb
atilde(i,k)=atilde(i,k)*dble(
* sy**(betD(i,iclpro,icltar)+epsG)
* *exp(-b2a/delD(i,iclpro,icltar)))
btildep(i,k)=dble(betDp(i,iclpro,icltar)-alppar+epsGp)
btildepp(i,k)=dble(betDpp(i,iclpro,icltar)-alppar+epsGt)
endif
if(btildep(i,k)+1d0.lt.-eps.or.btildepp(i,k)+1d0.lt.-eps)then
write(ifmt,*)' k,b,ip,it,gamb,alppar',k,bk(k),ip,it,gamb
* ,alppar
write(ifmt,*)' gammaV,epsGP1/2,epsGT1/2,epsGS1/2'
* ,gammaV(k),epsilongp(k,0),epsilongp(k,1)
* ,epsilongt(k,0),epsilongt(k,1),epsilongs(k,0),epsilongs(k,1)
write(ifmt,*)'*******************************************'
write(ifmt,*)" atilde,btildep,btildepp "
do ii=ntymin,ntymx
write(ifmt,*)ii,atilde(ii,k),btildep(ii,k),btildepp(ii,k)
enddo
call utstop('Error in epos-omg in GfunPark&')
endif
enddo
enddo
endif
!...........................................................................
if(ish.ge.10)then
do k=1,koll
ip=iproj(k)
it=itarg(k)
write(ifch,*)' k,b,ip,it,',k,bk(k),ip,it
write(ifch,*)' zparpro,zpartar,xzcutpar'
* ,zparpro(k),zpartar(k),xzcutpar(k)
write(ifch,*)' gammaV,epsilonGP1/2,epsilonGT1/2,epsilonGs1/2'
* ,gammaV(k),epsilongp(k,0),epsilongp(k,1)
* ,epsilongt(k,0),epsilongt(k,1),epsilongs(k,0),epsilongs(k,1)
write(ifch,*)'*******************************************'
write(ifch,*)" atilde,btildep,btildepp "
do i=ntymin,ntymx
write(ifch,*)i,atilde(i,k),btildep(i,k),btildepp(i,k)
enddo
write(ifch,*)" atildg,btildgp,btildgpp "
do i=ntymin,ntymx
write(ifch,*)i,atildg(i,k),btildgp(i,k),btildgpp(i,k)
enddo
call GfunPar(xpar7,xpar7,1,0,bk(k),sy,alp,bet,betp
& ,epsp,epst,epss,gamvv)
call GfunPar(xpar7,xpar7,1,1,bk(k),sy,alp,bet,betp
& ,epsp,epst,epss,gamvv)
enddo
endif
call utprjx('GfunParK ',ish,ishini,10)
return
end
c----------------------------------------------------------------------
subroutine GfunPar(zzip,zzit,m,i,b,spp,alp,bet,betp,epsp,epst
& ,epss,gamvv)
c----------------------------------------------------------------------
c calculates parameters alp,bet,betp of the G functions for pp (b-mode)
c and the screening exponents epsp,epst,epss,gamvv
c----------------------------------------------------------------------
c zzip:additional z component (nuclear effect projectile side)
c zzit:additional z component (nuclear effect target side)
c m=1: profile function Gpro, i=0: soft i=2: diff
c m=2: Gtilde, i=1: semi (no screening for diff)
c b: impact param, spp: pp energy squared
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incpar'
include 'epos.incems'
parameter(nbkbin=40)
common /kfitd/ xkappafit(nclegy,nclha,nclha,nbkbin),xkappa,bkbin
common /kwrite/ xkapZ
double precision utgam2,coefgdp,coefgdt,dalp,dbet,dbetp,eps
parameter(eps=1.d-20)
call utprj('GfunPar ',ish,ishini,10)
ee=sqrt(spp)
c bglaub2=FbGlaub2(ee)
rs=r2had(iclpro)+r2had(icltar)+slopom*log(spp)
bglaub2=4.*.0389*rs
if(ish.ge.10)write(ifch,*)'Gf in:',m,i,b,bglaub2,spp,zzip,zzit
& ,iclpro,icltar
b2=b*b
cfalpro=ucfpro
cfaltar=ucftar
gamb=gamD(i,iclpro,icltar)*b2
if(iscreen.ne.0)then
absb=max(1.e-9,abs(b))
b2a=absb*absb
b2x=2.*epscrp*bglaub2
zzp=epscrw*exp(-b2a/b2x)*fscra(ee/egyscr)
zzp=min(zzp,epscrx)+zzip !saturation
zzt=epscrw*exp(-b2a/b2x)*fscra(ee/egyscr)
zzt=min(zzt,epscrx)+zzit !saturation
x=zzp
epsG=0.
if(i.eq.1.and.spp.ge.sfshlim)then
epsGp=max(-betDp(i,iclpro,icltar)-0.95+alppar,
& epscrh*x)
c & min(epscrx,epscrh*x))
elseif(i.le.1)then
epsGp=max(-betDp(i,iclpro,icltar)-0.95+alppar,
& epscrs*x)
c & min(epscrx,epscrs*x))
else
if(epscrd.lt.0.)then
epsG=epsG+sign(abs(epscrd)*x,epscrd)
epsGp=0.
else
epsGp=sign(abs(epscrd)*x,epscrd)
epsG=0.
endif
endif
gamV=1.
x=zzt
if(i.eq.1.and.spp.ge.sfshlim)then
epsGt=max(-betDpp(i,iclpro,icltar)-0.95+alppar,
& epscrh*x)
c & min(epscrx,epscrh*x))
elseif(i.le.1)then
epsGt=max(-betDpp(i,iclpro,icltar)-0.95+alppar,
& epscrs*x)
c & min(epscrx,epscrs*x))
else
if(epscrd.lt.0.)then
epsG=epsG+sign(abs(epscrd)*x,epscrd)
epsGt=0.
else
epsGt=sign(abs(epscrd)*x,epscrd)
epsG=0.
endif
endif
c gamV=gamV
else
zzp=0.
zzt=0.
epsGp=0.
epsGt=0.
epsG=0.
gamV=1.
endif
gfactorp=1.!+(gfactor-1)*exp(-5*b/gwidth/bglaub)
gfactort=1.!+(gfactor-1)*exp(-5*b/gwidth/bglaub)
rho=betD(i,iclpro,icltar)+gamb+epsG
if(m.eq.1)then
dalp=dble(alpD(i,iclpro,icltar))
if(i.lt.2)then
dalp=dalp
* *exp(min(50d0,dble(rho*log(spp)-b2/delD(i,iclpro,icltar))))
dbet=dble(betDp(i,iclpro,icltar)
* +epsGp
* +gamb-alppar)
dbetp=dble(betDpp(i,iclpro,icltar)
* +epsGt
* +gamb-alppar)
else
absb=abs(b)-bmxdif(iclpro,icltar)
b2a=absb*absb
dalp=dalp
* *exp(min(50d0,dble((betD(i,iclpro,icltar)+epsG)*log(spp)
* -b2a/delD(i,iclpro,icltar))))
dbet=dble(betDp(i,iclpro,icltar)-alppar+epsGp)
dbetp=dble(betDpp(i,iclpro,icltar)-alppar+epsGt)
endif
if((dbet+1.d0).lt.-eps.or.(dbetp+1.d0).lt.-eps)then
write(*,*)'m,i,b,spp,alp,bet,betp',m,i,b,spp,dalp,dbet,dbetp
call utstop('Error : beta < -1 in Gfunpar in epos-omg&')
endif
elseif(m.eq.2)then
xkappa=1.
c if(i.eq.0.and.b.lt.(nbkbin-1)*bkbin)then
if(b.lt.(nbkbin-1)*bkbin)then
ibk=int(b/bkbin)+1
if(isetcs.gt.1.and.iclegy.lt.iclegy2)then
egy0=egylow*egyfac**float(iclegy-1)
xkappa1=xkappafit(iclegy,iclpro,icltar,ibk)
* +(b-bkbin*float(ibk-1))/bkbin
* *(xkappafit(iclegy,iclpro,icltar,ibk+1)
* -xkappafit(iclegy,iclpro,icltar,ibk))
xkappa2=xkappafit(iclegy+1,iclpro,icltar,ibk)
* +(b-bkbin*float(ibk-1))/bkbin
* *(xkappafit(iclegy+1,iclpro,icltar,ibk+1)
* -xkappafit(iclegy+1,iclpro,icltar,ibk))
xkappa=xkappa1+(xkappa2-xkappa1)/log(egyfac)
* *(log(ee)-log(egy0))
xkappa=facmc*xkappa
else
xkappa=xkappafit(iclegy,iclpro,icltar,ibk)
* +(b-bkbin*float(ibk-1))/bkbin
* *(xkappafit(iclegy,iclpro,icltar,ibk+1)
* -xkappafit(iclegy,iclpro,icltar,ibk))
xkappa=facmc*xkappa
endif
xkapZ=xkappa
endif
dalp=dble(alpD(i,iclpro,icltar)
* *cfalpro*cfaltar
* *gamV)
c if(i.eq.0)alp=alp
dalp=dalp
* *dble(xkappa)*dble(xkappa)
if(i.lt.2)then
dalp=dalp
* *exp(min(50d0,dble(rho*log(spp)-b2/delD(i,iclpro,icltar))))
* *dble(chad(iclpro)*chad(icltar)
* *gfactorp *gfactort)
dbet=dble(betDp(i,iclpro,icltar)
* +epsGp
* +gamb-alppar+1.)
dbetp=dble(betDpp(i,iclpro,icltar)
* +epsGt
* +gamb-alppar+1.)
else
absb=abs(b)-bmxdif(iclpro,icltar)
b2a=absb*absb
dalp=dalp
* *exp(min(50d0,dble((betD(i,iclpro,icltar)+epsG)*log(spp)
* -b2a/delD(i,iclpro,icltar))))
dbet=dble(betDp(i,iclpro,icltar)-alppar+1.+epsGp)
dbetp=dble(betDpp(i,iclpro,icltar)-alppar+1.+epsGt)
endif
coefgdp=utgam2(1.d0+dble(alplea(iclpro))+dbet)
coefgdt=utgam2(1.d0+dble(alplea(icltar))+dbetp)
dalp=dalp*utgam2(dbet)*utgam2(dbetp)/coefgdp/coefgdt
else
stop'GproPar: wrong m value. '
endif
alp=sngl(dalp)
bet=sngl(dbet)
betp=sngl(dbetp)
epsp=epsGp
epst=epsGt
epss=epsG
gamvv=gamV
alpUni(i,m)=dalp
betUni(i,m)=dbet
betpUni(i,m)=dbetp
if(ish.ge.10)write(ifch,*)' GfunPar :',alp,bet,betp,epsp,epst
& ,epss,gamvv
call utprjx('GfunPar ',ish,ishini,10)
end
c----------------------------------------------------------------------
subroutine GfomPar(b,spp)
c----------------------------------------------------------------------
c calculates parameters of the fom functions for pp (b-mode)
c----------------------------------------------------------------------
c b: impact param, spp: pp energy squared
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incpar'
include 'epos.incems'
call utprj('GfomPar ',ish,ishini,6)
ee=sqrt(spp)
rs=r2had(iclpro)+r2had(icltar)+slopom*log(spp)
bglaub2=4.*.0389*rs
if(iscreen.ne.0)then
absb=max(1.e-9,abs(b))
b2a=absb*absb
b2x=2.*epscrp*bglaub2
zzp=(epscrw*exp(-b2a/b2x))*fscra(ee/egyscr)
zzp=min(zzp,epscrx) !saturation
zzt=(epscrw*exp(-b2a/b2x))*fscra(ee/egyscr)
zzt=min(zzt,epscrx) !saturation
else
zzp=0.
zzt=0.
endif
z0=alpfomi!*epscrw*fscra(ee/egyscr)
if(z0.gt.0.)then
z1=zzp
zzpUni=dble(z1**gamfom/z0)*exp(-dble(b*b/delD(1,iclpro,icltar)))
c zzpUni=dble(4.*z0*(z1/z0)**1.5)
z1=zzt
zztUni=dble(z1**gamfom/z0)*exp(-dble(b*b/delD(1,iclpro,icltar)))
c zztUni=dble(4.*z0*(z1/z0)**1.5)
else
zzpUni=0d0
zztUni=0d0
endif
if(ish.ge.6)write(ifch,*)' GfomPar :',zzpUni,zztUni
call utprjx('GfomPar ',ish,ishini,6)
end
c----------------------------------------------------------------------
function fscra(x)
c----------------------------------------------------------------------
fscra=0
x2=x*x
if(x2.gt.1.)fscra=log(x2)!**2
end
c----------------------------------------------------------------------
function fscro(x,rho)
c----------------------------------------------------------------------
include 'epos.incpar'
fscro=znurho*rho
x2=x*x
c if(x2.gt.1.)fscro=sqrt(log(x2)**2+fscro**2)
if(x2.gt.1.)fscro=log(x2)*(1.+fscro)
end
c----------------------------------------------------------------------
function FbGlaub2(x)
c----------------------------------------------------------------------
c calculates (glauber radius)^2 from pp cross section (data fit)
c(estimated if not already calculated --> not any more to have smoother xs)
c----------------------------------------------------------------------
c x: pp energy
c----------------------------------------------------------------------
include 'epos.inc'
c if(sigine.eq.0.)then
if(iclpro+icltar.eq.3)then !pi+p
siginex=20.+0.08*log(x)**3.-0.004*log(x)**4.
elseif(iclpro+icltar.eq.5)then !K+p
siginex=16.+0.08*log(x)**3.-0.004*log(x)**4.
else
siginex=30.+0.095*log(x)**3.-0.004*log(x)**4.
endif
c else
c siginex=sigine
c endif
FbGlaub2=siginex/10./pi
return
end
c----------------------------------------------------------------------
subroutine recalcZPtn !???????????? not updated !!!
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incems'
stop'recalcZPtn not valid any more!!!!!!!'
c if(koll.eq.1.and.maproj.eq.1.and.matarg.eq.1)then
c npom=nprt(1)
c k=1
c ip=iproj(1)
c it=itarg(1)
c zparpro(k)=max(0,npom-1)*0.2
c zpartar(k)=0
c zpartar(k)=max(0,npom-1)*0.2
c ztav=zpartar(k)
c zpav=zparpro(k)
c zppevt=zpav
c zptevt=ztav
c endif
end
c----------------------------------------------------------------------
double precision function om1(xh,yp,b) !---test---
c----------------------------------------------------------------------
c om1 = G * C * gamd (C and gamd usually 1)
c xh - fraction of the energy squared s for the pomeron;
c b - impact parameter between the pomeron ends;
c yp - rapidity for the pomeron;
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incpar'
double precision Gf,xp,xm,xh,yp
Gf=0.d0
xp=sqrt(xh)*exp(yp)
xm=xh/xp
spp=engy**2
imax=idxD1
if(iomega.eq.2)imax=1
do i=idxDmin,imax
call Gfunpar(0.,0.,1,i,b,spp,alp,bet,betp,epsp,epst,epss,gamv)
Gf=Gf+dble(alp)*xp**dble(bet)*xm**dble(betp)
enddo
om1=Gf
* * dble(chad(iclpro)*chad(icltar))
end
c----------------------------------------------------------------------
double precision function om1intb(b) !---test---
c----------------------------------------------------------------------
c om1 integrated over xp and xm for given b
c Calculation by analytical integration
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incpar'
double precision cint,cint2,eps
parameter(eps=1.d-20)
spp=engy*engy
imax=idxD1
if(iomega.eq.2)imax=1
cint=0.d0
do i=idxDmin,imax
call Gfunpar(0.,0.,1,i,b,spp,alp,bet,betp,epsp,epst,epss,gamv)
cint2=dble(gamv*alp)
if((bet+1.0).gt.eps)then
cint2=cint2/dble(bet+1.0)
else
cint2=-cint2*log(xminDf)
endif
if((betp+1.0).gt.eps)then
cint2=cint2/dble(betp+1.0)
else
cint2=-cint2*log(xminDf)
endif
cint=cint+cint2
enddo
if(cint.lt.-eps)then
write(*,*) 'WARNING ! om1intb in epos-omg is <0 !!!!!'
write(*,*) 'WARNING ! => om1intb set to 1e-3 !!!!!'
write(*,*) 'WARNING ! => bmax=3.5 fm !!!!!'
cint=1.d-3
endif
om1intb=cint
* *dble(chad(iclpro)*chad(icltar))
return
end
c----------------------------------------------------------------------
double precision function om1intbk(k) !---MC---
c----------------------------------------------------------------------
c Diffractive part of om1 integrated over xp and xm for given pair k
c Calculation by analytical integration
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incems'
include 'epos.incpar'
double precision cint,cint2,eps,bet,betp
parameter(eps=1.d-20)
imax=idxD1
if(iomega.eq.2)imax=1
om1intbk=0d0
cint=0
do i=idxDmin,imax
bet=btildep(i,k)
betp=btildepp(i,k)
cint2=atilde(i,k)
if((bet+1.d0).gt.eps)then
cint2=cint2/(bet+1.d0)
else
cint2=-cint2*log(xminDf)
endif
if((betp+1.d0).gt.eps)then
cint2=cint2/(betp+1.d0)
else
cint2=-cint2*log(xminDf)
endif
cint=cint+cint2
enddo
if(cint.lt.-eps)then
write(*,*) 'WARNING ! om1intbk in epos-omg is <0 !!!!!'
write(*,*) 'WARNING ! => om1intbk set to 1e-3 !!!!!'
write(*,*) 'WARNING ! => bmax=3.5 fm !!!!!'
cint=1.d-3
endif
om1intbk=cint
* *dble(chad(iclpro)*chad(icltar))
return
end
c----------------------------------------------------------------------
double precision function om1intbi(b,iqq) !---MC---
c----------------------------------------------------------------------
c om1 integrated over xp and xm for given b
c Calculation by analytical integration of contribution iqq
c----------------------------------------------------------------------
include 'epos.inc'
include 'epos.incsem'
include 'epos.incpar'
double precision eps,cint
parameter(eps=1.d-20)
spp=engy*engy
call Gfunpar(0.,0.,1,iqq,b,spp,alp,bet,betp,epsp,epst,epss,gamv)
cint=dble(gamv*alp)
if(dble(bet+1.0).gt.eps)then
cint=cint/dble(bet+1.0)
else
cint=-cint*log(xminDf)
endif
if(dble(betp+1.0).gt.eps)then
cint=cint/dble(betp+1.0)
else
cint=-cint*log(xminDf)
endif
if(cint.lt.-eps)then