1 subroutine MRST2006evolve(x,Q,pdf)
3 c implicit real*8 (a-h,o-z)
4 include 'parmsetup.inc'
5 character*16 name(nmxset)
6 integer i,f,nhess,nx,nq,np,nqc0,nqb0,nqc,nqb,n,m,io
7 double precision x,q,xmin,xmax,qsqmin,qsqmax,emc2,emb2,eps,
8 & dummy,qsq,xlog,qsqlog,res,ExtrapLHAPDF,InterpLHAPDF
9 double precision emc2minuseps,emb2minuseps
10 double precision xsave,q2save
11 double precision pdf(-6:6),alfas,Qalfa,Q2fit,eorder
12 double precision upv,dnv,glu,usea,str,dsea,sbar,chm,bot
13 integer iset,nflav,mem,nset,nmem(nmxset),ndef(nmxset),mmem,imem
14 common/NAME/name,nmem,ndef,mmem
15 parameter(nx=64,nq=48,np=9,nqc0=4,nqb0=14,nqc=nq-nqc0,nqb=nq-nqb0)
16 parameter(xmin=1d-6,xmax=1d0,qsqmin=1d0,qsqmax=1d9,eps=1d-6)
17 parameter(nhess=30,emc2=1.43d0**2,emb2=4.3d0**2) ! MRST 2006 NNLO
18 parameter(emc2minuseps=emc2-eps,emb2minuseps=emb2-eps)
19 double precision f1(nx,nq),f2(nx,nq),f3(nx,nq),f4(nx,nq),
20 & f5(nx,nq),f6(nx,nq),f7(nx,nq),f8(nx,nq),fc(nx,nqc),
21 & fb(nx,nqb),f9(nx,nq)
22 double precision qq(nq),xx(nx),cc1(0:nhess,nx,nq,4,4),
23 & cc2(0:nhess,nx,nq,4,4),cc3(0:nhess,nx,nq,4,4),
24 & cc4(0:nhess,nx,nq,4,4),cc6(0:nhess,nx,nq,4,4),
25 & cc8(0:nhess,nx,nq,4,4),cc9(0:nhess,nx,nq,4,4),
26 & ccc(0:nhess,nx,nqc,4,4),ccb(0:nhess,nx,nqb,4,4)
27 double precision xxl(nx),qql(nq),qqlc(nqc),qqlb(nqb)
28 data xx/1d-6,2d-6,4d-6,6d-6,8d-6,
29 & 1d-5,2d-5,4d-5,6d-5,8d-5,
30 & 1d-4,2d-4,4d-4,6d-4,8d-4,
31 & 1d-3,2d-3,4d-3,6d-3,8d-3,
32 & 1d-2,1.4d-2,2d-2,3d-2,4d-2,6d-2,8d-2,
33 & .1d0,.125d0,.15d0,.175d0,.2d0,.225d0,.25d0,.275d0,
34 & .3d0,.325d0,.35d0,.375d0,.4d0,.425d0,.45d0,.475d0,
35 & .5d0,.525d0,.55d0,.575d0,.6d0,.625d0,.65d0,.675d0,
36 & .7d0,.725d0,.75d0,.775d0,.8d0,.825d0,.85d0,.875d0,
37 & .9d0,.925d0,.95d0,.975d0,1d0/
39 & 1.25d0,1.5d0,emc2minuseps,emc2,2.5d0,3.2d0,4d0,5d0,6.4d0,8d0,
40 & 1d1,1.2d1,emb2minuseps,emb2,2.6d1,4d1,6.4d1,1d2,
41 & 1.6d2,2.4d2,4d2,6.4d2,1d3,1.8d3,3.2d3,5.6d3,1d4,
42 & 1.8d4,3.2d4,5.6d4,1d5,1.8d5,3.2d5,5.6d5,1d6,
43 & 1.8d6,3.2d6,5.6d6,1d7,1.8d7,3.2d7,5.6d7,1d8,
44 & 1.8d8,3.2d8,5.6d8,1d9/
52 if (qsq.gt.qq(nqc0).and.qsq.lt.qq(nqc0+1)) qsq = qq(nqc0)
53 if (qsq.gt.qq(nqb0).and.qsq.lt.qq(nqb0+1)) qsq = qq(nqb0)
59 call getnmem(iset,imem)
71 if (x.le.0.d0.or.x.gt.xmax.or.qsq.lt.qsqmin) then
72 print *,"Error in GetOnePDF: x,qsq = ",x,qsq
73 else if (x.lt.xmin.or.qsq.gt.qsqmax) then ! extrapolate
74 print *, "Warning in GetOnePDF, extrapolating: f = ",f,
75 & ", x = ",x,", q = ",q
76 upv = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc1)
77 dnv = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc2)
78 glu = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc3)
79 usea = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc4)
80 str = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc6)
81 dsea = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc8)
82 sbar = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc9)
83 if (qsq.ge.emc2) then ! chm
84 chm = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nqc,xxl,qqlc,ccc)
86 if (qsq.ge.emb2) then ! bot
87 bot = ExtrapLHAPDF(imem,nhess,xlog,qsqlog,nx,nqb,xxl,qqlb,ccb)
90 upv = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc1)
91 dnv = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc2)
92 glu = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc3)
93 usea = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc4)
94 str = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc6)
95 dsea = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc8)
96 sbar = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nq,xxl,qql,cc9)
97 if (qsq.ge.emc2) then ! chm
98 chm = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nqc,xxl,qqlc,ccc)
100 if (qsq.ge.emb2) then ! bot
101 bot = InterpLHAPDF(imem,nhess,xlog,qsqlog,nx,nqb,xxl,qqlb,ccb)
125 C----------------------------------------------------------------------
127 entry MRST2006read(nset)
128 read(1,*)nmem(nset),ndef(nset)
137 read(1,50) f1(n,m),f2(n,m),f3(n,m),f4(n,m),
138 & f5(n,m),f7(n,m),f6(n,m),f8(n,m),f9(n,m)
139 c write(6,50) f1(n,m),f2(n,m),f3(n,m),f4(n,m),
140 c & f5(n,m),f7(n,m),f6(n,m),f8(n,m),f9(n,m)
141 C-- Notation:1=upv 2=dnv 3=glu 4=usea 5=chm 6=str 7=bot 8=dsea 9=sbar
144 c write(*,*)'PDF set ',i,' first element ',f1(1,1)
164 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f1,cc1)
165 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f2,cc2)
166 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f3,cc3)
167 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f4,cc4)
168 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f6,cc6)
169 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f8,cc8)
170 call InitialLHAPDF(i,nhess,nx,nq,nqc0,nqb0,xxl,qql,f9,cc9)
178 call InitialLHAPDF(i,nhess,nx,nqc,nqc0-nqc0,nqb0-nqc0,
187 call InitialLHAPDF(i,nhess,nx,nqb,nqc0-nqb0,nqb0-nqb0,
194 C----------------------------------------------------------------------
196 entry MRST2006alfa(nflav,alfas,Qalfa)
197 c print *,'MRST2006alfa',nflav,alfas,Qalfa
199 call alphamrs(nflav,alfas,Qalfa)
200 c print *,'MRST2006alfa',nflav,alfas,Qalfa
203 C----------------------------------------------------------------------
204 entry MRST2006init(Eorder,Q2fit)
206 C----------------------------------------------------------------------
208 entry MRST2006pdf(mem)
210 call setnmem(iset,mem)
215 C----------------------------------------------------------------------
217 subroutine InitialLHAPDF(i,nhess,nx,my,myc0,myb0,xx,yy,ff,cc)
219 integer nhess,i,nx,my,myc0,myb0,j,k,l,m,n
220 double precision xx(nx),yy(my),ff(nx,my),
221 & ff1(nx,my),ff2(nx,my),
222 & ff12(nx,my),yy0(4),yy1(4),yy2(4),yy12(4),z(16),
223 & cl(16),cc(0:nhess,nx,my,4,4),iwt(16,16),
224 & dx,dy,polderiv,d1,d2,d1d2,xxd
226 data iwt/1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,
227 & 0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,
228 & -3,0,0,3,0,0,0,0,-2,0,0,-1,0,0,0,0,
229 & 2,0,0,-2,0,0,0,0,1,0,0,1,0,0,0,0,
230 & 0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,
231 & 0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,
232 & 0,0,0,0,-3,0,0,3,0,0,0,0,-2,0,0,-1,
233 & 0,0,0,0,2,0,0,-2,0,0,0,0,1,0,0,1,
234 & -3,3,0,0,-2,-1,0,0,0,0,0,0,0,0,0,0,
235 & 0,0,0,0,0,0,0,0,-3,3,0,0,-2,-1,0,0,
236 & 9,-9,9,-9,6,3,-3,-6,6,-6,-3,3,4,2,1,2,
237 & -6,6,-6,6,-4,-2,2,4,-3,3,3,-3,-2,-1,-1,-2,
238 & 2,-2,0,0,1,1,0,0,0,0,0,0,0,0,0,0,
239 & 0,0,0,0,0,0,0,0,2,-2,0,0,1,1,0,0,
240 & -6,6,-6,6,-3,-3,3,3,-4,4,2,-2,-2,-2,-1,-1,
241 & 4,-4,4,-4,2,2,-2,-2,2,-2,-2,2,1,1,1,1/
245 ff1(1,m)=(ff(2,m)-ff(1,m))/dx
247 ff1(nx,m)=(ff(nx,m)-ff(nx-1,m))/dx
249 ff1(n,m)=polderiv(xx(n-1),xx(n),xx(n+1),ff(n-1,m),ff(n,m),
254 C-- Calculate the derivatives at qsq=emc2-eps,emc2,emb2-eps,emb2
255 C-- in a similar way as at the endpoints qsqmin and qsqmax.
258 if (m.eq.1.or.m.eq.myc0+1.or.m.eq.myb0+1) then
260 ff2(n,m)=(ff(n,m+1)-ff(n,m))/dy
261 else if (m.eq.my.or.m.eq.myc0.or.m.eq.myb0) then
263 ff2(n,m)=(ff(n,m)-ff(n,m-1))/dy
265 ff2(n,m)=polderiv(yy(m-1),yy(m),yy(m+1),ff(n,m-1),
273 ff12(1,m)=(ff2(2,m)-ff2(1,m))/dx
275 ff12(nx,m)=(ff2(nx,m)-ff2(nx-1,m))/dx
277 ff12(n,m)=polderiv(xx(n-1),xx(n),xx(n+1),ff2(n-1,m),
278 x ff2(n,m),ff2(n+1,m))
305 yy12(3)=ff12(n+1,m+1)
318 xxd=xxd+iwt(k,l)*z(k)
334 C----------------------------------------------------------------------
336 double precision function InterpLHAPDF(i,nhess,x,y,nx,my,xx,yy,
339 integer i,nx,my,nhess,locx,l,m,n
340 double precision xx(nx),yy(my),cc(0:nhess,nx,my,4,4),
346 t=(x-xx(n))/(xx(n+1)-xx(n))
347 u=(y-yy(m))/(yy(m+1)-yy(m))
351 z=t*z+((cc(i,n,m,l,4)*u+cc(i,n,m,l,3))*u
352 . +cc(i,n,m,l,2))*u+cc(i,n,m,l,1)
360 C----------------------------------------------------------------------
362 double precision function ExtrapLHAPDF(i,nhess,x,y,nx,my,xx,yy,
365 integer i,nx,my,nhess,locx,n,m
366 double precision xx(nx),yy(my),cc(0:nhess,nx,my,4,4),
367 & x,y,z,f0,f1,z0,z1,InterpLHAPDF
369 n=locx(xx,nx,x) ! 0: below xmin, nx: above xmax
370 m=locx(yy,my,y) ! 0: below qsqmin, my: above qsqmax
372 C-- If extrapolation in small x only:
373 if (n.eq.0.and.m.gt.0.and.m.lt.my) then
374 f0 = InterpLHAPDF(i,nhess,xx(1),y,nx,my,xx,yy,cc)
375 f1 = InterpLHAPDF(i,nhess,xx(2),y,nx,my,xx,yy,cc)
376 if (f0.gt.0.d0.and.f1.gt.0.d0) then
377 z = exp(log(f0)+(log(f1)-log(f0))/(xx(2)-xx(1))*(x-xx(1)))
379 z = f0+(f1-f0)/(xx(2)-xx(1))*(x-xx(1))
381 C-- If extrapolation into large q only:
382 else if (n.gt.0.and.m.eq.my) then
383 f0 = InterpLHAPDF(i,nhess,x,yy(my),nx,my,xx,yy,cc)
384 f1 = InterpLHAPDF(i,nhess,x,yy(my-1),nx,my,xx,yy,cc)
385 if (f0.gt.0.d0.and.f1.gt.0.d0) then
386 z = exp(log(f0)+(log(f0)-log(f1))/(yy(my)-yy(my-1))*
389 z = f0+(f0-f1)/(yy(my)-yy(my-1))*(y-yy(my))
391 C-- If extrapolation into large q AND small x:
392 else if (n.eq.0.and.m.eq.my) then
393 f0 = InterpLHAPDF(i,nhess,xx(1),yy(my),nx,my,xx,yy,cc)
394 f1 = InterpLHAPDF(i,nhess,xx(1),yy(my-1),nx,my,xx,yy,cc)
395 if (f0.gt.0.d0.and.f1.gt.0.d0) then
396 z0 = exp(log(f0)+(log(f0)-log(f1))/(yy(my)-yy(my-1))*
399 z0 = f0+(f0-f1)/(yy(my)-yy(my-1))*(y-yy(my))
401 f0 = InterpLHAPDF(i,nhess,xx(2),yy(my),nx,my,xx,yy,cc)
402 f1 = InterpLHAPDF(i,nhess,xx(2),yy(my-1),nx,my,xx,yy,cc)
403 if (f0.gt.0.d0.and.f1.gt.0.d0) then
404 z1 = exp(log(f0)+(log(f0)-log(f1))/(yy(my)-yy(my-1))*
407 z1 = f0+(f0-f1)/(yy(my)-yy(my-1))*(y-yy(my))
409 if (z0.gt.0.d0.and.z1.gt.0.d0) then
410 z = exp(log(z0)+(log(z1)-log(z0))/(xx(2)-xx(1))*(x-xx(1)))
412 z = z0+(z1-z0)/(xx(2)-xx(1))*(x-xx(1))
415 print *,"Error in ExtrapLHAPDF"
424 C----------------------------------------------------------------------