LHAPDF veraion 5.9.1
[u/mrichter/AliRoot.git] / LHAPDF / lhapdf-5.9.1 / src / wrapmrst.f
1 ! -*- F90 -*-
2
3
4       subroutine MRSTevolve(x,Q,pdf) 
5       implicit real*8(a-h,o-z) 
6       include 'parmsetup.inc' 
7       character*16 name(nmxset) 
8       integer nmem(nmxset),ndef(nmxset),mmem 
9 !      integer member(nmxset)                                           
10       integer nset,iset 
11       common/NAME/name,nmem,ndef,mmem 
12       parameter(nx=49,nq=37,np=8,nqc0=2,nqb0=11,nqc=35,nqb=26,          &
13      &nhess=30)                                                         
14       real*8 pdf(-6:6) 
15       double precision gridx(nmxgridx),gridq(nmxgridq)
16       integer ngridx,ngridq,jx,jq
17       real*8 f1(nx,nq)                                                  &
18      &,f2(nx,nq)                                                        &
19      &,f3(nx,nq)                                                        &
20      &,f4(nx,nq)                                                        &
21      &,f5(nx,nq)                                                        &
22      &,f6(nx,nq)                                                        &
23      &,f7(nx,nq)                                                        &
24      &,f8(nx,nq)                                                        &
25      &,fc(nx,nqc),fb(nx,nqb)                                            
26       real*8 qq(nq),xx(nx),                                             &
27      &cc1(0:nhess,nx,nq,4,4,nmxset),cc2(0:nhess,nx,nq,4,4,nmxset),      &
28      &cc3(0:nhess,nx,nq,4,4,nmxset),cc4(0:nhess,nx,nq,4,4,nmxset),      &
29      &cc6(0:nhess,nx,nq,4,4,nmxset),cc8(0:nhess,nx,nq,4,4,nmxset),      &
30      &ccc(0:nhess,nx,nqc,4,4,nmxset),ccb(0:nhess,nx,nqb,4,4,nmxset)     
31       real*8 xxl(nx),qql(nq),qqlc(nqc),qqlb(nqb) 
32       data xx/1d-5,2d-5,4d-5,6d-5,8d-5,                                 &
33      &              1d-4,2d-4,4d-4,6d-4,8d-4,                           &
34      &              1d-3,2d-3,4d-3,6d-3,8d-3,                           &
35      &              1d-2,1.4d-2,2d-2,3d-2,4d-2,6d-2,8d-2,               &
36      &           .1d0,.125d0,.15d0,.175d0,.2d0,.225d0,.25d0,.275d0,     &
37      &           .3d0,.325d0,.35d0,.375d0,.4d0,.425d0,.45d0,.475d0,     &
38      &           .5d0,.525d0,.55d0,.575d0,.6d0,.65d0,.7d0,.75d0,        &
39      &           .8d0,.9d0,1d0/                                         
40       data qq/1.25d0,1.5d0,2d0,2.5d0,3.2d0,4d0,5d0,6.4d0,8d0,1d1,       &
41      &        1.2d1,1.8d1,2.6d1,4d1,6.4d1,1d2,                          &
42      &        1.6d2,2.4d2,4d2,6.4d2,1d3,1.8d3,3.2d3,5.6d3,1d4,          &
43      &        1.8d4,3.2d4,5.6d4,1d5,1.8d5,3.2d5,5.6d5,1d6,              &
44      &        1.8d6,3.2d6,5.6d6,1d7/                                    
45       data xmin,xmax,qsqmin,qsqmax/1d-5,1d0,1.25d0,1d7/ 
46       save 
47                                                                         
48       xsave=x 
49       qsq = q*q 
50       q2save=qsq 
51                                                                         
52       xlog=dlog(x) 
53       qsqlog=dlog(qsq) 
54                                                                         
55       call getnset(iset) 
56 !      imem=member(iset)                                                
57       call getnmem(iset,imem) 
58                                                                         
59                                 ! G.W. 24/04/2008                       
60       if (qsq.lt.qsqmin) then 
61          qsqlog=dlog(1.01D0*qsqmin) 
62          call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                    &
63      &   cc1(0,1,1,1,1,iset),upv1)                                      
64          call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                    &
65      &         cc2(0,1,1,1,1,iset),dnv1)                                
66          call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                    &
67      &         cc3(0,1,1,1,1,iset),glu1)                                
68          call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                    &
69      &   cc4(0,1,1,1,1,iset),usea1)                                     
70          call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                    &
71      &   cc6(0,1,1,1,1,iset),str1)                                      
72          call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                    &
73      &   cc8(0,1,1,1,1,iset),dsea1)                                     
74          qsqlog=dlog(qsqmin) 
75       end if 
76       call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                       &
77      &cc1(0,1,1,1,1,iset),upv)                                          
78       call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                       &
79      &cc2(0,1,1,1,1,iset),dnv)                                          
80       call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                       &
81      &cc3(0,1,1,1,1,iset),glu)                                          
82       call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                       &
83      &cc4(0,1,1,1,1,iset),usea)                                         
84       call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                       &
85      &cc6(0,1,1,1,1,iset),str)                                          
86       call jeppe2(imem,xlog,qsqlog,nx,nq,xxl,qql,                       &
87      &cc8(0,1,1,1,1,iset),dsea)                                         
88                                 ! G.W. 24/04/2008                       
89       if (qsq.lt.qsqmin) then 
90 !--   Calculate the anomalous dimension, dlog(xf)/dlog(qsq),            
91 !--   evaluated at qsqmin.  Then extrapolate the PDFs to low            
92 !--   qsq < qsqmin by interpolating the anomalous dimenion between      
93 !--   the value at qsqmin and a value of 1 for qsq << qsqmin.           
94 !--   If value of PDF at qsqmin is very small, just set                 
95 !--   anomalous dimension to 1 to prevent rounding errors.              
96          if (abs(upv).lt.1.D-4)  then
97              anom = 1.D0
98          else  
99              anom = (upv1-upv)/upv/0.01D0 
100          endif   
101          upv = upv*(qsq/qsqmin)**(anom*qsq/qsqmin+1.D0-qsq/qsqmin)
102          if (abs(dnv).lt.1.D-4) then
103              anom = 1.D0
104          else 
105              anom = (dnv1-dnv)/dnv/0.01D0 
106          endif
107          dnv = dnv*(qsq/qsqmin)**(anom*qsq/qsqmin+1.D0-qsq/qsqmin)
108          if (abs(glu).lt.1.D-4) then
109              anom = 1.D0
110          else 
111              anom = (glu1-glu)/glu/0.01D0  
112          endif
113          glu = glu*(qsq/qsqmin)**(anom*qsq/qsqmin+1.D0-qsq/qsqmin) 
114          if (abs(usea).lt.1.D-4) then
115              anom = 1.D0 
116          else
117              anom = (usea1-usea)/usea/0.01D0
118          endif 
119          usea = usea*(qsq/qsqmin)**(anom*qsq/qsqmin+1.D0-qsq/qsqmin) 
120          if (abs(str).lt.1.D-4) then
121              anom = 1.D0 
122          else
123              anom = (str1-str)/str/0.01D0 
124          endif
125          str = str*(qsq/qsqmin)**(anom*qsq/qsqmin+1.D0-qsq/qsqmin) 
126          if (abs(dsea).lt.1.D-4) then
127              anom = 1.D0 
128          else
129              anom = (dsea1-dsea)/dsea/0.01D0 
130          endif
131          dsea = dsea*(qsq/qsqmin)**(anom*qsq/qsqmin+1.D0-qsq/qsqmin) 
132       end if 
133                                                                         
134                                                                         
135                                                                         
136       chm=0.d0 
137       if(qsq.gt.emc2) then 
138       call jeppe2(imem,xlog,qsqlog,nx,nqc,xxl,qqlc,                     &
139      &ccc(0,1,1,1,1,iset),chm)                                          
140       endif 
141                                                                         
142       bot=0.d0 
143       if(qsq.gt.emb2) then 
144       call jeppe2(imem,xlog,qsqlog,nx,nqb,xxl,qqlb,                     &
145      &ccb(0,1,1,1,1,iset),bot)                                          
146       endif 
147                                                                         
148       pdf(0)  = glu 
149       pdf(1)  = dnv+dsea 
150       pdf(-1) = dsea 
151       pdf(2)  = upv+usea 
152       pdf(-2) = usea 
153       pdf(3)  = str 
154       pdf(-3) = str 
155       pdf(4)  = chm 
156       pdf(-4) = chm 
157       pdf(5)  = bot 
158       pdf(-5) = bot 
159       pdf(6)  = 0.0d0 
160       pdf(-6) = 0.0d0 
161                                                                         
162       x=xsave 
163       qsq=q2save 
164       return 
165 !                                                                       
166       entry MRSTgetgrid(nset,ngridx,ngridq,gridx,gridq)
167       do jx=1,nx
168           gridx(jx)=xx(jx)
169       enddo
170       do jq=1,nq
171           gridq(jq)=qq(jq)
172       enddo
173       ngridx=nx
174       ngridq=nq        
175       return
176       
177       entry MRSTread(nset) 
178       read(1,*)nmem(nset),ndef(nset) 
179 !      print *,nmem(nset),ndef(nset)                                    
180 !      do nm = 0,nmem-1                                                 
181       do nm = 0,nmem(nset) 
182         do 20 n=1,nx-1 
183         do 20 m=1,nq 
184         read(1,50)f1(n,m),f2(n,m),f3(n,m),f4(n,m),                      &
185      &                  f5(n,m),f7(n,m),f6(n,m),f8(n,m)                 
186 ! notation: 1=uval 2=val 3=glue 4=usea 5=chm 6=str 7=btm 8=dsea         
187    20 continue 
188 !      write(*,*)'PDF set ',nm,' first element ',f1(1,1)                
189       do 40 m=1,nq 
190       f1(nx,m)=0.d0 
191       f2(nx,m)=0.d0 
192       f3(nx,m)=0.d0 
193       f4(nx,m)=0.d0 
194       f5(nx,m)=0.d0 
195       f6(nx,m)=0.d0 
196       f7(nx,m)=0.d0 
197       f8(nx,m)=0.d0 
198    40 continue 
199                                                                         
200       do n=1,nx 
201       xxl(n)=dlog(xx(n)) 
202       enddo 
203       do m=1,nq 
204       qql(m)=dlog(qq(m)) 
205       enddo 
206                                                                         
207       call jeppe1(nm,nx,nq,xxl,qql,f1,cc1(0,1,1,1,1,nset)) 
208       call jeppe1(nm,nx,nq,xxl,qql,f2,cc2(0,1,1,1,1,nset)) 
209       call jeppe1(nm,nx,nq,xxl,qql,f3,cc3(0,1,1,1,1,nset)) 
210       call jeppe1(nm,nx,nq,xxl,qql,f4,cc4(0,1,1,1,1,nset)) 
211       call jeppe1(nm,nx,nq,xxl,qql,f6,cc6(0,1,1,1,1,nset)) 
212       call jeppe1(nm,nx,nq,xxl,qql,f8,cc8(0,1,1,1,1,nset)) 
213                                                                         
214       emc2=2.045 
215       emb2=18.5 
216                                                                         
217       do 44 m=1,nqc 
218       qqlc(m)=qql(m+nqc0) 
219       do 44 n=1,nx 
220       fc(n,m)=f5(n,m+nqc0) 
221    44 continue 
222       qqlc(1)=dlog(emc2) 
223       call jeppe1(nm,nx,nqc,xxl,qqlc,fc,ccc(0,1,1,1,1,nset)) 
224                                                                         
225       do 45 m=1,nqb 
226       qqlb(m)=qql(m+nqb0) 
227       do 45 n=1,nx 
228       fb(n,m)=f7(n,m+nqb0) 
229    45 continue 
230       qqlb(1)=dlog(emb2) 
231       call jeppe1(nm,nx,nqb,xxl,qqlb,fb,ccb(0,1,1,1,1,nset)) 
232                                                                         
233                                                                         
234       enddo 
235    50 format(8f10.5) 
236       return 
237 !                                                                       
238                                                                         
239       entry MRSTalfa(nflav,alfas,Qalfa) 
240         call getnset(iset) 
241 !        mem = member(iset)                                             
242 !        call setnmem(member(iset))                                     
243         call alphamrs(nflav,alfas,Qalfa) 
244       return 
245 !                                                                       
246       entry MRSTinit(Eorder,Q2fit) 
247       return 
248 !                                                                       
249       entry MRSTpdf(mem) 
250 !      if(mem.eq.0) mem=ndef                                            
251 !      imem = mem                                                       
252         call getnset(iset) 
253 !        member(iset)=mem                                               
254         call setnmem(iset,mem) 
255                                                                         
256       return 
257 !                                                                       
258       END