MINT2
SBW.h
Go to the documentation of this file.
1 #ifndef SBW_LINESHAPE_HH
2 #define SBW_LINESHAPE_HH
3 // author: Philippe d'Argent (p.dargent@cern.ch)
4 // status: 19 March 2015 GMT
5 
6 #include <complex>
7 #include <string>
8 
9 #include "Mint/ILineshape.h"
10 #include "Mint/BW_BW.h"
11 #include "Mint/NamedParameter.h"
12 #include "Mint/FitParDependent.h"
15 
16 class SBW : public BW_BW, virtual public ILineshape{
17  public:
18 
19  SBW( const AssociatedDecayTree& tree, const std::string& namePrefix): BW_BW(tree, namePrefix){}
20 
21  virtual std::string name() const{
22  return "SBW("+prefix()+_theDecay.oneLiner() +")";
23  }
24 
25  virtual ~SBW(){}
26 
27  protected:
28 
29  virtual double GofM() {return mumsWidth();}
30  virtual std::complex<double> BreitWigner(){
31 
32  double mass = mumsMass();
33  double width = mumsWidth();
34 
35  double gamma = sqrt(mass*mass*(mass*mass+width*width));
36  double k = mass*width*gamma/sqrt(mass*mass+gamma);
37 
38  const double m2hh = mumsRecoMass2()/GeV/GeV;
39  //double p = twoBody_dgtPsq_in_MumsFrame(mumsRecoMass(), daughterPDGMass(0), daughterPDGMass(1));
40  //if(p <= 0) return 0.;
41 
42  std::complex<double> invBW(mumsRecoMass()-mass, - width/2.);
43  //return sqrt(k)/invBW;
44  return 1.*GeV/invBW;//*pow(pABSq(),GetAlpha());//*(1.+pABSq()/GeV*c1()+pow(pABSq()/GeV,2.)*c2()+pow(pABSq()/GeV,3.)*c3());
45  }
46 
47 // double GetAlpha() const{
48 // return _RPL->get(mumsPID())->alpha();
49 // }
50  double c1() const {
51  return _RPL->get(mumsPID())->c1();
52  }
53  double c2() const {
54  return _RPL->get(mumsPID())->c2();
55  }
56  double c3() const {
57  return _RPL->get(mumsPID())->c3();
58  }
59  double c4() const {
60  return _RPL->get(mumsPID())->c4();
61  }
62  double c5() const {
63  return _RPL->get(mumsPID())->c5();
64  }
65 
66 };
67 
68 #endif
69 //
ResonancePropertiesList * _RPL
Definition: BW_BW.h:97
SBW(const AssociatedDecayTree &tree, const std::string &namePrefix)
Definition: SBW.h:19
const AssociatedDecayTree & _theDecay
Definition: BW_BW.h:69
Definition: SBW.h:16
virtual double mumsMass() const
Definition: BW_BW.cpp:453
virtual std::string name() const
Definition: SBW.h:21
double c5() const
Definition: SBW.h:62
virtual int mumsPID() const
Definition: BW_BW.cpp:517
double c4() const
Definition: SBW.h:59
const ResonanceProperties * get(int i) const
Definition: BW_BW.h:30
virtual std::complex< double > BreitWigner()
Definition: SBW.h:30
double c3() const
Definition: SBW.h:56
double c1() const
Definition: SBW.h:50
virtual double mumsRecoMass2() const
Definition: BW_BW.cpp:548
double c2() const
Definition: SBW.h:53
static const double GeV
virtual double mumsRecoMass() const
Definition: BW_BW.cpp:631
virtual double mumsWidth() const
Definition: BW_BW.cpp:460
void oneLiner(std::stringstream &seam, int generation=0) const
Definition: DDTree.h:375
virtual double GofM()
Definition: SBW.h:29
const std::string & prefix() const
Definition: BW_BW.h:72
virtual ~SBW()
Definition: SBW.h:25