-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathinit.f90
3274 lines (2859 loc) · 117 KB
/
init.f90
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
! ****************************************************************************
! * *
! * ECHO-QGP *
! * *
! * Copyright (C) 2015-2019 The ECHO-QGP team *
! * *
! * File: init.f90 *
! * *
! * License: GPL version 2.0 (Please, read the file LICENSE.TXT) *
! * *
! * This program is free software; you can redistribute it and/or *
! * modify it under the terms of the GNU General Public License *
! * as published by the Free Software Foundation; either version 2 *
! * of the License, or (at your option) any later version. *
! * *
! * This program is distributed in the hope that it will be useful, *
! * but WITHOUT ANY WARRANTY; without even the implied warranty of *
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
! * GNU General Public License for more details. *
! * *
! * You should have received a copy of the GNU General Public License *
! * along with this program; if not, write to the Free Software *
! * Foundation, Inc., 51 Franklin Street, Fifth Floor, *
! * Boston, MA 02110-1301, USA. *
! * *
! * Authors: Gabriele Inghirami ([email protected]) *
! * Francesco Becattini ([email protected]) *
! * Vinod Chandra ([email protected]) *
! * Andrea Beraudo ([email protected]) *
! * Arturo De Pace ([email protected]) *
! * *
! * Contributors: Luca Del Zanna ([email protected]) *
! * *
! * Acknowledgments: P. Cea, J. Rizzo, L. Cosmai, G. Denicol *
! * *
! *****************************************************************************
module qgpinit
use eos
use parallel
use common
use viscous
use glaubermc
implicit none
real(8) :: xx,yy,eta
real(8),allocatable,dimension(:) :: thick !thickness vector
integer, parameter :: minus=-1, center=0, plus=1 !flag/parameters for changing the simmetry of energy distribution
contains
! ***************************************************************************
subroutine print_param()
implicit none
real(8) :: s_input, en_dens_output
integer ferror
character(len=120) :: outbuffer
010 format(a50,2x,i6)
011 format(a50,2(2x,i6))
012 format(a50,3(2x,i6))
013 format(a50,2x,e14.7)
014 format(a50,2(2x,e14.7))
015 format(a50,3(2x,e14.7))
016 format(a51,2x,e14.7)
if(ipe .eq. 0) then
write(*,*) ' Using as parameters file: '//param_file
write(*,*) ' Settings and parameters:'
write(*,*)
select case(init_type)
case(GLAUBER_GEOMETRIC)
write(*,"(2x,a4,1x,i2,2x,a44)") 'Test',GLAUBER_GEOMETRIC,'- optical-geometrical Glauber initialization'
case(SHOCK_TUBE_2D)
write(*,"(2x,a4,1x,i2,2x,a27)") 'Test',SHOCK_TUBE_2D,'- 2D+1 shock tube init_type'
case(SHEAR_VISCOUS_1D)
write(*,"(2x,a4,1x,i2,2x,a33)") 'Test',SHEAR_VISCOUS_1D,'- 1D viscous shear flow init_type'
case(GLAUBER_MONTECARLO)
if(avg_events .eq. 0) then
write(*,"(2x,a4,1x,i2,2x,a37)") 'Test',GLAUBER_MONTECARLO,'- Glauber-Monte Carlo initial profile'
else if(avg_events .eq. 1) then
write(*,"(2x,a4,1x,i2,2x,a54)") 'Test',GLAUBER_MONTECARLO,'- run with averaged Glauber-Monte Carlo in. conditions'
else
write(*,*) " Sorry, but it is not clear whether I should run event by event simulations"
write(*,*) " or average many initial conditions and perform just 1 simulation. I quit."
call exit(1)
end if
case(GUBSER_IDEAL)
write(*,"(2x,a4,1x,i2,2x,a31)") 'Test',GUBSER_IDEAL,'- ideal Gubser flow init_type'
case(GUBSER_VISCOUS)
write(*,"(2x,a4,1x,i2,2x,a31)") 'Test',GUBSER_VISCOUS,'- viscous Gubser flow init_type'
case(TABULATED_INIT)
write(*,"(2x,a4,1x,i2,2x,a53)") 'Test',TABULATED_INIT,'- tabulated initial entropy or energy density profile'
case(ALFVEN2D)
write(*,"(2x,a4,1x,i2,2x,a20)") 'Test',ALFVEN2D,'- alfven 2D mhd test'
case(BLAST3D)
write(*,"(2x,a4,1x,i2,2x,a23)") 'Test',BLAST3D,'- blast wave 3D mhd test'
case(MHD1DBJ)
write(*,"(2x,a4,1x,i2,2x,a28)") 'Test',MHD1DBJ,'- 0D+1 bjorken flow mhd test'
case(SODMHD)
write(*,"(2x,a4,1x,i2,2x,a18)") 'Test',SODMHD,'- 1 D Sod mhd test'
case(LYU)
write(*,"(2x,a4,1x,i2,2x,a24)") 'Test',LYU,'- 1 D Lyutikov rmhd test'
case(BERHAD)
write(*,"(2x,a4,1x,i2,2x,a31)") 'Test',BERHAD,'- 1 D entropy density rmhd test'
case(ROTOR)
write(*,"(2x,a4,1x,i2,2x,a36)") 'Test',ROTOR,'- 2D or 3D Rotor rmhd test'
case(OT)
write(*,"(2x,a4,1x,i2,2x,a42)") 'Test',OT,'- 2D Orszang-Tang test'
case(DBG)
write(*,"(2x,a4,1x,i2,2x,a42)") 'Test',DBG,'- 2D grid with B field on transverse plane'
case(EXT_CONS)
write(*,"(2x,a4,1x,i2,2x,a43)") 'Test',EXT_CONS,'- external file with conservative variables'
case(EXT_GLISSANDO)
write(*,"(2x,a4,1x,i2,2x,a35)") 'Test',EXT_GLISSANDO,'- external file with Glissando data'
end select
write(*,"(i1)",advance='no') coordinates_option
if(coordinates .eq. MINKOWSKI) then
write(*,*) ' Using Minkowski coordinates'
else
write(*,*)' Using Bjorken coordinates'
end if
end if
if(ipe .eq. 0) then
write(*,"(i1)",advance='no') boundary_cond_option
end if
if(boundary_cond .eq. 0) then
if(ipe .eq. 0) write(*,*) ' Using periodic boundary conditions'
ibc=0
else if (boundary_cond .eq. 2) then
if(ipe .eq. 0) write(*,*) ' Using outflow 0th order boundary conditions'
ibc=2
else
if(ipe .eq. 0) then
write(*,*) ' Sorry, but ECHO-QGP has been tested only with 0 (periodic) or 2 (outflow) boundary conditions'
write(*,*) ' I am not sure to provide correct results, so I quit.'
end if
call exit(2)
end if
!-- Redefine BCs if parallel
if (prl) then
if (ipe/= 0 .or. ibc(1,1,1)==0) ibc(:,1,1)=-2
if (ipe/=npe-1 .or. ibc(1,2,1)==0) ibc(:,2,1)=-2
end if
if(ipe .eq. 0) then
write(*,"(i1)",advance='no') mhd_flag_option
if(mhd) then
if(rmhd) then
write(*,*) ' This is a resistive MHD simulation'
else
write(*,*) ' This is an ideal MHD simulation'
end if
write(*,"(i1)",advance='no') divclean_flag_option
if(divclean) then
write(*,*) " with Dedner's method for divergence cleaning enabled"
write(*,"(i1)",advance='no') glm_alpha_option
write(*,013) " with a value of the dumping parameter glm_alpha: ",glm_alpha
else
write(*,*) " without any method to control the solenoidal condition"
end if
else
write(*,*) ' This is an HD simulation'
end if
write(*,"(i1)",advance='no') viscosity_flag_option
if(viscosity) then
write(*,*) ' This is a viscous simulation'
write(*,"(i1)",advance='no') bulkvis_flag_option
if(bulkvis) then
write(*,*) ' Bulk viscosity is taken into account'
else
write(*,*) ' Bulk viscosity is neglected'
end if
write(*,"(i1)",advance='no') obtained_option
if(obtained .eq. 'no') then
write(*,*) ' Evolved shear viscous tensor components: xx, yy, zz, xy, xz, zz'
write(*,*) ' tt, tx, ty and tz components are obtained imposing orthogonality'
else
write(*,*) ' Evolved shear viscous tensor components: xx, yy, zz, xy, xz'
write(*,*) ' tt, tx, ty, tz and zz are obtained imposing orth. and null trace'
end if
else
write(*,*) ' This is an unviscid simulation'
end if
if(prl) then
write(*,*) ' This simulation uses MPI and the number of processors is:',npe
else
write(*,*) ' This is a serial run (no MPI)'
end if
write(*,*)
write(*,*) ' Grid parameters:'
write(*,"(i1)",advance='no') nx_option
write(*,010) ' number of cells along x direction: ', nx
write(*,"(i1)",advance='no') ny_option
write(*,010) ' number of cells along x direction: ', ny
write(*,"(i1)",advance='no') nz_option
write(*,010) ' number of cells along z (or eta) direction: ', nz
write(*,"(i1)",advance='no') x1_option
write(*,013) ' x min (fm): ',x1
write(*,"(i1)",advance='no') x2_option
write(*,013) ' x max (fm): ',x2
write(*,"(i1)",advance='no') y1_option
write(*,013) ' y min (fm): ',y1
write(*,"(i1)",advance='no') y2_option
write(*,013) ' y max (fm): ',y2
write(*,"(i1)",advance='no') z1_option
write(*,013) ' z min (fm): ',z1
write(*,"(i1)",advance='no') z2_option
write(*,013) ' z max (fm): ',z2
write(*,*)
write(*,*) ' Time parameters:'
write(*,"(i1)",advance='no') tmin_option
write(*,013) ' starting time: ', tmin
write(*,"(i1)",advance='no') tmax_option
write(*,013) ' ending time: ', tmax
write(*,"(i1)",advance='no') temp_end_option
write(*,013) ' ending temeperature (MeV): ', temp_end*1000.
write(*,"(i1)",advance='no') maxdt_option
write(*,013) ' maximum timestep: ', maxdt
write(*,*)
if(viscosity) then
write(*,*) ' Viscous parameters:'
write(*,"(i1)",advance='no') eta_over_s_option
write(*,013) ' eta/s parameter for shear viscosity tensor: ', eta_over_s
write(*,"(i1)",advance='no') smooth_temp_option
if(enable_smooth_viscosity) write(*,013) ' Temperature limit for smoothing viscosity: ', smooth_temp
write(*,*)
end if
write(*,"(i1)",advance='no') eqstate_option
write(*,010) ' Equation of state: ', eqstate
if(eqstate .eq. 3) then
write(*,"(i1)",advance='no') eos_tab_name_option
write(*,*) ' File with tabulated eos: '//adjustl(eos_tab_name)
end if
write(*,"(i1)",advance='no') numerically_option
write(*,010) ' numerical derivatives with analytic eos: ', numerically
write(*,*)
if((init_type .eq. GLAUBER_GEOMETRIC) .or. (init_type .eq. GLAUBER_MONTECARLO) .or. (init_type .eq. TABULATED_INIT) &
& .or. (init_type .eq. DBG)) then
write(*,*) ' Nucleus parameters:'
write(*,"(i1)",advance='no') nuclei_data_file_option
write(*,*) ' file with nuclear parameters: '//adjustl(nuclei_data_file)
write(*,016) ' nucleus mass (GeV): ',projmass
write(*,016) ' nuclear radius (fm): ', radius
write(*,016) ' W.-S. width (fm): ',delta
write(*,"(i1)",advance='no') rads_option
write(*,013) ' sqrt(s) (GeV): ',rads
write(*,016) ' beam rapidity Y_b: ',yb
write(*,"(i1)",advance='no') sigma_in_option
write(*,013) ' cross section (mb): ',sigma_in
write(*,"(i1)",advance='no') bimpact_option
write(*,013) ' impact parameter (fm): ',bimpact
write(*,"(i1)",advance='no') ah_option
write(*,013) ' initial hardness parameter: ',ah
write(*,"(i1)",advance='no') ecenter_option
if(ienentr .eq. 0) then
write(*,013) ' central energy density (GeV/fm^3): ',ecenter
else if(ienentr .eq. 1) then
write(*,013) ' central entropy density (fm^-3): ',ecenter
if(eqstate .ne. 3) then
write(*,*) " Sorry, but currently you can initialize with entropy density only if you use a tabulated EoS..."
call exit(3)
end if
else
write(*,*) " Sorry, but I don't understand the value of the IENENTR parameter..."
call exit(3)
end if
write(*,"(i1)",advance='no') enezero_option
write(*,013) ' enezero (GeV/fm^3): ',enezero
write(*,016) ' przero (GeV/fm^3): ',przero
write(*,"(i1)",advance='no') rhocenter_option
write(*,013) ' central density: ',rhocenter
write(*,"(i1)",advance='no') deta_option
write(*,013) ' rapidity distrib. shift (deta or etaflat par.): ',deta
write(*,"(i1)",advance='no') sigeta_option
write(*,013) ' rapidity distrib. width (sigeta parameter): ',sigeta
write(*,*)
end if
if(init_type .eq. GLAUBER_GEOMETRIC) then
write(*,"(i1)",advance='no') ueta_coeff_option
write(*,013) ' ueta A coeff. so that u^eta=A*x: ', ueta_coeff
end if
write(*,*)
if(init_type .eq. EXT_CONS) then
write(*,"(i1)",advance='no') ext_cons_file_option
write(*,013) ' File containing initial conservative variables: '//adjustl(ext_cons_file)
end if
write(*,*)
if(init_type .eq. EXT_GLISSANDO) then
write(*,"(i1)",advance='no') ext_glissando_file_option
write(*,013) ' File containing initial variables from Glissando: '//adjustl(ext_glissando_file)
end if
write(*,*)
if(mhd) then
write(*,"(i1)",advance='no') input_B_field_option
write(*,*) ' File with initial magnetic field configuration: '//adjustl(input_B_field)
write(*,"(i1)",advance='no') B_amplify_factor_option
write(*,013) ' Initial B field amplification factor: ', B_amplify_factor
write(*,"(i1)",advance='no') dump_initial_B_flag_option
if(dump_initial_B) then
write(*,*) ' Initial B field dumping enabled'
else
write(*,*) ' Initial B field dumping disabled'
end if
if(init_type .eq. GLAUBER_MONTECARLO) then
write(*,"(i1)",advance='no') magfield_type_option
if(magfield_type .eq. 1 ) then
write(*,*) ' The initial B field will have both classical and chiral sources'
else if(magfield_type .eq. 2 ) then
write(*,*) ' The initial B field will have classical origin only'
else if(magfield_type .eq. 3 ) then
write(*,*) ' The initial B field will have chiral origin only'
else
write(*,*) ' Sorry, but I cannot understand which kind of initial magnetic field you want'
write(*,*) ' (i.e. classical, chiral or both)'
call exit(2)
end if
if((magfield_type .eq. 1) .or. (magfield_type .eq. 2)) then
write(*,"(i1)",advance='no') sigma_el_cond_option
write(*,013) ' Electrical cond. of the medium to compute in. B: ', sigma_el_cond
end if
if((magfield_type .eq. 1) .or. (magfield_type .eq. 3)) then
write(*,"(i1)",advance='no') sigma_chiral_cond_option
write(*,013) ' Chiral cond. of the medium to compute initial B: ', sigma_chiral_cond
end if
end if
write(*,"(i1)",advance='no') edens_for_B_dump_option
write(*,013) ' En. dens. limit for the B field supp.(GeV/fm^3):', edens_for_B_dump
write(*,"(i1)",advance='no') B_th_press_ratio_option
write(*,013) ' Max ratio between magnetic and thermal press.: ', B_th_press_ratio
if((out_sel%rc .eq. 1) .and. out_freeze) then
print_rho_comov=.true.
write(*,*) "Printing electric charge density in comoving frame on the freezeout hypersurface"
end if
end if
write(*,*)
if(init_type .eq. GLAUBER_MONTECARLO) then!for Glauber-MonteCarlo
write(*,*) ' Glauber-MonteCarlo specific parameters:'
write(*,"(i1)",advance='no') kind_of_collision_option
select case(kind_of_collision)
case(1)
write(*,*) ' kind of collision: nucleus-nucleus'
case(2)
write(*,*) ' kind of collision: deuton-nucleus'
case(3)
write(*,*) ' kind of collision: proton-nucleus'
end select
write(*,"(i1)",advance='no') nconf_option
write(*,010) ' number of nuclear configurations: ',nconf
if(fixed_b) then
write(*,"(i1)",advance='no') fixed_b_flag_option
write(*,*) ' Using a fixed impact parameter'
else
write(*,"(i1)",advance='no') fixed_b_flag_option
write(*,*) ' Using a randomly sampled impact parameter'
write(*,"(i1)",advance='no') nbcoll_option
write(*,010) ' number of impact parameters per configuration: ', nbcoll
end if
write(*,"(i1)",advance='no') events_start_option
write(*,010) ' id number of starting event: ', events_start
write(*,"(i1)",advance='no') events_stop_option
write(*,010) ' id number of ending event: ', events_stop
write(*,"(i1)",advance='no') kappa_option
write(*,013) ' k model parameters: ',kappa
write(*,"(i1)",advance='no') sig_mc_option
write(*,013) ' sigma model smearing parameter: ',sig_mc
write(*,"(i1)",advance='no') min_participants_option
write(*,010) ' min number of participants to accept the event: ',min_participants
if(kind_of_collision .eq. 3) then
write(*,"(i1)",advance='no') eprot_option
write(*,013) ' proton beam energy: ',eprot
end if
end if
write(*,*)
if(out_freeze) then !for freeze-out hypersurface computation
write(*,"(i1)",advance='no') out_freeze_flag_option
write(*,*) ' Computing freeze-out hypersurface'
if(freeze_type .eq. 0) then
write(*,"(i1)",advance='no') freeze_type_option
write(*,*) ' hypersurface computation based on temperature'
write(*,"(i1)",advance='no') freeze_value_option
write(*,013) ' freezeout threshold (MeV): ', freeze_value*1000.
else
write(*,"(i1)",advance='no') freeze_type_option
write(*,*) ' hypersurface computation based on energy density'
write(*,"(i1)",advance='no') freeze_value_option
write(*,013) ' freezeout threshold (GeV/f^3): ', freeze_value
end if
write(*,"(i1)",advance='no') freeze_time_interval_option
write(*,013) ' time interval between hypersurf. computations: ', freeze_time_interval
end if
if(print_rho_comov) then
write(*,"(i1)",advance='no') print_rho_comov_flag_option
write(*,*) ' Computing the electric charge density in the comoving frame on f.o. hypersurface'
end if
write(*,*)
write(*,*) ' Other numerical parameters:'
if(mhd) then
write(*,"(i1)",advance='no') algo_solver_option
write(*,010) ' Inv. algorithm for obtaining prim. variables: ',algo_solver
end if
write(*,"(i1)",advance='no') cfl_option
write(*,013) ' Courant-Fr.-Lew. condition parameter: ',cfl
write(*,"(i1)",advance='no') nrk_option
if(nrk .eq. 2) then
write(*,*) ' Time integration algorithm: RK2'
else if(nrk .eq. 3) then
write(*,*) ' Time integration algorithm: RK3'
else
write(*,*) ' Time integration algorithm: IMEX-SSP3(4,3,3)'
end if
write(*,"(i1)",advance='no') recal_option
write(*,*) ' Reconstruction algorithm: ',recal
write(*,"(i1)",advance='no') maxspeed_option
write(*,013) ' If the velocity is > c, the speed is reduced to:',maxspeed
write(*,*)
write(*,*) ' Output parameters:'
write(*,"(i1)",advance='no') dtlog_option
write(*,013) ' interval between log updating: ', dtlog
write(*,"(i1)",advance='no') dtout_option
write(*,013) ' interval between output printing: ', dtout
write(*,"(i1)",advance='no') output_precision_option
if(output_precision .eq. 8) then
write(*,*) ' output precision:', 'double - 8 bytes'
else
write(*,*) ' output precision:', 'float - 4 bytes'
end if
write(*,*)
write(*,*) ' Variables printed in the output files:'
write(*,"(i1)",advance='no') density_option
if(out_sel%density .eq. 1) write(*,*) ' density'
write(*,"(i1)",advance='no') vx_option
if(out_sel%vx .eq. 1) write(*,*) ' vx'
write(*,"(i1)",advance='no') vy_option
if(out_sel%vy .eq. 1) write(*,*) ' vy'
write(*,"(i1)",advance='no') vz_option
if(out_sel%vz .eq. 1) write(*,*) ' vz'
write(*,"(i1)",advance='no') pressure_option
if(out_sel%pressure .eq. 1) write(*,*) ' pressure'
write(*,"(i1)",advance='no') energy_density_option
if(out_sel%energy_density .eq. 1) write(*,*) ' energy density'
write(*,"(i1)",advance='no') temperature_option
if(out_sel%temperature .eq. 1) write(*,*) ' temperature'
write(*,"(i1)",advance='no') entropy_density_option
if(out_sel%entropy_density .eq. 1) write(*,*) ' entropy density'
if(viscosity) then
write(*,"(i1)",advance='no') bulk_option
if(out_sel%bulk .eq. 1) then
write(*,*) ' bulk viscosity'
else
write(*,*) ' no bulk viscosity'
end if
write(*,"(i1)",advance='no') pitt_option
if(out_sel%pitt .eq. 1) then
write(*,*) ' pi^tt'
else
write(*,*) ' no pi^tt'
end if
write(*,"(i1)",advance='no') pitx_option
if(out_sel%pitx .eq. 1) then
write(*,*) ' pi^tx'
else
write(*,*) ' no pi^tx'
end if
write(*,"(i1)",advance='no') pity_option
if(out_sel%pity .eq. 1) then
write(*,*) ' pi^ty'
else
write(*,*) ' no pi^ty'
end if
write(*,"(i1)",advance='no') pitz_option
if(out_sel%pitz .eq. 1) then
write(*,*) ' pi^tz'
else
write(*,*) ' no pi^tz'
end if
write(*,"(i1)",advance='no') pixy_option
if(out_sel%pixy .eq. 1) then
write(*,*) ' pi^xy'
else
write(*,*) ' no pi^xy'
end if
write(*,"(i1)",advance='no') pixz_option
if(out_sel%pixz .eq. 1) then
write(*,*) ' pi^xz'
else
write(*,*) ' no pi^xz'
end if
write(*,"(i1)",advance='no') piyz_option
if(out_sel%piyz .eq. 1) then
write(*,*) ' pi^yz'
else
write(*,*) ' no pi^yz'
end if
write(*,"(i1)",advance='no') pixx_option
if(out_sel%pixx .eq. 1) then
write(*,*) ' pi^xx'
else
write(*,*) ' no pi^xx'
end if
write(*,"(i1)",advance='no') piyy_option
if(out_sel%piyy .eq. 1) then
write(*,*) ' pi^yy'
else
write(*,*) ' no pi^yy'
end if
write(*,"(i1)",advance='no') pizz_option
if(out_sel%pizz .eq. 1) then
write(*,*) ' pi^zz'
else
write(*,*) ' no pi^zz'
end if
end if
write(*,"(i1)",advance='no') v0_option
if(out_sel%v0 .eq. 1) write(*,*) ' u0 or gamma Lorentz factor'
if(mhd) then
write(*,"(i1)",advance='no') bx_option
if(out_sel%bx .eq. 1) write(*,*) ' bx'
write(*,"(i1)",advance='no') by_option
if(out_sel%by .eq. 1) write(*,*) ' by'
write(*,"(i1)",advance='no') bz_option
if(out_sel%bz .eq. 1) write(*,*) ' bz'
write(*,"(i1)",advance='no') ex_option
if(out_sel%ex .eq. 1) write(*,*) ' ex'
write(*,"(i1)",advance='no') ey_option
if(out_sel%ey .eq. 1) write(*,*) ' ey'
write(*,"(i1)",advance='no') ez_option
if(out_sel%ez .eq. 1) write(*,*) ' ez'
write(*,"(i1)",advance='no') glm_option
if(out_sel%glm .eq. 1) write(*,*) ' glm'
write(*,"(i1)",advance='no') rc_option
if(out_sel%rc .eq. 1) write(*,*) ' rc'
end if
if(derivatives_out) then
write(*,"(i1)",advance='no') derivatives_flag_option
write(*,*) ' derivatives will also be printed into separated output files'
end if
if(flows_out) then
write(*,"(i1)",advance='no') flows_flag_option
write(*,*) ' flows will also be printed into separated output files'
end if
end if !it closes (ipe .eq. 0)
end subroutine print_param
!---------------------------------------------------------------
subroutine print_config_summary
implicit none
integer ferror
character(len=120) :: outbuffer
010 format(a50,2x,i6)
011 format(a50,2(2x,i6))
012 format(a50,3(2x,i6))
013 format(a50,2x,e14.7)
014 format(a50,2(2x,e14.7))
015 format(a50,3(2x,e14.7))
open(unit=25,file="config_summary.dat",status='replace',IOSTAT=ferror)
if(ferror .ne. 0) then
write(*,*) "Sorry, but I cannot open the file config_summary.dat"
write(*,*) "Although this is not an essential file, nevertheless this means that there are problems in writing to HD"
write(*,*) "and therefore I prefer to quit now and let you solve this issue. Bye!"
call exit(2)
end if
write (outbuffer, '(i2)') init_type
write(25,"(a10,a2)") "INIT_TYPE=",adjustl(outbuffer)
write(25,"(a10,i1)") "COORD....=",coordinates
if(viscosity) then
write(25,"(a11)") "VISCOUS..=1"
else
write(25,"(a11)") "VISCOUS..=0"
end if
if(bulkvis) then
write(25,"(a11)") "BULK.....=1"
else
write(25,"(a11)") "BULK.....=0"
end if
if(mhd) then
write(25,"(a11)") "MHD......=1"
else
write(25,"(a11)") "MHD......=0"
end if
if(divclean) then
write(25,"(a11)") "DIVCLEAN.=1"
else
write(25,"(a11)") "DIVCLEAN.=0"
end if
write (outbuffer, '(f14.4)') glm_alpha
write(25,"(a10,a16)") "GLM_PARAM=",adjustl(outbuffer)
if(dump_initial_B) then
write(25,"(a11)") "DUMP_IN_B=1"
else
write(25,"(a11)") "DUMP_IN_B=0"
end if
write (outbuffer, '(f14.8)') edens_for_B_dump
write(25,"(a10,a16)") "B_DUM_EN.=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') B_th_press_ratio
write(25,"(a10,a16)") "Bp_ov_Tp.=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') smooth_temp
write(25,"(a10,a16)") "CUT_TEMP.=",adjustl(outbuffer)
write (outbuffer, '(i6)') nx
write(25,"(a10,a6)") "NX.......=",adjustl(outbuffer)
write (outbuffer, '(i6)') ny
write(25,"(a10,a6)") "NY.......=",adjustl(outbuffer)
write (outbuffer, '(i6)') nz
write(25,"(a10,a6)") "NZ.......=",adjustl(outbuffer)
write(*,*) " config_summary.dat written"
write (outbuffer, '(f14.8)') x1
write(25,"(a10,a16)") "XMIN.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') x2
write(25,"(a10,a16)") "XMAX.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') y1
write(25,"(a10,a16)") "YMIN.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') y2
write(25,"(a10,a16)") "YMAX.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') z1
write(25,"(a10,a16)") "ZMIN.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') z2
write(25,"(a10,a16)") "ZMAX.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') tmin
write(25,"(a10,a16)") "TSTART...=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') tmax
write(25,"(a10,a16)") "TSTOP....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') temp_end
write(25,"(a10,a16)") "TEMP_END.=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') dtlog
write(25,"(a10,a16)") "DTLOG....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') dtout
write(25,"(a10,a16)") "DTOUT....=",adjustl(outbuffer)
write (outbuffer, '(i1)') output_precision
write(25,"(a10,a1)") "OUTP_PREC=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') maxdt
write(25,"(a10,a16)") "MAXDT....=",adjustl(outbuffer)
write (outbuffer, '(i1)') restart_type
write(25,"(a10,a1)") "RESTART..=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') cfl
write(25,"(a10,a16)") "CFL......=",adjustl(outbuffer)
write(25,"(a10,a5)") "REC_ALGO.=",adjustl(recal)
write (outbuffer, '(i1)') algo_solver
write(25,"(a10,a1)") "ALG_SOLV.=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') maxspeed
write(25,"(a10,a16)") "MAXSPEED.=",adjustl(outbuffer)
write (outbuffer, '(a60)') nuclei_data_file
write(25,"(a10,a60)") "NUC_DATAF=",adjustl(outbuffer)
write(25,"(a10,a5)") "NUCLEUS..=",adjustl(nucleus)
write (outbuffer, '(f14.8)') rads
write(25,"(a10,a16)") "RADS.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') sigma_in*10.
! sigma has been transformed into fm^2, but it was given in mb = 1/10 fm^2 -> we multiply by 10 to get mb again
write(25,"(a10,a16)") "SIGMA_IN.=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') bimpact
write(25,"(a10,a16)") "B........=",adjustl(outbuffer)
write (outbuffer, '(i1)') ienentr
write(25,"(a10,a1)") "IENENTR..=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') ah
write(25,"(a10,a16)") "AH.......=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') ecenter
write(25,"(a10,a16)") "ECENTER..=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') enezero
write(25,"(a10,a16)") "ENEZERO..=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') rhocenter
write(25,"(a10,a16)") "RHOCENTER=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') deta
write(25,"(a10,a16)") "DETA.....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') sigeta
write(25,"(a10,a16)") "SIGETA...=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') eta_over_s
write(25,"(a10,a16)") "ETA_S....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') tau_pi_coeff
write(25,"(a10,a16)") "TAU_PI_C.=",adjustl(outbuffer)
write(25,"(a10,a2)") "TRACE_IMP=",adjustl(obtained)
write (outbuffer, '(i1)') eqstate
write(25,"(a10,a1)") "EOS......=",adjustl(outbuffer)
write (outbuffer, '(a60)') eos_tab_name
write(25,"(a10,a60)") "EOS_FILE.=",adjustl(outbuffer)
write (outbuffer, '(i1)') numerically
write(25,"(a10,a1)") "NUM_DER..=",adjustl(outbuffer)
write (outbuffer, '(i6)') nconf
write(25,"(a10,a6)") "NCONF....=",adjustl(outbuffer)
write (outbuffer, '(i6)') nbcoll
write(25,"(a10,a6)") "NBCOLL...=",adjustl(outbuffer)
if(fixed_b) then
write(25,"(a11)") "FIXED_B..=1"
else
write(25,"(a11)") "FIXED_B..=0"
end if
write (outbuffer, '(i6)') events_start
write(25,"(a10,a6)") "EV_START.=",adjustl(outbuffer)
write (outbuffer, '(i6)') events_stop
write(25,"(a10,a6)") "EV_STOP..=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') kappa
write(25,"(a10,a16)") "KAPPA....=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') sig_mc
write(25,"(a10,a16)") "SIG......=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') eprot
write(25,"(a10,a16)") "EPROT....=",adjustl(outbuffer)
write (outbuffer, '(i1)') avg_events
write(25,"(a10,a1)") "AVG_IC...=",adjustl(outbuffer)
write (outbuffer, '(i3)') min_participants
write(25,"(a10,a6)") "MIN_PART.=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') pw
write(25,"(a10,a16)") "PW.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') kind_of_collision
write(25,"(a10,a1)") "COLLISION=",adjustl(outbuffer)
if(out_freeze) then
write(25,"(a11)") "HYP_COMPU=1"
else
write(25,"(a11)") "HYP_COMPU=0"
end if
write (outbuffer, '(i1)') freeze_type
write(25,"(a10,a1)") "FREEZKIND=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') freeze_value
write(25,"(a10,a16)") "FREEZEVAL=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') freeze_time_interval
write(25,"(a10,a16)") "HYPSURFTI=",adjustl(outbuffer)
write (outbuffer, '(a60)') input_edf
write(25,"(a10,a60)") "IN_D_FILE=",adjustl(outbuffer)
write (outbuffer, '(a60)') ext_cons_file
write(25,"(a10,a60)") "EX_C_FILE=",adjustl(outbuffer)
write (outbuffer, '(a60)') ext_glissando_file
write(25,"(a10,a60)") "EX_G_FILE=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') etam
write(25,"(a10,a16)") "ETAM_TILT=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') ueta_coeff
write(25,"(a10,a16)") "UETA_COEF=",adjustl(outbuffer)
write (outbuffer, '(a60)') input_B_field
write(25,"(a10,a60)") "B_in_FILE=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') sigma_el_cond
write(25,"(a10,a16)") "EL_COND..=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') sigma_chiral_cond
write(25,"(a10,a16)") "CHIR_COND=",adjustl(outbuffer)
write (outbuffer, '(i1)') magfield_type
write(25,"(a10,a1)") "MAGFIELDT=",adjustl(outbuffer)
write (outbuffer, '(f14.8)') B_amplify_factor
write(25,"(a10,a16)") "B_amp_fac=",adjustl(outbuffer)
write (outbuffer, '(i1)') boundary_cond
write(25,"(a10,a1)") "bound_con=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%density
write(25,"(a10,a1)") "density..=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%vx
write(25,"(a10,a1)") "vx.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%vy
write(25,"(a10,a1)") "vy.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%vz
write(25,"(a10,a1)") "vz.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pressure
write(25,"(a10,a1)") "pressure.=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%energy_density
write(25,"(a10,a1)") "ene_dens.=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%temperature
write(25,"(a10,a1)") "temper...=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%entropy_density
write(25,"(a10,a1)") "entr_dens=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%bulk
write(25,"(a10,a1)") "bulk_visc=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pitt
write(25,"(a10,a1)") "pi^tt....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pitx
write(25,"(a10,a1)") "pi^tx....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pity
write(25,"(a10,a1)") "pi^ty....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pitz
write(25,"(a10,a1)") "pi^tz....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pixy
write(25,"(a10,a1)") "pi^xy....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pixz
write(25,"(a10,a1)") "pi^xz....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%piyz
write(25,"(a10,a1)") "pi^yz....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pixx
write(25,"(a10,a1)") "pi^xx....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%piyy
write(25,"(a10,a1)") "pi^yy....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%pizz
write(25,"(a10,a1)") "pi^zz....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%v0
write(25,"(a10,a1)") "gamma....=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%bx
write(25,"(a10,a1)") "bx.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%by
write(25,"(a10,a1)") "by.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%bz
write(25,"(a10,a1)") "bz.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%ex
write(25,"(a10,a1)") "ex.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%ey
write(25,"(a10,a1)") "ey.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%ez
write(25,"(a10,a1)") "ez.......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%glm
write(25,"(a10,a1)") "glm......=",adjustl(outbuffer)
write (outbuffer, '(i1)') out_sel%rc
write(25,"(a10,a1)") "rc.......=",adjustl(outbuffer)
if(derivatives_out) then
write(25,"(a11)") "derivativ=1"
else
write(25,"(a11)") "derivativ=0"
end if
if(flows_out) then
write(25,"(a11)") "flows....=1"
else
write(25,"(a11)") "flows....=0"
end if
close(25)
end subroutine print_config_summary
! ***************************************************************************
subroutine initio (nu,xini,vv,error_flag)
implicit none
integer, intent(in):: nu
real(8),intent(in) :: xini(3)
real(8) :: vv(nu)
real(8) ymax,tr_weight
real(8) rho,edens,press,dpr1,dpr2
real(8) velo(3)
integer error_flag
error_flag=0
xx =xini(1)
yy =xini(2)
eta=xini(3)
!call calc_enedens(edens,errcode) !it matches the version of calc_enedens provided by A.DP
call calc_enedens(edens)
call calc_chargedens(rho)
if(init_type .eq. GLAUBER_GEOMETRIC) call calc_vel_long_nz_ueta(velo,xx,yy,eta)
if(ienentr .eq. 0) then
call eos_pressure(rho,edens,press,dpr1,dpr2,error_flag)
else
call eos_pressure_from_s(edens,press,error_flag)
end if
press=press+przero
if(error_flag .gt. 0) then
write(*,*) "An error occurred into module init, subroutine initio, when computing the pressure profile..."
errcode=error_flag
write(*,*) "Error code:",errcode
write(*,*) "Position: x=",xx,"y=",yy,"z=",eta
run_crashed=.true.
return
end if
vv(1)=rho
if(.not. mhd) then
vv(2)=velo(1)
vv(3)=velo(2)
vv(4)=velo(3)
end if
vv(5)=press
return
end subroutine initio
! ***************************************************************************
subroutine calc_enedens(edens)
implicit none
real(8) edens,tt1,tt2,H, HP, HM, thplus, thminus,sigeta2_inv
sigeta2_inv=1./sigeta**2
if (nz .eq. 1) then ! 2+1 with Glauber (ah*N_coll+(1-ah)*N_wound)
call wound(amassa,tt1,tt2,thplus,thminus)
edens=ecenter*(ah*bincoll()+(1.-ah)*(tt1+tt2))/((ah*weightbin)+((1-ah)*weightwound))
else
! 3+1 used both by Hirano and Shenke.
! The reference here is hep-ph 1004.1408v1 (Shenke, Jeon, Gale) Formulas 51,52
call wound(amassa,tt1,tt2, thplus, thminus)
H=H_of_eta(eta, deta, sigeta2_inv, thplus, thminus,center)
HP=H_of_eta(eta, deta, sigeta2_inv, thplus, thminus,plus)
HM=H_of_eta(eta, deta, sigeta2_inv, thplus, thminus,minus)
if(tilting) then !initial energy density tilting applied
edens=ecenter*(ah*bincoll()*H+2.*(1.-ah)*(tt1*HM+tt2*HP))/((ah*weightbin)&
&+((1-ah)*weightwound))
!write(*,*) 'edens:',edens
else
edens=H*ecenter*(ah*bincoll()+(1.-ah)*(tt1+tt2))/((ah*weightbin)+((1-ah)*weightwound))
endif
endif
end subroutine calc_enedens
! ***************************************************************************
function H_of_eta(etac, etaflat, sigmaetasq_inv, thplus, thminus, simmetry_flag)
! The reference here is hep-ph 1004.1408v1 (Shenke, Jeon, Gale) Formulas 51,52
implicit none
real(8), intent(in) :: thplus, thminus
real(8) H_of_eta
real(8) etac, etaflat, sigmaetasq_inv, eta0
integer simmetry_flag
if(tilting) then !in the case of an initial energy density tilting,eta0=0
eta0=0.
else
eta0=0.5*log(((thminus+thplus)*gamma_col+(thminus-thplus)*gamma_col*beta_col)/&
&((thminus+thplus)*gamma_col-(thminus-thplus)*gamma_col*beta_col))
end if
if (abs(etac-eta0) .gt. (etaflat*0.5)) then
H_of_eta=exp(-1.0* (( abs(etac-eta0)-(0.5*etaflat) )**2) * 0.5* sigmaetasq_inv )
else
H_of_eta=1.0
endif
if(simmetry_flag .eq. plus) then
if(etac-eta0 .lt. -etam) then
H_of_eta=0.
else if(etac-eta0 .le. etam) then
H_of_eta=H_of_eta*(etac-eta0+etam)/(2.*etam)
end if
else if(simmetry_flag .eq. minus) then
if(etac-eta0 .gt. etam) then
H_of_eta=0.
else if(etac-eta0 .ge. -etam) then
H_of_eta=H_of_eta*(-etac+eta0+etam)/(2.*etam)
end if
end if
end function H_of_eta
! ***************************************************************************
subroutine calc_chargedens(rho)
implicit none
real(8) rho
rho=rhocenter
return
end subroutine calc_chargedens
! ***************************************************************************
subroutine calc_vel_long_nz_ueta(velo,xx,yy,eta)
implicit none
real(8) velo(3)
real(8) :: vsound = 1/sqrt(3.)-1.e-6
real(8) :: glg_ueta
real(8) :: ueta
real(8), intent(in) :: xx, yy, eta
integer thickarray_index
real(8) thick_plus, thick_minus, thick_zero, thick_limit, cut_factor, rplus, rminus
velo(1)=0.
velo(2)=0.
velo(3)=0.
! here we define u^eta - yb is the beam rapidity, xx and eta the x and eta coordinates, t is the time
! ueta_coeff is chosen inside the parameter file (usually param.dat )