-
Notifications
You must be signed in to change notification settings - Fork 3
/
medium-simple.f
818 lines (727 loc) · 23.8 KB
/
medium-simple.f
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
C+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
C++ Copyright (C) 2017 Korinna C. Zapp [[email protected]] ++
C++ ++
C++ This file is part of JEWEL 2.2.0 ++
C++ ++
C++ The JEWEL homepage is jewel.hepforge.org ++
C++ ++
C++ The medium model was partly implemented by Jochen Klein. ++
C++ Raghav Kunnawalkam Elayavalli helped with the implementation ++
C++ of the V+jet processes. ++
C++ ++
C++ Please follow the MCnet GUIDELINES and cite Eur.Phys.J. C74 ++
C++ (2014) no.2, 2762 [arXiv:1311.0048] for the code and ++
C++ JHEP 1303 (2013) 080 [arXiv:1212.1599] and ++
C++ optionally EPJC 60 (2009) 617 [arXiv:0804.3568] for the ++
C++ physics. The reference for V+jet processes is EPJC 76 (2016) ++
C++ no.12 695 [arXiv:1608.03099] and for recoil effects it is ++
C++ arXiv:1707.01539.
C++ ++
C++ JEWEL relies heavily on PYTHIA 6 for the event generation. The ++
C++ modified version of PYTHIA 6.4.25 that is distributed with ++
C++ JEWEL is, however, not an official PYTHIA release and must not ++
C++ be used for anything else. Please refer to results as ++
C++ "JEWEL+PYTHIA". ++
C++ ++
C++ JEWEL also uses code provided by S. Zhang and J. M. Jing ++
C++ (Computation of Special Functions, John Wiley & Sons, New York, ++
C++ 1996 and http://jin.ece.illinois.edu) for computing the ++
C++ exponential integral Ei(x). ++
C++ ++
C++ ++
C++ JEWEL is free software; you can redistribute it and/or ++
C++ modify it under the terms of the GNU General Public License ++
C++ as published by the Free Software Foundation; either version 2 ++
C++ of the License, or (at your option) any later version. ++
C++ ++
C++ JEWEL is distributed in the hope that it will be useful, ++
C++ but WITHOUT ANY WARRANTY; without even the implied warranty of ++
C++ MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the ++
C++ GNU General Public License for more details. ++
C++ ++
C++ You should have received a copy of the GNU General Public ++
C++ License along with this program; if not, write to the Free ++
C++ Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, ++
C++ MA 02110-1301 USA ++
C++ ++
C++ Linking JEWEL statically or dynamically with other modules is ++
C++ making a combined work based on JEWEL. Thus, the terms and ++
C++ conditions of the GNU General Public License cover the whole ++
C++ combination. ++
C++ ++
C++ In addition, as a special exception, I give you permission to ++
C++ combine JEWEL with the code for the computation of special ++
C++ functions provided by S. Zhang and J. M. Jing. You may copy and ++
C++ distribute such a system following the terms of the GNU GPL for ++
C++ JEWEL and the licenses of the other code concerned, provided ++
C++ that you include the source code of that other code when and as ++
C++ the GNU GPL requires distribution of source code. ++
C+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
SUBROUTINE MEDINIT(FILE,id,etam,mass)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--longitudinal boost of momentum distribution
common/boostmed/boost
logical boost
C--factor to vary Debye mass
COMMON/MDFAC/MDFACTOR,MDSCALEFAC
DOUBLE PRECISION MDFACTOR,MDSCALEFAC
C--nuclear thickness function
COMMON /THICKFNC/ RMAX,TA(100,2)
DOUBLE PRECISION RMAX,TA
C--geometrical cross section
COMMON /CROSSSEC/ IMPMAX,CROSS(200,3)
DOUBLE PRECISION IMPMAX,CROSS
C--identifier of log file
common/logfile/logfid
integer logfid
DATA RAU/10./
DATA D3/0.9d0/
DATA ZETA3/1.2d0/
C--local variables
INTEGER I,LUN,POS,IOS,id,mass
double precision etam
CHARACTER*100 BUFFER,LABEL,tempbuf
CHARACTER*80 FILE
character firstchar
logical fileexist
etamax2 = etam
logfid = id
IOS=0
LUN=77
C--default settings
TAUI=0.6d0
TI=0.36d0
TC=0.17d0
WOODSSAXON=.TRUE.
CENTRMIN=0.d0
CENTRMAX=10.d0
NF=3
A=mass
N0=0.17d0
D=0.54d0
SIGMANN=6.2
MDFACTOR=0.45d0
MDSCALEFAC=0.9d0
boost = .true.
C--read settings from file
write(logfid,*)
inquire(file=FILE,exist=fileexist)
if(fileexist)then
write(logfid,*)'Reading medium parameters from ',FILE
OPEN(unit=LUN,file=FILE,status='old',err=10)
do 20 i=1,1000
READ(LUN, '(A)', iostat=ios) BUFFER
if (ios.ne.0) goto 30
firstchar = buffer(1:1)
if (firstchar.eq.'#') goto 20
POS=SCAN(BUFFER,' ')
LABEL=BUFFER(1:POS)
BUFFER=BUFFER(POS+1:)
IF (LABEL=="TAUI")THEN
READ(BUFFER,*,IOSTAT=IOS) TAUI
ELSE IF (LABEL=="TI") THEN
READ(BUFFER,*,IOSTAT=IOS) TI
ELSE IF (LABEL=="TC") THEN
READ(BUFFER,*,IOSTAT=IOS) TC
ELSE IF (LABEL=="WOODSSAXON") THEN
READ(BUFFER,*,IOSTAT=IOS) WOODSSAXON
ELSE IF (LABEL=="CENTRMIN") THEN
READ(BUFFER,*,IOSTAT=IOS) CENTRMIN
ELSE IF (LABEL=="CENTRMAX") THEN
READ(BUFFER,*,IOSTAT=IOS) CENTRMAX
ELSE IF (LABEL=="NF") THEN
READ(BUFFER,*,IOSTAT=IOS) NF
ELSE IF (LABEL=="N0") THEN
READ(BUFFER,*,IOSTAT=IOS) N0
ELSE IF (LABEL=="D") THEN
READ(BUFFER,*,IOSTAT=IOS) D
ELSE IF (LABEL=="SIGMANN") THEN
READ(BUFFER,*,IOSTAT=IOS) SIGMANN
ELSE IF (LABEL=="MDFACTOR") THEN
READ(BUFFER,*,IOSTAT=IOS) MDFACTOR
ELSE IF (LABEL=="MDSCALEFAC") THEN
READ(BUFFER,*,IOSTAT=IOS) MDSCALEFAC
else
write(logfid,*)'unknown label ',label
endif
20 continue
30 close(LUN,status='keep')
write(logfid,*)'...done'
goto 40
10 write(logfid,*)'Could not open medium parameter file, '//
& 'will run with default settings.'
else
write(logfid,*)'No medium parameter file found, '//
& 'will run with default settings.'
endif
40 write(logfid,*)'using parameters:'
write(logfid,*)'TAUI =',TAUI
write(logfid,*)'TI =',TI
write(logfid,*)'TC =',TC
write(logfid,*)'WOODSSAXON =',WOODSSAXON
write(logfid,*)'CENTRMIN =',CENTRMIN
write(logfid,*)'CENTRMAX =',CENTRMAX
write(logfid,*)'NF =',NF
write(logfid,*)'A =',A
write(logfid,*)'N0 =',N0
write(logfid,*)'D =',D
write(logfid,*)'SIGMANN =',SIGMANN
write(logfid,*)'MDFACTOR =',MDFACTOR
write(logfid,*)'MDSCALEFAC =',MDSCALEFAC
write(logfid,*)
write(logfid,*)
write(logfid,*)
C--calculate T_A(x,y)
CALL CALCTA
C--calculate geometrical cross section
CALL CALCXSECTION
END
SUBROUTINE MEDNEXTEVT
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--geometrical cross section
COMMON /CROSSSEC/ IMPMAX,CROSS(200,3)
DOUBLE PRECISION IMPMAX,CROSS
C--local variables
integer i,j
DOUBLE PRECISION PYR,R,b1,b2,gettemp
C--pick an impact parameter
r=(pyr(0)*(centrmax-centrmin)+centrmin)/100.
i=0
do 130 j=1,200
if ((r-cross(j,3)/cross(200,3)).ge.0.) then
i=i+1
else
goto 132
endif
130 continue
132 continue
b1 = (i-1)*0.1d0
b2 = i*0.1d0
breal = (b2*(cross(i,3)/cross(200,3)-r)
& +b1*(r-cross(i+1,3)/cross(200,3)))/
& (cross(i,3)/cross(200,3)-cross(i+1,3)/cross(200,3))
centr = r;
END
double precision function getcentrality()
implicit none
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
getcentrality=centr
end
SUBROUTINE PICKVTX(X,Y)
IMPLICIT NONE
DOUBLE PRECISION X,Y
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
C--local variables
DOUBLE PRECISION X1,X2,Y1,Y2,Z,XVAL,YVAL,ZVAL,NTHICK,PYR
X1=BREAL/2.-RAU
X2=RAU-BREAL/2.
Y1=-SQRT(4*RAU**2-BREAL**2)/2.
Y2=SQRT(4*RAU**2-BREAL**2)/2.
131 XVAL=PYR(0)*(X2-X1)+X1
YVAL=PYR(0)*(Y2-Y1)+Y1
IF((NTHICK(XVAL-BREAL/2.,YVAL).EQ.0.d0).OR.
& NTHICK(XVAL+BREAL/2.,YVAL).EQ.0.d0) GOTO 131
ZVAL=PYR(0)*NTHICK(-BREAL/2.,0d0)*NTHICK(BREAL/2.,0d0)
Z=NTHICK(XVAL-BREAL/2.,YVAL)*NTHICK(XVAL+BREAL/2.,YVAL)
IF(ZVAL.GT.Z) GOTO 131
X=XVAL
Y=YVAL
END
SUBROUTINE SETB(BVAL)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
DOUBLE PRECISION BVAL
BREAL=BVAL
END
SUBROUTINE GETSCATTERER(X,Y,Z,T,TYPE,PX,PY,PZ,E,MS)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
C--internal medium parameters
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--longitudinal boost of momentum distribution
common/boostmed/boost
logical boost
C--function calls
DOUBLE PRECISION GETTEMP,GETMD,GETMOM,GETMS
C--identifier of log file
common/logfile/logfid
integer logfid
C--local variables
DOUBLE PRECISION X,Y,Z,T,MS,PX,PY,PZ,E,MD,TEMP
INTEGER TYPE
DOUBLE PRECISION R,PYR,pmax,wt,tau,theta,phi,pi,p,ys,pz2,e2
DATA PI/3.141592653589793d0/
R=PYR(0)
IF(R.LT.(2.*12.*NF*D3/3.)/(2.*12.*NF*D3/3.+3.*16.*ZETA3/2.))THEN
TYPE=2
ELSE
TYPE=21
ENDIF
MS=GETMS(X,Y,Z,T)
MD=GETMD(X,Y,Z,T)
TEMP=GETTEMP(X,Y,Z,T)
tau=sqrt(t**2-z**2)
if (boost) then
ys = 0.5*log((t+z)/(t-z))
else
ys = 0.d0
endif
pmax = 10.*temp
IF(TEMP.LT.1.D-2)THEN
write(logfid,*)'asking for a scattering centre without medium:'
write(logfid,*)'at (x,y,z,t)=',X,Y,Z,T
write(logfid,*)'making one up to continue but '//
& 'something is wrong!'
TYPE=21
PX=0.d0
PY=0.d0
PZ=0.d0
MS=GETMS(0.d0,0.d0,0.d0,0.d0)
MD=GETMD(0.d0,0.d0,0.d0,0.d0)
E=SQRT(PX**2+PY**2+PZ**2+MS**2)
RETURN
ENDIF
10 p = pyr(0)**0.3333333*pmax
E2 = sqrt(p**2+ms**2)
if (type.eq.2) then
wt = (exp(ms/temp)-1.)/(exp(E2/temp)-1.)
else
wt = (exp(ms/temp)+1.)/(exp(E2/temp)+1.)
endif
if (wt.gt.1.) write(logfid,*)'Error in getscatterer: weight = ',wt
if (wt.lt.0.) write(logfid,*)'Error in getscatterer: weight = ',wt
if (pyr(0).gt.wt) goto 10
phi = pyr(0)*2.*pi
theta = -acos(2.*pyr(0)-1.)+pi
px = p*sin(theta)*cos(phi)
py = p*sin(theta)*sin(phi)
pz2 = p*cos(theta)
E = cosh(ys)*E2 + sinh(ys)*pz2
pz = sinh(ys)*E2 + cosh(ys)*pz2
END
SUBROUTINE AVSCATCEN(X,Y,Z,T,PX,PY,PZ,E,m)
IMPLICIT NONE
C--longitudinal boost of momentum distribution
common/boostmed/boost
logical boost
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--local variables
double precision x,y,z,t,px,py,pz,e,getms,m,ys
if (boost) then
ys = 0.5*log((t+z)/(t-z))
if ((z.eq.0.d0).and.(t.eq.0.d0)) ys =0.d0
if (ys.gt.etamax2) ys=etamax2
if (ys.lt.-etamax2) ys=-etamax2
else
ys = 0.d0
endif
m = getms(x,y,z,t)
e = m*cosh(ys)
px = 0.d0
py = 0.d0
pz = m*sinh(ys)
end
SUBROUTINE maxscatcen(PX,PY,PZ,E,m)
IMPLICIT NONE
C--longitudinal boost of momentum distribution
common/boostmed/boost
logical boost
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--local variables
double precision px,py,pz,e,getmsmax,m,ys
if (boost) then
ys = etamax2
else
ys = 0.d0
endif
m = getmsmax()
e = m*cosh(ys)
px = 0.d0
py = 0.d0
pz = m*sinh(ys)
end
DOUBLE PRECISION FUNCTION GETMD(X1,Y1,Z1,T1)
IMPLICIT NONE
C--factor to vary Debye mass
COMMON/MDFAC/MDFACTOR,MDSCALEFAC
DOUBLE PRECISION MDFACTOR,MDSCALEFAC
DOUBLE PRECISION X1,Y1,Z1,T1,GETTEMP
GETMD=MDSCALEFAC*3.*GETTEMP(X1,Y1,Z1,T1)
GETMD=MAX(GETMD,MDFACTOR)
END
DOUBLE PRECISION FUNCTION GETMS(X2,Y2,Z2,T2)
IMPLICIT NONE
DOUBLE PRECISION X2,Y2,Z2,T2,GETMD
GETMS=GETMD(X2,Y2,Z2,T2)/SQRT(2.)
END
DOUBLE PRECISION FUNCTION GETNEFF(X3,Y3,Z3,T3)
IMPLICIT NONE
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C-- local variables
DOUBLE PRECISION X3,Y3,Z3,T3,PI,GETTEMP,tau,cosheta
DATA PI/3.141592653589793d0/
tau = sqrt(t3**2-z3**2)
cosheta = t3/tau
GETNEFF=(2.*6.*NF*D3*2./3. + 16.*ZETA3*3./2.)
& *GETTEMP(X3,Y3,Z3,T3)**3/PI**2
getneff = getneff/cosheta
END
DOUBLE PRECISION FUNCTION GETTEMP(X4,Y4,Z4,T4)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--local variables
DOUBLE PRECISION X4,Y4,Z4,T4,TAU,NPART,EPS0,EPSIN,TEMPIN,PI,
&NTHICK,ys
DATA PI/3.141592653589793d0/
GETTEMP=0.D0
IF(ABS(Z4).GT.T4)RETURN
TAU=SQRT(T4**2-Z4**2)
C--check for overlap region
IF((NTHICK(X4-BREAL/2.,Y4).EQ.0.d0).OR.
&NTHICK(X4+BREAL/2.,Y4).EQ.0.d0) RETURN
ys = 0.5*log((t4+z4)/(t4-z4))
if (abs(ys).gt.etamax2) return
C--determine initial temperature at transverse position
IF(WOODSSAXON)THEN
EPS0=(16.*8.+7.*2.*6.*NF)*PI**2*TI**4/240.
EPSIN=EPS0*NPART(X4-BREAL/2.,Y4,X4+BREAL/2.,Y4)
& *PI*RAU**2/(2.*A)
TEMPIN=(EPSIN*240./(PI**2*(16.*8.+7.*2.*6.*NF)))**0.25
ELSE
TEMPIN=TI
ENDIF
C--calculate temperature if before initial time
IF(TAU.LE.TAUI)THEN
GETTEMP=TEMPIN*TAU/TAUI
ELSE
C--evolve temperature
GETTEMP=TEMPIN*(TAUI/TAU)**0.3333
ENDIF
IF(GETTEMP.LT.TC) GETTEMP=0.d0
END
DOUBLE PRECISION FUNCTION GETTEMPMAX()
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--function call
DOUBLE PRECISION GETTEMP
GETTEMPMAX=GETTEMP(0.D0,0.D0,0.D0,TAUI)
END
DOUBLE PRECISION FUNCTION GETMDMAX()
IMPLICIT NONE
C--factor to vary Debye mass
COMMON/MDFAC/MDFACTOR,MDSCALEFAC
DOUBLE PRECISION MDFACTOR,MDSCALEFAC
DOUBLE PRECISION GETTEMPMAX
GETMDMAX=MDSCALEFAC*3.*GETTEMPMAX()
GETMDMAX=MAX(GETMDMAX,MDFACTOR)
END
DOUBLE PRECISION FUNCTION GETMDMIN()
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--factor to vary Debye mass
COMMON/MDFAC/MDFACTOR,MDSCALEFAC
DOUBLE PRECISION MDFACTOR,MDSCALEFAC
DOUBLE PRECISION GETTEMPMAX
GETMDMIN=MDSCALEFAC*3.*TC
GETMDMIN=MAX(GETMDMIN,MDFACTOR)
END
DOUBLE PRECISION FUNCTION GETMSMAX()
IMPLICIT NONE
DOUBLE PRECISION GETMDMAX,SQRT
GETMSMAX=GETMDMAX()/SQRT(2.D0)
END
DOUBLE PRECISION FUNCTION GETNATMDMIN()
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--factor to vary Debye mass
COMMON/MDFAC/MDFACTOR,MDSCALEFAC
DOUBLE PRECISION MDFACTOR,MDSCALEFAC,PI
DATA PI/3.141592653589793d0/
C--local variables
DOUBLE PRECISION T,GETMDMIN
T=GETMDMIN()/(MDSCALEFAC*3.)
GETNATMDMIN=(2.*6.*NF*D3*2./3. + 16.*ZETA3*3./2.)
& *T**3/PI**2
END
DOUBLE PRECISION FUNCTION GETLTIMEMAX()
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--function call
DOUBLE PRECISION GETTEMPMAX
GETLTIMEMAX=TAUI*(GETTEMPMAX()/TC)**3*cosh(etamax2)
END
DOUBLE PRECISION FUNCTION GETNEFFMAX()
IMPLICIT NONE
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C-- local variables
DOUBLE PRECISION PI,GETTEMPMAX
DATA PI/3.141592653589793d0/
GETNEFFMAX=(2.*6.*NF*D3*2./3. + 16.*ZETA3*3./2.)
& *GETTEMPMAX()**3/PI**2
END
DOUBLE PRECISION FUNCTION NPART(XX1,YY1,XX2,YY2)
IMPLICIT NONE
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--local variables
DOUBLE PRECISION XX1,YY1,XX2,YY2,NTHICK
NPART = NTHICK(XX1,YY1)*(1.-EXP(-SIGMANN*NTHICK(XX2,YY2))) +
& NTHICK(XX2,YY2)*(1.-EXP(-SIGMANN*NTHICK(XX1,YY1)))
END
DOUBLE PRECISION FUNCTION NTHICK(X1,Y1)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--identifier of log file
common/logfile/logfid
integer logfid
C--nuclear thickness function
COMMON /THICKFNC/ RMAX,TA(100,2)
DOUBLE PRECISION RMAX,TA
INTEGER LINE,LMIN,LMAX,I
DOUBLE PRECISION X1,Y1,XA(4),YA(4),Y,DY,R,C,B,DELTA
R=SQRT(X1**2+Y1**2)
IF(R.GT.TA(100,1))THEN
NTHICK=0.
ELSE
LINE=INT(R*99.d0/TA(100,1)+1)
LMIN=MAX(LINE,1)
LMIN=MIN(LMIN,99)
IF((R.LT.TA(LMIN,1)).OR.(R.GT.TA(LMIN+1,1)))
& write(logfid,*)LINE,LMIN,R,TA(LMIN,1),TA(LMIN+1,1)
XA(1)=TA(LMIN,1)
XA(2)=TA(LMIN+1,1)
YA(1)=TA(LMIN,2)
YA(2)=TA(LMIN+1,2)
C=(YA(2)-YA(1))/(XA(2)-XA(1))
B=YA(1)-C*XA(1)
NTHICK=C*R+B
ENDIF
END
SUBROUTINE CALCTA()
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C-- nuclear thickness function
COMMON /THICKFNC/ RMAX,TA(100,2)
DOUBLE PRECISION RMAX,TA
C--variables for integration
COMMON/INTEG/B,R
DOUBLE PRECISION B,R
C--local variables
INTEGER NSTEPS,I
DOUBLE PRECISION EPS,HFIRST,Y
NSTEPS=100
EPS=1.E-4
HFIRST=0.1D0
R=1.12*A**(0.33333)-0.86*A**(-0.33333)
RMAX=2.*R
DO 10 I=1,NSTEPS
C--set transverse position
B=(I-1)*2.D0*R/NSTEPS
Y=0.D0
C--integrate along longitudinal line
CALL ODEINT(Y,-2*R,2*R,EPS,HFIRST,0.d0,101)
TA(I,1)=B
TA(I,2)=Y
10 CONTINUE
END
SUBROUTINE CALCXSECTION()
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,NF
INTEGER NF
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C-- geometrical cross section
COMMON /CROSSSEC/ IMPMAX,CROSS(200,3)
DOUBLE PRECISION IMPMAX,CROSS
C--local variables
INTEGER IX,IY,IB
DOUBLE PRECISION B,P,PROD,X,Y,NTHICK,NPART,pprev
pprev=0.
DO 30 IB=1,200
B=0.1d0*IB
PROD=1.d0
DO 10 IX=1,100
DO 20 IY=1,100
X=-20.d0+IX*0.4d0
Y=-20.d0+IY*0.4d0
PROD=PROD*
&EXP(-NTHICK(X+B/2.D0,Y)*SIGMANN)**(0.16d0*NTHICK(X-B/2.D0,Y))
20 CONTINUE
10 CONTINUE
P=(1.D0-PROD)*8.8D0/14.D0*B
CROSS(IB,1)=B
CROSS(IB,2)=P
if (ib.eq.1) then
cross(ib,3)=0.
else
cross(ib,3)=cross(ib-1,3)+(p+pprev)/2.*0.1
endif
pprev=p
30 CONTINUE
IMPMAX=19.95
END
DOUBLE PRECISION FUNCTION MEDDERIV(XVAL,W)
IMPLICIT NONE
DOUBLE PRECISION XVAL
INTEGER W
C--medium parameters
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON
C--variables for integration
COMMON/INTEG/B,R
DOUBLE PRECISION B,R
IF (W.EQ.1) THEN
C--XVAL corresponds to z-coordinate
MEDDERIV=N0/(1+EXP((SQRT(B**2+XVAL**2)-R)/D))
ELSE
MEDDERIV=0.D0
ENDIF
END