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

Public Types

enum  LauResonanceModel {
  BW, RelBW, GS, Flatte,
  Sigma, Kappa, Dabba, LASS,
  LASS_BW, LASS_NR, EFKLLM, KMatrix,
  FlatNR, NRModel, BelleNR, PowerLawNR,
  BelleSymNR, BelleSymNRNoInter, TaylorNR, PolNR,
  MIPW_MagPhase, MIPW_RealImag, GaussIncoh, RhoOmegaMix_GS,
  RhoOmegaMix_RBW, RhoOmegaMix_GS_1, RhoOmegaMix_RBW_1
}
 Define the allowed resonance types. More...
 
enum  LauSpinType { Zemach_P, Zemach_Pstar, Covariant, Legendre }
 Define the allowed spin formalisms. 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...
 
LauSpinType getSpinType () const
 Get the spin 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...
 
virtual Bool_t preSymmetrised () const
 Is the amplitude pre-symmetrised? 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 momenta flag. More...
 
Bool_t ignoreSpin () const
 Get the ignore spin flag. More...
 
void ignoreSpin (const Bool_t boolean)
 Set the ignore spin flag. More...
 
Bool_t ignoreBarrierScaling () const
 Get the ignore barrier factor scaling flag. More...
 
void ignoreBarrierScaling (const Bool_t boolean)
 Set the ignore barrier factor scaling 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 setSpinType (const LauSpinType spinType)
 Set the spin formalism to be used. 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...
 

Static Public Member Functions

static bool isIncoherentModel (LauResonanceModel model)
 Is the resonance model incoherent? 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...
 
Double_t getCovFactor () const
 Get the current value of the full spin-dependent covariant factor. 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...
 
Double_t calcZemachSpinFactor (const Double_t pProd) const
 Calculate the amplitude spin term using the Zemach tensor formalism. More...
 
Double_t calcCovSpinFactor (const Double_t pProd)
 Calculate the amplitude spin term using the covariant tensor formalism. More...
 
void calcCovFactor (const Double_t erm)
 Calculate the spin-dependent covariant factor. More...
 
Double_t calcLegendrePoly () const
 Calculate the Legendre polynomial for the spin factor. More...
 
Double_t calcLegendrePoly (const Double_t cosHel)
 Calculate the Legendre polynomial for the spin factor (specifying the cosHel value) 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...
 
LauSpinType spinType_
 Spin formalism. More...
 
Bool_t flipHelicity_
 Boolean to flip helicity. More...
 
Bool_t ignoreMomenta_
 Boolean to ignore the momentum factors in both the spin factor and the mass-dependent width. More...
 
Bool_t ignoreSpin_
 Boolean to set the spinTerm to unity always. More...
 
Bool_t ignoreBarrierScaling_
 Boolean to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width. More...
 
Double_t mass_
 Invariant mass. More...
 
Double_t cosHel_
 Helicity angle cosine. 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...
 
Double_t erm_
 Covariant factor. More...
 
Double_t covFactor_
 Covariant factor (full spin-dependent expression) 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

EFKLLM 

a form-factor-based description of the Kpi S-wave

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

BelleSymNRNoInter 

an empirical exponential nonresonant amplitude for symmetrised DPs without interference

TaylorNR 

an empirical Taylor expansion nonresonant amplitude for symmetrised DPs

PolNR 

an empirical polynomial nonresonant amplitude

MIPW_MagPhase 

a model independent partial wave - magnitude and phase representation

MIPW_RealImag 

a model independent partial wave - real and imaginary part representation

GaussIncoh 

an incoherent Gaussian shape

RhoOmegaMix_GS 

mass mixing model using GS for res 1 and RBW for res 2

RhoOmegaMix_RBW 

mass mixing model using two RBWs

RhoOmegaMix_GS_1 

mass mixing model using GS for res 1 and RBW for res 2, with denominator factor = 1

RhoOmegaMix_RBW_1 

mass mixing model using two RBWs, with denominator factor = 1

Definition at line 41 of file LauAbsResonance.hh.

Define the allowed spin formalisms.

Enumerator
Zemach_P 

Zemach tensor formalism, bachelor momentum in resonance rest frame

Zemach_Pstar 

Zemach tensor formalism, bachelor momentum in parent rest frame

Covariant 

Covariant tensor formalism

Legendre 

Legendre polynomials only

Definition at line 72 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 64 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 124 of file LauAbsResonance.cc.

LauAbsResonance::~LauAbsResonance ( )
virtual

Destructor.

Definition at line 179 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 429 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 LauRhoOmegaMix, LauKMatrixProdPole, LauKMatrixProdSVP, LauBelleSymNR, LauNRAmplitude, and LauFlatNR.

Definition at line 183 of file LauAbsResonance.cc.

void LauAbsResonance::calcCovFactor ( const Double_t  erm)
protected

Calculate the spin-dependent covariant factor.

Parameters
[in]ermE_ij in the parent rest-frame divided by m_ij (equivalent to sqrt(1 + p^2/mParent^2))

Definition at line 273 of file LauAbsResonance.cc.

Double_t LauAbsResonance::calcCovSpinFactor ( const Double_t  pProd)
protected

Calculate the amplitude spin term using the covariant tensor formalism.

Parameters
[in]pProdthe momentum factor (q * pstar)

Definition at line 292 of file LauAbsResonance.cc.

Double_t LauAbsResonance::calcLegendrePoly ( ) const
protected

Calculate the Legendre polynomial for the spin factor.

Uses the current-event value of cosHel_

Definition at line 346 of file LauAbsResonance.cc.

Double_t LauAbsResonance::calcLegendrePoly ( const Double_t  cosHel)
protected

Calculate the Legendre polynomial for the spin factor (specifying the cosHel value)

Parameters
[in]cosHelthe cosine of the helicity angle

Definition at line 340 of file LauAbsResonance.cc.

Double_t LauAbsResonance::calcZemachSpinFactor ( const Double_t  pProd) const
protected

Calculate the amplitude spin term using the Zemach tensor formalism.

Parameters
[in]pProdthe momentum factor (either q * p or q * pstar)

Definition at line 314 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 392 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 372 of file LauAbsResonance.cc.

void LauAbsResonance::clearFloatingParameters ( )
inlineprotected

Clear list of floating parameters.

Definition at line 436 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 442 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 298 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 310 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::fixParRadius ( ) const

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

Definition at line 472 of file LauAbsResonance.cc.

Bool_t LauAbsResonance::fixResRadius ( ) const

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

Definition at line 461 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 304 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 316 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::flipHelicity ( ) const
inline

Get the helicity flip flag.

Returns
the flip helicity flag

Definition at line 199 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 205 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 LauRhoOmegaMix, LauBelleSymNR, LauBelleNR, LauEFKLLMRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 416 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getCharge ( ) const
inline

Get the charge of the resonance.

Returns
the resonance charge

Definition at line 155 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getChargeBachelor ( ) const
protected

Get the charge of the bachelor daughter.

Definition at line 619 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeDaug1 ( ) const
protected

Get the charge of daughter 1.

Definition at line 583 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeDaug2 ( ) const
protected

Get the charge of daughter 2.

Definition at line 601 of file LauAbsResonance.cc.

Int_t LauAbsResonance::getChargeParent ( ) const
protected

Get the Charge of the parent particle.

Definition at line 571 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getCovFactor ( ) const
inlineprotected

Get the current value of the full spin-dependent covariant factor.

Definition at line 383 of file LauAbsResonance.hh.

const LauDaughters* LauAbsResonance::getDaughters ( ) const
inlineprotected

Access the daughters object.

Definition at line 396 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 LauRhoOmegaMix, LauAbsModIndPartWave, LauBelleSymNR, LauBelleNR, LauEFKLLMRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, LauSigmaRes, LauGaussIncohRes, LauBreitWignerRes, LauGounarisSakuraiRes, and LauRelBreitWignerRes.

Definition at line 185 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMass ( ) const
inline

Get the mass of the resonance.

Returns
the resonance mass

Definition at line 161 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMassBachelor ( ) const
protected

Get the mass of the bachelor daughter.

Definition at line 553 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassDaug1 ( ) const
protected

Get the mass of daughter 1.

Definition at line 517 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getMassDaug2 ( ) const
protected

Get the mass of daughter 2.

Definition at line 535 of file LauAbsResonance.cc.

LauParameter* LauAbsResonance::getMassPar ( )
inline

Get the mass parameter of the resonance.

Returns
the resonance mass parameter

Definition at line 173 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getMassParent ( ) const
protected

Get the parent particle mass.

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

TString LauAbsResonance::getNameDaug1 ( ) const
protected

Get the name of the first daughter of the resonance.

Definition at line 649 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameDaug2 ( ) const
protected

Get the name of the second daughter of the resonance.

Definition at line 667 of file LauAbsResonance.cc.

TString LauAbsResonance::getNameParent ( ) const
protected

Get the name of the parent particle.

Definition at line 637 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 379 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 143 of file LauAbsResonance.hh.

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

Access the list of floating parameters.

Definition at line 445 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::getParBWFactor ( )
inlineprotected

Get the centrifugal barrier for the parent decay.

Definition at line 386 of file LauAbsResonance.hh.

const LauBlattWeisskopfFactor* LauAbsResonance::getParBWFactor ( ) const
inlineprotected

Definition at line 387 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getParRadius ( ) const

Get the radius of the parent barrier factor.

Definition at line 494 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 381 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 377 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::getResBWFactor ( )
inlineprotected

Get the centrifugal barrier for the resonance decay.

Definition at line 389 of file LauAbsResonance.hh.

const LauBlattWeisskopfFactor* LauAbsResonance::getResBWFactor ( ) const
inlineprotected

Definition at line 390 of file LauAbsResonance.hh.

LauResonanceInfo* LauAbsResonance::getResInfo ( ) const
inlineprotected

Access the resonance info object.

Definition at line 393 of file LauAbsResonance.hh.

const TString& LauAbsResonance::getResonanceName ( ) const
inline

Get the name of the resonance.

Returns
the resonance name

Definition at line 131 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 LauRhoOmegaMix, LauBelleSymNR, LauBelleNR, LauEFKLLMRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 422 of file LauAbsResonance.cc.

Double_t LauAbsResonance::getResRadius ( ) const

Get the radius of the resonance barrier factor.

Definition at line 483 of file LauAbsResonance.cc.

const TString& LauAbsResonance::getSanitisedName ( ) const
inline

Get the name of the resonance.

Returns
the resonance name

Definition at line 137 of file LauAbsResonance.hh.

Int_t LauAbsResonance::getSpin ( ) const
inline

Get the spin of the resonance.

Returns
the resonance spin

Definition at line 149 of file LauAbsResonance.hh.

LauSpinType LauAbsResonance::getSpinType ( ) const
inline

Get the spin type.

Returns
the spin formalism

Definition at line 125 of file LauAbsResonance.hh.

Double_t LauAbsResonance::getWidth ( ) const
inline

Get the width of the resonance.

Returns
the resonance width

Definition at line 167 of file LauAbsResonance.hh.

LauParameter* LauAbsResonance::getWidthPar ( )
inline

Get the width parameter of the resonance.

Returns
the resonance width parameter

Definition at line 179 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreBarrierScaling ( ) const
inline

Get the ignore barrier factor scaling flag.

Whether to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width

Returns
the ignore barrier amplitude scaling flag

Definition at line 242 of file LauAbsResonance.hh.

void LauAbsResonance::ignoreBarrierScaling ( const Bool_t  boolean)
inline

Set the ignore barrier factor scaling flag.

Whether to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width

Parameters
[in]booleanthe ignore barrier factor scaling status

Definition at line 249 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreMomenta ( ) const
inline

Get the ignore momenta flag.

Whether to ignore the momentum factors in both the spin factor and the mass-dependent width

Returns
the ignore momenta flag

Definition at line 213 of file LauAbsResonance.hh.

void LauAbsResonance::ignoreMomenta ( const Bool_t  boolean)
inline

Set the ignore momenta flag.

Whether to ignore the momentum factors in both the spin factor and the mass-dependent width

Parameters
[in]booleanthe ignore momenta status

Definition at line 221 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreSpin ( ) const
inline

Get the ignore spin flag.

Whether to set the spinTerm to unity always

Returns
the ignore spin flag

Definition at line 228 of file LauAbsResonance.hh.

void LauAbsResonance::ignoreSpin ( const Bool_t  boolean)
inline

Set the ignore spin flag.

Whether to set the spinTerm to unity always

Parameters
[in]booleanthe ignore spin status

Definition at line 235 of file LauAbsResonance.hh.

bool LauAbsResonance::isIncoherentModel ( LauResonanceModel  model)
static

Is the resonance model incoherent?

Parameters
[in]modelthe resonance model
Returns
true if the model is incoherent

Definition at line 28 of file LauAbsResonance.cc.

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

Copy assignment operator (not implemented)

virtual Bool_t LauAbsResonance::preSymmetrised ( ) const
inlinevirtual

Is the amplitude pre-symmetrised?

The default value is kFALSE, so pre-symmetrised lineshapes should override this.

Returns
whether the amplitude is already symmetrised

Reimplemented in LauBelleSymNR.

Definition at line 193 of file LauAbsResonance.hh.

virtual LauComplex LauAbsResonance::resAmp ( Double_t  mass,
Double_t  spinTerm 
)
protectedpure virtual
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 329 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 LauRhoOmegaMix, LauBelleSymNR, LauBelleNR, LauEFKLLMRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 410 of file LauAbsResonance.cc.

void LauAbsResonance::setSpinType ( const LauSpinType  spinType)
inline

Set the spin formalism to be used.

Parameters
[in]spinTypethe spin formalism

Definition at line 322 of file LauAbsResonance.hh.

Member Data Documentation

Int_t LauAbsResonance::chargeBachelor_
private

Bachelor charge.

Definition at line 476 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeDaug1_
private

Daughter 1 charge.

Definition at line 472 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeDaug2_
private

Daughter 2 charge.

Definition at line 474 of file LauAbsResonance.hh.

Int_t LauAbsResonance::chargeParent_
private

Parent charge.

Definition at line 470 of file LauAbsResonance.hh.

Double_t LauAbsResonance::cosHel_
private

Helicity angle cosine.

Definition at line 529 of file LauAbsResonance.hh.

Double_t LauAbsResonance::covFactor_
private

Covariant factor (full spin-dependent expression)

Definition at line 551 of file LauAbsResonance.hh.

const LauDaughters* LauAbsResonance::daughters_
private

Information on the particles.

Definition at line 458 of file LauAbsResonance.hh.

Double_t LauAbsResonance::erm_
private

Covariant factor.

sqrt(1 + z*z), where z = p / mParent

Can also be expressed as E_ij in the parent rest-frame divided by m_ij - indeed this is how LauKinematics calculates it.

See Also
LauKinematics::getcov12
LauKinematics::getcov13
LauKinematics::getcov23

Definition at line 548 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::flipHelicity_
private

Boolean to flip helicity.

Definition at line 516 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreBarrierScaling_
private

Boolean to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width.

Definition at line 522 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreMomenta_
private

Boolean to ignore the momentum factors in both the spin factor and the mass-dependent width.

Definition at line 518 of file LauAbsResonance.hh.

Bool_t LauAbsResonance::ignoreSpin_
private

Boolean to set the spinTerm to unity always.

Definition at line 520 of file LauAbsResonance.hh.

Double_t LauAbsResonance::mass_
private

Invariant mass.

Definition at line 527 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massBachelor_
private

Definition at line 485 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massDaug1_
private

Daughter 1 mass.

Definition at line 481 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massDaug2_
private

Daughter 2 mass.

Definition at line 483 of file LauAbsResonance.hh.

Double_t LauAbsResonance::massParent_
private

Parent mass.

Definition at line 479 of file LauAbsResonance.hh.

TString LauAbsResonance::nameBachelor_
private

Bachelor name.

Definition at line 467 of file LauAbsResonance.hh.

TString LauAbsResonance::nameDaug1_
private

Daughter 1 name.

Definition at line 463 of file LauAbsResonance.hh.

TString LauAbsResonance::nameDaug2_
private

Daughter 2 name.

Definition at line 465 of file LauAbsResonance.hh.

TString LauAbsResonance::nameParent_
private

Parent name.

Definition at line 461 of file LauAbsResonance.hh.

Double_t LauAbsResonance::p_
private

Bachelor momentum in resonance rest frame.

Definition at line 534 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::parBWFactor_
private

Blatt Weisskopf barrier for parent decay.

Definition at line 508 of file LauAbsResonance.hh.

Double_t LauAbsResonance::pstar_
private

Bachelor momentum in parent rest frame.

Definition at line 536 of file LauAbsResonance.hh.

Double_t LauAbsResonance::q_
private

Daughter momentum in resonance rest frame.

Definition at line 532 of file LauAbsResonance.hh.

LauBlattWeisskopfFactor* LauAbsResonance::resBWFactor_
private

Blatt Weisskopf barrier for resonance decay.

Definition at line 510 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resCharge_
private

Resonance charge.

Definition at line 504 of file LauAbsResonance.hh.

LauResonanceInfo* LauAbsResonance::resInfo_
private

Information on the resonance.

Definition at line 455 of file LauAbsResonance.hh.

LauParameter* LauAbsResonance::resMass_
private

Resonance mass.

Definition at line 494 of file LauAbsResonance.hh.

TString LauAbsResonance::resName_
private

Resonance name.

Definition at line 488 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resPairAmpInt_
private

DP axis identifier.

Definition at line 506 of file LauAbsResonance.hh.

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

All parameters of the resonance.

Definition at line 499 of file LauAbsResonance.hh.

Int_t LauAbsResonance::resSpin_
private

Resonance spin.

Definition at line 502 of file LauAbsResonance.hh.

LauParameter* LauAbsResonance::resWidth_
private

Resonance width.

Definition at line 496 of file LauAbsResonance.hh.

TString LauAbsResonance::sanitisedName_
private

Resonance name with illegal characters removed.

Definition at line 491 of file LauAbsResonance.hh.

LauSpinType LauAbsResonance::spinType_
private

Spin formalism.

Definition at line 513 of file LauAbsResonance.hh.


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