1 | // HelicityMatrixElements.cc is a part of the PYTHIA event generator. |
---|
2 | // Copyright (C) 2012 Philip Ilten, Torbjorn Sjostrand. |
---|
3 | // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details. |
---|
4 | // Please respect the MCnet Guidelines, see GUIDELINES for details. |
---|
5 | |
---|
6 | // Function definitions (not found in the header) for physics classes |
---|
7 | // used in tau decays. |
---|
8 | |
---|
9 | #include "HelicityMatrixElements.h" |
---|
10 | |
---|
11 | namespace Pythia8 { |
---|
12 | |
---|
13 | //========================================================================== |
---|
14 | |
---|
15 | // The HelicityMatrixElements class. |
---|
16 | |
---|
17 | //-------------------------------------------------------------------------- |
---|
18 | |
---|
19 | // Initialize the helicity matrix element. |
---|
20 | |
---|
21 | void HelicityMatrixElement::initPointers(ParticleData* particleDataPtrIn, |
---|
22 | Couplings* couplingsPtrIn) { |
---|
23 | |
---|
24 | particleDataPtr = particleDataPtrIn; |
---|
25 | couplingsPtr = couplingsPtrIn; |
---|
26 | for(int i = 0; i <= 5; i++) |
---|
27 | gamma.push_back(GammaMatrix(i)); |
---|
28 | |
---|
29 | } |
---|
30 | |
---|
31 | //-------------------------------------------------------------------------- |
---|
32 | |
---|
33 | // Initialize the channel for the helicity matrix element. |
---|
34 | |
---|
35 | HelicityMatrixElement* HelicityMatrixElement::initChannel( |
---|
36 | vector<HelicityParticle>& p) { |
---|
37 | |
---|
38 | pID.clear(); |
---|
39 | pM.clear(); |
---|
40 | for(int i = 0; i < static_cast<int>(p.size()); i++) { |
---|
41 | pID.push_back(p[i].id()); |
---|
42 | pM.push_back(p[i].m()); |
---|
43 | } |
---|
44 | initConstants(); |
---|
45 | return this; |
---|
46 | |
---|
47 | } |
---|
48 | |
---|
49 | //-------------------------------------------------------------------------- |
---|
50 | |
---|
51 | // Calculate a particle's decay matrix. |
---|
52 | |
---|
53 | void HelicityMatrixElement::calculateD(vector<HelicityParticle>& p) { |
---|
54 | |
---|
55 | // Reset the D matrix to zero. |
---|
56 | for (int i = 0; i < p[0].spinStates(); i++) { |
---|
57 | for (int j = 0; j < p[0].spinStates(); j++) { |
---|
58 | p[0].D[i][j] = 0; |
---|
59 | } |
---|
60 | } |
---|
61 | |
---|
62 | // Initialize the wave functions. |
---|
63 | initWaves(p); |
---|
64 | |
---|
65 | // Create the helicity vectors. |
---|
66 | vector<int> h1(p.size(),0); |
---|
67 | vector<int> h2(p.size(),0); |
---|
68 | |
---|
69 | // Call the recursive sub-method. |
---|
70 | calculateD(p, h1, h2, 0); |
---|
71 | |
---|
72 | // Normalize the decay matrix. |
---|
73 | p[0].normalize(p[0].D); |
---|
74 | |
---|
75 | } |
---|
76 | |
---|
77 | //-------------------------------------------------------------------------- |
---|
78 | |
---|
79 | // Recursive sub-method for calculating a particle's decay matrix. |
---|
80 | |
---|
81 | void HelicityMatrixElement::calculateD(vector<HelicityParticle>& p, |
---|
82 | vector<int>& h1, vector<int>& h2, unsigned int i) { |
---|
83 | |
---|
84 | if (i < p.size()) { |
---|
85 | for (h1[i] = 0; h1[i] < p[i].spinStates(); h1[i]++) { |
---|
86 | for (h2[i] = 0; h2[i] < p[i].spinStates(); h2[i]++) { |
---|
87 | calculateD(p, h1, h2, i+1); |
---|
88 | } |
---|
89 | } |
---|
90 | } |
---|
91 | else { |
---|
92 | p[0].D[h1[0]][h2[0]] += calculateME(h1) * conj(calculateME(h2)) * |
---|
93 | calculateProductD(p, h1, h2); |
---|
94 | } |
---|
95 | |
---|
96 | } |
---|
97 | |
---|
98 | //-------------------------------------------------------------------------- |
---|
99 | |
---|
100 | // Calculate a particle's helicity density matrix. |
---|
101 | |
---|
102 | void HelicityMatrixElement::calculateRho(unsigned int idx, |
---|
103 | vector<HelicityParticle>& p) { |
---|
104 | |
---|
105 | // Reset the rho matrix to zero. |
---|
106 | for (int i = 0; i < p[idx].spinStates(); i++) { |
---|
107 | for (int j = 0; j < p[idx].spinStates(); j++) { |
---|
108 | p[idx].rho[i][j] = 0; |
---|
109 | } |
---|
110 | } |
---|
111 | |
---|
112 | // Initialize the wave functions. |
---|
113 | initWaves(p); |
---|
114 | |
---|
115 | // Create the helicity vectors. |
---|
116 | vector<int> h1(p.size(),0); |
---|
117 | vector<int> h2(p.size(),0); |
---|
118 | |
---|
119 | // Call the recursive sub-method. |
---|
120 | calculateRho(idx, p, h1, h2, 0); |
---|
121 | |
---|
122 | // Normalize the density matrix. |
---|
123 | p[idx].normalize(p[idx].rho); |
---|
124 | |
---|
125 | } |
---|
126 | |
---|
127 | //-------------------------------------------------------------------------- |
---|
128 | |
---|
129 | // Recursive sub-method for calculating a particle's helicity density matrix. |
---|
130 | |
---|
131 | void HelicityMatrixElement::calculateRho(unsigned int idx, |
---|
132 | vector<HelicityParticle>& p, vector<int>& h1, vector<int>& h2, |
---|
133 | unsigned int i) { |
---|
134 | |
---|
135 | if (i < p.size()) { |
---|
136 | for (h1[i] = 0; h1[i] < p[i].spinStates(); h1[i]++) { |
---|
137 | for (h2[i] = 0; h2[i] < p[i].spinStates(); h2[i]++) { |
---|
138 | calculateRho(idx, p, h1, h2, i+1); |
---|
139 | } |
---|
140 | } |
---|
141 | } |
---|
142 | else { |
---|
143 | // Calculate rho from a hard process. |
---|
144 | if (p[1].direction < 0) |
---|
145 | p[idx].rho[h1[idx]][h2[idx]] += p[0].rho[h1[0]][h2[0]] * |
---|
146 | p[1].rho[h1[1]][h2[1]] * calculateME(h1)*conj(calculateME(h2)) * |
---|
147 | calculateProductD(idx, 2, p, h1, h2); |
---|
148 | // Calculate rho from a decay. |
---|
149 | else |
---|
150 | p[idx].rho[h1[idx]][h2[idx]] += p[0].rho[h1[0]][h2[0]] * |
---|
151 | calculateME(h1)*conj(calculateME(h2)) * |
---|
152 | calculateProductD(idx, 1, p, h1, h2); |
---|
153 | return; |
---|
154 | } |
---|
155 | |
---|
156 | } |
---|
157 | |
---|
158 | //-------------------------------------------------------------------------- |
---|
159 | |
---|
160 | // Calculate a decay's weight. |
---|
161 | |
---|
162 | double HelicityMatrixElement::decayWeight(vector<HelicityParticle>& p) { |
---|
163 | |
---|
164 | complex weight = complex(0,0); |
---|
165 | |
---|
166 | // Initialize the wave functions. |
---|
167 | initWaves(p); |
---|
168 | |
---|
169 | // Create the helicity vectors. |
---|
170 | vector<int> h1(p.size(),0); |
---|
171 | vector<int> h2(p.size(),0); |
---|
172 | |
---|
173 | // Call the recursive sub-method. |
---|
174 | decayWeight(p, h1, h2, weight, 0); |
---|
175 | |
---|
176 | return real(weight); |
---|
177 | |
---|
178 | } |
---|
179 | |
---|
180 | //-------------------------------------------------------------------------- |
---|
181 | |
---|
182 | // Recursive sub-method for calculating a decay's weight. |
---|
183 | |
---|
184 | void HelicityMatrixElement::decayWeight(vector<HelicityParticle>& p, |
---|
185 | vector<int>& h1, vector<int>& h2, complex& weight, unsigned int i) { |
---|
186 | |
---|
187 | if (i < p.size()) { |
---|
188 | for (h1[i] = 0; h1[i] < p[i].spinStates(); h1[i]++) { |
---|
189 | for (h2[i] = 0; h2[i] < p[i].spinStates(); h2[i]++) { |
---|
190 | decayWeight(p, h1, h2, weight, i+1); |
---|
191 | } |
---|
192 | } |
---|
193 | } |
---|
194 | else { |
---|
195 | weight += p[0].rho[h1[0]][h2[0]] * calculateME(h1) * |
---|
196 | conj(calculateME(h2)) * calculateProductD(p, h1, h2); |
---|
197 | } |
---|
198 | |
---|
199 | } |
---|
200 | |
---|
201 | //-------------------------------------------------------------------------- |
---|
202 | |
---|
203 | // Calculate the product of the decay matrices (hard process). |
---|
204 | |
---|
205 | complex HelicityMatrixElement::calculateProductD(unsigned int idx, |
---|
206 | unsigned int start, vector<HelicityParticle>& p, |
---|
207 | vector<int>& h1, vector<int>& h2) { |
---|
208 | |
---|
209 | complex answer(1,0); |
---|
210 | for (unsigned int i = start; i < p.size(); i++) { |
---|
211 | if (i != idx) { |
---|
212 | answer *= p[i].D[h1[i]][h2[i]]; |
---|
213 | } |
---|
214 | } |
---|
215 | return answer; |
---|
216 | |
---|
217 | } |
---|
218 | |
---|
219 | //-------------------------------------------------------------------------- |
---|
220 | |
---|
221 | // Calculate the product of the decay matrices (decay process). |
---|
222 | |
---|
223 | complex HelicityMatrixElement::calculateProductD( |
---|
224 | vector<HelicityParticle>& p, vector<int>& h1, vector<int>& h2) { |
---|
225 | |
---|
226 | complex answer(1,0); |
---|
227 | for (unsigned int i = 1; i < p.size(); i++) { |
---|
228 | answer *= p[i].D[h1[i]][h2[i]]; |
---|
229 | } |
---|
230 | return answer; |
---|
231 | |
---|
232 | } |
---|
233 | |
---|
234 | //-------------------------------------------------------------------------- |
---|
235 | |
---|
236 | // Initialize a fermion line. |
---|
237 | |
---|
238 | void HelicityMatrixElement::setFermionLine(int position, |
---|
239 | HelicityParticle& p0, HelicityParticle& p1) { |
---|
240 | |
---|
241 | vector< Wave4 > u0, u1; |
---|
242 | |
---|
243 | // First particle is incoming and particle, or outgoing and anti-particle. |
---|
244 | if (p0.id()*p0.direction < 0) { |
---|
245 | pMap[position] = position; pMap[position+1] = position+1; |
---|
246 | for (int h = 0; h < p0.spinStates(); h++) u0.push_back(p0.wave(h)); |
---|
247 | for (int h = 0; h < p1.spinStates(); h++) u1.push_back(p1.waveBar(h)); |
---|
248 | } |
---|
249 | // First particle is outgoing and particle, or incoming and anti-particle. |
---|
250 | else { |
---|
251 | pMap[position] = position+1; pMap[position+1] = position; |
---|
252 | for (int h = 0; h < p0.spinStates(); h++) u1.push_back(p0.waveBar(h)); |
---|
253 | for (int h = 0; h < p1.spinStates(); h++) u0.push_back(p1.wave(h)); |
---|
254 | } |
---|
255 | u.push_back(u0); u.push_back(u1); |
---|
256 | |
---|
257 | } |
---|
258 | |
---|
259 | //-------------------------------------------------------------------------- |
---|
260 | |
---|
261 | // Return a fixed width Breit-Wigner. |
---|
262 | |
---|
263 | complex HelicityMatrixElement::breitWigner(double s, double M, double G) { |
---|
264 | |
---|
265 | return (-M * M + complex(0, 1) * M * G) / (s - M * M + complex(0, 1) * M * G); |
---|
266 | |
---|
267 | } |
---|
268 | |
---|
269 | //-------------------------------------------------------------------------- |
---|
270 | |
---|
271 | // Return an s-wave BreitWigner. |
---|
272 | |
---|
273 | complex HelicityMatrixElement::sBreitWigner(double m0, double m1, double s, |
---|
274 | double M, double G) { |
---|
275 | |
---|
276 | double gs = sqrtpos((s - pow2(m0+m1)) * (s - pow2(m0-m1))) / (2*sqrtpos(s)); |
---|
277 | double gM = sqrtpos((M*M - pow2(m0+m1)) * (M*M - pow2(m0-m1))) / (2*M); |
---|
278 | return M*M / (M*M - s - complex(0,1)*G*M*M/sqrtpos(s)*(gs/gM)); |
---|
279 | |
---|
280 | } |
---|
281 | |
---|
282 | //-------------------------------------------------------------------------- |
---|
283 | |
---|
284 | // Return a p-wave BreitWigner. |
---|
285 | |
---|
286 | complex HelicityMatrixElement::pBreitWigner(double m0, double m1, double s, |
---|
287 | double M, double G) { |
---|
288 | |
---|
289 | double gs = sqrtpos((s - pow2(m0+m1)) * (s - pow2(m0-m1))) / (2*sqrtpos(s)); |
---|
290 | double gM = sqrtpos((M*M - pow2(m0+m1)) * (M*M - pow2(m0-m1))) / (2*M); |
---|
291 | return M*M / (M*M - s - complex(0,1)*G*M*M/sqrtpos(s)*pow3(gs/gM)); |
---|
292 | |
---|
293 | } |
---|
294 | |
---|
295 | //-------------------------------------------------------------------------- |
---|
296 | |
---|
297 | // Return a d-wave BreitWigner. |
---|
298 | |
---|
299 | complex HelicityMatrixElement::dBreitWigner(double m0, double m1, double s, |
---|
300 | double M, double G) { |
---|
301 | |
---|
302 | double gs = sqrtpos((s - pow2(m0+m1)) * (s - pow2(m0-m1))) / (2*sqrtpos(s)); |
---|
303 | double gM = sqrtpos((M*M - pow2(m0+m1)) * (M*M - pow2(m0-m1))) / (2*M); |
---|
304 | return M*M / (M*M - s - complex(0,1)*G*M*M/sqrtpos(s)*pow5(gs/gM)); |
---|
305 | |
---|
306 | } |
---|
307 | |
---|
308 | //========================================================================== |
---|
309 | |
---|
310 | // Helicity matrix element for two fermions -> W -> two fermions. This matrix |
---|
311 | // element handles s-channel hard processes in addition to t-channel, assuming |
---|
312 | // the first two particles are a fermion line and the second two particles |
---|
313 | // are a fermion line. This matrix element is not scaled with respect to W |
---|
314 | // propagator energy as currently this matrix element is used only for |
---|
315 | // calculating helicity density matrices. |
---|
316 | |
---|
317 | //-------------------------------------------------------------------------- |
---|
318 | |
---|
319 | // Initialize spinors for the helicity matrix element. |
---|
320 | |
---|
321 | void HMETwoFermions2W2TwoFermions::initWaves(vector<HelicityParticle>& p) { |
---|
322 | |
---|
323 | u.clear(); |
---|
324 | pMap.resize(4); |
---|
325 | setFermionLine(0,p[0],p[1]); |
---|
326 | setFermionLine(2,p[2],p[3]); |
---|
327 | |
---|
328 | } |
---|
329 | |
---|
330 | //-------------------------------------------------------------------------- |
---|
331 | |
---|
332 | // Return element for the helicity matrix element. |
---|
333 | |
---|
334 | complex HMETwoFermions2W2TwoFermions::calculateME(vector<int> h) { |
---|
335 | |
---|
336 | complex answer(0,0); |
---|
337 | for (int mu = 0; mu <= 3; mu++) { |
---|
338 | answer += (u[1][h[pMap[1]]] * gamma[mu] * (1 - gamma[5]) |
---|
339 | * u[0][h[pMap[0]]]) * gamma[4](mu,mu) * (u[3][h[pMap[3]]] |
---|
340 | * gamma[mu] * (1 - gamma[5]) * u[2][h[pMap[2]]]); |
---|
341 | } |
---|
342 | return answer; |
---|
343 | |
---|
344 | } |
---|
345 | |
---|
346 | //========================================================================== |
---|
347 | |
---|
348 | // Helicity matrix element for two fermions -> photon -> two fermions. This |
---|
349 | // matrix element can be combined with the Z matrix element to provide full |
---|
350 | // interference effects. |
---|
351 | |
---|
352 | // p0Q: charge of the incoming fermion line |
---|
353 | // p2Q: charge of the outgoing fermion line |
---|
354 | // s: center of mass energy |
---|
355 | |
---|
356 | //-------------------------------------------------------------------------- |
---|
357 | |
---|
358 | // Initialize wave functions for the helicity matrix element. |
---|
359 | |
---|
360 | void HMETwoFermions2Gamma2TwoFermions::initWaves( |
---|
361 | vector<HelicityParticle>& p) { |
---|
362 | |
---|
363 | u.clear(); |
---|
364 | pMap.resize(4); |
---|
365 | setFermionLine(0, p[0], p[1]); |
---|
366 | setFermionLine(2, p[2], p[3]); |
---|
367 | s = max( 1., pow2(p[4].m())); |
---|
368 | p0Q = p[0].charge(); p2Q = p[2].charge(); |
---|
369 | |
---|
370 | } |
---|
371 | |
---|
372 | //-------------------------------------------------------------------------- |
---|
373 | |
---|
374 | // Return element for the helicity matrix element. |
---|
375 | |
---|
376 | |
---|
377 | complex HMETwoFermions2Gamma2TwoFermions::calculateME(vector<int> h) { |
---|
378 | |
---|
379 | complex answer(0,0); |
---|
380 | for (int mu = 0; mu <= 3; mu++) { |
---|
381 | answer += (u[1][h[pMap[1]]] * gamma[mu] * u[0][h[pMap[0]]]) |
---|
382 | * gamma[4](mu,mu) * (u[3][h[pMap[3]]] * gamma[mu] * u[2][h[pMap[2]]]); |
---|
383 | } |
---|
384 | return p0Q*p2Q * answer / s; |
---|
385 | |
---|
386 | } |
---|
387 | |
---|
388 | //========================================================================== |
---|
389 | |
---|
390 | // Helicity matrix element for two fermions -> Z -> two fermions. This matrix |
---|
391 | // element can be combined with the photon matrix element to provide full |
---|
392 | // interference effects. |
---|
393 | |
---|
394 | // Note that there is a double contraction in the Z matrix element, which can |
---|
395 | // be very time consuming. If the two incoming fermions are oriented along |
---|
396 | // the z-axis, their helicities must be opposite for a non-zero matrix element |
---|
397 | // term. Consequently, this check is made to help speed up the matrix element. |
---|
398 | |
---|
399 | // sin2W: sine of the Weinberg angle |
---|
400 | // cos2W: cosine of the Weinberg angle |
---|
401 | // zM: on-shell mass of the Z |
---|
402 | // zG: on-shell width of the Z |
---|
403 | // p0CA: axial coupling of particle 0 to the Z |
---|
404 | // p2CA: axial coupling of particle 2 to the Z |
---|
405 | // p0CV: vector coupling of particle 0 to the Z |
---|
406 | // p2CV: vector coupling of particle 2 to the Z |
---|
407 | // zaxis: true if the incoming fermions are oriented along the z-axis |
---|
408 | |
---|
409 | //-------------------------------------------------------------------------- |
---|
410 | |
---|
411 | // Initialize the constant for the helicity matrix element. |
---|
412 | |
---|
413 | void HMETwoFermions2Z2TwoFermions::initConstants() { |
---|
414 | |
---|
415 | // Set the Weinberg angle. |
---|
416 | sin2W = couplingsPtr->sin2thetaW(); |
---|
417 | cos2W = couplingsPtr->cos2thetaW(); |
---|
418 | // Set the on-shell Z mass and width. |
---|
419 | zG = particleDataPtr->mWidth(23); |
---|
420 | zM = particleDataPtr->m0(23); |
---|
421 | // Set the vector and axial couplings to the fermions. |
---|
422 | p0CA = couplingsPtr->af(abs(pID[0])); |
---|
423 | p2CA = couplingsPtr->af(abs(pID[2])); |
---|
424 | p0CV = couplingsPtr->vf(abs(pID[0])); |
---|
425 | p2CV = couplingsPtr->vf(abs(pID[2])); |
---|
426 | |
---|
427 | } |
---|
428 | |
---|
429 | //-------------------------------------------------------------------------- |
---|
430 | |
---|
431 | // Initialize wave functions for the helicity matrix element. |
---|
432 | |
---|
433 | void HMETwoFermions2Z2TwoFermions::initWaves(vector<HelicityParticle>& p) { |
---|
434 | |
---|
435 | vector< Wave4 > u4; |
---|
436 | u.clear(); |
---|
437 | pMap.resize(4); |
---|
438 | setFermionLine(0, p[0], p[1]); |
---|
439 | setFermionLine(2, p[2], p[3]); |
---|
440 | u4.push_back(Wave4(p[2].p() + p[3].p())); |
---|
441 | u.push_back(u4); |
---|
442 | // Center of mass energy. |
---|
443 | s = max( 1., pow2(p[4].m())); |
---|
444 | // Check if incoming fermions are oriented along z-axis. |
---|
445 | zaxis = (p[0].pAbs() == fabs(p[0].pz())) && |
---|
446 | (p[1].pAbs() == fabs(p[1].pz())); |
---|
447 | |
---|
448 | } |
---|
449 | |
---|
450 | //-------------------------------------------------------------------------- |
---|
451 | |
---|
452 | // Return element for helicity matrix element. |
---|
453 | |
---|
454 | complex HMETwoFermions2Z2TwoFermions::calculateME(vector<int> h) { |
---|
455 | |
---|
456 | complex answer(0,0); |
---|
457 | // Return zero if correct helicity conditions. |
---|
458 | if (h[0] == h[1] && zaxis) return answer; |
---|
459 | for (int mu = 0; mu <= 3; mu++) { |
---|
460 | for (int nu = 0; nu <= 3; nu++) { |
---|
461 | answer += |
---|
462 | (u[1][h[pMap[1]]] * gamma[mu] * (p0CV - p0CA * gamma[5]) * |
---|
463 | u[0][h[pMap[0]]]) * |
---|
464 | (gamma[4](mu,nu) - gamma[4](mu,mu)*u[4][0](mu) * |
---|
465 | gamma[4](nu,nu) * u[4][0](nu) / (zM*zM)) * |
---|
466 | (u[3][h[pMap[3]]] * gamma[nu] * (p2CV - p2CA * gamma[5]) * |
---|
467 | u[2][h[pMap[2]]]); |
---|
468 | } |
---|
469 | } |
---|
470 | return answer / (16 * pow2(sin2W * cos2W) * |
---|
471 | (s - zM*zM + complex(0, s*zG/zM))); |
---|
472 | |
---|
473 | } |
---|
474 | |
---|
475 | //========================================================================== |
---|
476 | |
---|
477 | // Helicity matrix element for two fermions -> photon/Z -> two fermions. Full |
---|
478 | // interference is obtained by combining the photon and Z helicity matrix |
---|
479 | // elements. |
---|
480 | |
---|
481 | // In general the initPointers and initChannel methods should not be |
---|
482 | // redeclared. |
---|
483 | |
---|
484 | //-------------------------------------------------------------------------- |
---|
485 | |
---|
486 | // Initialize the matrix element. |
---|
487 | |
---|
488 | void HMETwoFermions2GammaZ2TwoFermions::initPointers( |
---|
489 | ParticleData* particleDataPtrIn, Couplings* couplingsPtrIn) { |
---|
490 | |
---|
491 | zHME.initPointers(particleDataPtrIn, couplingsPtrIn); |
---|
492 | gHME.initPointers(particleDataPtrIn, couplingsPtrIn); |
---|
493 | |
---|
494 | } |
---|
495 | |
---|
496 | //-------------------------------------------------------------------------- |
---|
497 | |
---|
498 | // Initialize the channel for the helicity matrix element. |
---|
499 | |
---|
500 | HelicityMatrixElement* HMETwoFermions2GammaZ2TwoFermions::initChannel( |
---|
501 | vector<HelicityParticle>& p) { |
---|
502 | |
---|
503 | zHME.initChannel(p); |
---|
504 | zHME.initChannel(p); |
---|
505 | return this; |
---|
506 | |
---|
507 | } |
---|
508 | |
---|
509 | //-------------------------------------------------------------------------- |
---|
510 | |
---|
511 | // Initialize wave functions for the helicity matrix element. |
---|
512 | |
---|
513 | void HMETwoFermions2GammaZ2TwoFermions::initWaves( |
---|
514 | vector<HelicityParticle>& p) { |
---|
515 | |
---|
516 | zHME.initWaves(p); |
---|
517 | gHME.initWaves(p); |
---|
518 | |
---|
519 | } |
---|
520 | |
---|
521 | //-------------------------------------------------------------------------- |
---|
522 | |
---|
523 | // Return element for the helicity matrix element. |
---|
524 | |
---|
525 | complex HMETwoFermions2GammaZ2TwoFermions::calculateME(vector<int> h) { |
---|
526 | |
---|
527 | return zHME.calculateME(h) + gHME.calculateME(h); |
---|
528 | |
---|
529 | } |
---|
530 | |
---|
531 | //========================================================================== |
---|
532 | |
---|
533 | // Helicity matrix element for Z -> two fermions. |
---|
534 | |
---|
535 | // Helicity matrix element for Z -> two fermions. This matrix element is used |
---|
536 | // when the production of the Z is from an unknown process. |
---|
537 | |
---|
538 | // p2CA: axial coupling of particle 2 to the Z |
---|
539 | // p2CV: vector coupling of particle 2 to the Z |
---|
540 | |
---|
541 | //-------------------------------------------------------------------------- |
---|
542 | |
---|
543 | // Initialize the constant for the helicity matrix element. |
---|
544 | |
---|
545 | void HMEZ2TwoFermions::initConstants() { |
---|
546 | |
---|
547 | // Set the vector and axial couplings to the fermions. |
---|
548 | p2CA = couplingsPtr->af(abs(pID[2])); |
---|
549 | p2CV = couplingsPtr->vf(abs(pID[2])); |
---|
550 | |
---|
551 | } |
---|
552 | |
---|
553 | //-------------------------------------------------------------------------- |
---|
554 | |
---|
555 | // Initialize wave functions for the helicity matrix element. |
---|
556 | |
---|
557 | void HMEZ2TwoFermions::initWaves(vector<HelicityParticle>& p) { |
---|
558 | |
---|
559 | u.clear(); |
---|
560 | pMap.resize(4); |
---|
561 | // Initialize Z wave function. |
---|
562 | vector< Wave4 > u1; |
---|
563 | pMap[1] = 1; |
---|
564 | for (int h = 0; h < p[pMap[1]].spinStates(); h++) |
---|
565 | u1.push_back(p[pMap[1]].wave(h)); |
---|
566 | u.push_back(u1); |
---|
567 | // Initialize fermion wave functions. |
---|
568 | setFermionLine(2, p[2], p[3]); |
---|
569 | |
---|
570 | } |
---|
571 | |
---|
572 | //-------------------------------------------------------------------------- |
---|
573 | |
---|
574 | // Return element for helicity matrix element. |
---|
575 | |
---|
576 | complex HMEZ2TwoFermions::calculateME(vector<int> h) { |
---|
577 | |
---|
578 | complex answer(0,0); |
---|
579 | for (int mu = 0; mu <= 3; mu++) { |
---|
580 | answer += |
---|
581 | u[0][h[pMap[1]]](mu) * (u[2][h[pMap[3]]] * gamma[mu] |
---|
582 | * (p2CV - p2CA * gamma[5]) * u[1][h[pMap[2]]]); |
---|
583 | } |
---|
584 | return answer; |
---|
585 | } |
---|
586 | |
---|
587 | //========================================================================== |
---|
588 | |
---|
589 | // Helicity matrix element for the decay of a CP even Higgs to two fermions. |
---|
590 | // All SM and MSSM Higgses couple to fermions with a vertex factor of |
---|
591 | // (pfCV - pfCA * gamma[5]) where pf indicates the type of fermion line. For |
---|
592 | // simplicity for the SM and MSSM CP even Higgses pfCV is set to one, and |
---|
593 | // pfCA to zero, as this matrix element is used only for calculating helicity |
---|
594 | // density matrices. |
---|
595 | |
---|
596 | // p2CA: in the SM and MSSM this coupling is zero |
---|
597 | // p2CV: in the SM and MSSM this coupling is given by: |
---|
598 | // i * g_w * m_f / (2 * m_W) |
---|
599 | // * -1 for the SM H |
---|
600 | // * -sin(alpha) / sin(beta) for H^0 u-type |
---|
601 | // * -cos(alpha) / cos(beta) for H^0 d-type |
---|
602 | // * -cos(alpha) / sin(beta) for h^0 u-type |
---|
603 | // * sin(alpha) / cos(beta) for h^0 d-type |
---|
604 | |
---|
605 | //-------------------------------------------------------------------------- |
---|
606 | |
---|
607 | // Initialize wave functions for the helicity matrix element. |
---|
608 | |
---|
609 | void HMEHiggsEven2TwoFermions::initWaves(vector<HelicityParticle>& p) { |
---|
610 | |
---|
611 | u.clear(); |
---|
612 | pMap.resize(4); |
---|
613 | p2CA = 0; p2CV = 1; |
---|
614 | setFermionLine(2, p[2], p[3]); |
---|
615 | |
---|
616 | } |
---|
617 | |
---|
618 | //-------------------------------------------------------------------------- |
---|
619 | |
---|
620 | // Return element for the helicity matrix element. |
---|
621 | |
---|
622 | complex HMEHiggsEven2TwoFermions::calculateME(vector<int> h) { |
---|
623 | |
---|
624 | return (u[1][h[pMap[3]]] * (p2CV - p2CA * gamma[5]) * u[0][h[pMap[2]]]); |
---|
625 | |
---|
626 | } |
---|
627 | |
---|
628 | //========================================================================== |
---|
629 | |
---|
630 | // Helicity matrix element for the decay of a CP odd Higgs to two fermions. |
---|
631 | // See HMEHiggsEven2TwoFermions for more details. For the MSSM CP odd Higgs |
---|
632 | // pfCA is set to one and pfCV is set to zero. |
---|
633 | |
---|
634 | // p2CA: in the MSSM this coupling is given by: |
---|
635 | // -g_w * m_f / (2 * m_W) |
---|
636 | // * cot(beta) for A^0 u-type |
---|
637 | // * tan(beta) for A^0 d-type |
---|
638 | // p2CV: in the MSSM this coupling is zero |
---|
639 | |
---|
640 | //-------------------------------------------------------------------------- |
---|
641 | |
---|
642 | // Initialize wave functions for the helicity matrix element. |
---|
643 | |
---|
644 | void HMEHiggsOdd2TwoFermions::initWaves(vector<HelicityParticle>& p) { |
---|
645 | |
---|
646 | u.clear(); |
---|
647 | pMap.resize(4); |
---|
648 | p2CA = 1; p2CV = 0; |
---|
649 | setFermionLine(2, p[2], p[3]); |
---|
650 | |
---|
651 | } |
---|
652 | |
---|
653 | //-------------------------------------------------------------------------- |
---|
654 | |
---|
655 | // Return element for the helicity matrix element. |
---|
656 | |
---|
657 | complex HMEHiggsOdd2TwoFermions::calculateME(vector<int> h) { |
---|
658 | |
---|
659 | return (u[1][h[pMap[3]]] * (p2CV - p2CA * gamma[5]) * u[0][h[pMap[2]]]); |
---|
660 | |
---|
661 | } |
---|
662 | |
---|
663 | //========================================================================== |
---|
664 | |
---|
665 | // Helicity matrix element for the decay of a charged Higgs to two fermions. |
---|
666 | // See HMEHiggsEven2TwoFermions for more details. For the MSSM charged Higgs |
---|
667 | // pfCA is set to +/- one given an H^+/- and pfCV is set to one. |
---|
668 | |
---|
669 | // p2CA: in the MSSM this coupling is given by: |
---|
670 | // i * g / (sqrt(8) * m_W) * (m_d * tan(beta) + m_u * cot(beta)) |
---|
671 | // p2CV: in the MSSM this coupling is given by: |
---|
672 | // +/- i * g / (sqrt(8) * m_W) * (m_d * tan(beta) - m_u * cot(beta)) |
---|
673 | |
---|
674 | //-------------------------------------------------------------------------- |
---|
675 | |
---|
676 | // Initialize wave functions for the helicity matrix element. |
---|
677 | |
---|
678 | void HMEHiggsCharged2TwoFermions::initWaves(vector<HelicityParticle>& p) { |
---|
679 | |
---|
680 | u.clear(); |
---|
681 | pMap.resize(4); |
---|
682 | p2CV = 1; |
---|
683 | if (pID[3] == 15 || pID[3] == -16) p2CA = 1; |
---|
684 | else p2CA = -1; |
---|
685 | setFermionLine(2, p[2], p[3]); |
---|
686 | |
---|
687 | } |
---|
688 | |
---|
689 | //-------------------------------------------------------------------------- |
---|
690 | |
---|
691 | // Return element for the helicity matrix element. |
---|
692 | |
---|
693 | complex HMEHiggsCharged2TwoFermions::calculateME(vector<int> h) { |
---|
694 | |
---|
695 | return (u[1][h[pMap[3]]] * (p2CV - p2CA * gamma[5]) * u[0][h[pMap[2]]]); |
---|
696 | |
---|
697 | } |
---|
698 | |
---|
699 | //========================================================================== |
---|
700 | |
---|
701 | // Helicity matrix element which provides an unpolarized helicity |
---|
702 | // density matrix. This matrix element is used for unkown hard processes. |
---|
703 | |
---|
704 | // Note that calculateRho is redefined for this special case, but that in |
---|
705 | // general calculateRho should not be redefined. |
---|
706 | |
---|
707 | //-------------------------------------------------------------------------- |
---|
708 | |
---|
709 | // Calculate a particle's helicity density matrix. |
---|
710 | |
---|
711 | void HMEUnpolarized::calculateRho(unsigned int idx, |
---|
712 | vector<HelicityParticle>& p) { |
---|
713 | |
---|
714 | for (int i = 0; i < p[idx].spinStates(); i++ ) { |
---|
715 | for (int j = 1; j < p[idx].spinStates(); j++) { |
---|
716 | if (i == j) p[idx].rho[i][j] = 1.0 / |
---|
717 | static_cast<double>(p[idx].spinStates()); |
---|
718 | else p[idx].rho[i][j] = 0; |
---|
719 | } |
---|
720 | } |
---|
721 | |
---|
722 | } |
---|
723 | |
---|
724 | //========================================================================== |
---|
725 | |
---|
726 | // Base class for all tau decay matrix elements. This class derives from |
---|
727 | // the HelicityMatrixElement class and redefines some of the virtual functions. |
---|
728 | |
---|
729 | // One new method, initHadronicCurrent is defined which initializes the |
---|
730 | // hadronic current in the initWaves method. For each tau decay matrix element |
---|
731 | // the hadronic current method must be redefined accordingly, but initWaves |
---|
732 | // should not be redefined. |
---|
733 | |
---|
734 | //-------------------------------------------------------------------------- |
---|
735 | |
---|
736 | // Initialize wave functions for the helicity matrix element. |
---|
737 | void HMETauDecay::initWaves(vector<HelicityParticle>& p) { |
---|
738 | |
---|
739 | u.clear(); |
---|
740 | pMap.resize(p.size()); |
---|
741 | setFermionLine(0, p[0], p[1]); |
---|
742 | initHadronicCurrent(p); |
---|
743 | |
---|
744 | } |
---|
745 | |
---|
746 | //-------------------------------------------------------------------------- |
---|
747 | |
---|
748 | // Return element for the helicity matrix element. |
---|
749 | complex HMETauDecay::calculateME(vector<int> h) { |
---|
750 | |
---|
751 | complex answer(0,0); |
---|
752 | for (int mu = 0; mu <= 3; mu++) { |
---|
753 | answer += |
---|
754 | (u[1][h[pMap[1]]] * gamma[mu] * (1 - gamma[5]) * u[0][h[pMap[0]]]) |
---|
755 | * gamma[4](mu,mu) * u[2][0](mu); |
---|
756 | } |
---|
757 | return answer; |
---|
758 | |
---|
759 | } |
---|
760 | |
---|
761 | //-------------------------------------------------------------------------- |
---|
762 | |
---|
763 | // Return the maximum decay weight for the helicity matrix element. |
---|
764 | |
---|
765 | double HMETauDecay::decayWeightMax(vector<HelicityParticle>& p) { |
---|
766 | |
---|
767 | // Determine the maximum on-diagonal element of rho. |
---|
768 | double on = real(p[0].rho[0][0]) > real(p[0].rho[1][1]) ? |
---|
769 | real(p[0].rho[0][0]) : real(p[0].rho[1][1]); |
---|
770 | // Determine the maximum off-diagonal element of rho. |
---|
771 | double off = fabs(real(p[0].rho[0][1])) + fabs(imag(p[0].rho[0][1])); |
---|
772 | return DECAYWEIGHTMAX * (on + off); |
---|
773 | |
---|
774 | } |
---|
775 | |
---|
776 | //-------------------------------------------------------------------------- |
---|
777 | |
---|
778 | // Calculate complex resonance weights given a phase and amplitude vector. |
---|
779 | |
---|
780 | void HMETauDecay::calculateResonanceWeights(vector<double>& phase, |
---|
781 | vector<double>& amplitude, vector<complex>& weight) { |
---|
782 | |
---|
783 | for (unsigned int i = 0; i < phase.size(); i++) |
---|
784 | weight.push_back(amplitude[i] * (cos(phase[i]) + |
---|
785 | complex(0,1) * sin(phase[i]))); |
---|
786 | |
---|
787 | } |
---|
788 | |
---|
789 | //========================================================================== |
---|
790 | |
---|
791 | // Tau decay matrix element for tau decay into a single scalar meson. |
---|
792 | |
---|
793 | // The maximum decay weight for this matrix element can be found analytically |
---|
794 | // to be 4 * m_tau^2 * (m_tau^2 - m_meson^2). However, because m_tau >> m_meson |
---|
795 | // for the relevant tau decay channels, this expression is approximated by |
---|
796 | // m_tau^4. |
---|
797 | |
---|
798 | //-------------------------------------------------------------------------- |
---|
799 | |
---|
800 | // Initialize constants for the helicity matrix element. |
---|
801 | |
---|
802 | void HMETau2Meson::initConstants() { |
---|
803 | |
---|
804 | DECAYWEIGHTMAX = 4*pow4(pM[0]); |
---|
805 | |
---|
806 | } |
---|
807 | |
---|
808 | //-------------------------------------------------------------------------- |
---|
809 | |
---|
810 | // Initialize the hadronic current for the helicity matrix element. |
---|
811 | |
---|
812 | void HMETau2Meson::initHadronicCurrent(vector<HelicityParticle>& p) { |
---|
813 | |
---|
814 | vector< Wave4 > u2; |
---|
815 | pMap[2] = 2; |
---|
816 | u2.push_back(Wave4(p[2].p())); |
---|
817 | u.push_back(u2); |
---|
818 | |
---|
819 | } |
---|
820 | |
---|
821 | //========================================================================== |
---|
822 | |
---|
823 | // Tau decay matrix element for tau decay into two leptons. Because there is |
---|
824 | // no hadronic current, but rather a leptonic current, the calculateME and |
---|
825 | // initWaves methods must be redefined. |
---|
826 | |
---|
827 | //-------------------------------------------------------------------------- |
---|
828 | |
---|
829 | // Initialize constants for the helicity matrix element. |
---|
830 | |
---|
831 | void HMETau2TwoLeptons::initConstants() { |
---|
832 | |
---|
833 | DECAYWEIGHTMAX = 16*pow4(pM[0]); |
---|
834 | |
---|
835 | } |
---|
836 | |
---|
837 | //-------------------------------------------------------------------------- |
---|
838 | |
---|
839 | // Initialize spinors for the helicity matrix element. |
---|
840 | |
---|
841 | void HMETau2TwoLeptons::initWaves(vector<HelicityParticle>& p) { |
---|
842 | |
---|
843 | u.clear(); |
---|
844 | pMap.resize(4); |
---|
845 | setFermionLine(0,p[0],p[1]); |
---|
846 | setFermionLine(2,p[2],p[3]); |
---|
847 | |
---|
848 | } |
---|
849 | |
---|
850 | //-------------------------------------------------------------------------- |
---|
851 | |
---|
852 | // Return element for the helicity matrix element. |
---|
853 | |
---|
854 | complex HMETau2TwoLeptons::calculateME(vector<int> h) { |
---|
855 | |
---|
856 | complex answer(0,0); |
---|
857 | for (int mu = 0; mu <= 3; mu++) { |
---|
858 | answer += (u[1][h[pMap[1]]] * gamma[mu] * (1 - gamma[5]) |
---|
859 | * u[0][h[pMap[0]]]) * gamma[4](mu,mu) * (u[3][h[pMap[3]]] |
---|
860 | * gamma[mu] * (1 - gamma[5]) * u[2][h[pMap[2]]]); |
---|
861 | } |
---|
862 | return answer; |
---|
863 | |
---|
864 | } |
---|
865 | |
---|
866 | //========================================================================== |
---|
867 | |
---|
868 | // Tau decay matrix element for tau decay into two mesons through an |
---|
869 | // intermediate vector meson. This matrix element is used for pi^0 + pi^- |
---|
870 | // decays (rho resonances), K^0 + K^- decays (rho resonances), and eta + K^- |
---|
871 | // decays (K^* resonances). Note that for the rho resonances the pi^0 + pi^- |
---|
872 | // running width dominates while for the K^* resonances the pi^- + K^0 running |
---|
873 | // width dominates. |
---|
874 | |
---|
875 | // vecM: on-shell masses for the vector resonances |
---|
876 | // vecG: on-shell widths for the vector resonances |
---|
877 | // vecP: phases used to calculate vector resonance weights |
---|
878 | // vecA: amplitudes used to calculate vector resonance weights |
---|
879 | // vecW: vector resonance weights |
---|
880 | |
---|
881 | //-------------------------------------------------------------------------- |
---|
882 | |
---|
883 | // Initialize constants for the helicity matrix element. |
---|
884 | |
---|
885 | void HMETau2TwoMesonsViaVector::initConstants() { |
---|
886 | |
---|
887 | // Clear the vectors from previous decays. |
---|
888 | vecM.clear(); vecG.clear(); vecP.clear(); vecA.clear(); vecW.clear(); |
---|
889 | |
---|
890 | // Decay through K^* resonances (eta + K^- decay). |
---|
891 | if (abs(pID[2]) == 221) { |
---|
892 | DECAYWEIGHTMAX = 10; |
---|
893 | pM[2] = particleDataPtr->m0(211); pM[3] = particleDataPtr->m0(311); |
---|
894 | vecM.push_back(0.8921); vecM.push_back(1.700); |
---|
895 | vecG.push_back(0.0513); vecG.push_back(0.235); |
---|
896 | vecP.push_back(0); vecP.push_back(M_PI); |
---|
897 | vecA.push_back(1); vecA.push_back(0.038); |
---|
898 | } |
---|
899 | |
---|
900 | // Decay through rho resonances (pi^0 + pi^- and K^0 + K^- decays). |
---|
901 | else { |
---|
902 | if (abs(pID[2]) == 111) DECAYWEIGHTMAX = 800; |
---|
903 | else if (abs(pID[2]) == 311) DECAYWEIGHTMAX = 6; |
---|
904 | pM[2] = particleDataPtr->m0(111); pM[3] = particleDataPtr->m0(211); |
---|
905 | vecM.push_back(0.7746); vecM.push_back(1.4080); vecM.push_back(1.700); |
---|
906 | vecG.push_back(0.1490); vecG.push_back(0.5020); vecG.push_back(0.235); |
---|
907 | vecP.push_back(0); vecP.push_back(M_PI); vecP.push_back(0); |
---|
908 | vecA.push_back(1.0); vecA.push_back(0.167); vecA.push_back(0.050); |
---|
909 | } |
---|
910 | calculateResonanceWeights(vecP, vecA, vecW); |
---|
911 | |
---|
912 | } |
---|
913 | |
---|
914 | //-------------------------------------------------------------------------- |
---|
915 | |
---|
916 | // Initialize the hadronic current for the helicity matrix element. |
---|
917 | |
---|
918 | void HMETau2TwoMesonsViaVector::initHadronicCurrent( |
---|
919 | vector<HelicityParticle>& p) { |
---|
920 | |
---|
921 | vector< Wave4 > u2; |
---|
922 | Wave4 u3(p[3].p() - p[2].p()); |
---|
923 | Wave4 u4(p[2].p() + p[3].p()); |
---|
924 | double s1 = m2(u3, u4); |
---|
925 | double s2 = m2(u4); |
---|
926 | complex sumBW = 0; |
---|
927 | for (unsigned int i = 0; i < vecW.size(); i++) |
---|
928 | sumBW += vecW[i] * pBreitWigner(pM[2], pM[3], s2, vecM[i], vecG[i]); |
---|
929 | u2.push_back((u3 - s1 / s2 * u4) * sumBW); |
---|
930 | u.push_back(u2); |
---|
931 | |
---|
932 | } |
---|
933 | |
---|
934 | //========================================================================== |
---|
935 | |
---|
936 | // Tau decay matrix element for tau decay into two mesons through both |
---|
937 | // intermediate vector and scalar mesons. |
---|
938 | |
---|
939 | // scaC: scalar coupling constant |
---|
940 | // scaM: on-shell masses for the scalar resonances |
---|
941 | // scaG: on-shell widths for the scalar resonances |
---|
942 | // scaP: phases used to calculate scalar resonance weights |
---|
943 | // scaA: amplitudes used to calculate scalar resonance weights |
---|
944 | // scaW: scalar resonance weights |
---|
945 | // vecC: scalar coupling constant |
---|
946 | // vecM: on-shell masses for the vector resonances |
---|
947 | // vecG: on-shell widths for the vector resonances |
---|
948 | // vecP: phases used to calculate vector resonance weights |
---|
949 | // vecA: amplitudes used to calculate vector resonance weights |
---|
950 | // vecW: vector resonance weights |
---|
951 | |
---|
952 | //-------------------------------------------------------------------------- |
---|
953 | |
---|
954 | // Initialize constants for the helicity matrix element. |
---|
955 | |
---|
956 | void HMETau2TwoMesonsViaVectorScalar::initConstants() { |
---|
957 | |
---|
958 | DECAYWEIGHTMAX = 5400; |
---|
959 | // Clear the vectors from previous decays. |
---|
960 | scaM.clear(); scaG.clear(); scaP.clear(); scaA.clear(); scaW.clear(); |
---|
961 | vecM.clear(); vecG.clear(); vecP.clear(); vecA.clear(); vecW.clear(); |
---|
962 | // Scalar resonance parameters. |
---|
963 | scaC = 0.465; |
---|
964 | scaM.push_back(0.878); |
---|
965 | scaG.push_back(0.499); |
---|
966 | scaP.push_back(0); |
---|
967 | scaA.push_back(1); |
---|
968 | calculateResonanceWeights(scaP, scaA, scaW); |
---|
969 | // Vector resonance parameters. |
---|
970 | vecC = 1; |
---|
971 | vecM.push_back(0.89547); vecM.push_back(1.414); |
---|
972 | vecG.push_back(0.04619); vecG.push_back(0.232); |
---|
973 | vecP.push_back(0); vecP.push_back(1.4399); |
---|
974 | vecA.push_back(1); vecA.push_back(0.075); |
---|
975 | calculateResonanceWeights(vecP, vecA, vecW); |
---|
976 | |
---|
977 | } |
---|
978 | |
---|
979 | //-------------------------------------------------------------------------- |
---|
980 | |
---|
981 | // Initialize the hadronic current for the helicity matrix element. |
---|
982 | |
---|
983 | void HMETau2TwoMesonsViaVectorScalar::initHadronicCurrent( |
---|
984 | vector<HelicityParticle>& p) { |
---|
985 | |
---|
986 | vector< Wave4 > u2; |
---|
987 | Wave4 u3(p[3].p() - p[2].p()); |
---|
988 | Wave4 u4(p[2].p() + p[3].p()); |
---|
989 | double s1 = m2(u3,u4); |
---|
990 | double s2 = m2(u4); |
---|
991 | complex scaSumBW = 0; complex scaSumW = 0; |
---|
992 | complex vecSumBW = 0; complex vecSumW = 0; complex vecSumBWM = 0; |
---|
993 | for (unsigned int i = 0; i < scaW.size(); i++) { |
---|
994 | scaSumBW += scaW[i] * sBreitWigner(pM[2], pM[3], s2, scaM[i], scaG[i]); |
---|
995 | scaSumW += scaW[i]; |
---|
996 | } |
---|
997 | for (unsigned int i = 0; i < vecW.size(); i++) { |
---|
998 | vecSumBW += vecW[i] * pBreitWigner(pM[2], pM[3], s2, vecM[i], vecG[i]); |
---|
999 | vecSumBWM += vecW[i] * pBreitWigner(pM[2], pM[3], s2, vecM[i], vecG[i]) / |
---|
1000 | pow2(vecM[i]); |
---|
1001 | vecSumW += vecW[i]; |
---|
1002 | } |
---|
1003 | u2.push_back(vecC * (vecSumBW * u3 - s1 * vecSumBWM * u4) / vecSumW + |
---|
1004 | scaC * u4 * scaSumBW / scaSumW); |
---|
1005 | u.push_back(u2); |
---|
1006 | |
---|
1007 | } |
---|
1008 | |
---|
1009 | //========================================================================== |
---|
1010 | |
---|
1011 | // Tau decay matrix element for tau decay into three mesons. This matrix |
---|
1012 | // element provides a base class for all implemented three meson decays. |
---|
1013 | |
---|
1014 | // mode: three meson decay mode of the tau |
---|
1015 | // initMode(): initialize the decay mode |
---|
1016 | // initResonances(): initialize the resonance constants |
---|
1017 | // s1, s2, s3, s4: center-of-mass energies |
---|
1018 | // q, q2, q3, q4: summed and individual hadronic momentum four-vectors |
---|
1019 | // a1BW: stored value of a1BreitWigner for speed |
---|
1020 | // a1PhaseSpace(s): phase space factor for the a1 |
---|
1021 | // a1BreitWigner(s): Breit-Wigner for the a1 |
---|
1022 | // T(m0, m1, s, M, G, W): sum weighted p-wave Breit-Wigners |
---|
1023 | // T(s, M, G, W): sum weighted fixed width Breit-Wigners |
---|
1024 | // F1(), F2(), F3(), F4(): sub-current form factors |
---|
1025 | |
---|
1026 | //-------------------------------------------------------------------------- |
---|
1027 | |
---|
1028 | // Initialize constants for the helicity matrix element. |
---|
1029 | |
---|
1030 | void HMETau2ThreeMesons::initConstants() { |
---|
1031 | |
---|
1032 | initMode(); |
---|
1033 | initResonances(); |
---|
1034 | |
---|
1035 | } |
---|
1036 | |
---|
1037 | //-------------------------------------------------------------------------- |
---|
1038 | |
---|
1039 | // Initialize the hadronic current for the helicity matrix element. |
---|
1040 | |
---|
1041 | void HMETau2ThreeMesons::initHadronicCurrent(vector<HelicityParticle>& p) { |
---|
1042 | |
---|
1043 | vector< Wave4 > u2; |
---|
1044 | |
---|
1045 | // Initialize the momenta. |
---|
1046 | initMomenta(p); |
---|
1047 | |
---|
1048 | // Calculate the center of mass energies. |
---|
1049 | s1 = m2(q); |
---|
1050 | s2 = m2(q3 + q4); |
---|
1051 | s3 = m2(q2 + q4); |
---|
1052 | s4 = m2(q2 + q3); |
---|
1053 | |
---|
1054 | // Calculate the form factors. |
---|
1055 | a1BW = a1BreitWigner(s1); |
---|
1056 | complex f1 = F1(); |
---|
1057 | complex f2 = F2(); |
---|
1058 | complex f3 = F3(); |
---|
1059 | complex f4 = F4(); |
---|
1060 | |
---|
1061 | // Calculate the hadronic current. |
---|
1062 | Wave4 u3 = (f3 - f2) * q2 + (f1 - f3) * q3 + (f2 - f1) * q4; |
---|
1063 | u3 = u3 - (u3 * gamma[4] * q / s1) * q; |
---|
1064 | if (f4 != complex(0, 0)) |
---|
1065 | u3 = u3 + complex(0, 1) * f4 * epsilon(q2, q3, q4); |
---|
1066 | u2.push_back(u3); |
---|
1067 | u.push_back(u2); |
---|
1068 | |
---|
1069 | } |
---|
1070 | |
---|
1071 | //-------------------------------------------------------------------------- |
---|
1072 | |
---|
1073 | // Initialize the tau decay mode. |
---|
1074 | |
---|
1075 | void HMETau2ThreeMesons::initMode() { |
---|
1076 | |
---|
1077 | if (abs(pID[2]) == 111 && abs(pID[3]) == 111 && abs(pID[4]) == 211) |
---|
1078 | mode = Pi0Pi0Pim; |
---|
1079 | else if (abs(pID[2]) == 211 && abs(pID[3]) == 211 && abs(pID[4]) == 211) |
---|
1080 | mode = PimPimPip; |
---|
1081 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 211 && abs(pID[4]) == 311) |
---|
1082 | mode = Pi0PimK0b; |
---|
1083 | else if (abs(pID[2]) == 211 && abs(pID[3]) == 211 && abs(pID[4]) == 321) |
---|
1084 | mode = PimPipKm; |
---|
1085 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 211 && abs(pID[4]) == 221) |
---|
1086 | mode = Pi0PimEta; |
---|
1087 | else if (abs(pID[2]) == 211 && abs(pID[3]) == 321 && abs(pID[4]) == 321) |
---|
1088 | mode = PimKmKp; |
---|
1089 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 311 && abs(pID[4]) == 321) |
---|
1090 | mode = Pi0K0Km; |
---|
1091 | else if (abs(pID[2]) == 130 && abs(pID[3]) == 211 && abs(pID[4]) == 310) |
---|
1092 | mode = KlPimKs; |
---|
1093 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 111 && abs(pID[4]) == 321) |
---|
1094 | mode = Pi0Pi0Km; |
---|
1095 | else if (abs(pID[2]) == 130 && abs(pID[3]) == 130 && abs(pID[4]) == 211) |
---|
1096 | mode = KlKlPim; |
---|
1097 | else if (abs(pID[2]) == 211 && abs(pID[3]) == 310 && abs(pID[4]) == 310) |
---|
1098 | mode = PimKsKs; |
---|
1099 | else if (abs(pID[2]) == 211 && abs(pID[3]) == 311 && abs(pID[4]) == 311) |
---|
1100 | mode = PimK0bK0; |
---|
1101 | else |
---|
1102 | mode = Uknown; |
---|
1103 | } |
---|
1104 | |
---|
1105 | //-------------------------------------------------------------------------- |
---|
1106 | |
---|
1107 | // Initialize the momenta for the helicity matrix element. |
---|
1108 | |
---|
1109 | void HMETau2ThreeMesons::initMomenta(vector<HelicityParticle>& p) { |
---|
1110 | |
---|
1111 | q = Wave4(p[2].p() + p[3].p() + p[4].p()); |
---|
1112 | // pi-, pi-, pi+ decay and pi0, pi0, pi- decay. |
---|
1113 | if (mode == PimPimPip || mode == Pi0Pi0Pim) { |
---|
1114 | q2 = Wave4(p[2].p()); q3 = Wave4(p[3].p()); q4 = Wave4(p[4].p()); |
---|
1115 | // K-, pi-, K+ decay. |
---|
1116 | } else if (mode == PimKmKp) { |
---|
1117 | q2 = Wave4(p[3].p()); q3 = Wave4(p[2].p()); q4 = Wave4(p[4].p()); |
---|
1118 | // K0, pi-, Kbar0 decay. |
---|
1119 | } else if (mode == PimK0bK0) { |
---|
1120 | q2 = Wave4(p[3].p()); q3 = Wave4(p[2].p()); q4 = Wave4(p[4].p()); |
---|
1121 | // K_S0, pi-, K_S0 decay. |
---|
1122 | } else if (mode == PimKsKs) { |
---|
1123 | q2 = Wave4(p[3].p()); q3 = Wave4(p[2].p()); q4 = Wave4(p[4].p()); |
---|
1124 | // K_L0, pi-, K_L0 decay. |
---|
1125 | } else if (mode == KlKlPim) { |
---|
1126 | q2 = Wave4(p[2].p()); q3 = Wave4(p[4].p()); q4 = Wave4(p[3].p()); |
---|
1127 | // K_S0, pi-, K_L0 decay. |
---|
1128 | } else if (mode == KlPimKs) { |
---|
1129 | q2 = Wave4(p[4].p()); q3 = Wave4(p[3].p()); q4 = Wave4(p[2].p()); |
---|
1130 | } // K-, pi0, K0 decay. |
---|
1131 | else if (mode == Pi0K0Km) { |
---|
1132 | q2 = Wave4(p[4].p()); q3 = Wave4(p[2].p()); q4 = Wave4(p[3].p()); |
---|
1133 | } // pi0, pi0, K- decay. |
---|
1134 | else if (mode == Pi0Pi0Km) { |
---|
1135 | q2 = Wave4(p[2].p()); q3 = Wave4(p[3].p()); q4 = Wave4(p[4].p()); |
---|
1136 | } // K-, pi-, pi+ decay. |
---|
1137 | else if (mode == PimPipKm) { |
---|
1138 | q2 = Wave4(p[4].p()); q3 = Wave4(p[2].p()); q4 = Wave4(p[3].p()); |
---|
1139 | } // pi-, Kbar0, pi0 decay. |
---|
1140 | else if (mode == Pi0PimK0b) { |
---|
1141 | q2 = Wave4(p[3].p()); q3 = Wave4(p[4].p()); q4 = Wave4(p[2].p()); |
---|
1142 | } // pi-, pi0, eta decay. |
---|
1143 | else if (mode == Pi0PimEta) { |
---|
1144 | q2 = Wave4(p[3].p()); q3 = Wave4(p[2].p()); q4 = Wave4(p[4].p()); |
---|
1145 | } |
---|
1146 | |
---|
1147 | } |
---|
1148 | |
---|
1149 | //-------------------------------------------------------------------------- |
---|
1150 | |
---|
1151 | // Return the phase space factor for the a1. Implements equation 3.16 of Z. |
---|
1152 | // Phys. C48 (1990) 445-452. |
---|
1153 | |
---|
1154 | double HMETau2ThreeMesons::a1PhaseSpace(double s) { |
---|
1155 | |
---|
1156 | double piM = 0.13957; // Mass of the charged pion. |
---|
1157 | double rhoM = 0.773; // Mass of the rho. |
---|
1158 | if (s < pow2(3 * piM)) |
---|
1159 | return 0; |
---|
1160 | else if (s < pow2(rhoM + piM)) { |
---|
1161 | double sum = (s - 9 * piM * piM); |
---|
1162 | return 4.1 * sum * sum * sum * (1 - 3.3 * sum + 5.8 * sum * sum); |
---|
1163 | } |
---|
1164 | else |
---|
1165 | return s * (1.623 + 10.38 / s - 9.32 / (s * s) + 0.65 / (s * s * s)); |
---|
1166 | |
---|
1167 | } |
---|
1168 | |
---|
1169 | //-------------------------------------------------------------------------- |
---|
1170 | |
---|
1171 | // Return the Breit-Wigner for the a1. Implements equation 3.18 of Z. Phys. C48 |
---|
1172 | // (1990) 445-452. |
---|
1173 | |
---|
1174 | complex HMETau2ThreeMesons::a1BreitWigner(double s) { |
---|
1175 | |
---|
1176 | |
---|
1177 | double a1M = 1.251; // Mass of the a1. |
---|
1178 | double a1G = 0.475; // Width of the a1. |
---|
1179 | return a1M * a1M / (a1M * a1M - s - complex(0,1) * a1M * a1G |
---|
1180 | * a1PhaseSpace(s) / a1PhaseSpace(a1M * a1M)); |
---|
1181 | |
---|
1182 | } |
---|
1183 | |
---|
1184 | //-------------------------------------------------------------------------- |
---|
1185 | |
---|
1186 | // Return summed weighted running p Breit-Wigner resonances. |
---|
1187 | |
---|
1188 | complex HMETau2ThreeMesons::T(double m0, double m1, double s, vector<double> &M, |
---|
1189 | vector<double> &G, vector<double> &W) { |
---|
1190 | |
---|
1191 | complex num(0, 0); |
---|
1192 | double den(0); |
---|
1193 | for (unsigned int i = 0; i < M.size(); i++) { |
---|
1194 | num += W[i] * pBreitWigner(m0, m1, s, M[i], G[i]); |
---|
1195 | den += W[i]; |
---|
1196 | } |
---|
1197 | return num / den; |
---|
1198 | |
---|
1199 | } |
---|
1200 | |
---|
1201 | //-------------------------------------------------------------------------- |
---|
1202 | |
---|
1203 | // Return summed weighted fixed width Breit-Wigner resonances. |
---|
1204 | |
---|
1205 | complex HMETau2ThreeMesons::T(double s, vector<double> &M, |
---|
1206 | vector<double> &G, vector<double> &W) { |
---|
1207 | |
---|
1208 | complex num(0, 0); |
---|
1209 | double den(0); |
---|
1210 | for (unsigned int i = 0; i < M.size(); i++) { |
---|
1211 | num += W[i] * breitWigner(s, M[i], G[i]); |
---|
1212 | den += W[i]; |
---|
1213 | } |
---|
1214 | return num / den; |
---|
1215 | |
---|
1216 | } |
---|
1217 | |
---|
1218 | //========================================================================== |
---|
1219 | |
---|
1220 | // Tau decay matrix element for tau decay into three pions. This matrix element |
---|
1221 | // is taken from the Herwig++ implementation based on the CLEO fits. |
---|
1222 | |
---|
1223 | // rhoM: on-shell masses for the rho resonances |
---|
1224 | // rhoG: on-shell widths for the rho resonances |
---|
1225 | // rhoPp: p-wave phase for the rho coupling to the a1 |
---|
1226 | // rhoAp: p-wave amplitude for the rho coupling to the a1 |
---|
1227 | // rhoPd: d-wave phase for the rho coupling to the a1 |
---|
1228 | // rhoAd: d-wave amplitude for the rho coupling to the a1 |
---|
1229 | // f0M: f0 on-shell mass |
---|
1230 | // f0G: f0 on-shell width |
---|
1231 | // f0P: phase for the coupling of the f0 to the a1 |
---|
1232 | // f0A: amplitude for the coupling of the f0 to the a1 |
---|
1233 | // f2M: f2 on-shell mass |
---|
1234 | // f2G: f2 on-shell width |
---|
1235 | // f2P: phase for the coupling of the f2 to the a1 |
---|
1236 | // f2P: amplitude for the coupling of the f2 to the a1 |
---|
1237 | // sigM: sigma on-shell mass |
---|
1238 | // sigG: sigma on-shell width |
---|
1239 | // sigP: phase for the coupling of the sigma to the a1 |
---|
1240 | // sigA: amplitude for the coupling of the sigma to the a1 |
---|
1241 | |
---|
1242 | //-------------------------------------------------------------------------- |
---|
1243 | |
---|
1244 | // Initialize resonance constants for the helicity matrix element. |
---|
1245 | |
---|
1246 | void HMETau2ThreePions::initResonances() { |
---|
1247 | |
---|
1248 | // Three charged pion decay. |
---|
1249 | if (mode == PimPimPip) DECAYWEIGHTMAX = 6000; |
---|
1250 | |
---|
1251 | // Two neutral and one charged pion decay. |
---|
1252 | else DECAYWEIGHTMAX = 3000; |
---|
1253 | |
---|
1254 | // Clear the vectors from previous decays. |
---|
1255 | rhoM.clear(); rhoG.clear(); |
---|
1256 | rhoPp.clear(); rhoAp.clear(); rhoWp.clear(); |
---|
1257 | rhoPd.clear(); rhoAd.clear(); rhoWd.clear(); |
---|
1258 | |
---|
1259 | // Rho parameters. |
---|
1260 | rhoM.push_back(.7743); rhoM.push_back(1.370); rhoM.push_back(1.720); |
---|
1261 | rhoG.push_back(.1491); rhoG.push_back(.386); rhoG.push_back(.250); |
---|
1262 | rhoPp.push_back(0); rhoPp.push_back(3.11018); rhoPp.push_back(0); |
---|
1263 | rhoAp.push_back(1); rhoAp.push_back(0.12); rhoAp.push_back(0); |
---|
1264 | rhoPd.push_back(-0.471239); rhoPd.push_back(1.66504); rhoPd.push_back(0); |
---|
1265 | rhoAd.push_back(3.7e-07); rhoAd.push_back(8.7e-07); rhoAd.push_back(0); |
---|
1266 | |
---|
1267 | // Scalar and tensor parameters. |
---|
1268 | f0M = 1.186; f2M = 1.275; sigM = 0.860; |
---|
1269 | f0G = 0.350; f2G = 0.185; sigG = 0.880; |
---|
1270 | f0P = -1.69646; f2P = 1.75929; sigP = 0.722566; |
---|
1271 | f0A = 0.77; f2A = 7.1e-07; sigA = 2.1; |
---|
1272 | |
---|
1273 | // Calculate the weights from the phases and amplitudes. |
---|
1274 | calculateResonanceWeights(rhoPp, rhoAp, rhoWp); |
---|
1275 | calculateResonanceWeights(rhoPd, rhoAd, rhoWd); |
---|
1276 | f0W = f0A * (cos(f0P) + complex(0,1) * sin(f0P)); |
---|
1277 | f2W = f2A * (cos(f2P) + complex(0,1) * sin(f2P)); |
---|
1278 | sigW = sigA * (cos(sigP) + complex(0,1) * sin(sigP)); |
---|
1279 | |
---|
1280 | } |
---|
1281 | |
---|
1282 | //-------------------------------------------------------------------------- |
---|
1283 | |
---|
1284 | // Return the first form factor. |
---|
1285 | |
---|
1286 | complex HMETau2ThreePions::F1() { |
---|
1287 | |
---|
1288 | complex answer(0,0); |
---|
1289 | |
---|
1290 | // Three charged pion decay. |
---|
1291 | if (mode == PimPimPip) { |
---|
1292 | for (unsigned int i = 0; i < rhoM.size(); i++) { |
---|
1293 | answer += - rhoWp[i] * pBreitWigner(pM[3], pM[4], s2, rhoM[i], rhoG[i]) |
---|
1294 | - rhoWd[i] / 3.0 * pBreitWigner(pM[2], pM[4], s3, rhoM[i], rhoG[i]) |
---|
1295 | * (s2 - s4); |
---|
1296 | } |
---|
1297 | answer += -2.0 / 3.0 * (sigW * sBreitWigner(pM[2], pM[4], s3, sigM, sigG) |
---|
1298 | + f0W * sBreitWigner(pM[2], pM[4], s3, f0M, f0G)); |
---|
1299 | answer += f2W * (0.5 * (s4 - s3) |
---|
1300 | * dBreitWigner(pM[3], pM[4], s2, f2M, f2G) |
---|
1301 | - 1.0 / (18 * s3) * (4 * pow2(pM[2]) - s3) |
---|
1302 | * (s1 + s3 - pow2(pM[2])) |
---|
1303 | * dBreitWigner(pM[2], pM[4], s3, f2M, f2G)); |
---|
1304 | } |
---|
1305 | |
---|
1306 | // Two neutral and one charged pion decay. |
---|
1307 | else { |
---|
1308 | for (unsigned int i = 0; i < rhoM.size(); i++) { |
---|
1309 | answer += rhoWp[i] * pBreitWigner(pM[3], pM[4], s2, rhoM[i], rhoG[i]) |
---|
1310 | - rhoWd[i] / 3.0 * pBreitWigner(pM[2], pM[4], s3, rhoM[i], rhoG[i]) |
---|
1311 | * (s4 - s2 - pow2(pM[4]) + pow2(pM[2])); |
---|
1312 | } |
---|
1313 | answer += 2.0 / 3.0 * (sigW * sBreitWigner(pM[2], pM[3], s4, sigM, sigG) |
---|
1314 | + f0W * sBreitWigner(pM[2], pM[3], s4, f0M, f0G)); |
---|
1315 | answer += f2W / (18 * s4) * (s1 - pow2(pM[4]) + s4) |
---|
1316 | * (4 * pow2(pM[2]) - s4) * dBreitWigner(pM[2], pM[3], s4, f2M, f2G); |
---|
1317 | } |
---|
1318 | return a1BW * answer; |
---|
1319 | |
---|
1320 | } |
---|
1321 | |
---|
1322 | //-------------------------------------------------------------------------- |
---|
1323 | |
---|
1324 | // Return the second form factor. |
---|
1325 | |
---|
1326 | complex HMETau2ThreePions::F2() { |
---|
1327 | |
---|
1328 | complex answer(0,0); |
---|
1329 | |
---|
1330 | // Three charged pion decay. |
---|
1331 | if (mode == PimPimPip) { |
---|
1332 | for (unsigned int i = 0; i < rhoM.size(); i++) { |
---|
1333 | answer += -rhoWp[i] * pBreitWigner(pM[2], pM[4], s3, rhoM[i], rhoG[i]) |
---|
1334 | - rhoWd[i] / 3.0 * pBreitWigner(pM[3], pM[4], s2, rhoM[i], rhoG[i]) |
---|
1335 | * (s3 - s4); |
---|
1336 | } |
---|
1337 | answer += -2.0 / 3.0 * (sigW * sBreitWigner(pM[3], pM[4], s2, sigM, sigG) |
---|
1338 | + f0W * sBreitWigner(pM[3], pM[4], s2, f0M, f0G)); |
---|
1339 | answer += f2W * (0.5 * (s4 - s2) |
---|
1340 | * dBreitWigner(pM[2], pM[4], s3, f2M, f2G) |
---|
1341 | - 1.0 / (18 * s2) * (4 * pow2(pM[2]) - s2) * (s1 + s2 - pow2(pM[2])) |
---|
1342 | * dBreitWigner(pM[3], pM[4], s2, f2M, f2G)); |
---|
1343 | } |
---|
1344 | |
---|
1345 | // Two neutral and one charged pion decay. |
---|
1346 | else { |
---|
1347 | for (unsigned int i = 0; i < rhoM.size(); i++) { |
---|
1348 | answer += -rhoWp[i] / 3.0 * |
---|
1349 | pBreitWigner(pM[2], pM[4], s3, rhoM[i], rhoG[i]) - |
---|
1350 | rhoWd[i] * pBreitWigner(pM[3], pM[4], s2, rhoM[i], rhoG[i]) * |
---|
1351 | (s4 - s3 - pow2(pM[4]) + pow2(pM[3])); |
---|
1352 | } |
---|
1353 | answer += 2.0 / 3.0 * (sigW * sBreitWigner(pM[2], pM[3], s4, sigM, sigG) |
---|
1354 | + f0W * sBreitWigner(pM[2], pM[3], s4, f0M, f0G)); |
---|
1355 | answer += f2W / (18 * s4) * (s1 - pow2(pM[4]) + s4) * |
---|
1356 | (4 * pow2(pM[2]) - s4) * dBreitWigner(pM[2], pM[3], s4, f2M, f2G); |
---|
1357 | } |
---|
1358 | return -a1BW * answer; |
---|
1359 | |
---|
1360 | } |
---|
1361 | |
---|
1362 | //-------------------------------------------------------------------------- |
---|
1363 | |
---|
1364 | // Return the third form factor. |
---|
1365 | |
---|
1366 | complex HMETau2ThreePions::F3() { |
---|
1367 | |
---|
1368 | complex answer(0,0); |
---|
1369 | |
---|
1370 | // Three charged pion decay. |
---|
1371 | if (mode == PimPimPip) { |
---|
1372 | for (unsigned int i = 0; i < rhoM.size(); i++) { |
---|
1373 | answer += -rhoWd[i] * (1.0 / 3.0 * (s3 - s4) * |
---|
1374 | pBreitWigner(pM[3], pM[4], s2, rhoM[i], rhoG[i]) |
---|
1375 | - 1.0 / 3.0 * (s2 - s4) * |
---|
1376 | pBreitWigner(pM[2], pM[4], s3, rhoM[i], |
---|
1377 | rhoG[i])); |
---|
1378 | } |
---|
1379 | answer += -2.0 / 3.0 * (sigW * sBreitWigner(pM[3], pM[4], s2, sigM, sigG) |
---|
1380 | + f0W * sBreitWigner(pM[3], pM[4], s2, f0M, f0G)); |
---|
1381 | answer += 2.0 / 3.0 * (sigW * sBreitWigner(pM[2], pM[4], s3, sigM, sigG) |
---|
1382 | + f0W * sBreitWigner(pM[2], pM[4], s3, f0M, f0G)); |
---|
1383 | answer += f2W * (-1.0 / (18 * s2) * (4 * pow2(pM[2]) - s2) * |
---|
1384 | (s1 + s2 - pow2(pM[2])) * |
---|
1385 | dBreitWigner(pM[3], pM[4], s2, f2M, f2G) + |
---|
1386 | 1.0 / (18 * s3) * (4 * pow2(pM[2]) - s3) * |
---|
1387 | (s1 + s3 - pow2(pM[2])) * |
---|
1388 | dBreitWigner(pM[2], pM[4], s3, f2M, f2G)); |
---|
1389 | } |
---|
1390 | |
---|
1391 | // Two neutral and one charged pion decay. |
---|
1392 | else { |
---|
1393 | for (unsigned int i = 0; i < rhoM.size(); i++) { |
---|
1394 | answer += rhoWd[i] * (-1.0 / 3.0 * |
---|
1395 | (s4 - s3 - pow2(pM[4]) + pow2(pM[3])) * |
---|
1396 | pBreitWigner(pM[3], pM[4], s2, rhoM[i], rhoG[i]) + |
---|
1397 | 1.0 / 3.0 * (s4 - s2 - pow2(pM[4]) + pow2(pM[2])) |
---|
1398 | * pBreitWigner(pM[2], pM[4], s3, rhoM[i], |
---|
1399 | rhoG[i])); |
---|
1400 | } |
---|
1401 | answer += -f2W / 2.0 * (s2 - s3) * |
---|
1402 | dBreitWigner(pM[2], pM[3], s4, f2M, f2G); |
---|
1403 | } |
---|
1404 | return a1BW * answer; |
---|
1405 | |
---|
1406 | } |
---|
1407 | |
---|
1408 | //-------------------------------------------------------------------------- |
---|
1409 | |
---|
1410 | // Return the running width for the a1 (multiplied by a factor of a1M). |
---|
1411 | |
---|
1412 | double HMETau2ThreePions::a1PhaseSpace(double s) { |
---|
1413 | |
---|
1414 | double picM = 0.1753; // (m_pi^- + m_pi^- + m_pi^+)^2 |
---|
1415 | double pinM = 0.1676; // (m_pi^0 + m_pi^0 + m_pi^-)^2 |
---|
1416 | double kM = 0.496; // K mass. |
---|
1417 | double ksM = 0.894; // K^* mass. |
---|
1418 | double picG = 0; // Width contribution from three charged pions. |
---|
1419 | double pinG = 0; // Width contribution from two neutral one charged. |
---|
1420 | double kG = 0; // Width contributions from s-wave K K^*. |
---|
1421 | double piW = pow2(0.2384)/1.0252088; // Overall weight factor. |
---|
1422 | double kW = pow2(4.7621); // K K^* width weight factor. |
---|
1423 | |
---|
1424 | // Three charged pion width contribution. |
---|
1425 | if (s < picM) |
---|
1426 | picG = 0; |
---|
1427 | else if (s < 0.823) |
---|
1428 | picG = 5.80900 * pow3(s - picM) * (1.0 - 3.00980 * (s - picM) + |
---|
1429 | 4.5792 * pow2(s - picM)); |
---|
1430 | else |
---|
1431 | picG = -13.91400 + 27.67900 * s - 13.39300 * pow2(s) + 3.19240 * pow3(s) |
---|
1432 | - 0.10487 * pow4(s); |
---|
1433 | |
---|
1434 | // Two neutral and one charged pion width contribution. |
---|
1435 | if (s < pinM) |
---|
1436 | pinG = 0; |
---|
1437 | else if (s < 0.823) |
---|
1438 | pinG = 6.28450 * pow3(s - pinM) * (1.0 - 2.95950 * (s - pinM) + |
---|
1439 | 4.33550 * pow2(s - pinM)); |
---|
1440 | else |
---|
1441 | pinG = -15.41100 + 32.08800 * s - 17.66600 * pow2(s) + 4.93550 * pow3(s) |
---|
1442 | - 0.37498 * pow4(s); |
---|
1443 | |
---|
1444 | // K and K^* width contribution. |
---|
1445 | if (s > pow2(ksM + kM)) |
---|
1446 | kG = 0.5 * sqrt((s - pow2(ksM + kM)) * (s - pow2(ksM - kM))) / s; |
---|
1447 | return piW*(picG + pinG + kW*kG); |
---|
1448 | |
---|
1449 | } |
---|
1450 | |
---|
1451 | //-------------------------------------------------------------------------- |
---|
1452 | |
---|
1453 | // Return the Breit-Wigner for the a1. |
---|
1454 | |
---|
1455 | complex HMETau2ThreePions::a1BreitWigner(double s) { |
---|
1456 | |
---|
1457 | double a1M = 1.331; // Mass of the a1. |
---|
1458 | return a1M*a1M/(a1M*a1M - s - complex(0,1)*a1PhaseSpace(s)); |
---|
1459 | |
---|
1460 | } |
---|
1461 | |
---|
1462 | //========================================================================== |
---|
1463 | |
---|
1464 | // Tau decay matrix element for tau decay into three mesons with kaons. The form |
---|
1465 | // factors are taken from hep-ph/9503474. |
---|
1466 | |
---|
1467 | // rhoMa(v): on-shell masses for the axial (vector) rho resonances |
---|
1468 | // rhoGa(v): widths for the axial (vector) rho resonances |
---|
1469 | // rhoWa(v): weights for the axial (vector) rho resonances |
---|
1470 | // kstarXa(v): on-shell masses, widths, and weights for the K* resonances |
---|
1471 | // k1Xa(b): on-shell masses, width, and weight for the K1 resonances, |
---|
1472 | // the a constants are for K1 -> K* pi, K* -> K pi |
---|
1473 | // the b constants are for K1 -> rho K, rho -> pi pi |
---|
1474 | // omegaX: on-shell masses, width, and weights for the omega reosnances |
---|
1475 | // kM: kaon mass |
---|
1476 | // piM: charged pion mass |
---|
1477 | // piW: pion coupling, f_W |
---|
1478 | |
---|
1479 | //-------------------------------------------------------------------------- |
---|
1480 | |
---|
1481 | // Initialize resonance constants for the helicity matrix element. |
---|
1482 | |
---|
1483 | void HMETau2ThreeMesonsWithKaons::initResonances() { |
---|
1484 | |
---|
1485 | // K-, pi-, K+ decay. |
---|
1486 | if (mode == PimKmKp) DECAYWEIGHTMAX = 130; |
---|
1487 | // K0, pi-, Kbar0 decay. |
---|
1488 | else if (mode == PimK0bK0) DECAYWEIGHTMAX = 115; |
---|
1489 | // K_S0, pi-, K_S0 and K_L0, pi-, K_L0 decay. |
---|
1490 | else if (mode == PimKsKs || mode == KlKlPim) DECAYWEIGHTMAX = 230; |
---|
1491 | // K_S0, pi-, K_L0 decay. |
---|
1492 | else if (mode == KlPimKs) DECAYWEIGHTMAX = 230; |
---|
1493 | // K-, pi0, K0 decay. |
---|
1494 | else if (mode == Pi0K0Km) DECAYWEIGHTMAX = 125; |
---|
1495 | // pi0, pi0, K- decay. |
---|
1496 | else if (mode == Pi0Pi0Km) DECAYWEIGHTMAX = 2.5e4; |
---|
1497 | // K-, pi-, pi+ decay. |
---|
1498 | else if (mode == PimPipKm) DECAYWEIGHTMAX = 1.8e4; |
---|
1499 | // pi-, Kbar0, pi0 decay. |
---|
1500 | else if (mode == Pi0PimK0b) DECAYWEIGHTMAX = 3.9e4; |
---|
1501 | |
---|
1502 | // Clear the vectors from previous decays. |
---|
1503 | rhoMa.clear(); rhoGa.clear(); rhoWa.clear(); |
---|
1504 | rhoMv.clear(); rhoGv.clear(); rhoWv.clear(); |
---|
1505 | kstarMa.clear(); kstarGa.clear(); kstarWa.clear(); |
---|
1506 | kstarMv.clear(); kstarGv.clear(); kstarWv.clear(); |
---|
1507 | k1Ma.clear(); k1Ga.clear(); k1Wa.clear(); |
---|
1508 | k1Mb.clear(); k1Gb.clear(); k1Wb.clear(); |
---|
1509 | omegaM.clear(); omegaG.clear(); omegaW.clear(); |
---|
1510 | |
---|
1511 | // Rho parameters. |
---|
1512 | rhoMa.push_back(0.773); rhoGa.push_back(0.145); rhoWa.push_back(1); |
---|
1513 | rhoMa.push_back(1.370); rhoGa.push_back(0.510); rhoWa.push_back(-0.145); |
---|
1514 | rhoMv.push_back(0.773); rhoGv.push_back(0.145); rhoWv.push_back(1); |
---|
1515 | rhoMv.push_back(1.500); rhoGv.push_back(0.220); rhoWv.push_back(-6.5 / 26.0); |
---|
1516 | rhoMv.push_back(1.750); rhoGv.push_back(0.120); rhoWv.push_back(-1.0 / 26.0); |
---|
1517 | |
---|
1518 | // Kstar parameters. |
---|
1519 | kstarMa.push_back(0.892); kstarGa.push_back(0.050); |
---|
1520 | kstarMa.push_back(1.412); kstarGa.push_back(0.227); |
---|
1521 | kstarWa.push_back(1); |
---|
1522 | kstarWa.push_back(-0.135); |
---|
1523 | kstarMv.push_back(0.892); kstarGv.push_back(0.050); |
---|
1524 | kstarMv.push_back(1.412); kstarGv.push_back(0.227); |
---|
1525 | kstarMv.push_back(1.714); kstarGv.push_back(0.323); |
---|
1526 | kstarWv.push_back(1); |
---|
1527 | kstarWv.push_back(-6.5 / 26.0); |
---|
1528 | kstarWv.push_back(-1.0 / 26.0); |
---|
1529 | |
---|
1530 | // K1 parameters. |
---|
1531 | k1Ma.push_back(1.270); k1Ga.push_back(0.090); k1Wa.push_back(0.33); |
---|
1532 | k1Ma.push_back(1.402); k1Ga.push_back(0.174); k1Wa.push_back(1); |
---|
1533 | k1Mb.push_back(1.270); k1Gb.push_back(0.090); k1Wb.push_back(1); |
---|
1534 | |
---|
1535 | // Omega and phi parameters. |
---|
1536 | omegaM.push_back(0.782); omegaG.push_back(0.00843); omegaW.push_back(1); |
---|
1537 | omegaM.push_back(1.020); omegaG.push_back(0.00443); omegaW.push_back(0.05); |
---|
1538 | |
---|
1539 | // Kaon and pion parameters |
---|
1540 | kM = 0.49765; piM = 0.13957; piW = 0.0942; |
---|
1541 | |
---|
1542 | } |
---|
1543 | |
---|
1544 | //-------------------------------------------------------------------------- |
---|
1545 | |
---|
1546 | // Return the first form factor. |
---|
1547 | |
---|
1548 | complex HMETau2ThreeMesonsWithKaons::F1() { |
---|
1549 | |
---|
1550 | complex answer; |
---|
1551 | // K-, pi-, K+ decay. |
---|
1552 | if (mode == PimKmKp) |
---|
1553 | answer = a1BW * T(piM, kM, s2, kstarMa, kstarGa, kstarWa) / 2.0; |
---|
1554 | // K0, pi-, Kbar0 decay. |
---|
1555 | else if (mode == PimK0bK0) |
---|
1556 | answer = a1BW * T(piM, kM, s2, kstarMa, kstarGa, kstarWa) / 2.0; |
---|
1557 | // K_S0, pi-, K_S0 decay and K_L0, pi-, K_L0 decay. |
---|
1558 | else if (mode == PimKsKs || mode == KlKlPim) |
---|
1559 | answer = -a1BW * (T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1560 | + T(piM, kM, s4, kstarMa, kstarGa, kstarWa)) / 2.0; |
---|
1561 | // K_S0, pi-, K_L0 decay. |
---|
1562 | else if (mode == KlPimKs) |
---|
1563 | answer = a1BW * (T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1564 | - T(piM, kM, s4, kstarMa, kstarGa, kstarWa)) / 2.0; |
---|
1565 | // K-, pi0, K0 decay. |
---|
1566 | else if (mode == Pi0K0Km) |
---|
1567 | answer = a1BW * (T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1568 | - T(piM, kM, s4, kstarMa, kstarGa, kstarWa)) / 2.0; |
---|
1569 | // pi0, pi0, K- decay. |
---|
1570 | else if (mode == Pi0Pi0Km) |
---|
1571 | answer = T(s1, k1Ma, k1Ga, k1Wa) |
---|
1572 | * T(piM, kM, s2, kstarMa, kstarGa, kstarWa); |
---|
1573 | // K-, pi-, pi+ decay. |
---|
1574 | else if (mode == PimPipKm) |
---|
1575 | answer = T(s1, k1Mb, k1Gb, k1Wb) |
---|
1576 | * T(piM, piM, s2, rhoMa, rhoGa, rhoWa); |
---|
1577 | // pi-, Kbar0, pi0 decay. |
---|
1578 | else if (mode == Pi0PimK0b) |
---|
1579 | answer = T(s1, k1Ma, k1Ga, k1Wa) |
---|
1580 | * (T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1581 | - T(piM, kM, s4, kstarMa, kstarGa, kstarWa)); |
---|
1582 | return -1.0 / 3.0 * answer; |
---|
1583 | } |
---|
1584 | |
---|
1585 | //-------------------------------------------------------------------------- |
---|
1586 | |
---|
1587 | // Return the second form factor. |
---|
1588 | |
---|
1589 | complex HMETau2ThreeMesonsWithKaons::F2() { |
---|
1590 | |
---|
1591 | complex answer; |
---|
1592 | // K-, pi-, K+ decay. |
---|
1593 | if (mode == PimKmKp) |
---|
1594 | answer = a1BW * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) / 2.0; |
---|
1595 | // K0, pi-, Kbar0 decay. |
---|
1596 | else if (mode == PimK0bK0) |
---|
1597 | answer = a1BW * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) / 2.0; |
---|
1598 | // K_S0, pi-, K_S0 decay and K_L0, pi-, K_L0 decay. |
---|
1599 | else if (mode == PimKsKs || mode == KlKlPim) |
---|
1600 | answer = a1BW * T(piM, kM, s4, kstarMa, kstarGa, kstarWa) / 2.0; |
---|
1601 | // K_S0, pi-, K_L0 decay. |
---|
1602 | else if (mode == KlPimKs) |
---|
1603 | answer = a1BW * (2.0 * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1604 | + T(piM, kM, s4, kstarMa, kstarGa, kstarWa)) / 2.0; |
---|
1605 | // K-, pi0, K0 decay. |
---|
1606 | else if (mode == Pi0K0Km) |
---|
1607 | answer = a1BW * (2.0 * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1608 | + T(piM, kM, s4, kstarMa, kstarGa, kstarWa)) / 2.0; |
---|
1609 | // pi0, pi0, K- decay. |
---|
1610 | else if (mode == Pi0Pi0Km) |
---|
1611 | answer = T(s1, k1Ma, k1Ga, k1Wa) |
---|
1612 | * T(piM, kM, s3, kstarMa, kstarGa, kstarWa); |
---|
1613 | // K-, pi-, pi+ decay. |
---|
1614 | else if (mode == PimPipKm) |
---|
1615 | answer = T(s1, k1Ma, k1Ga, k1Wa) |
---|
1616 | * T(piM, kM, s3, kstarMa, kstarGa, kstarWa); |
---|
1617 | // pi-, Kbar0, pi0 decay. |
---|
1618 | else if (mode == Pi0PimK0b) |
---|
1619 | answer = 2.0 * T(s1, k1Mb, k1Gb, k1Wb) |
---|
1620 | * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1621 | + T(s1, k1Ma, k1Ga, k1Wa) * T(piM, kM, s4, kstarMa, kstarGa, kstarWa); |
---|
1622 | return 1.0 / 3.0 * answer; |
---|
1623 | |
---|
1624 | } |
---|
1625 | |
---|
1626 | //-------------------------------------------------------------------------- |
---|
1627 | |
---|
1628 | // Return the fourth form factor. |
---|
1629 | |
---|
1630 | complex HMETau2ThreeMesonsWithKaons::F4() { |
---|
1631 | |
---|
1632 | complex answer; |
---|
1633 | // K-, pi-, K+ decay. |
---|
1634 | if (mode == PimKmKp) |
---|
1635 | answer = (sqrt(2) - 1) * T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1636 | * (sqrt(2) * T(s3, omegaM, omegaG, omegaW) |
---|
1637 | + T(piM, kM, s2, kstarMa, kstarGa, kstarWa)); |
---|
1638 | // K0, pi-, Kbar0 decay. |
---|
1639 | else if (mode == PimK0bK0) |
---|
1640 | answer = -(sqrt(2) - 1) * T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1641 | * (sqrt(2) * T(s3, omegaM, omegaG, omegaW) |
---|
1642 | + T(piM, kM, s2, kstarMa, kstarGa, kstarWa)); |
---|
1643 | // K_S0, pi-, K_S0 decay and K_L0, pi-, K_L0 decay. |
---|
1644 | else if (mode == PimKsKs || mode == KlKlPim) |
---|
1645 | answer = (sqrt(2) - 1) * T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1646 | * (T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1647 | - T(piM, kM, s4, kstarMa, kstarGa, kstarWa)); |
---|
1648 | // K_S0, pi-, K_L0 decay. |
---|
1649 | else if (mode == KlPimKs) |
---|
1650 | answer = -(sqrt(2) - 1) * T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1651 | * (2 * sqrt(2) * T(s3, omegaM, omegaG, omegaW) |
---|
1652 | + T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1653 | + T(piM, kM, s4, kstarMa, kstarGa, kstarWa)); |
---|
1654 | // K-, pi0, K0 decay. |
---|
1655 | else if (mode == Pi0K0Km) |
---|
1656 | answer = -(sqrt(2) - 1) * T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1657 | * (T(piM, kM, s4, kstarMa, kstarGa, kstarWa) |
---|
1658 | - T(piM, kM, s2, kstarMa, kstarGa, kstarWa)); |
---|
1659 | // pi0, pi0, K- decay. |
---|
1660 | else if (mode == Pi0Pi0Km) |
---|
1661 | answer = T(piM, kM, s1, kstarMv, kstarGv, kstarWv) |
---|
1662 | * (T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1663 | - T(piM, kM, s3, kstarMa, kstarGa, kstarWa)); |
---|
1664 | // K-, pi-, pi+ decay. |
---|
1665 | else if (mode == PimPipKm) |
---|
1666 | answer = -T(piM, kM, s1, kstarMv, kstarGv, kstarWv) |
---|
1667 | * (T(piM, piM, s2, rhoMa, rhoGa, rhoWa) |
---|
1668 | + T(piM, kM, s3, kstarMa, kstarGa, kstarWa)); |
---|
1669 | // pi-, Kbar0, pi0 decay. |
---|
1670 | else if (mode == Pi0PimK0b) |
---|
1671 | answer = T(piM, kM, s1, kstarMv, kstarGv, kstarWv) |
---|
1672 | * (2.0 * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1673 | + T(piM, kM, s2, kstarMa, kstarGa, kstarWa) |
---|
1674 | + T(piM, kM, s4, kstarMa, kstarGa, kstarWa)); |
---|
1675 | return 1.0 / (8.0 * M_PI * M_PI * piW * piW) * answer; |
---|
1676 | |
---|
1677 | } |
---|
1678 | |
---|
1679 | //========================================================================== |
---|
1680 | |
---|
1681 | // Tau decay matrix element for tau decay into three mesons. The form |
---|
1682 | // factors are taken from Comp. Phys. Com. 76 (1993) 361-380. |
---|
1683 | |
---|
1684 | // rhoMa(v): on-shell masses for the axial (vector) rho resonances |
---|
1685 | // rhoGa(v): widths for the axial (vector) rho resonances |
---|
1686 | // rhoWa(v): weights for the axial (vector) rho resonances |
---|
1687 | // kstarX: on-shell masses, widths, and weights for the K* resonances |
---|
1688 | // k1X: on-shell masses, width, and weight for the K1 resonances |
---|
1689 | // kM: kaon mass |
---|
1690 | // piM: charged pion mass |
---|
1691 | // piW: pion coupling, f_W |
---|
1692 | |
---|
1693 | //-------------------------------------------------------------------------- |
---|
1694 | |
---|
1695 | // Initialize resonances for the helicity matrix element. |
---|
1696 | |
---|
1697 | void HMETau2ThreeMesonsGeneric::initResonances() { |
---|
1698 | |
---|
1699 | // pi-, pi-, pi+ decay and pi0, pi0, pi- decay. |
---|
1700 | if (mode == PimPimPip || mode == Pi0Pi0Pim) DECAYWEIGHTMAX = 1.3e4; |
---|
1701 | // K-, pi-, K+ decay. |
---|
1702 | else if (mode == PimKmKp) DECAYWEIGHTMAX = 330; |
---|
1703 | // K0, pi-, Kbar0 decay. |
---|
1704 | else if (mode == PimK0bK0) DECAYWEIGHTMAX = 300; |
---|
1705 | // K-, pi0, K0 decay. |
---|
1706 | else if (mode == Pi0K0Km) DECAYWEIGHTMAX = 40; |
---|
1707 | // pi0, pi0, K- decay. |
---|
1708 | else if (mode == Pi0Pi0Km) DECAYWEIGHTMAX = 9.4e4; |
---|
1709 | // K-, pi-, pi+ decay. |
---|
1710 | else if (mode == PimPipKm) DECAYWEIGHTMAX = 9.0e3; |
---|
1711 | // pi-, Kbar0, pi0 decay. |
---|
1712 | else if (mode == Pi0PimK0b) DECAYWEIGHTMAX = 1.2e4; |
---|
1713 | // pi-, pi0, eta decay. |
---|
1714 | else if (mode == Pi0PimEta) DECAYWEIGHTMAX = 360; |
---|
1715 | |
---|
1716 | // Clear the vectors from previous decays. |
---|
1717 | rhoMa.clear(); rhoGa.clear(); rhoWa.clear(); |
---|
1718 | rhoMv.clear(); rhoGv.clear(); rhoWv.clear(); |
---|
1719 | kstarM.clear(); kstarG.clear(); kstarW.clear(); |
---|
1720 | k1M.clear(); k1G.clear(); k1W.clear(); |
---|
1721 | |
---|
1722 | // Rho parameters. |
---|
1723 | rhoMa.push_back(0.773); rhoGa.push_back(0.145); rhoWa.push_back(1); |
---|
1724 | rhoMa.push_back(1.370); rhoGa.push_back(0.510); rhoWa.push_back(-0.145); |
---|
1725 | rhoMv.push_back(0.773); rhoGv.push_back(0.145); rhoWv.push_back(-26); |
---|
1726 | rhoMv.push_back(1.5); rhoGv.push_back(0.220); rhoWv.push_back(6.5); |
---|
1727 | rhoMv.push_back(1.75); rhoGv.push_back(0.120); rhoWv.push_back(1); |
---|
1728 | |
---|
1729 | // Kaon parameters. |
---|
1730 | kstarM.push_back(0.892); kstarG.push_back(0.0513); kstarW.push_back(1); |
---|
1731 | k1M.push_back(1.402); k1G.push_back(0.174); k1W.push_back(1); |
---|
1732 | |
---|
1733 | // Kaon and pion parameters |
---|
1734 | kM = 0.49765; piM = 0.13957; piW = 0.0942; |
---|
1735 | |
---|
1736 | } |
---|
1737 | |
---|
1738 | //-------------------------------------------------------------------------- |
---|
1739 | |
---|
1740 | // Return the first form factor. |
---|
1741 | |
---|
1742 | complex HMETau2ThreeMesonsGeneric::F1() { |
---|
1743 | |
---|
1744 | complex answer; |
---|
1745 | // pi-, pi-, pi+ decay and pi0, pi0, pi- decay. |
---|
1746 | if (mode == PimPimPip || mode == Pi0Pi0Pim) |
---|
1747 | answer = a1BW * T(piM, piM, s2, rhoMa, rhoGa, rhoWa); |
---|
1748 | // K-, pi-, K+ decay. |
---|
1749 | else if (mode == PimKmKp) |
---|
1750 | answer = -a1BW * T(piM, kM, s2, kstarM, kstarG, kstarW) / 3.0; |
---|
1751 | // K0, pi-, Kbar0 decay. |
---|
1752 | else if (mode == PimK0bK0) |
---|
1753 | answer = -a1BW * T(piM, kM, s2, kstarM, kstarG, kstarW) / 3.0; |
---|
1754 | // K-, pi0, K0 decay. |
---|
1755 | else if (mode == Pi0K0Km) |
---|
1756 | answer = 0; |
---|
1757 | // pi0, pi0, K- decay. |
---|
1758 | else if (mode == Pi0Pi0Km) |
---|
1759 | answer = T(s1, k1M, k1G, k1W) * T(piM, kM, s2, kstarM, kstarG, kstarW); |
---|
1760 | // K-, pi-, pi+ decay. |
---|
1761 | else if (mode == PimPipKm) |
---|
1762 | answer = -T(s1, k1M, k1G, k1W) * T(piM, piM, s2, rhoMa, rhoGa, rhoWa) / 3.0; |
---|
1763 | // pi-, Kbar0, pi0 decay. |
---|
1764 | else if (mode == Pi0PimK0b) |
---|
1765 | answer = 0; |
---|
1766 | // pi-, pi0, eta decay. |
---|
1767 | else if (mode == Pi0PimEta) |
---|
1768 | answer = 0; |
---|
1769 | return answer; |
---|
1770 | |
---|
1771 | } |
---|
1772 | |
---|
1773 | //-------------------------------------------------------------------------- |
---|
1774 | |
---|
1775 | // Return the second form factor. |
---|
1776 | |
---|
1777 | complex HMETau2ThreeMesonsGeneric::F2() { |
---|
1778 | |
---|
1779 | complex answer; |
---|
1780 | // pi-, pi-, pi+ decay and pi0, pi0, pi- decay. |
---|
1781 | if (mode == PimPimPip || mode == Pi0Pi0Pim) |
---|
1782 | answer = -a1BW * T(piM, piM, s3, rhoMa, rhoGa, rhoWa); |
---|
1783 | // K-, pi-, K+ decay. |
---|
1784 | else if (mode == PimKmKp) |
---|
1785 | answer = a1BW * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) / 3.0; |
---|
1786 | // K0, pi-, Kbar0 decay. |
---|
1787 | else if (mode == PimK0bK0) |
---|
1788 | answer = a1BW * T(piM, piM, s3, rhoMa, rhoGa, rhoWa) / 3.0; |
---|
1789 | // K-, pi0, K0 decay. |
---|
1790 | else if (mode == Pi0K0Km) |
---|
1791 | answer = a1BW * T(piM, piM, s3, rhoMa, rhoGa, rhoWa); |
---|
1792 | // pi0, pi0, K- decay. |
---|
1793 | else if (mode == Pi0Pi0Km) |
---|
1794 | answer = -T(s1, k1M, k1G, k1W) * T(piM, kM, s3, kstarM, kstarG, kstarW); |
---|
1795 | // K-, pi-, pi+ decay. |
---|
1796 | else if (mode == PimPipKm) |
---|
1797 | answer = T(s1, k1M, k1G, k1W) |
---|
1798 | * T(piM, kM, s3, kstarM, kstarG, kstarW) / 3.0; |
---|
1799 | // pi-, Kbar0, pi0 decay. |
---|
1800 | else if (mode == Pi0PimK0b) |
---|
1801 | answer = T(s1, k1M, k1G, k1W) * T(piM, piM, s3, rhoMa, rhoGa, rhoWa); |
---|
1802 | // pi-, pi0, eta decay. |
---|
1803 | else if (mode == Pi0PimEta) |
---|
1804 | answer = 0; |
---|
1805 | return answer; |
---|
1806 | |
---|
1807 | } |
---|
1808 | |
---|
1809 | //-------------------------------------------------------------------------- |
---|
1810 | |
---|
1811 | // Return the fourth form factor. |
---|
1812 | |
---|
1813 | complex HMETau2ThreeMesonsGeneric::F4() { |
---|
1814 | |
---|
1815 | complex answer; |
---|
1816 | // pi-, pi-, pi+ decay and pi0, pi0, pi- decay. |
---|
1817 | if (mode == PimPimPip || mode == Pi0Pi0Pim) |
---|
1818 | answer = 0; |
---|
1819 | // K-, pi-, K+ decay. |
---|
1820 | else if (mode == PimKmKp) |
---|
1821 | answer = T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1822 | * (T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1823 | - 0.2 * T(piM, kM, s2, kstarM, kstarG, kstarW)) * (1.25); |
---|
1824 | // K0, pi-, Kbar0 decay. |
---|
1825 | else if (mode == PimK0bK0) |
---|
1826 | answer = -T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1827 | * (T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1828 | - 0.2 * T(piM, kM, s2, kstarM, kstarG, kstarW)) * (1.25); |
---|
1829 | // K-, pi0, K0 decay. |
---|
1830 | else if (mode == Pi0K0Km) |
---|
1831 | answer = 0; |
---|
1832 | // pi0, pi0, K- decay. |
---|
1833 | else if (mode == Pi0Pi0Km) |
---|
1834 | answer = 0; |
---|
1835 | // K-, pi-, pi+ decay. |
---|
1836 | else if (mode == PimPipKm) |
---|
1837 | answer = -T(piM, kM, s1, kstarM, kstarG, kstarW) |
---|
1838 | * (T(piM, piM, s2, rhoMa, rhoGa, rhoWa) |
---|
1839 | - 0.2 * T(piM, kM, s3, kstarM, kstarG, kstarW)) * (1.25); |
---|
1840 | // pi-, Kbar0, pi0 decay. |
---|
1841 | else if (mode == Pi0PimK0b) |
---|
1842 | answer = 2.0 * T(piM, kM, s1, kstarM, kstarG, kstarW) |
---|
1843 | * (T(piM, piM, s3, rhoMa, rhoGa, rhoWa) |
---|
1844 | - 0.2 * T(piM, kM, s2, kstarM, kstarG, kstarW)) * (1.25); |
---|
1845 | // pi-, pi0, eta decay. |
---|
1846 | else if (mode == Pi0PimEta) |
---|
1847 | answer = T(piM, piM, s1, rhoMv, rhoGv, rhoWv) |
---|
1848 | * T(piM, piM, s4, rhoMa, rhoGa, rhoWa); |
---|
1849 | return 1.0 / (4.0 * M_PI * M_PI * piW * piW) * answer; |
---|
1850 | |
---|
1851 | } |
---|
1852 | |
---|
1853 | //========================================================================== |
---|
1854 | |
---|
1855 | // Tau decay matrix element for tau decay into two pions with a photons taken |
---|
1856 | // from Comp. Phys. Com. 76 (1993) 361-380. Because a photon is in the final |
---|
1857 | // state the matrix element is reimplented to handle the two spin states. |
---|
1858 | |
---|
1859 | // F(s, M, G, W): form factor for resonance |
---|
1860 | // rhoM: on-shell mass of the rho resonances |
---|
1861 | // rhoG: width of the rho resonances |
---|
1862 | // rhoW: weight of the rho resonances |
---|
1863 | // omegaX: on-shell mass, width, and weight of the omega resonances |
---|
1864 | // piM: charged pion mass |
---|
1865 | |
---|
1866 | //-------------------------------------------------------------------------- |
---|
1867 | |
---|
1868 | // Initialize constants for the helicity matrix element. |
---|
1869 | |
---|
1870 | void HMETau2TwoPionsGamma::initConstants() { |
---|
1871 | |
---|
1872 | DECAYWEIGHTMAX = 4e4; |
---|
1873 | |
---|
1874 | // Clear the vectors from previous decays. |
---|
1875 | rhoM.clear(); rhoG.clear(); rhoW.clear(); |
---|
1876 | omegaM.clear(); omegaG.clear(); omegaW.clear(); |
---|
1877 | |
---|
1878 | // Set parameters. |
---|
1879 | rhoM.push_back(0.773); rhoG.push_back(0.145); rhoW.push_back(1); |
---|
1880 | rhoM.push_back(1.7); rhoG.push_back(0.26); rhoW.push_back(-0.1); |
---|
1881 | omegaM.push_back(0.782); omegaG.push_back(0.0085); omegaW.push_back(1); |
---|
1882 | piM = 0.13957; |
---|
1883 | |
---|
1884 | } |
---|
1885 | |
---|
1886 | //-------------------------------------------------------------------------- |
---|
1887 | |
---|
1888 | // Initialize wave functions for the helicity matrix element. |
---|
1889 | void HMETau2TwoPionsGamma::initWaves(vector<HelicityParticle>& p) { |
---|
1890 | |
---|
1891 | // Calculate the hadronic current. |
---|
1892 | u.clear(); |
---|
1893 | pMap.resize(p.size()); |
---|
1894 | setFermionLine(0, p[0], p[1]); |
---|
1895 | |
---|
1896 | // Calculate the hadronic current. |
---|
1897 | vector< Wave4 > u2; |
---|
1898 | Wave4 q(p[2].p() + p[3].p() + p[4].p()); |
---|
1899 | Wave4 q2(p[2].p()), q3(p[3].p()), q4(p[4].p()); |
---|
1900 | double s1 = m2(q); |
---|
1901 | double s2 = m2(q3 + q2); |
---|
1902 | complex f = F(s1, rhoM, rhoG, rhoW) * F(0, rhoM, rhoG, rhoW) |
---|
1903 | * F(s2, omegaM, omegaG, omegaW); |
---|
1904 | double q4q2 = m2(q4, q2); |
---|
1905 | double q4q3 = m2(q4, q3); |
---|
1906 | double q3q2 = m2(q3, q2); |
---|
1907 | for (int h = 0; h < 2; h++) { |
---|
1908 | Wave4 e = p[2].wave(h); |
---|
1909 | complex q4e = q4*gamma[4]*e; |
---|
1910 | complex q3e = q3*gamma[4]*e; |
---|
1911 | u2.push_back(f * (e * (piM*piM*q4q2 - q3q2*(q4q3 - q4q2)) |
---|
1912 | - q3 * (q3e*q4q2 - q4e*q3q2) |
---|
1913 | + q2 * (q3e*q4q3 - q4e*(piM*piM + q3q2)))); |
---|
1914 | } |
---|
1915 | u.push_back(u2); |
---|
1916 | |
---|
1917 | } |
---|
1918 | |
---|
1919 | //-------------------------------------------------------------------------- |
---|
1920 | |
---|
1921 | // Return element for the helicity matrix element. |
---|
1922 | complex HMETau2TwoPionsGamma::calculateME(vector<int> h) { |
---|
1923 | |
---|
1924 | complex answer(0,0); |
---|
1925 | for (int mu = 0; mu <= 3; mu++) { |
---|
1926 | answer += |
---|
1927 | (u[1][h[pMap[1]]] * gamma[mu] * (1 - gamma[5]) * u[0][h[pMap[0]]]) |
---|
1928 | * gamma[4](mu,mu) * u[2][h[2]](mu); |
---|
1929 | } |
---|
1930 | return answer; |
---|
1931 | |
---|
1932 | } |
---|
1933 | |
---|
1934 | //-------------------------------------------------------------------------- |
---|
1935 | |
---|
1936 | // Return the form factor. |
---|
1937 | complex HMETau2TwoPionsGamma::F(double s, vector<double> M, vector<double> G, |
---|
1938 | vector<double> W) { |
---|
1939 | |
---|
1940 | complex answer(0, 0); |
---|
1941 | for (unsigned int i = 0; i < M.size(); i++) |
---|
1942 | answer += W[i] / (M[i]*M[i] - s - complex(0, 1) * M[i] * G[i]); |
---|
1943 | return answer; |
---|
1944 | |
---|
1945 | } |
---|
1946 | |
---|
1947 | //========================================================================== |
---|
1948 | |
---|
1949 | // Tau decay matrix element for tau decay into pions using the Novosibirsk |
---|
1950 | // model of Comp. Phys. Com. 174 (2006) 818-835. |
---|
1951 | |
---|
1952 | // G(i,s): G-factor for index 1, 2, or 3 |
---|
1953 | // tX(q,q1,q2,q3,q4): combined resonance current |
---|
1954 | // a1D(s): a1 Breit-Wigner denominator |
---|
1955 | // rhoD(s): rho Breit-Wigner denominator |
---|
1956 | // sigD(s): sigma Breit-Wigner denominator |
---|
1957 | // omeD(s): omega Breit-Wigner denominator |
---|
1958 | // a1FormFactor(s): form factor for the a1 resonance |
---|
1959 | // rhoFormFactor1(2)(s): form factor for the rho resonance |
---|
1960 | // omeFormFactor(s): form factor for the omega resonance |
---|
1961 | // sigM: on-shell mass of the sigma resonance |
---|
1962 | // sigG: width of the sigma resonance |
---|
1963 | // sigA: amplitude of the sigma resonance |
---|
1964 | // sigP: phase of the sigma resonance |
---|
1965 | // sigW: weight of the sigma resonance (from amplitude and weight) |
---|
1966 | // omeX: mass, width, amplitude, phase, and weight of the omega resonance |
---|
1967 | // a1X: mass and width of the a1 resonance |
---|
1968 | // rhoX: mass and width of the rho resonance |
---|
1969 | // picM: charge pion mass |
---|
1970 | // pinM: neutral pion mass |
---|
1971 | // lambda2: a1 form factor cut-off |
---|
1972 | |
---|
1973 | //-------------------------------------------------------------------------- |
---|
1974 | |
---|
1975 | // Initialize constants for the helicity matrix element. |
---|
1976 | |
---|
1977 | void HMETau2FourPions::initConstants() { |
---|
1978 | |
---|
1979 | if (abs(pID[3]) == 111) DECAYWEIGHTMAX = 5e8; |
---|
1980 | else DECAYWEIGHTMAX = 5e9; |
---|
1981 | pinM = particleDataPtr->m0(111); |
---|
1982 | picM = particleDataPtr->m0(211); |
---|
1983 | sigM = 0.8; omeM = 0.782; a1M = 1.23; rhoM = 0.7761; |
---|
1984 | sigG = 0.8; omeG = 0.00841; a1G = 0.45; rhoG = 0.1445; |
---|
1985 | sigP = 0.43585; omeP = 0.0; |
---|
1986 | sigA = 1.39987; omeA = 1.0; |
---|
1987 | sigW = sigA*(cos(sigP)+complex(0,1)*sin(sigP)); |
---|
1988 | omeW = omeA*(cos(omeP)+complex(0,1)*sin(omeP)); |
---|
1989 | lambda2 = 1.2; |
---|
1990 | |
---|
1991 | } |
---|
1992 | |
---|
1993 | //-------------------------------------------------------------------------- |
---|
1994 | |
---|
1995 | // Initialize the hadronic current for the helicity matrix element. |
---|
1996 | |
---|
1997 | void HMETau2FourPions::initHadronicCurrent(vector<HelicityParticle>& p) { |
---|
1998 | |
---|
1999 | vector< Wave4 > u2; |
---|
2000 | |
---|
2001 | // Create pion momenta. |
---|
2002 | Wave4 q(p[2].p() + p[3].p() + p[4].p()+ p[5].p()); |
---|
2003 | Wave4 q2(p[2].p()), q3(p[3].p()), q4(p[4].p()), q5(p[5].p()); |
---|
2004 | |
---|
2005 | // Calculate the four pion system energy. |
---|
2006 | double s = m2(q); |
---|
2007 | |
---|
2008 | // Create the hadronic current for the 3 neutral pion channel. |
---|
2009 | if (abs(pID[3]) == 111) |
---|
2010 | u2.push_back(G(1,s)*(t1(q,q3,q4,q5,q2) + t1(q,q3,q2,q5,q4) + |
---|
2011 | t1(q,q4,q3,q5,q2) + t1(q,q4,q2,q5,q3) + |
---|
2012 | t1(q,q2,q3,q5,q4) + t1(q,q2,q4,q5,q3) + |
---|
2013 | t2(q,q3,q5,q4,q2) + t2(q,q4,q5,q3,q2) + |
---|
2014 | t2(q,q2,q5,q4,q3) - t2(q,q5,q3,q4,q2) - |
---|
2015 | t2(q,q5,q4,q3,q2) - t2(q,q5,q2,q4,q3))); |
---|
2016 | |
---|
2017 | // Create the hadronic current for the 3 charged pion channel. |
---|
2018 | else if (abs(pID[3]) == 211) |
---|
2019 | u2.push_back(G(2,s)*(t1(q,q3,q5,q4,q2) + t1(q,q4,q5,q3,q2) + |
---|
2020 | t1(q,q3,q4,q5,q2) + t1(q,q4,q3,q5,q2) + |
---|
2021 | t1(q,q2,q4,q3,q5) + t1(q,q2,q3,q4,q5) + |
---|
2022 | t2(q,q2,q4,q3,q5) + t2(q,q2,q3,q4,q5) - |
---|
2023 | t2(q,q3,q2,q4,q5) - t2(q,q4,q2,q3,q5)) + |
---|
2024 | G(3,s)*(t3(q,q3,q5,q4,q2) + t3(q,q4,q5,q3,q2) - |
---|
2025 | t3(q,q3,q4,q5,q2) - t3(q,q4,q3,q5,q2) - |
---|
2026 | t3(q,q3,q2,q4,q5) - t3(q,q4,q2,q3,q5))); |
---|
2027 | u.push_back(u2); |
---|
2028 | |
---|
2029 | } |
---|
2030 | |
---|
2031 | //-------------------------------------------------------------------------- |
---|
2032 | |
---|
2033 | // Return the first t-vector. |
---|
2034 | |
---|
2035 | Wave4 HMETau2FourPions::t1(Wave4 &q, Wave4 &q1, Wave4 &q2, |
---|
2036 | Wave4 &q3, Wave4 &q4) { |
---|
2037 | |
---|
2038 | Wave4 a1Q(q2 + q3 + q4); |
---|
2039 | Wave4 rhoQ(q3 + q4); |
---|
2040 | double a1S = m2(a1Q); |
---|
2041 | double rhoS = m2(rhoQ); |
---|
2042 | |
---|
2043 | // Needed to match Herwig++. |
---|
2044 | double gM = sqrtpos(rhoM*rhoM - 4*picM*picM) * (rhoM*rhoM - 4*picM*picM) |
---|
2045 | / rhoM; |
---|
2046 | double dm = (rhoFormFactor1(0) - rhoFormFactor1(rhoM*rhoM) + |
---|
2047 | rhoM*rhoM * rhoFormFactor2(rhoM*rhoM)) / gM; |
---|
2048 | return - a1FormFactor(a1S) / (a1D(a1S) * rhoD(rhoS)) * pow2(a1M) * |
---|
2049 | (rhoM*rhoM + rhoM*rhoG*dm) * |
---|
2050 | (m2(q,a1Q) * (m2(q3,a1Q) * q4 - m2(q4,a1Q) * q3) + |
---|
2051 | (m2(q,q4) * m2(q1,q3) - m2(q,q3) * m2(q1,q4)) * a1Q); |
---|
2052 | |
---|
2053 | } |
---|
2054 | |
---|
2055 | //-------------------------------------------------------------------------- |
---|
2056 | |
---|
2057 | // Return the second t-vector. |
---|
2058 | |
---|
2059 | Wave4 HMETau2FourPions::t2(Wave4 &q, Wave4 &/*q1*/, Wave4 &q2, |
---|
2060 | Wave4 &q3, Wave4 &q4) { |
---|
2061 | |
---|
2062 | Wave4 a1Q(q2 + q3 + q4); |
---|
2063 | Wave4 sigQ(q3 + q4); |
---|
2064 | double a1S = m2(a1Q); |
---|
2065 | double sigS = m2(sigQ); |
---|
2066 | return sigW * a1FormFactor(a1S) / (a1D(a1S) * sigD(sigS)) * |
---|
2067 | pow2(a1M) * pow2(sigM) * (m2(q,a1Q) * a1S * q2 - m2(q,q2) * a1S * a1Q); |
---|
2068 | |
---|
2069 | } |
---|
2070 | |
---|
2071 | //-------------------------------------------------------------------------- |
---|
2072 | |
---|
2073 | // Return the third t-vector. |
---|
2074 | |
---|
2075 | Wave4 HMETau2FourPions::t3(Wave4 &q, Wave4 &q1, Wave4 &q2, |
---|
2076 | Wave4 &q3, Wave4 &q4) { |
---|
2077 | Wave4 omeQ(q2 + q3 + q4); |
---|
2078 | Wave4 rhoQ(q3 + q4); |
---|
2079 | double omeS = m2(omeQ); |
---|
2080 | double rhoS = m2(rhoQ); |
---|
2081 | |
---|
2082 | // Needed to match Herwig++. |
---|
2083 | double gM = sqrtpos(rhoM*rhoM - 4*picM*picM) * (rhoM*rhoM - 4*picM*picM) |
---|
2084 | / rhoM; |
---|
2085 | double dm = (rhoFormFactor1(0) - rhoFormFactor1(rhoM*rhoM) + |
---|
2086 | rhoM*rhoM * rhoFormFactor2(rhoM*rhoM)) / gM; |
---|
2087 | return omeW * omeFormFactor(omeS) / (omeD(omeS) * rhoD(rhoS)) * |
---|
2088 | pow2(omeM) * (rhoM*rhoM + rhoM*rhoG*dm) * |
---|
2089 | ((m2(q,q3) * m2(q1,q4) - m2(q,q4) * m2(q1,q3)) * q2 + |
---|
2090 | (m2(q,q4) * m2(q1,q2) - m2(q,q2) * m2(q1,q4)) * q3 + |
---|
2091 | (m2(q,q2) * m2(q1,q3) - m2(q,q3) * m2(q1,q2)) * q4); |
---|
2092 | |
---|
2093 | } |
---|
2094 | |
---|
2095 | //-------------------------------------------------------------------------- |
---|
2096 | |
---|
2097 | // Return the D function for the a1(1260). |
---|
2098 | |
---|
2099 | complex HMETau2FourPions::a1D(double s) { |
---|
2100 | |
---|
2101 | // rG is defined as the running width. |
---|
2102 | double rG = 0; |
---|
2103 | |
---|
2104 | // The rho and pion cut off thresholds defined in the fit. |
---|
2105 | double piM = 0.16960; |
---|
2106 | double rM = 0.83425; |
---|
2107 | |
---|
2108 | // Fit of width below three pion mass threshold. |
---|
2109 | if (s < piM) |
---|
2110 | rG = 0; |
---|
2111 | |
---|
2112 | // Fit of width below pion and rho mass threshold. |
---|
2113 | else if (s < rM) |
---|
2114 | rG = 0.003052*pow3(s - piM)*(1.0 + 151.088*(s - piM) + |
---|
2115 | 174.495*pow2(s - piM)); |
---|
2116 | |
---|
2117 | // Fit of width above pion and rho mass threshold. |
---|
2118 | else |
---|
2119 | rG = 2.60817 - 2.47790*s + 0.66539*pow2(s) - 0.0678183*pow3(s) + |
---|
2120 | 1.66577*(s-1.23701)/s; |
---|
2121 | return s - a1M*a1M + complex(0,1) * sqrtpos(s) * rG; |
---|
2122 | |
---|
2123 | } |
---|
2124 | |
---|
2125 | //-------------------------------------------------------------------------- |
---|
2126 | |
---|
2127 | // Return the D function for the rho(770). |
---|
2128 | |
---|
2129 | complex HMETau2FourPions::rhoD(double s) { |
---|
2130 | |
---|
2131 | double gQ = sqrtpos(s - 4*picM*picM) * (s - 4*picM*picM) / sqrtpos(s); |
---|
2132 | double gM = sqrtpos(rhoM*rhoM - 4*picM*picM) * (rhoM*rhoM - 4*picM*picM) |
---|
2133 | / rhoM; |
---|
2134 | double dm = (rhoFormFactor1(s) - rhoFormFactor1(rhoM*rhoM) - |
---|
2135 | (s - rhoM*rhoM) * rhoFormFactor2(rhoM*rhoM)) / gM; |
---|
2136 | |
---|
2137 | // Ensure complex part is zero below available channel. |
---|
2138 | if (s < 4*picM*picM) gQ = 0; |
---|
2139 | return s - rhoM*rhoM - rhoM*rhoG*dm + complex(0,1)*rhoM*rhoG*(gQ/gM); |
---|
2140 | |
---|
2141 | } |
---|
2142 | |
---|
2143 | //-------------------------------------------------------------------------- |
---|
2144 | |
---|
2145 | // Return the D function for the sigma(800) (just s-wave running width). |
---|
2146 | |
---|
2147 | complex HMETau2FourPions::sigD(double s) { |
---|
2148 | |
---|
2149 | // Sigma decay to two neutral pions for three neutral pion channel. |
---|
2150 | double piM = abs(pID[3]) == 111 ? pinM : picM; |
---|
2151 | double gQ = sqrtpos(1.0 - 4*piM*piM/s); |
---|
2152 | double gM = sqrtpos(1.0 - 4*piM*piM/(sigM*sigM)); |
---|
2153 | return s - sigM*sigM + complex(0,1)*sigM*sigG*gQ/gM; |
---|
2154 | |
---|
2155 | } |
---|
2156 | |
---|
2157 | //-------------------------------------------------------------------------- |
---|
2158 | |
---|
2159 | // Return the D function for the omega(782). |
---|
2160 | |
---|
2161 | complex HMETau2FourPions::omeD(double s) { |
---|
2162 | |
---|
2163 | double g = 0; |
---|
2164 | double q = sqrtpos(s); |
---|
2165 | double x = q - omeM; |
---|
2166 | |
---|
2167 | // Fit of width given in TAUOLA. |
---|
2168 | if (s < 1) |
---|
2169 | g = 1 + 17.560*x + 141.110*pow2(x) + 894.884*pow3(x) + 4977.35*pow4(x) + |
---|
2170 | 7610.66*pow5(x) - 42524.4*pow6(x); |
---|
2171 | else |
---|
2172 | g = -1333.26 + 4860*q - 6000.81*pow2(q) + 2504.97*pow3(q); |
---|
2173 | if (g < 0) g = 0; |
---|
2174 | return s - omeM*omeM + complex(0,1)*omeM*omeG*g; |
---|
2175 | |
---|
2176 | } |
---|
2177 | |
---|
2178 | //-------------------------------------------------------------------------- |
---|
2179 | |
---|
2180 | // Return the form factor for the a1. |
---|
2181 | |
---|
2182 | double HMETau2FourPions::a1FormFactor(double s) { |
---|
2183 | |
---|
2184 | return pow2((1.0 + a1M*a1M/lambda2) / (1.0 + s/lambda2)); |
---|
2185 | |
---|
2186 | } |
---|
2187 | |
---|
2188 | //-------------------------------------------------------------------------- |
---|
2189 | |
---|
2190 | // Return the form factor for the rho(770) (equivalent to h(s) in TAUOLA). |
---|
2191 | |
---|
2192 | double HMETau2FourPions::rhoFormFactor1(double s) { |
---|
2193 | |
---|
2194 | double f = sqrtpos(1 - 4*picM*picM/s); |
---|
2195 | if (s > 4*picM*picM) |
---|
2196 | f = f * log((1 + f) / (1 - f)) * (s - 4*picM*picM) / M_PI; |
---|
2197 | else if (s < 0.0000001) |
---|
2198 | f = -8 * picM*picM / M_PI; |
---|
2199 | else |
---|
2200 | f = 0; |
---|
2201 | return f; |
---|
2202 | |
---|
2203 | } |
---|
2204 | |
---|
2205 | //-------------------------------------------------------------------------- |
---|
2206 | |
---|
2207 | // Return the form factor for the rho(770) (equivalent to h(s) derivative). |
---|
2208 | |
---|
2209 | double HMETau2FourPions::rhoFormFactor2(double s) { |
---|
2210 | |
---|
2211 | double f = sqrtpos(1 - 4*picM*picM/s); |
---|
2212 | if (s > 4*picM*picM) |
---|
2213 | f = f / (M_PI * s) * (s*f + (2*picM*picM + s)*log((1 + f) / (1 - f))); |
---|
2214 | else |
---|
2215 | f = 0; |
---|
2216 | return f; |
---|
2217 | |
---|
2218 | } |
---|
2219 | |
---|
2220 | //-------------------------------------------------------------------------- |
---|
2221 | |
---|
2222 | // Return the form factor for the omega(782). |
---|
2223 | |
---|
2224 | double HMETau2FourPions::omeFormFactor(double /*s*/) { |
---|
2225 | |
---|
2226 | return 1.0; |
---|
2227 | |
---|
2228 | } |
---|
2229 | |
---|
2230 | //-------------------------------------------------------------------------- |
---|
2231 | |
---|
2232 | // Return the G-functions given in TAUOLA using a piece-wise fit. |
---|
2233 | |
---|
2234 | double HMETau2FourPions::G(int i, double s) { |
---|
2235 | |
---|
2236 | // Break points for the fits. |
---|
2237 | double s0(0), s1(0), s2(0), s3(0), s4(0), s5(0); |
---|
2238 | |
---|
2239 | // Parameters for the fits. |
---|
2240 | double a1(0), b1(0); |
---|
2241 | double a2(0), b2(0), c2(0), d2(0), e2(0); |
---|
2242 | double a3(0), b3(0), c3(0), d3(0), e3(0); |
---|
2243 | double a4(0), b4(0); |
---|
2244 | double a5(0), b5(0); |
---|
2245 | |
---|
2246 | // Three neutral pion parameters. |
---|
2247 | if (i == 1) { |
---|
2248 | s0 = 0.614403; s1 = 0.656264; s2 = 1.57896; |
---|
2249 | s3 = 3.08198; s4 = 3.12825; s5 = 3.17488; |
---|
2250 | a1 = -23383.7; b1 = 38059.2; |
---|
2251 | a2 = 230.368; b2 = -4.39368; c2 = 687.002; |
---|
2252 | d2 = -732.581; e2 = 207.087; |
---|
2253 | a3 = 1633.92; b3 = -2596.21; c3 = 1703.08; |
---|
2254 | d3 = -501.407; e3 = 54.5919; |
---|
2255 | a4 = -2982.44; b4 = 986.009; |
---|
2256 | a5 = 6948.99; b5 = -2188.74; |
---|
2257 | } |
---|
2258 | |
---|
2259 | // Three charged pion parameters. |
---|
2260 | else if (i == 2) { |
---|
2261 | s0 = 0.614403; s1 = 0.635161; s2 = 2.30794; |
---|
2262 | s3 = 3.08198; s4 = 3.12825; s5 = 3.17488; |
---|
2263 | a1 = -54171.5; b1 = 88169.3; |
---|
2264 | a2 = 454.638; b2 = -3.07152; c2 = -48.7086; |
---|
2265 | d2 = 81.9702; e2 = -24.0564; |
---|
2266 | a3 = -162.421; b3 = 308.977; c3 = -27.7887; |
---|
2267 | d3 = -48.5957; e3 = 10.6168; |
---|
2268 | a4 = -2650.29; b4 = 879.776; |
---|
2269 | a5 = 6936.99; b5 = -2184.97; |
---|
2270 | } |
---|
2271 | |
---|
2272 | // Omega mediated three charged pion parameters. |
---|
2273 | else if (i == 3) { |
---|
2274 | s0 = 0.81364; s1 = 0.861709; s2 = 1.92621; |
---|
2275 | s3 = 3.08198; s4 = 3.12825; s5 = 3.17488; |
---|
2276 | a1 = -84888.9; b1 = 104332; |
---|
2277 | a2 = 2698.15; b2 = -3.08302; c2 = 1936.11; |
---|
2278 | d2 = -1254.59; e2 = 201.291; |
---|
2279 | a3 = 7171.65; b3 = -6387.9; c3 = 3056.27; |
---|
2280 | d3 = -888.63; e3 = 108.632; |
---|
2281 | a4 = -5607.48; b4 = 1917.27; |
---|
2282 | a5 = 26573; b5 = -8369.76; |
---|
2283 | } |
---|
2284 | |
---|
2285 | // Return the appropriate fit. |
---|
2286 | if (s < s0) |
---|
2287 | return 0.0; |
---|
2288 | else if (s < s1) |
---|
2289 | return a1 + b1*s; |
---|
2290 | else if (s < s2) |
---|
2291 | return a2*pow(s,b2) + c2*pow2(s) + d2*pow3(s) + e2*pow4(s); |
---|
2292 | else if (s < s3) |
---|
2293 | return a3 + b3*s + c3*pow2(s) + d3*pow3(s) + e3*pow4(s); |
---|
2294 | else if (s < s4) |
---|
2295 | return a4 + b4*s; |
---|
2296 | else if (s < s5) |
---|
2297 | return a5 + b5*s; |
---|
2298 | else |
---|
2299 | return 0.0; |
---|
2300 | |
---|
2301 | } |
---|
2302 | |
---|
2303 | //========================================================================== |
---|
2304 | |
---|
2305 | // Tau decay matrix element for tau decay into five pions using the model given |
---|
2306 | // in hep-ph/0602162v1. |
---|
2307 | |
---|
2308 | // Ja(q,q1,q2,q3,q4,q5): current through rho and omega resonances |
---|
2309 | // Jb(q,q1,q2,q3,q4,q5): current through a1 and sigma resonances |
---|
2310 | // breitWigner(s,M,G): s-wave Breit-Wigner assuming massless products |
---|
2311 | // a1M: on-shell mass of the a1 resonance |
---|
2312 | // a1G: width of the a1 resonance |
---|
2313 | // rhoX: mass and width of the rho resonance |
---|
2314 | // omegaX: mass, width, and weight of the omega resonance |
---|
2315 | // sigmaX: mass, width, and weight of the sigma resonance |
---|
2316 | |
---|
2317 | //-------------------------------------------------------------------------- |
---|
2318 | |
---|
2319 | // Initialize constants for the helicity matrix element. |
---|
2320 | |
---|
2321 | void HMETau2FivePions::initConstants() { |
---|
2322 | |
---|
2323 | // pi-, pi-, pi+, pi+, pi- decay. |
---|
2324 | if (abs(pID[2]) == 211 && abs(pID[3]) == 211 && abs(pID[4]) == 211 && |
---|
2325 | abs(pID[5]) == 211 && abs(pID[6]) == 211) |
---|
2326 | DECAYWEIGHTMAX = 4e4; |
---|
2327 | // pi+, pi-, pi0, pi-, pi0 decay. |
---|
2328 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 111 && abs(pID[4]) == 211 && |
---|
2329 | abs(pID[5]) == 211 && abs(pID[6]) == 211) |
---|
2330 | DECAYWEIGHTMAX = 1e7; |
---|
2331 | // pi0, pi0, pi-, pi0, pi0 decay. |
---|
2332 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 111 && abs(pID[4]) == 111 && |
---|
2333 | abs(pID[5]) == 111 && abs(pID[6]) == 211) |
---|
2334 | DECAYWEIGHTMAX = 1e5; |
---|
2335 | |
---|
2336 | // Set resonances. |
---|
2337 | a1M = 1.260; a1G = 0.400; |
---|
2338 | rhoM = 0.776; rhoG = 0.150; |
---|
2339 | omegaM = 0.782; omegaG = 0.0085; omegaW = 11.5; |
---|
2340 | sigmaM = 0.800; sigmaG = 0.600; sigmaW = 1; |
---|
2341 | |
---|
2342 | } |
---|
2343 | |
---|
2344 | //-------------------------------------------------------------------------- |
---|
2345 | |
---|
2346 | // Initialize the hadronic current for the helicity matrix element. |
---|
2347 | |
---|
2348 | void HMETau2FivePions::initHadronicCurrent(vector<HelicityParticle>& p) { |
---|
2349 | |
---|
2350 | vector< Wave4 > u2; |
---|
2351 | |
---|
2352 | Wave4 q(p[2].p() + p[3].p() + p[4].p() + p[5].p() + p[6].p()); |
---|
2353 | Wave4 q2(p[2].p()), q3(p[3].p()), q4(p[4].p()), q5(p[5].p()), q6(p[6].p()); |
---|
2354 | // pi-, pi-, pi+, pi+, pi- decay. |
---|
2355 | if (abs(pID[2]) == 211 && abs(pID[3]) == 211 && abs(pID[4]) == 211 && |
---|
2356 | abs(pID[5]) == 211 && abs(pID[6]) == 211) |
---|
2357 | u2.push_back(Jb(q, q2, q3, q5, q6, q4) + Jb(q, q4, q3, q5, q6, q2) |
---|
2358 | + Jb(q, q2, q4, q5, q6, q3) + Jb(q, q2, q3, q6, q5, q4) |
---|
2359 | + Jb(q, q4, q3, q6, q5, q2) + Jb(q, q2, q4, q6, q5, q3)); |
---|
2360 | // pi+, pi-, pi0, pi-, pi0 decay. |
---|
2361 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 111 && abs(pID[4]) == 211 && |
---|
2362 | abs(pID[5]) == 211 && abs(pID[6]) == 211) |
---|
2363 | u2.push_back(Ja(q, q6, q4, q2, q5, q3) + Ja(q, q6, q5, q2, q4, q3) |
---|
2364 | + Ja(q, q6, q4, q3, q5, q2) + Ja(q, q6, q5, q3, q4, q2) |
---|
2365 | + Jb(q, q4, q5, q6, q2, q3) + Jb(q, q2, q3, q4, q6, q5) |
---|
2366 | + Jb(q, q2, q3, q5, q6, q4)); |
---|
2367 | // pi0, pi0, pi-, pi0, pi0 decay. |
---|
2368 | else if (abs(pID[2]) == 111 && abs(pID[3]) == 111 && abs(pID[4]) == 111 && |
---|
2369 | abs(pID[5]) == 111 && abs(pID[6]) == 211) |
---|
2370 | u2.push_back(Jb(q, q2, q3, q6, q4, q5) + Jb(q, q5, q3, q6, q4, q2) |
---|
2371 | + Jb(q, q3, q4, q6, q2, q5) + Jb(q, q2, q4, q6, q3, q5) |
---|
2372 | + Jb(q, q2, q5, q6, q4, q3) + Jb(q, q4, q5, q6, q2, q3)); |
---|
2373 | |
---|
2374 | u.push_back(u2); |
---|
2375 | |
---|
2376 | } |
---|
2377 | |
---|
2378 | //-------------------------------------------------------------------------- |
---|
2379 | |
---|
2380 | // Return the omega-rho hadronic current. |
---|
2381 | |
---|
2382 | Wave4 HMETau2FivePions::Ja(Wave4 &q, Wave4 &q1, Wave4 &q2, |
---|
2383 | Wave4 &q3, Wave4 &q4, Wave4 &q5) { |
---|
2384 | |
---|
2385 | Wave4 j = epsilon(q1, q2, q3); |
---|
2386 | return omegaW * (breitWigner(m2(q), a1M, a1G) |
---|
2387 | * breitWigner(m2(q1 + q2 + q3), omegaM, omegaG) |
---|
2388 | * breitWigner(m2(q4 + q5), rhoM, rhoG) |
---|
2389 | * epsilon(q4 - q5, j, q) |
---|
2390 | * (breitWigner(m2(q2 + q3), rhoM, rhoG) |
---|
2391 | + breitWigner(m2(q1 + q3), rhoM, rhoG) |
---|
2392 | + breitWigner(m2(q1 + q2), rhoM, rhoG))); |
---|
2393 | |
---|
2394 | } |
---|
2395 | |
---|
2396 | //-------------------------------------------------------------------------- |
---|
2397 | |
---|
2398 | // Return the a1-sigma hadronic current. |
---|
2399 | |
---|
2400 | Wave4 HMETau2FivePions::Jb(Wave4 &q, Wave4 &q1, Wave4 &q2, |
---|
2401 | Wave4 &q3, Wave4 &q4, Wave4 &q5) { |
---|
2402 | |
---|
2403 | double s = m2(q); |
---|
2404 | Wave4 a1Q = q1 + q2 + q3; |
---|
2405 | double a1S = m2(a1Q); |
---|
2406 | Wave4 j = (m2(q2, q1 - q3) / a1S * a1Q - q1 + q3) |
---|
2407 | * breitWigner(m2(q1 + q3), rhoM, rhoG) |
---|
2408 | + (m2(q1, q2 - q3) / a1S * a1Q - q2 + q3) |
---|
2409 | * breitWigner(m2(q2 + q3), rhoM, rhoG); |
---|
2410 | j = (j * gamma[4] * q / s) * q - j; |
---|
2411 | return sigmaW * (breitWigner(s, a1M, a1G) * breitWigner(a1S, a1M, a1G) |
---|
2412 | * breitWigner(m2(q4 + q5), sigmaM, sigmaG) * j); |
---|
2413 | |
---|
2414 | } |
---|
2415 | |
---|
2416 | complex HMETau2FivePions::breitWigner(double s, double M, double G) { |
---|
2417 | |
---|
2418 | return M * M / (M * M - s - complex(0, 1) * M * G); |
---|
2419 | |
---|
2420 | } |
---|
2421 | |
---|
2422 | //========================================================================== |
---|
2423 | |
---|
2424 | } // end namespace Pythia8 |
---|