Package aloha :: Package template_files :: Module wavefunctions
[hide private]
[frames] | no frames]

Source Code for Module aloha.template_files.wavefunctions

  1  from __future__ import division  
  2  import math 
  3  from math import sqrt, pow 
  4  from itertools import product 
  5   
6 -class WaveFunction(list):
7 """a objet for a WaveFunction""" 8 9 spin_to_size={0:1, 10 1:3, 11 2:6, 12 3:6, 13 4:18, 14 5:18} 15
16 - def __init__(self, spin= None, size=None):
17 """Init the list with zero value""" 18 19 if spin: 20 size = self.spin_to_size[spin] 21 list.__init__(self, [0]*size)
22 23
24 -def ixxxxx(p,fmass,nhel,nsf):
25 """Defines an inflow fermion.""" 26 27 fi = WaveFunction(2) 28 29 fi[0] = complex(-p[0]*nsf,-p[3]*nsf) 30 fi[1] = complex(-p[1]*nsf,-p[2]*nsf) 31 nh = nhel*nsf 32 if (fmass != 0.): 33 pp = min(p[0],sqrt(p[1]**2 + p[2]**2 + p[3]**2 )) 34 if (pp == 0.): 35 sqm = [sqrt(abs(fmass))] 36 sqm.append(sign(sqm[0],fmass)) 37 ip = (1+nh)//2 38 im = (1-nh)//2 39 40 fi[2] = ip*sqm[ip] 41 fi[3] = im*nsf*sqm[ip] 42 fi[4] = ip*nsf*sqm[im] 43 fi[5] = im*sqm[im] 44 45 else: 46 sf = [(1+nsf+(1-nsf)*nh)*0.5,(1+nsf-(1-nsf)*nh)*0.5] 47 omega = [sqrt(p[0]+pp),fmass/(sqrt(p[0]+pp))] 48 ip = (1+nh)//2 49 im = (1-nh)//2 50 sfomeg = [sf[0]*omega[ip],sf[1]*omega[im]] 51 pp3 = max(pp+p[3],0.) 52 if (pp3 == 0.): 53 chi1 = complex(-nh,0.) 54 else: 55 chi1 = complex(nh*p[1]/sqrt(2.*pp*pp3),\ 56 p[2]/sqrt(2.*pp*pp3)) 57 chi = [complex(sqrt(pp3*0.5/pp)),chi1] 58 59 fi[2] = sfomeg[0]*chi[im] 60 fi[3] = sfomeg[0]*chi[ip] 61 fi[4] = sfomeg[1]*chi[im] 62 fi[5] = sfomeg[1]*chi[ip] 63 64 else: 65 sqp0p3 = sqrt(max(p[0]+p[3],0.))*nsf 66 if (sqp0p3 == 0.): 67 chi1 = complex(-nhel*sqrt(2.*p[0]),0.) 68 else: 69 chi1 = complex(nh*p[1]/sqp0p3,p[2]/sqp0p3) 70 chi = [complex(sqp0p3,0.),chi1] 71 if (nh == 1): 72 fi[2] = complex(0.,0.) 73 fi[3] = complex(0.,0.) 74 fi[4] = chi[0] 75 fi[5] = chi[1] 76 else: 77 fi[2] = chi[1] 78 fi[3] = chi[0] 79 fi[4] = complex(0.,0.) 80 fi[5] = complex(0.,0.) 81 82 return fi
83
84 -def oxxxxx(p,fmass,nhel,nsf):
85 """ initialize an outgoing fermion""" 86 87 fo = WaveFunction(2) 88 fo[0] = complex(p[0]*nsf,p[3]*nsf) 89 fo[1] = complex(p[1]*nsf,p[2]*nsf) 90 nh = nhel*nsf 91 if (fmass != 0.): 92 pp = min(p[0],sqrt(p[1]**2 + p[2]**2 + p[3]**2 )) 93 if (pp == 0.): 94 sqm = [sqrt(abs(fmass))] 95 sqm.append( sign(sqm[0], fmass)) 96 ip = int(-((1-nh)//2) * nhel) 97 im = int((1+nh)//2 * nhel) 98 99 fo[2] = im*sqm[abs(im)] 100 fo[3] = ip*nsf*sqm[abs(im)] 101 fo[4] = im*nsf*sqm[abs(ip)] 102 fo[5] = ip*sqm[abs(ip)] 103 104 else: 105 sf = [(1+nsf+(1-nsf)*nh)*0.5,(1+nsf-(1-nsf)*nh)*0.5] 106 omega = [sqrt(p[0]+pp),fmass/(sqrt(p[0]+pp))] 107 ip = (1+nh)//2 108 im = (1-nh)//2 109 sfomeg = [sf[0]*omega[ip],sf[1]*omega[im]] 110 pp3 = max(pp+p[3],0.) 111 if (pp3 == 0.): 112 chi1 = complex(-nh,0.) 113 else: 114 chi1 = complex(nh*p[1]/sqrt(2.*pp*pp3),\ 115 -p[2]/sqrt(2.*pp*pp3)) 116 chi = [complex(sqrt(pp3*0.5/pp)),chi1] 117 118 fo[2] = sfomeg[1]*chi[im] 119 fo[3] = sfomeg[1]*chi[ip] 120 fo[4] = sfomeg[0]*chi[im] 121 fo[5] = sfomeg[0]*chi[ip] 122 123 else: 124 sqp0p3 = sqrt(max(p[0]+p[3],0.))*nsf 125 if (sqp0p3 == 0.): 126 chi1 = complex(-nhel*sqrt(2.*p[0]),0.) 127 else: 128 chi1 = complex(nh*p[1]/sqp0p3,-p[2]/sqp0p3) 129 chi = [complex(sqp0p3,0.),chi1] 130 if (nh == 1): 131 fo[2] = chi[0] 132 fo[3] = chi[1] 133 fo[4] = complex(0.,0.) 134 fo[5] = complex(0.,0.) 135 else: 136 fo[2] = complex(0.,0.) 137 fo[3] = complex(0.,0.) 138 fo[4] = chi[1] 139 fo[5] = chi[0] 140 141 return fo
142
143 -def vxxxxx(p,vmass,nhel,nsv):
144 """ initialize a vector wavefunction. nhel=4 is for checking BRST""" 145 146 vc = WaveFunction(3) 147 148 sqh = sqrt(0.5) 149 nsvahl = nsv*abs(nhel) 150 pt2 = p[1]**2 + p[2]**2 151 pp = min(p[0],sqrt(pt2 + p[3]**2)) 152 pt = min(pp,sqrt(pt2)) 153 154 vc[0] = complex(p[0]*nsv,p[3]*nsv) 155 vc[1] = complex(p[1]*nsv,p[2]*nsv) 156 157 if (nhel == 4): 158 if (vmass == 0.): 159 vc[2] = 1. 160 vc[3]=p[1]/p[0] 161 vc[4]=p[2]/p[0] 162 vc[5]=p[3]/p[0] 163 else: 164 vc[2] = p[0]/vmass 165 vc[3] = p[1]/vmass 166 vc[4] = p[2]/vmass 167 vc[5] = p[3]/vmass 168 169 return vc 170 171 if (vmass != 0.): 172 hel0 = 1.-abs(nhel) 173 174 if (pp == 0.): 175 vc[2] = complex(0.,0.) 176 vc[3] = complex(-nhel*sqh,0.) 177 vc[4] = complex(0.,nsvahl*sqh) 178 vc[5] = complex(hel0,0.) 179 180 else: 181 emp = p[0]/(vmass*pp) 182 vc[2] = complex(hel0*pp/vmass,0.) 183 vc[5] = complex(hel0*p[3]*emp+nhel*pt/pp*sqh) 184 if (pt != 0.): 185 pzpt = p[3]/(pp*pt)*sqh*nhel 186 vc[3] = complex(hel0*p[1]*emp-p[1]*pzpt, \ 187 -nsvahl*p[2]/pt*sqh) 188 vc[4] = complex(hel0*p[2]*emp-p[2]*pzpt, \ 189 nsvahl*p[1]/pt*sqh) 190 else: 191 vc[3] = complex(-nhel*sqh,0.) 192 vc[4] = complex(0.,nsvahl*sign(sqh,p[3])) 193 else: 194 pp = p[0] 195 pt = sqrt(p[1]**2 + p[2]**2) 196 vc[2] = complex(0.,0.) 197 vc[5] = complex(nhel*pt/pp*sqh) 198 if (pt != 0.): 199 pzpt = p[3]/(pp*pt)*sqh*nhel 200 vc[3] = complex(-p[1]*pzpt,-nsv*p[2]/pt*sqh) 201 vc[4] = complex(-p[2]*pzpt,nsv*p[1]/pt*sqh) 202 else: 203 vc[3] = complex(-nhel*sqh,0.) 204 vc[4] = complex(0.,nsv*sign(sqh,p[3])) 205 206 return vc
207
208 -def sign(x,y):
209 """Fortran's sign transfer function""" 210 try: 211 cmp = (y < 0.) 212 except TypeError: 213 # y should be complex 214 if abs(y.imag) < 1e-6 * abs(y.real): 215 y = y.real 216 else: 217 raise 218 finally: 219 if (y < 0.): 220 return -abs(x) 221 else: 222 return abs(x)
223
224 -def sxxxxx(p,nss):
225 """initialize a scalar wavefunction""" 226 227 sc = WaveFunction(1) 228 229 sc[2] = complex(1.,0.) 230 sc[0] = complex(p[0]*nss,p[3]*nss) 231 sc[1] = complex(p[1]*nss,p[2]*nss) 232 return sc
233 234
235 -def txxxxx(p, tmass, nhel, nst):
236 """ initialize a tensor wavefunction""" 237 238 tc = WaveFunction(5) 239 240 sqh = sqrt(0.5) 241 sqs = sqrt(1/6) 242 243 pt2 = p[1]**2 + p[2]**2 244 pp = min(p[0],sqrt(pt2+p[3]**2)) 245 pt = min(pp,sqrt(pt2)) 246 247 ft = {} 248 ft[(4,0)] = complex(p[0], p[3]) * nst 249 ft[(5,0)] = complex(p[1], p[2]) * nst 250 251 if ( nhel >= 0 ): 252 #construct eps+ 253 ep = [0] * 4 254 255 if ( pp == 0 ): 256 #ep[0] = 0 257 ep[1] = -sqh 258 ep[2] = complex(0, nst*sqh) 259 #ep[3] = 0 260 else: 261 #ep[0] = 0 262 ep[3] = pt/pp*sqh 263 if (pt != 0): 264 pzpt = p[3]/(pp*pt)*sqh 265 ep[1] = complex( -p[1]*pzpt , -nst*p[2]/pt*sqh ) 266 ep[2] = complex( -p[2]*pzpt , nst*p[1]/pt*sqh ) 267 else: 268 ep[1] = -sqh 269 ep[2] = complex( 0 , nst*sign(sqh,p[3]) ) 270 271 272 273 if ( nhel <= 0 ): 274 #construct eps- 275 em = [0] * 4 276 if ( pp == 0 ): 277 #em[0] = 0 278 em[1] = sqh 279 em[2] = complex( 0 , nst*sqh ) 280 #em[3] = 0 281 else: 282 #em[0] = 0 283 em[3] = -pt/pp*sqh 284 if pt: 285 pzpt = -p[3]/(pp*pt)*sqh 286 em[1] = complex( -p[1]*pzpt , -nst*p[2]/pt*sqh ) 287 em[2] = complex( -p[2]*pzpt , nst*p[1]/pt*sqh ) 288 else: 289 em[1] = sqh 290 em[2] = complex( 0 , nst*sign(sqh,p[3]) ) 291 292 293 if ( abs(nhel) <= 1 ): 294 #construct eps0 295 e0 = [0] * 4 296 if ( pp == 0 ): 297 #e0[0] = complex( rZero ) 298 #e0[1] = complex( rZero ) 299 #e0[2] = complex( rZero ) 300 e0[3] = 1 301 else: 302 emp = p[0]/(tmass*pp) 303 e0[0] = pp/tmass 304 e0[3] = p[3]*emp 305 if pt: 306 e0[1] = p[1]*emp 307 e0[2] = p[2]*emp 308 #else: 309 # e0[1] = complex( rZero ) 310 # e0[2] = complex( rZero ) 311 312 if nhel == 2: 313 for j in range(4): 314 for i in range(4): 315 ft[(i,j)] = ep[i]*ep[j] 316 elif nhel == -2: 317 for j in range(4): 318 for i in range(4): 319 ft[(i,j)] = em[i]*em[j] 320 elif tmass == 0: 321 for j in range(4): 322 for i in range(4): 323 ft[(i,j)] = 0 324 elif nhel == 1: 325 for j in range(4): 326 for i in range(4): 327 ft[(i,j)] = sqh*( ep[i]*e0[j] + e0[i]*ep[j] ) 328 elif nhel == 0: 329 for j in range(4): 330 for i in range(4): 331 ft[(i,j)] = sqs*( ep[i]*em[j] + em[i]*ep[j] + 2 *e0[i]*e0[j] ) 332 elif nhel == -1: 333 for j in range(4): 334 for i in range(4): 335 ft[(i,j)] = sqh*( em[i]*e0[j] + e0[i]*em[j] ) 336 337 else: 338 raise Exception, 'invalid helicity TXXXXXX' 339 340 341 342 tc[2] = ft[(0,0)] 343 tc[3] = ft[(0,1)] 344 tc[4] = ft[(0,2)] 345 tc[5] = ft[(0,3)] 346 tc[6] = ft[(1,0)] 347 tc[7] = ft[(1,1)] 348 tc[8] = ft[(1,2)] 349 tc[9] = ft[(1,3)] 350 tc[10] = ft[(2,0)] 351 tc[11] = ft[(2,1)] 352 tc[12] = ft[(2,2)] 353 tc[13] = ft[(2,3)] 354 tc[14] = ft[(3,0)] 355 tc[15] = ft[(3,1)] 356 tc[16] = ft[(3,2)] 357 tc[17] = ft[(3,3)] 358 tc[0] = ft[(4,0)] 359 tc[1] = ft[(5,0)] 360 361 return tc
362
363 -def irxxxx(p, mass, nhel, nsr):
364 """ initialize a incoming spin 3/2 wavefunction.""" 365 366 # This routines is a python translation of a routine written by 367 # K.Mawatari in fortran dated from the 2008/02/26 368 369 ri = WaveFunction(4) 370 371 sqh = sqrt(0.5) 372 sq2 = sqrt(2) 373 sq3 = sqrt(3) 374 375 pt2 = p[1]**2 + p[2]**2 376 pp = min([p[0], sqrt(pt2+p[3]**2)]) 377 pt = min([pp,sqrt(pt2)]) 378 379 rc = {} 380 rc[(4,0)] = -1*complex(p[0],p[3])*nsr 381 rc[(5,0)] = -1*complex(p[1],p[2])*nsr 382 383 384 nsv = -nsr # nsv=+1 for final, -1 for initial 385 386 # Construct eps+ 387 if nhel > 0: 388 ep = [0] * 4 389 if pp: 390 #ep[0] = 0 391 ep[3] = pt/pp*sqh 392 if pt: 393 pzpt = p[3]/(pp*pt)*sqh 394 ep[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 395 ep[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 396 else: 397 ep[1] = -sqh 398 ep[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 399 else: 400 #ep[0] = 0d0 401 ep[1] = -sqh 402 ep[2] = complex(0, nsv * sqh) 403 #ep[3] = 0d0 404 405 406 if ( nhel < 0 ): 407 #construct eps- 408 em = [0] * 4 409 if ( pp == 0 ): 410 #em[0] = 0 411 em[1] = sqh 412 em[2] = complex( 0 , nsv*sqh ) 413 #em[3] = 0 414 else: 415 #em[0] = 0 416 em[3] = -pt/pp*sqh 417 if pt: 418 pzpt = -p[3]/(pp*pt)*sqh 419 em[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 420 em[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 421 else: 422 em[1] = sqh 423 em[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 424 425 if ( abs(nhel) <= 1 ): 426 #construct eps0 427 e0 = [0] * 4 428 if ( pp == 0 ): 429 #e0[0] = complex( rZero ) 430 #e0[1] = complex( rZero ) 431 #e0[2] = complex( rZero ) 432 e0[3] = 1 433 else: 434 emp = p[0]/(mass*pp) 435 e0[0] = pp/mass 436 e0[3] = p[3]*emp 437 if pt: 438 e0[1] = p[1]*emp 439 e0[2] = p[2]*emp 440 #else: 441 # e0[1] = complex( rZero ) 442 # e0[2] = complex( rZero ) 443 444 445 446 if ( nhel >= -1 ): 447 # constract spinor+ 448 fip = [0] * 4 449 sf, omega, sfomeg, chi = [0, 0], [0,0], [0,0], [0,0] 450 nh = nsr 451 if mass: 452 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 453 if pp == 0: 454 sqm = sqrt(mass) 455 ip = (1+nh)//2 456 im = (1-nh)//2 457 fip[0] = ip * sqm 458 fip[1] = im * nsr * sqm 459 fip[2] = ip * nsr * sqm 460 fip[3] = im * sqm 461 else: 462 sf[0] = float(1+nsr+(1-nsr)*nh)*0.5 463 sf[1] = float(1+nsr-(1-nsr)*nh)*0.5 464 omega[0] = sqrt(p[0]+pp) 465 omega[1] = mass/omega[0] 466 ip = ((3+nh)//2) -1 # -1 since they are index 467 im = ((3-nh)//2) -1 # -1 since they are index 468 sfomeg[0] = sf[0]*omega[ip] 469 sfomeg[1] = sf[1]*omega[im] 470 pp3 = max([pp+p[3],0]) 471 chi[0] = sqrt(pp3*0.5/pp) 472 if pp3 ==0: 473 chi[1] = -nh 474 else: 475 chi[1] = complex( nh*p[1] , p[2] )/sqrt(2*pp*pp3) 476 477 fip[0] = sfomeg[0]*chi[im] 478 fip[1] = sfomeg[0]*chi[ip] 479 fip[2] = sfomeg[1]*chi[im] 480 fip[3] = sfomeg[1]*chi[ip] 481 else: 482 sqp0p3 = sqrt(max([p[0]+p[3],0])) * nsr 483 chi[0] = sqp0p3 484 if sqp0p3 == 0: 485 chi[1] = -nhel * sqrt(2*p[0]) 486 else: 487 chi[1] = complex( nh*p[1], p[2] )/sqp0p3 488 if nh == 1: 489 #fip[0] = complex( rZero ) 490 #fip[1] = complex( rZero ) 491 fip[2] = chi[0] 492 fip[3] = chi[1] 493 else: 494 fip[0] = chi[1] 495 fip[1] = chi[0] 496 #fip(3) = complex( rZero ) 497 #fip(4) = complex( rZero ) 498 499 if ( nhel <= 1 ): 500 # constract spinor- 501 fim = [0] * 4 502 sf, omega, sfomeg, chi = [0, 0], [0,0], [0,0], [0,0] 503 nh = -nsr 504 if mass: 505 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 506 if pp == 0: 507 sqm = sqrt(mass) 508 ip = (1+nh)/2 509 im = (1-nh)/2 510 fim[0] = ip * sqm 511 fim[1] = im * nsr * sqm 512 fim[2] = ip * nsr * sqm 513 fim[3] = im * sqm 514 else: 515 sf[0] = float(1+nsr+(1-nsr)*nh)*0.5 516 sf[1] = float(1+nsr-(1-nsr)*nh)*0.5 517 omega[0] = sqrt(p[0]+pp) 518 omega[1] = mass/omega[0] 519 ip = (3+nh)//2 -1 520 im = (3-nh)//2 -1 521 sfomeg[0] = sf[0]*omega[ip] 522 sfomeg[1] = sf[1]*omega[im] 523 pp3 = max([pp+p[3],0]) 524 chi[0] = sqrt(pp3*0.5/pp) 525 if pp3 ==0: 526 chi[1] = -nh 527 else: 528 chi[1] = complex( nh*p[1] , p[2] )/sqrt(2*pp*pp3) 529 530 fim[0] = sfomeg[0]*chi[im] 531 fim[1] = sfomeg[0]*chi[ip] 532 fim[2] = sfomeg[1]*chi[im] 533 fim[3] = sfomeg[1]*chi[ip] 534 else: 535 sqp0p3 = sqrt(max([p[0]+p[3],0])) * nsr 536 chi[0] = sqp0p3 537 if sqp0p3 == 0: 538 chi[1] = -nhel * sqrt(2*p[0]) 539 else: 540 chi[1] = complex( nh*p[1], p[2] )/sqp0p3 541 if nh == 1: 542 #fip[0] = complex( rZero ) 543 #fip[1] = complex( rZero ) 544 fim[2] = chi[0] 545 fim[3] = chi[1] 546 else: 547 fim[0] = chi[1] 548 fim[1] = chi[0] 549 #fip(3) = complex( rZero ) 550 #fip(4) = complex( rZero ) 551 552 553 554 # recurent relation put her for optimization 555 cond1 = (pt == 0 and p[3] >= 0) 556 cond2 = (pt == 0 and p[3] < 0) 557 558 # spin-3/2 fermion wavefunction 559 if nhel == 3: 560 for i,j in product(range(4), range(4)): 561 rc[(i, j)] = ep[i] *fip[j] 562 563 elif nhel == 1: 564 for i,j in product(range(4), range(4)): 565 if cond1: 566 rc[(i,j)] = sq2/sq3*e0[i]*fip[j] +1/sq3*ep[i]*fim[j] 567 elif cond2: 568 rc[(i,j)] = sq2/sq3*e0[i]*fip[j] -1/sq3*ep[i]*fim[j] 569 else: 570 rc[(i,j)] = sq2/sq3*e0[i]*fip[j] + \ 571 1/sq3*ep[i]*fim[j] *complex(p[1],nsr*p[2])/pt 572 elif nhel == -1: 573 for i,j in product(range(4), range(4)): 574 if cond1: 575 rc[(i,j)] = 1/sq3*em[i]*fip[j] +sq2/sq3*e0[i]*fim[j] 576 elif cond2: 577 rc[(i,j)] = 1/sq3*em[i]*fip[j] -sq2/sq3*e0[i]*fim[j] 578 else: 579 rc[(i,j)] = 1/sq3*em[i]*fip[j] + \ 580 sq2/sq3*e0[i]*fim[j] *complex(p[1],nsr*p[2])/pt 581 else: 582 for i,j in product(range(4), range(4)): 583 if cond1: 584 rc[(i, j)] = em[i] *fim[j] 585 elif cond2: 586 rc[(i, j)] = -em[i] *fim[j] 587 else: 588 rc[(i, j)] = em[i]*fim[j] *complex(p[1],nsr*p[2])/pt 589 590 ri[2] = rc[(0,0)] 591 ri[3] = rc[(0,1)] 592 ri[4] = rc[(0,2)] 593 ri[5] = rc[(0,3)] 594 ri[6] = rc[(1,0)] 595 ri[7] = rc[(1,1)] 596 ri[8] = rc[(1,2)] 597 ri[9] = rc[(1,3)] 598 ri[10] = rc[(2,0)] 599 ri[11] = rc[(2,1)] 600 ri[12] = rc[(2,2)] 601 ri[13] = rc[(2,3)] 602 ri[14] = rc[(3,0)] 603 ri[15] = rc[(3,1)] 604 ri[16] = rc[(3,2)] 605 ri[17] = rc[(3,3)] 606 ri[0] = rc[(4,0)] 607 ri[1] = rc[(5,0)] 608 609 return ri
610
611 -def orxxxx(p, mass, nhel, nsr):
612 """ initialize a incoming spin 3/2 wavefunction.""" 613 614 # This routines is a python translation of a routine written by 615 # K.Mawatari in fortran dated from the 2008/02/26 616 617 618 ro = WaveFunction(spin=4) 619 620 sqh = sqrt(0.5) 621 sq2 = sqrt(2) 622 sq3 = sqrt(3) 623 624 pt2 = p[1]**2 + p[2]**2 625 pp = min([p[0], sqrt(pt2+p[3]**2)]) 626 pt = min([pp,sqrt(pt2)]) 627 rc = {} 628 rc[(4,0)] = complex(p[0],p[3])*nsr 629 rc[(5,0)] = complex(p[1],p[2])*nsr 630 631 632 nsv = nsr # nsv=+1 for final, -1 for initial 633 634 # Construct eps+ 635 if nhel > 0: 636 ep = [0] * 4 637 if pp: 638 #ep[0] = 0 639 ep[3] = pt/pp*sqh 640 if pt: 641 pzpt = p[3]/(pp*pt)*sqh 642 ep[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 643 ep[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 644 else: 645 ep[1] = -sqh 646 ep[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 647 else: 648 #ep[0] = 0d0 649 ep[1] = -sqh 650 ep[2] = complex(0, nsv * sqh) 651 #ep[3] = 0d0 652 653 if ( nhel < 0 ): 654 #construct eps- 655 em = [0] * 4 656 if ( pp == 0 ): 657 #em[0] = 0 658 em[1] = sqh 659 em[2] = complex( 0 , nsv*sqh ) 660 #em[3] = 0 661 else: 662 #em[0] = 0 663 em[3] = -pt/pp*sqh 664 if pt: 665 pzpt = -p[3]/(pp*pt)*sqh 666 em[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 667 em[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 668 else: 669 em[1] = sqh 670 em[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 671 672 if ( abs(nhel) <= 1 ): 673 #construct eps0 674 e0 = [0] * 4 675 if ( pp == 0 ): 676 #e0[0] = complex( rZero ) 677 #e0[1] = complex( rZero ) 678 #e0[2] = complex( rZero ) 679 e0[3] = 1 680 else: 681 emp = p[0]/(mass*pp) 682 e0[0] = pp/mass 683 e0[3] = p[3]*emp 684 if pt: 685 e0[1] = p[1]*emp 686 e0[2] = p[2]*emp 687 #else: 688 # e0[1] = complex( rZero ) 689 # e0[2] = complex( rZero ) 690 691 if nhel >= -1: 692 #constract spinor+ 693 nh = nsr 694 sqm, fop, omega, sf, sfomeg = [0]*2,[0]*4,[0]*2,[0]*2,[0]*2 695 chi = [0]*2 696 if mass: 697 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 698 if ( pp == 0): 699 sqm[0] = sqrt(abs(mass)) # possibility of negative fermion masses 700 sqm[1] = sign(sqm[0],mass) # possibility of negative fermion masses 701 ip = -((1+nh)/2) 702 im = (1-nh)/2 703 fop[0] = im * sqm[im] 704 fop[1] = ip*nsr * sqm[im] 705 fop[2] = im*nsr * sqm[-ip] 706 fop[3] = ip * sqm[-ip] 707 else: 708 pp = min(p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)) 709 sf[0] = (1+nsr+(1-nsr)*nh)*0.5 710 sf[1] = (1+nsr-(1-nsr)*nh)*0.5 711 omega[0] = sqrt(p[0]+pp) 712 omega[1] = mass/omega[0] 713 ip = (3+nh)//2 -1 # -1 since this is index 714 im = (3-nh)//2 -1 # -1 since this is index 715 sfomeg[0] = sf[0]*omega[ip] 716 sfomeg[1] = sf[1]*omega[im] 717 pp3 = max([pp+p[3],0]) 718 chi[0] = sqrt(pp3*0.5/pp) 719 if pp3 == 0: 720 chi[1] = -nh 721 else: 722 chi[1] = complex( nh*p[1] , -p[2] )/sqrt(2*pp*pp3) 723 724 725 fop[0] = sfomeg[1]*chi[im] 726 fop[1] = sfomeg[1]*chi[ip] 727 fop[2] = sfomeg[0]*chi[im] 728 fop[3] = sfomeg[0]*chi[ip] 729 730 else: 731 if(p[1] == 0 and p[2] == 0 and p[3] < 0): 732 sqp0p3 = 0 733 else: 734 sqp0p3 = sqrt(max(p[0]+p[3], 0))*nsr 735 736 chi[0] = sqp0p3 737 if ( sqp0p3 == 0 ): 738 chi[1] = complex(-nhel )*sqrt(2*p[0]) 739 else: 740 chi[1] = complex( nh*p[1], -p[2] )/sqp0p3 741 742 if ( nh == 1 ): 743 fop[0] = chi[0] 744 fop[1] = chi[1] 745 #fop[2] = 0 746 #fop[3] = 0 747 else: 748 #fop[0] = 0 749 #fop[1] = 0 750 fop[2] = chi[1] 751 fop[3] = chi[0] 752 753 754 if ( nhel < 2 ): 755 # constract spinor+ 756 sqm, fom, omega, sf, sfomeg = [0]*2,[0]*4,[0]*2,[0]*2,[0]*2 757 chi = [0]*2 758 759 760 nh = -nsr 761 if mass: 762 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 763 if ( pp == 0): 764 sqm[0] = sqrt(abs(mass)) # possibility of negative fermion masses 765 sqm[1] = sign(sqm[0],mass) # possibility of negative fermion masses 766 ip = -((1+nh)/2) 767 im = (1-nh)/2 768 769 fom[0] = im * sqm[im] 770 fom[1] = ip*nsr * sqm[im] 771 fom[2] = im*nsr * sqm[-ip] 772 fom[3] = ip * sqm[-ip] 773 774 else: 775 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 776 sf[0] = (1+nsr+(1-nsr)*nh)*0.5 777 sf[1] = (1+nsr-(1-nsr)*nh)*0.5 778 omega[0] = sqrt(p[0]+pp) 779 omega[1] = mass/omega[0] 780 ip = (3+nh)//2 -1 #-1 since ip is an index 781 im = (3-nh)//2 -1 782 sfomeg[0] = sf[0]*omega[ip] 783 sfomeg[1] = sf[1]*omega[im] 784 pp3 = max([pp+p[3], 0]) 785 chi[0] = sqrt(pp3*0.5/pp) 786 if ( pp3 == 0): 787 chi[1] = -nh 788 else: 789 chi[1] = complex( nh*p[1] , -p[2] )/sqrt(2*pp*pp3) 790 791 792 fom[0] = sfomeg[1]*chi[im] 793 fom[1] = sfomeg[1]*chi[ip] 794 fom[2] = sfomeg[0]*chi[im] 795 fom[3] = sfomeg[0]*chi[ip] 796 else: 797 if(p[1] == 0 == p[2] and p[3] < 0): 798 sqp0p3 = 0 799 else: 800 sqp0p3 = sqrt(max([p[0]+p[3],0]))*nsr 801 chi[0] = sqp0p3 802 if ( sqp0p3 == 0): 803 chi[1] = complex(-nhel )*sqrt(2*p[0]) 804 else: 805 chi[1] = complex( nh*p[1], -p[2] )/sqp0p3 806 if ( nh == 1 ): 807 fom[0] = chi[0] 808 fom[1] = chi[1] 809 #fom[2] = 0 810 #fom[3] = 0 811 else: 812 #fom[1] = 0 813 #fom[2] = 0 814 fom[2] = chi[1] 815 fom[3] = chi[0] 816 817 cond1 = ( pt==0 and p[3]>=0) 818 cond2= (pt==0 and p[3]<0) 819 820 821 # spin-3/2 fermion wavefunction 822 if nhel == 3: 823 for i,j in product(range(4), range(4)): 824 rc[(i, j)] = ep[i] *fop[j] 825 826 827 elif nhel == 1: 828 for i,j in product(range(4), range(4)): 829 if cond1: 830 rc[(i,j)] = sq2/sq3*e0[i]*fop[j] + 1/sq3*ep[i]*fom[j] 831 elif cond2: 832 rc[(i,j)] = sq2/sq3*e0[i]*fop[j] - 1/sq3*ep[i]*fom[j] 833 else: 834 rc[(i,j)] = sq2/sq3*e0[i]*fop[j] + 1/sq3*ep[i]*fom[j] * \ 835 complex(p[1],-nsr*p[2])/pt 836 837 elif nhel == -1: 838 for i,j in product(range(4), range(4)): 839 if cond1: 840 rc[(i,j)] = 1/sq3*em[i]*fop[j]+sq2/sq3*e0[i]*fom[j] 841 elif cond2: 842 rc[(i,j)] =1/sq3*em[i]*fop[j]-sq2/sq3*e0[i]*fom[j] 843 else: 844 rc[(i,j)] = 1/sq3*em[i]*fop[j] + sq2/sq3*e0[i]*fom[j] *\ 845 complex(p[1],-nsr*p[2])/pt 846 else: 847 for i,j in product(range(4), range(4)): 848 if cond1: 849 rc[(i,j)] = em[i] * fom[j] 850 elif cond2: 851 rc[(i,j)] = - em[i] * fom[j] 852 else: 853 rc[(i,j)] = em[i] * fom[j] * complex(p[1],-nsr*p[2])/pt 854 855 856 857 ro[2] = rc[(0,0)] 858 ro[3] = rc[(0,1)] 859 ro[4] = rc[(0,2)] 860 ro[5] = rc[(0,3)] 861 ro[6] = rc[(1,0)] 862 ro[7] = rc[(1,1)] 863 ro[8] = rc[(1,2)] 864 ro[9] = rc[(1,3)] 865 ro[10] = rc[(2,0)] 866 ro[11] = rc[(2,1)] 867 ro[12] = rc[(2,2)] 868 ro[13] = rc[(2,3)] 869 ro[14] = rc[(3,0)] 870 ro[15] = rc[(3,1)] 871 ro[16] = rc[(3,2)] 872 ro[17] = rc[(3,3)] 873 ro[0] = rc[(4,0)] 874 ro[1] = rc[(5,0)] 875 876 return ro
877