1 C*********************************************************************
3 C... LHAGLUE Interface to LHAPDF library of modern parton
4 C... density functions (PDF) with uncertainties
6 C...Authors for v4: Dimitri Bourilkov, Craig Group, Mike Whalley
8 C...Authors for v3: Dimitri Bourilkov, Craig Group, Mike Whalley
10 C...Author for v1 and v2: Dimitri Bourilkov bourilkov@mailaps.org
11 C... University of Florida
13 C...HERWIG interface by Dimitri Bourilkov and Craig Group
15 C...New numbering scheme and upgrade for LHAPDF v2.1
16 C...by Dimitri Bourilkov and Mike Whalley
18 C...For more information, or when you cite this interface, currently
19 C...the official reference is:
20 C...D.Bourilkov, "Study of Parton Density Function Uncertainties with
21 C...LHAPDF and PYTHIA at LHC", hep-ph/0305126.
23 C...The official LHAPDF page is:
25 C... http://durpdg.dur.ac.uk/lhapdf/index.html
27 C...The interface contains four subroutines (similar to PDFLIB).
28 C...It can be used seamlessly by Monte Carlo generators
29 C...interfaced to PDFLIB or in stand-alone mode.
31 C... For initialization (called once)
33 C... PDFSET(PARM,VALUE)
35 C... For the proton/pion structure functions
37 C... STRUCTM(X,Q,UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU)
39 C... For the photon structure functions
41 C... STRUCTP(X,Q2,P2,IP2,UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU)
43 C... For statistics ON structure functions (under/over-flows)
47 C...This interface can be invoked in 3 ways depending
48 C...on the value of PARM(1) provided by the user when
49 C...calling PDFSET(PARM,VALUE):
51 C... For PYTHIA: PARM(1).EQ.'NPTYPE'
52 C... (this is set automatically by PYTHIA)
54 C... For HERWIG: PARM(1).EQ.'HWLHAPDF'
55 C... (set by the USER e.g. in the main program like this:
56 C... AUTPDF(1) = 'HWLHAPDF'
57 C... AUTPDF(2) = 'HWLHAPDF' )
59 C... For Stand-alone: PARM(1).EQ.'DEFAULT'
60 C... (can be used for PDF studies or when interfacing
63 C...The LHAPDF set/member is selected depending on the value of:
65 C... PYTHIA: ABS(MSTP(51)) - proton
66 C... ABS(MSTP(53)) - pion
67 C... ABS(MSTP(55)) - photon
69 C... HERWIG: ABS(INT(VALUE(1)))
71 C... STAND-ALONE: ABS(INT(VALUE(1)))
74 C... CONTROL switches:
75 C... ==================
77 C... THE LOCATION OF THE LHAPDF LIBRARY HAS TO BE SPECIFIED
78 C... AS DESCRIBED BELOW (the rest is optional)
80 C... if the user does nothing, sensible defaults
81 C... are active; to change the behaviour, the corresponding
82 C... values of LHAPARM() should be set to the values given below
84 C... Location of the LHAPDF library of PDFs (pathname):
85 C... uses common block /LHAPDFC/
87 C... If the user does nothing => default = subdir PDFsets of the
88 C... current directory (can be real subdir
89 C... OR a soft link to the real location)
90 C... If the user sets LHAPATH => supplied by the USER who defines the
91 C... path in common block COMMON/LHAPDFC/LHAPATH
92 C... BEFORE calling PDFSET
96 C... use common block /LHACONTROL/
98 C... Collect statistics on under/over-flow requests for PDFs
99 C... outside their validity ranges in X and Q**2
100 C... (call PDFSTA at end of run to print it out)
102 C... LHAPARM(16).EQ.'NOSTAT' => No statistics (faster)
103 C... LHAPARM(16).NE.'NOSTAT' => Default: collect statistics
105 C... Option to use the values for the strong coupling alpha_s
106 C... as computed in LHAPDF in the MC generator
107 C... (to ensure uniformity between the MC generator and the PDF set)
108 C... WARNING: implemented ONLY for PYTHIA in LHAPDFv4
110 C... LHAPARM(17).EQ.'LHAPDF' => Use alpha_s from LHAPDF
111 C... LHAPARM(17).NE.'LHAPDF' => Default (same as LHAPDF v1/v3)
113 C... Extrapolation of PDFs outside LHAPDF validity range given by
114 C... [Xmin,Xmax] and [Q2min,Q2max]; DEFAULT => PDFs "frozen" at the
117 C... LHAPARM(18).EQ.'EXTRAPOLATE' => Extrapolate PDFs on OWN RISK
118 C... WARNING: Crazy values can be returned
120 C... Printout of initialization information in PDFSET (by default)
122 C... LHAPARM(19).EQ.'SILENT' => No printout (silent mode)
123 C... LHAPARM(19).EQ.'LOWKEY' => Print 5 times (almost silent mode)
126 C...v5.0 06-Oct-2005 Major change to allow multiset-initializations
127 C...v4.0 28-Apr-2005 PDFSTA routine; option to use Alfa_s from LHAPDF
128 C...v4.0 21-Mar-2005 Photon/pion/new p PDFs, updated for LHAPDF v4
129 C...v3.1 26-Apr-2004 New numbering scheme, updated for LHAPDF v2/v3
130 C...v3.0 23-Jan-2004 HERWIG interface added
131 C...v2.0 20-Sep-2003 PDFLIB style adopted
132 C...v1.0 05-Mar-2003 First working version from PYTHIA to LHAPDF v1
134 C...interface to LHAPDF library
136 C*********************************************************************
138 C...Initialization for use of parton distributions
139 C... according to the LHAPDF interface.
141 C...v4.0 28-Apr-2005 Option to use Alfa_s from LHAPDF
142 C...v4.0 21-Mar-2005 Photon/pion/new p PDFs, updated for LHAPDF v4
143 C...v3.1 26-Apr-2004 New numbering scheme
144 C...v3.0 23-Jan-2004 HERWIG interface added
146 C...interface to LHAPDF library
148 SUBROUTINE PDFSET(PARM,VALUE,
149 > MSTU11,MSTP51,MSTP53,MSTP55,
151 > AXMIN,AXMAX,AQ2MIN,AQ2MAX)
152 C...Double precision and integer declarations.
153 IMPLICIT DOUBLE PRECISION(A-H, O-Z)
154 IMPLICIT INTEGER(I-N)
155 c...additions for multiset use
156 include 'parmsetup.inc'
157 include 'pathsetup.inc'
158 c character*172 LHANAMES(nmxset)
159 integer LHAMEMBERS(nmxset),LHANUMBERS(nmxset)
160 common/LHASETS/LHANAMES,LHANUMBERS,LHAMEMBERS,nsets
161 real*8 xxmin(nmxset),xxmax(nmxset),qq2min(nmxset),qq2max(nmxset)
162 save xxmin,xxmax,qq2min,qq2max
163 C...Interface to LHAPDFLIB.
164 c CHARACTER*172 LHANAME
165 INTEGER LHASET, LHAMEMB
166 COMMON/LHAPDF/LHANAME, LHASET, LHAMEMB
168 DOUBLE PRECISION QCDLHA4, QCDLHA5
170 COMMON/LHAPDFR/QCDLHA4, QCDLHA5, NFLLHA
172 c CHARACTER*132 LHAPATH
173 COMMON/LHAPDFC/LHAPATH
175 CHARACTER*20 LHAPARM(20)
176 DOUBLE PRECISION LHAVALUE(20)
177 COMMON/LHACONTROL/LHAPARM,LHAVALUE
180 COMMON/LHAPDFE/LHAEXTRP
183 COMMON/LHASILENT/LHASILENT
185 DOUBLE PRECISION XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
186 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
187 COMMON/LHAGLSTA/ XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
188 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
190 C...Interface to PDFLIB.
191 COMMON/W50511/ NPTYPEPDFL,NGROUPPDFL,NSETPDFL,MODEPDFL,
192 > NFLPDFL,LOPDFL,TMASPDFL
194 DOUBLE PRECISION TMASPDFL
195 C...Interface to PDFLIB.
196 COMMON/W50513/XMIN,XMAX,Q2MIN,Q2MAX
198 DOUBLE PRECISION XMIN,XMAX,Q2MIN,Q2MAX
199 C...Local arrays and character variables (NOT USED here DB)
200 CHARACTER*20 PARM(20)
201 DOUBLE PRECISION VALUE(20)
203 DOUBLE PRECISION QCDL4,QCDL5
204 DOUBLE PRECISION AXMIN,AXMAX,AQ2MIN,AQ2MAX
217 CHARACTER*1000 CHROOT
223 if(first .AND. (LHAPARM(20).NE.'LHAPATH')) then
224 c...overide the default PDFsets path
225 c ... check first if the environmental variable LHAPATH is set
226 call getenv('LHAPATH',lhapath)
227 if(lhapath.eq.'') then
228 C The environment variable LHAPATH is not set.
229 C Take the data from $ALICE_ROOT/LHAPDF/PDFsets
230 CALL GETENV('ALICE_ROOT',CHROOT)
231 LNROOT = LNBLNK(CHROOT)
233 LHAPATH='PDFsets' ! Default value
235 LHAPATH=CHROOT(1:LNROOT)//'/LHAPDF/PDFsets'
244 IF(LHAPARM(18).EQ.'EXTRAPOLATE')
245 > THEN ! Extrapolate PDFs on own risk
249 IF(LHAPARM(19).EQ.'SILENT') THEN ! No printout (silent MODE)
251 ELSEIF(LHAPARM(19).EQ.'LOWKEY') THEN ! Print 5 times (lowkey MODE)
252 IF(LHAFIVE .LT. 6) THEN
253 LHAFIVE = LHAFIVE + 1
258 IF(PARM(1).EQ.'NPTYPE') THEN ! PYTHIA
260 IF(VALUE(1) .EQ. 1) THEN ! nucleon
261 LHAINPUT = ABS(MSTP51)
262 ELSEIF(VALUE(1) .EQ. 2) THEN ! pion
263 LHAINPUT = ABS(MSTP53)
264 ELSEIF(VALUE(1) .EQ. 3) THEN ! photon
265 LHAINPUT = ABS(MSTP55)
268 > PRINT *,'==== PYTHIA WILL USE LHAPDF ===='
269 ELSEIF(PARM(1).EQ.'HWLHAPDF') THEN ! HERWIG
270 LHAINPUT = ABS(INT(VALUE(1)))
271 IF(LHAONCE.EQ.LHAINPUT) RETURN
273 > PRINT *,'==== HERWIG WILL USE LHAPDF ===='
276 ELSEIF(PARM(1).EQ.'DEFAULT') THEN ! Stand-alone
277 LHAINPUT = ABS(INT(VALUE(1)))
278 IF(LHAONCE.EQ.LHAINPUT) RETURN
280 > PRINT *,'==== STAND-ALONE LHAGLUE MODE TO USE LHAPDF ===='
284 PRINT *,'== UNKNOWN LHAPDF INTERFACE CALL! STOP EXECUTION! =='
287 C...Initialize parton distributions: LHAPDFLIB.
288 LHAPATHLEN=INDEX(LHAPATH,' ')-1
290 XMIN = 1.0D-6 ! X_min for current PDF set
291 XMAX = 1.0D0 ! X_max for current PDF set
292 Q2MIN = 1.0D0**2 ! Q**2_min scale for current PDF set [GeV]
293 Q2MAX = 1.0D5**2 ! Q**2_max scale for current PDF set [GeV]
298 IF((LHAINPUT .GE. 10000) .AND. (LHAINPUT .LE. 19999)) THEN
300 IF((LHAINPUT .GE. 10000) .AND. (LHAINPUT .LE. 10040)) THEN
302 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq6.LHpdf'
304 ELSEIF((LHAINPUT .GE. 10041) .AND. (LHAINPUT .LE. 10041)) THEN
306 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq6l.LHpdf'
308 ELSEIF((LHAINPUT .GE. 10042) .AND. (LHAINPUT .LE. 10042)) THEN
310 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq6ll.LHpdf'
312 ELSEIF((LHAINPUT .GE. 10050) .AND. (LHAINPUT .LE. 10090)) THEN
314 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq6mE.LHgrid'
316 ELSEIF((LHAINPUT .GE. 10100) .AND. (LHAINPUT .LE. 10140)) THEN
318 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq61.LHpdf'
320 ELSEIF((LHAINPUT .GE. 10150) .AND. (LHAINPUT .LE. 10190)) THEN
322 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq61.LHgrid'
324 ELSEIF((LHAINPUT .GE. 10250) .AND. (LHAINPUT .LE. 10269)) THEN
326 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq6AB.LHgrid'
328 ELSEIF((LHAINPUT .GE. 10350) .AND. (LHAINPUT .LE. 10390)) THEN
330 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq65.LHgrid'
334 ELSEIF((LHAINPUT .GE. 10450) .AND. (LHAINPUT .LE. 10456)) THEN
336 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq65c.LHgrid'
340 ELSEIF((LHAINPUT .GE. 19050) .AND. (LHAINPUT .LE. 19050)) THEN
342 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq5m.LHgrid'
344 ELSEIF((LHAINPUT .GE. 19051) .AND. (LHAINPUT .LE. 19051)) THEN
346 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq5m1.LHgrid'
348 ELSEIF((LHAINPUT .GE. 19053) .AND. (LHAINPUT .LE. 19053)) THEN
350 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq5f3.LHgrid'
352 ELSEIF((LHAINPUT .GE. 19054) .AND. (LHAINPUT .LE. 19054)) THEN
354 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq5f4.LHgrid'
356 ELSEIF((LHAINPUT .GE. 19060) .AND. (LHAINPUT .LE. 19060)) THEN
358 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq5d.LHgrid'
360 ELSEIF((LHAINPUT .GE. 19070) .AND. (LHAINPUT .LE. 19070)) THEN
362 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq5l.LHgrid'
364 ELSEIF((LHAINPUT .GE. 19150) .AND. (LHAINPUT .LE. 19150)) THEN
366 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq4m.LHgrid'
369 ELSEIF((LHAINPUT .GE. 19160) .AND. (LHAINPUT .LE. 19160)) THEN
371 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq4d.LHgrid'
374 ELSEIF((LHAINPUT .GE. 19170) .AND. (LHAINPUT .LE. 19170)) THEN
376 LHANAME=LHAPATH(1:LHAPATHLEN)//'/cteq4l.LHgrid'
380 WRITE(LHAPRINT,5150) LHASET
384 ELSEIF((LHAINPUT .GE. 20000) .AND. (LHAINPUT .LE. 29999)) THEN
388 IF((LHAINPUT .GE. 20000) .AND. (LHAINPUT .LE. 20004)) THEN
390 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2001nlo.LHpdf'
391 ELSEIF((LHAINPUT .GE. 20050) .AND. (LHAINPUT .LE. 20054)) THEN
393 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2001nlo.LHgrid'
394 ELSEIF((LHAINPUT .GE. 20060) .AND. (LHAINPUT .LE. 20061)) THEN
396 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2001lo.LHgrid'
397 ELSEIF((LHAINPUT .GE. 20070) .AND. (LHAINPUT .LE. 20074)) THEN
399 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2001nnlo.LHgrid'
400 ELSEIF((LHAINPUT .GE. 20100) .AND. (LHAINPUT .LE. 20130)) THEN
402 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2001E.LHpdf'
403 ELSEIF((LHAINPUT .GE. 20150) .AND. (LHAINPUT .LE. 20180)) THEN
405 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2001E.LHgrid'
406 ELSEIF((LHAINPUT .GE. 20200) .AND. (LHAINPUT .LE. 20201)) THEN
408 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2002nlo.LHpdf'
409 ELSEIF((LHAINPUT .GE. 20250) .AND. (LHAINPUT .LE. 20251)) THEN
411 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2002nlo.LHgrid'
412 ELSEIF((LHAINPUT .GE. 20270) .AND. (LHAINPUT .LE. 20271)) THEN
414 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2002nnlo.LHgrid'
415 ELSEIF((LHAINPUT .GE. 20300) .AND. (LHAINPUT .LE. 20301)) THEN
417 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2003cnlo.LHpdf'
420 ELSEIF((LHAINPUT .GE. 20350) .AND. (LHAINPUT .LE. 20351)) THEN
422 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2003cnlo.LHgrid'
425 ELSEIF((LHAINPUT .GE. 20370) .AND. (LHAINPUT .LE. 20371)) THEN
427 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2003cnnlo.LHgrid'
430 ELSEIF((LHAINPUT .GE. 20400) .AND. (LHAINPUT .LE. 20401)) THEN
432 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004nlo.LHpdf'
433 ELSEIF((LHAINPUT .GE. 20406) .AND. (LHAINPUT .LE. 20407)) THEN
435 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004FF3nlo.LHpdf'
436 ELSEIF((LHAINPUT .GE. 20408) .AND. (LHAINPUT .LE. 20409)) THEN
438 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004FF4nlo.LHpdf'
439 ELSEIF((LHAINPUT .GE. 20450) .AND. (LHAINPUT .LE. 20451)) THEN
441 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004nlo.LHgrid'
442 ELSEIF((LHAINPUT .GE. 20452) .AND. (LHAINPUT .LE. 20453)) THEN
444 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004FF3lo.LHgrid'
445 ELSEIF((LHAINPUT .GE. 20454) .AND. (LHAINPUT .LE. 20455)) THEN
447 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004FF4lo.LHgrid'
448 ELSEIF((LHAINPUT .GE. 20456) .AND. (LHAINPUT .LE. 20457)) THEN
450 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004FF3nlo.LHgrid'
451 ELSEIF((LHAINPUT .GE. 20458) .AND. (LHAINPUT .LE. 20459)) THEN
453 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004FF4nlo.LHgrid'
454 ELSEIF((LHAINPUT .GE. 20460) .AND. (LHAINPUT .LE. 20462)) THEN
456 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004qed.LHgrid'
457 ELSEIF((LHAINPUT .GE. 20470) .AND. (LHAINPUT .LE. 20471)) THEN
459 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2004nnlo.LHgrid'
460 ELSEIF((LHAINPUT .GE. 20550) .AND. (LHAINPUT .LE. 20580)) THEN
462 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST2006nnlo.LHgrid'
466 ELSEIF((LHAINPUT .GE. 29000) .AND. (LHAINPUT .LE. 29003)) THEN
468 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST98.LHpdf'
469 ELSEIF((LHAINPUT .GE. 29040) .AND. (LHAINPUT .LE. 29045)) THEN
471 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST98lo.LHgrid'
472 ELSEIF((LHAINPUT .GE. 29050) .AND. (LHAINPUT .LE. 29055)) THEN
474 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST98nlo.LHgrid'
475 ELSEIF((LHAINPUT .GE. 29060) .AND. (LHAINPUT .LE. 29065)) THEN
477 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST98dis.LHgrid'
478 ELSEIF((LHAINPUT .GE. 29070) .AND. (LHAINPUT .LE. 29071)) THEN
480 LHANAME=LHAPATH(1:LHAPATHLEN)//'/MRST98ht.LHgrid'
482 WRITE(LHAPRINT,5150) LHASET
486 ELSEIF((LHAINPUT .GE. 30000) .AND. (LHAINPUT .LE. 39999)) THEN
487 IF((LHAINPUT .GE. 30100) .AND. (LHAINPUT .LE. 30200)) THEN
489 LHANAME=LHAPATH(1:LHAPATHLEN)//'/Fermi2002_100.LHpdf'
490 ELSEIF((LHAINPUT .GE. 31000) .AND. (LHAINPUT .LE. 32000)) THEN
492 LHANAME=LHAPATH(1:LHAPATHLEN)//'/Fermi2002_1000.LHpdf'
494 WRITE(LHAPRINT,5150) LHASET
498 ELSEIF((LHAINPUT .GE. 40000) .AND. (LHAINPUT .LE. 49999)) THEN
499 IF((LHAINPUT .GE. 40100) .AND. (LHAINPUT .LE. 40200)) THEN
501 LHANAME=LHAPATH(1:LHAPATHLEN)//'/Alekhin_100.LHpdf'
502 ELSEIF((LHAINPUT .GE. 41000) .AND. (LHAINPUT .LE. 41999)) THEN
504 LHANAME=LHAPATH(1:LHAPATHLEN)//'/Alekhin_1000.LHpdf'
505 ELSEIF((LHAINPUT .GE. 40350) .AND. (LHAINPUT .LE. 40367)) THEN
507 LHANAME=LHAPATH(1:LHAPATHLEN)//'/a02m_lo.LHgrid'
511 ELSEIF((LHAINPUT .GE. 40450) .AND. (LHAINPUT .LE. 40467)) THEN
513 LHANAME=LHAPATH(1:LHAPATHLEN)//'/a02m_nlo.LHgrid'
517 ELSEIF((LHAINPUT .GE. 40550) .AND. (LHAINPUT .LE. 40567)) THEN
519 LHANAME=LHAPATH(1:LHAPATHLEN)//'/a02m_nnlo.LHgrid'
524 WRITE(LHAPRINT,5150) LHASET
528 ELSEIF((LHAINPUT .GE. 50000) .AND. (LHAINPUT .LE. 59999)) THEN
529 IF((LHAINPUT .GE. 50100) .AND. (LHAINPUT .LE. 50200)) THEN
531 LHANAME=LHAPATH(1:LHAPATHLEN)//'/Botje_100.LHpdf'
532 ELSEIF((LHAINPUT .GE. 51000) .AND. (LHAINPUT .LE. 51999)) THEN
534 LHANAME=LHAPATH(1:LHAPATHLEN)//'/Botje_1000.LHpdf'
536 WRITE(LHAPRINT,5150) LHASET
540 ELSEIF((LHAINPUT .GE. 60000) .AND. (LHAINPUT .LE. 69999)) THEN
543 IF((LHAINPUT .GE. 60000) .AND. (LHAINPUT .LE. 60022)) THEN
545 LHANAME=LHAPATH(1:LHAPATHLEN)//'/ZEUS2002_TR.LHpdf'
546 ELSEIF((LHAINPUT .GE. 60100) .AND. (LHAINPUT .LE. 60122)) THEN
548 LHANAME=LHAPATH(1:LHAPATHLEN)//'/ZEUS2002_ZM.LHpdf'
549 ELSEIF((LHAINPUT .GE. 60200) .AND. (LHAINPUT .LE. 60222)) THEN
551 LHANAME=LHAPATH(1:LHAPATHLEN)//'/ZEUS2002_FF.LHpdf'
552 ELSEIF((LHAINPUT .GE. 60300) .AND. (LHAINPUT .LE. 60322)) THEN
554 LHANAME=LHAPATH(1:LHAPATHLEN)//'/ZEUS2005_ZJ.LHpdf'
556 WRITE(LHAPRINT,5150) LHASET
560 ELSEIF((LHAINPUT .GE. 70000) .AND. (LHAINPUT .LE. 79999)) THEN
564 IF((LHAINPUT .GE. 70050) .AND. (LHAINPUT .LE. 70050)) THEN
566 LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000ms.LHgrid'
567 ELSEIF((LHAINPUT .GE. 70051) .AND. (LHAINPUT .LE. 70070)) THEN
569 LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000msE.LHgrid'
570 ELSEIF((LHAINPUT .GE. 70150) .AND. (LHAINPUT .LE. 70150)) THEN
572 LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000dis.LHgrid'
573 ELSEIF((LHAINPUT .GE. 70151) .AND. (LHAINPUT .LE. 70170)) THEN
575 LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000disE.LHgrid'
576 ELSEIF((LHAINPUT .GE. 70250) .AND. (LHAINPUT .LE. 70250)) THEN
578 LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000lo.LHgrid'
579 ELSEIF((LHAINPUT .GE. 70251) .AND. (LHAINPUT .LE. 70270)) THEN
581 LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000loE.LHgrid'
582 c tempoararily removed on returning to original H!2000 files
583 c ELSEIF((LHAINPUT .GE. 70350) .AND. (LHAINPUT .LE. 70350)) THEN
585 c LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000lo2.LHgrid'
586 c ELSEIF((LHAINPUT .GE. 70351) .AND. (LHAINPUT .LE. 70370)) THEN
588 c LHANAME=LHAPATH(1:LHAPATHLEN)//'/H12000lo2E.LHgrid'
590 WRITE(LHAPRINT,5150) LHASET
594 ELSEIF((LHAINPUT .GE. 80000) .AND. (LHAINPUT .LE. 89999)) THEN
598 IF((LHAINPUT .GE. 80050) .AND. (LHAINPUT .LE. 80051)) THEN
600 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRV98nlo.LHgrid'
601 ELSEIF((LHAINPUT .GE. 80060) .AND. (LHAINPUT .LE. 80060)) THEN
603 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRV98lo.LHgrid'
605 WRITE(LHAPRINT,5150) LHASET
612 ELSEIF((LHAINPUT .GE. 210) .AND. (LHAINPUT .LE. 212)) THEN
617 IF((LHAINPUT .GE. 210) .AND. (LHAINPUT .LE. 212)) THEN
619 LHANAME=LHAPATH(1:LHAPATHLEN)//'/OWPI.LHgrid'
621 WRITE(LHAPRINT,5150) LHASET
625 ELSEIF((LHAINPUT .GE. 230) .AND. (LHAINPUT .LE. 233)) THEN
630 IF((LHAINPUT .GE. 230) .AND. (LHAINPUT .LE. 233)) THEN
632 LHANAME=LHAPATH(1:LHAPATHLEN)//'/SMRSPI.LHgrid'
634 WRITE(LHAPRINT,5150) LHASET
638 ELSEIF((LHAINPUT .GE. 250) .AND. (LHAINPUT .LE. 252)) THEN
642 IF((LHAINPUT .GE. 250) .AND. (LHAINPUT .LE. 251)) THEN
645 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRVPI1.LHgrid'
646 ELSEIF((LHAINPUT .GE. 252) .AND. (LHAINPUT .LE. 252)) THEN
649 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRVPI0.LHgrid'
651 WRITE(LHAPRINT,5150) LHASET
655 ELSEIF((LHAINPUT .GE. 260) .AND. (LHAINPUT .LE. 263)) THEN
660 IF((LHAINPUT .GE. 260) .AND. (LHAINPUT .LE. 263)) THEN
662 LHANAME=LHAPATH(1:LHAPATHLEN)//'/ABFKWPI.LHgrid'
664 WRITE(LHAPRINT,5150) LHASET
671 ELSEIF((LHAINPUT .GE. 310) .AND. (LHAINPUT .LE. 312)) THEN
676 IF((LHAINPUT .GE. 310) .AND. (LHAINPUT .LE. 311)) THEN
678 LHANAME=LHAPATH(1:LHAPATHLEN)//'/DOG0.LHgrid'
679 ELSEIF((LHAINPUT .GE. 312) .AND. (LHAINPUT .LE. 312)) THEN
681 LHANAME=LHAPATH(1:LHAPATHLEN)//'/DOG1.LHgrid'
683 WRITE(LHAPRINT,5150) LHASET
687 ELSEIF((LHAINPUT .GE. 320) .AND. (LHAINPUT .LE. 324)) THEN
691 IF((LHAINPUT .GE. 320) .AND. (LHAINPUT .LE. 321)) THEN
695 LHANAME=LHAPATH(1:LHAPATHLEN)//'/DGG.LHgrid'
696 ELSEIF((LHAINPUT .GE. 322) .AND. (LHAINPUT .LE. 322)) THEN
700 LHANAME=LHAPATH(1:LHAPATHLEN)//'/DGG.LHgrid'
701 ELSEIF((LHAINPUT .GE. 323) .AND. (LHAINPUT .LE. 323)) THEN
705 LHANAME=LHAPATH(1:LHAPATHLEN)//'/DGG.LHgrid'
706 ELSEIF((LHAINPUT .GE. 324) .AND. (LHAINPUT .LE. 324)) THEN
710 LHANAME=LHAPATH(1:LHAPATHLEN)//'/DGG.LHgrid'
712 WRITE(LHAPRINT,5150) LHASET
716 ELSEIF((LHAINPUT .GE. 330) .AND. (LHAINPUT .LE. 334)) THEN
722 IF((LHAINPUT .GE. 330) .AND. (LHAINPUT .LE. 332)) THEN
724 LHANAME=LHAPATH(1:LHAPATHLEN)//'/LACG.LHgrid'
725 ELSEIF((LHAINPUT .GE. 333) .AND. (LHAINPUT .LE. 333)) THEN
728 LHANAME=LHAPATH(1:LHAPATHLEN)//'/LACG.LHgrid'
729 ELSEIF((LHAINPUT .GE. 334) .AND. (LHAINPUT .LE. 334)) THEN
732 LHANAME=LHAPATH(1:LHAPATHLEN)//'/LACG.LHgrid'
734 WRITE(LHAPRINT,5150) LHASET
737 C...GSG/GSG96-G Family
738 ELSEIF((LHAINPUT .GE. 340) .AND. (LHAINPUT .LE. 345)) THEN
743 IF((LHAINPUT .GE. 340) .AND. (LHAINPUT .LE. 341)) THEN
745 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GSG1.LHgrid'
746 ELSEIF((LHAINPUT .GE. 342) .AND. (LHAINPUT .LE. 343)) THEN
748 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GSG0.LHgrid'
749 ELSEIF((LHAINPUT .GE. 344) .AND. (LHAINPUT .LE. 344)) THEN
751 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GSG961.LHgrid'
752 ELSEIF((LHAINPUT .GE. 345) .AND. (LHAINPUT .LE. 345)) THEN
754 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GSG960.LHgrid'
756 WRITE(LHAPRINT,5150) LHASET
760 ELSEIF((LHAINPUT .GE. 350) .AND. (LHAINPUT .LE. 354)) THEN
765 IF((LHAINPUT .GE. 350) .AND. (LHAINPUT .LE. 352)) THEN
767 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRVG1.LHgrid'
768 ELSEIF((LHAINPUT .GE. 353) .AND. (LHAINPUT .LE. 353)) THEN
771 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRVG0.LHgrid'
772 ELSEIF((LHAINPUT .GE. 354) .AND. (LHAINPUT .LE. 354)) THEN
776 LHANAME=LHAPATH(1:LHAPATHLEN)//'/GRVG0.LHgrid'
778 WRITE(LHAPRINT,5150) LHASET
782 ELSEIF((LHAINPUT .GE. 360) .AND. (LHAINPUT .LE. 363)) THEN
787 IF((LHAINPUT .GE. 360) .AND. (LHAINPUT .LE. 363)) THEN
789 LHANAME=LHAPATH(1:LHAPATHLEN)//'/ACFGPG.LHgrid'
791 WRITE(LHAPRINT,5150) LHASET
795 ELSEIF((LHAINPUT .GE. 380) .AND. (LHAINPUT .LE. 386)) THEN
800 IF((LHAINPUT .GE. 380) .AND. (LHAINPUT .LE. 386)) THEN
802 LHANAME=LHAPATH(1:LHAPATHLEN)//'/WHITG.LHgrid'
804 WRITE(LHAPRINT,5150) LHASET
808 ELSEIF((LHAINPUT .GE. 390) .AND. (LHAINPUT .LE. 398)) THEN
813 IF((LHAINPUT .GE. 390) .AND. (LHAINPUT .LE. 392)) THEN
816 LHANAME=LHAPATH(1:LHAPATHLEN)//'/SASG.LHgrid'
817 ELSEIF((LHAINPUT .GE. 393) .AND. (LHAINPUT .LE. 394)) THEN
820 LHANAME=LHAPATH(1:LHAPATHLEN)//'/SASG.LHgrid'
821 ELSEIF((LHAINPUT .GE. 395) .AND. (LHAINPUT .LE. 396)) THEN
824 LHANAME=LHAPATH(1:LHAPATHLEN)//'/SASG.LHgrid'
825 ELSEIF((LHAINPUT .GE. 397) .AND. (LHAINPUT .LE. 398)) THEN
828 LHANAME=LHAPATH(1:LHAPATHLEN)//'/SASG.LHgrid'
830 WRITE(LHAPRINT,5150) LHASET
833 C...Unknown Family ?! Giving up
835 WRITE(LHAPRINT,5150) LHASET
839 LHAMEMB=LHAINPUT-LHASET
840 c....Now work out if we have already called this set/member
843 if(lhaname.eq.lhanames(j).and.
844 + lhamemb.eq.lhamembers(j)) then
850 if(nsets.gt.nmxset) then
851 if(LHASILENT.ne.1) then
852 print *,'WARNING:too many sets initialised'
853 print *,'overwriting from set 1 again'
859 lhanames(iset)=lhaname
860 lhanumbers(iset)=lhainput
861 lhamembers(iset)=lhamemb
866 CALL INITPDFSETM(iset,LHANAME)
867 CALL NUMBERPDFM(iset,LHAALLMEM)
868 IF(LHASILENT .NE. 1) THEN
870 WRITE(LHAPRINT,5152) LHANAME
871 WRITE(LHAPRINT,5153) LHAALLMEM
874 IF ((LHAMEMB.LT.0) .OR. (LHAMEMB.GT.LHAALLMEM)) THEN
875 WRITE(LHAPRINT,5155) LHAMEMB
876 WRITE(LHAPRINT,5156) LHAALLMEM
880 c print *,'calling initpdf',lhamemb
881 c print *,'calling initpdfm ',iset,lhaname,lhamemb
882 c print *,'LHAGLUE .... initializing set,member ',iset,lhamemb
883 CALL INITPDFM(iset,LHAMEMB)
885 c... the rest is done every time pdfset is called
886 c print *,'setting nset to:',iset
888 call setnmem(iset,lhamemb)
893 call GetLam4M(iset,LHAMEMB,qcdl4)
894 call GetLam5M(iset,LHAMEMB,qcdl5)
897 alphasLHA = alphasPDFM(iset,QMZ)
899 > WRITE(LHAPRINT,5158) alphasLHA
901 IF(LHAPARM(17).EQ.'LHAPDF') THEN
902 NPTYPEPDFL = 1 ! Proton PDFs
907 > WRITE(LHAPRINT,5159) QCDL4, QCDL5
909 NPTYPEPDFL = 1 ! Proton PDFs
914 IF(PARM(1).EQ.'NPTYPE') THEN ! PYTHIA
925 C...Formats for initialization information.
926 5150 FORMAT(1X,'WRONG LHAPDF set number =',I12,' given! STOP EXE!')
927 5151 FORMAT(1X,'==============================================')
928 5152 FORMAT(1X,'PDFset name ',A80)
929 5153 FORMAT(1X,'with ',I10,' members')
930 5154 FORMAT(1X,'==== initialized. ===========================')
931 5155 FORMAT(1X,'LHAPDF problem => YOU asked for member = ',I10)
932 5156 FORMAT(1X,'Valid range is: 0 - ',I10,' Execution stopped.')
933 5157 FORMAT(1X,'Number of flavors for PDF is:',I4)
934 5158 FORMAT(1X,'Strong coupling at Mz for PDF is:',F9.5)
935 5159 FORMAT(1X,'Will use for PYTHIA QCDL4, QCDL5:',2F9.5)
940 c********************************************************************
942 c -- copy of PDFLIB to use the eks98 nuclear correction factors
944 SUBROUTINE STRUCTA(X,Q,A,UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU)
945 IMPLICIT DOUBLE PRECISION (A-H,O-Z)
946 CALL EKS98(X,Q,A,RUV,RDV,RU,RD,RS,RC,RB,RT,RG)
947 CALL STRUCTM(X,Q,UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU)
960 C*********************************************************************
963 C...Gives parton distributions according to the LHAPDF interface.
964 C...Two evolution codes used:
965 C... EVLCTEQ for CTEQ PDF sets
966 C... QCDNUM for Other PDF sets
968 C...Author: Dimitri Bourilkov bourilkov@mailaps.org
970 C...v4.0 21-Mar-2005 Photon/pion/new p PDFs, updated for LHAPDF v4
973 C...interface to LHAPDF library
975 SUBROUTINE STRUCTM(DX,DQ,UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU)
977 C...Double precision and integer declarations.
978 IMPLICIT DOUBLE PRECISION(A-H, O-Z)
979 IMPLICIT INTEGER(I-N)
980 include 'parmsetup.inc'
981 C...Interface to LHAPDFLIB.
982 include 'pathsetup.inc'
983 c CHARACTER*172 LHANAME
984 INTEGER LHASET, LHAMEMB
985 COMMON/LHAPDF/LHANAME, LHASET, LHAMEMB
987 c...added next 2 lines for structp fix
988 integer LHAMEMBERS(nmxset),LHANUMBERS(nmxset)
989 common/LHASETS/LHANAMES,LHANUMBERS,LHAMEMBERS,nsets
991 DOUBLE PRECISION QCDLHA4, QCDLHA5
993 COMMON/LHAPDFR/QCDLHA4, QCDLHA5, NFLLHA
995 CHARACTER*20 LHAPARM(20)
996 DOUBLE PRECISION LHAVALUE(20)
997 COMMON/LHACONTROL/LHAPARM,LHAVALUE
1000 COMMON/LHAPDFE/LHAEXTRP
1002 DOUBLE PRECISION XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
1003 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
1004 COMMON/LHAGLSTA/ XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
1005 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
1007 C...Interface to PDFLIB.
1008 COMMON/W50513/XMIN,XMAX,Q2MIN,Q2MAX
1010 DOUBLE PRECISION XMIN,XMAX,Q2MIN,Q2MAX
1012 DOUBLE PRECISION UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU
1013 DOUBLE PRECISION DX,DQ,X,Q,F(-6:6)
1019 IF(LHAPARM(16).NE.'NOSTAT') THEN
1020 TOTNUM = TOTNUM+1.D0
1021 IF(X .LT. XMIN) XMINNUM = XMINNUM+1.D0
1022 IF(X .GT. XMAX) XMAXNUM = XMAXNUM+1.D0
1023 IF(Q2 .LT. Q2MIN) Q2MINNUM = Q2MINNUM+1.D0
1024 IF(Q2 .GT. Q2MAX) Q2MAXNUM = Q2MAXNUM+1.D0
1027 C...Range of validity e.g. 10^-6 < x < 1, Q2MIN < Q^2 extended by
1028 C...freezing x*f(x,Q2) at borders.
1029 IF(LHAEXTRP .NE. 1) THEN ! safe mode == "freeze"
1030 XIN=MAX(XMIN,MIN(XMAX,X))
1031 Q=SQRT(MAX(0D0,Q2MIN,MIN(Q2MAX,Q2)))
1032 ELSE ! adventurous mode == OWN RISK !
1037 c print *,'calling evolvepdfm:',iset
1039 C...fix to allow STRUCTM to work for photon PDFs (Herwig does this)
1040 C...set P2 = 0.0d0 and IP2 = 0
1041 if(LHANUMBERS(iset).ge.300.and.LHANUMBERS(iset).le.399) then
1044 CALL EVOLVEPDFPM(iset,XIN,Q,P2,IP2,F)
1046 CALL EVOLVEPDFM(iset,XIN,Q,F)
1061 C*********************************************************************
1064 C...Gives parton distributions according to the LHAPDF interface.
1065 C...Used for photons.
1067 C...v4.0 21-Mar-2005 Photon/pion/new p PDFs, updated for LHAPDF v4
1069 C...interface to LHAPDF library
1072 > (DX,DQ2,P2,IP2,UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU)
1074 C...Double precision and integer declarations.
1075 IMPLICIT DOUBLE PRECISION(A-H, O-Z)
1076 IMPLICIT INTEGER(I-N)
1077 include 'parmsetup.inc'
1078 C...Interface to LHAPDFLIB.
1079 include 'pathsetup.inc'
1080 c CHARACTER*172 LHANAME
1081 INTEGER LHASET, LHAMEMB
1082 COMMON/LHAPDF/LHANAME, LHASET, LHAMEMB
1084 DOUBLE PRECISION QCDLHA4, QCDLHA5
1086 COMMON/LHAPDFR/QCDLHA4, QCDLHA5, NFLLHA
1088 CHARACTER*20 LHAPARM(20)
1089 DOUBLE PRECISION LHAVALUE(20)
1090 COMMON/LHACONTROL/LHAPARM,LHAVALUE
1093 COMMON/LHAPDFE/LHAEXTRP
1095 DOUBLE PRECISION XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
1096 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
1097 COMMON/LHAGLSTA/ XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
1098 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
1100 C...Interface to PDFLIB.
1101 COMMON/W50513/XMIN,XMAX,Q2MIN,Q2MAX
1103 DOUBLE PRECISION XMIN,XMAX,Q2MIN,Q2MAX
1105 DOUBLE PRECISION UPV,DNV,USEA,DSEA,STR,CHM,BOT,TOP,GLU
1106 DOUBLE PRECISION DX,DQ2,Q2,X,Q,F(-6:6)
1111 IF(LHAPARM(16).NE.'NOSTAT') THEN
1112 TOTNUP = TOTNUP+1.D0
1113 IF(X .LT. XMIN) XMINNUP = XMINNUP+1.D0
1114 IF(X .GT. XMAX) XMAXNUP = XMAXNUP+1.D0
1115 IF(Q2 .LT. Q2MIN) Q2MINNUP = Q2MINNUP+1.D0
1116 IF(Q2 .GT. Q2MAX) Q2MAXNUP = Q2MAXNUP+1.D0
1119 C...Range of validity e.g. 10^-6 < x < 1, Q2MIN < Q^2 extended by
1120 C...freezing x*f(x,Q2) at borders.
1122 IF(LHAEXTRP .NE. 1) THEN ! safe mode == "freeze"
1123 XIN=MAX(XMIN,MIN(XMAX,X))
1124 Q=SQRT(MAX(0D0,Q2MIN,MIN(Q2MAX,Q2)))
1125 ELSE ! adventurous mode == OWN RISK !
1129 CALL EVOLVEPDFPM(iset,XIN,Q,P2,IP2,F)
1144 C*********************************************************************
1147 C...For statistics ON structure functions (under/over-flows)
1149 C...Author: Dimitri Bourilkov bourilkov@mailaps.org
1152 C...first introduced in v4.0 28-Apr-2005
1157 C...Double precision and integer declarations.
1158 IMPLICIT DOUBLE PRECISION(A-H, O-Z)
1159 IMPLICIT INTEGER(I-N)
1160 C...Interface to LHAPDFLIB.
1161 DOUBLE PRECISION XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
1162 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
1163 COMMON/LHAGLSTA/ XMINNUM,XMAXNUM,Q2MINNUM,Q2MAXNUM,TOTNUM,
1164 > XMINNUP,XMAXNUP,Q2MINNUP,Q2MAXNUP,TOTNUP
1168 PRINT *,'===== PDFSTA statistics for PDF under/over-flows ===='
1170 PRINT *,'====== STRUCTM statistics for nucleon/pion PDFs ====='
1172 PRINT *,' total # of calls ',TOTNUM
1173 IF(TOTNUM .GT. 0.D0) THEN
1174 PERCBELOW = 100.D0*XMINNUM/TOTNUM
1175 PERCABOVE = 100.D0*XMAXNUM/TOTNUM
1176 PRINT *,' X below PDF min ',XMINNUM,' or ',PERCBELOW, ' %'
1177 PRINT *,' X above PDF max ',XMAXNUM,' or ',PERCABOVE, ' %'
1178 PERCBELOW = 100.D0*Q2MINNUM/TOTNUM
1179 PERCABOVE = 100.D0*Q2MAXNUM/TOTNUM
1180 PRINT *,' Q2 below PDF min ',Q2MINNUM,' or ',PERCBELOW, ' %'
1181 PRINT *,' Q2 above PDF max ',Q2MAXNUM,' or ',PERCABOVE, ' %'
1184 PRINT *,'========= STRUCTP statistics for photon PDFs ========'
1186 PRINT *,' total # of calls ',TOTNUP
1187 IF(TOTNUP .GT. 0.D0) THEN
1188 PERCBELOW = 100.D0*XMINNUP/TOTNUP
1189 PERCABOVE = 100.D0*XMAXNUP/TOTNUP
1190 PRINT *,' X below PDF min ',XMINNUP,' or ',PERCBELOW, ' %'
1191 PRINT *,' X above PDF max ',XMAXNUP,' or ',PERCABOVE, ' %'
1192 PERCBELOW = 100.D0*Q2MINNUP/TOTNUP
1193 PERCABOVE = 100.D0*Q2MAXNUP/TOTNUP
1194 PRINT *,' Q2 below PDF min ',Q2MINNUP,' or ',PERCBELOW, ' %'
1195 PRINT *,' Q2 above PDF max ',Q2MAXNUP,' or ',PERCABOVE, ' %'
1201 **********************************************************************
1203 * $Id: lhaglue.f 209 2007-11-16 15:14:45Z whalley $
1206 * Revision 1.7 2005/12/02 14:50:54 whalley
1207 * Changes for new CTEQ code/AB sets
1209 * Revision 1.6 2005/10/18 15:35:48 whalley
1210 * fix to allow LHAPATH to be user defined as well as lhapdf-config
1212 * Revision 1.5 2005/10/18 11:47:48 whalley
1213 * Change to only set LHAPATH once per run
1215 * Revision 1.1.1.2 1996/10/30 08:29:06 cernlib
1218 * Revision 1.1.1.1 1996/04/12 15:29:26 plothow
1222 SUBROUTINE PFTOPDG(DX,DSCALE,DXPDF)
1224 Cinclude "pdf/expdp.inc"
1226 + DX,DSCALE,DUPV,DDNV,DUSEA,DDSEA,DSTR,DCHM,DBOT,DTOP,DGL,
1228 C... call STRUCTM in PDFLIB to get flavour content
1229 CALL STRUCTM(DX,DSCALE,
1230 + DUPV,DDNV,DUSEA,DDSEA,DSTR,DCHM,DBOT,DTOP,DGL)
1231 C... convert flavour convention of PDFLIB to PDG convention
1233 DXPDF(1) = DDNV + DDSEA
1234 DXPDF(2) = DUPV + DUSEA
1248 ****************************************************************************
1249 subroutine setPDFpath(pathname)
1250 implicit real*8 (A-H,O-Z)
1251 include 'parmsetup.inc'
1252 character*(*) pathname
1253 include 'pathsetup.inc'
1254 c character*132 lhapath
1255 common/LHAPDFC/lhapath
1256 character*20 lhaparm(20)
1258 common/LHACONTROL/lhaparm,lhavalue
1259 lhaparm(20) = 'LHAPATH'
1260 c do j=1,lnblnk(lhapath)
1261 do j=1,LEN_TRIM(lhapath)
1267 ***********************************************************************
1268 subroutine lhaset(lhaparm2,lhavalue2)
1269 implicit real*8 (a-h,o-z)
1270 character*20 lhaparm(20),lhaparm2(20)
1271 real*8 lhavalue(20),lhavalue2(20)
1272 common/LHACONTROL/lhaparm,lhavalue
1274 lhaparm(j)=lhaparm2(j)
1275 lhavalue(j)=lhavalue2(j)
1279 ******************************************************************
1280 subroutine setlhaparm(lparm)
1281 implicit real*8 (a-h,o-z)
1282 c character*(*) lparm,parm
1284 character*20 lhaparm(20)
1286 common/LHACONTROL/lhaparm,lhavalue
1287 if(lparm.eq.'NOSTAT') then
1288 lhaparm(16)='NOSTAT'
1289 else if (lparm.eq.'16') then
1291 else if (lparm.eq.'LHAPDF') then
1292 lhaparm(17)='LHAPDF'
1293 else if (lparm.eq.'17') then
1295 else if (lparm.eq.'EXTRAPOLATE') then
1296 lhaparm(18)='EXTRAPOLATE'
1297 else if (lparm.eq.'18') then
1299 else if (lparm.eq.'SILENT') then
1300 lhaparm(19)='SILENT'
1301 else if (lparm.eq.'LOWKEY') then
1302 lhaparm(19)='LOWKEY'
1303 else if (lparm.eq.'19') then
1306 print *,'WARNING from SetLHAPARM - value',lparm,'
1311 entry getlhaparm(nparm,lparm)
1312 lparm = lhaparm(nparm)
1316 ***************************************************************