~maddevelopers/mg5amcnlo/3.0.2-alpha0

« back to all changes in this revision

Viewing changes to HELAS/fvodmx.F

Added Template and HELAS into bzr

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
      subroutine fvodmx(fo,vc,gc,fmass,fwidth , fvo)
 
2
c
 
3
c This subroutine computes a dipole moment off-shell fermion
 
4
c wavefunction from a flowing-OUT external fermion and a vector boson.
 
5
c
 
6
c input:
 
7
c       complex fo(6)          : flow-out fermion                   <fo|
 
8
c       complex vc(6)          : input    vector                      v
 
9
c       complex gc(2,2)        : coupling constants                  gvf
 
10
c                              : first index is L,R as normal
 
11
c                              : second index is EDM,-MDM
 
12
c       real    fmass          : mass  of output fermion f'
 
13
c       real    fwidth         : width of output fermion f'
 
14
c
 
15
c output:
 
16
c       complex fvo(6)         : off-shell fermion             <fo,v,f'|
 
17
c
 
18
      implicit none
 
19
      double complex fo(6), vc(6), fvo(6), sl1, sl2, sr1, sr2, d
 
20
      double complex gc(2,2), gL, gR
 
21
      double precision  pf(0:3), fmass, fwidth, pf2
 
22
 
 
23
      double complex kvc21, kvc31, kvc41, kvc32, kvc42, kvc43
 
24
      double precision k(1:4)
 
25
      double precision rZero, rOne
 
26
      parameter( rZero = 0.0d0, rOne = 1.0d0 )
 
27
      double complex cImag, cZero
 
28
      parameter( cImag = ( 0.0d0, 1.0d0 ), cZero = ( 0.0d0, 0.0d0 ) )
 
29
 
 
30
#ifdef HELAS_CHECK
 
31
      integer stdo
 
32
      parameter( stdo = 6 )
 
33
#endif
 
34
c
 
35
#ifdef HELAS_CHECK
 
36
      if ( abs(fo(1))+abs(fo(2))+abs(fo(3))+abs(fo(4)).eq.rZero ) then
 
37
         write(stdo,*) ' helas-warn  : fo in fvodmx is zero spinor'
 
38
      endif
 
39
      if ( abs(fo(5))+abs(fo(6)).eq.rZero ) then
 
40
         write(stdo,*)
 
41
     &        ' helas-error : fo in fvodmx has zero momentum'
 
42
      endif
 
43
      if ( abs(vc(1))+abs(vc(2))+abs(vc(3))+abs(vc(4)).eq.rZero ) then
 
44
         write(stdo,*) ' helas-warn  : vc in fvodmx is zero vector'
 
45
      endif
 
46
      if ( abs(vc(5))+abs(vc(6)).eq.rZero ) then
 
47
         write(stdo,*)
 
48
     &        ' helas-error : vc in fvodmx has zero momentum'
 
49
      endif
 
50
      if ( gc(1,1).eq.cZero .and. gc(2,1).eq.cZero .and.
 
51
     &     gc(1,2).eq.cZero .and. gc(2,2).eq.cZero      ) then
 
52
         write(stdo,*)
 
53
     &        ' helas-error : gc in fvodmx is zero coupling'
 
54
      endif
 
55
      if ( fmass.lt.rZero ) then
 
56
         write(stdo,*) ' helas-error : fmass in fvodmx is negative'
 
57
         write(stdo,*) '             : fmass = ',fmass
 
58
      endif
 
59
      if ( fwidth.lt.rZero ) then
 
60
         write(stdo,*) ' helas-error : fwidth in fvodmx is negative'
 
61
         write(stdo,*) '             : fwidth = ',fwidth
 
62
      endif
 
63
#endif
 
64
 
 
65
      gL = -gc(1,1) + cImag*gc(1,2)
 
66
      gR =  gc(2,1) + cImag*gc(2,2)
 
67
 
 
68
c k in vertex formula defined as (pi - po)
 
69
      k(1) = dble( vc(5))
 
70
      k(2) = dble( vc(6))
 
71
      k(3) = dimag(vc(6))
 
72
      k(4) = dimag(vc(5))
 
73
 
 
74
      kvc21 = (k(2)*vc(1) - k(1)*vc(2))*cImag
 
75
      kvc31 =  k(3)*vc(1) - k(1)*vc(3)
 
76
      kvc41 = (k(4)*vc(1) - k(1)*vc(4))*cImag
 
77
      kvc32 =  k(3)*vc(2) - k(2)*vc(3)
 
78
      kvc42 = (k(4)*vc(2) - k(2)*vc(4))*cImag
 
79
      kvc43 =  k(4)*vc(3) - k(3)*vc(4)
 
80
 
 
81
      fvo(5) = fo(5) + vc(5)
 
82
      fvo(6) = fo(6) + vc(6)
 
83
 
 
84
      pf(0) = dble( fvo(5))
 
85
      pf(1) = dble( fvo(6))
 
86
      pf(2) = dimag(fvo(6))
 
87
      pf(3) = dimag(fvo(5))
 
88
      pf2 = pf(0)**2-(pf(1)**2+pf(2)**2+pf(3)**2)
 
89
 
 
90
#ifdef HELAS_CHECK
 
91
      if ( abs(fvo(5))+abs(fvo(6)).eq.rZero ) then
 
92
         write(stdo,*)
 
93
     &        ' helas-error : fvo in fvodmx has zero momentum'
 
94
      endif
 
95
      if ( fwidth.eq.rZero .and. pf2.eq.fmass**2 ) then
 
96
         write(stdo,*)
 
97
     &         ' helas-error : fvo in fvodmx is on fmass pole'
 
98
         write(stdo,*)
 
99
     &         '             : p     = ',pf(0),pf(1),pf(2),pf(3)
 
100
         write(stdo,*)
 
101
     &         '             : abs(p)= ',sqrt(abs(pf2))
 
102
         fvo(1) = cZero
 
103
         fvo(2) = cZero
 
104
         fvo(3) = cZero
 
105
         fvo(4) = cZero
 
106
         return
 
107
      endif
 
108
#endif
 
109
 
 
110
      d = -rOne/dcmplx( pf2-fmass**2, fmass*fwidth )
 
111
 
 
112
      sl1 = gL*(  fo(1)*(kvc41 + kvc32)
 
113
     &          - fo(2)*(kvc42 - kvc21 - kvc43 + kvc31) )
 
114
      sl2 = gL*(  fo(1)*(kvc42 + kvc21 + kvc43 + kvc31)
 
115
     &          - fo(2)*(kvc41 + kvc32)                 )
 
116
 
 
117
      if ( gc(2,1).ne.cZero .or.
 
118
     &     gc(2,2).ne.cZero      ) then
 
119
         sr1 = gR*(- fo(3)*(kvc41 - kvc32)
 
120
     &             - fo(4)*(kvc42 + kvc21 - kvc43 - kvc31) )
 
121
         sr2 = gR*(  fo(3)*(kvc42 - kvc21 + kvc43 - kvc31)
 
122
     &             + fo(4)*(kvc41 - kvc32)                 )
 
123
 
 
124
         fvo(1) = (  (pf(0)+pf(3))*sr1       + fvo(6)*sr2
 
125
     &             + fmass*sl1                            )*d
 
126
         fvo(2) = ( dconjg(fvo(6))*sr1 +(pf(0)-pf(3))*sr2
 
127
     &             + fmass*sl2                            )*d
 
128
         fvo(3) = (  (pf(0)-pf(3))*sl1       - fvo(6)*sl2
 
129
     &             + fmass*sr1                            )*d
 
130
         fvo(4) = (-dconjg(fvo(6))*sl1 +(pf(0)+pf(3))*sl2
 
131
     &             + fmass*sr2                            )*d
 
132
 
 
133
      else
 
134
         fvo(1) = fmass*sl1*d
 
135
         fvo(2) = fmass*sl2*d
 
136
         fvo(3) = (  (pf(0)-pf(3))*sl1        - fvo(6)*sl2)*d
 
137
         fvo(4) = (-dconjg(fvo(6))*sl1 + (pf(0)+pf(3))*sl2)*d
 
138
      end if
 
139
c
 
140
      return
 
141
      end