1
subroutine fvixxx(fi,vc,gc,fmass,fwidth , fvi)
3
c This subroutine computes an off-shell fermion wavefunction from a
4
c flowing-IN external fermion and a vector boson.
7
c complex fi(6) : flow-in fermion |fi>
8
c complex vc(6) : input vector v
9
c complex gc(2) : coupling constants gvf
10
c real fmass : mass of output fermion f'
11
c real fwidth : width of output fermion f'
14
c complex fvi(6) : off-shell fermion |f',v,fi>
17
double complex fi(6),vc(6),gc(2),fvi(6),sl1,sl2,sr1,sr2,d
18
double precision pf(0:3),fmass,fwidth,pf2
20
double precision rZero, rOne
21
parameter( rZero = 0.0d0, rOne = 1.0d0 )
22
double complex cImag, cZero
23
parameter( cImag = ( 0.0d0, 1.0d0 ), cZero = ( 0.0d0, 0.0d0 ) )
31
if ( abs(fi(1))+abs(fi(2))+abs(fi(3))+abs(fi(4)).eq.rZero ) then
32
write(stdo,*) ' helas-warn : fi in fvixxx is zero spinor'
34
if ( abs(fi(5))+abs(fi(6)).eq.rZero ) then
36
& ' helas-error : fi in fvixxx has zero momentum'
38
if ( abs(vc(1))+abs(vc(2))+abs(vc(3))+abs(vc(4)).eq.rZero ) then
39
write(stdo,*) ' helas-warn : vc in fvixxx is zero vector'
41
if ( abs(vc(5))+abs(vc(6)).eq.rZero ) then
43
& ' helas-error : vc in fvixxx has zero momentum'
45
if ( gc(1).eq.cZero .and. gc(2).eq.cZero ) then
47
& ' helas-error : gc in fvixxx is zero coupling'
49
if ( fmass.lt.rZero ) then
50
write(stdo,*) ' helas-error : fmass in fvixxx is negative'
51
write(stdo,*) ' : fmass = ',fmass
53
if ( fwidth.lt.rZero ) then
54
write(stdo,*) ' helas-error : fwidth in fvixxx is negative'
55
write(stdo,*) ' : fwidth = ',fwidth
66
pf2 = pf(0)**2-(pf(1)**2+pf(2)**2+pf(3)**2)
69
if ( abs(fvi(5))+abs(fvi(6)).eq.rZero ) then
71
& ' helas-error : fvi in fvixxx has zero momentum'
73
if ( fwidth.eq.rZero .and. pf2.eq.fmass**2 ) then
75
& ' helas-error : fvi in fvixxx is on fmass pole'
77
& ' : p = ',pf(0),pf(1),pf(2),pf(3)
79
& ' : abs(p)= ',dsqrt(dabs(pf2))
88
d = -rOne/dcmplx( pf2-fmass**2, fmass*fwidth )
89
sl1 = (vc(1)+ vc(4))*fi(1)
90
& + (vc(2)-cImag*vc(3))*fi(2)
91
sl2 = (vc(2)+cImag*vc(3))*fi(1)
92
& + (vc(1)- vc(4))*fi(2)
94
if ( gc(2).ne.cZero ) then
95
sr1 = (vc(1)- vc(4))*fi(3)
96
& - (vc(2)-cImag*vc(3))*fi(4)
97
sr2 = - (vc(2)+cImag*vc(3))*fi(3)
98
& + (vc(1)+ vc(4))*fi(4)
100
fvi(1) = ( gc(1)*((pf(0)-pf(3))*sl1 - dconjg(fvi(6))*sl2)
101
& +gc(2)*fmass*sr1 )*d
102
fvi(2) = ( gc(1)*( -fvi(6)*sl1 + (pf(0)+pf(3))*sl2)
103
& +gc(2)*fmass*sr2 )*d
104
fvi(3) = ( gc(2)*((pf(0)+pf(3))*sr1 + dconjg(fvi(6))*sr2)
105
& +gc(1)*fmass*sl1 )*d
106
fvi(4) = ( gc(2)*( fvi(6)*sr1 + (pf(0)-pf(3))*sr2)
107
& +gc(1)*fmass*sl2 )*d
111
fvi(1) = ((pf(0)-pf(3))*sl1 - dconjg(fvi(6))*sl2)*d
112
fvi(2) = ( -fvi(6)*sl1 + (pf(0)+pf(3))*sl2)*d