-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathmedium-hydro.f
1644 lines (1434 loc) · 54.1 KB
/
medium-hydro.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
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++ JEWEL Add-On to Read v-USPhydro 2+1D Profiles ++
C++ ++
C++ This file is an unofficial modification of the original ++
C++ JEWEL program (copyright header below). The original JEWEL team ++
C++ does not hold any responsibility regarding it. ++
C++ ++
C++ The program is part of the developed interface between the ++
C++ parton propagation of JEWEL and an external hydrodynamic 2+1D ++
C++ medium profile, intended for v-USPhydro. ++
C++ ++
C++ ++
C++ Created by: ++
C++ - Fabio M. Canedo [fabio.canedo@usp.br] ++
C++ - Leonardo Barreto [leonardo.barreto.campos@usp.br] ++
C++ Instituto de Fisica, Universidade de Sao Paulo, Brazil ++
C++ 2019 ++
C+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
C+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
C++ Copyright (C) 2017 Korinna C. Zapp [Korinna.Zapp@cern.ch] ++
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,
& NX,NY,NT,NF,DX,DT,XMAX,XMIN,TMAX
INTEGER NX,NY,NT,NF
DOUBLE PRECISION DX,DT,XMAX,XMIN,TMAX
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
common/grid/timesteps(60),tprofile(834,834,60),vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
common/gridvel/ux(834,834,60),uy(834,834,60)
double precision ux,uy
COMMON/MEDFILEC/MEDFILE,NLIST,endoff
CHARACTER*200 MEDFILE
INTEGER NLIST
logical endoff
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON,MODMED,MEDFILELIST
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED,MEDFILELIST
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
common/temperature/tempfac
double precision tempfac
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--hydrodynamic quantities
COMMON /HYDROLIM/ MIDRAPLIM, TMAXLIM, TVELMAXLIM, BOOSTTR,
&GLOBALLIMS, PRETAUHYDRO
DOUBLE PRECISION MIDRAPLIM, TMAXLIM, TVELMAXLIM
LOGICAL BOOSTTR, GLOBALLIMS, PRETAUHYDRO
C--hydro auxiliary files
COMMON /HYDROF/ INITVTXF, IDEALN0, CUSTOMN0F
CHARACTER*200 INITVTXF, CUSTOMN0F
LOGICAL IDEALN0
C--identifier of log file
common/logfile/logfid
integer logfid
C--grid parameters
common/gridpar/ gdt,gdx,gxmax,gxmin,gnx,gny,gnt
double precision gdt,gdx,gxmax,gxmin
integer gnx,gny,gnt
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*100 FILE
character firstchar
logical fileexist
! Use MIDRAPLIM for limits and
! etamax2 for the simulation
etamax2 = etam
logfid = id
IOS=0
LUN=77
NLIST=1
endoff=.false.
C--default settings
TAUI=0.6d0 ! Not used for hydro
TI=0.36d0 ! Not used for hydro
TC=0.17d0
WOODSSAXON=.TRUE. ! Not used for hydro
CENTRMIN=0.d0 ! Not used for hydro
CENTRMAX=10.d0 ! Not used for hydro
NF=3
A=mass
N0=0.17d0
D=0.54d0
SIGMANN=6.2
MDFACTOR=0.45d0
MDSCALEFAC=0.9d0
tempfac=1.0d0
boost =.true.
breal=0.d0
C--hydro settings
MODMED=.TRUE. ! Deprecated
MEDFILELIST=.FALSE. ! Deprecated
GLOBALLIMS = .true. ! Global limits vs JEWEL default impl.
! Extrapolation of soft mid rap lim, +0.1 then tuning
! (see nucl-ex/1612.08966)
MIDRAPLIM=3.3d0
TMAXLIM=0.675d0 ! Approx vUSP+TRENTo temp lim (PbPb 5TeV)
BOOSTTR=.true. ! Logical for transverse velocity (F => u=0)
TVELMAXLIM=0.927d0 ! Approx vUSP+TRENTo uT lim (PbPb 5TeV)
! Temperature assumption before TAU0, i.e. F => T(tau < tau0) = 0
PRETAUHYDRO=.false.
!GRIDN=834 ! Number of points in grid (per dimension)
INITVTXF='initvertexmap.dat' ! Initial vertex map file
! Logical for using n0 propto T ** 3, ignores CUSTOMN0F if true
IDEALN0=.false.
CUSTOMN0F='n0proptoentropy_derv21_Tmax675.dat' ! Custom n0 function file
! Print program info
call printhydroheader(logfid)
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=="MODMED") THEN
READ(BUFFER,*,IOSTAT=IOS) MODMED
ELSE IF (LABEL=="MEDFILE") THEN
READ(BUFFER,'(50A)',IOSTAT=IOS) MEDFILE
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 IF (LABEL=="GLOBALLIMS") THEN
READ(BUFFER,*,IOSTAT=IOS) GLOBALLIMS
ELSE IF (LABEL=="MIDRAPLIM") THEN
READ(BUFFER,*,IOSTAT=IOS) MIDRAPLIM
ELSE IF (LABEL=="TMAXLIM") THEN
READ(BUFFER,*,IOSTAT=IOS) TMAXLIM
ELSE IF (LABEL=="BOOSTTR") THEN
READ(BUFFER,*,IOSTAT=IOS) BOOSTTR
ELSE IF (LABEL=="TVELMAXLIM") THEN
READ(BUFFER,*,IOSTAT=IOS) TVELMAXLIM
ELSE IF (LABEL=="PRETAUHYDRO") THEN
READ(BUFFER,*,IOSTAT=IOS) PRETAUHYDRO
!ELSE IF (LABEL=="GRIDN") THEN
!READ(BUFFER,*,IOSTAT=IOS) GRIDN
ELSE IF (LABEL=="INITVTXF") THEN
READ(BUFFER,*,IOSTAT=IOS) INITVTXF
ELSE IF (LABEL=="IDEALN0") THEN
READ(BUFFER,*,IOSTAT=IOS) IDEALN0
ELSE IF (LABEL=="CUSTOMN0F") THEN
READ(BUFFER,*,IOSTAT=IOS) CUSTOMN0F
! JEWEL original boost (longitudinal) as user option
ELSE IF (LABEL=="BOOSTZ") THEN
READ(BUFFER,*,IOSTAT=IOS) boost
ELSE IF (LABEL=="BREAL") THEN
READ(BUFFER,*,IOSTAT=IOS) breal
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,*)'MODMED = ',MODMED
write(logfid,*)'MEDFILELIST = ',MEDFILELIST
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,*)'BREAL = ',breal
write(logfid,*)'GLOBALLIMS = ',GLOBALLIMS
write(logfid,*)'MIDRAPLIM = ',MIDRAPLIM
write(logfid,*)'TMAXLIM = ',TMAXLIM
write(logfid,*)'BOOSTTR = ',BOOSTTR
write(logfid,*)'TVELMAXLIM = ',TVELMAXLIM
write(logfid,*)'PRETAUHYDRO = ',PRETAUHYDRO
!write(logfid,*)'GRIDN = ',GRIDN
write(logfid,*)'INITVTXF = ',INITVTXF
write(logfid,*)'IDEALN0 = ',IDEALN0
write(logfid,*)'CUSTOMN0F = ',CUSTOMN0F
write(logfid,*)'BOOSTZ = ',boost
write(logfid,*)
write(logfid,*)
write(logfid,*)
if (.not. boost) then
write(logfid,*) 'No longitudinal boost => MIDRAPLIM = 0'
MIDRAPLIM = 0.0
else if (MIDRAPLIM.lt.etamax2) then
MIDRAPLIM = etamax2
write(logfid,*) 'ETAMAX > MIDRAPLIM'
write(logfid,*) 'Extrapolating mid-rapidity limit'
write(logfid,*)
end if
if (.not. BOOSTTR) then
TVELMAXLIM = 0.d0
write(logfid,*) 'No transverse u => TVELMAXLIM = 0'
write(logfid,*)
end if
C--Call the modified medium setup
CALL MYMED()
write(logfid,*) 'Hydrodynamic profile loaded: ', MEDFILE
write(logfid,*)
END
SUBROUTINE MEDNEXTEVT
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,
& NX,NY,NT,NF,DX,DT,XMAX,XMIN,TMAX
INTEGER NX,NY,NT,NF
DOUBLE PRECISION DX,DT,XMAX,XMIN,TMAX
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
common/grid/timesteps(60),tprofile(834,834,60),vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
common/gridvel/ux(834,834,60),uy(834,834,60)
double precision ux,uy
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON,MODMED
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED
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
! Dummy function
END
SUBROUTINE PICKVTX(x, y)
IMPLICIT NONE
DOUBLE PRECISION x, y
C--medium parameters
common/grid/timesteps(60), tprofile(834,834,60),vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
integer tries, maxpos, minpos, irand, jrand, counter
double precision pyr, zval
common/logfile/logfid
integer logfid
do counter = 1, 5
zval = pyr(0) * maxval(vtxmap)
! Consider only a 10 fm by 10 fm (k = 250, 583)
! since anywhere else wont have entropy (PbPb)
maxpos = 583
minpos = 250
do tries = 1, 1000000
irand = int(pyr(0) * (maxpos - minpos) + minpos)
jrand = int(pyr(0) * (maxpos - minpos) + minpos)
if (vtxmap(irand, jrand) .gt. zval) then
x = -25.d0 + (irand - 1) * (50.d0 / 833.d0)
y = -25.d0 + (jrand - 1) * (50.d0 / 833.d0)
return
end if
end do
write(*,*) 'Failed to find vtx, restarting selection.'
write(*,*) 'This should not happen often.'
end do
! Only do process 5 times, kill simulation otherwise
write(*,*) 'No initial vertex found. Check initialvtx table.'
write(logfid,*) 'No initial vertex found. Check initialvtx table.'
stop
END SUBROUTINE
SUBROUTINE GETSCATTERER(X,Y,Z,T,TYPE,PX,PY,PZ,E,MS)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,
& NX,NY,NT,NF,DX,DT,XMAX,XMIN,TMAX
INTEGER NX,NY,NT,NF
DOUBLE PRECISION DX,DT,XMAX,XMIN,TMAX
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
common/grid/timesteps(60),tprofile(834,834,60),vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
C--internal medium parameters
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON,MODMED
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED
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
double precision u,ux,uy,px2,py2,px3,py3,e3
double precision getu,getutheta,eta,vel
INTEGER TYPE
DOUBLE PRECISION R,PYR,pmax,wt,tau,theta,phi,pi,p,ys,pz2,e2
double precision e4
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
! JEWEL original
!MS=GETMS(X,Y,Z,T)
!MD=GETMD(X,Y,Z,T)
MD=GETMD(X,Y,Z,T)
MS = MD / sqrt(2.)
TEMP=GETTEMP(X,Y,Z,T)
tau=sqrt(t**2-z**2)
if (boost) then
ys = 0.5d0*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
! JEWEL original
!MS=GETMS(0.d0,0.d0,0.d0,0.d0)
!MD=GETMD(0.d0,0.d0,0.d0,0.d0)
MD=GETMD(0.d0,0.d0,0.d0,0.d0)
MS=MD / sqrt(2.)
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)
!JEWEL original longitudinal boost, note that boost is using -ys
!E = cosh(ys)*E2 + sinh(ys)*pz2
!pz = sinh(ys)*E2 + cosh(ys)*pz2
!Perform boost
E = E2
pz = pz2
call LorentzBoostToLabFrame(E,px,py,pz,x,y,z,t)
IF (E.lt.0.d0) THEN
write(logfid,*) 'Negative energy (', E, ') in GETSCATTERER'
END IF
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
double precision getu,getutheta,theta,u,gettemp
C--local variables
double precision x,y,z,t,px,py,pz,e,getms,m,ys,temp
!Original JEWEL implementation
!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)
m = getms(x,y,z,t)
e = m
px = 0.d0
py = 0.d0
pz = 0.d0
!Apply boost
call LorentzBoostToLabFrame(e,px,py,pz,x,y,z,t)
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
double precision getu,getutheta,eta,u,theta
COMMON/GAMMAMAX/GAMMAMAXIMUM,VELMAXIMUM,DENSITYMAXIMUM,
&DENSITYMINIMUM, TAUMIN
DOUBLE PRECISION GAMMAMAXIMUM,VELMAXIMUM,DENSITYMAXIMUM,
&DENSITYMINIMUM, TAUMIN
double precision pi
DATA PI/3.141592653589793d0/
double precision gettempmax
C--hydrodynamic limits
COMMON /HYDROLIM/ MIDRAPLIM, TMAXLIM, TVELMAXLIM, BOOSTTR,
&GLOBALLIMS, PRETAUHYDRO
DOUBLE PRECISION MIDRAPLIM, TMAXLIM, TVELMAXLIM
LOGICAL BOOSTTR, GLOBALLIMS, PRETAUHYDRO
C--local variables
double precision px,py,pz,e,getmsmax,m,ys
!Original JEWEL implementation
!if (boost) then
! ys = etamax2
!else
! ys = 0.d0
!endif
! GETSSCAT limits are boost invariant
! thus we keep the original code without
! a transverse boost of the 4-momentum
! This must be checked.
if (boost) then
ys = MIDRAPLIM
else
ys = 0.d0
endif
m = getmsmax()
e = m * cosh(ys)
px = 0
py = 0
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
double precision getmdmin
GETMD=MDSCALEFAC*3.*GETTEMP(X1,Y1,Z1,T1)
!GETMD=MAX(GETMD,MDFACTOR)
GETMD=MAX(GETMD,GETMDMIN())
END
DOUBLE PRECISION FUNCTION GETMS(X2,Y2,Z2,T2)
IMPLICIT NONE
DOUBLE PRECISION X2,Y2,Z2,T2,GETMD
! Deprecated
GETMS=GETMD(X2,Y2,Z2,T2)/SQRT(2.)
END
DOUBLE PRECISION FUNCTION GETN0(X3,Y3,Z3,T3)
IMPLICIT NONE
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON,MODMED
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED
C--identifier of log file
common/logfile/logfid
integer logfid
C--hydro auxiliary files
COMMON /HYDROF/ INITVTXF, IDEALN0, CUSTOMN0F
CHARACTER*200 INITVTXF, CUSTOMN0F
LOGICAL IDEALN0
C--number density parameters
common/n0par/ densconst, n0array(2000), tempn0array(2000),
&temparraymaxpos, temparrayminpos
double precision densconst, n0array, tempn0array
integer temparraymaxpos, temparrayminpos
C-- local variables
DOUBLE PRECISION X3,Y3,Z3,T3,PI,GETTEMP,tau,cosheta
double precision localtemp, interpolaten0
DATA PI/3.141592653589793d0/
getn0 = 0.d0
localtemp = GETTEMP(X3,Y3,Z3,T3)
IF ((ABS(Z3).gt.T3) .OR. (localtemp.le.TC)) THEN
RETURN
END IF
if (IDEALN0) then
getn0 = densconst * localtemp ** 3
else
! Interpolate
getn0 = interpolaten0(localtemp)
end if
END
DOUBLE PRECISION FUNCTION GETNEFF(X3,Y3,Z3,T3,P0,P1,P2,P3)
IMPLICIT NONE
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,
& NX,NY,NT,NF,DX,DT,XMAX,XMIN,TMAX
INTEGER NX,NY,NT,NF
DOUBLE PRECISION DX,DT,XMAX,XMIN,TMAX
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
common/grid/timesteps(60),tprofile(834,834,60),vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
common/gridvel/ux(834,834,60),uy(834,834,60)
double precision ux,uy
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON,MODMED
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED
common/temperature/tempfac
double precision tempfac
common/boostmed/boost
logical boost
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
C--number density parameters
common/n0par/ densconst, n0array(2000), tempn0array(2000),
&temparraymaxpos, temparrayminpos
double precision densconst, n0array, tempn0array
integer temparraymaxpos, temparrayminpos
C-- local variables
DOUBLE PRECISION X3,Y3,Z3,T3,PI,GETTEMP,tau,cosheta
double precision getu,getutheta
double precision umx, umy, umz, umr, umtheta, gam, vp, gamp,
&localtemp
double precision J0, J1, J2, J3, P0, P1, P2, P3, ys
double precision getn0
DATA PI/3.141592653589793d0/
! getneff=0.d0
getneff=getn0(X3, Y3, Z3, T3)
if (getneff .eq. 0.d0) then
return
end if
! localtemp = GETTEMP(X3,Y3,Z3,T3)
! IF ((ABS(Z3).gt.T3) .OR. (localtemp.le.TC)) THEN
! RETURN
! END IF
tau = sqrt(t3**2-z3**2)
if (boost) then
umz = z3 / tau
else
umz = 0.d0
end if
!Medium 4-velocity
umr = getu(x3, y3, z3, t3, localtemp)
umtheta = getutheta(x3, y3, z3, t3, localtemp)
umx = umr * cos(umtheta)
umy = umr * sin(umtheta)
gam = sqrt(1 + umr ** 2 + umz ** 2)
!cosheta = t3/tau
! GETNEFF=(2.*6.*NF*D3*2./3. + 16.*ZETA3*3./2.)
! & *localtemp**3/PI**2
!write(*,*) "Original neff = ", GETNEFF/cosheta
!Scattering center 4-current
J0 = getneff * gam
J1 = getneff * umx
J2 = getneff * umy
J3 = getneff * umz
!Effective density transform as
!J0 = p_mu J^mu / p0, check nucl-th/0612068
getneff = (P0 * J0 - P1 * J1 - P2 * J2 - P3 * J3) / P0
END
DOUBLE PRECISION FUNCTION GETTEMP(X4,Y4,Z4,T4)
IMPLICIT NONE
C--medium parameters
COMMON/MEDPARAM/CENTRMIN,CENTRMAX,BREAL,CENTR,RAU,
& NX,NY,NT,NF,DX,DT,XMAX,XMIN,TMAX
INTEGER NX,NY,NT,NF
DOUBLE PRECISION DX,DT,XMAX,XMIN,TMAX
DOUBLE PRECISION CENTRMIN,CENTRMAX,BREAL,CENTR,RAU
common/grid/timesteps(60),tprofile(834,834,60),vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
common/gridvel/u(834,834,60),utheta(834,834,60)
double precision u,utheta
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
&N0,SIGMANN,A,WOODSSAXON,MODMED
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
&SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED
COMMON/GAMMAMAX/GAMMAMAXIMUM,VELMAXIMUM,DENSITYMAXIMUM,
&DENSITYMINIMUM, TAUMIN
DOUBLE PRECISION GAMMAMAXIMUM,VELMAXIMUM,DENSITYMAXIMUM,
&DENSITYMINIMUM, TAUMIN
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
common/temperature/tempfac
double precision tempfac
C--hydrodynamic quantities
COMMON /HYDROLIM/ MIDRAPLIM, TMAXLIM, TVELMAXLIM, BOOSTTR,
& GLOBALLIMS, PRETAUHYDRO
DOUBLE PRECISION MIDRAPLIM, TMAXLIM, TVELMAXLIM
LOGICAL BOOSTTR, GLOBALLIMS, PRETAUHYDRO
C--local variables
DOUBLE PRECISION X4,Y4,Z4,T4,TAU,NPART,EPS0,EPSIN,TEMPIN,PI,
&NTHICK,ys,MEDPART,interpolate
double precision gettempmax
DATA PI/3.141592653589793d0/
GETTEMP=0.D0
IF (ABS(Z4).gt.T4) RETURN
TAU=SQRT(T4**2-Z4**2)
! NO ASSUMPTION BEFORE TAU0
if ((tau.lt.TAUMIN-0.001) .and. (.not.PRETAUHYDRO)) return
! Consider only relevant regions for calculation
! otherwise temp = 0 => prob of interaction = 0
ys = 0.5d0 * log((T4 + Z4) / (T4 - Z4))
if (abs(ys).gt.etamax2) return
! Grid values are between -25 and 25
if ((abs(X4).gt.25.d0) .or. (abs(Y4).gt.25.d0)) return
GETTEMP=tempfac*interpolate(X4,Y4,tau,1)
if(gettemp.lt.tc) gettemp=0.0d0
if(gettemp.ge.gettempmax()) gettemp=gettempmax()
RETURN
END
double precision function getu(x,y,z,t,localtemperature)
implicit none
integer np
double precision x,y,z,t,tau,localtemperature
double precision interpol
common/grid/timesteps(60),tprofile(834,834,60),
& vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
common/gridvel/u(834,834,60),utheta(834,834,60)
double precision u,utheta
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
& N0,SIGMANN,A,WOODSSAXON,MODMED,MEDFILELIST
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
& SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED,MEDFILELIST
C--hydrodynamic quantities
COMMON /HYDROLIM/ MIDRAPLIM, TMAXLIM, TVELMAXLIM, BOOSTTR,
& GLOBALLIMS, PRETAUHYDRO
DOUBLE PRECISION MIDRAPLIM, TMAXLIM, TVELMAXLIM
LOGICAL BOOSTTR, GLOBALLIMS, PRETAUHYDRO
getu = 0.d0
IF (BOOSTTR .eqv. .false.) then
RETURN
END IF
tau = sqrt(t**2 - z**2)
IF ((tau.le.TAUI) .OR. (localtemperature.le.TC)) THEN
RETURN
END IF
getu=interpol(tau,x,y,834,timesteps,u,.true.)
return
end function
double precision function getutheta(x,y,z,t,localtemperature)
implicit none
integer np
double precision x,y,z,t,tau,localtemperature
double precision interpol
common/grid/timesteps(60),tprofile(834,834,60),
& vtxmap(834,834)
double precision timesteps,tprofile,vtxmap
common/gridvel/u(834,834,60),utheta(834,834,60)
double precision u,utheta
COMMON/MEDPARAMINT/TAUI,TI,TC,D3,ZETA3,D,
& N0,SIGMANN,A,WOODSSAXON,MODMED,MEDFILELIST
DOUBLE PRECISION TAUI,TI,TC,ALPHA,BETA,GAMMA,D3,ZETA3,D,N0,
& SIGMANN
INTEGER A
LOGICAL WOODSSAXON,MODMED,MEDFILELIST
C--hydrodynamic quantities
COMMON /HYDROLIM/ MIDRAPLIM, TMAXLIM, TVELMAXLIM, BOOSTTR,
& GLOBALLIMS, PRETAUHYDRO
DOUBLE PRECISION MIDRAPLIM, TMAXLIM, TVELMAXLIM
LOGICAL BOOSTTR, GLOBALLIMS, PRETAUHYDRO
getutheta = 0.d0
IF (BOOSTTR .eqv. .false.) then
RETURN
END IF
tau = sqrt(t**2 - z**2)
IF ((tau.le.TAUI) .OR. (localtemperature.le.TC)) THEN
RETURN
END IF
getutheta=interpol(tau,x,y,834,timesteps,utheta,.false.)
return
end function
subroutine LorentzBoostToLabFrame(e, px, py, pz, x, y, z, t)
implicit none
C--longitudinal boost of momentum distribution
common/boostmed/boost
logical boost
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
double precision fourmom(4,1), fourmomboost(4,1),
& boostm(4,4)
double precision ux, uy, uz, gam, u2, ux2, uy2, uz2,
& unorm, uangle, tau
double precision getu,getutheta
C--local variables
double precision e, px, py, pz, x, y, z, t
! Use temp as 10.d0 to always calculate (T_C < 10)
unorm = getu(x, y, z, t, 10.d0)
uangle = getutheta(x, y, z, t, 10.d0)
tau = sqrt(t ** 2 - z ** 2)
! Note that all velocity components must be flipped
! since they are defined in the lab frame, thus the frame
! must be boosted as -u so the scattering centers are
! boosted as u.
if (tau.gt.0.d0 .and. boost) then
uz = -z / tau
else
uz = 0.d0
end if
! Respect simulation limit in uz
if (abs(uz).gt.sinh(etamax2)) then
if (z.gt.0d0) then
uz = -sinh(etamax2)
else
uz = sinh(etamax2)
end if
end if
! Only transform if there is velocity (careful with sign)
if (unorm.gt.0.d0 .or. abs(uz).gt.0.d0) then
ux = -unorm * cos(uangle)
uy = -unorm * sin(uangle)
uz2 = uz ** 2
u2 = unorm ** 2 + uz2
ux2 = ux ** 2
uy2 = uy ** 2
gam = sqrt(1 + u2)
! Define initial four-momentum
fourmom(1,1) = e
fourmom(2,1) = px
fourmom(3,1) = py
fourmom(4,1) = pz
! Define boost matrix
boostm(1,1) = gam
boostm(1,2) = -ux
boostm(1,3) = -uy
boostm(1,4) = -uz
boostm(2,1) = -ux
boostm(2,2) = 1 + (gam - 1) * ux2 / u2
boostm(2,3) = (gam - 1) * ux * uy / u2
boostm(2,4) = (gam - 1) * ux * uz / u2
boostm(3,1) = -uy
boostm(3,2) = (gam - 1) * uy * ux / u2
boostm(3,3) = 1 + (gam - 1) * uy2 / u2
boostm(3,4) = (gam - 1) * uy * uz / u2
boostm(4,1) = -uz
boostm(4,2) = (gam - 1) * ux * uz / u2
boostm(4,3) = (gam - 1) * uy * uz / u2
boostm(4,4) = 1 + (gam - 1) * uz2 / u2
! Calculate new four-momentum
fourmomboost = matmul(boostm, fourmom)
e = fourmomboost(1,1)
px = fourmomboost(2,1)
py = fourmomboost(3,1)
pz = fourmomboost(4,1)
end if
end subroutine
subroutine LorentzBoostToFluidFrame(e, px, py, pz, x, y, z, t)
implicit none
C--longitudinal boost of momentum distribution
common/boostmed/boost
logical boost
C--max rapidity
common/rapmax2/etamax2
double precision etamax2
double precision fourmom(4,1), fourmomboost(4,1),
& boostm(4,4)
double precision ux, uy, uz, gam, u2, ux2, uy2, uz2,
& unorm, uangle, tau
double precision getu,getutheta
C--local variables
double precision e, px, py, pz, x, y, z, t
! Use temp as 10.d0 to always calculate (T_C < 10)
unorm = getu(x, y, z, t, 10.d0)