BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
D0Topippim2pi0 Class Reference

#include <D0Topippim2pi0.h>

Public Member Functions

 D0Topippim2pi0 ()
 
virtual ~D0Topippim2pi0 ()
 
void init ()
 
complex< double > Amp (vector< double > Pip, vector< double > Pim, vector< double > Pi01, vector< double > Pi02)
 

Detailed Description

Definition at line 9 of file D0Topippim2pi0.h.

Constructor & Destructor Documentation

◆ D0Topippim2pi0()

D0Topippim2pi0::D0Topippim2pi0 ( )
inline

Definition at line 13 of file D0Topippim2pi0.h.

13{}

◆ ~D0Topippim2pi0()

D0Topippim2pi0::~D0Topippim2pi0 ( )
virtual

Definition at line 11 of file D0Topippim2pi0.cxx.

11{}

Member Function Documentation

◆ Amp()

complex< double > D0Topippim2pi0::Amp ( vector< double >  Pip,
vector< double >  Pim,
vector< double >  Pi01,
vector< double >  Pi02 
)

Definition at line 1384 of file D0Topippim2pi0.cxx.

1385{
1386
1387 vector<double> PipPim; PipPim.clear();
1388 vector<double> PipPi01; PipPi01.clear();
1389 vector<double> PipPi02; PipPi02.clear();
1390 vector<double> PimPi01; PimPi01.clear();
1391 vector<double> PimPi02; PimPi02.clear();
1392 vector<double> Pi01Pi02; Pi01Pi02.clear();
1393
1394 PipPim = sum_tensor(Pip, Pim);
1395 PipPi01 = sum_tensor(Pip, Pi01);
1396 PipPi02 = sum_tensor(Pip, Pi02);
1397 PimPi01 = sum_tensor(Pim, Pi01);
1398 PimPi02 = sum_tensor(Pim, Pi02);
1399 Pi01Pi02 = sum_tensor(Pi01, Pi02);
1400
1401 vector<double> PipPimPi01; PipPimPi01.clear();
1402 vector<double> PipPimPi02; PipPimPi02.clear();
1403 vector<double> PipPi01Pi02; PipPi01Pi02.clear();
1404 vector<double> PimPi01Pi02; PimPi01Pi02.clear();
1405
1406 PipPimPi01 = sum_tensor(PipPim, Pi01);
1407 PipPimPi02 = sum_tensor(PipPim, Pi02);
1408 PipPi01Pi02 = sum_tensor(PipPi01, Pi02);
1409 PimPi01Pi02 = sum_tensor(PimPi01, Pi02);
1410
1411 vector<double> D0; D0.clear();
1412 D0 = sum_tensor(PipPimPi01, Pi02);
1413
1414 double M2_PipPim = contract_11_0(PipPim, PipPim);
1415 double M2_PipPi01 = contract_11_0(PipPi01, PipPi01);
1416 double M2_PipPi02 = contract_11_0(PipPi02, PipPi02);
1417 double M2_PimPi01 = contract_11_0(PimPi01, PimPi01);
1418 double M2_PimPi02 = contract_11_0(PimPi02, PimPi02);
1419 double M2_Pi01Pi02 = contract_11_0(Pi01Pi02, Pi01Pi02);
1420
1421 double M2_PipPimPi01 = contract_11_0(PipPimPi01, PipPimPi01);
1422 double M2_PipPimPi02 = contract_11_0(PipPimPi02, PipPimPi02);
1423 double M2_PipPi01Pi02 = contract_11_0(PipPi01Pi02, PipPi01Pi02);
1424 double M2_PimPi01Pi02 = contract_11_0(PimPi01Pi02, PimPi01Pi02);
1425 double M2_D0 = contract_11_0(D0, D0);
1426
1427 complex<double> GS_rho770_pm = GS(M2_PipPim, m0_rho7700, w0_rho7700, m2_Pi, m2_Pi, rRes, 1);
1428 complex<double> GS_rho770_p1 = GS(M2_PipPi01, m0_rho770p, w0_rho770p, m2_Pi, m2_Pi0, rRes, 1);
1429 complex<double> GS_rho770_p2 = GS(M2_PipPi02, m0_rho770p, w0_rho770p, m2_Pi, m2_Pi0, rRes, 1);
1430 complex<double> GS_rho770_m1 = GS(M2_PimPi01, m0_rho770p, w0_rho770p, m2_Pi, m2_Pi0, rRes, 1);
1431 complex<double> GS_rho770_m2 = GS(M2_PimPi02, m0_rho770p, w0_rho770p, m2_Pi, m2_Pi0, rRes, 1);
1432
1433 complex<double> RBW_f21270_pm = RBW(M2_PipPim, m0_f21270, w0_f21270, m2_Pi, m2_Pi, rRes, 2);
1434 complex<double> RBW_f21270_00 = RBW(M2_Pi01Pi02, m0_f21270, w0_f21270, m2_Pi0, m2_Pi0, rRes, 2);
1435
1436 complex<double> PiPiS_pm_0 = Fvector(M2_PipPim, s0_prod, 0);
1437 complex<double> PiPiS_00_0 = Fvector(M2_Pi01Pi02, s0_prod, 0);
1438
1439 complex<double> PiPiS_pm_1 = Fvector(M2_PipPim, s0_prod, 1);
1440 complex<double> PiPiS_00_1 = Fvector(M2_Pi01Pi02, s0_prod, 1);
1441
1442 complex<double> PiPiS_pm_5 = Fvector(M2_PipPim, s0_prod, 5);
1443 complex<double> PiPiS_00_5 = Fvector(M2_Pi01Pi02, s0_prod, 5);
1444
1445 complex<double> PiPiS_pm_6 = Fvector(M2_PipPim, s0_prod, 6);
1446 complex<double> PiPiS_00_6 = Fvector(M2_Pi01Pi02, s0_prod, 6);
1447
1448 complex<double> RBW_a11260_p = RBWa1260(M2_PipPi01Pi02, m0_a11260, g1_a11260, g2_a11260);
1449 complex<double> RBW_a11260_m = RBWa1260(M2_PimPi01Pi02, m0_a11260, g1_a11260, g2_a11260);
1450 complex<double> RBW_a11260_01 = RBWa1260(M2_PipPimPi01, m0_a11260, g1_a11260, g2_a11260);
1451 complex<double> RBW_a11260_02 = RBWa1260(M2_PipPimPi02, m0_a11260, g1_a11260, g2_a11260);
1452
1453 complex<double> RBW_a11420_p = RBW(M2_PipPi01Pi02, m0_a11420, w0_a11420,-1,-1,-1,-1);
1454 complex<double> RBW_a11420_m = RBW(M2_PimPi01Pi02, m0_a11420, w0_a11420,-1,-1,-1,-1);
1455 complex<double> RBW_a11420_01 = RBW(M2_PipPimPi01, m0_a11420, w0_a11420,-1,-1,-1,-1);
1456 complex<double> RBW_a11420_02 = RBW(M2_PipPimPi02, m0_a11420, w0_a11420,-1,-1,-1,-1);
1457
1458 complex<double> RBW_omega_01 = RBW(M2_PipPimPi01, m0_omega, w0_omega,-1,-1,-1,-1);
1459 complex<double> RBW_omega_02 = RBW(M2_PipPimPi02, m0_omega, w0_omega,-1,-1,-1,-1);
1460
1461 complex<double> RBW_phi_01 = RBW(M2_PipPimPi01, m0_phi, w0_phi,-1,-1,-1,-1);
1462 complex<double> RBW_phi_02 = RBW(M2_PipPimPi02, m0_phi, w0_phi,-1,-1,-1,-1);
1463
1464 complex<double> RBW_a21320_p = RBW(M2_PipPi01Pi02, m0_a21320, w0_a21320,-1,-1,-1,-1);
1465 complex<double> RBW_a21320_m = RBW(M2_PimPi01Pi02, m0_a21320, w0_a21320,-1,-1,-1,-1);
1466
1467 complex<double> RBW_pi1300_p = RBWpi1300(M2_PipPi01Pi02, m0_pi1300, w0_pi1300);
1468 complex<double> RBW_pi1300_m = RBWpi1300(M2_PimPi01Pi02, m0_pi1300, w0_pi1300);
1469 complex<double> RBW_pi1300_01 = RBWpi1300(M2_PipPimPi01, m0_pi1300, w0_pi1300);
1470 complex<double> RBW_pi1300_02 = RBWpi1300(M2_PipPimPi02, m0_pi1300, w0_pi1300);
1471
1472 complex<double> RBW_h11170_01 = RBW(M2_PipPimPi01, m0_h11170, w0_h11170,-1,-1,-1,-1);
1473 complex<double> RBW_h11170_02 = RBW(M2_PipPimPi02, m0_h11170, w0_h11170,-1,-1,-1,-1);
1474
1475 complex<double> RBW_pi21670_01 = RBW(M2_PipPimPi01, m0_pi21670, w0_pi21670,-1,-1,-1,-1);
1476 complex<double> RBW_pi21670_02 = RBW(M2_PipPimPi02, m0_pi21670, w0_pi21670,-1,-1,-1,-1);
1477
1478 // D->XX Projection
1479 vector<double> Proj1_3p; Proj1_3p.clear();
1480 vector<double> Proj1_3m; Proj1_3m.clear();
1481 vector<double> Proj1_3z1; Proj1_3z1.clear();
1482 vector<double> Proj1_3z2; Proj1_3z2.clear();
1483
1484 Proj1_3p = ProjectionTensors(PipPi01Pi02,1);
1485 Proj1_3m = ProjectionTensors(PimPi01Pi02,1);
1486 Proj1_3z1 = ProjectionTensors(PipPimPi01,1);
1487 Proj1_3z2 = ProjectionTensors(PipPimPi02,1);
1488
1489 vector<double> Proj2_3p; Proj2_3p.clear();
1490 vector<double> Proj2_3m; Proj2_3m.clear();
1491 vector<double> Proj2_3z1; Proj2_3z1.clear();
1492 vector<double> Proj2_3z2; Proj2_3z2.clear();
1493
1494 Proj2_3p = ProjectionTensors(PipPi01Pi02,2);
1495 Proj2_3m = ProjectionTensors(PimPi01Pi02,2);
1496 Proj2_3z1 = ProjectionTensors(PipPimPi01,2);
1497 Proj2_3z2 = ProjectionTensors(PipPimPi02,2);
1498
1499 // X->PP Orbital
1500 vector<double> T1_PipPim; T1_PipPim.clear();
1501 vector<double> T1_PipPi01; T1_PipPi01.clear();
1502 vector<double> T1_PipPi02; T1_PipPi02.clear();
1503 vector<double> T1_PimPi01; T1_PimPi01.clear();
1504 vector<double> T1_PimPi02; T1_PimPi02.clear();
1505 vector<double> T1_Pi01Pi02; T1_Pi01Pi02.clear();
1506
1507 T1_PipPim = OrbitalTensors(PipPim, Pip, Pim, rRes, 1);
1508 T1_PipPi01 = OrbitalTensors(PipPi01, Pip, Pi01, rRes, 1);
1509 T1_PipPi02 = OrbitalTensors(PipPi02, Pip, Pi02, rRes, 1);
1510 T1_PimPi01 = OrbitalTensors(PimPi01, Pim, Pi01, rRes, 1);
1511 T1_PimPi02 = OrbitalTensors(PimPi02, Pim, Pi02, rRes, 1);
1512 T1_Pi01Pi02 = OrbitalTensors(Pi01Pi02, Pi01, Pi02, rRes, 1);
1513
1514 vector<double> T2_PipPim; T2_PipPim.clear();
1515 vector<double> T2_Pi01Pi02; T2_Pi01Pi02.clear();
1516
1517 T2_PipPim = OrbitalTensors(PipPim, Pip, Pim, rRes, 2);
1518 T2_Pi01Pi02 = OrbitalTensors(Pi01Pi02, Pi01, Pi02, rRes, 2);
1519
1520 // X->YP Orbital
1521 vector<double> T1_PipPimPi01; T1_PipPimPi01.clear();
1522 vector<double> T1_PipPimPi02; T1_PipPimPi02.clear();
1523 vector<double> T1_PipPi01Pi02; T1_PipPi01Pi02.clear();
1524 vector<double> T1_PipPi02Pi01; T1_PipPi02Pi01.clear();
1525 vector<double> T1_PimPi01Pi02; T1_PimPi01Pi02.clear();
1526 vector<double> T1_PimPi02Pi01; T1_PimPi02Pi01.clear();
1527 vector<double> T1_PipPi01Pim; T1_PipPi01Pim.clear();
1528 vector<double> T1_PipPi02Pim; T1_PipPi02Pim.clear();
1529 vector<double> T1_PimPi01Pip; T1_PimPi01Pip.clear();
1530 vector<double> T1_PimPi02Pip; T1_PimPi02Pip.clear();
1531 vector<double> T1_Pi01Pi02Pip; T1_Pi01Pi02Pip.clear();
1532 vector<double> T1_Pi01Pi02Pim; T1_Pi01Pi02Pim.clear();
1533
1534 T1_PipPimPi01 = OrbitalTensors(PipPimPi01, PipPim, Pi01, rRes, 1);
1535 T1_PipPimPi02 = OrbitalTensors(PipPimPi02, PipPim, Pi02, rRes, 1);
1536 T1_PipPi01Pi02 = OrbitalTensors(PipPi01Pi02, PipPi01, Pi02, rRes, 1);
1537 T1_PipPi02Pi01 = OrbitalTensors(PipPi01Pi02, PipPi02, Pi01, rRes, 1);
1538 T1_PimPi01Pi02 = OrbitalTensors(PimPi01Pi02, PimPi01, Pi02, rRes, 1);
1539 T1_PimPi02Pi01 = OrbitalTensors(PimPi01Pi02, PimPi02, Pi01, rRes, 1);
1540 T1_PipPi01Pim = OrbitalTensors(PipPimPi01, PipPi01, Pim, rRes, 1);
1541 T1_PipPi02Pim = OrbitalTensors(PipPimPi02, PipPi02, Pim, rRes, 1);
1542 T1_PimPi01Pip = OrbitalTensors(PipPimPi01, PimPi01, Pip, rRes, 1);
1543 T1_PimPi02Pip = OrbitalTensors(PipPimPi02, PimPi02, Pip, rRes, 1);
1544 T1_Pi01Pi02Pip = OrbitalTensors(PipPi01Pi02, Pi01Pi02, Pip, rRes, 1);
1545 T1_Pi01Pi02Pim = OrbitalTensors(PimPi01Pi02, Pi01Pi02, Pim, rRes, 1);
1546
1547 vector<double> T2_PipPimPi01; T2_PipPimPi01.clear();
1548 vector<double> T2_PipPimPi02; T2_PipPimPi02.clear();
1549 vector<double> T2_PipPi01Pi02; T2_PipPi01Pi02.clear();
1550 vector<double> T2_PipPi02Pi01; T2_PipPi02Pi01.clear();
1551 vector<double> T2_PimPi01Pi02; T2_PimPi01Pi02.clear();
1552 vector<double> T2_PimPi02Pi01; T2_PimPi02Pi01.clear();
1553 vector<double> T2_PipPi01Pim; T2_PipPi01Pim.clear();
1554 vector<double> T2_PipPi02Pim; T2_PipPi02Pim.clear();
1555 vector<double> T2_PimPi01Pip; T2_PimPi01Pip.clear();
1556 vector<double> T2_PimPi02Pip; T2_PimPi02Pip.clear();
1557 vector<double> T2_Pi01Pi02Pip; T2_Pi01Pi02Pip.clear();
1558 vector<double> T2_Pi01Pi02Pim; T2_Pi01Pi02Pim.clear();
1559
1560 T2_PipPimPi01 = OrbitalTensors(PipPimPi01, PipPim, Pi01, rRes, 2);
1561 T2_PipPimPi02 = OrbitalTensors(PipPimPi02, PipPim, Pi02, rRes, 2);
1562 T2_PipPi01Pi02 = OrbitalTensors(PipPi01Pi02, PipPi01, Pi02, rRes, 2);
1563 T2_PipPi02Pi01 = OrbitalTensors(PipPi01Pi02, PipPi02, Pi01, rRes, 2);
1564 T2_PimPi01Pi02 = OrbitalTensors(PimPi01Pi02, PimPi01, Pi02, rRes, 2);
1565 T2_PimPi02Pi01 = OrbitalTensors(PimPi01Pi02, PimPi02, Pi01, rRes, 2);
1566 T2_PipPi01Pim = OrbitalTensors(PipPimPi01, PipPi01, Pim, rRes, 2);
1567 T2_PipPi02Pim = OrbitalTensors(PipPimPi02, PipPi02, Pim, rRes, 2);
1568 T2_PimPi01Pip = OrbitalTensors(PipPimPi01, PimPi01, Pip, rRes, 2);
1569 T2_PimPi02Pip = OrbitalTensors(PipPimPi02, PimPi02, Pip, rRes, 2);
1570 T2_Pi01Pi02Pip = OrbitalTensors(PipPi01Pi02, Pi01Pi02, Pip, rRes, 2);
1571 T2_Pi01Pi02Pim = OrbitalTensors(PimPi01Pi02, Pi01Pi02, Pim, rRes, 2);
1572
1573 // D->XX Orbital
1574 vector<double> T1_2pm12; T1_2pm12.clear();
1575 vector<double> T1_2p1m2; T1_2p1m2.clear();
1576 vector<double> T1_2p2m1; T1_2p2m1.clear();
1577
1578 T1_2pm12 = OrbitalTensors(D0, PipPim, Pi01Pi02, rD, 1);
1579 T1_2p1m2 = OrbitalTensors(D0, PipPi01, PimPi02, rD, 1);
1580 T1_2p2m1 = OrbitalTensors(D0, PipPi02, PimPi01, rD, 1);
1581
1582 vector<double> T2_2pm12; T2_2pm12.clear();
1583 vector<double> T2_2p1m2; T2_2p1m2.clear();
1584 vector<double> T2_2p2m1; T2_2p2m1.clear();
1585
1586 T2_2pm12 = OrbitalTensors(D0, PipPim, Pi01Pi02, rD, 2);
1587 T2_2p1m2 = OrbitalTensors(D0, PipPi01, PimPi02, rD, 2);
1588 T2_2p2m1 = OrbitalTensors(D0, PipPi02, PimPi01, rD, 2);
1589
1590 // D->XP Orbital
1591 vector<double> T1_3pm; T1_3pm.clear();
1592 vector<double> T1_3mp; T1_3mp.clear();
1593 vector<double> T1_3z12; T1_3z12.clear();
1594 vector<double> T1_3z21; T1_3z21.clear();
1595
1596 T1_3pm = OrbitalTensors(D0, PipPi01Pi02, Pim, rD, 1);
1597 T1_3mp = OrbitalTensors(D0, PimPi01Pi02, Pip, rD, 1);
1598 T1_3z12 = OrbitalTensors(D0, PipPimPi01, Pi02, rD, 1);
1599 T1_3z21 = OrbitalTensors(D0, PipPimPi02, Pi01, rD, 1);
1600
1601 vector<double> T2_3pm; T2_3pm.clear();
1602 vector<double> T2_3mp; T2_3mp.clear();
1603 vector<double> T2_3z12; T2_3z12.clear();
1604 vector<double> T2_3z21; T2_3z21.clear();
1605
1606 T2_3pm = OrbitalTensors(D0, PipPi01Pi02, Pim, rD, 2);
1607 T2_3mp = OrbitalTensors(D0, PimPi01Pi02, Pip, rD, 2);
1608 T2_3z12 = OrbitalTensors(D0, PipPimPi01, Pi02, rD, 2);
1609 T2_3z21 = OrbitalTensors(D0, PipPimPi02, Pi01, rD, 2);
1610
1611 complex<double> amplitude(0,0);
1612
1613 // D0 -> a1(1260)+ {rho(770)+ pi0[S]} pi-
1614 double SF_Ap_S_Vp1P = contract_11_0(contract_21_1(Proj1_3p, T1_PipPi01), T1_3pm);
1615 double SF_Ap_S_Vp2P = contract_11_0(contract_21_1(Proj1_3p, T1_PipPi02), T1_3pm);
1616
1617 amplitude += fitpara[0]*(SF_Ap_S_Vp1P*RBW_a11260_p*GS_rho770_p1 + SF_Ap_S_Vp2P*RBW_a11260_p*GS_rho770_p2);
1618
1619 // D0 -> a1(1260)+ {rho(770)+ pi0[D]} pi-
1620 double SF_Ap_D_Vp1P = contract_11_0(contract_21_1(T2_PipPi01Pi02, T1_PipPi01), T1_3pm);
1621 double SF_Ap_D_Vp2P = contract_11_0(contract_21_1(T2_PipPi02Pi01, T1_PipPi02), T1_3pm);
1622
1623// cout<<SF_Ap_D_VP_1<<","<<SF_Ap_D_VP_2<<","<<SF_Ap_D_VP_3<<","<<SF_Ap_D_VP_4<<","<<endl;
1624// cout<<"-------"<<endl;
1625 amplitude += fitpara[1]*(SF_Ap_D_Vp1P*RBW_a11260_p*GS_rho770_p1 + SF_Ap_D_Vp2P*RBW_a11260_p*GS_rho770_p2);
1626
1627 // D0 -> a1(1260)+ {f2(1270) pi+ [P]} pi-
1628 double SF_Ap_P_TP = contract_11_0(contract_21_1(contract_42_2(Proj2_3p, T2_Pi01Pi02), T1_Pi01Pi02Pip), T1_3pm);
1629
1630 amplitude += fitpara[2]*(SF_Ap_P_TP*RBW_a11260_p*RBW_f21270_00);
1631
1632 // D0 -> a1(1260)+ {f0 pi+ [P]} pi-
1633 double SF_Ap_P_SP = contract_11_0(T1_3pm, T1_Pi01Pi02Pip);
1634
1635 amplitude += fitpara[3]*(SF_Ap_P_SP*RBW_a11260_p*PiPiS_00_0);
1636 amplitude += fitpara[4]*(SF_Ap_P_SP*RBW_a11260_p*PiPiS_00_1);
1637 amplitude += fitpara[5]*(SF_Ap_P_SP*RBW_a11260_p*PiPiS_00_5);
1638
1639 // D0 -> a1(1260)- { rho(770)- pi0 [S]} pi+
1640 double SF_Am_S_Vm1P = contract_11_0(contract_21_1(Proj1_3m, T1_PimPi01), T1_3mp);
1641 double SF_Am_S_Vm2P = contract_11_0(contract_21_1(Proj1_3m, T1_PimPi02), T1_3mp);
1642
1643 amplitude += fitpara[6]*fitpara[0]*(SF_Am_S_Vm1P*RBW_a11260_m*GS_rho770_m1 + SF_Am_S_Vm2P*RBW_a11260_m*GS_rho770_m2);
1644
1645 // D0 -> a1(1260)- {rho(770)- pi0[D]} pi+
1646 double SF_Am_D_Vm1P = contract_11_0(contract_21_1(T2_PimPi01Pi02, T1_PimPi01), T1_3mp);
1647 double SF_Am_D_Vm2P = contract_11_0(contract_21_1(T2_PimPi02Pi01, T1_PimPi02), T1_3mp);
1648
1649 amplitude += fitpara[6]*fitpara[1]*(SF_Am_D_Vm1P*RBW_a11260_m*GS_rho770_m1 + SF_Am_D_Vm2P*RBW_a11260_m*GS_rho770_m2);
1650
1651 // D0 -> a1(1260)- {f2(1270) pi- [P]} pi+
1652 double SF_Am_P_TP = contract_11_0(contract_21_1(contract_42_2(Proj2_3m, T2_Pi01Pi02), T1_Pi01Pi02Pim), T1_3mp);
1653
1654 amplitude += fitpara[6]*fitpara[2]*(SF_Am_P_TP*RBW_a11260_m*RBW_f21270_00);
1655
1656 // D0 -> a1(1260)- {f0 pi- [P]} pi+
1657 double SF_Am_P_SP = contract_11_0(T1_3mp, T1_Pi01Pi02Pim);
1658
1659 amplitude += fitpara[6]*fitpara[3]*(SF_Am_P_SP*RBW_a11260_m*PiPiS_00_0);
1660 amplitude += fitpara[6]*fitpara[4]*(SF_Am_P_SP*RBW_a11260_m*PiPiS_00_1);
1661 amplitude += fitpara[6]*fitpara[5]*(SF_Am_P_SP*RBW_a11260_m*PiPiS_00_5);
1662
1663 // D -> a1(1260)0 pi0
1664 double SF_A01_S_Vp1P = contract_11_0(contract_21_1(Proj1_3z1, T1_PipPi01), T1_3z12);
1665 double SF_A02_S_Vp2P = contract_11_0(contract_21_1(Proj1_3z2, T1_PipPi02), T1_3z21);
1666 double SF_A01_S_Vm1P = contract_11_0(contract_21_1(Proj1_3z1, T1_PimPi01), T1_3z12);
1667 double SF_A02_S_Vm2P = contract_11_0(contract_21_1(Proj1_3z2, T1_PimPi02), T1_3z21);
1668 double SF_A01_S_VzP = contract_11_0(contract_21_1(Proj1_3z1, T1_PipPim), T1_3z12);
1669 double SF_A02_S_VzP = contract_11_0(contract_21_1(Proj1_3z2, T1_PipPim), T1_3z21);
1670
1671 amplitude += fitpara[7]*fitpara[0]*(SF_A01_S_Vp1P*RBW_a11260_01*GS_rho770_p1 + SF_A02_S_Vp2P*RBW_a11260_02*GS_rho770_p2 + SF_A01_S_Vm1P*RBW_a11260_01*GS_rho770_m1 + SF_A02_S_Vm2P*RBW_a11260_02*GS_rho770_m2);
1672
1673 double SF_A01_D_Vp1P = contract_11_0(contract_21_1(T2_PipPi01Pim, T1_PipPi01), T1_3z12);
1674 double SF_A02_D_Vp2P = contract_11_0(contract_21_1(T2_PipPi02Pim, T1_PipPi02), T1_3z21);
1675 double SF_A01_D_Vm1P = contract_11_0(contract_21_1(T2_PimPi01Pip, T1_PimPi01), T1_3z12);
1676 double SF_A02_D_Vm2P = contract_11_0(contract_21_1(T2_PimPi02Pip, T1_PimPi02), T1_3z21);
1677
1678 amplitude += fitpara[7]*fitpara[1]*(SF_A01_D_Vp1P*RBW_a11260_01*GS_rho770_p1 + SF_A02_D_Vp2P*RBW_a11260_02*GS_rho770_p2 + SF_A01_D_Vm1P*RBW_a11260_01*GS_rho770_m1 + SF_A02_D_Vm2P*RBW_a11260_02*GS_rho770_m2);
1679
1680 double SF_A01_P_TP = contract_11_0(contract_21_1(contract_42_2(Proj2_3z1,T2_PipPim), T1_PipPimPi01), T1_3z12);
1681 double SF_A02_P_TP = contract_11_0(contract_21_1(contract_42_2(Proj2_3z2,T2_PipPim), T1_PipPimPi02), T1_3z21);
1682
1683 amplitude += fitpara[7]*fitpara[2]*(-1.0)*(SF_A01_P_TP*RBW_a11260_01*RBW_f21270_pm + SF_A02_P_TP*RBW_a11260_02*RBW_f21270_pm);
1684
1685 double SF_A01_P_SP = contract_11_0(T1_3z12, T1_PipPimPi01);
1686 double SF_A02_P_SP = contract_11_0(T1_3z21, T1_PipPimPi02);
1687
1688 amplitude += fitpara[7]*fitpara[3]*(-1.0)*(SF_A01_P_SP*RBW_a11260_01*PiPiS_pm_0 + SF_A02_P_SP*RBW_a11260_02*PiPiS_pm_0);
1689 amplitude += fitpara[7]*fitpara[4]*(-1.0)*(SF_A01_P_SP*RBW_a11260_01*PiPiS_pm_1 + SF_A02_P_SP*RBW_a11260_02*PiPiS_pm_1);
1690 amplitude += fitpara[7]*fitpara[5]*(-1.0)*(SF_A01_P_SP*RBW_a11260_01*PiPiS_pm_5 + SF_A02_P_SP*RBW_a11260_02*PiPiS_pm_5);
1691
1692 // D0 -> a1(1420)+ {rho(770)+ pi0[S]} pi-
1693// amplitude += fitpara[8]*(SF_Ap_S_Vp1P*RBW_a11420_p*GS_rho770_p1 + SF_Ap_S_Vp2P*RBW_a11420_p*GS_rho770_p2);
1694// amplitude += fitpara[9]*(SF_A01_S_Vp1P*RBW_a11420_01*GS_rho770_p1 + SF_A02_S_Vp2P*RBW_a11420_02*GS_rho770_p2 + SF_A01_S_Vm1P*RBW_a11420_01*GS_rho770_m1 + SF_A02_S_Vm2P*RBW_a11420_02*GS_rho770_m2);
1695
1696 // D0 -> a1(1420)+ {f0 pi0+[S]} pi-
1697 amplitude += fitpara[8]*(SF_Ap_P_SP*RBW_a11420_p*PiPiS_00_5);
1698 amplitude += fitpara[9]*(SF_Ap_P_SP*RBW_a11420_p*PiPiS_00_6);
1699
1700 // D0 -> a2(1320)+ {rho+ pi0} pi-
1701 double SF_Tp_D_Vp1P = contract_22_0(contract_22_2(contract_31_2(contract_41_3(epsilon_uvmn, contract_21_1(Proj1_3p, T1_PipPi01)), PipPi01Pi02), contract_42_2(Proj2_3p, T2_3pm)), T2_PipPi01Pi02);
1702 double SF_Tp_D_Vp2P = contract_22_0(contract_22_2(contract_31_2(contract_41_3(epsilon_uvmn, contract_21_1(Proj1_3p, T1_PipPi02)), PipPi01Pi02), contract_42_2(Proj2_3p, T2_3pm)), T2_PipPi02Pi01);
1703
1704 amplitude += fitpara[10]*(SF_Tp_D_Vp1P*GS_rho770_p1*RBW_a21320_p + SF_Tp_D_Vp2P*GS_rho770_p2*RBW_a21320_p);
1705
1706 // D0 -> a2(1320)- {rho- pi0} pi+
1707 double SF_Tm_D_Vm1P = contract_22_0(contract_22_2(contract_31_2(contract_41_3(epsilon_uvmn, contract_21_1(Proj1_3m, T1_PimPi01)), PimPi01Pi02), contract_42_2(Proj2_3m, T2_3mp)), T2_PimPi01Pi02);
1708 double SF_Tm_D_Vm2P = contract_22_0(contract_22_2(contract_31_2(contract_41_3(epsilon_uvmn, contract_21_1(Proj1_3m, T1_PimPi02)), PimPi01Pi02), contract_42_2(Proj2_3m, T2_3mp)), T2_PimPi02Pi01);
1709 amplitude += fitpara[11]*(SF_Tm_D_Vm1P*GS_rho770_m1*RBW_a21320_m + SF_Tm_D_Vm2P*GS_rho770_m2*RBW_a21320_m);
1710
1711 // D0 -> h1(1170)0 {rho pi0} pi0
1712 amplitude += fitpara[12]*(SF_A01_S_Vp1P*RBW_h11170_01*GS_rho770_p1 + SF_A02_S_Vp2P*RBW_h11170_02*GS_rho770_p2 - SF_A01_S_Vm1P*RBW_h11170_01*GS_rho770_m1 - SF_A02_S_Vm2P*RBW_h11170_02*GS_rho770_m2 - SF_A01_S_VzP*RBW_h11170_01*GS_rho770_pm - SF_A02_S_VzP*RBW_h11170_02*GS_rho770_pm);
1713
1714 // D0 -> pi(1300)- {rho(770)- pi0} pi+
1715 double SF_Pm_P_Vm1P = contract_11_0(T1_PimPi01,T1_PimPi01Pi02);
1716 double SF_Pm_P_Vm2P = contract_11_0(T1_PimPi02,T1_PimPi02Pi01);
1717
1718 amplitude += fitpara[13]*(SF_Pm_P_Vm1P*GS_rho770_m1*RBW_pi1300_m + SF_Pm_P_Vm2P*GS_rho770_m2*RBW_pi1300_m);
1719
1720 // D0 -> pi(1300)- {f2 pi-} pi+
1721// double SF_Pm_D_TP = contract_22_0(T2_Pi01Pi02,T2_Pi01Pi02Pim);
1722// amplitude += fitpara[14]*fitpara[13]*(SF_Pm_D_TP*RBW_f21270_00*RBW_pi1300_m);
1723
1724 // D0 -> pi(1300)- {f0 pi-} pi+
1725 amplitude += fitpara[14]*fitpara[13]*(RBW_pi1300_m*PiPiS_00_0);
1726// amplitude += fitpara[15]*fitpara[13]*(RBW_pi1300_m*PiPiS_00_5);
1727 amplitude += fitpara[15]*fitpara[13]*(RBW_pi1300_m*PiPiS_00_6);
1728
1729 // D0 -> pi(1300)+ {rho(770)+ pi0} pi-
1730 double SF_Pp_P_Vp1P = contract_11_0(T1_PipPi01,T1_PipPi01Pi02);
1731 double SF_Pp_P_Vp2P = contract_11_0(T1_PipPi02,T1_PipPi02Pi01);
1732
1733 amplitude += fitpara[16]*(SF_Pp_P_Vp1P*GS_rho770_p1*RBW_pi1300_p + SF_Pp_P_Vp2P*GS_rho770_p2*RBW_pi1300_p);
1734
1735 // D0 -> pi(1300)+ {f2 pi-} pi+
1736// double SF_Pp_D_TP = contract_22_0(T2_Pi01Pi02,T2_Pi01Pi02Pip);
1737// amplitude += fitpara[14]*fitpara[17]*(SF_Pp_D_TP*RBW_f21270_00*RBW_pi1300_p);
1738
1739 // D0 -> pi(1300)+ {f0 pi+} pi-
1740 amplitude += fitpara[14]*fitpara[16]*(RBW_pi1300_p*PiPiS_00_0);
1741// amplitude += fitpara[15]*fitpara[17]*(RBW_pi1300_p*PiPiS_00_5);
1742 amplitude += fitpara[15]*fitpara[16]*(RBW_pi1300_p*PiPiS_00_6);
1743
1744 // D0 -> pi(1300)0 {rho pi} pi0
1745 double SF_P01_P_Vp1P = contract_11_0(T1_PipPi01,T1_PipPi01Pim);
1746 double SF_P02_P_Vp2P = contract_11_0(T1_PipPi02,T1_PipPi02Pim);
1747 double SF_P01_P_Vm1P = contract_11_0(T1_PimPi01,T1_PimPi01Pip);
1748 double SF_P02_P_Vm2P = contract_11_0(T1_PimPi02,T1_PimPi02Pip);
1749
1750 amplitude += fitpara[17]*(SF_P01_P_Vp1P*RBW_pi1300_01*GS_rho770_p1 + SF_P02_P_Vp2P*RBW_pi1300_02*GS_rho770_p2 + SF_P01_P_Vm1P*RBW_pi1300_01*GS_rho770_m1 + SF_P02_P_Vm2P*RBW_pi1300_02*GS_rho770_m2);
1751
1752
1753// double SF_P01_D_TP = contract_22_0(T2_PipPim,T2_PipPimPi01);
1754// double SF_P02_D_TP = contract_22_0(T2_PipPim,T2_PipPimPi02);
1755// amplitude += fitpara[14]*fitpara[18]*(-1.0)*(SF_P01_D_TP*RBW_f21270_pm*RBW_pi1300_01 + SF_P02_D_TP*RBW_f21270_pm*RBW_pi1300_02);
1756
1757 amplitude += fitpara[14]*fitpara[17]*(-1.0)*(RBW_pi1300_01*PiPiS_pm_0 + RBW_pi1300_02*PiPiS_pm_0);
1758// amplitude += fitpara[15]*fitpara[18]*(-1.0)*(RBW_pi1300_01*PiPiS_pm_5 + RBW_pi1300_02*PiPiS_pm_5);
1759 amplitude += fitpara[15]*fitpara[17]*(-1.0)*(RBW_pi1300_01*PiPiS_pm_6 + RBW_pi1300_02*PiPiS_pm_6);
1760
1761 // D0 -> rho+ rho- [S]
1762 double SF_Vp1Vm2_S = contract_11_0(T1_PipPi01, T1_PimPi02);
1763 double SF_Vp2Vm1_S = contract_11_0(T1_PipPi02, T1_PimPi01);
1764
1765 amplitude += fitpara[18]*(SF_Vp1Vm2_S*GS_rho770_p1*GS_rho770_m2 + SF_Vp2Vm1_S*GS_rho770_p2*GS_rho770_m1);
1766
1767 // D0 -> rho+ rho- [P]
1768 double SF_Vp1Vm2_P = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, T1_PipPi01),T1_PimPi02),T1_2p1m2), D0);
1769 double SF_Vp2Vm1_P = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, T1_PipPi02),T1_PimPi01),T1_2p2m1), D0);
1770
1771 amplitude += fitpara[19]*(SF_Vp1Vm2_P*GS_rho770_p1*GS_rho770_m2 + SF_Vp2Vm1_P*GS_rho770_p2*GS_rho770_m1);
1772
1773 // D0 -> rho+ rho- [D]
1774 double SF_Vp1Vm2_D = contract_11_0(contract_21_1(T2_2p1m2,T1_PipPi01), T1_PimPi02);
1775 double SF_Vp2Vm1_D = contract_11_0(contract_21_1(T2_2p2m1,T1_PipPi02), T1_PimPi01);
1776 amplitude += fitpara[20]*(SF_Vp1Vm2_D*GS_rho770_p1*GS_rho770_m2 + SF_Vp2Vm1_D*GS_rho770_p2*GS_rho770_m1);
1777
1778 // D0 -> rho0 (Pi0 Pi0)S
1779 double SF_VpmS12_P = contract_11_0(T1_PipPim,T1_2pm12);
1780
1781 amplitude += fitpara[21]*(SF_VpmS12_P*GS_rho770_pm*PiPiS_00_0);
1782 amplitude += fitpara[22]*(SF_VpmS12_P*GS_rho770_pm*PiPiS_00_5);
1783 amplitude += fitpara[23]*(SF_VpmS12_P*GS_rho770_pm*PiPiS_00_6);
1784
1785 // D0 -> f0f0
1786 //00
1787 amplitude += fitpara[24]*(PiPiS_pm_0*PiPiS_00_0 + PiPiS_00_0*PiPiS_pm_0);
1788 amplitude += fitpara[25]*(PiPiS_pm_0*PiPiS_00_1 + PiPiS_00_0*PiPiS_pm_1);
1789 amplitude += fitpara[26]*(PiPiS_pm_1*PiPiS_00_1 + PiPiS_00_1*PiPiS_pm_1);
1790 amplitude += fitpara[27]*(PiPiS_pm_1*PiPiS_00_5 + PiPiS_00_1*PiPiS_pm_5);
1791 amplitude += fitpara[28]*(PiPiS_pm_5*PiPiS_00_5 + PiPiS_00_5*PiPiS_pm_5);
1792 amplitude += fitpara[29]*(PiPiS_pm_5*PiPiS_00_6 + PiPiS_00_5*PiPiS_pm_6);
1793
1794 // D0 -> f2(1270) f0
1795 double SF_TpmS00_D = contract_22_0(T2_PipPim, T2_2pm12);
1796 double SF_T00Spm_D = contract_22_0(T2_Pi01Pi02, T2_2pm12);
1797
1798 amplitude += fitpara[30]*(SF_TpmS00_D*RBW_f21270_pm*PiPiS_00_5 + SF_T00Spm_D*RBW_f21270_00*PiPiS_pm_5);
1799 amplitude += fitpara[31]*(SF_TpmS00_D*RBW_f21270_pm*PiPiS_00_6 + SF_T00Spm_D*RBW_f21270_00*PiPiS_pm_6);
1800
1801 // D0 -> pi2(1670)0[f2(1270) pi0] pi0
1802 double SF_PT01_S_TP = contract_22_0(contract_42_2(Proj2_3z1, T2_PipPim), T2_3z12);
1803 double SF_PT02_S_TP = contract_22_0(contract_42_2(Proj2_3z2, T2_PipPim), T2_3z21);
1804
1805 amplitude += fitpara[32]*(-1.0)*(SF_PT01_S_TP*RBW_f21270_pm*RBW_pi21670_01 + SF_PT02_S_TP*RBW_f21270_pm*RBW_pi21670_02);
1806
1807 // D -> 1--(rho pi) pi0
1808 double SF_V1_Vz = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, PipPimPi01), T1_PipPimPi01), T1_PipPim), contract_21_1(Proj1_3z1, T1_3z12));
1809 double SF_V1_Vp1 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, PipPimPi01), T1_PipPi01Pim), T1_PipPi01), contract_21_1(Proj1_3z1, T1_3z12));
1810 double SF_V1_Vm1 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, PipPimPi01), T1_PimPi01Pip), T1_PimPi01), contract_21_1(Proj1_3z1, T1_3z12));
1811
1812 double SF_V2_Vz = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, PipPimPi02), T1_PipPimPi02), T1_PipPim), contract_21_1(Proj1_3z2, T1_3z21));
1813 double SF_V2_Vp2 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, PipPimPi02), T1_PipPi02Pim), T1_PipPi02), contract_21_1(Proj1_3z2, T1_3z21));
1814 double SF_V1_Vm2 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(epsilon_uvmn, PipPimPi02), T1_PimPi02Pip), T1_PimPi02), contract_21_1(Proj1_3z2, T1_3z21));
1815
1816
1817 // D0 -> omega(rho pi) pi0
1818 amplitude += (-1.0) * fitpara[33]*(SF_V1_Vp1*RBW_omega_01*GS_rho770_p1 - SF_V1_Vz*RBW_omega_01*GS_rho770_pm - SF_V1_Vm1*RBW_omega_01*GS_rho770_m1 + SF_V2_Vp2*RBW_omega_02*GS_rho770_p2 - SF_V2_Vz*RBW_omega_02*GS_rho770_pm - SF_V1_Vm2*RBW_omega_02*GS_rho770_m2);
1819
1820
1821 // D0 -> phi(rho pi) pi0
1822 amplitude += (-1.0) * fitpara[34]*(SF_V1_Vp1*RBW_phi_01*GS_rho770_p1 - SF_V1_Vz*RBW_phi_01*GS_rho770_pm - SF_V1_Vm1*RBW_phi_01*GS_rho770_m1 + SF_V2_Vp2*RBW_phi_02*GS_rho770_p2 - SF_V2_Vz*RBW_phi_02*GS_rho770_pm - SF_V1_Vm2*RBW_phi_02*GS_rho770_m2);
1823
1824 return amplitude;
1825
1826}

Referenced by QCMCFilter::findD0Decay().

◆ init()

void D0Topippim2pi0::init ( )

Definition at line 13 of file D0Topippim2pi0.cxx.

13 {
14
15 //std::cout << "D0Topippim2pi0 ==> Initialization !" << std::endl;
16
17 double mag[35], pha[35];
18 mag[0]= 100.0; pha[0]= 0.0;
19 mag[1]= 7.95507; pha[1]= -0.0687407;
20 mag[2]= 37.5559; pha[2]= -1.74946;
21 mag[3]= 61.2172; pha[3]= 2.98079;
22 mag[4]= 187.79; pha[4]= 2.64471;
23 mag[5]= 385.474; pha[5]= -0.137107;
24 mag[6]= 0.330788; pha[6]= 0.268133;
25 mag[7]= 0.584175; pha[7]= -2.89693;
26 mag[8]= 127.158; pha[8]= -2.47773;
27 mag[9]= 339.914; pha[9]= 2.22856;
28 mag[10]=0.320888; pha[10]=-2.6194;
29 mag[11]=0.366283; pha[11]=-0.26867;
30 mag[12]=14.1344; pha[12]=-0.41164;
31 mag[13]=86.0865; pha[13]=-2.49649;
32 mag[14]=6.1541; pha[14]=-1.18299;
33 mag[15]=56.6067; pha[15]=0.142977;
34 mag[16]=92.3073; pha[16]=-2.15881;
35 mag[17]=80.9453; pha[17]=0.825815;
36 mag[18]=16.9555; pha[18]=-2.98994;
37 mag[19]=9.72524; pha[19]=-1.39929;
38 mag[20]=5.71448; pha[20]=0.271902;
39 mag[21]=21.4195; pha[21]=-1.23701;
40 mag[22]=56.8867; pha[22]=-0.385837;
41 mag[23]=231.626; pha[23]=2.14842;
42 mag[24]=2938.45; pha[24]=-0.693491;
43 mag[25]=7252.7; pha[25]=2.23659;
44 mag[26]=5165.87; pha[26]=0.913557;
45 mag[27]=11508.6; pha[27]=-1.07187;
46 mag[28]=2461.86; pha[28]=1.8709;
47 mag[29]=8757.75; pha[29]=2.40756;
48 mag[30]=19.7413; pha[30]=-1.0753;
49 mag[31]=66.3826; pha[31]=2.34666;
50 mag[32]=11.2904; pha[32]=-0.822345;
51 mag[33]=2.04576; pha[33]=-0.281429;
52 mag[34]=0.57927; pha[34]=2.7182;
53
54 fitpara.clear();
55 for(int i=0; i<35; i++){
56 complex<double> ctemp(mag[i]*cos(pha[i]),mag[i]*sin(pha[i]));
57 fitpara.push_back(ctemp);
58 }
59
60 g_uv.clear();
61 for(int i=0; i<4; i++){
62 for(int j=0; j<4; j++){
63 if(i!=j){
64 g_uv.push_back(0.0);
65 }else if(i<3){
66 g_uv.push_back(-1.0);
67 }else if(i==3){
68 g_uv.push_back(1.0);
69 }
70 }
71 }
72
73 epsilon_uvmn.clear();
74 for(int i=0; i<4; i++){
75 for(int j=0; j<4; j++){
76 for(int k=0; k<4; k++){
77 for(int l=0; l<4; l++){
78 if(i==j || i==k || i==l || j==k || j==l || k==l){
79 epsilon_uvmn.push_back(0.0);
80 }else{
81 if(i==0 && j==1 && k==2 && l==3) epsilon_uvmn.push_back(1.0);
82 if(i==0 && j==1 && k==3 && l==2) epsilon_uvmn.push_back(-1.0);
83 if(i==0 && j==2 && k==1 && l==3) epsilon_uvmn.push_back(-1.0);
84 if(i==0 && j==2 && k==3 && l==1) epsilon_uvmn.push_back(1.0);
85 if(i==0 && j==3 && k==1 && l==2) epsilon_uvmn.push_back(1.0);
86 if(i==0 && j==3 && k==2 && l==1) epsilon_uvmn.push_back(-1.0);
87
88 if(i==1 && j==0 && k==2 && l==3) epsilon_uvmn.push_back(-1.0);
89 if(i==1 && j==0 && k==3 && l==2) epsilon_uvmn.push_back(1.0);
90 if(i==1 && j==2 && k==0 && l==3) epsilon_uvmn.push_back(1.0);
91 if(i==1 && j==2 && k==3 && l==0) epsilon_uvmn.push_back(-1.0);
92 if(i==1 && j==3 && k==0 && l==2) epsilon_uvmn.push_back(-1.0);
93 if(i==1 && j==3 && k==2 && l==0) epsilon_uvmn.push_back(1.0);
94
95 if(i==2 && j==0 && k==1 && l==3) epsilon_uvmn.push_back(1.0);
96 if(i==2 && j==0 && k==3 && l==1) epsilon_uvmn.push_back(-1.0);
97 if(i==2 && j==1 && k==0 && l==3) epsilon_uvmn.push_back(-1.0);
98 if(i==2 && j==1 && k==3 && l==0) epsilon_uvmn.push_back(1.0);
99 if(i==2 && j==3 && k==0 && l==1) epsilon_uvmn.push_back(1.0);
100 if(i==2 && j==3 && k==1 && l==0) epsilon_uvmn.push_back(-1.0);
101
102 if(i==3 && j==0 && k==1 && l==2) epsilon_uvmn.push_back(-1.0);
103 if(i==3 && j==0 && k==2 && l==1) epsilon_uvmn.push_back(1.0);
104 if(i==3 && j==1 && k==0 && l==2) epsilon_uvmn.push_back(1.0);
105 if(i==3 && j==1 && k==2 && l==0) epsilon_uvmn.push_back(-1.0);
106 if(i==3 && j==2 && k==0 && l==1) epsilon_uvmn.push_back(-1.0);
107 if(i==3 && j==2 && k==1 && l==0) epsilon_uvmn.push_back(1.0);
108
109 }
110 }
111 }
112 }
113 }
114
115 _nd = 4;
116 math_pi = 3.1415926f;
117 mass_Pion = 0.13957f;
118
119 rRes = 3.0*0.197321;
120 rD = 5.0*0.197321;
121 m_Pi = 0.139570;
122 m_Pi0 = 0.134977;
123 m2_Pi = m_Pi*m_Pi;
124 m2_Pi0 = m_Pi0*m_Pi0;
125
126 m0_rho7700 = 0.77526;
127 w0_rho7700 = 0.1478;
128
129 m0_rho770p = 0.77511;
130 w0_rho770p = 0.1491;
131
132 m0_rho1450 = 1.465;
133 w0_rho1450 = 0.400;
134
135 m0_f21270 = 1.2755;
136 w0_f21270 = 0.1867;
137
138 m0_a11260 = 1.1337;
139 g1_a11260 = 0.00335;
140 g2_a11260 = 0.0;
141
142 m0_pi1300 = 1.498;
143 w0_pi1300 = 0.590;
144
145 m0_a11420 = 1.411;
146 w0_a11420 = 0.161;
147
148 m0_a11640 = 1.655;
149 w0_a11640 = 0.254;
150
151 m0_a21320 = 1.3186;
152 w0_a21320 = 0.105;
153
154 m0_pi21670 = 1.6706;
155 w0_pi21670 = 0.258;
156
157 m0_pi11400 = 1.354;
158 w0_pi11400 = 0.330;
159
160 m0_h11170 = 1.166;
161 w0_h11170 = 0.375;
162
163 m0_omega = 0.78265;
164 w0_omega = 0.00849;
165
166 m0_phi = 1.019461;
167 w0_phi = 0.004249;
168
169 s0_prod = -5.0;
170
171}
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213

Referenced by QCMCFilter::findD0Decay().


The documentation for this class was generated from the following files: