laura is hosted by Hepforge, IPPP Durham
Laura++  v2r1
A maximum likelihood fitting package for performing Dalitz-plot analysis.

Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc.) More...

#include <LauAbsResonance.hh>

Inheritance diagram for LauAbsResonance:
LauBelleNR LauBelleSymNR LauBreitWignerRes LauDabbaRes LauFlatNR LauFlatteRes LauGounarisSakuraiRes LauKappaRes LauKMatrixProdPole LauKMatrixProdSVP LauLASSBWRes LauLASSNRRes LauLASSRes LauNRAmplitude LauPolNR LauRelBreitWignerRes LauSigmaRes

Public Types

enum  LauResonanceModel {
  BW, RelBW, GS, Flatte,
  Sigma, Kappa, Dabba, LASS,
  LASS_BW, LASS_NR, KMatrix, FlatNR,
  NRModel, BelleNR, PolNR
}
 Define the allowed resonance types. More...
 
enum  BarrierType { BWBarrier, BWPrimeBarrier, ExpBarrier }
 Define the allowed types of barrier factors. More...
 

Public Member Functions

 LauAbsResonance (const TString &resName, Double_t resMass, Double_t resWidth, Int_t resSpin, Int_t resCharge, Int_t resPairAmpInt, const LauDaughters *daughters)
 Constructor. More...
 
virtual ~LauAbsResonance ()
 Destructor. More...
 
virtual void initialise ()=0
 Initialise the model. More...
 
virtual LauComplex amplitude (const LauKinematics *kinematics)
 Calculate the complex amplitude. More...
 
virtual LauResonanceModel getResonanceModel () const =0
 Get the resonance model type. More...
 
const TString & getResonanceName () const
 Get the name of the resonance. More...
 
Int_t getPairInt () const
 Get the integer to identify which DP axis the resonance belongs to. More...
 
Int_t getSpin () const
 Get the spin of the resonance. More...
 
Int_t getCharge () const
 Get the charge of the resonance. More...
 
Double_t getMass () const
 Get the mass of the resonance. More...
 
Double_t getWidth () const
 Get the width of the resonance. More...
 
Bool_t flipHelicity () const
 Get the helicity flip flag. More...
 
void flipHelicity (Bool_t boolean)
 Set the helicity flip flag. More...
 
Bool_t ignoreMomenta () const
 Get the ignore momenta flag. More...
 
void ignoreMomenta (Bool_t boolean)
 Set the ignore p_ and q_ flag. More...
 
void changeResonance (Double_t newMass, Double_t newWidth, Int_t newSpin)
 Allow the mass, width and spin of the resonance to be changed. More...
 
virtual void setResonanceParameter (Double_t value, const TString &name)
 Set the updated parameters from changeResonance. More...
 

Protected Member Functions

TString getNameParent () const
 Get the name of the parent particle. More...
 
TString getNameDaug1 () const
 Get the name of the first daughter of the resonance. More...
 
TString getNameDaug2 () const
 Get the name of the second daughter of the resonance. More...
 
TString getNameBachelor () const
 Get the name of the daughter that does not originate form the resonance. More...
 
Double_t getMassParent () const
 Get the parent particle mass. More...
 
Double_t getMassDaug1 () const
 Get the mass of daughter 1. More...
 
Double_t getMassDaug2 () const
 Get the mass of daughter 2. More...
 
Double_t getMassBachelor () const
 Get the mass of the bachelor daughter. More...
 
Int_t getChargeParent () const
 Get the Charge of the parent particle. More...
 
Int_t getChargeDaug1 () const
 Get the charge of daughter 1. More...
 
Int_t getChargeDaug2 () const
 Get the charge of daughter 2. More...
 
Int_t getChargeBachelor () const
 Get the charge of the bachelor daughter. More...
 
Double_t getQ () const
 Get the current value of the daughter momentum in the resonance rest frame. More...
 
Double_t getP () const
 Get the current value of the bachelor momentum in the resonance rest frame. More...
 
Double_t getPstar () const
 Get the current value of the bachelor momentum in the parent rest frame. More...
 
const LauDaughtersgetDaughters () const
 Access the daughters object. More...
 
virtual LauComplex resAmp (Double_t mass, Double_t spinTerm)=0
 Complex resonant amplitude. More...
 

Private Attributes

const LauDaughtersdaughters_
 Information on the particles. More...
 
TString nameParent_
 Parent name. More...
 
TString nameDaug1_
 Daughter 1 name. More...
 
TString nameDaug2_
 Daughter 2 name. More...
 
TString nameBachelor_
 Bachelor name. More...
 
Int_t chargeParent_
 Parent charge. More...
 
Int_t chargeDaug1_
 Daughter 1 charge. More...
 
Int_t chargeDaug2_
 Daughter 2 charge. More...
 
Int_t chargeBachelor_
 Bachelor charge. More...
 
Double_t massParent_
 Parent mass. More...
 
Double_t massDaug1_
 Daughter 1 mass. More...
 
Double_t massDaug2_
 Daughter 2 mass. More...
 
Double_t massBachelor_
 
TString resName_
 Resonance name. More...
 
Double_t resMass_
 Resonance mass. More...
 
Double_t resWidth_
 Resonance width. More...
 
Int_t resSpin_
 Resonance spin. More...
 
Int_t resCharge_
 Resonance charge. More...
 
Int_t resPairAmpInt_
 DP axis identifier. More...
 
Bool_t flipHelicity_
 Boolean to flip helicity. More...
 
Bool_t ignoreMomenta_
 Boolean to ignore q_ and p_ in spinTerm. More...
 
Double_t q_
 Daughter momentum in resonance rest frame. More...
 
Double_t p_
 Bachelor momentum in resonance rest frame. More...
 
Double_t pstar_
 Bachelor momentum in parent rest frame. More...
 

Detailed Description

Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc.)

Abstract Class for defining the type for all classes used to model resonances in the Dalitz plot, such as Breit-Wigner functions. In addition, some common functionality is implemented, including data such as the mass and width of the desired state.

Definition at line 34 of file LauAbsResonance.hh.

Member Enumeration Documentation

Define the allowed types of barrier factors.

Enumerator
BWBarrier 

Blatt-Weisskopf barrier factor (for use when momentum terms not used in angular term)

BWPrimeBarrier 

Blatt-Weisskopf barrier factor (for use when momentum terms are used in angular term) - the default

ExpBarrier 

expoential barrier factor (mostly used for virtual contributions)

Definition at line 57 of file LauAbsResonance.hh.

Define the allowed resonance types.

Enumerator
BW 

simple Breit-Wigner

RelBW 

relativistic Breit-Wigner

GS 

a modified Breit-Wigner from Gounaris-Sakurai

Flatte 

Flatte or coupled-channel Breit-Wigner

Sigma 

special shape for the sigma or f_0(600)

Kappa 

special shape for the kappa, a low-mass Kpi scalar

Dabba 

special shape for the dabba, a low-mass Dpi scalar

LASS 

the LASS amplitude to describe the Kpi S-wave

LASS_BW 

the resonant part of the LASS amplitude

LASS_NR 

the nonresonant part of the LASS amplitude

KMatrix 

S-wave description using K-matrix and P-vector

FlatNR 

a uniform nonresonant amplitude

NRModel 

a theoretical model nonresonant amplitude

BelleNR 

an empirical exponential nonresonant amplitude

PolNR 

an empirical polynomial nonresonant amplitude

Definition at line 38 of file LauAbsResonance.hh.

Constructor & Destructor Documentation

LauAbsResonance::LauAbsResonance ( const TString &  resName,
Double_t  resMass,
Double_t  resWidth,
Int_t  resSpin,
Int_t  resCharge,
Int_t  resPairAmpInt,
const LauDaughters daughters 
)

Constructor.

Parameters
[in]resNamethe name of the resonance
[in]resMassthe mass of the resonance
[in]resWidththe width of the resonance
[in]resSpinthe spin of the resonance
[in]resChargethe charge of the resonance
[in]resPairAmpIntthe number of the daughter not produced by the resonance
[in]daughtersthe daughter particles

Definition at line 38 of file LauAbsResonance.cc.

LauAbsResonance::~LauAbsResonance ( )
virtual

Destructor.

Definition at line 83 of file LauAbsResonance.cc.

Member Function Documentation

LauComplex LauAbsResonance::amplitude ( const LauKinematics kinematics)
virtual

Calculate the complex amplitude.

Parameters
[in]kinematicsthe kinematic variables of the current event
Returns
the complex amplitude

Reimplemented in LauBelleSymNR, LauNRAmplitude, LauFlatNR, LauKMatrixProdPole, and LauKMatrixProdSVP.

Definition at line 87 of file LauAbsResonance.cc.

void LauAbsResonance::changeResonance ( Double_t  newMass,
Double_t  newWidth,
Int_t  newSpin 
)

Allow the mass, width and spin of the resonance to be changed.

Parameters
[in]newMassnew value of the resonance mass
[in]newWidthnew value of the resonance width
[in]newSpinnew value of the resonance spin

Definition at line 154 of file LauAbsResonance.cc.

Bool_t LauAbsResonance::flipHelicity ( ) const
inline

Get the helicity flip flag.

Returns
the flip helicity flag

Definition at line 136 of file LauAbsResonance.hh.

void LauAbsResonance::flipHelicity ( Bool_t  boolean)
inline

Set the helicity flip flag.

Parameters
[in]booleanthe helicity flip status

Definition at line 142 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getCharge ( ) const
inline

Get the charge of the resonance.

Returns
the resonance charge

Definition at line 118 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getChargeBachelor ( ) const
protected

Get the charge of the bachelor daughter.

Definition at line 291 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeDaug1 ( ) const
protected

Get the charge of daughter 1.

Definition at line 255 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeDaug2 ( ) const
protected

Get the charge of daughter 2.

Definition at line 273 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeParent ( ) const
protected

Get the Charge of the parent particle.

Definition at line 243 of file LauAbsResonance.cc.

const LauDaughters* LauAbsResonance::getDaughters ( ) const
inlineprotected

Access the daughters object.

Definition at line 205 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMass ( ) const
inline

Get the mass of the resonance.

Returns
the resonance mass

Definition at line 124 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMassBachelor ( ) const
protected

Get the mass of the bachelor daughter.

Definition at line 225 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassDaug1 ( ) const
protected

Get the mass of daughter 1.

Definition at line 189 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassDaug2 ( ) const
protected

Get the mass of daughter 2.

Definition at line 207 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassParent ( ) const
protected

Get the parent particle mass.

Definition at line 177 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameBachelor ( ) const
protected

Get the name of the daughter that does not originate form the resonance.

Definition at line 357 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameDaug1 ( ) const
protected

Get the name of the first daughter of the resonance.

Definition at line 321 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameDaug2 ( ) const
protected

Get the name of the second daughter of the resonance.

Definition at line 339 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameParent ( ) const
protected

Get the name of the parent particle.

Definition at line 309 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getP ( ) const
inlineprotected

Get the current value of the bachelor momentum in the resonance rest frame.

Definition at line 200 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getPairInt ( ) const
inline

Get the integer to identify which DP axis the resonance belongs to.

Returns
the DP axis identification number, the ID of the bachelor

Definition at line 106 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getPstar ( ) const
inlineprotected

Get the current value of the bachelor momentum in the parent rest frame.

Definition at line 202 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getQ ( ) const
inlineprotected

Get the current value of the daughter momentum in the resonance rest frame.

Definition at line 198 of file LauAbsResonance.hh.

virtual LauResonanceModel LauAbsResonance::getResonanceModel ( ) const
pure virtual
const TString& LauAbsResonance::getResonanceName ( ) const
inline

Get the name of the resonance.

Returns
the resonance name

Definition at line 100 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getSpin ( ) const
inline

Get the spin of the resonance.

Returns
the resonance spin

Definition at line 112 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getWidth ( ) const
inline

Get the width of the resonance.

Returns
the resonance width

Definition at line 130 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreMomenta ( ) const
inline

Get the ignore momenta flag.

Returns
the ignore momenta flag

Definition at line 148 of file LauAbsResonance.hh.

void LauAbsResonance::ignoreMomenta ( Bool_t  boolean)
inline

Set the ignore p_ and q_ flag.

Parameters
[in]booleanthe ignore momenta status

Definition at line 154 of file LauAbsResonance.hh.

virtual LauComplex LauAbsResonance::resAmp ( Double_t  mass,
Double_t  spinTerm 
)
protectedpure virtual

Complex resonant amplitude.

Parameters
[in]massappropriate invariant mass for the resonance
[in]spinTermZemach spin term

Implemented in LauLASSBWRes, LauLASSNRRes, LauLASSRes, LauFlatteRes, LauBelleSymNR, LauBelleNR, LauPolNR, LauKappaRes, LauSigmaRes, LauDabbaRes, LauGounarisSakuraiRes, LauRelBreitWignerRes, LauNRAmplitude, LauFlatNR, LauBreitWignerRes, LauKMatrixProdPole, and LauKMatrixProdSVP.

void LauAbsResonance::setResonanceParameter ( Double_t  value,
const TString &  name 
)
virtual

Set the updated parameters from changeResonance.

Parameters
[in]valuethe updated parameter value
[in]namethe updated parameter name

Reimplemented in LauLASSBWRes, LauLASSNRRes, LauLASSRes, and LauFlatteRes.

Definition at line 171 of file LauAbsResonance.cc.

Member Data Documentation

Int_t LauAbsResonance::chargeBachelor_
private

Bachelor charge.

Definition at line 234 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeDaug1_
private

Daughter 1 charge.

Definition at line 230 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeDaug2_
private

Daughter 2 charge.

Definition at line 232 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeParent_
private

Parent charge.

Definition at line 228 of file LauAbsResonance.hh.

const LauDaughters* LauAbsResonance::daughters_
private

Information on the particles.

Definition at line 216 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::flipHelicity_
private

Boolean to flip helicity.

Definition at line 259 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreMomenta_
private

Boolean to ignore q_ and p_ in spinTerm.

Definition at line 261 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massBachelor_
private

Definition at line 243 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massDaug1_
private

Daughter 1 mass.

Definition at line 239 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massDaug2_
private

Daughter 2 mass.

Definition at line 241 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massParent_
private

Parent mass.

Definition at line 237 of file LauAbsResonance.hh.

TString LauAbsResonance::nameBachelor_
private

Bachelor name.

Definition at line 225 of file LauAbsResonance.hh.

TString LauAbsResonance::nameDaug1_
private

Daughter 1 name.

Definition at line 221 of file LauAbsResonance.hh.

TString LauAbsResonance::nameDaug2_
private

Daughter 2 name.

Definition at line 223 of file LauAbsResonance.hh.

TString LauAbsResonance::nameParent_
private

Parent name.

Definition at line 219 of file LauAbsResonance.hh.

Double_t LauAbsResonance::p_
private

Bachelor momentum in resonance rest frame.

Definition at line 266 of file LauAbsResonance.hh.

Double_t LauAbsResonance::pstar_
private

Bachelor momentum in parent rest frame.

Definition at line 268 of file LauAbsResonance.hh.

Double_t LauAbsResonance::q_
private

Daughter momentum in resonance rest frame.

Definition at line 264 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resCharge_
private

Resonance charge.

Definition at line 254 of file LauAbsResonance.hh.

Double_t LauAbsResonance::resMass_
private

Resonance mass.

Definition at line 248 of file LauAbsResonance.hh.

TString LauAbsResonance::resName_
private

Resonance name.

Definition at line 246 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resPairAmpInt_
private

DP axis identifier.

Definition at line 256 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resSpin_
private

Resonance spin.

Definition at line 252 of file LauAbsResonance.hh.

Double_t LauAbsResonance::resWidth_
private

Resonance width.

Definition at line 250 of file LauAbsResonance.hh.


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