EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
pepsiMaineRHIC_noradcorr.v1.F
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file pepsiMaineRHIC_noradcorr.v1.F
1  program pepsimainerhic
2 
3  implicit none
4  include "leptou.inc"
5  include "lujets.inc"
6  include "ludat1.inc"
7  include "mcRadCor.inc"
8  include "phiout.inc"
9  include "mc_set.inc"
10 
11  integer nevgen, i, j, lastgenevent, genevent
12  integer iev, nev, nprt
13  integer idum1, idum2, initseed
14  integer ltype, tracknr, t1, t2
15  integer*4 today(3), now(3)
16  real rval, alpha
17  real mcfixedweight, mcextraweight
18  double precision dilut, f1, f2, r, da1, da2
19  double precision gendilut, genf1, genf2, genr, gena1, gena2
20  DOUBLE PRECISION radgamenucl
21  DOUBLE PRECISION pbeame, pbeta, pgamma, ebeame, epznucl
22  DOUBLE PRECISION depol, gamma2, d, eta, eps, chi
23  DOUBLE PRECISION gendepol, gend, geneta, geneps, genchi
24  DOUBLE PRECISION dbrack, dxsec
25 
26  CHARACTER param*20
27 
28 c---------------------------------------------------------------------
29 c ascii output file
30 c ---------------------------------------------------------------------
31  integer asciilun
32  parameter(asciilun=29)
33  CHARACTER*256 outputfilename
34  CHARACTER*256 outname
35 
36  alpha = 1./137.
38  qedrad = 0
39  nevgen = 0
40  genevent = 0
41  lastgenevent = 0.
42  dilut= 1.
43  r = 0.
44  depol =0.
45  da1 = 0.
46  da2 = 0.
47  f1 = 0.
48  f2 = 0.
49  gamma2= 1.
50  d = 0.
51  eta = 0.
52  eps = 0.
53  chi = 0.
54  gendilut= 1.
55  genr = 0.
56  gendepol =0.
57  gena1 = 0.
58  gena2 = 0.
59  genf1 = 0.
60  genf2 = 0.
61  gend = 0.
62  geneta = 0.
63  geneps = 0.
64  genchi = 0.
65  radgamenucl = 0.
66  mcfixedweight=1.
67  mcextraweight=1.
68  weight=1.
69 c init all the parameters and read them in
70  ltype = 11
71  pbeam = 100.
72  ebeam = 4.
73 
74 c...Read parameters for linit call
75 c...Read output file name
76  READ(*,*) outname
77 c...Read beam target particle type
78  READ(*,*) ltype
79 c...Read beam and target particle energy
80  READ(*,*) pbeam, ebeam
81 c...Read number of events to generate, and to print.
82  READ(*,*) nev,nprt
83 c...Read min/max x of radgen lookup table
84  READ(*,*) mcset_xmin, mcset_xmax
85 c...Read min/max y of generation range
86  READ(*,*) mcset_ymin, mcset_ymax
87 c...Read min/max q2 of generation range
88  READ(*,*) mcset_q2min, mcset_q2max
89 c...Read information for cross section used in radgen
90  READ(*,*) genset_fstruct, genset_r
91 c...Read switch for rad-corrections
92  READ(*,*) qedrad
93 c...Read min/max w2 of generation range
94  READ(*,*) mcset_w2min,mcset_w2max
95 c...Read min/max nu of generation range
96  READ(*,*) mcset_numin,mcset_numax
97 c...Read beam and Target polarisation direction
98  READ(*,*) mcset_pbvalue, mcset_ptvalue
99 c...Read beam and Target polarisation
100  READ(*,*) mcset_pbeam, mcset_ptarget
101 c...Read all the pepsi/lepto parameters
102  100 READ(*,*, end=200) param, j, rval
103  if (param.eq.'LST') then
104  lst(j)=rval
105  write(*,*) 'LST(',j,')=',rval
106  endif
107  if (param.eq.'PARL') then
108  parl(j)=rval
109  write(*,*) 'PARL(',j,')=',rval
110  endif
111  if (param.eq.'PARJ') then
112  parj(j)=rval
113  write(*,*) 'PARJ(',j,')=',rval
114  endif
115  if (param.eq.'MSTJ') then
116  mstj(j)=rval
117  write(*,*) 'MSTJ(',j,')=',rval
118  endif
119  if (param.eq.'MSTU') then
120  mstu(j)=rval
121  write(*,*) 'MSTU(',j,')=',rval
122  endif
123  if (param.eq.'PARU') then
124  paru(j)=rval
125  write(*,*) 'PARU(',j,')=',rval
126  endif
127  goto 100
128  200 write(*,*) '*********************************************'
129  write(*,*) 'NOW all parameters are read by PEPSI'
130  write(*,*) '*********************************************'
131 
132 c ... pepsi polarization
133  if (mcset_pbvalue*mcset_ptvalue.lt.0.) then
134  lst(40) = 1.
135  elseif (mcset_pbvalue*mcset_ptvalue.gt.0.) then
136  lst(40) = -1.
137  else
138  lst(40) = 0.
139  endif
140 
141 c ... polarization for helium 3 !
142 c use lst(39)=1 for ideal helium 3 (no s prime and d state)
143 c use lst(39)=3 or (02 should work too ?) for a real helium 3.
144 c lst(39) = 0
145  if (lst(40).ne.0.and.parl(1).eq.3.and.parl(2).eq.2) then
146  lst(39) = 3
147  lst(36) = 1
148  endif
149 
150 c ... xmin,xmax
151  cut(1)=mcset_xmin
152  cut(2)=mcset_xmax
153 c ... ymin,ymax
154  cut(3)=mcset_ymin
155  cut(4)=mcset_ymax
156 c ... q2min, q2max
157  cut(5)=mcset_q2min
158  cut(6)=mcset_q2max
159 c ... w2min, w2max
160  cut(7)=mcset_w2min
161 c cut(8)=mcset_w2max
162 c ... numin, numax
163 c cut(9)=mcset_numin
164 c cut(10)=mcset_numax
165 
166  mcset_tara=parl(1)
167  mcset_tarz=parl(2)
168 
169 c getting the date and time of the event generation
170  call idate(today) ! today(1)=day, (2)=month, (3)=year
171  call itime(now) ! now(1)=hour, (2)=minute, (3)=second
172 c
173 c take date as the seed for the random number generation
174  initseed = today(1) + 10*today(2) + today(3) + now(1) + 5*now(3)
175  write(6,*) 'SEED = ', initseed
176  call rndmq(idum1,idum2,initseed,' ')
177 
178  if (abs(ltype).eq.11) then
179  masse=0.000510998910
180  else
181  write(*,*) "Only positron and electron beam supported"
182  endif
183 
184  if ((parl(1).eq.1).and.(parl(2).eq.1)) then
185  massp=0.938272013
186  elseif ((parl(1).eq.1).and.(parl(2).eq.0)) then
187  massp=0.93956536
188  elseif ((parl(2).eq.2).and.(parl(2).eq.1)) then
189  massp=1.875612793
190  else
191  write(*,*) "Currently only p,n,d as target supported"
192  endif
193 
194  pbeame=sqrt(pbeam**2+massp**2)
195  pbeta=pbeam/pbeame
196  pgamma=pbeame/massp
197  ebeame=sqrt(ebeam**2+masse**2)
198  ebeamenucl=pgamma*ebeame-pgamma*pbeta*(-ebeam)
199  epznucl=-pgamma*pbeta*(ebeame)+pgamma*(-ebeam)
200  write(*,*) ebeamenucl, ebeame, epznucl, -ebeam
201  mcset_enebeam=sngl(ebeamenucl)
202 
203 c initialize pepsi
204  CALL timex(t1)
205  call linit(0,ltype,(-ebeam),pbeam,lst(23))
206  CALL timex(t2)
207 
208 c ---------------------------------------------------------------------
209 c Open ascii output file
210 c ---------------------------------------------------------------------
211  outputfilename=outname
212  open(asciilun, file=outputfilename)
213  write(*,*) 'the outputfile will be named: ', outname
214 
215 c the total integrated xsec, in pbarns
216  mcextraweight = parl(23)
217 
218 ! ============== Generate until happy (past cuts) =================
219 
220 c this is what we write in the ascii-file
221 
222  write(29,*)' PEPSI EVENT FILE '
223  write(29,*)'============================================'
224  write(29,30)
225 30 format('I, ievent, genevent, process, subprocess, nucleon,
226  & struckparton, partontrck, trueY, trueQ2, trueX, trueW2,
227  & trueNu, FixedWeight, weight, dxsec, Extraweight,
228  & dilut, F1, F2, A1, A2, R, Depol, d, eta, eps, chi,
229  & gendilut, genF1, genF2, genA1, genA2, genR, genDepol, gend,
230  & geneta, geneps, genchi, Sigcor, radgamEnucl, nrTracks')
231 
232  write(29,*)'============================================'
233  write(29,*)' I K(I,1) K(I,2) K(I,3) K(I,4) K(I,5)
234  & P(I,1) P(I,2) P(I,3) P(I,4) P(I,5) V(I,1) V(I,2) V(I,3)'
235  write(29,*)'============================================'
236 
237  DO 300 iev=1,nev
238 
239 c count generated events
240 10 nevgen = nevgen+1
241 
242  if (qedrad.eq.0) then
243  If (iev.eq.1) then
244  write(*,*) 'Generating full event with LEPTO/PEPSI.'
245  endif
246  DO i=1,mstu(4)
247  DO j=1,5
248  k(i,j)=0.
249  p(i,j)=0.
250  v(i,j)=0.
251  enddo
252  enddo
253  lst(7) = 1
254  call lepto
255  if (lst(21).ne.0) then
256  write(*,*) 'PEPSI failed, trying again,
257  & no rad-corrections'
258  goto 10
259  endif
260  da2=0.
261  da1=0.
262  f1=0.
263  f2=0.
264  dilut=1.
265  r=0.
266  depol=0.
267  gamma2=0.
268  d=0.
269  eta=0.
270  eps=0.
271  chi=0.
272  genevent=nevgen-lastgenevent
273 
274  endif
275 
276  DO i=1,mstu(4)
277  if ((k(i,1).eq.0).and.(k(i,2).eq.0).and.(k(i,3).eq.0).and.
278  & (k(i,4).eq.0).and.(k(i,5).eq.0)) then
279  tracknr=i-1
280  goto 400
281  endif
282  ENDDO
283 400 continue
284 
285 c get some info for the event
286  call mkasym(dble(q2), dble(x), mcset_tara, mcset_tarz, da1, da2)
287  call mkf2(dble(q2), dble(x), mcset_tara, mcset_tarz, f2, f1)
288  call mkr(dble(q2), dble(x), r)
289  call fdilut(dble(q2), dble(x), mcset_tara, dilut)
290 
291 c... calculate some kinematic variables
292  gamma2=q2/u**2
293  depol=(y*((1.+gamma2*y/2.)*(2.-y)-2*y**2*masse**2/q2))/
294  & (y**2*(1.-2.*masse**2/q2)*(1.+gamma2)+2.*(1.+r)*
295  & (1.-y-gamma2*(y**2)/4.))
296  d=depol*((sqrt(1.-y-gamma2*(y**2)/4.)*(1+gamma2*y/2.)) /
297  & (1.-y/2.)*(1+gamma2*y/2.)-y**2*masse**2/q2)
298  eta=(sqrt(gamma2)*(1.-y-gamma2*y**2/4.-y**2*masse**2/q2)) /
299  & ((1.+gamma2**2*y/2.)*(1-y/2.)-y**2*masse**2/q2)
300  chi=(sqrt(gamma2)*(1.-y/2.-y**2.*masse**2/q2))/(1.+gamma2*y/2.)
301  eps=1./(1.+1./2.*(1.-2.*masse**2/q2)*((y**2+gamma2*y**2)/
302  & (1.-y-gamma2*(y**2)/4.)))
303 
304  dbrack = (4.*paru(1)*alpha**2)/q2**2
305 c cross section in hbar^2c^2/gev^2/dq^2/dx
306  dxsec = (1.-y-(massp*x*y/2./ebeame))*f2/x+y**2*f1
307 c cross sectin in mbarn/dq^2/dx
308  dxsec = dbrack * dxsec * 0.389379292
309 
310  write(29,32) 0, iev, genevent, lst(23), lst(24), lst(22),
311  & lst(25), lst(26), y, q2, x, w2, u, mcfixedweight, weight,
312  & dxsec, mcextraweight, real(dilut), real(F1),
313  & real(F2), real(dA1), real(dA2), real(R), real(Depol), real(d),
314  & real(eta), real(eps), real(chi), real(gendilut), real(genF1),
315  & real(genF2), real(genA1), real(genA2), real(genR),
316  & real(genDepol), real(gend), real(geneta), real(geneps),
317  & real(genchi), mcradcor_sigcor, radgamenucl, tracknr
318  32 format((i4,1x,$),(i10,1x,$),6(i6,1x,$),
319  & 33(e18.10,1x,$),i12,/)
320  write(29,*)'============================================'
321 
322  DO i=1,tracknr
323  write(29,34) i,k(i,1),k(i,2),k(i,3),k(i,4),k(i,5),
324  & p(i,1),p(i,2),p(i,3),p(i,4),p(i,5),
325  & v(i,1),v(i,2),v(i,3)
326  ENDDO
327  34 format(2(i6,1x,$),i10,1x,$,3(i6,1x,$),8(f15.6,1x,$),/)
328  write(29,*)'=============== Event finished ==============='
329 
330  lastgenevent=nevgen
331 
332  300 CONTINUE
333 
334 c...print the pythia cross section which is needed to get an absolut
335 c normalisation the number is in pbarns
336  write(*,*)'==================================================='
337  write(*,*)'Pepsi total cross section in pb from numerical',
338  & ' integration', parl(23)
339  write(*,*)'Pepsi total cross section in pb from MC simulation',
340  & parl(24)
341  write(*,*)'Total Number of generated events', nev
342  write(*,*)'Total Number of trials', nevgen
343  write(*,*)'==================================================='
344  close(29)
345  END