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