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 if (y < 0.): 211 return -abs(x) 212 else: 213 return abs(x)
214 215
216 -def sxxxxx(p,nss):
217 """initialize a scalar wavefunction""" 218 219 sc = WaveFunction(1) 220 221 sc[2] = complex(1.,0.) 222 sc[0] = complex(p[0]*nss,p[3]*nss) 223 sc[1] = complex(p[1]*nss,p[2]*nss) 224 return sc
225 226
227 -def txxxxx(p, tmass, nhel, nst):
228 """ initialize a tensor wavefunction""" 229 230 tc = WaveFunction(5) 231 232 sqh = sqrt(0.5) 233 sqs = sqrt(1/6) 234 235 pt2 = p[1]**2 + p[2]**2 236 pp = min(p[0],sqrt(pt2+p[3]**2)) 237 pt = min(pp,sqrt(pt2)) 238 239 ft = {} 240 ft[(4,0)] = complex(p[0], p[3]) * nst 241 ft[(5,0)] = complex(p[1], p[2]) * nst 242 243 if ( nhel >= 0 ): 244 #construct eps+ 245 ep = [0] * 4 246 247 if ( pp == 0 ): 248 #ep[0] = 0 249 ep[1] = -sqh 250 ep[2] = complex(0, nst*sqh) 251 #ep[3] = 0 252 else: 253 #ep[0] = 0 254 ep[3] = pt/pp*sqh 255 if (pt != 0): 256 pzpt = p[3]/(pp*pt)*sqh 257 ep[1] = complex( -p[1]*pzpt , -nst*p[2]/pt*sqh ) 258 ep[2] = complex( -p[2]*pzpt , nst*p[1]/pt*sqh ) 259 else: 260 ep[1] = -sqh 261 ep[2] = complex( 0 , nst*sign(sqh,p[3]) ) 262 263 264 265 if ( nhel <= 0 ): 266 #construct eps- 267 em = [0] * 4 268 if ( pp == 0 ): 269 #em[0] = 0 270 em[1] = sqh 271 em[2] = complex( 0 , nst*sqh ) 272 #em[3] = 0 273 else: 274 #em[0] = 0 275 em[3] = -pt/pp*sqh 276 if pt: 277 pzpt = -p[3]/(pp*pt)*sqh 278 em[1] = complex( -p[1]*pzpt , -nst*p[2]/pt*sqh ) 279 em[2] = complex( -p[2]*pzpt , nst*p[1]/pt*sqh ) 280 else: 281 em[1] = sqh 282 em[2] = complex( 0 , nst*sign(sqh,p[3]) ) 283 284 285 if ( abs(nhel) <= 1 ): 286 #construct eps0 287 e0 = [0] * 4 288 if ( pp == 0 ): 289 #e0[0] = complex( rZero ) 290 #e0[1] = complex( rZero ) 291 #e0[2] = complex( rZero ) 292 e0[3] = 1 293 else: 294 emp = p[0]/(tmass*pp) 295 e0[0] = pp/tmass 296 e0[3] = p[3]*emp 297 if pt: 298 e0[1] = p[1]*emp 299 e0[2] = p[2]*emp 300 #else: 301 # e0[1] = complex( rZero ) 302 # e0[2] = complex( rZero ) 303 304 if nhel == 2: 305 for j in range(4): 306 for i in range(4): 307 ft[(i,j)] = ep[i]*ep[j] 308 elif nhel == -2: 309 for j in range(4): 310 for i in range(4): 311 ft[(i,j)] = em[i]*em[j] 312 elif tmass == 0: 313 for j in range(4): 314 for i in range(4): 315 ft[(i,j)] = 0 316 elif nhel == 1: 317 for j in range(4): 318 for i in range(4): 319 ft[(i,j)] = sqh*( ep[i]*e0[j] + e0[i]*ep[j] ) 320 elif nhel == 0: 321 for j in range(4): 322 for i in range(4): 323 ft[(i,j)] = sqs*( ep[i]*em[j] + em[i]*ep[j] + 2 *e0[i]*e0[j] ) 324 elif nhel == -1: 325 for j in range(4): 326 for i in range(4): 327 ft[(i,j)] = sqh*( em[i]*e0[j] + e0[i]*em[j] ) 328 329 else: 330 raise Exception, 'invalid helicity TXXXXXX' 331 332 333 334 tc[2] = ft[(0,0)] 335 tc[3] = ft[(0,1)] 336 tc[4] = ft[(0,2)] 337 tc[5] = ft[(0,3)] 338 tc[6] = ft[(1,0)] 339 tc[7] = ft[(1,1)] 340 tc[8] = ft[(1,2)] 341 tc[9] = ft[(1,3)] 342 tc[10] = ft[(2,0)] 343 tc[11] = ft[(2,1)] 344 tc[12] = ft[(2,2)] 345 tc[13] = ft[(2,3)] 346 tc[14] = ft[(3,0)] 347 tc[15] = ft[(3,1)] 348 tc[16] = ft[(3,2)] 349 tc[17] = ft[(3,3)] 350 tc[0] = ft[(4,0)] 351 tc[1] = ft[(5,0)] 352 353 return tc
354
355 -def irxxxx(p, mass, nhel, nsr):
356 """ initialize a incoming spin 3/2 wavefunction.""" 357 358 # This routines is a python translation of a routine written by 359 # K.Mawatari in fortran dated from the 2008/02/26 360 361 ri = WaveFunction(4) 362 363 sqh = sqrt(0.5) 364 sq2 = sqrt(2) 365 sq3 = sqrt(3) 366 367 pt2 = p[1]**2 + p[2]**2 368 pp = min([p[0], sqrt(pt2+p[3]**2)]) 369 pt = min([pp,sqrt(pt2)]) 370 371 rc = {} 372 rc[(4,0)] = -1*complex(p[0],p[3])*nsr 373 rc[(5,0)] = -1*complex(p[1],p[2])*nsr 374 375 376 nsv = -nsr # nsv=+1 for final, -1 for initial 377 378 # Construct eps+ 379 if nhel > 0: 380 ep = [0] * 4 381 if pp: 382 #ep[0] = 0 383 ep[3] = pt/pp*sqh 384 if pt: 385 pzpt = p[3]/(pp*pt)*sqh 386 ep[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 387 ep[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 388 else: 389 ep[1] = -sqh 390 ep[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 391 else: 392 #ep[0] = 0d0 393 ep[1] = -sqh 394 ep[2] = complex(0, nsv * sqh) 395 #ep[3] = 0d0 396 397 398 if ( nhel < 0 ): 399 #construct eps- 400 em = [0] * 4 401 if ( pp == 0 ): 402 #em[0] = 0 403 em[1] = sqh 404 em[2] = complex( 0 , nsv*sqh ) 405 #em[3] = 0 406 else: 407 #em[0] = 0 408 em[3] = -pt/pp*sqh 409 if pt: 410 pzpt = -p[3]/(pp*pt)*sqh 411 em[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 412 em[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 413 else: 414 em[1] = sqh 415 em[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 416 417 if ( abs(nhel) <= 1 ): 418 #construct eps0 419 e0 = [0] * 4 420 if ( pp == 0 ): 421 #e0[0] = complex( rZero ) 422 #e0[1] = complex( rZero ) 423 #e0[2] = complex( rZero ) 424 e0[3] = 1 425 else: 426 emp = p[0]/(mass*pp) 427 e0[0] = pp/mass 428 e0[3] = p[3]*emp 429 if pt: 430 e0[1] = p[1]*emp 431 e0[2] = p[2]*emp 432 #else: 433 # e0[1] = complex( rZero ) 434 # e0[2] = complex( rZero ) 435 436 437 438 if ( nhel >= -1 ): 439 # constract spinor+ 440 fip = [0] * 4 441 sf, omega, sfomeg, chi = [0, 0], [0,0], [0,0], [0,0] 442 nh = nsr 443 if mass: 444 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 445 if pp == 0: 446 sqm = sqrt(mass) 447 ip = (1+nh)//2 448 im = (1-nh)//2 449 fip[0] = ip * sqm 450 fip[1] = im * nsr * sqm 451 fip[2] = ip * nsr * sqm 452 fip[3] = im * sqm 453 else: 454 sf[0] = float(1+nsr+(1-nsr)*nh)*0.5 455 sf[1] = float(1+nsr-(1-nsr)*nh)*0.5 456 omega[0] = sqrt(p[0]+pp) 457 omega[1] = mass/omega[0] 458 ip = ((3+nh)//2) -1 # -1 since they are index 459 im = ((3-nh)//2) -1 # -1 since they are index 460 sfomeg[0] = sf[0]*omega[ip] 461 sfomeg[1] = sf[1]*omega[im] 462 pp3 = max([pp+p[3],0]) 463 chi[0] = sqrt(pp3*0.5/pp) 464 if pp3 ==0: 465 chi[1] = -nh 466 else: 467 chi[1] = complex( nh*p[1] , p[2] )/sqrt(2*pp*pp3) 468 469 fip[0] = sfomeg[0]*chi[im] 470 fip[1] = sfomeg[0]*chi[ip] 471 fip[2] = sfomeg[1]*chi[im] 472 fip[3] = sfomeg[1]*chi[ip] 473 else: 474 sqp0p3 = sqrt(max([p[0]+p[3],0])) * nsr 475 chi[0] = sqp0p3 476 if sqp0p3 == 0: 477 chi[1] = -nhel * sqrt(2*p[0]) 478 else: 479 chi[1] = complex( nh*p[1], p[2] )/sqp0p3 480 if nh == 1: 481 #fip[0] = complex( rZero ) 482 #fip[1] = complex( rZero ) 483 fip[2] = chi[0] 484 fip[3] = chi[1] 485 else: 486 fip[0] = chi[1] 487 fip[1] = chi[0] 488 #fip(3) = complex( rZero ) 489 #fip(4) = complex( rZero ) 490 491 if ( nhel <= 1 ): 492 # constract spinor- 493 fim = [0] * 4 494 sf, omega, sfomeg, chi = [0, 0], [0,0], [0,0], [0,0] 495 nh = -nsr 496 if mass: 497 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 498 if pp == 0: 499 sqm = sqrt(mass) 500 ip = (1+nh)/2 501 im = (1-nh)/2 502 fim[0] = ip * sqm 503 fim[1] = im * nsr * sqm 504 fim[2] = ip * nsr * sqm 505 fim[3] = im * sqm 506 else: 507 sf[0] = float(1+nsr+(1-nsr)*nh)*0.5 508 sf[1] = float(1+nsr-(1-nsr)*nh)*0.5 509 omega[0] = sqrt(p[0]+pp) 510 omega[1] = mass/omega[0] 511 ip = (3+nh)//2 -1 512 im = (3-nh)//2 -1 513 sfomeg[0] = sf[0]*omega[ip] 514 sfomeg[1] = sf[1]*omega[im] 515 pp3 = max([pp+p[3],0]) 516 chi[0] = sqrt(pp3*0.5/pp) 517 if pp3 ==0: 518 chi[1] = -nh 519 else: 520 chi[1] = complex( nh*p[1] , p[2] )/sqrt(2*pp*pp3) 521 522 fim[0] = sfomeg[0]*chi[im] 523 fim[1] = sfomeg[0]*chi[ip] 524 fim[2] = sfomeg[1]*chi[im] 525 fim[3] = sfomeg[1]*chi[ip] 526 else: 527 sqp0p3 = sqrt(max([p[0]+p[3],0])) * nsr 528 chi[0] = sqp0p3 529 if sqp0p3 == 0: 530 chi[1] = -nhel * sqrt(2*p[0]) 531 else: 532 chi[1] = complex( nh*p[1], p[2] )/sqp0p3 533 if nh == 1: 534 #fip[0] = complex( rZero ) 535 #fip[1] = complex( rZero ) 536 fim[2] = chi[0] 537 fim[3] = chi[1] 538 else: 539 fim[0] = chi[1] 540 fim[1] = chi[0] 541 #fip(3) = complex( rZero ) 542 #fip(4) = complex( rZero ) 543 544 545 546 # recurent relation put her for optimization 547 cond1 = (pt == 0 and p[3] >= 0) 548 cond2 = (pt == 0 and p[3] < 0) 549 550 # spin-3/2 fermion wavefunction 551 if nhel == 3: 552 for i,j in product(range(4), range(4)): 553 rc[(i, j)] = ep[i] *fip[j] 554 555 elif nhel == 1: 556 for i,j in product(range(4), range(4)): 557 if cond1: 558 rc[(i,j)] = sq2/sq3*e0[i]*fip[j] +1/sq3*ep[i]*fim[j] 559 elif cond2: 560 rc[(i,j)] = sq2/sq3*e0[i]*fip[j] -1/sq3*ep[i]*fim[j] 561 else: 562 rc[(i,j)] = sq2/sq3*e0[i]*fip[j] + \ 563 1/sq3*ep[i]*fim[j] *complex(p[1],nsr*p[2])/pt 564 elif nhel == -1: 565 for i,j in product(range(4), range(4)): 566 if cond1: 567 rc[(i,j)] = 1/sq3*em[i]*fip[j] +sq2/sq3*e0[i]*fim[j] 568 elif cond2: 569 rc[(i,j)] = 1/sq3*em[i]*fip[j] -sq2/sq3*e0[i]*fim[j] 570 else: 571 rc[(i,j)] = 1/sq3*em[i]*fip[j] + \ 572 sq2/sq3*e0[i]*fim[j] *complex(p[1],nsr*p[2])/pt 573 else: 574 for i,j in product(range(4), range(4)): 575 if cond1: 576 rc[(i, j)] = em[i] *fim[j] 577 elif cond2: 578 rc[(i, j)] = -em[i] *fim[j] 579 else: 580 rc[(i, j)] = em[i]*fim[j] *complex(p[1],nsr*p[2])/pt 581 582 ri[2] = rc[(0,0)] 583 ri[3] = rc[(0,1)] 584 ri[4] = rc[(0,2)] 585 ri[5] = rc[(0,3)] 586 ri[6] = rc[(1,0)] 587 ri[7] = rc[(1,1)] 588 ri[8] = rc[(1,2)] 589 ri[9] = rc[(1,3)] 590 ri[10] = rc[(2,0)] 591 ri[11] = rc[(2,1)] 592 ri[12] = rc[(2,2)] 593 ri[13] = rc[(2,3)] 594 ri[14] = rc[(3,0)] 595 ri[15] = rc[(3,1)] 596 ri[16] = rc[(3,2)] 597 ri[17] = rc[(3,3)] 598 ri[0] = rc[(4,0)] 599 ri[1] = rc[(5,0)] 600 601 return ri
602
603 -def orxxxx(p, mass, nhel, nsr):
604 """ initialize a incoming spin 3/2 wavefunction.""" 605 606 # This routines is a python translation of a routine written by 607 # K.Mawatari in fortran dated from the 2008/02/26 608 609 610 ro = WaveFunction(spin=4) 611 612 sqh = sqrt(0.5) 613 sq2 = sqrt(2) 614 sq3 = sqrt(3) 615 616 pt2 = p[1]**2 + p[2]**2 617 pp = min([p[0], sqrt(pt2+p[3]**2)]) 618 pt = min([pp,sqrt(pt2)]) 619 rc = {} 620 rc[(4,0)] = complex(p[0],p[3])*nsr 621 rc[(5,0)] = complex(p[1],p[2])*nsr 622 623 624 nsv = nsr # nsv=+1 for final, -1 for initial 625 626 # Construct eps+ 627 if nhel > 0: 628 ep = [0] * 4 629 if pp: 630 #ep[0] = 0 631 ep[3] = pt/pp*sqh 632 if pt: 633 pzpt = p[3]/(pp*pt)*sqh 634 ep[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 635 ep[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 636 else: 637 ep[1] = -sqh 638 ep[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 639 else: 640 #ep[0] = 0d0 641 ep[1] = -sqh 642 ep[2] = complex(0, nsv * sqh) 643 #ep[3] = 0d0 644 645 if ( nhel < 0 ): 646 #construct eps- 647 em = [0] * 4 648 if ( pp == 0 ): 649 #em[0] = 0 650 em[1] = sqh 651 em[2] = complex( 0 , nsv*sqh ) 652 #em[3] = 0 653 else: 654 #em[0] = 0 655 em[3] = -pt/pp*sqh 656 if pt: 657 pzpt = -p[3]/(pp*pt)*sqh 658 em[1] = complex( -p[1]*pzpt , -nsv*p[2]/pt*sqh ) 659 em[2] = complex( -p[2]*pzpt , nsv*p[1]/pt*sqh ) 660 else: 661 em[1] = sqh 662 em[2] = complex( 0 , nsv*sign(sqh,p[3]) ) 663 664 if ( abs(nhel) <= 1 ): 665 #construct eps0 666 e0 = [0] * 4 667 if ( pp == 0 ): 668 #e0[0] = complex( rZero ) 669 #e0[1] = complex( rZero ) 670 #e0[2] = complex( rZero ) 671 e0[3] = 1 672 else: 673 emp = p[0]/(mass*pp) 674 e0[0] = pp/mass 675 e0[3] = p[3]*emp 676 if pt: 677 e0[1] = p[1]*emp 678 e0[2] = p[2]*emp 679 #else: 680 # e0[1] = complex( rZero ) 681 # e0[2] = complex( rZero ) 682 683 if nhel >= -1: 684 #constract spinor+ 685 nh = nsr 686 sqm, fop, omega, sf, sfomeg = [0]*2,[0]*4,[0]*2,[0]*2,[0]*2 687 chi = [0]*2 688 if mass: 689 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 690 if ( pp == 0): 691 sqm[0] = sqrt(abs(mass)) # possibility of negative fermion masses 692 sqm[1] = sign(sqm[0],mass) # possibility of negative fermion masses 693 ip = -((1+nh)/2) 694 im = (1-nh)/2 695 fop[0] = im * sqm[im] 696 fop[1] = ip*nsr * sqm[im] 697 fop[2] = im*nsr * sqm[-ip] 698 fop[3] = ip * sqm[-ip] 699 else: 700 pp = min(p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)) 701 sf[0] = (1+nsr+(1-nsr)*nh)*0.5 702 sf[1] = (1+nsr-(1-nsr)*nh)*0.5 703 omega[0] = sqrt(p[0]+pp) 704 omega[1] = mass/omega[0] 705 ip = (3+nh)//2 -1 # -1 since this is index 706 im = (3-nh)//2 -1 # -1 since this is index 707 sfomeg[0] = sf[0]*omega[ip] 708 sfomeg[1] = sf[1]*omega[im] 709 pp3 = max([pp+p[3],0]) 710 chi[0] = sqrt(pp3*0.5/pp) 711 if pp3 == 0: 712 chi[1] = -nh 713 else: 714 chi[1] = complex( nh*p[1] , -p[2] )/sqrt(2*pp*pp3) 715 716 717 fop[0] = sfomeg[1]*chi[im] 718 fop[1] = sfomeg[1]*chi[ip] 719 fop[2] = sfomeg[0]*chi[im] 720 fop[3] = sfomeg[0]*chi[ip] 721 722 else: 723 if(p[1] == 0 and p[2] == 0 and p[3] < 0): 724 sqp0p3 = 0 725 else: 726 sqp0p3 = sqrt(max(p[0]+p[3], 0))*nsr 727 728 chi[0] = sqp0p3 729 if ( sqp0p3 == 0 ): 730 chi[1] = complex(-nhel )*sqrt(2*p[0]) 731 else: 732 chi[1] = complex( nh*p[1], -p[2] )/sqp0p3 733 734 if ( nh == 1 ): 735 fop[0] = chi[0] 736 fop[1] = chi[1] 737 #fop[2] = 0 738 #fop[3] = 0 739 else: 740 #fop[0] = 0 741 #fop[1] = 0 742 fop[2] = chi[1] 743 fop[3] = chi[0] 744 745 746 if ( nhel < 2 ): 747 # constract spinor+ 748 sqm, fom, omega, sf, sfomeg = [0]*2,[0]*4,[0]*2,[0]*2,[0]*2 749 chi = [0]*2 750 751 752 nh = -nsr 753 if mass: 754 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 755 if ( pp == 0): 756 sqm[0] = sqrt(abs(mass)) # possibility of negative fermion masses 757 sqm[1] = sign(sqm[0],mass) # possibility of negative fermion masses 758 ip = -((1+nh)/2) 759 im = (1-nh)/2 760 761 fom[0] = im * sqm[im] 762 fom[1] = ip*nsr * sqm[im] 763 fom[2] = im*nsr * sqm[-ip] 764 fom[3] = ip * sqm[-ip] 765 766 else: 767 pp = min([p[0],sqrt(p[1]**2+p[2]**2+p[3]**2)]) 768 sf[0] = (1+nsr+(1-nsr)*nh)*0.5 769 sf[1] = (1+nsr-(1-nsr)*nh)*0.5 770 omega[0] = sqrt(p[0]+pp) 771 omega[1] = mass/omega[0] 772 ip = (3+nh)//2 -1 #-1 since ip is an index 773 im = (3-nh)//2 -1 774 sfomeg[0] = sf[0]*omega[ip] 775 sfomeg[1] = sf[1]*omega[im] 776 pp3 = max([pp+p[3], 0]) 777 chi[0] = sqrt(pp3*0.5/pp) 778 if ( pp3 == 0): 779 chi[1] = -nh 780 else: 781 chi[1] = complex( nh*p[1] , -p[2] )/sqrt(2*pp*pp3) 782 783 784 fom[0] = sfomeg[1]*chi[im] 785 fom[1] = sfomeg[1]*chi[ip] 786 fom[2] = sfomeg[0]*chi[im] 787 fom[3] = sfomeg[0]*chi[ip] 788 else: 789 if(p[1] == 0 == p[2] and p[3] < 0): 790 sqp0p3 = 0 791 else: 792 sqp0p3 = sqrt(max([p[0]+p[3],0]))*nsr 793 chi[0] = sqp0p3 794 if ( sqp0p3 == 0): 795 chi[1] = complex(-nhel )*sqrt(2*p[0]) 796 else: 797 chi[1] = complex( nh*p[1], -p[2] )/sqp0p3 798 if ( nh == 1 ): 799 fom[0] = chi[0] 800 fom[1] = chi[1] 801 #fom[2] = 0 802 #fom[3] = 0 803 else: 804 #fom[1] = 0 805 #fom[2] = 0 806 fom[2] = chi[1] 807 fom[3] = chi[0] 808 809 cond1 = ( pt==0 and p[3]>=0) 810 cond2= (pt==0 and p[3]<0) 811 812 813 # spin-3/2 fermion wavefunction 814 if nhel == 3: 815 for i,j in product(range(4), range(4)): 816 rc[(i, j)] = ep[i] *fop[j] 817 818 819 elif nhel == 1: 820 for i,j in product(range(4), range(4)): 821 if cond1: 822 rc[(i,j)] = sq2/sq3*e0[i]*fop[j] + 1/sq3*ep[i]*fom[j] 823 elif cond2: 824 rc[(i,j)] = sq2/sq3*e0[i]*fop[j] - 1/sq3*ep[i]*fom[j] 825 else: 826 rc[(i,j)] = sq2/sq3*e0[i]*fop[j] + 1/sq3*ep[i]*fom[j] * \ 827 complex(p[1],-nsr*p[2])/pt 828 829 elif nhel == -1: 830 for i,j in product(range(4), range(4)): 831 if cond1: 832 rc[(i,j)] = 1/sq3*em[i]*fop[j]+sq2/sq3*e0[i]*fom[j] 833 elif cond2: 834 rc[(i,j)] =1/sq3*em[i]*fop[j]-sq2/sq3*e0[i]*fom[j] 835 else: 836 rc[(i,j)] = 1/sq3*em[i]*fop[j] + sq2/sq3*e0[i]*fom[j] *\ 837 complex(p[1],-nsr*p[2])/pt 838 else: 839 for i,j in product(range(4), range(4)): 840 if cond1: 841 rc[(i,j)] = em[i] * fom[j] 842 elif cond2: 843 rc[(i,j)] = - em[i] * fom[j] 844 else: 845 rc[(i,j)] = em[i] * fom[j] * complex(p[1],-nsr*p[2])/pt 846 847 848 849 ro[2] = rc[(0,0)] 850 ro[3] = rc[(0,1)] 851 ro[4] = rc[(0,2)] 852 ro[5] = rc[(0,3)] 853 ro[6] = rc[(1,0)] 854 ro[7] = rc[(1,1)] 855 ro[8] = rc[(1,2)] 856 ro[9] = rc[(1,3)] 857 ro[10] = rc[(2,0)] 858 ro[11] = rc[(2,1)] 859 ro[12] = rc[(2,2)] 860 ro[13] = rc[(2,3)] 861 ro[14] = rc[(3,0)] 862 ro[15] = rc[(3,1)] 863 ro[16] = rc[(3,2)] 864 ro[17] = rc[(3,3)] 865 ro[0] = rc[(4,0)] 866 ro[1] = rc[(5,0)] 867 868 return ro
869