laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
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 LauPolarFormFactorNR LauPolarFormFactorSymNR LauPoleRes LauPolNR LauRelBreitWignerRes LauRescattering2Res LauRescatteringRes 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,
  Pole, PolarFFNR, PolarFFSymNR, PolarFFSymNRNoInter,
  Rescattering, Rescattering2, RescatteringNoInter, 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, const Int_t resSpin)
 Constructor (for use by K-matrix components) More...
 
virtual ~LauAbsResonance ()
 Destructor.
 
virtual void initialise ()=0
 Initialise the model.
 
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.
 
Bool_t fixResRadius () const
 Get the status of resonance barrier radius (fixed or released)
 
Bool_t fixParRadius () const
 Get the status of parent barrier radius (fixed or released)
 
Double_t getResRadius () const
 Get the radius of the resonance barrier factor.
 
Double_t getParRadius () const
 Get the radius of the parent barrier factor.
 

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.
 
TString getNameDaug1 () const
 Get the name of the first daughter of the resonance.
 
TString getNameDaug2 () const
 Get the name of the second daughter of the resonance.
 
TString getNameBachelor () const
 Get the name of the daughter that does not originate form the resonance.
 
Double_t getMassParent () const
 Get the parent particle mass.
 
Double_t getMassDaug1 () const
 Get the mass of daughter 1.
 
Double_t getMassDaug2 () const
 Get the mass of daughter 2.
 
Double_t getMassBachelor () const
 Get the mass of the bachelor daughter.
 
Int_t getChargeParent () const
 Get the Charge of the parent particle.
 
Int_t getChargeDaug1 () const
 Get the charge of daughter 1.
 
Int_t getChargeDaug2 () const
 Get the charge of daughter 2.
 
Int_t getChargeBachelor () const
 Get the charge of the bachelor daughter.
 
Double_t getQ () const
 Get the current value of the daughter momentum in the resonance rest frame.
 
Double_t getP () const
 Get the current value of the bachelor momentum in the resonance rest frame.
 
Double_t getPstar () const
 Get the current value of the bachelor momentum in the parent rest frame.
 
Double_t getCovFactor () const
 Get the current value of the full spin-dependent covariant factor.
 
LauBlattWeisskopfFactorgetParBWFactor ()
 Get the centrifugal barrier for the parent decay.
 
const LauBlattWeisskopfFactorgetParBWFactor () const
 Get the centrifugal barrier for the parent decay.
 
LauBlattWeisskopfFactorgetResBWFactor ()
 Get the centrifugal barrier for the resonance decay.
 
const LauBlattWeisskopfFactorgetResBWFactor () const
 Get the centrifugal barrier for the resonance decay.
 
LauResonanceInfogetResInfo () const
 Access the resonance info object.
 
const LauDaughtersgetDaughters () const
 Access the daughters object.
 
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.
 
void addFloatingParameter (LauParameter *param)
 Add parameter to the list of floating parameters. More...
 
std::vector< LauParameter * > & getParameters ()
 Access the list of floating parameters.
 

Private Member Functions

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

Private Attributes

LauResonanceInforesInfo_ { 0 }
 Information on the resonance.
 
const LauDaughtersdaughters_
 Information on the particles.
 
TString nameParent_ { "" }
 Parent name.
 
TString nameDaug1_ { "" }
 Daughter 1 name.
 
TString nameDaug2_ { "" }
 Daughter 2 name.
 
TString nameBachelor_ { "" }
 Bachelor name.
 
Int_t chargeParent_ { 0 }
 Parent charge.
 
Int_t chargeDaug1_ { 0 }
 Daughter 1 charge.
 
Int_t chargeDaug2_ { 0 }
 Daughter 2 charge.
 
Int_t chargeBachelor_ { 0 }
 Bachelor charge.
 
Double_t massParent_ { 0.0 }
 Parent mass.
 
Double_t massDaug1_ { 0.0 }
 Daughter 1 mass.
 
Double_t massDaug2_ { 0.0 }
 Daughter 2 mass.
 
Double_t massBachelor_ { 0.0 }
 Bachelor mass.
 
TString resName_
 Resonance name.
 
TString sanitisedName_
 Resonance name with illegal characters removed.
 
LauParameterresMass_ { 0 }
 Resonance mass.
 
LauParameterresWidth_ { 0 }
 Resonance width.
 
std::vector< LauParameter * > resParameters_
 All parameters of the resonance.
 
Int_t resSpin_
 Resonance spin.
 
Int_t resCharge_ { 0 }
 Resonance charge.
 
Int_t resPairAmpInt_
 DP axis identifier.
 
LauBlattWeisskopfFactorparBWFactor_ { 0 }
 Blatt Weisskopf barrier for parent decay.
 
LauBlattWeisskopfFactorresBWFactor_ { 0 }
 Blatt Weisskopf barrier for resonance decay.
 
LauSpinType spinType_ { Zemach_P }
 Spin formalism.
 
Bool_t flipHelicity_ { kFALSE }
 Boolean to flip helicity.
 
Bool_t ignoreMomenta_ { kFALSE }
 Boolean to ignore the momentum factors in both the spin factor and the mass-dependent width.
 
Bool_t ignoreSpin_ { kFALSE }
 Boolean to set the spinTerm to unity always.
 
Bool_t ignoreBarrierScaling_ { kFALSE }
 Boolean to ignore barrier factor scaling in the amplitude numerator, they are still used for the mass-dependent width.
 
Double_t mass_ { 0.0 }
 Invariant mass.
 
Double_t cosHel_ { 0.0 }
 Helicity angle cosine.
 
Double_t q_ { 0.0 }
 Daughter momentum in resonance rest frame.
 
Double_t p_ { 0.0 }
 Bachelor momentum in resonance rest frame.
 
Double_t pstar_ { 0.0 }
 Bachelor momentum in parent rest frame.
 
Double_t erm_ { 1.0 }
 Covariant factor. More...
 
Double_t covFactor_ { 1.0 }
 Covariant factor (full spin-dependent expression)
 

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 51 of file LauAbsResonance.hh.

Member Enumeration Documentation

◆ LauResonanceModel

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 

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

Pole 

scalar Pole lineshape

PolarFFNR 

Polar Form Factor nonresonant amplitude

PolarFFSymNR 

Polar Form Factor nonresonant amplitude for symmetrised DPs

PolarFFSymNRNoInter 

Polar Form Factor nonresonant amplitude for symmetrised DPs without interference

Rescattering 

KK-PiPi inelastic scattering amplitude

Rescattering2 

KK-PiPi inelastic scattering amplitude

RescatteringNoInter 

KK-PiPi inelastic scattering 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 55 of file LauAbsResonance.hh.

◆ LauSpinType

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 93 of file LauAbsResonance.hh.

Constructor & Destructor Documentation

◆ LauAbsResonance() [1/2]

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

◆ LauAbsResonance() [2/2]

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

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
[in]resSpinthe spin of the final channel into which the K-matrix scatters

Definition at line 133 of file LauAbsResonance.cc.

Member Function Documentation

◆ addFloatingParameter()

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

◆ amplitude()

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, LauPolarFormFactorSymNR, LauRescatteringRes, LauBelleSymNR, LauNRAmplitude, and LauFlatNR.

Definition at line 172 of file LauAbsResonance.cc.

◆ calcCovFactor()

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

◆ calcCovSpinFactor()

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

◆ calcLegendrePoly() [1/2]

Double_t LauAbsResonance::calcLegendrePoly ( ) const
protected

Calculate the Legendre polynomial for the spin factor.

Uses the current-event value of cosHel_

Definition at line 342 of file LauAbsResonance.cc.

◆ calcLegendrePoly() [2/2]

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

◆ calcZemachSpinFactor()

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

◆ changeBWBarrierRadii()

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

◆ changeResonance()

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

◆ fixMass() [1/2]

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 354 of file LauAbsResonance.hh.

◆ fixMass() [2/2]

void LauAbsResonance::fixMass ( const Bool_t  parFixed)
inline

Fix or release the resonance mass.

Parameters
[in]parFixednew status of mass

Definition at line 332 of file LauAbsResonance.hh.

◆ fixWidth() [1/2]

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 360 of file LauAbsResonance.hh.

◆ fixWidth() [2/2]

void LauAbsResonance::fixWidth ( const Bool_t  parFixed)
inline

Fix or release the resonance width.

Parameters
[in]parFixednew status of width

Definition at line 343 of file LauAbsResonance.hh.

◆ flipHelicity() [1/2]

Bool_t LauAbsResonance::flipHelicity ( ) const
inline

Get the helicity flip flag.

Returns
the flip helicity flag

Definition at line 229 of file LauAbsResonance.hh.

◆ flipHelicity() [2/2]

void LauAbsResonance::flipHelicity ( const Bool_t  boolean)
inline

Set the helicity flip flag.

Parameters
[in]booleanthe helicity flip status

Definition at line 235 of file LauAbsResonance.hh.

◆ floatResonanceParameter()

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, LauBelleNR, LauBelleSymNR, LauEFKLLMRes, LauPolarFormFactorSymNR, LauRescatteringRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauPolarFormFactorNR, LauRescattering2Res, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 426 of file LauAbsResonance.cc.

◆ getCharge()

Int_t LauAbsResonance::getCharge ( ) const
inline

Get the charge of the resonance.

Returns
the resonance charge

Definition at line 182 of file LauAbsResonance.hh.

◆ getFloatingParameters()

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

◆ getMass()

Double_t LauAbsResonance::getMass ( ) const
inline

Get the mass of the resonance.

Returns
the resonance mass

Definition at line 188 of file LauAbsResonance.hh.

◆ getMassPar()

LauParameter* LauAbsResonance::getMassPar ( )
inline

Get the mass parameter of the resonance.

Returns
the resonance mass parameter

Definition at line 200 of file LauAbsResonance.hh.

◆ getPairInt()

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 170 of file LauAbsResonance.hh.

◆ getResonanceModel()

◆ getResonanceName()

const TString& LauAbsResonance::getResonanceName ( ) const
inline

Get the name of the resonance.

Returns
the resonance name

Definition at line 158 of file LauAbsResonance.hh.

◆ getResonanceParameter()

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, LauBelleNR, LauBelleSymNR, LauEFKLLMRes, LauPolarFormFactorSymNR, LauRescatteringRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauPolarFormFactorNR, LauRescattering2Res, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 433 of file LauAbsResonance.cc.

◆ getSanitisedName()

const TString& LauAbsResonance::getSanitisedName ( ) const
inline

Get the name of the resonance.

Returns
the resonance name

Definition at line 164 of file LauAbsResonance.hh.

◆ getSpin()

Int_t LauAbsResonance::getSpin ( ) const
inline

Get the spin of the resonance.

Returns
the resonance spin

Definition at line 176 of file LauAbsResonance.hh.

◆ getSpinType()

LauSpinType LauAbsResonance::getSpinType ( ) const
inline

Get the spin type.

Returns
the spin formalism

Definition at line 152 of file LauAbsResonance.hh.

◆ getWidth()

Double_t LauAbsResonance::getWidth ( ) const
inline

Get the width of the resonance.

Returns
the resonance width

Definition at line 194 of file LauAbsResonance.hh.

◆ getWidthPar()

LauParameter* LauAbsResonance::getWidthPar ( )
inline

Get the width parameter of the resonance.

Returns
the resonance width parameter

Definition at line 206 of file LauAbsResonance.hh.

◆ ignoreBarrierScaling() [1/2]

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 275 of file LauAbsResonance.hh.

◆ ignoreBarrierScaling() [2/2]

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 283 of file LauAbsResonance.hh.

◆ ignoreMomenta() [1/2]

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 243 of file LauAbsResonance.hh.

◆ ignoreMomenta() [2/2]

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 251 of file LauAbsResonance.hh.

◆ ignoreSpin() [1/2]

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 259 of file LauAbsResonance.hh.

◆ ignoreSpin() [2/2]

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 267 of file LauAbsResonance.hh.

◆ isIncoherentModel()

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

◆ preSymmetrised()

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 223 of file LauAbsResonance.hh.

◆ resAmp()

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

◆ setBarrierRadii()

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 373 of file LauAbsResonance.hh.

◆ setResonanceParameter()

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, LauBelleNR, LauBelleSymNR, LauEFKLLMRes, LauPolarFormFactorSymNR, LauRescatteringRes, LauLASSNRRes, LauLASSRes, LauNRAmplitude, LauPolarFormFactorNR, LauRescattering2Res, LauKappaRes, LauDabbaRes, LauFlatteRes, LauLASSBWRes, and LauSigmaRes.

Definition at line 419 of file LauAbsResonance.cc.

◆ setSpinType()

void LauAbsResonance::setSpinType ( const LauSpinType  spinType)
inline

Set the spin formalism to be used.

Parameters
[in]spinTypethe spin formalism

Definition at line 366 of file LauAbsResonance.hh.

Member Data Documentation

◆ erm_

Double_t LauAbsResonance::erm_ { 1.0 }
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 600 of file LauAbsResonance.hh.


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