5 #ifndef COMPWA_PHYSICS_DYNAMICS_FLATTE_HPP_ 6 #define COMPWA_PHYSICS_DYNAMICS_FLATTE_HPP_ 19 std::shared_ptr<ComPWA::FunctionTree::FitParameter>
G;
25 inline std::complex<double>
27 double massB,
unsigned int J,
double mesonRadius,
29 auto qR = std::abs(
qValueAC(mR, massA, massB));
31 auto ffR = FormFactorFunctor->operator()(qR *qR, J, mesonRadius);
32 auto qS = std::abs(
qValueAC(sqrtS, massA, massB));
33 auto barrierA = FormFactorFunctor->operator()(qS *qS, J, mesonRadius) / ffR;
36 std::complex<double> vtxA(1, 0);
37 if (J > 0 || FormFactorFunctor->getName() ==
"CrystalBarrel") {
38 vtxA = ffR * std::pow(qR, J);
47 return (width * barrierA * barrierA);
60 inline std::complex<double>
62 std::complex<double> termB,
63 std::complex<double> termC = std::complex<double>(0, 0)) {
64 std::complex<double> i(0, 1);
65 double sqrtS = sqrt(mSq);
67 std::complex<double> denom = std::complex<double>(mR * mR - mSq, 0) +
68 (-1.0) * i * sqrtS * (termA + termB + termC);
70 std::complex<double> result = std::complex<double>(gA, 0) / denom;
73 if (std::isnan(result.real()) || std::isnan(result.imag())) {
74 std::cout <<
"AmpFlatteRes::dynamicalFunction() | " << mR <<
" " << mSq
75 <<
" " << termA <<
" " << termB <<
" " << termC << std::endl;
101 inline std::complex<double>
103 double gA,
double massB1,
double massB2,
double couplingB,
104 double massC1,
double massC2,
double couplingC,
105 unsigned int L,
double mesonRadius,
107 double sqrtS = sqrt(mSq);
114 mesonRadius, FormFactorFunctor);
117 std::complex<double> termC;
118 if (couplingC != 0.0) {
120 mesonRadius, FormFactorFunctor);
139 std::shared_ptr<ComPWA::FunctionTree::Parameter> &out);
std::complex< double > phspFactorAC(double sqrtS, double ma, double mb)
Two body phsp factor.
std::shared_ptr< FormFactor > FormFactorFunctor
std::complex< double > dynamicalFunction(double mSq, double mR, double gA, std::complex< double > termA, std::complex< double > termB, std::complex< double > termC=std::complex< double >(0, 0))
Dynamical function for two coupled channel approach.
FlatteStrategy(std::shared_ptr< FormFactor > FF)
std::shared_ptr< ComPWA::FunctionTree::TreeNode > createFunctionTree(InputInfo Params, std::shared_ptr< ComPWA::FunctionTree::Value< std::vector< double >>> InvMassSquared)
std::complex< double > couplingToWidth(double mR, double g, std::complex< double > gamma, std::complex< double > phspFactor)
Convert width to complex coupling.
std::complex< double > qValueAC(double sqrtS, double ma, double mb)
Calculate Break-up momentum.
std::complex< double > flatteCouplingTerm(double sqrtS, double mR, double coupling, double massA, double massB, unsigned int J, double mesonRadius, std::shared_ptr< FormFactor > FormFactorFunctor)
Helper function to calculate the coupling terms for the Flatte formular.
ParType
Enums for the type of the parameter, should be extended if an new parameter type is added...
Virtual base class for operations of FunctionTree nodes.
This class provides a list of parameters and values of different types.