BOSS 7.1.1
BESIII Offline Software System
Loading...
Searching...
No Matches
D0To2pip2pim Class Reference

#include <D0To2pip2pim.h>

Public Member Functions

 D0To2pip2pim ()
 
virtual ~D0To2pip2pim ()
 
void init ()
 
complex< double > Amp (vector< double > Pip1, vector< double > Pim1, vector< double > Pip2, vector< double > Pim2)
 

Detailed Description

Definition at line 9 of file D0To2pip2pim.h.

Constructor & Destructor Documentation

◆ D0To2pip2pim()

D0To2pip2pim::D0To2pip2pim ( )
inline

Definition at line 13 of file D0To2pip2pim.h.

13{}

◆ ~D0To2pip2pim()

D0To2pip2pim::~D0To2pip2pim ( )
virtual

Definition at line 11 of file D0To2pip2pim.cxx.

11{}

Member Function Documentation

◆ Amp()

complex< double > D0To2pip2pim::Amp ( vector< double > Pip1,
vector< double > Pim1,
vector< double > Pip2,
vector< double > Pim2 )

Definition at line 1426 of file D0To2pip2pim.cxx.

1427{
1428
1429 vector<double> Pip1Pim1; Pip1Pim1.clear();
1430 vector<double> Pip1Pim2; Pip1Pim2.clear();
1431 vector<double> Pip2Pim1; Pip2Pim1.clear();
1432 vector<double> Pip2Pim2; Pip2Pim2.clear();
1433
1434 Pip1Pim1 = sum_tensor(Pip1, Pim1);
1435 Pip1Pim2 = sum_tensor(Pip1, Pim2);
1436 Pip2Pim1 = sum_tensor(Pip2, Pim1);
1437 Pip2Pim2 = sum_tensor(Pip2, Pim2);
1438
1439 vector<double> Pip1Pip2Pim1; Pip1Pip2Pim1.clear();
1440 vector<double> Pip1Pip2Pim2; Pip1Pip2Pim2.clear();
1441 vector<double> Pim1Pim2Pip1; Pim1Pim2Pip1.clear();
1442 vector<double> Pim1Pim2Pip2; Pim1Pim2Pip2.clear();
1443
1444 Pip1Pip2Pim1 = sum_tensor(Pip1Pim1, Pip2);
1445 Pip1Pip2Pim2 = sum_tensor(Pip1Pim2, Pip2);
1446 Pim1Pim2Pip1 = sum_tensor(Pip1Pim1, Pim2);
1447 Pim1Pim2Pip2 = sum_tensor(Pip2Pim1, Pim2);
1448
1449 vector<double> D0; D0.clear();
1450 D0 = sum_tensor(Pip1Pip2Pim1, Pim2);
1451
1452 double M2_Pip1Pim1 = contract_11_0(Pip1Pim1, Pip1Pim1);
1453 double M2_Pip1Pim2 = contract_11_0(Pip1Pim2, Pip1Pim2);
1454 double M2_Pip2Pim1 = contract_11_0(Pip2Pim1, Pip2Pim1);
1455 double M2_Pip2Pim2 = contract_11_0(Pip2Pim2, Pip2Pim2);
1456
1457 double M2_Pip1Pip2Pim1 = contract_11_0(Pip1Pip2Pim1, Pip1Pip2Pim1);
1458 double M2_Pip1Pip2Pim2 = contract_11_0(Pip1Pip2Pim2, Pip1Pip2Pim2);
1459 double M2_Pim1Pim2Pip1 = contract_11_0(Pim1Pim2Pip1, Pim1Pim2Pip1);
1460 double M2_Pim1Pim2Pip2 = contract_11_0(Pim1Pim2Pip2, Pim1Pim2Pip2);
1461 double M2_D0 = contract_11_0(D0, D0);
1462
1463 complex<double> GS_rho770_11 = GS(M2_Pip1Pim1, m0_rho770, w0_rho770, m2_Pi, m2_Pi, rRes, 1);
1464 complex<double> GS_rho770_12 = GS(M2_Pip1Pim2, m0_rho770, w0_rho770, m2_Pi, m2_Pi, rRes, 1);
1465 complex<double> GS_rho770_21 = GS(M2_Pip2Pim1, m0_rho770, w0_rho770, m2_Pi, m2_Pi, rRes, 1);
1466 complex<double> GS_rho770_22 = GS(M2_Pip2Pim2, m0_rho770, w0_rho770, m2_Pi, m2_Pi, rRes, 1);
1467
1468 complex<double> GS_rho1450_11 = GS(M2_Pip1Pim1, m0_rho1450, w0_rho1450, m2_Pi, m2_Pi, rRes, 1);
1469 complex<double> GS_rho1450_12 = GS(M2_Pip1Pim2, m0_rho1450, w0_rho1450, m2_Pi, m2_Pi, rRes, 1);
1470 complex<double> GS_rho1450_21 = GS(M2_Pip2Pim1, m0_rho1450, w0_rho1450, m2_Pi, m2_Pi, rRes, 1);
1471 complex<double> GS_rho1450_22 = GS(M2_Pip2Pim2, m0_rho1450, w0_rho1450, m2_Pi, m2_Pi, rRes, 1);
1472
1473 complex<double> RBW_f21270_11 = RBW(M2_Pip1Pim1, m0_f21270, w0_f21270, m2_Pi, m2_Pi, rRes, 2);
1474 complex<double> RBW_f21270_12 = RBW(M2_Pip1Pim2, m0_f21270, w0_f21270, m2_Pi, m2_Pi, rRes, 2);
1475 complex<double> RBW_f21270_21 = RBW(M2_Pip2Pim1, m0_f21270, w0_f21270, m2_Pi, m2_Pi, rRes, 2);
1476 complex<double> RBW_f21270_22 = RBW(M2_Pip2Pim2, m0_f21270, w0_f21270, m2_Pi, m2_Pi, rRes, 2);
1477
1478 complex<double> PiPiS_11_0 = Fvector(M2_Pip1Pim1, s0_prod, 0);
1479 complex<double> PiPiS_12_0 = Fvector(M2_Pip1Pim2, s0_prod, 0);
1480 complex<double> PiPiS_21_0 = Fvector(M2_Pip2Pim1, s0_prod, 0);
1481 complex<double> PiPiS_22_0 = Fvector(M2_Pip2Pim2, s0_prod, 0);
1482
1483 complex<double> PiPiS_11_1 = Fvector(M2_Pip1Pim1, s0_prod, 1);
1484 complex<double> PiPiS_12_1 = Fvector(M2_Pip1Pim2, s0_prod, 1);
1485 complex<double> PiPiS_21_1 = Fvector(M2_Pip2Pim1, s0_prod, 1);
1486 complex<double> PiPiS_22_1 = Fvector(M2_Pip2Pim2, s0_prod, 1);
1487
1488 complex<double> PiPiS_11_5 = Fvector(M2_Pip1Pim1, s0_prod, 5);
1489 complex<double> PiPiS_12_5 = Fvector(M2_Pip1Pim2, s0_prod, 5);
1490 complex<double> PiPiS_21_5 = Fvector(M2_Pip2Pim1, s0_prod, 5);
1491 complex<double> PiPiS_22_5 = Fvector(M2_Pip2Pim2, s0_prod, 5);
1492
1493 complex<double> PiPiS_11_6 = Fvector(M2_Pip1Pim1, s0_prod, 6);
1494 complex<double> PiPiS_12_6 = Fvector(M2_Pip1Pim2, s0_prod, 6);
1495 complex<double> PiPiS_21_6 = Fvector(M2_Pip2Pim1, s0_prod, 6);
1496 complex<double> PiPiS_22_6 = Fvector(M2_Pip2Pim2, s0_prod, 6);
1497
1498 complex<double> RBW_a11260p_1 = RBWa1260(M2_Pip1Pip2Pim1, m0_a11260, g1_a11260, g2_a11260);
1499 complex<double> RBW_a11260p_2 = RBWa1260(M2_Pip1Pip2Pim2, m0_a11260, g1_a11260, g2_a11260);
1500 complex<double> RBW_a11260m_1 = RBWa1260(M2_Pim1Pim2Pip1, m0_a11260, g1_a11260, g2_a11260);
1501 complex<double> RBW_a11260m_2 = RBWa1260(M2_Pim1Pim2Pip2, m0_a11260, g1_a11260, g2_a11260);
1502
1503 complex<double> RBW_a21320p_1 = RBW(M2_Pip1Pip2Pim1, m0_a21320, w0_a21320, -1.0, -1.0, -1, -1);
1504 complex<double> RBW_a21320p_2 = RBW(M2_Pip1Pip2Pim2, m0_a21320, w0_a21320, -1.0, -1.0, -1, -1);
1505 complex<double> RBW_a21320m_1 = RBW(M2_Pim1Pim2Pip1, m0_a21320, w0_a21320, -1.0, -1.0, -1, -1);
1506 complex<double> RBW_a21320m_2 = RBW(M2_Pim1Pim2Pip2, m0_a21320, w0_a21320, -1.0, -1.0, -1, -1);
1507
1508 complex<double> RBW_pi1300p_1 = RBWpi1300(M2_Pip1Pip2Pim1, m0_pi1300, w0_pi1300);
1509 complex<double> RBW_pi1300p_2 = RBWpi1300(M2_Pip1Pip2Pim2, m0_pi1300, w0_pi1300);
1510 complex<double> RBW_pi1300m_1 = RBWpi1300(M2_Pim1Pim2Pip1, m0_pi1300, w0_pi1300);
1511 complex<double> RBW_pi1300m_2 = RBWpi1300(M2_Pim1Pim2Pip2, m0_pi1300, w0_pi1300);
1512
1513 complex<double> RBW_a11420p_1 = RBW(M2_Pip1Pip2Pim1, m0_a11420, w0_a11420, -1.0, -1.0, -1, -1);
1514 complex<double> RBW_a11420p_2 = RBW(M2_Pip1Pip2Pim2, m0_a11420, w0_a11420, -1.0, -1.0, -1, -1);
1515 complex<double> RBW_a11420m_1 = RBW(M2_Pim1Pim2Pip1, m0_a11420, w0_a11420, -1.0, -1.0, -1, -1);
1516 complex<double> RBW_a11420m_2 = RBW(M2_Pim1Pim2Pip2, m0_a11420, w0_a11420, -1.0, -1.0, -1, -1);
1517
1518 // D->XP 1-rank Projection
1519 vector<double> Proj1_3p1; Proj1_3p1.clear();
1520 vector<double> Proj1_3p2; Proj1_3p2.clear();
1521 vector<double> Proj1_3m1; Proj1_3m1.clear();
1522 vector<double> Proj1_3m2; Proj1_3m2.clear();
1523
1524 Proj1_3p1 = ProjectionTensors(Pip1Pip2Pim1,1);
1525 Proj1_3p2 = ProjectionTensors(Pip1Pip2Pim2,1);
1526 Proj1_3m1 = ProjectionTensors(Pim1Pim2Pip1,1);
1527 Proj1_3m2 = ProjectionTensors(Pim1Pim2Pip2,1);
1528
1529 // D->XP 2-rank Projection
1530 vector<double> Proj2_3p1; Proj2_3p1.clear();
1531 vector<double> Proj2_3p2; Proj2_3p2.clear();
1532 vector<double> Proj2_3m1; Proj2_3m1.clear();
1533 vector<double> Proj2_3m2; Proj2_3m2.clear();
1534
1535 Proj2_3p1 = ProjectionTensors(Pip1Pip2Pim1,2);
1536 Proj2_3p2 = ProjectionTensors(Pip1Pip2Pim2,2);
1537 Proj2_3m1 = ProjectionTensors(Pim1Pim2Pip1,2);
1538 Proj2_3m2 = ProjectionTensors(Pim1Pim2Pip2,2);
1539
1540 // X->PP Orbital
1541 vector<double> T1_Pip1Pim1; T1_Pip1Pim1.clear();
1542 vector<double> T1_Pip1Pim2; T1_Pip1Pim2.clear();
1543 vector<double> T1_Pip2Pim1; T1_Pip2Pim1.clear();
1544 vector<double> T1_Pip2Pim2; T1_Pip2Pim2.clear();
1545
1546 T1_Pip1Pim1 = OrbitalTensors(Pip1Pim1, Pip1, Pim1, rRes, 1);
1547 T1_Pip1Pim2 = OrbitalTensors(Pip1Pim2, Pip1, Pim2, rRes, 1);
1548 T1_Pip2Pim1 = OrbitalTensors(Pip2Pim1, Pip2, Pim1, rRes, 1);
1549 T1_Pip2Pim2 = OrbitalTensors(Pip2Pim2, Pip2, Pim2, rRes, 1);
1550
1551 vector<double> T1_Pim1Pip1; T1_Pim1Pip1.clear();
1552 vector<double> T1_Pim1Pip2; T1_Pim1Pip2.clear();
1553 vector<double> T1_Pim2Pip1; T1_Pim2Pip1.clear();
1554 vector<double> T1_Pim2Pip2; T1_Pim2Pip2.clear();
1555
1556 T1_Pim1Pip1 = OrbitalTensors(Pip1Pim1, Pim1, Pip1, rRes, 1);
1557 T1_Pim1Pip2 = OrbitalTensors(Pip2Pim1, Pim1, Pip2, rRes, 1);
1558 T1_Pim2Pip1 = OrbitalTensors(Pip1Pim2, Pim2, Pip1, rRes, 1);
1559 T1_Pim2Pip2 = OrbitalTensors(Pip2Pim2, Pim2, Pip2, rRes, 1);
1560
1561 vector<double> T2_Pip1Pim1; T2_Pip1Pim1.clear();
1562 vector<double> T2_Pip1Pim2; T2_Pip1Pim2.clear();
1563 vector<double> T2_Pip2Pim1; T2_Pip2Pim1.clear();
1564 vector<double> T2_Pip2Pim2; T2_Pip2Pim1.clear();
1565
1566 T2_Pip1Pim1 = OrbitalTensors(Pip1Pim1, Pip1, Pim1, rRes, 2);
1567 T2_Pip1Pim2 = OrbitalTensors(Pip1Pim2, Pip1, Pim2, rRes, 2);
1568 T2_Pip2Pim1 = OrbitalTensors(Pip2Pim1, Pip2, Pim1, rRes, 2);
1569 T2_Pip2Pim2 = OrbitalTensors(Pip2Pim2, Pip2, Pim2, rRes, 2);
1570
1571 // X->YP Orbital
1572 vector<double> T1_Pip1Pim1Pip2; T1_Pip1Pim1Pip2.clear();
1573 vector<double> T1_Pip2Pim1Pip1; T1_Pip2Pim1Pip1.clear();
1574 vector<double> T1_Pip1Pim2Pip2; T1_Pip1Pim2Pip2.clear();
1575 vector<double> T1_Pip2Pim2Pip1; T1_Pip2Pim2Pip1.clear();
1576 vector<double> T1_Pip1Pim1Pim2; T1_Pip1Pim1Pim2.clear();
1577 vector<double> T1_Pip1Pim2Pim1; T1_Pip1Pim2Pim1.clear();
1578 vector<double> T1_Pip2Pim1Pim2; T1_Pip2Pim1Pim2.clear();
1579 vector<double> T1_Pip2Pim2Pim1; T1_Pip2Pim2Pim1.clear();
1580
1581 T1_Pip1Pim1Pip2 = OrbitalTensors(Pip1Pip2Pim1, Pip1Pim1, Pip2, rRes, 1);
1582 T1_Pip2Pim1Pip1 = OrbitalTensors(Pip1Pip2Pim1, Pip2Pim1, Pip1, rRes, 1);
1583 T1_Pip1Pim2Pip2 = OrbitalTensors(Pip1Pip2Pim2, Pip1Pim2, Pip2, rRes, 1);
1584 T1_Pip2Pim2Pip1 = OrbitalTensors(Pip1Pip2Pim2, Pip2Pim2, Pip1, rRes, 1);
1585 T1_Pip1Pim1Pim2 = OrbitalTensors(Pim1Pim2Pip1, Pip1Pim1, Pim2, rRes, 1);
1586 T1_Pip2Pim1Pim2 = OrbitalTensors(Pim1Pim2Pip2, Pip2Pim1, Pim2, rRes, 1);
1587 T1_Pip1Pim2Pim1 = OrbitalTensors(Pim1Pim2Pip1, Pip1Pim2, Pim1, rRes, 1);
1588 T1_Pip2Pim2Pim1 = OrbitalTensors(Pim1Pim2Pip2, Pip2Pim2, Pim1, rRes, 1);
1589
1590 vector<double> T2_Pip1Pim1Pip2; T2_Pip1Pim1Pip2.clear();
1591 vector<double> T2_Pip2Pim1Pip1; T2_Pip2Pim1Pip1.clear();
1592 vector<double> T2_Pip1Pim2Pip2; T2_Pip1Pim2Pip2.clear();
1593 vector<double> T2_Pip2Pim2Pip1; T2_Pip2Pim2Pip1.clear();
1594 vector<double> T2_Pip1Pim1Pim2; T2_Pip1Pim1Pim2.clear();
1595 vector<double> T2_Pip2Pim1Pim2; T2_Pip2Pim1Pim2.clear();
1596 vector<double> T2_Pip1Pim2Pim1; T2_Pip1Pim2Pim1.clear();
1597 vector<double> T2_Pip2Pim2Pim1; T2_Pip2Pim2Pim1.clear();
1598
1599 T2_Pip1Pim1Pip2 = OrbitalTensors(Pip1Pip2Pim1, Pip1Pim1, Pip2, rRes, 2);
1600 T2_Pip2Pim1Pip1 = OrbitalTensors(Pip1Pip2Pim1, Pip2Pim1, Pip1, rRes, 2);
1601 T2_Pip1Pim2Pip2 = OrbitalTensors(Pip1Pip2Pim2, Pip1Pim2, Pip2, rRes, 2);
1602 T2_Pip2Pim2Pip1 = OrbitalTensors(Pip1Pip2Pim2, Pip2Pim2, Pip1, rRes, 2);
1603 T2_Pip1Pim1Pim2 = OrbitalTensors(Pim1Pim2Pip1, Pip1Pim1, Pim2, rRes, 2);
1604 T2_Pip2Pim1Pim2 = OrbitalTensors(Pim1Pim2Pip2, Pip2Pim1, Pim2, rRes, 2);
1605 T2_Pip1Pim2Pim1 = OrbitalTensors(Pim1Pim2Pip1, Pip1Pim2, Pim1, rRes, 2);
1606 T2_Pip2Pim2Pim1 = OrbitalTensors(Pim1Pim2Pip2, Pip2Pim2, Pim1, rRes, 2);
1607
1608 // D->XX Orbital
1609 vector<double> T1_2z11; T1_2z11.clear();
1610 vector<double> T1_2z12; T1_2z12.clear();
1611 vector<double> T1_2z21; T1_2z21.clear();
1612 vector<double> T1_2z22; T1_2z22.clear();
1613
1614 T1_2z11 = OrbitalTensors(D0, Pip1Pim1, Pip2Pim2, rD, 1);
1615 T1_2z12 = OrbitalTensors(D0, Pip2Pim2, Pip1Pim1, rD, 1);
1616 T1_2z21 = OrbitalTensors(D0, Pip1Pim2, Pip2Pim1, rD, 1);
1617 T1_2z22 = OrbitalTensors(D0, Pip2Pim1, Pip1Pim2, rD, 1);
1618
1619 vector<double> T2_2z11; T2_2z11.clear();
1620 vector<double> T2_2z12; T2_2z12.clear();
1621 vector<double> T2_2z21; T2_2z21.clear();
1622 vector<double> T2_2z22; T2_2z22.clear();
1623
1624 T2_2z11 = OrbitalTensors(D0, Pip1Pim1, Pip2Pim2, rD, 2);
1625 T2_2z12 = OrbitalTensors(D0, Pip2Pim2, Pip1Pim1, rD, 2);
1626 T2_2z21 = OrbitalTensors(D0, Pip1Pim2, Pip2Pim1, rD, 2);
1627 T2_2z22 = OrbitalTensors(D0, Pip2Pim1, Pip1Pim2, rD, 2);
1628
1629 // D->XP Orbital
1630 vector<double> T1_3p1; T1_3p1.clear();
1631 vector<double> T1_3p2; T1_3p2.clear();
1632 vector<double> T1_3m1; T1_3m1.clear();
1633 vector<double> T1_3m2; T1_3m2.clear();
1634
1635 T1_3p1 = OrbitalTensors(D0, Pip1Pip2Pim1, Pim2, rD, 1);
1636 T1_3p2 = OrbitalTensors(D0, Pip1Pip2Pim2, Pim1, rD, 1);
1637 T1_3m1 = OrbitalTensors(D0, Pim1Pim2Pip1, Pip2, rD, 1);
1638 T1_3m2 = OrbitalTensors(D0, Pim1Pim2Pip2, Pip1, rD, 1);
1639
1640 vector<double> T2_3p1; T2_3p1.clear();
1641 vector<double> T2_3p2; T2_3p2.clear();
1642 vector<double> T2_3m1; T2_3m1.clear();
1643 vector<double> T2_3m2; T2_3m2.clear();
1644
1645 T2_3p1 = OrbitalTensors(D0, Pip1Pip2Pim1, Pim2, rD, 2);
1646 T2_3p2 = OrbitalTensors(D0, Pip1Pip2Pim2, Pim1, rD, 2);
1647 T2_3m1 = OrbitalTensors(D0, Pim1Pim2Pip1, Pip2, rD, 2);
1648 T2_3m2 = OrbitalTensors(D0, Pim1Pim2Pip2, Pip1, rD, 2);
1649
1650 complex<double> amplitude(0,0);
1651 vector< complex<double> > g_fitpara;g_fitpara.clear();
1652 g_fitpara.push_back(complex<double>(100.0*cos(0.0), 100.0*sin(0.0)));
1653 g_fitpara.push_back(complex<double>(7.95507 *cos(-0.0687407), 7.95507 *sin(-0.0687407)));
1654 g_fitpara.push_back(complex<double>(37.5559 *cos(-1.74946 ), 37.5559 *sin(-1.74946 )));
1655 g_fitpara.push_back(complex<double>(61.2172 *cos(2.98079 ), 61.2172 *sin(2.98079 )));
1656 g_fitpara.push_back(complex<double>(187.79 *cos(2.64471 ), 187.79 *sin(2.64471 )));
1657 g_fitpara.push_back(complex<double>(385.474 *cos(-0.137107 ), 385.474 *sin(-0.137107 )));
1658 g_fitpara.push_back(complex<double>(0.330788*cos(0.268133 ), 0.330788*sin(0.268133 )));
1659 g_fitpara.push_back(complex<double>(127.158 *cos(-2.47773 ), 127.158 *sin(-2.47773 )));
1660 g_fitpara.push_back(complex<double>(339.914 *cos(2.22856 ), 339.914 *sin(2.22856 )));
1661 g_fitpara.push_back(complex<double>(0.320888*cos(-2.6194 ), 0.320888*sin(-2.6194 )));
1662 g_fitpara.push_back(complex<double>(0.366283*cos(-0.26867 ), 0.366283*sin(-0.26867 )));
1663 g_fitpara.push_back(complex<double>(86.0865 *cos(-2.49649 ), 86.0865 *sin(-2.49649 )));
1664 g_fitpara.push_back(complex<double>(6.1541 *cos(-1.18299 ), 6.1541 *sin(-1.18299 )));
1665 g_fitpara.push_back(complex<double>(56.6067 *cos(0.142977 ), 56.6067 *sin(0.142977 )));
1666 g_fitpara.push_back(complex<double>(92.3073 *cos(-2.15881 ), 92.3073 *sin(-2.15881 )));
1667 g_fitpara.push_back(complex<double>(10.5885 *cos(-3.03166 ), 10.5885 *sin(-3.03166 )));
1668 g_fitpara.push_back(complex<double>(8.36765 *cos(1.8417 ), 8.36765 *sin(1.8417 )));
1669 g_fitpara.push_back(complex<double>(6.56437 *cos(-2.93087 ), 6.56437 *sin(-2.93087 )));
1670 g_fitpara.push_back(complex<double>(15.7197 *cos(0.96925 ), 15.7197 *sin(0.96925 )));
1671 g_fitpara.push_back(complex<double>(21.4195 *cos(-1.23701 ), 21.4195 *sin(-1.23701 )));
1672 g_fitpara.push_back(complex<double>(56.8867 *cos(-0.385837 ), 56.8867 *sin(-0.385837 )));
1673 g_fitpara.push_back(complex<double>(231.626 *cos(2.14842 ), 231.626 *sin(2.14842 )));
1674 g_fitpara.push_back(complex<double>(2938.45 *cos(-0.693491 ), 2938.45 *sin(-0.693491 )));
1675 g_fitpara.push_back(complex<double>(7252.7 *cos(2.23659 ), 7252.7 *sin(2.23659 )));
1676 g_fitpara.push_back(complex<double>(5165.87 *cos(0.913557 ), 5165.87 *sin(0.913557 )));
1677 g_fitpara.push_back(complex<double>(11508.6 *cos(-1.07187 ), 11508.6 *sin(-1.07187 )));
1678 g_fitpara.push_back(complex<double>(2461.86 *cos(1.8709 ), 2461.86 *sin(1.8709 )));
1679 g_fitpara.push_back(complex<double>(8757.75 *cos(2.40756 ), 8757.75 *sin(2.40756 )));
1680 g_fitpara.push_back(complex<double>(19.7413 *cos(-1.0753 ), 19.7413 *sin(-1.0753 )));
1681 g_fitpara.push_back(complex<double>(66.3826 *cos(2.34666 ), 66.3826 *sin(2.34666 )));
1682
1683 // D0 -> a1(1260)+ {rho(770)0 pi+[S]} pi-
1684 double SF_Ap_S_VP_1 = contract_11_0(contract_21_1(Proj1_3p1, T1_Pip2Pim1), T1_3p1);
1685 double SF_Ap_S_VP_2 = contract_11_0(contract_21_1(Proj1_3p1, T1_Pip1Pim1), T1_3p1);
1686 double SF_Ap_S_VP_3 = contract_11_0(contract_21_1(Proj1_3p2, T1_Pip2Pim2), T1_3p2);
1687 double SF_Ap_S_VP_4 = contract_11_0(contract_21_1(Proj1_3p2, T1_Pip1Pim2), T1_3p2);
1688 amplitude += g_fitpara[0]*(SF_Ap_S_VP_1*RBW_a11260p_1*GS_rho770_21 + SF_Ap_S_VP_2*RBW_a11260p_1*GS_rho770_11 + SF_Ap_S_VP_3*RBW_a11260p_2*GS_rho770_22 + SF_Ap_S_VP_4*RBW_a11260p_2*GS_rho770_12);
1689
1690 // D0 -> a1(1260)+ {rho(770)0 pi+[D]} pi-
1691 double SF_Ap_D_VP_1 = contract_11_0(contract_21_1(T2_Pip2Pim1Pip1, T1_Pip2Pim1), T1_3p1);
1692 double SF_Ap_D_VP_2 = contract_11_0(contract_21_1(T2_Pip1Pim1Pip2, T1_Pip1Pim1), T1_3p1);
1693 double SF_Ap_D_VP_3 = contract_11_0(contract_21_1(T2_Pip2Pim2Pip1, T1_Pip2Pim2), T1_3p2);
1694 double SF_Ap_D_VP_4 = contract_11_0(contract_21_1(T2_Pip1Pim2Pip2, T1_Pip1Pim2), T1_3p2);
1695
1696 amplitude += g_fitpara[1]*(SF_Ap_D_VP_1*RBW_a11260p_1*GS_rho770_21 + SF_Ap_D_VP_2*RBW_a11260p_1*GS_rho770_11 + SF_Ap_D_VP_3*RBW_a11260p_2*GS_rho770_22 + SF_Ap_D_VP_4*RBW_a11260p_2*GS_rho770_12);
1697
1698 // D0 -> a1(1260)+ {f2(1270) pi+ [P]} pi-
1699 double SF_Ap_P_TP_1 = contract_11_0(contract_21_1(contract_42_2(Proj2_3p1, T2_Pip2Pim1), T1_Pip2Pim1Pip1), T1_3p1);
1700 double SF_Ap_P_TP_2 = contract_11_0(contract_21_1(contract_42_2(Proj2_3p1, T2_Pip1Pim1), T1_Pip1Pim1Pip2), T1_3p1);
1701 double SF_Ap_P_TP_3 = contract_11_0(contract_21_1(contract_42_2(Proj2_3p2, T2_Pip2Pim2), T1_Pip2Pim2Pip1), T1_3p2);
1702 double SF_Ap_P_TP_4 = contract_11_0(contract_21_1(contract_42_2(Proj2_3p2, T2_Pip1Pim2), T1_Pip1Pim2Pip2), T1_3p2);
1703
1704 amplitude += g_fitpara[2]*(SF_Ap_P_TP_1*RBW_a11260p_1*RBW_f21270_21 + SF_Ap_P_TP_2*RBW_a11260p_1*RBW_f21270_11 + SF_Ap_P_TP_3*RBW_a11260p_2*RBW_f21270_22 + SF_Ap_P_TP_4*RBW_a11260p_2*RBW_f21270_12);
1705
1706 // D0 -> a1(1260)+ {f0 pi+ [P]} pi-
1707 double SF_Ap_P_SP_1 = contract_11_0(T1_3p1, T1_Pip2Pim1Pip1);
1708 double SF_Ap_P_SP_2 = contract_11_0(T1_3p1, T1_Pip1Pim1Pip2);
1709 double SF_Ap_P_SP_3 = contract_11_0(T1_3p2, T1_Pip2Pim2Pip1);
1710 double SF_Ap_P_SP_4 = contract_11_0(T1_3p2, T1_Pip1Pim2Pip2);
1711
1712 amplitude += g_fitpara[3]*(SF_Ap_P_SP_1*RBW_a11260p_1*PiPiS_21_0 + SF_Ap_P_SP_2*RBW_a11260p_1*PiPiS_11_0 + SF_Ap_P_SP_3*RBW_a11260p_2*PiPiS_22_0 + SF_Ap_P_SP_4*RBW_a11260p_2*PiPiS_12_0);
1713 amplitude += g_fitpara[4]*(SF_Ap_P_SP_1*RBW_a11260p_1*PiPiS_21_1 + SF_Ap_P_SP_2*RBW_a11260p_1*PiPiS_11_1 + SF_Ap_P_SP_3*RBW_a11260p_2*PiPiS_22_1 + SF_Ap_P_SP_4*RBW_a11260p_2*PiPiS_12_1);
1714 amplitude += g_fitpara[5]*(SF_Ap_P_SP_1*RBW_a11260p_1*PiPiS_21_5 + SF_Ap_P_SP_2*RBW_a11260p_1*PiPiS_11_5 + SF_Ap_P_SP_3*RBW_a11260p_2*PiPiS_22_5 + SF_Ap_P_SP_4*RBW_a11260p_2*PiPiS_12_5);
1715
1716 // D0->a1(1260)- {rho(770)0 pi- [S]} pi+
1717 double SF_Am_S_VP_1 = contract_11_0(contract_21_1(Proj1_3m1, T1_Pim2Pip1), T1_3m1);
1718 double SF_Am_S_VP_2 = contract_11_0(contract_21_1(Proj1_3m1, T1_Pim1Pip1), T1_3m1);
1719 double SF_Am_S_VP_3 = contract_11_0(contract_21_1(Proj1_3m2, T1_Pim2Pip2), T1_3m2);
1720 double SF_Am_S_VP_4 = contract_11_0(contract_21_1(Proj1_3m2, T1_Pim1Pip2), T1_3m2);
1721
1722 amplitude += g_fitpara[0]*g_fitpara[6]*(SF_Am_S_VP_1*RBW_a11260m_1*GS_rho770_12 + SF_Am_S_VP_2*RBW_a11260m_1*GS_rho770_11 + SF_Am_S_VP_3*RBW_a11260m_2*GS_rho770_22 + SF_Am_S_VP_4*RBW_a11260m_2*GS_rho770_21);
1723
1724 // D0 -> a1(1260)- {rho(770)0 pi-[D]} pi+
1725 double SF_Am_D_VP_1 = contract_11_0(contract_21_1(T2_Pip1Pim2Pim1, T1_Pim2Pip1), T1_3m1);
1726 double SF_Am_D_VP_2 = contract_11_0(contract_21_1(T2_Pip1Pim1Pim2, T1_Pim1Pip1), T1_3m1);
1727 double SF_Am_D_VP_3 = contract_11_0(contract_21_1(T2_Pip2Pim2Pim1, T1_Pim2Pip2), T1_3m2);
1728 double SF_Am_D_VP_4 = contract_11_0(contract_21_1(T2_Pip2Pim1Pim2, T1_Pim1Pip2), T1_3m2);
1729
1730 amplitude += g_fitpara[1]*g_fitpara[6]*(SF_Am_D_VP_1*RBW_a11260m_1*GS_rho770_12 + SF_Am_D_VP_2*RBW_a11260m_1*GS_rho770_11 + SF_Am_D_VP_3*RBW_a11260m_2*GS_rho770_22 + SF_Am_D_VP_4*RBW_a11260m_2*GS_rho770_21);
1731
1732 // D0 -> a1(1260)- {f2(1270) pi- [P]} pi+
1733 double SF_Am_P_TP_1 = contract_11_0(contract_21_1(contract_42_2(Proj2_3m1, T2_Pip1Pim2), T1_Pip1Pim2Pim1), T1_3m1);
1734 double SF_Am_P_TP_2 = contract_11_0(contract_21_1(contract_42_2(Proj2_3m1, T2_Pip1Pim1), T1_Pip1Pim1Pim2), T1_3m1);
1735 double SF_Am_P_TP_3 = contract_11_0(contract_21_1(contract_42_2(Proj2_3m2, T2_Pip2Pim2), T1_Pip2Pim2Pim1), T1_3m2);
1736 double SF_Am_P_TP_4 = contract_11_0(contract_21_1(contract_42_2(Proj2_3m2, T2_Pip2Pim1), T1_Pip2Pim1Pim2), T1_3m2);
1737
1738 amplitude += g_fitpara[2]*g_fitpara[6]*(SF_Am_P_TP_1*RBW_a11260m_1*RBW_f21270_12 + SF_Am_P_TP_2*RBW_a11260m_1*RBW_f21270_11 + SF_Am_P_TP_3*RBW_a11260m_2*RBW_f21270_22 + SF_Am_P_TP_4*RBW_a11260m_2*RBW_f21270_21);
1739
1740 // D0 -> a1(1260)- {f0 pi- [P]} pi+
1741 double SF_Am_P_SP_1 = contract_11_0(T1_3m1, T1_Pip1Pim2Pim1);
1742 double SF_Am_P_SP_2 = contract_11_0(T1_3m1, T1_Pip1Pim1Pim2);
1743 double SF_Am_P_SP_3 = contract_11_0(T1_3m2, T1_Pip2Pim2Pim1);
1744 double SF_Am_P_SP_4 = contract_11_0(T1_3m2, T1_Pip2Pim1Pim2);
1745
1746 amplitude += g_fitpara[3]*g_fitpara[6]*(SF_Am_P_SP_1*RBW_a11260m_1*PiPiS_12_0 + SF_Am_P_SP_2*RBW_a11260m_1*PiPiS_11_0 + SF_Am_P_SP_3*RBW_a11260m_2*PiPiS_22_0 + SF_Am_P_SP_4*RBW_a11260m_2*PiPiS_21_0);
1747 amplitude += g_fitpara[4]*g_fitpara[6]*(SF_Am_P_SP_1*RBW_a11260m_1*PiPiS_12_1 + SF_Am_P_SP_2*RBW_a11260m_1*PiPiS_11_1 + SF_Am_P_SP_3*RBW_a11260m_2*PiPiS_22_1 + SF_Am_P_SP_4*RBW_a11260m_2*PiPiS_21_1);
1748 amplitude += g_fitpara[5]*g_fitpara[6]*(SF_Am_P_SP_1*RBW_a11260m_1*PiPiS_12_5 + SF_Am_P_SP_2*RBW_a11260m_1*PiPiS_11_5 + SF_Am_P_SP_3*RBW_a11260m_2*PiPiS_22_5 + SF_Am_P_SP_4*RBW_a11260m_2*PiPiS_21_5);
1749
1750 // D0 -> a1(1420)+ {rho(770)0 pi+[S]} pi-
1751// amplitude += g_fitpara[7]*(SF_Ap_S_VP_1*RBW_a11420p_1*GS_rho770_21 + SF_Ap_S_VP_2*RBW_a11420p_1*GS_rho770_11 + SF_Ap_S_VP_3*RBW_a11420p_2*GS_rho770_22 + SF_Ap_S_VP_4*RBW_a11420p_2*GS_rho770_12);
1752
1753 // D0 -> a1(1420)+ {f0 pi+[P]} pi-
1754 amplitude += g_fitpara[7]*(SF_Ap_P_SP_1*RBW_a11420p_1*PiPiS_21_5 + SF_Ap_P_SP_2*RBW_a11420p_1*PiPiS_11_5 + SF_Ap_P_SP_3*RBW_a11420p_2*PiPiS_22_5 + SF_Ap_P_SP_4*RBW_a11420p_2*PiPiS_12_5);
1755 amplitude += g_fitpara[8]*(SF_Ap_P_SP_1*RBW_a11420p_1*PiPiS_21_6 + SF_Ap_P_SP_2*RBW_a11420p_1*PiPiS_11_6 + SF_Ap_P_SP_3*RBW_a11420p_2*PiPiS_22_6 + SF_Ap_P_SP_4*RBW_a11420p_2*PiPiS_12_6);
1756
1757 vector<double> m_epsilon_uvmn;
1758 m_epsilon_uvmn.clear();
1759 for(int i=0; i<4; i++){
1760 for(int j=0; j<4; j++){
1761 for(int k=0; k<4; k++){
1762 for(int l=0; l<4; l++){
1763 if(i==j || i==k || i==l || j==k || j==l || k==l){
1764 m_epsilon_uvmn.push_back(0.0);
1765 }else{
1766 if(i==0 && j==1 && k==2 && l==3) m_epsilon_uvmn.push_back(1.0);
1767 if(i==0 && j==1 && k==3 && l==2) m_epsilon_uvmn.push_back(-1.0);
1768 if(i==0 && j==2 && k==1 && l==3) m_epsilon_uvmn.push_back(-1.0);
1769 if(i==0 && j==2 && k==3 && l==1) m_epsilon_uvmn.push_back(1.0);
1770 if(i==0 && j==3 && k==1 && l==2) m_epsilon_uvmn.push_back(1.0);
1771 if(i==0 && j==3 && k==2 && l==1) m_epsilon_uvmn.push_back(-1.0);
1772
1773 if(i==1 && j==0 && k==2 && l==3) m_epsilon_uvmn.push_back(-1.0);
1774 if(i==1 && j==0 && k==3 && l==2) m_epsilon_uvmn.push_back(1.0);
1775 if(i==1 && j==2 && k==0 && l==3) m_epsilon_uvmn.push_back(1.0);
1776 if(i==1 && j==2 && k==3 && l==0) m_epsilon_uvmn.push_back(-1.0);
1777 if(i==1 && j==3 && k==0 && l==2) m_epsilon_uvmn.push_back(-1.0);
1778 if(i==1 && j==3 && k==2 && l==0) m_epsilon_uvmn.push_back(1.0);
1779
1780 if(i==2 && j==0 && k==1 && l==3) m_epsilon_uvmn.push_back(1.0);
1781 if(i==2 && j==0 && k==3 && l==1) m_epsilon_uvmn.push_back(-1.0);
1782 if(i==2 && j==1 && k==0 && l==3) m_epsilon_uvmn.push_back(-1.0);
1783 if(i==2 && j==1 && k==3 && l==0) m_epsilon_uvmn.push_back(1.0);
1784 if(i==2 && j==3 && k==0 && l==1) m_epsilon_uvmn.push_back(1.0);
1785 if(i==2 && j==3 && k==1 && l==0) m_epsilon_uvmn.push_back(-1.0);
1786
1787 if(i==3 && j==0 && k==1 && l==2) m_epsilon_uvmn.push_back(-1.0);
1788 if(i==3 && j==0 && k==2 && l==1) m_epsilon_uvmn.push_back(1.0);
1789 if(i==3 && j==1 && k==0 && l==2) m_epsilon_uvmn.push_back(1.0);
1790 if(i==3 && j==1 && k==2 && l==0) m_epsilon_uvmn.push_back(-1.0);
1791 if(i==3 && j==2 && k==0 && l==1) m_epsilon_uvmn.push_back(-1.0);
1792 if(i==3 && j==2 && k==1 && l==0) m_epsilon_uvmn.push_back(1.0);
1793
1794 }
1795 }
1796 }
1797 }
1798 }
1799 // D0 -> a2(1320)+ {rho(770)0 pi+[D]} pi-
1800 double SF_Tp_D_VP_1 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3p1, T1_Pip2Pim1)), Pip1Pip2Pim1), contract_42_2(Proj2_3p1, T2_3p1)), T2_Pip2Pim1Pip1);
1801 double SF_Tp_D_VP_2 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3p1, T1_Pip1Pim1)), Pip1Pip2Pim1), contract_42_2(Proj2_3p1, T2_3p1)), T2_Pip1Pim1Pip2);
1802 double SF_Tp_D_VP_3 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3p2, T1_Pip2Pim2)), Pip1Pip2Pim2), contract_42_2(Proj2_3p2, T2_3p2)), T2_Pip2Pim2Pip1);
1803 double SF_Tp_D_VP_4 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3p2, T1_Pip1Pim2)), Pip1Pip2Pim2), contract_42_2(Proj2_3p2, T2_3p2)), T2_Pip1Pim2Pip2);
1804
1805 amplitude += g_fitpara[9]*(SF_Tp_D_VP_1*RBW_a21320p_1*GS_rho770_21 + SF_Tp_D_VP_2*RBW_a21320p_1*GS_rho770_11 + SF_Tp_D_VP_3*RBW_a21320p_2*GS_rho770_22 + SF_Tp_D_VP_4*RBW_a21320p_2*GS_rho770_12);
1806
1807 // D0 -> a2(1320)- {rho(770)0 pi-[D]} pi+
1808 double SF_Tm_D_VP_1 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3m1, T1_Pim2Pip1)), Pim1Pim2Pip1), contract_42_2(Proj2_3m1, T2_3m1)), T2_Pip1Pim2Pim1);
1809 double SF_Tm_D_VP_2 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3m1, T1_Pim1Pip1)), Pim1Pim2Pip1), contract_42_2(Proj2_3m1, T2_3m1)), T2_Pip1Pim1Pim2);
1810 double SF_Tm_D_VP_3 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3m2, T1_Pim2Pip2)), Pim1Pim2Pip2), contract_42_2(Proj2_3m2, T2_3m2)), T2_Pip2Pim2Pim1);
1811 double SF_Tm_D_VP_4 = contract_22_0(contract_22_2(contract_31_2(contract_41_3(m_epsilon_uvmn, contract_21_1(Proj1_3m2, T1_Pim1Pip2)), Pim1Pim2Pip2), contract_42_2(Proj2_3m2, T2_3m2)), T2_Pip2Pim1Pim2);
1812
1813 amplitude += g_fitpara[10]*(SF_Tm_D_VP_1*RBW_a21320m_1*GS_rho770_12 + SF_Tm_D_VP_2*RBW_a21320m_1*GS_rho770_11 + SF_Tm_D_VP_3*RBW_a21320m_2*GS_rho770_22 + SF_Tm_D_VP_4*RBW_a21320m_2*GS_rho770_21);
1814
1815 // D0 -> pi(1300)- {rho(770)0 pi-} pi+
1816 double SF_Pm_P_VP_1 = contract_11_0(T1_Pim2Pip1,T1_Pip1Pim2Pim1);
1817 double SF_Pm_P_VP_2 = contract_11_0(T1_Pim1Pip1,T1_Pip1Pim1Pim2);
1818 double SF_Pm_P_VP_3 = contract_11_0(T1_Pim2Pip2,T1_Pip2Pim2Pim1);
1819 double SF_Pm_P_VP_4 = contract_11_0(T1_Pim1Pip2,T1_Pip2Pim1Pim2);
1820
1821 amplitude += g_fitpara[11]*(SF_Pm_P_VP_1*GS_rho770_12*RBW_pi1300m_1 + SF_Pm_P_VP_2*GS_rho770_11*RBW_pi1300m_1 + SF_Pm_P_VP_3*GS_rho770_22*RBW_pi1300m_2 + SF_Pm_P_VP_4*GS_rho770_21*RBW_pi1300m_2);
1822/*
1823 // D0 -> pi(1300)- {f2(1270)0 pi-} pi+
1824 double SF_Pm_D_TP_1 = contract_22_0(T2_Pip1Pim2,T2_Pip1Pim2Pim1);
1825 double SF_Pm_D_TP_2 = contract_22_0(T2_Pip1Pim1,T2_Pip1Pim1Pim2);
1826 double SF_Pm_D_TP_3 = contract_22_0(T2_Pip2Pim2,T2_Pip2Pim2Pim1);
1827 double SF_Pm_D_TP_4 = contract_22_0(T2_Pip2Pim1,T2_Pip2Pim1Pim2);
1828
1829 amplitude += g_fitpara[12]*g_fitpara[11]*(SF_Pm_D_TP_1*RBW_f21270_12*RBW_pi1300m_1 + SF_Pm_D_TP_2*RBW_f21270_11*RBW_pi1300m_1 + SF_Pm_D_TP_3*RBW_f21270_22*RBW_pi1300m_2 + SF_Pm_D_TP_4*RBW_f21270_21*RBW_pi1300m_2);
1830*/
1831 // D0 -> pi(1300)- {f0 pi-} pi+
1832 amplitude += g_fitpara[12]*g_fitpara[11]*(RBW_pi1300m_1*PiPiS_12_0 + RBW_pi1300m_1*PiPiS_11_0 + RBW_pi1300m_2*PiPiS_22_0 + RBW_pi1300m_2*PiPiS_21_0);
1833// amplitude += g_fitpara[13]*g_fitpara[11]*(RBW_pi1300m_1*PiPiS_12_5 + RBW_pi1300m_1*PiPiS_11_5 + RBW_pi1300m_2*PiPiS_22_5 + RBW_pi1300m_2*PiPiS_21_5);
1834
1835 amplitude += g_fitpara[13]*g_fitpara[11]*(RBW_pi1300m_1*PiPiS_12_6 + RBW_pi1300m_1*PiPiS_11_6 + RBW_pi1300m_2*PiPiS_22_6 + RBW_pi1300m_2*PiPiS_21_6);
1836
1837 // D0 -> pi(1300)+ {rho(770)0 pi+} pi-
1838 double SF_Pp_P_VP_1 = contract_11_0(T1_Pip2Pim1,T1_Pip2Pim1Pip1);
1839 double SF_Pp_P_VP_2 = contract_11_0(T1_Pip1Pim1,T1_Pip1Pim1Pip2);
1840 double SF_Pp_P_VP_3 = contract_11_0(T1_Pip2Pim2,T1_Pip2Pim2Pip1);
1841 double SF_Pp_P_VP_4 = contract_11_0(T1_Pip1Pim2,T1_Pip1Pim2Pip2);
1842
1843 amplitude += g_fitpara[14]*(SF_Pp_P_VP_1*GS_rho770_21*RBW_pi1300p_1 + SF_Pp_P_VP_2*GS_rho770_11*RBW_pi1300p_1 + SF_Pp_P_VP_3*GS_rho770_22*RBW_pi1300p_2 + SF_Pp_P_VP_4*GS_rho770_12*RBW_pi1300p_2);
1844/*
1845 // D0 -> pi(1300)+ {f2(1270)0 pi+} pi-
1846 double SF_Pp_D_TP_1 = contract_22_0(T2_Pip2Pim1,T2_Pip2Pim1Pip1);
1847 double SF_Pp_D_TP_2 = contract_22_0(T2_Pip1Pim1,T2_Pip1Pim1Pip2);
1848 double SF_Pp_D_TP_3 = contract_22_0(T2_Pip2Pim2,T2_Pip2Pim2Pip1);
1849 double SF_Pp_D_TP_4 = contract_22_0(T2_Pip1Pim2,T2_Pip1Pim2Pip2);
1850
1851 amplitude += g_fitpara[12]*g_fitpara[15]*(SF_Pp_D_TP_1*RBW_f21270_21*RBW_pi1300p_1 + SF_Pp_D_TP_2*RBW_f21270_11*RBW_pi1300p_1 + SF_Pp_D_TP_3*RBW_f21270_22*RBW_pi1300p_2 + SF_Pp_D_TP_4*RBW_f21270_12*RBW_pi1300p_2);
1852*/
1853 // D0 -> pi(1300)+ {f0 pi+} pi-
1854 amplitude += g_fitpara[12]*g_fitpara[14]*(RBW_pi1300p_1*PiPiS_21_0 + RBW_pi1300p_1*PiPiS_11_0 + RBW_pi1300p_2*PiPiS_22_0 + RBW_pi1300p_2*PiPiS_12_0);
1855// amplitude += g_fitpara[13]*g_fitpara[15]*(RBW_pi1300p_1*PiPiS_21_5 + RBW_pi1300p_1*PiPiS_11_5 + RBW_pi1300p_2*PiPiS_22_5 + RBW_pi1300p_2*PiPiS_12_5);
1856
1857 amplitude += g_fitpara[13]*g_fitpara[14]*(RBW_pi1300p_1*PiPiS_21_6 + RBW_pi1300p_1*PiPiS_11_6 + RBW_pi1300p_2*PiPiS_22_6 + RBW_pi1300p_2*PiPiS_12_6);
1858
1859 // D0 -> rho rho [S]
1860 double SF_VV_S_1 = contract_11_0(T1_Pip1Pim1, T1_Pip2Pim2);
1861 double SF_VV_S_3 = contract_11_0(T1_Pip1Pim2, T1_Pip2Pim1);
1862
1863 amplitude += g_fitpara[15]*(SF_VV_S_1*GS_rho770_11*GS_rho770_22+SF_VV_S_3*GS_rho770_12*GS_rho770_21);
1864
1865 // D0 -> rho rho [P]
1866 double SF_VV_P_1 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(m_epsilon_uvmn, T1_Pip1Pim1),T1_Pip2Pim2),T1_2z11), D0);
1867 double SF_VV_P_2 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(m_epsilon_uvmn, T1_Pip2Pim2),T1_Pip1Pim1),T1_2z12), D0);
1868 double SF_VV_P_3 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(m_epsilon_uvmn, T1_Pip1Pim2),T1_Pip2Pim1),T1_2z21), D0);
1869 double SF_VV_P_4 = contract_11_0(contract_21_1(contract_31_2(contract_41_3(m_epsilon_uvmn, T1_Pip2Pim1),T1_Pip1Pim2),T1_2z22), D0);
1870
1871 amplitude += g_fitpara[16]*(SF_VV_P_1*GS_rho770_11*GS_rho770_22+SF_VV_P_3*GS_rho770_12*GS_rho770_21);
1872
1873 // D0 -> rho rho [D]
1874 double SF_VV_D_1 = contract_11_0(contract_21_1(T2_2z11,T1_Pip2Pim2), T1_Pip1Pim1);
1875 double SF_VV_D_3 = contract_11_0(contract_21_1(T2_2z21,T1_Pip2Pim1), T1_Pip1Pim2);
1876
1877 amplitude += g_fitpara[17]*(SF_VV_D_1*GS_rho770_11*GS_rho770_22+SF_VV_D_3*GS_rho770_12*GS_rho770_21);
1878
1879 // D0 -> rho rho' [P]
1880 amplitude += g_fitpara[18]*(SF_VV_P_1*GS_rho770_11*GS_rho1450_22+SF_VV_P_2*GS_rho770_22*GS_rho1450_11 + SF_VV_P_3*GS_rho770_12*GS_rho1450_21 + SF_VV_P_3*GS_rho770_21*GS_rho1450_12);
1881
1882 // D0 ->rho f0
1883 double SF_VS_P_1 = contract_11_0(T1_Pip1Pim1,T1_2z11);
1884 double SF_VS_P_2 = contract_11_0(T1_Pip2Pim2,T1_2z12);
1885 double SF_VS_P_3 = contract_11_0(T1_Pip1Pim2,T1_2z21);
1886 double SF_VS_P_4 = contract_11_0(T1_Pip2Pim1,T1_2z22);
1887
1888 amplitude += g_fitpara[19]*(SF_VS_P_1*GS_rho770_11*PiPiS_22_0 + SF_VS_P_2*GS_rho770_22*PiPiS_11_0 + SF_VS_P_3*GS_rho770_12*PiPiS_21_0 + SF_VS_P_4*GS_rho770_21*PiPiS_12_0);
1889 amplitude += g_fitpara[20]*(SF_VS_P_1*GS_rho770_11*PiPiS_22_5 + SF_VS_P_2*GS_rho770_22*PiPiS_11_5 + SF_VS_P_3*GS_rho770_12*PiPiS_21_5 + SF_VS_P_4*GS_rho770_21*PiPiS_12_5);
1890 amplitude += g_fitpara[21]*(SF_VS_P_1*GS_rho770_11*PiPiS_22_6 + SF_VS_P_2*GS_rho770_22*PiPiS_11_6 + SF_VS_P_3*GS_rho770_12*PiPiS_21_6 + SF_VS_P_4*GS_rho770_21*PiPiS_12_6);
1891
1892 // D0 -> f0f0
1893 //00
1894 amplitude += g_fitpara[22]*(PiPiS_11_0*PiPiS_22_0 + PiPiS_12_0*PiPiS_21_0 + PiPiS_22_0*PiPiS_11_0 + PiPiS_21_0*PiPiS_12_0);
1895 amplitude += g_fitpara[23]*(PiPiS_11_0*PiPiS_22_1 + PiPiS_12_0*PiPiS_21_1 + PiPiS_22_0*PiPiS_11_1 + PiPiS_21_0*PiPiS_12_1);
1896 amplitude += g_fitpara[24]*(PiPiS_11_1*PiPiS_22_1 + PiPiS_12_1*PiPiS_21_1 + PiPiS_22_1*PiPiS_11_1 + PiPiS_21_1*PiPiS_12_1);
1897 amplitude += g_fitpara[25]*(PiPiS_11_1*PiPiS_22_5 + PiPiS_12_1*PiPiS_21_5 + PiPiS_22_1*PiPiS_11_5 + PiPiS_21_1*PiPiS_12_5);
1898 amplitude += g_fitpara[26]*(PiPiS_11_5*PiPiS_22_5 + PiPiS_12_5*PiPiS_21_5 + PiPiS_22_5*PiPiS_11_5 + PiPiS_21_5*PiPiS_12_5);
1899 amplitude += g_fitpara[27]*(PiPiS_11_5*PiPiS_22_6 + PiPiS_12_5*PiPiS_21_6 + PiPiS_22_5*PiPiS_11_6 + PiPiS_21_5*PiPiS_12_6);
1900
1901 // D0 -> f2 f0
1902 double SF_TS_D_1 = contract_22_0(T2_Pip1Pim1,T2_2z11);
1903 double SF_TS_D_2 = contract_22_0(T2_Pip2Pim2,T2_2z12);
1904 double SF_TS_D_3 = contract_22_0(T2_Pip1Pim2,T2_2z21);
1905 double SF_TS_D_4 = contract_22_0(T2_Pip2Pim1,T2_2z22);
1906
1907 amplitude += g_fitpara[28]*(SF_TS_D_1*RBW_f21270_11*PiPiS_22_5 + SF_TS_D_2*RBW_f21270_22*PiPiS_11_5 + SF_TS_D_3*RBW_f21270_12*PiPiS_21_5 + SF_TS_D_4*RBW_f21270_21*PiPiS_12_5);
1908 amplitude += g_fitpara[29]*(SF_TS_D_1*RBW_f21270_11*PiPiS_22_6 + SF_TS_D_2*RBW_f21270_22*PiPiS_11_6 + SF_TS_D_3*RBW_f21270_12*PiPiS_21_6 + SF_TS_D_4*RBW_f21270_21*PiPiS_12_6);
1909
1910 return amplitude;
1911
1912}
double sin(const BesAngle a)
Definition BesAngle.h:210
double cos(const BesAngle a)
Definition BesAngle.h:213

Referenced by QCMCFilter::findD0Decay().

◆ init()

void D0To2pip2pim::init ( )

Definition at line 13 of file D0To2pip2pim.cxx.

13 {
14
15 //std::cout << "D0To2pip2pim ==> Initialization !" << std::endl;
16
17 g_mag[0]= 100.0; g_pha[0]= 0.0;
18 g_mag[1]= 7.95507 ; g_pha[1]= -0.0687407;
19 g_mag[2]= 37.5559 ; g_pha[2]= -1.74946 ;
20 g_mag[3]= 61.2172 ; g_pha[3]= 2.98079 ;
21 g_mag[4]= 187.79 ; g_pha[4]= 2.64471 ;
22 g_mag[5]= 385.474 ; g_pha[5]= -0.137107 ;
23 g_mag[6]= 0.330788; g_pha[6]= 0.268133 ;
24 g_mag[7]= 127.158 ; g_pha[7]= -2.47773 ;
25 g_mag[8]= 339.914 ; g_pha[8]= 2.22856 ;
26 g_mag[9]= 0.320888; g_pha[9]= -2.6194 ;
27 g_mag[10]=0.366283; g_pha[10]=-0.26867 ;
28 g_mag[11]=86.0865 ; g_pha[11]=-2.49649 ;
29 g_mag[12]=6.1541 ; g_pha[12]=-1.18299 ;
30 g_mag[13]=56.6067 ; g_pha[13]=0.142977 ;
31 g_mag[14]=92.3073 ; g_pha[14]=-2.15881 ;
32 g_mag[15]=10.5885 ; g_pha[15]=-3.03166 ;
33 g_mag[16]=8.36765 ; g_pha[16]=1.8417 ;
34 g_mag[17]=6.56437 ; g_pha[17]=-2.93087 ;
35 g_mag[18]=15.7197 ; g_pha[18]=0.96925 ;
36 g_mag[19]=21.4195 ; g_pha[19]=-1.23701 ;
37 g_mag[20]=56.8867 ; g_pha[20]=-0.385837 ;
38 g_mag[21]=231.626 ; g_pha[21]=2.14842 ;
39 g_mag[22]=2938.45 ; g_pha[22]=-0.693491 ;
40 g_mag[23]=7252.7 ; g_pha[23]=2.23659 ;
41 g_mag[24]=5165.87 ; g_pha[24]=0.913557 ;
42 g_mag[25]=11508.6 ; g_pha[25]=-1.07187 ;
43 g_mag[26]=2461.86 ; g_pha[26]=1.8709 ;
44 g_mag[27]=8757.75 ; g_pha[27]=2.40756 ;
45 g_mag[28]=19.7413 ; g_pha[28]=-1.0753 ;
46 g_mag[29]=66.3826 ; g_pha[29]=2.34666 ;
47 fitpara.clear();
48 for(int i=0; i<30; i++){
49 complex<double> ctemp(g_mag[i]*cos(g_pha[i]),g_mag[i]*sin(g_pha[i]));
50 fitpara.push_back(ctemp);
51 }
52
53 //g_uv.clear();
54 //for(int i=0; i<4; i++){
55 // for(int j=0; j<4; j++){
56 // if(i!=j){
57 // g_uv.push_back(0.0);
58 // }else if(i<3){
59 // g_uv.push_back(-1.0);
60 // }else if(i==3){
61 // g_uv.push_back(1.0);
62 // }
63 // }
64 //}
65
66 epsilon_uvmn.clear();
67 for(int i=0; i<4; i++){
68 for(int j=0; j<4; j++){
69 for(int k=0; k<4; k++){
70 for(int l=0; l<4; l++){
71 if(i==j || i==k || i==l || j==k || j==l || k==l){
72 epsilon_uvmn.push_back(0.0);
73 }else{
74 if(i==0 && j==1 && k==2 && l==3) epsilon_uvmn.push_back(1.0);
75 if(i==0 && j==1 && k==3 && l==2) epsilon_uvmn.push_back(-1.0);
76 if(i==0 && j==2 && k==1 && l==3) epsilon_uvmn.push_back(-1.0);
77 if(i==0 && j==2 && k==3 && l==1) epsilon_uvmn.push_back(1.0);
78 if(i==0 && j==3 && k==1 && l==2) epsilon_uvmn.push_back(1.0);
79 if(i==0 && j==3 && k==2 && l==1) epsilon_uvmn.push_back(-1.0);
80
81 if(i==1 && j==0 && k==2 && l==3) epsilon_uvmn.push_back(-1.0);
82 if(i==1 && j==0 && k==3 && l==2) epsilon_uvmn.push_back(1.0);
83 if(i==1 && j==2 && k==0 && l==3) epsilon_uvmn.push_back(1.0);
84 if(i==1 && j==2 && k==3 && l==0) epsilon_uvmn.push_back(-1.0);
85 if(i==1 && j==3 && k==0 && l==2) epsilon_uvmn.push_back(-1.0);
86 if(i==1 && j==3 && k==2 && l==0) epsilon_uvmn.push_back(1.0);
87
88 if(i==2 && j==0 && k==1 && l==3) epsilon_uvmn.push_back(1.0);
89 if(i==2 && j==0 && k==3 && l==1) epsilon_uvmn.push_back(-1.0);
90 if(i==2 && j==1 && k==0 && l==3) epsilon_uvmn.push_back(-1.0);
91 if(i==2 && j==1 && k==3 && l==0) epsilon_uvmn.push_back(1.0);
92 if(i==2 && j==3 && k==0 && l==1) epsilon_uvmn.push_back(1.0);
93 if(i==2 && j==3 && k==1 && l==0) epsilon_uvmn.push_back(-1.0);
94
95 if(i==3 && j==0 && k==1 && l==2) epsilon_uvmn.push_back(-1.0);
96 if(i==3 && j==0 && k==2 && l==1) epsilon_uvmn.push_back(1.0);
97 if(i==3 && j==1 && k==0 && l==2) epsilon_uvmn.push_back(1.0);
98 if(i==3 && j==1 && k==2 && l==0) epsilon_uvmn.push_back(-1.0);
99 if(i==3 && j==2 && k==0 && l==1) epsilon_uvmn.push_back(-1.0);
100 if(i==3 && j==2 && k==1 && l==0) epsilon_uvmn.push_back(1.0);
101
102 }
103 }
104 }
105 }
106 }
107
108 _nd = 4;
109 math_pi = 3.1415926;
110 mass_Pion = 0.13957;
111
112 rRes = 3.0*0.197321;
113 rD = 5.0*0.197321;
114 m_Pi = 0.139570;
115 m2_Pi = m_Pi*m_Pi;
116
117 m0_rho770 = 0.77526;
118 w0_rho770 = 0.1478;
119
120 m0_rho1450 = 1.465;
121 w0_rho1450 = 0.400;
122
123 m0_f21270 = 1.2755;
124 w0_f21270 = 0.1867;
125
126 m0_a11260 = 1.1337;
127 g1_a11260 = 0.00335;
128 g2_a11260 = 0.0;
129
130 m0_pi1300 = 1.498;
131 w0_pi1300 = 0.590;
132
133 m0_a11420 = 1.411;
134 w0_a11420 = 0.161;
135
136 m0_a11640 = 1.655;
137 w0_a11640 = 0.254;
138
139 m0_a21320 = 1.3186;
140 w0_a21320 = 0.105;
141
142 m0_pi11400 = 1.354;
143 w0_pi11400 = 0.330;
144
145 s0_prod = -5.0;
146
147 return;
148}

Referenced by QCMCFilter::findD0Decay().


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