MINT2
Public Member Functions | Protected Member Functions | List of all members
GLass Class Reference

#include <GLass.h>

Inheritance diagram for GLass:
BW_BW ILineshape ILineshape MINT::FitParDependent MINT::IFitParDependent MINT::IFitParDependent MINT::IFitParRegister MINT::PolymorphVector< FitParRef > MINT::IFitParDependent

Public Member Functions

 GLass (const AssociatedDecayTree &tree, const std::string &namePrefix="")
 
virtual std::string name () const
 
virtual ~GLass ()
 
- Public Member Functions inherited from BW_BW
virtual double prSq () const
 
virtual double prSqForGofM () const
 
virtual double pABSq ()
 
virtual double prSqMax () const
 
 BW_BW (const AssociatedDecayTree &decay, const std::string &lineshapePrefix="", MINT::MinuitParameterSet *mps=0)
 
 BW_BW (const BW_BW &other)
 
virtual ~BW_BW ()
 
virtual std::complex< double > getVal (IDalitzEvent &evt)
 
virtual void print (IDalitzEvent &evt, std::ostream &out=std::cout)
 
virtual void print (std::ostream &out=std::cout) const
 
virtual DalitzCoordinate getDalitzCoordinate (double nSigma=3) const
 
virtual MINT::counted_ptr< IGenFctgeneratingFunction () const
 
virtual void setGenerationLimits (double mi, double ma)
 
virtual int twoLPlusOne () const
 
virtual std::complex< double > EvtGenValue (IDalitzEvent &evt)
 
- Public Member Functions inherited from ILineshape
virtual ~ILineshape ()
 
- Public Member Functions inherited from MINT::FitParDependent
virtual unsigned int size () const
 
virtual const FitParRefoperator[] (unsigned int i) const
 
virtual FitParRefoperator[] (unsigned int i)
 
virtual bool changedSinceLastCall () const
 
virtual void rememberFitParValues ()
 
virtual bool registerFitParDependence (const IFitParDependent &fpd)
 
bool registerFitParDependence (const FitParRef &fpr)
 
void removeAllFitParDependencies ()
 
 FitParDependent (IFitParRegister *daddy=0)
 
 FitParDependent (const FitParDependent &other, IFitParRegister *newDaddy=0)
 
void listFitParDependencies (std::ostream &os=std::cout) const
 
- Public Member Functions inherited from MINT::PolymorphVector< FitParRef >
 PolymorphVector ()
 
 PolymorphVector (unsigned int N)
 
 PolymorphVector (unsigned int N, const FitParRef &c)
 
 PolymorphVector (const PolymorphVector &other)
 
 PolymorphVector (const typename std::vector< FitParRef > &other)
 
virtual ~PolymorphVector ()
 
std::vector< FitParRef > & theVector ()
 
const std::vector< FitParRef > & theVector () const
 
FitParRefoperator[] (unsigned int i)
 
const FitParRefoperator[] (unsigned int i) const
 
FitParRefat (unsigned int i)
 
const FitParRefat (unsigned int i) const
 
std::vector< FitParRef >::iterator begin ()
 
std::vector< FitParRef >::const_iterator begin () const
 
std::vector< FitParRef >::iterator end ()
 
std::vector< FitParRef >::const_iterator end () const
 
std::vector< FitParRef >::iterator find (const FitParRef &c)
 
std::vector< FitParRef >::const_iterator find (const FitParRef &c) const
 
FitParReffront ()
 
const FitParReffront () const
 
FitParRefback ()
 
const FitParRefback () const
 
unsigned int size () const
 
bool empty () const
 
void push_back (const FitParRef &c)
 
void pop_back ()
 
void erase (typename std::vector< FitParRef >::iterator pos)
 
void erase (typename std::vector< FitParRef >::iterator first, typename std::vector< FitParRef >::iterator last)
 
PolymorphVector< FitParRef > & operator= (const PolymorphVector< FitParRef > &other)
 
void clear ()
 
void resize (unsigned int N)
 
void resize (unsigned int N, const FitParRef &c)
 
 operator const typename std::vector< FitParRef > & () const
 
 operator typename std::vector< FitParRef > & ()
 
bool operator== (const MINT::PolymorphVector< FitParRef > &v2) const
 
bool operator!= (const MINT::PolymorphVector< FitParRef > &v2) const
 
bool operator< (const MINT::PolymorphVector< FitParRef > &v2) const
 
bool operator> (const MINT::PolymorphVector< FitParRef > &v2) const
 

Protected Member Functions

virtual std::complex< double > BreitWigner ()
 
double a () const
 
double r () const
 
double R () const
 
double phiR () const
 
double F () const
 
double phiF () const
 
double alpha1 () const
 
double alpha2 () const
 
double alpha3 () const
 
- Protected Member Functions inherited from BW_BW
const MINT::MinuitParameterSetgetMinuitParameterSet () const
 
MINT::MinuitParameterSetgetMinuitParameterSet ()
 
const std::string & prefix () const
 
bool substitutePDGForReco () const
 
virtual double twoBody_dgtPsq_in_MumsFrame (double mumsMass, double mA, double mB) const
 
virtual double twoBody_dgtPsq_in_MumsPDGFrame () const
 
virtual double twoBody_recodgtPsq_in_MumsPDGFrame () const
 
virtual double twoBody_dgtPsq_in_MumsRecoFrame ()
 
bool setEventPtr (IDalitzEvent &evt) const
 
IDalitzEventgetEvent () const
 
bool setAllFitParameters ()
 
ResonancePropertiesListresonancePropertiesList () const
 
const ResonancePropertiesresonanceProperties () const
 
virtual const ParticlePropertiesmumsProperties () const
 
virtual ResonancePropertiesFitRefmumsFittableProperties () const
 
virtual int mumsPID () const
 
virtual double mumsMass () const
 
virtual double mumsWidth () const
 
virtual double mumsRadius () const
 
virtual double globalRadius () const
 
virtual double Radius () const
 
double GetAlpha () const
 
virtual int lowestPossibleTwoLPlusOne () const
 
virtual int maxDaughterSpinSum () const
 
virtual int minDaughterSpinSum () const
 
virtual int minDaughterSpinSum2 () const
 
virtual int minDaughterSpinSum3 () const
 
virtual int minDaughterSpinSum4 () const
 
virtual int maxDaughterPairSpinSum (int i, int j) const
 
virtual int minDaughterPairSpinSum (int i, int j) const
 
virtual std::string mumsSpin () const
 
virtual int mumsSpinValue () const
 
virtual int mumsParity () const
 
virtual bool startOfDecayChain () const
 
virtual bool parityConservingL (int L) const
 
virtual double mumsRecoMass2 () const
 
virtual double mumsRecoMass () const
 
virtual MultiQuarkContent mumsQuarkContent () const
 
bool nonResonant () const
 
virtual int numDaughters () const
 
virtual TLorentzVector daughterP4 (int i) const
 
virtual double daughterPDGMass (const int &i) const
 
virtual double daughterWidth (int i) const
 
virtual double daughterRecoMass2 (int i) const
 
virtual double daughterRecoMass (int i) const
 
virtual std::string daughterSpin (int i) const
 
virtual int daughterSpinValue (int i) const
 
virtual MultiQuarkContent daughterQuarkContent (int i) const
 
virtual int daughterP (int i) const
 
virtual int dgtrsInternalParity () const
 
virtual MultiQuarkContent dgtrsQuarkContent () const
 
virtual bool isWeakDecay () const
 
virtual double Fr_PDG_BL ()
 
virtual double Fr_BELLE (double prSquared)
 
virtual double Fr_BELLE_Max ()
 
virtual double Fr ()
 
virtual double FrForGofM ()
 
virtual double FrMax ()
 
virtual double GofM ()
 
virtual void resetInternals ()
 
virtual void resetPDG ()
 

Additional Inherited Members

- Protected Attributes inherited from BW_BW
MINT::MinuitParameterSet_mps
 
std::string _prefix
 
MINT::NamedParameter< int > _normBF
 
bool _useGlobalRadius
 
const AssociatedDecayTree_theDecay
 
double _gen_s_mi
 
double _gen_s_ma
 
ResonancePropertiesList_RPL
 
ResonancePropertiesFitRef_fittableResonancePropertiesPtr
 
MINT::FitParRef_fittableGlobalRadiusPtr
 
- Protected Attributes inherited from MINT::PolymorphVector< FitParRef >
std::vector< FitParRef_vec
 

Detailed Description

Definition at line 23 of file GLass.h.

Constructor & Destructor Documentation

◆ GLass()

GLass::GLass ( const AssociatedDecayTree tree,
const std::string &  namePrefix = "" 
)
inline

Definition at line 40 of file GLass.h.

40  :
41  BW_BW(tree, namePrefix)
42  {}
BW_BW(const AssociatedDecayTree &decay, const std::string &lineshapePrefix="", MINT::MinuitParameterSet *mps=0)
Definition: BW_BW.cpp:26

◆ ~GLass()

virtual GLass::~GLass ( )
inlinevirtual

Definition at line 47 of file GLass.h.

47 {}

Member Function Documentation

◆ a()

double GLass::a ( ) const
inlineprotected

Definition at line 28 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414
static const double GeV

◆ alpha1()

double GLass::alpha1 ( ) const
inlineprotected

Definition at line 34 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ alpha2()

double GLass::alpha2 ( ) const
inlineprotected

Definition at line 35 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ alpha3()

double GLass::alpha3 ( ) const
inlineprotected

Definition at line 36 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ BreitWigner()

std::complex< double > GLass::BreitWigner ( )
protectedvirtual

Reimplemented from BW_BW.

Definition at line 9 of file GLass.cpp.

9  {
10 
11  double q2 = pABSq();
12  double q = sqrt( q2 );
13  double gamma = GofM();
14  double scatteringPhase = phiF() + atan( 2. * a() * q / ( 2. + a() * r() * q2 ) );
15  double resonancePhase = phiR() + atan( mumsMass() * gamma / ( mumsMass() * mumsMass() - mumsRecoMass2() ) );
16  double rho = sqrt( q2 / s );
17  complex<double> returnValue = ( F() * sin( scatteringPhase ) * std::polar(1.,scatteringPhase) +
18  R() * sin( resonancePhase ) * std::polar(1., resonancePhase + 2. * scatteringPhase ) ) / rho;
19 
20  std::vector<int> asi = _theDecay.getVal().asi();
21  double min = getEvent()->eventPattern().sijMin(asi);
22  double max = getEvent()->eventPattern().sijMax(asi);
23  double x = 2.* (mumsRecoMass2() - min)/(max - min) - 1.;
24 
25  returnValue *= exp(alpha1()*x) + exp(alpha2()*x*x) + exp(alpha3()*x*x*x);
26 
27  return returnValue;
28 }
double sijMax(const MINT::PolymorphVector< int > &indices) const
const AssociatedDecayTree & _theDecay
Definition: BW_BW.h:69
virtual double mumsMass() const
Definition: BW_BW.cpp:453
static const double s
const ValueType & getVal() const
Definition: DDTree.h:102
double alpha1() const
Definition: GLass.h:34
double F() const
Definition: GLass.h:32
double r() const
Definition: GLass.h:29
IDalitzEvent * getEvent() const
Definition: BW_BW.cpp:147
double sijMin(const MINT::PolymorphVector< int > &indices) const
virtual double mumsRecoMass2() const
Definition: BW_BW.cpp:548
virtual const DalitzEventPattern & eventPattern() const =0
virtual double pABSq()
Definition: BW_BW.cpp:1014
virtual double GofM()
Definition: BW_BW.cpp:1035
double phiR() const
Definition: GLass.h:31
double alpha3() const
Definition: GLass.h:36
const std::vector< int > & asi() const
double alpha2() const
Definition: GLass.h:35
double a() const
Definition: GLass.h:28
double R() const
Definition: GLass.h:30
double phiF() const
Definition: GLass.h:33

◆ F()

double GLass::F ( ) const
inlineprotected

Definition at line 32 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ name()

virtual std::string GLass::name ( ) const
inlinevirtual

Reimplemented from BW_BW.

Definition at line 44 of file GLass.h.

44  {
45  return "GLASS("+prefix()+_theDecay.oneLiner() +")";
46  }
const AssociatedDecayTree & _theDecay
Definition: BW_BW.h:69
void oneLiner(std::stringstream &seam, int generation=0) const
Definition: DDTree.h:375
const std::string & prefix() const
Definition: BW_BW.h:72

◆ phiF()

double GLass::phiF ( ) const
inlineprotected

Definition at line 33 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ phiR()

double GLass::phiR ( ) const
inlineprotected

Definition at line 31 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ r()

double GLass::r ( ) const
inlineprotected

Definition at line 29 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414
static const double GeV

◆ R()

double GLass::R ( ) const
inlineprotected

Definition at line 30 of file GLass.h.

virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

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