laura is hosted by Hepforge, IPPP Durham
Laura++  v3r0
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:
LauAbsIncohRes LauBelleNR LauBelleSymNR LauBreitWignerRes LauDabbaRes LauFlatNR LauFlatteRes LauGounarisSakuraiRes LauKappaRes LauKMatrixProdPole LauKMatrixProdSVP LauLASSBWRes LauLASSNRRes LauLASSRes LauModIndPartWave LauNRAmplitude LauPolNR LauRelBreitWignerRes LauSigmaRes

Public Types

enum  LauResonanceModel {
  BW, RelBW, GS, Flatte,
  Sigma, Kappa, Dabba, LASS,
  LASS_BW, LASS_NR, KMatrix, FlatNR,
  NRModel, BelleNR, PowerLawNR, BelleSymNR,
  TaylorNR, PolNR, MIPW, GaussIncoh
}
 Define the allowed resonance types. More...
 

Public Member Functions

 LauAbsResonance (LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters)
 Constructor (for use by standard resonances) More...
 
 LauAbsResonance (const TString &resName, const Int_t resPairAmpInt, const LauDaughters *daughters)
 Constructor (for use by K-matrix components) 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...
 
const TString & getSanitisedName () 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...
 
LauParametergetMassPar ()
 Get the mass parameter of the resonance. More...
 
LauParametergetWidthPar ()
 Get the width parameter of the resonance. More...
 
virtual const std::vector
< LauParameter * > & 
getFloatingParameters ()
 Retrieve the resonance parameters, e.g. so that they can be loaded into a fit. More...
 
Bool_t flipHelicity () const
 Get the helicity flip flag. More...
 
void flipHelicity (const Bool_t boolean)
 Set the helicity flip flag. More...
 
Bool_t ignoreMomenta () const
 Get the ignore momenta flag. More...
 
void ignoreMomenta (const Bool_t boolean)
 Set the ignore p_ and q_ flag. More...
 
void changeResonance (const Double_t newMass, const Double_t newWidth, const Int_t newSpin)
 Allow the mass, width and spin of the resonance to be changed. More...
 
void changeBWBarrierRadii (const Double_t resRadius, const Double_t parRadius)
 Allow the Blatt-Weisskopf radius for the resonance and parent factors to be changed. More...
 
virtual void setResonanceParameter (const TString &name, const Double_t value)
 Set value of the various parameters. More...
 
virtual void floatResonanceParameter (const TString &name)
 Allow the various parameters to float in the fit. More...
 
virtual LauParametergetResonanceParameter (const TString &name)
 Access the given resonance parameter. More...
 
void fixMass (const Bool_t parFixed)
 Fix or release the resonance mass. More...
 
void fixWidth (const Bool_t parFixed)
 Fix or release the resonance width. More...
 
Bool_t fixMass () const
 Get the status of resonance mass (fixed or released) More...
 
Bool_t fixWidth () const
 Get the status of resonance width (fixed or released) More...
 
void setBarrierRadii (LauBlattWeisskopfFactor *resFactor, LauBlattWeisskopfFactor *parFactor)
 Set the form factor model and parameters. More...
 
void fixBarrierRadii (const Bool_t fixResRadius, const Bool_t fixParRadius)
 Fix or release the Blatt-Weisskopf barrier radii. More...
 
Bool_t fixResRadius () const
 Get the status of resonance barrier radius (fixed or released) More...
 
Bool_t fixParRadius () const
 Get the status of parent barrier radius (fixed or released) More...
 
Double_t getResRadius () const
 Get the radius of the resonance barrier factor. More...
 
Double_t getParRadius () const
 Get the radius of the parent barrier factor. 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...
 
LauBlattWeisskopfFactorgetParBWFactor ()
 Get the centrifugal barrier for the parent decay. More...
 
const LauBlattWeisskopfFactorgetParBWFactor () const
 
LauBlattWeisskopfFactorgetResBWFactor ()
 Get the centrifugal barrier for the resonance decay. More...
 
const LauBlattWeisskopfFactorgetResBWFactor () const
 
LauResonanceInfogetResInfo () const
 Access the resonance info object. More...
 
const LauDaughtersgetDaughters () const
 Access the daughters object. More...
 
virtual LauComplex resAmp (Double_t mass, Double_t spinTerm)=0
 Complex resonant amplitude. More...
 
void clearFloatingParameters ()
 Clear list of floating parameters. More...
 
void addFloatingParameter (LauParameter *param)
 Add parameter to the list of floating parameters. More...
 
std::vector< LauParameter * > & getParameters ()
 Access the list of floating parameters. More...
 

Private Member Functions

 LauAbsResonance (const LauAbsResonance &rhs)
 Copy constructor (not implemented) More...
 
LauAbsResonanceoperator= (const LauAbsResonance &rhs)
 Copy assignment operator (not implemented) More...
 

Private Attributes

LauResonanceInforesInfo_
 Information on the resonance. More...
 
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...
 
TString sanitisedName_
 Resonance name with illegal characters removed. More...
 
LauParameterresMass_
 Resonance mass. More...
 
LauParameterresWidth_
 Resonance width. More...
 
std::vector< LauParameter * > resParameters_
 All parameters of the resonance. More...
 
Int_t resSpin_
 Resonance spin. More...
 
Int_t resCharge_
 Resonance charge. More...
 
Int_t resPairAmpInt_
 DP axis identifier. More...
 
LauBlattWeisskopfFactorparBWFactor_
 Blatt Weisskopf barrier for parent decay. More...
 
LauBlattWeisskopfFactorresBWFactor_
 Blatt Weisskopf barrier for resonance decay. 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 37 of file LauAbsResonance.hh.

Member Enumeration Documentation

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

PowerLawNR 

an empirical power law nonresonant amplitude

BelleSymNR 

an empirical exponential nonresonant amplitude for symmetrised DPs

TaylorNR 

an empirical Taylor expansion nonresonant amplitude for symmetrised DPs

PolNR 

an empirical polynomial nonresonant amplitude

MIPW 

a model independent partial wave

GaussIncoh 

an incoherent Gaussian shape

Definition at line 41 of file LauAbsResonance.hh.

Constructor & Destructor Documentation

LauAbsResonance::LauAbsResonance ( LauResonanceInfo resInfo,
const Int_t  resPairAmpInt,
const LauDaughters daughters 
)

Constructor (for use by standard resonances)

Parameters
[in]resInfothe object containing information on the resonance name, mass, width, spin, charge, etc.
[in]resPairAmpIntthe number of the daughter not produced by the resonance
[in]daughtersthe daughter particles

Definition at line 30 of file LauAbsResonance.cc.

LauAbsResonance::LauAbsResonance ( const TString &  resName,
const Int_t  resPairAmpInt,
const LauDaughters daughters 
)

Constructor (for use by K-matrix components)

Parameters
[in]resNamethe name of the component
[in]resPairAmpIntthe number of the daughter not produced by the resonance
[in]daughtersthe daughter particles

Definition at line 83 of file LauAbsResonance.cc.

LauAbsResonance::~LauAbsResonance ( )
virtual

Destructor.

Definition at line 131 of file LauAbsResonance.cc.

LauAbsResonance::LauAbsResonance ( const LauAbsResonance rhs)
private

Copy constructor (not implemented)

Member Function Documentation

void LauAbsResonance::addFloatingParameter ( LauParameter param)
protected

Add parameter to the list of floating parameters.

Parameters
[in]paramthe parameter to be added to the list

Definition at line 278 of file LauAbsResonance.cc.

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, LauKMatrixProdPole, LauKMatrixProdSVP, LauNRAmplitude, and LauFlatNR.

Definition at line 135 of file LauAbsResonance.cc.

void LauAbsResonance::changeBWBarrierRadii ( const Double_t  resRadius,
const Double_t  parRadius 
)

Allow the Blatt-Weisskopf radius for the resonance and parent factors to be changed.

Negative values wil be ignored, so if, for example, you want to only change the parent radius you can provide a negative value for the resonance radius

Parameters
[in]resRadiusnew value of the resonance radius
[in]parRadiusnew value of the parent radius

Definition at line 241 of file LauAbsResonance.cc.

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

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

Negative values wil be ignored, so if, for example, you want to only change the spin you can provide negative values for the mass and width

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 221 of file LauAbsResonance.cc.

void LauAbsResonance::clearFloatingParameters ( )
inlineprotected

Clear list of floating parameters.

Definition at line 330 of file LauAbsResonance.hh.

void LauAbsResonance::fixBarrierRadii ( const Bool_t  fixResRadius,
const Bool_t  fixParRadius 
)

Fix or release the Blatt-Weisskopf barrier radii.

Definition at line 291 of file LauAbsResonance.cc.

void LauAbsResonance::fixMass ( const Bool_t  parFixed)
inline

Fix or release the resonance mass.

Parameters
[in]parFixednew status of mass

Definition at line 230 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::fixMass ( ) const
inline

Get the status of resonance mass (fixed or released)

Returns
the status of resonance mass (fixed or released)

Definition at line 242 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::fixParRadius ( ) const

Get the status of parent barrier radius (fixed or released)

Definition at line 321 of file LauAbsResonance.cc.

Bool_t LauAbsResonance::fixResRadius ( ) const

Get the status of resonance barrier radius (fixed or released)

Definition at line 310 of file LauAbsResonance.cc.

void LauAbsResonance::fixWidth ( const Bool_t  parFixed)
inline

Fix or release the resonance width.

Parameters
[in]parFixednew status of width

Definition at line 236 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::fixWidth ( ) const
inline

Get the status of resonance width (fixed or released)

Returns
the status of resonance width (fixed or released)

Definition at line 248 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::flipHelicity ( ) const
inline

Get the helicity flip flag.

Returns
the flip helicity flag

Definition at line 163 of file LauAbsResonance.hh.

void LauAbsResonance::flipHelicity ( const Bool_t  boolean)
inline

Set the helicity flip flag.

Parameters
[in]booleanthe helicity flip status

Definition at line 169 of file LauAbsResonance.hh.

void LauAbsResonance::floatResonanceParameter ( const TString &  name)
virtual

Allow the various parameters to float in the fit.

Parameters
[in]namethe name of the parameter to be floated

Reimplemented in LauLASSNRRes, LauLASSRes, LauBelleSymNR, LauNRAmplitude, LauBelleNR, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 265 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getCharge ( ) const
inline

Get the charge of the resonance.

Returns
the resonance charge

Definition at line 127 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getChargeBachelor ( ) const
protected

Get the charge of the bachelor daughter.

Definition at line 468 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeDaug1 ( ) const
protected

Get the charge of daughter 1.

Definition at line 432 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeDaug2 ( ) const
protected

Get the charge of daughter 2.

Definition at line 450 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeParent ( ) const
protected

Get the Charge of the parent particle.

Definition at line 420 of file LauAbsResonance.cc.

const LauDaughters* LauAbsResonance::getDaughters ( ) const
inlineprotected

Access the daughters object.

Definition at line 320 of file LauAbsResonance.hh.

virtual const std::vector<LauParameter*>& LauAbsResonance::getFloatingParameters ( )
inlinevirtual

Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.

Returns
floating parameters of the resonance

Reimplemented in LauLASSNRRes, LauLASSRes, LauBelleSymNR, LauNRAmplitude, LauBelleNR, LauModIndPartWave, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, LauSigmaRes, LauGaussIncohRes, LauBreitWignerRes, LauGounarisSakuraiRes, and LauRelBreitWignerRes.

Definition at line 157 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMass ( ) const
inline

Get the mass of the resonance.

Returns
the resonance mass

Definition at line 133 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMassBachelor ( ) const
protected

Get the mass of the bachelor daughter.

Definition at line 402 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassDaug1 ( ) const
protected

Get the mass of daughter 1.

Definition at line 366 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassDaug2 ( ) const
protected

Get the mass of daughter 2.

Definition at line 384 of file LauAbsResonance.cc.

LauParameter* LauAbsResonance::getMassPar ( )
inline

Get the mass parameter of the resonance.

Returns
the resonance mass parameter

Definition at line 145 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMassParent ( ) const
protected

Get the parent particle mass.

Definition at line 354 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 534 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameDaug1 ( ) const
protected

Get the name of the first daughter of the resonance.

Definition at line 498 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameDaug2 ( ) const
protected

Get the name of the second daughter of the resonance.

Definition at line 516 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameParent ( ) const
protected

Get the name of the parent particle.

Definition at line 486 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 305 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 115 of file LauAbsResonance.hh.

std::vector<LauParameter*>& LauAbsResonance::getParameters ( )
inlineprotected

Access the list of floating parameters.

Definition at line 339 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::getParBWFactor ( )
inlineprotected

Get the centrifugal barrier for the parent decay.

Definition at line 310 of file LauAbsResonance.hh.

const LauBlattWeisskopfFactor* LauAbsResonance::getParBWFactor ( ) const
inlineprotected

Definition at line 311 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getParRadius ( ) const

Get the radius of the parent barrier factor.

Definition at line 343 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getPstar ( ) const
inlineprotected

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

Definition at line 307 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 303 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::getResBWFactor ( )
inlineprotected

Get the centrifugal barrier for the resonance decay.

Definition at line 313 of file LauAbsResonance.hh.

const LauBlattWeisskopfFactor* LauAbsResonance::getResBWFactor ( ) const
inlineprotected

Definition at line 314 of file LauAbsResonance.hh.

LauResonanceInfo* LauAbsResonance::getResInfo ( ) const
inlineprotected

Access the resonance info object.

Definition at line 317 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 103 of file LauAbsResonance.hh.

LauParameter * LauAbsResonance::getResonanceParameter ( const TString &  name)
virtual

Access the given resonance parameter.

Parameters
[in]namethe name of the parameter
Returns
the corresponding parameter

Reimplemented in LauLASSNRRes, LauLASSRes, LauBelleSymNR, LauNRAmplitude, LauBelleNR, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 271 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getResRadius ( ) const

Get the radius of the resonance barrier factor.

Definition at line 332 of file LauAbsResonance.cc.

const TString& LauAbsResonance::getSanitisedName ( ) const
inline

Get the name of the resonance.

Returns
the resonance name

Definition at line 109 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getSpin ( ) const
inline

Get the spin of the resonance.

Returns
the resonance spin

Definition at line 121 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getWidth ( ) const
inline

Get the width of the resonance.

Returns
the resonance width

Definition at line 139 of file LauAbsResonance.hh.

LauParameter* LauAbsResonance::getWidthPar ( )
inline

Get the width parameter of the resonance.

Returns
the resonance width parameter

Definition at line 151 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreMomenta ( ) const
inline

Get the ignore momenta flag.

Returns
the ignore momenta flag

Definition at line 175 of file LauAbsResonance.hh.

void LauAbsResonance::ignoreMomenta ( const Bool_t  boolean)
inline

Set the ignore p_ and q_ flag.

Parameters
[in]booleanthe ignore momenta status

Definition at line 181 of file LauAbsResonance.hh.

LauAbsResonance& LauAbsResonance::operator= ( const LauAbsResonance rhs)
private

Copy assignment operator (not implemented)

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 LauKappaRes, LauSigmaRes, LauDabbaRes, LauLASSRes, LauLASSNRRes, LauLASSBWRes, LauFlatteRes, LauBelleNR, LauBelleSymNR, LauNRAmplitude, LauModIndPartWave, LauPolNR, LauBreitWignerRes, LauGounarisSakuraiRes, LauKMatrixProdPole, LauKMatrixProdSVP, LauRelBreitWignerRes, LauFlatNR, and LauAbsIncohRes.

void LauAbsResonance::setBarrierRadii ( LauBlattWeisskopfFactor resFactor,
LauBlattWeisskopfFactor parFactor 
)
inline

Set the form factor model and parameters.

Parameters
[in]resFactorthe barrier factor for the resonance decay
[in]parFactorthe barrier factor for the parent decay

Definition at line 255 of file LauAbsResonance.hh.

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

Set value of the various parameters.

Parameters
[in]namethe name of the parameter to be changed
[in]valuethe new parameter value

Reimplemented in LauLASSNRRes, LauLASSRes, LauBelleSymNR, LauNRAmplitude, LauBelleNR, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 259 of file LauAbsResonance.cc.

Member Data Documentation

Int_t LauAbsResonance::chargeBachelor_
private

Bachelor charge.

Definition at line 371 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeDaug1_
private

Daughter 1 charge.

Definition at line 367 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeDaug2_
private

Daughter 2 charge.

Definition at line 369 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeParent_
private

Parent charge.

Definition at line 365 of file LauAbsResonance.hh.

const LauDaughters* LauAbsResonance::daughters_
private

Information on the particles.

Definition at line 353 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::flipHelicity_
private

Boolean to flip helicity.

Definition at line 408 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreMomenta_
private

Boolean to ignore q_ and p_ in spinTerm.

Definition at line 410 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massBachelor_
private

Definition at line 380 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massDaug1_
private

Daughter 1 mass.

Definition at line 376 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massDaug2_
private

Daughter 2 mass.

Definition at line 378 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massParent_
private

Parent mass.

Definition at line 374 of file LauAbsResonance.hh.

TString LauAbsResonance::nameBachelor_
private

Bachelor name.

Definition at line 362 of file LauAbsResonance.hh.

TString LauAbsResonance::nameDaug1_
private

Daughter 1 name.

Definition at line 358 of file LauAbsResonance.hh.

TString LauAbsResonance::nameDaug2_
private

Daughter 2 name.

Definition at line 360 of file LauAbsResonance.hh.

TString LauAbsResonance::nameParent_
private

Parent name.

Definition at line 356 of file LauAbsResonance.hh.

Double_t LauAbsResonance::p_
private

Bachelor momentum in resonance rest frame.

Definition at line 415 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::parBWFactor_
private

Blatt Weisskopf barrier for parent decay.

Definition at line 403 of file LauAbsResonance.hh.

Double_t LauAbsResonance::pstar_
private

Bachelor momentum in parent rest frame.

Definition at line 417 of file LauAbsResonance.hh.

Double_t LauAbsResonance::q_
private

Daughter momentum in resonance rest frame.

Definition at line 413 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::resBWFactor_
private

Blatt Weisskopf barrier for resonance decay.

Definition at line 405 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resCharge_
private

Resonance charge.

Definition at line 399 of file LauAbsResonance.hh.

LauResonanceInfo* LauAbsResonance::resInfo_
private

Information on the resonance.

Definition at line 350 of file LauAbsResonance.hh.

LauParameter* LauAbsResonance::resMass_
private

Resonance mass.

Definition at line 389 of file LauAbsResonance.hh.

TString LauAbsResonance::resName_
private

Resonance name.

Definition at line 383 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resPairAmpInt_
private

DP axis identifier.

Definition at line 401 of file LauAbsResonance.hh.

std::vector<LauParameter*> LauAbsResonance::resParameters_
private

All parameters of the resonance.

Definition at line 394 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resSpin_
private

Resonance spin.

Definition at line 397 of file LauAbsResonance.hh.

LauParameter* LauAbsResonance::resWidth_
private

Resonance width.

Definition at line 391 of file LauAbsResonance.hh.

TString LauAbsResonance::sanitisedName_
private

Resonance name with illegal characters removed.

Definition at line 386 of file LauAbsResonance.hh.


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