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

Class for defining the EFKLLM K-pi S-wave model. More...

#include <LauEFKLLMRes.hh>

Inheritance diagram for LauEFKLLMRes:
LauAbsResonance

Public Member Functions

 LauEFKLLMRes (LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters)
 Constructor. More...
 
virtual ~LauEFKLLMRes ()
 Destructor. More...
 
virtual void initialise ()
 Initialise the model. More...
 
virtual
LauAbsResonance::LauResonanceModel 
getResonanceModel () const
 Get the resonance model type. More...
 
virtual void setResonanceParameter (const TString &name, const Double_t value)
 Set value of a resonance parameter. 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...
 
virtual const std::vector
< LauParameter * > & 
getFloatingParameters ()
 Retrieve the resonance parameters, e.g. so that they can be loaded into a fit. More...
 
- Public Member Functions inherited from LauAbsResonance
 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 LauComplex amplitude (const LauKinematics *kinematics)
 Calculate the complex amplitude. 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 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...
 
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 void setupFormFactor (const TString &inputFile)
 Read the form factor information from text file. More...
 
- Static Public Member Functions inherited from LauAbsResonance
static bool isIncoherentModel (LauResonanceModel model)
 Is the resonance model incoherent? More...
 

Protected Member Functions

void setMassFactor (const Double_t massFactor)
 Set the power of the mass dependence. More...
 
Double_t getMassFactor () const
 Get the power of the mass dependence. More...
 
Bool_t fixMassFactor () const
 See if the mass factor parameter is fixed or floating. More...
 
virtual LauComplex resAmp (Double_t mass, Double_t spinTerm)
 Complex resonant amplitude. More...
 
- Protected Member Functions inherited from LauAbsResonance
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...
 
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 Attributes

LauParametermassFactor_
 The power of the mass dependence. More...
 

Static Private Attributes

static Lau1DCubicSplinemagSpline_ = 0
 Spline describing the magnitude variation of the form-factor. More...
 
static Lau1DCubicSplinephaseSpline_ = 0
 Spline describing the phase variation of the form-factor. More...
 

Additional Inherited Members

- Public Types inherited from LauAbsResonance
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...
 

Detailed Description

Class for defining the EFKLLM K-pi S-wave model.

Class for defining the EFKLLM form-factor model for the K-pi S-wave. The model consists of a tabulated form-factor, which is interpolated using cubic splines (one for magnitude values and one for phase values), multiplied by a mass-dependence (e.g. constant, 1/m^2, etc.). The massFactor resonance parameter is the power of the mass dependence - defaults to zero, i.e. constant. For more details see B. El-Bennich et al. Phys. Rev. D 79, 094005 (2009), arXiv:0902.3645 [hep-ph]. (The acronym EFKLLM is constructed from the surnames of the authors of the above paper.)

Definition at line 50 of file LauEFKLLMRes.hh.

Constructor & Destructor Documentation

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

Constructor.

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 45 of file LauEFKLLMRes.cc.

LauEFKLLMRes::~LauEFKLLMRes ( )
virtual

Destructor.

Definition at line 63 of file LauEFKLLMRes.cc.

Member Function Documentation

Bool_t LauEFKLLMRes::fixMassFactor ( ) const
inlineprotected

See if the mass factor parameter is fixed or floating.

Returns
kTRUE if the mass factor parameter is fixed, kFALSE otherwise

Definition at line 125 of file LauEFKLLMRes.hh.

void LauEFKLLMRes::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 from LauAbsResonance.

Definition at line 81 of file LauEFKLLMRes.cc.

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

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

Returns
floating parameters of the resonance

Reimplemented from LauAbsResonance.

Definition at line 105 of file LauEFKLLMRes.cc.

Double_t LauEFKLLMRes::getMassFactor ( ) const
inlineprotected

Get the power of the mass dependence.

Returns
the power of the mass dependence

Definition at line 119 of file LauEFKLLMRes.hh.

virtual LauAbsResonance::LauResonanceModel LauEFKLLMRes::getResonanceModel ( ) const
inlinevirtual

Get the resonance model type.

Returns
the resonance model type

Implements LauAbsResonance.

Definition at line 80 of file LauEFKLLMRes.hh.

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

Access the given resonance parameter.

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

Reimplemented from LauAbsResonance.

Definition at line 95 of file LauEFKLLMRes.cc.

void LauEFKLLMRes::initialise ( )
virtual

Initialise the model.

Implements LauAbsResonance.

Definition at line 67 of file LauEFKLLMRes.cc.

LauComplex LauEFKLLMRes::resAmp ( Double_t  mass,
Double_t  spinTerm 
)
protectedvirtual

Complex resonant amplitude.

Parameters
[in]massappropriate invariant mass for the resonance
[in]spinTermspin term

Implements LauAbsResonance.

Definition at line 116 of file LauEFKLLMRes.cc.

void LauEFKLLMRes::setMassFactor ( const Double_t  massFactor)
protected

Set the power of the mass dependence.

Parameters
[in]massFactorthe new power of the mass dependence

Definition at line 169 of file LauEFKLLMRes.cc.

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

Set value of a resonance parameter.

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

Reimplemented from LauAbsResonance.

Definition at line 71 of file LauEFKLLMRes.cc.

void LauEFKLLMRes::setupFormFactor ( const TString &  inputFile)
static

Read the form factor information from text file.

Creates the splines from the tabulated form factor data. These are shared between all instances of this class.

Parameters
[in]inputFilethe name of the file to be read

Definition at line 135 of file LauEFKLLMRes.cc.

Member Data Documentation

Lau1DCubicSpline * LauEFKLLMRes::magSpline_ = 0
staticprivate

Spline describing the magnitude variation of the form-factor.

Definition at line 136 of file LauEFKLLMRes.hh.

LauParameter* LauEFKLLMRes::massFactor_
private

The power of the mass dependence.

Definition at line 141 of file LauEFKLLMRes.hh.

Lau1DCubicSpline * LauEFKLLMRes::phaseSpline_ = 0
staticprivate

Spline describing the phase variation of the form-factor.

Definition at line 138 of file LauEFKLLMRes.hh.


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