]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PDF/spdf/sfwhi4.F
dbcf134f98b88bf305d29e23a485472b00a2a1aa
[u/mrichter/AliRoot.git] / PDF / spdf / sfwhi4.F
1 *
2 * $Id$
3 *
4 * $Log$
5 * Revision 1.1.1.2  1996/10/30 08:30:45  cernlib
6 * Version 7.04
7 *
8 * Revision 1.1.1.1  1996/04/12 15:29:46  plothow
9 * Version 7.01
10 *
11 *
12 #include "pdf/pilot.h"
13 c-------------------------------------------------------
14       subroutine SFWHI4(ZX,ZQ,ZUV,ZDV,ZUB,ZDB,ZSB,ZCB,ZGL)
15 c-------------------------------------------------------
16 c     WHIT4 parton distribution in the photon
17 c
18 c     INPUT:  integer ic  : if ic=0 then qc=0
19 c                           else qc is calculated
20 c             DOUBLE PRECISION  Q2  : energy scale Q^2 (GeV^2)
21 c             DOUBLE PRECISION  x   : energy fraction
22 c
23 c     OUTPUT: DOUBLE PRECISION  qu  : up-quark dist.
24 c             DOUBLE PRECISION  qd  : down- or strange-quark dist.
25 c             DOUBLE PRECISION  qc  : charm-quark dist.
26 c             DOUBLE PRECISION  g   : gluon dist.
27 c-------------------------------------------------------
28 c     Modified by M.Tanaka on July 22, 1994.
29 c     The bug pointed out by M.Drees is fixed.
30 c-------------------------------------------------------
31 c     Modified by I.Watanabe on July 22, 1994.
32 c-------------------------------------------------------
33       implicit none
34       external WHIT4G
35 #include "pdf/expdp.inc"
36      +       ZX,ZQ,ZUV,ZDV,ZUB,ZDB,ZSB,ZCB,ZGL
37 c arg
38       integer ic
39       DOUBLE PRECISION Q2,x
40       DOUBLE PRECISION qu,qd,qc,g
41 c const
42       DOUBLE PRECISION q42it,q52it,lam42,lam52
43       DOUBLE PRECISION alinv,mc,PI
44 c local
45       DOUBLE PRECISION qv,qsea,cv,cs,dcv,dcs
46       DOUBLE PRECISION A0val,A1val,A2val,Bval,Cval,
47      $       A0sea,B0sea,BB0sea,C0sea
48       DOUBLE PRECISION A0dcv,A1dcv,A2dcv,A3dcv,Bdcv,Cdcv
49       DOUBLE PRECISION Adcs, B0dcs, B1dcs, Cdcs
50       DOUBLE PRECISION x1,x2,mc2q2
51       DOUBLE PRECISION s,s2,s3,s4,prsccf,alstpi
52       DOUBLE PRECISION WHIT4G
53 c parameters
54       parameter(lam42=0.16d0, lam52=0.091411319d0)
55       parameter(Q42IT=4.0d0, Q52IT=100.0d0)
56       parameter(alinv=137.036d0, mc=1.5d0)
57       parameter(pi=3.14159265358979323846d0)
58       common /scale/ s,s2,s3,s4,prsccf
59 c
60 c begin
61       x=ZX
62       Q2=ZQ*ZQ
63       ic=1
64 c
65       x1=1.0d0-x
66       x2=x**2
67       mc2q2=mc**2/Q2
68 c
69       if(Q2.lt.100.0d0) then
70 c  under 100 GeV^2
71 c
72 c  set  scale s
73          if(Q2.lt.4.0d0) then
74 cccc  for under 4GeV^2 prescription
75             s=  0.0d0
76             prsccf =  log(Q2/LAM42)/ log(Q42IT/LAM42)
77             alstpi = 6.0d0/25.0d0/ log(Q42IT/LAM42)
78          else
79             s=   log(  log(Q2/LAM42)/ log(Q42IT/LAM42))
80             prsccf = 1.0d0
81             alstpi = 6.0d0/25.0d0/ log(Q2/LAM42)
82          endif
83             s2=s**2
84             s3=s2*s
85             s4=s2**2
86 c
87 cccccc   WHIT4 quark (U100)
88 c
89       A0val= 2.540000d+00+s*( 2.000000d+00)+s2*( 7.180000d-01)
90       A1val= 6.230000d-02+s*(-7.010000d+00)+s2*( 1.251000d-01)
91       A2val=-1.642000d-01+s*(-4.360000d-01)+s2*( 1.048000d+01)
92      $           +s3*(-5.200000d+00)
93       Bval = 6.990000d-01+s*(-2.796000d-02)+s2*(-3.650000d-03)
94       Cval = 4.420000d-01+s*(-1.255000d+00)+s2*( 1.941000d+00)
95      $           +s3*(-9.950000d-01)
96       A0sea= 1.308000d+00+s*( 2.315000d+00)+s2*(-7.880000d+00)
97      $           +s3*( 8.260000d+00)+s4*(-3.004000d+00)
98       B0sea=-3.730000d-02+s*( 5.630000d-02)+s2*(-1.133000d+00)
99      $           +s3*( 1.185000d+00)+s4*(-4.180000d-01)
100       BB0sea=2.103000d+00+s*( 4.850000d+00)+s2*(-1.781000d+01)
101      $           +s3*( 2.062000d+01)+s4*(-7.940000d+00)
102       C0sea= 7.000000d+00+s*(-1.017000d+01)+s2*( 2.600000d+01)
103      $           +s3*(-2.960000d+01)+s4*( 1.227000d+01)
104 c
105          qv  = prsccf/alinv/x*
106      $         (A0val+A1val*x+A2val*x2) * x**Bval * x1**Cval
107          qsea= prsccf/alinv/x*
108      $         A0sea * x**(B0sea+BB0sea*x) * x1**C0sea
109 c
110          qu  = qv/3.0d0  + qsea/6.0d0
111          qu  = qu*x
112          ZUV=qu
113          ZUB=qu
114          qd  = qv/12.0d0 + qsea/6.0d0
115          qd  = qd*x
116          ZDV=qd
117          ZDB=qd
118          ZSB=qd
119 c
120          if((ic.ne.0) .and. (x*(1.0d0+4.0d0*mc2q2).lt.1.0d0)) then
121             call WHIT4Q(x,mc2q2,cv,cs)
122             qc = cv/alinv/2.0d0/PI + cs*alstpi
123             qc  = qc*x
124             ZCB=qc
125          else
126             qc = 0.0d0
127             ZCB=qc
128          endif
129 c
130          g   = WHIT4G(x,Q2)
131          g   = g*x
132          ZGL=g
133 c
134       else
135 c over 100 GeV^2
136 c
137 c set scale s
138          s=   log(  log(Q2/LAM52)/ log(Q52IT/LAM52))
139          prsccf = 1.0d0
140          alstpi = 6.0d0/23.0d0/ log(Q2/LAM52)
141             s2=s**2
142             s3=s2*s
143             s4=s2**2
144 c
145 cccccc   WHIT4 quark (O100)
146 c
147       A0val= 4.270000d+00+s*( 3.096000d+00)+s2*( 1.619000d+00)
148       A1val=-4.740000d+00+s*(-6.900000d+00)+s2*(-2.430000d+00)
149       A2val= 2.837000d+00+s*( 6.470000d+00)+s2*( 4.090000d+00)
150       Bval = 6.780000d-01+s*(-3.940000d-02)+s2*( 1.756000d-02)
151       Cval = 1.728000d-01+s*(-2.479000d-02)+s2*( 1.446000d-01)
152       A0sea= 1.188000d+00+s*(-1.396000d+00)+s2*( 8.710000d+00)
153      $           +s3*(-2.542000d+01)+s4*( 2.492000d+01)
154       B0sea=-2.448000d-01+s*(-4.190000d-01)+s2*( 1.007000d+00)
155      $           +s3*(-2.689000d+00)+s4*( 2.517000d+00)
156       BB0sea=1.942000d+00+s*(-6.040000d+00)+s2*( 5.030000d+01)
157      $           +s3*(-1.478000d+02)+s4*( 1.481000d+02)
158       C0sea= 5.420000d+00+s*( 6.110000d+00)+s2*(-5.380000d+01)
159      $           +s3*( 1.632000d+02)+s4*(-1.716000d+02)
160 c
161          qv  = 1.0d0/alinv/x*
162      $         (A0val+A1val*x+A2val*x2) * x**Bval * x1**Cval
163          qsea= 1.0d0/alinv/x*
164      $         A0sea * x**(B0sea+BB0sea*x) * x1**C0sea
165 c
166          qu  = qv/3.0d0  + qsea/6.0d0
167          qu  = qu*x
168          ZUV=qu
169          ZUB=qu
170          qd  = qv/12.0d0 + qsea/6.0d0
171          qd  = qd*x
172          ZDV=qd
173          ZDB=qd
174          ZSB=qd
175          g   = WHIT4G(x,Q2)
176          g   = g*x
177          ZGL=g
178 c
179          if((ic.ne.0) .and. (x*(1.0d0+4.0d0*mc2q2).lt.1.0d0)) then
180       A0dcv=              s*( 1.219000d-01)+s2*( 6.200000d+00)
181      $           +s3*(-2.504000d+01)+s4*( 3.098000d+01)
182       A1dcv=              s*( 1.913000d+00)+s2*(-7.690000d+01)
183      $           +s3*( 3.180000d+02)+s4*(-3.920000d+02)
184       A2dcv=              s*(-7.160000d+00)+s2*( 2.503000d+02)
185      $           +s3*(-1.062000d+03)+s4*( 1.308000d+03)
186       A3dcv=              s*( 3.190000d+00)+s2*(-2.301000d+02)
187      $           +s3*( 1.012000d+03)+s4*(-1.250000d+03)
188       Bdcv = 4.990000d-01+s*( 3.470000d+00)+s2*(-1.526000d+01)
189      $           +s3*( 1.967000d+01)
190       Cdcv = 3.290000d-01+s*( 8.240000d+00)+s2*(-3.800000d+01)
191      $           +s3*( 4.630000d+01)
192       Adcs =              s*(-2.821000d-02)+s2*(-2.649000d-04)
193      $           +s3*( 7.040000d-03)
194       B0dcs=-3.270000d-01+s*(-2.298000d-01)+s2*( 3.500000d-02)
195       B1dcs= 1.254000d+00+s*( 8.780000d-01)+s2*( 2.086000d-01)
196       Cdcs = 4.170000d+00+s*( 6.400000d-01)+s2*(-7.630000d+00)
197      $           +s3*( 7.170000d+00)
198 c
199          dcv = 1.0d0/alinv/x*
200      $         (A0dcv+x*A1dcv+x2*A2dcv+x2*x*A3dcv) * x**Bdcv * x1**Cdcv
201          dcs = 1.0d0/alinv/x*
202      $         Adcs * x**(B0dcs+B1dcs*x) * x1**Cdcs
203 c
204            call WHIT4Q(x,mc*mc/Q2,cv,cs)
205            qc = cv/alinv/2.0d0/PI + cs*alstpi + dcs + dcv
206            qc  = qc*x
207            ZCB=qc
208          else
209            qc = 0.0d0
210            ZCB=qc
211          endif
212       endif
213 c
214       return
215       end