Skip to content
Snippets Groups Projects
Commit 2f857715 authored by Julian Pychy's avatar Julian Pychy
Browse files

fixed phasespace factor for pipi-S wave

parent 06f4ef78
No related branches found
No related tags found
No related merge requests found
......@@ -24,7 +24,7 @@
#include "PwaDynamics/KMatrixPiPiS.hh"
#include "PwaDynamics/KPole.hh"
#include "PwaDynamics/AbsPhaseSpace.hh"
#include "PwaDynamics/PhaseSpaceIsobar.hh"
#include "PwaDynamics/PhaseSpaceIsobarAS.hh"
#include "PwaDynamics/PhaseSpace4Pi.hh"
#include "qft++/relativistic-quantum-mechanics/Utils.hh"
#include "qft++/matrix/IdentityMatrix.hh"
......@@ -88,11 +88,11 @@ void KMatrixPiPiS::initASParam1900(){
gFactorsMap[4].push_back(-0.00355); //eta eta
gFactorsMap[4].push_back(0.22358); //eta eta'
std::shared_ptr<AbsPhaseSpace> pipiPhp(new PhaseSpaceIsobar(piMass, piMass));
std::shared_ptr<AbsPhaseSpace> kkPhp(new PhaseSpaceIsobar(KplusMass, K0Mass));
std::shared_ptr<AbsPhaseSpace> pipiPhp(new PhaseSpaceIsobarAS(piMass, piMass));
std::shared_ptr<AbsPhaseSpace> kkPhp(new PhaseSpaceIsobarAS(KplusMass, K0Mass));
std::shared_ptr<AbsPhaseSpace> pipipipiPhp(new PhaseSpace4Pi());
std::shared_ptr<AbsPhaseSpace> etaetaPhp(new PhaseSpaceIsobar(etaMass, etaMass));
std::shared_ptr<AbsPhaseSpace> etaetapPhp(new PhaseSpaceIsobar(etaMass, etaprimeMass));
std::shared_ptr<AbsPhaseSpace> etaetaPhp(new PhaseSpaceIsobarAS(etaMass, etaMass));
std::shared_ptr<AbsPhaseSpace> etaetapPhp(new PhaseSpaceIsobarAS(etaMass, etaprimeMass));
_phpVecs.push_back(pipiPhp);
_phpVecs.push_back(kkPhp);
......
//************************************************************************//
// //
// Copyright 2014 Bertram Kopf (bertram@ep1.rub.de) //
// Julian Pychy (julian@ep1.rub.de) //
// - Ruhr-Universität Bochum //
// //
// This file is part of Pawian. //
// //
// Pawian is free software: you can redistribute it and/or modify //
// it under the terms of the GNU General Public License as published by //
// the Free Software Foundation, either version 3 of the License, or //
// (at your option) any later version. //
// //
// Pawian is distributed in the hope that it will be useful, //
// but WITHOUT ANY WARRANTY; without even the implied warranty of //
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
// GNU General Public License for more details. //
// //
// You should have received a copy of the GNU General Public License //
// along with Pawian. If not, see <http://www.gnu.org/licenses/>. //
// //
//************************************************************************//
#include "PwaDynamics/PhaseSpaceIsobarAS.hh"
#include "qft++/relativistic-quantum-mechanics/Utils.hh"
PhaseSpaceIsobarAS::PhaseSpaceIsobarAS(double mass1, double mass2):
AbsPhaseSpace()
, _mass1(mass1)
,_mass2(mass2)
{
}
PhaseSpaceIsobarAS::~PhaseSpaceIsobarAS(){
}
complex<double> PhaseSpaceIsobarAS::factor(const double mass){
return phaseSpaceFacAS(mass,_mass1, _mass2);
}
complex<double> PhaseSpaceIsobarAS::breakUpMom(const double mass){
return breakupMomQAS(mass,_mass1, _mass2);
}
complex<double> PhaseSpaceIsobarAS::factor(const complex<double> mass){
// Calc from the breakup momentum to account for chosen sign
complex<double> q = breakupMomQAS(mass,_mass1, _mass2);
CorrectForChosenSign(q);
return q * 2. / mass;
}
complex<double> PhaseSpaceIsobarAS::breakUpMom(const complex<double> mass){
complex<double> q = breakupMomQAS(mass,_mass1, _mass2);
CorrectForChosenSign(q);
return q;
}
//************************************************************************//
// //
// Copyright 2014 Bertram Kopf (bertram@ep1.rub.de) //
// Julian Pychy (julian@ep1.rub.de) //
// - Ruhr-Universität Bochum //
// //
// This file is part of Pawian. //
// //
// Pawian is free software: you can redistribute it and/or modify //
// it under the terms of the GNU General Public License as published by //
// the Free Software Foundation, either version 3 of the License, or //
// (at your option) any later version. //
// //
// Pawian is distributed in the hope that it will be useful, //
// but WITHOUT ANY WARRANTY; without even the implied warranty of //
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
// GNU General Public License for more details. //
// //
// You should have received a copy of the GNU General Public License //
// along with Pawian. If not, see <http://www.gnu.org/licenses/>. //
// //
//************************************************************************//
// PhaseSpaceIsobarAS class definition file. -*- C++ -*-
// Copyright 2014 Julian Pychy
#pragma once
//_____________________________________________________________________________
// @file PhaseSpaceIsobarAS.h
//_____________________________________________________________________________
#include <iostream>
#include <vector>
#include "PwaDynamics/AbsPhaseSpace.hh"
using namespace std;
//_____________________________________________________________________________
//_____________________________________________________________________________
class PhaseSpaceIsobarAS: public AbsPhaseSpace {
public:
/// Constructor
PhaseSpaceIsobarAS(double mass1, double mass2);
/// Destructor
virtual ~PhaseSpaceIsobarAS();
// operators:
// functions:
virtual complex<double> factor(const double mass);
virtual complex<double> breakUpMom(const double mass);
virtual complex<double> factor(const complex<double> mass);
virtual complex<double> breakUpMom(const complex<double> mass);
protected:
private:
double _mass1;
double _mass2;
};
//_____________________________________________________________________________
......@@ -579,6 +579,21 @@ complex<double> phaseSpaceFac(MassType mass, double massDec1, double massDec2){
return result;
}
template<typename MassType>
complex<double> phaseSpaceFacAS(MassType mass, double massDec1, double massDec2){
complex<double> result(0.,0.);
if(fabs(mass) < 1e-8) {
std::cout << "mass " << mass << " very close to 0; not possible to calculate phasespace factor " <<std::endl;
}
MassType termPlus=(massDec1+massDec2)/mass;
MassType tmpVal= 1. - termPlus*termPlus;
result = std::sqrt(std::complex<double>(tmpVal));
return result;
}
template<typename T>
complex<double> breakupMomQ(T mass, double massDec1, double massDec2){
......@@ -587,9 +602,22 @@ complex<double> breakupMomQ(T mass, double massDec1, double massDec2){
return result;
}
template<typename T>
complex<double> breakupMomQAS(T mass, double massDec1, double massDec2){
complex<double> result=phaseSpaceFacAS(mass, massDec1, massDec2)*mass/2.;
return result;
}
template complex<double> phaseSpaceFac(double, double, double);
template complex<double> phaseSpaceFac(complex<double> , double, double);
template complex<double> phaseSpaceFacAS(double, double, double);
template complex<double> phaseSpaceFacAS(complex<double> , double, double);
template complex<double> breakupMomQ(double, double, double);
template complex<double> breakupMomQ(complex<double>, double, double);
template complex<double> breakupMomQAS(double, double, double);
template complex<double> breakupMomQAS(complex<double>, double, double);
......@@ -158,5 +158,9 @@ complex<double> phaseSpaceFac(MassType mass, double massDec1, double massDec2);
template<typename MassType>
complex<double> breakupMomQ(MassType mass, double massDec1, double massDec2);
template<typename MassType>
complex<double> phaseSpaceFacAS(MassType mass, double massDec1, double massDec2);
template<typename MassType>
complex<double> breakupMomQAS(MassType mass, double massDec1, double massDec2);
#endif /* _Utils_H */
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment