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

Class for defining signal dynamics using the isobar model. More...

#include <LauIsobarDynamics.hh>

Inheritance diagram for LauIsobarDynamics:
LauAbsDPDynamics

Public Member Functions

 LauIsobarDynamics (LauDaughters *daughters, LauEffModel *effModel, LauEffModel *scfFractionModel=0)
 Constructor. More...
 
 LauIsobarDynamics (LauDaughters *daughters, LauEffModel *effModel, LauTagCatScfFractionModelMap scfFractionModel)
 Constructor. More...
 
virtual ~LauIsobarDynamics ()
 Destructor. More...
 
virtual void initialise (const std::vector< LauComplex > &coeffs)
 Initialise the Dalitz plot dynamics. More...
 
void setIntFileName (const TString &fileName)
 Set the name of the file to which to save the results of the integrals. More...
 
virtual void setIntegralBinWidths (Double_t m13BinWidth, Double_t m23BinWidth)
 Set the widths of the bins to use when integrating across the Dalitz plot. More...
 
virtual void addResonance (const TString &resName, Int_t resPairAmpInt, const TString &resType, Double_t newMass=-1.0, Double_t newWidth=-1.0, Int_t newSpin=-1)
 Add a resonance to the Dalitz plot. More...
 
virtual void defineKMatrixPropagator (const TString &propName, const TString &paramFileName, Int_t resPairAmpInt, Int_t nChannels, Int_t nPoles, Int_t rowIndex=1)
 Define a new K-matrix Propagator. More...
 
virtual void addKMatrixProdPole (const TString &poleName, const TString &propName, Int_t poleIndex)
 Add a K-matrix production pole term to the model. More...
 
virtual void addKMatrixProdSVP (const TString &SVPName, const TString &propName, Int_t channelIndex)
 Add a K-matrix slowly-varying part (SVP) term to the model. More...
 
virtual void changeResonance (const TString &resName, Double_t newMass=-1.0, Double_t newWidth=-1.0, Int_t newSpin=-1)
 Change the properties of a resonance particle within this model. More...
 
virtual Bool_t generate ()
 Generate a toy MC signal event. More...
 
virtual ToyMCStatus checkToyMC (Bool_t printErrorMessages=kTRUE, Bool_t printInfoMessages=kFALSE)
 Check the status of the toy MC generation. More...
 
Int_t maxGenIterations () const
 Retrieve the maximum number of iterations allowed when generating an event. More...
 
virtual void calcLikelihoodInfo (UInt_t iEvt)
 Calculate the likelihood (and all associated information) for the given event number. More...
 
virtual void calcLikelihoodInfo (Double_t m13Sq, Double_t m23Sq)
 Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates. More...
 
virtual void calcLikelihoodInfo (Double_t m13Sq, Double_t m23Sq, Int_t tagCat)
 Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates and the tagging category. More...
 
virtual const LauComplexgetEvtDPAmp () const
 Retrieve the total amplitude for the current event. More...
 
virtual Double_t getEvtm13Sq () const
 Retrieve the invariant mass squared of the first and third daughters in the current event. More...
 
virtual Double_t getEvtm23Sq () const
 Retrieve the invariant mass squared of the second and third daughters in the current event. More...
 
virtual Double_t getEvtmPrime () const
 Retrieve the square Dalitz plot coordinate, m', for the current event. More...
 
virtual Double_t getEvtthPrime () const
 Retrieve the square Dalitz plot coordinate, theta', for the current event. More...
 
virtual Double_t getEvtEff () const
 Retrieve the efficiency for the current event. More...
 
virtual Double_t getEvtScfFraction () const
 Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the current event. More...
 
virtual Double_t getEvtJacobian () const
 Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event. More...
 
virtual Double_t getEvtLikelihood () const
 Retrieve the likelihood for the current event. More...
 
virtual void calcExtraInfo (Bool_t init=kFALSE)
 Calculate the fit fractions, mean efficiency and total DP rate. More...
 
virtual Bool_t gotReweightedEvent ()
 Calculates whether an event with the current kinematics should be accepted in order to produce a distribution of events that matches the model e.g. when reweighting embedded data. More...
 
virtual Double_t getEventWeight ()
 Calculate the acceptance rate, for events with the current kinematics, when generating events according to the model. More...
 
void setASqMaxValue (Double_t value)
 Set the maximum value of A squared to be used in the accept/reject. More...
 
Double_t getASqMaxSetValue () const
 Retrieve the maximum value of A squared to be used in the accept/reject. More...
 
Double_t getASqMaxVarValue () const
 Retrieve the maximum of A squared that has been found while generating. More...
 
virtual LauComplex getDynamicAmp (Int_t resID) const
 Retrieve the normalised dynamic part of the amplitude of the given amplitude component at the current point in the Dalitz plot. More...
 
const std::vector< std::vector
< LauComplex > > & 
getFiFjSum () const
 Retrieve the event-by-event running totals of amplitude cross terms for all pairs of amplitude components. More...
 
const std::vector< std::vector
< LauComplex > > & 
getFiFjEffSum () const
 Retrieve the event-by-event running totals of efficiency corrected amplitude cross terms for all pairs of amplitude components. More...
 
const std::vector< Double_t > & getFNorm () const
 Retrieve the normalisation factors for the dynamic parts of the amplitudes for all of the amplitude components. More...
 
virtual void fillDataTree (const LauFitDataTree &fitDataTree)
 Fill the internal data structure that caches the resonance dynamics. More...
 
virtual void setDataEventNo (UInt_t iEvt)
 Load the data for a given event. More...
 
void setBelleNRAlpha (Double_t alpha)
 Set the alpha parameter for new Belle non-resonant components. More...
 
Double_t getBelleNRAlpha () const
 Retreive the alpha parameter for new Belle non-resonant components. More...
 
void setLASSParameters (Double_t a, Double_t r, Double_t R=1.0, Double_t phiR=0.0, Double_t B=1.0, Double_t phiB=0.0, Double_t cutOff=1.8)
 Set the parameters for new LASS components. More...
 
void setLASSScatteringLength (Double_t a)
 Set the scattering length for new LASS components. More...
 
void setLASSEffectiveRange (Double_t r)
 Set the effective range for new LASS components. More...
 
void setLASSResonanceMag (Double_t R)
 Set the resonance magnitude for new LASS components. More...
 
void setLASSResonancePhase (Double_t phiR)
 Set the resonance phase for new LASS components. More...
 
void setLASSBackgroundMag (Double_t B)
 Set the background magnitude for new LASS components. More...
 
void setLASSBackgroundPhase (Double_t phiB)
 Set the background phase for new LASS components. More...
 
void setLASSCutOff (Double_t cutOff)
 Set the cutoff value for new LASS components. More...
 
Double_t getLASSScatteringLength () const
 Retrieve the scattering length for new LASS components. More...
 
Double_t getLASSEffectiveRange () const
 Retrieve the effective range for new LASS components. More...
 
Double_t getLASSResonanceMag () const
 Retrieve the resonance magnitude for new LASS components. More...
 
Double_t getLASSResonancePhase () const
 Retrieve the resonance phase for new LASS components. More...
 
Double_t getLASSBackgroundMag () const
 Retrieve the background magnitude for new LASS components. More...
 
Double_t getLASSBackgroundPhase () const
 Retrieve the background phase for new LASS components. More...
 
Double_t getLASSCutOff () const
 Retrieve the cutoff value for new LASS components. More...
 
void setFlatteParameters (Double_t g1, Double_t g2)
 Set the parameters for new Flatte components. More...
 
void setFlatteg1 (Double_t g1)
 Set the g1 parameter for new Flatte components. More...
 
void setFlatteg2 (Double_t g2)
 Set the g2 parameter for new Flatte components. More...
 
Double_t getFlatteg1 () const
 Retrieve the g1 parameter for new Flatte components. More...
 
Double_t getFlatteg2 () const
 Retrieve the g2 parameter for new Flatte components. More...
 
void setBarrierRadii (Double_t resRadius, Double_t parRadius=4.0, LauAbsResonance::BarrierType type=LauAbsResonance::BWPrimeBarrier)
 Set the parameters for the barrier factors for new resonances. More...
 
void setResBarrierRadius (Double_t radius)
 Set the radius of the barrier factor due to the resonance to use for new amplitude components. More...
 
void setParBarrierRadius (Double_t radius)
 Set the radius of the barrier factor due to the parent to use for new amplitude components. More...
 
void setBarrierType (LauAbsResonance::BarrierType type)
 Set the type of barrier factor to use for new amplitude components. More...
 
Double_t getResBarrierRadius () const
 Retrieve the radius of the barrier factor due to the resonance to use for new amplitude components. More...
 
Double_t getParBarrierRadius () const
 Retrieve the radius of the barrier factor due to the parent to use for new amplitude components. More...
 
void flipHelicityForCPEigenstates (Bool_t boolean)
 Set the helicity flip flag for new amplitude components. More...
 
- Public Member Functions inherited from LauAbsDPDynamics
 LauAbsDPDynamics (LauDaughters *daughters, LauEffModel *effModel, LauEffModel *scfFractionModel=0)
 Constructor. More...
 
 LauAbsDPDynamics (LauDaughters *daughters, LauEffModel *effModel, const LauTagCatScfFractionModelMap &scfFractionModel)
 Constructor. More...
 
virtual ~LauAbsDPDynamics ()
 Destructor. More...
 
virtual Double_t retrieveEfficiency ()
 Obtain the efficiency of the current event from the model. More...
 
virtual Double_t retrieveScfFraction (Int_t tagCat)
 Obtain the self cross feed fraction of the current event from the model. More...
 
virtual Bool_t hasResonance (const TString &resName) const
 Check whether this model includes a named resonance. More...
 
virtual TString getConjResName (const TString &resName) const
 Retrieve the name of the charge conjugate of a named resonance. More...
 
virtual void updateCoeffs (const std::vector< LauComplex > &coeffs)
 Update the complex coefficients for the resonances. More...
 
LauParameter getMeanEff () const
 Retrieve the mean efficiency across the Dalitz plot. More...
 
LauParameter getDPRate () const
 Retrieve the overall Dalitz plot rate. More...
 
const LauParArraygetFitFractions () const
 Retrieve the fit fractions for the amplitude components. More...
 
UInt_t getnDefinedResonances () const
 Retrieve the number of defined resonances in the resonance maker. More...
 
UInt_t getnAmp () const
 Retrieve the number of amplitude components. More...
 
Double_t getDPNorm () const
 Retrieve the normalisation factor for the log-likelihood function. More...
 
UInt_t nData () const
 Retrieve the number of cached events. More...
 
const std::vector
< LauCacheData * > & 
getCacheData () const
 Retrieve the cached data. More...
 
LauDaughtersgetDaughters ()
 Retrieve the daughters. More...
 
LauResonanceMakergetResonanceMaker ()
 Retrieve the resonance maker object. More...
 
LauKinematicsgetKinematics ()
 Retrieve the Dalitz plot kinematics. More...
 
LauEffModelgetEffModel ()
 Retrieve the model for the efficiency across the Dalitz plot. More...
 
LauEffModelgetScfFractionModel ()
 Retrieve the model for the fraction of events that are poorly reconstructed (the self cross feed fraction) in each Dalitz plot bin for the first (or only) tagging category. More...
 
std::map< Int_t, LauEffModel * > getScfFractionModels ()
 Retrieve the model for the fraction of events that are poorly reconstructed (the self cross feed fraction) in each Dalitz plot bin for all tagging categories. More...
 
Bool_t usingScfModel ()
 Check whether a self cross feed fraction model is being used. More...
 
std::vector< LauParametergetExtraParameters ()
 Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions) More...
 

Protected Member Functions

virtual void initSummary ()
 Print a summary of the model to be used. More...
 
virtual void initialiseVectors ()
 Initialise the internal storage for this model. More...
 
virtual void calcDPNormalisation ()
 Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot. More...
 
virtual void calcDPPartialIntegral (Double_t minm13, Double_t maxm13, Double_t minm23, Double_t maxm23, Double_t m13BinWidth, Double_t m23BinWidth)
 Calculate the Dalitz plot normalisation integrals over a given range. More...
 
virtual void writeIntegralsFile ()
 Write the results of the integrals (and related information) to a file. More...
 
virtual void setFFTerm (UInt_t index, Double_t realPart, Double_t imagPart)
 Set the dynamic part of the amplitude for a given amplitude component at the current point in the Dalitz plot. More...
 
virtual void dynamics (Bool_t cacheResData=kTRUE, Double_t weight=1.0, Bool_t useEff=kTRUE)
 Calculate the total Dalitz plot amplitude at the current point in the Dalitz plot. More...
 
void setASqMaxVarValue (Double_t value)
 Set the maximum of A squared that has been found. More...
 
virtual Double_t calcSigDPNorm ()
 Calculate the normalisation factor for the log-likelihood function. More...
 
virtual LauComplex resAmp (Int_t index)
 Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz plot. More...
 
virtual LauAbsResonancefindResonance (const TString &name)
 Retrieve the named resonance. More...
 
virtual const LauAbsResonancefindResonance (const TString &name) const
 Retrieve the named resonance. More...
 
virtual void removeCharge (TString &string) const
 Remove the charge from the given particle name. More...
 
Bool_t gotKMatrixMatch (UInt_t resAmpInt, const TString &propName) const
 Check whether a resonance is a K-matrix component of a given propagator. More...
 

Private Types

typedef std::map< TString,
LauKMatrixPropagator * > 
KMPropMap
 The type used for containing the K-matrix propagators. More...
 
typedef std::map< TString,
TString > 
KMStringMap
 The type used for mapping K-matrix components to their propagators. More...
 

Private Attributes

std::vector< LauAbsResonance * > sigResonances_
 The resonances in the model. More...
 
KMPropMap kMatrixPropagators_
 The K-matrix propagators. More...
 
KMStringMap kMatrixPropSet_
 The names of the M-matrix components in the model mapped to their propagators. More...
 
std::vector< TString > resTypAmp_
 The resonance types of all of the amplitude components. More...
 
std::vector< Int_t > resIntAmp_
 The index within the resonance maker for each amplitude component. More...
 
std::vector< Int_t > resPairAmp_
 The index of the daughter not produced by the resonance for each amplitude component. More...
 
std::vector< Int_t > typDaug_
 The PDG codes of the daughters. More...
 
Bool_t symmetricalDP_
 Whether the Dalitz plot is symmetrical. More...
 
Bool_t integralsDone_
 Whether the integrals have been performed. More...
 
TString intFileName_
 The name of the file to save integrals to. More...
 
Double_t m13BinWidth_
 The bin width to use when integrating over m13. More...
 
Double_t m23BinWidth_
 The bin width to use when integrating over m23. More...
 
Double_t m13Sq_
 The invariant mass squared of the first and third daughters. More...
 
Double_t m23Sq_
 The invariant mass squared of the second and third daughters. More...
 
Double_t mPrime_
 The square Dalitz plot coordinate, m'. More...
 
Double_t thPrime_
 The square Dalitz plot coordinate theta'. More...
 
Double_t eff_
 The efficiency at the current point in the Dalitz plot. More...
 
Double_t scfFraction_
 The fraction of events that are poorly reconstructed (the self cross feed fraction) at the current point in the Dalitz plot. More...
 
Double_t jacobian_
 The Jacobian, for the transformation into square DP coordinates at the current point in the Dalitz plot. More...
 
Double_t ASq_
 The value of A squared for the current event. More...
 
Double_t evtLike_
 The normalised likelihood for the current event. More...
 
LauComplex totAmp_
 The total amplitude for the current event. More...
 
std::vector< std::vector
< LauComplex > > 
fifjEffSum_
 The event-by-event running total of efficiency corrected amplitude cross terms for each pair of amplitude components. More...
 
std::vector< std::vector
< LauComplex > > 
fifjSum_
 The event-by-event running total of the amplitude cross terms for each pair of amplitude components. More...
 
std::vector< LauComplexff_
 The dynamic part of the amplitude for each amplitude component at the current point in the Dalitz plot. More...
 
std::vector< Double_t > fSqSum_
 The event-by-event running total of the dynamical amplitude squared for each amplitude component. More...
 
std::vector< Double_t > fNorm_
 The normalisation factors for the dynamic parts of the amplitude for each amplitude component. More...
 
Int_t iterationsMax_
 The maximum allowed number of attempts when generating an event. More...
 
Int_t nSigGenLoop_
 The number of unsucessful attempts to generate an event so far. More...
 
Double_t aSqMaxSet_
 The maximum allowed value of A squared. More...
 
Double_t aSqMaxVar_
 The maximum value of A squared that has been seen so far while generating. More...
 
Double_t BelleNRAlpha_
 The alpha parameter for new Belle non-resonant components. More...
 
Double_t LASSScatteringLength_
 The scattering length for new LASS components. More...
 
Double_t LASSEffectiveRange_
 The effective range for new LASS components. More...
 
Double_t LASSResonanceMag_
 The resonance magnitude for new LASS components. More...
 
Double_t LASSResonancePhase_
 The resonance phase for new LASS components. More...
 
Double_t LASSBackgroundMag_
 The background magnitude for new LASS components. More...
 
Double_t LASSBackgroundPhase_
 The background phase for new LASS components. More...
 
Double_t LASSCutOff_
 The cutoff value for new LASS components. More...
 
Bool_t changeLASSScatteringLength_
 Whether the default value of the LASS scattering length has been changed. More...
 
Bool_t changeLASSEffectiveRange_
 Whether the default value of the LASS effective range has been changed. More...
 
Bool_t changeLASSResonanceMag_
 Whether the default value of the LASS resonance magnitude has been changed. More...
 
Bool_t changeLASSResonancePhase_
 Whether the default value of the LASS resonance phase has been changed. More...
 
Bool_t changeLASSBackgroundMag_
 Whether the default value of the LASS background magnitude has been changed. More...
 
Bool_t changeLASSBackgroundPhase_
 Whether the default value of the LASS background phase has been changed. More...
 
Bool_t changeLASSCutOff_
 Whether the default value of the LASS cutoff has been changed. More...
 
Double_t FlatteParameterg1_
 The constant parameter g1 for new Flatte components. More...
 
Double_t FlatteParameterg2_
 The constant parameter g2 for new Flatte components. More...
 
Bool_t changeFlatteParameterg1_
 Whether the default value of the Flatte parameter g1 has been changed. More...
 
Bool_t changeFlatteParameterg2_
 Whether the default value of the Flatte parameter g2 has been changed. More...
 
Double_t resBarrierRadius_
 The radius of the resonance barrier factor for new amplitude components. More...
 
Double_t parBarrierRadius_
 The radius of the parent barrier factor for new amplitude components. More...
 
LauAbsResonance::BarrierType barrierType_
 The type of the barrier factor for new amplitude components. More...
 
Bool_t flipHelicity_
 The helicity flip flag for new amplitude components. More...
 

Additional Inherited Members

- Public Types inherited from LauAbsDPDynamics
enum  ToyMCStatus { GenOK, MaxIterError, ASqMaxError }
 The possible statuses for toy MC generation. More...
 
typedef std::map< Int_t,
LauEffModel * > 
LauTagCatScfFractionModelMap
 The type used for containing multiple self cross feed fraction models for different categories (e.g. tagging categories) More...
 
- Protected Attributes inherited from LauAbsDPDynamics
LauDaughtersdaughters_
 The daughters of the decay. More...
 
LauResonanceMakerresonanceMaker_
 Object to create resonances. More...
 
LauKinematicskinematics_
 The kinematics of the decay. More...
 
LauEffModeleffModel_
 The efficiency model across the Dalitz plot. More...
 
LauTagCatScfFractionModelMap scfFractionModel_
 The self cross feed fraction models across the Dalitz plot. More...
 
UInt_t nAmp_
 The number of amplitude components. More...
 
UInt_t nResDefMax_
 The number of resonances defined in the resonance maker. More...
 
std::vector< LauComplexAmp_
 The complex coefficients for the amplitude components. More...
 
Double_t DPNorm_
 The normalisation factor for the log-likelihood function. More...
 
LauParArray fitFrac_
 The fit fractions for the amplitude components. More...
 
LauParameter DPRate_
 The overall Dalitz plot rate. More...
 
LauParameter meanDPEff_
 The mean efficiency across the Dalitz plot. More...
 
std::vector< LauCacheData * > data_
 The cached data for all events. More...
 
LauCacheDatacurrentEvent_
 The cached data for the current event. More...
 
std::vector< LauParameterextraParameters_
 any extra parameters/quantities (e.g. K-matrix total fit fractions) More...
 

Detailed Description

Class for defining signal dynamics using the isobar model.

Definition at line 35 of file LauIsobarDynamics.hh.

Member Typedef Documentation

typedef std::map<TString, LauKMatrixPropagator*> LauIsobarDynamics::KMPropMap
private

The type used for containing the K-matrix propagators.

Definition at line 586 of file LauIsobarDynamics.hh.

typedef std::map<TString, TString> LauIsobarDynamics::KMStringMap
private

The type used for mapping K-matrix components to their propagators.

Definition at line 589 of file LauIsobarDynamics.hh.

Constructor & Destructor Documentation

LauIsobarDynamics::LauIsobarDynamics ( LauDaughters daughters,
LauEffModel effModel,
LauEffModel scfFractionModel = 0 
)

Constructor.

Parameters
[in]daughtersthe daughters of the decay
[in]effModelthe model to describe the efficiency across the Dalitz plot
[in]scfFractionModelthe model to describe the fraction of poorly constructed events (the self cross feed fraction) across the Dalitz plot

Definition at line 49 of file LauIsobarDynamics.cc.

LauIsobarDynamics::LauIsobarDynamics ( LauDaughters daughters,
LauEffModel effModel,
LauTagCatScfFractionModelMap  scfFractionModel 
)

Constructor.

Parameters
[in]daughtersthe daughters of the decay
[in]effModelthe model to describe the efficiency across the Dalitz plot
[in]scfFractionModelthe models to describe the fraction of poorly constructed events (the self cross feed fraction) across the Dalitz plot for various tagging categories

Definition at line 107 of file LauIsobarDynamics.cc.

LauIsobarDynamics::~LauIsobarDynamics ( )
virtual

Destructor.

Definition at line 163 of file LauIsobarDynamics.cc.

Member Function Documentation

void LauIsobarDynamics::addKMatrixProdPole ( const TString &  poleName,
const TString &  propName,
Int_t  poleIndex 
)
virtual

Add a K-matrix production pole term to the model.

Parameters
[in]poleNamethe name of the pole
[in]propNamethe name of the propagator to use
[in]poleIndexthe index of the pole within the propagator

Definition at line 557 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::addKMatrixProdSVP ( const TString &  SVPName,
const TString &  propName,
Int_t  channelIndex 
)
virtual

Add a K-matrix slowly-varying part (SVP) term to the model.

Parameters
[in]SVPNamethe name of the term
[in]propNamethe name of the propagator to use
[in]channelIndexthe index of the channel within the propagator

Definition at line 605 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::addResonance ( const TString &  resName,
Int_t  resPairAmpInt,
const TString &  resType,
Double_t  newMass = -1.0,
Double_t  newWidth = -1.0,
Int_t  newSpin = -1 
)
virtual

Add a resonance to the Dalitz plot.

Parameters
[in]resNamethe name of the resonant particle
[in]resPairAmpIntthe index of the daughter not produced by the resonance
[in]resTypethe type of the resonance. Allowed types are: flatte, relbw, dabba, kappa, sigma, lass-bw, lass-nr, lass, gs, nrmodel, bellesymnr and bellenr
[in]newMassset a custom mass for the resonance
[in]newWidthset a custom width for the resonance
[in]newSpinset a custom spin for the resonance

Implements LauAbsDPDynamics.

Definition at line 403 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::calcDPNormalisation ( )
protectedvirtual

Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.

Definition at line 721 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::calcDPPartialIntegral ( Double_t  minm13,
Double_t  maxm13,
Double_t  minm23,
Double_t  maxm23,
Double_t  m13BinWidth,
Double_t  m23BinWidth 
)
protectedvirtual

Calculate the Dalitz plot normalisation integrals over a given range.

Parameters
[in]minm13the minimum value of m13 in the integration range
[in]maxm13the maximum value of m13 in the integration range
[in]minm23the minimum value of m23 in the integration range
[in]maxm23the maximum value of m23 in the integration range
[in]m13BinWidththe bin width in m13
[in]m23BinWidththe bin width in m23

Definition at line 1165 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::calcExtraInfo ( Bool_t  init = kFALSE)
virtual

Calculate the fit fractions, mean efficiency and total DP rate.

Parameters
[in]initwhether the calculated values should be stored as the initial/generated values or the fitted values

Implements LauAbsDPDynamics.

Definition at line 1422 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::calcLikelihoodInfo ( UInt_t  iEvt)
virtual

Calculate the likelihood (and all associated information) for the given event number.

Parameters
[in]iEvtthe event number

Implements LauAbsDPDynamics.

Definition at line 1718 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::calcLikelihoodInfo ( Double_t  m13Sq,
Double_t  m23Sq 
)
virtual

Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates.

Parameters
[in]m13Sqthe invariant mass squared of the first and third daughters
[in]m23Sqthe invariant mass squared of the second and third daughters

Implements LauAbsDPDynamics.

Definition at line 1743 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::calcLikelihoodInfo ( Double_t  m13Sq,
Double_t  m23Sq,
Int_t  tagCat 
)
virtual

Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates and the tagging category.

Also obtain the self cross feed fraction to cache with the rest of the Dalitz plot quantities.

Parameters
[in]m13Sqthe invariant mass squared of the first and third daughters
[in]m23Sqthe invariant mass squared of the second and third daughters
[in]tagCatthe tagging category

Definition at line 1748 of file LauIsobarDynamics.cc.

Double_t LauIsobarDynamics::calcSigDPNorm ( )
protectedvirtual

Calculate the normalisation factor for the log-likelihood function.

Returns
the normalisation factor

Implements LauAbsDPDynamics.

Definition at line 1581 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::changeResonance ( const TString &  resName,
Double_t  newMass = -1.0,
Double_t  newWidth = -1.0,
Int_t  newSpin = -1 
)
virtual

Change the properties of a resonance particle within this model.

Note that parameters set to -1 are ignored.

Parameters
[in]resNamethe name of the resonance to modify
[in]newMassthe new mass for this resonance
[in]newWidththe new width for this resonance
[in]newSpinthe new spin for this resonance

Definition at line 705 of file LauIsobarDynamics.cc.

LauAbsDPDynamics::ToyMCStatus LauIsobarDynamics::checkToyMC ( Bool_t  printErrorMessages = kTRUE,
Bool_t  printInfoMessages = kFALSE 
)
virtual

Check the status of the toy MC generation.

Parameters
[in]printErrorMessageswhether error messages should be printed
[in]printInfoMessageswhether info messages should be printed
Returns
the status of the toy MC generation

Implements LauAbsDPDynamics.

Definition at line 1674 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::defineKMatrixPropagator ( const TString &  propName,
const TString &  paramFileName,
Int_t  resPairAmpInt,
Int_t  nChannels,
Int_t  nPoles,
Int_t  rowIndex = 1 
)
virtual

Define a new K-matrix Propagator.

Parameters
[in]propNamethe name of the propagator
[in]paramFileNamethe file that defines the propagator
[in]resPairAmpIntthe index of the bachelor
[in]nChannelsthe number of channels
[in]nPolesthe number of poles
[in]rowIndexthe index of the row to be used when summing over all amplitude channels: S-wave corresponds to rowIndex = 1.

Definition at line 534 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::dynamics ( Bool_t  cacheResData = kTRUE,
Double_t  weight = 1.0,
Bool_t  useEff = kTRUE 
)
protectedvirtual

Calculate the total Dalitz plot amplitude at the current point in the Dalitz plot.

Parameters
[in]cacheResDatawhether the amplitudes have already been cached
[in]weightthe weight to apply (used when calculating the integrals)
[in]useEffwhether to apply efficiency corrections

Definition at line 1281 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::fillDataTree ( const LauFitDataTree fitDataTree)
virtual

Fill the internal data structure that caches the resonance dynamics.

Parameters
[in]fitDataTreethe data source

Implements LauAbsDPDynamics.

Definition at line 1777 of file LauIsobarDynamics.cc.

LauAbsResonance * LauIsobarDynamics::findResonance ( const TString &  name)
protectedvirtual

Retrieve the named resonance.

Parameters
[in]namethe name of the resonance to retrieve
Returns
the named resonance

Implements LauAbsDPDynamics.

Definition at line 651 of file LauIsobarDynamics.cc.

const LauAbsResonance * LauIsobarDynamics::findResonance ( const TString &  name) const
protectedvirtual

Retrieve the named resonance.

Parameters
[in]namethe name of the resonance to retrieve
Returns
the named resonance

Implements LauAbsDPDynamics.

Definition at line 672 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::flipHelicityForCPEigenstates ( Bool_t  boolean)
inline

Set the helicity flip flag for new amplitude components.

Parameters
[in]booleanthe helicity flip flag

Definition at line 495 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::generate ( )
virtual

Generate a toy MC signal event.

Returns
kTRUE if the event is successfully generated, kFALSE otherwise

Implements LauAbsDPDynamics.

Definition at line 1615 of file LauIsobarDynamics.cc.

Double_t LauIsobarDynamics::getASqMaxSetValue ( ) const
inline

Retrieve the maximum value of A squared to be used in the accept/reject.

Returns
the maximum value of A squared

Definition at line 254 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getASqMaxVarValue ( ) const
inline

Retrieve the maximum of A squared that has been found while generating.

Returns
the maximum of A squared that has been found

Definition at line 260 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getBelleNRAlpha ( ) const
inline

Retreive the alpha parameter for new Belle non-resonant components.

Returns
the alpha parameter

Definition at line 309 of file LauIsobarDynamics.hh.

virtual LauComplex LauIsobarDynamics::getDynamicAmp ( Int_t  resID) const
inlinevirtual

Retrieve the normalised dynamic part of the amplitude of the given amplitude component at the current point in the Dalitz plot.

Parameters
[in]resIDthe index of the amplitude component within the model
Returns
the amplitude of the given amplitude component

Implements LauAbsDPDynamics.

Definition at line 267 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getEventWeight ( )
virtual

Calculate the acceptance rate, for events with the current kinematics, when generating events according to the model.

Returns
the weight for the current kinematics

Implements LauAbsDPDynamics.

Definition at line 1874 of file LauIsobarDynamics.cc.

virtual const LauComplex& LauIsobarDynamics::getEvtDPAmp ( ) const
inlinevirtual

Retrieve the total amplitude for the current event.

Returns
the total amplitude

Implements LauAbsDPDynamics.

Definition at line 173 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtEff ( ) const
inlinevirtual

Retrieve the efficiency for the current event.

Returns
the efficiency for the current event

Implements LauAbsDPDynamics.

Definition at line 203 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtJacobian ( ) const
inlinevirtual

Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event.

Returns
the Jacobian for the current event

Implements LauAbsDPDynamics.

Definition at line 215 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtLikelihood ( ) const
inlinevirtual

Retrieve the likelihood for the current event.

evtLike_ = totAmp_.abs2()*eff_/DPNorm_

Returns
the likelihood for the current event

Implements LauAbsDPDynamics.

Definition at line 223 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtm13Sq ( ) const
inlinevirtual

Retrieve the invariant mass squared of the first and third daughters in the current event.

Returns
the invariant mass squared of the first and third daughters in the current event

Implements LauAbsDPDynamics.

Definition at line 179 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtm23Sq ( ) const
inlinevirtual

Retrieve the invariant mass squared of the second and third daughters in the current event.

Returns
the invariant mass squared of the second and third daughters in the current event

Implements LauAbsDPDynamics.

Definition at line 185 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtmPrime ( ) const
inlinevirtual

Retrieve the square Dalitz plot coordinate, m', for the current event.

Returns
the square Dalitz plot coordinate, m', for the current event

Implements LauAbsDPDynamics.

Definition at line 191 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtScfFraction ( ) const
inlinevirtual

Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the current event.

Returns
the self cross feed fraction for the current event

Implements LauAbsDPDynamics.

Definition at line 209 of file LauIsobarDynamics.hh.

virtual Double_t LauIsobarDynamics::getEvtthPrime ( ) const
inlinevirtual

Retrieve the square Dalitz plot coordinate, theta', for the current event.

Returns
the square Dalitz plot coordinate, theta', for the current event

Implements LauAbsDPDynamics.

Definition at line 197 of file LauIsobarDynamics.hh.

const std::vector< std::vector<LauComplex> >& LauIsobarDynamics::getFiFjEffSum ( ) const
inline

Retrieve the event-by-event running totals of efficiency corrected amplitude cross terms for all pairs of amplitude components.

Returns
the event-by-event running totals of amplitude cross terms with efficiency corrections applied

Definition at line 279 of file LauIsobarDynamics.hh.

const std::vector< std::vector<LauComplex> >& LauIsobarDynamics::getFiFjSum ( ) const
inline

Retrieve the event-by-event running totals of amplitude cross terms for all pairs of amplitude components.

Returns
the event-by-event running totals of amplitude cross terms

Definition at line 273 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getFlatteg1 ( ) const
inline

Retrieve the g1 parameter for new Flatte components.

Returns
the constant parameter g1

Definition at line 439 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getFlatteg2 ( ) const
inline

Retrieve the g2 parameter for new Flatte components.

Returns
the constant parameter g2

Definition at line 445 of file LauIsobarDynamics.hh.

const std::vector<Double_t>& LauIsobarDynamics::getFNorm ( ) const
inline

Retrieve the normalisation factors for the dynamic parts of the amplitudes for all of the amplitude components.

Returns
the normalisation factors

Definition at line 285 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSBackgroundMag ( ) const
inline

Retrieve the background magnitude for new LASS components.

Returns
the background magnitude

Definition at line 402 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSBackgroundPhase ( ) const
inline

Retrieve the background phase for new LASS components.

Returns
the background phase

Definition at line 408 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSCutOff ( ) const
inline

Retrieve the cutoff value for new LASS components.

Returns
the cutoff value

Definition at line 414 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSEffectiveRange ( ) const
inline

Retrieve the effective range for new LASS components.

Returns
the effective range

Definition at line 384 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSResonanceMag ( ) const
inline

Retrieve the resonance magnitude for new LASS components.

Returns
the resonance magnitude

Definition at line 390 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSResonancePhase ( ) const
inline

Retrieve the resonance phase for new LASS components.

Returns
the resonance phase

Definition at line 396 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getLASSScatteringLength ( ) const
inline

Retrieve the scattering length for new LASS components.

Returns
the scattering length

Definition at line 378 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getParBarrierRadius ( ) const
inline

Retrieve the radius of the barrier factor due to the parent to use for new amplitude components.

Returns
the barrier radius

Definition at line 489 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::getResBarrierRadius ( ) const
inline

Retrieve the radius of the barrier factor due to the resonance to use for new amplitude components.

Returns
the barrier radius

Definition at line 483 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::gotKMatrixMatch ( UInt_t  resAmpInt,
const TString &  propName 
) const
protected

Check whether a resonance is a K-matrix component of a given propagator.

Parameters
[in]resAmpIntthe index of the resonance within the model
[in]propNamethe name of the K-matrix propagator
Returns
true if the resonance is a component of the given propagator, otherwise return false

Definition at line 1551 of file LauIsobarDynamics.cc.

Bool_t LauIsobarDynamics::gotReweightedEvent ( )
virtual

Calculates whether an event with the current kinematics should be accepted in order to produce a distribution of events that matches the model e.g. when reweighting embedded data.

Uses the accept/reject method.

Returns
kTRUE if the event has been accepted, kFALSE otherwise

Implements LauAbsDPDynamics.

Definition at line 1856 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::initialise ( const std::vector< LauComplex > &  coeffs)
virtual

Initialise the Dalitz plot dynamics.

Parameters
[in]coeffsthe complex coefficients for the resonances

Implements LauAbsDPDynamics.

Definition at line 167 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::initialiseVectors ( )
protectedvirtual

Initialise the internal storage for this model.

Definition at line 294 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::initSummary ( )
protectedvirtual

Print a summary of the model to be used.

Definition at line 245 of file LauIsobarDynamics.cc.

Int_t LauIsobarDynamics::maxGenIterations ( ) const
inline

Retrieve the maximum number of iterations allowed when generating an event.

Returns
the maximum number of iterations allowed

Definition at line 145 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::removeCharge ( TString &  string) const
protectedvirtual

Remove the charge from the given particle name.

Parameters
[in,out]stringthe particle name

Definition at line 693 of file LauIsobarDynamics.cc.

LauComplex LauIsobarDynamics::resAmp ( Int_t  index)
protectedvirtual

Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz plot.

Parameters
[in]indexthe index of the amplitude component within the model

Definition at line 1378 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::setASqMaxValue ( Double_t  value)
inline

Set the maximum value of A squared to be used in the accept/reject.

Parameters
[in]valuethe new value

Definition at line 248 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setASqMaxVarValue ( Double_t  value)
inlineprotected

Set the maximum of A squared that has been found.

Parameters
[in]valuethe new value

Definition at line 542 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setBarrierRadii ( Double_t  resRadius,
Double_t  parRadius = 4.0,
LauAbsResonance::BarrierType  type = LauAbsResonance::BWPrimeBarrier 
)
inline

Set the parameters for the barrier factors for new resonances.

Parameters
[in]resRadiusthe radius due to the resonance
[in]parRadiusthe radius due to the parent
[in]typethe type of the barrier factor

Definition at line 453 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setBarrierType ( LauAbsResonance::BarrierType  type)
inline

Set the type of barrier factor to use for new amplitude components.

Parameters
[in]typethe type of barrier factor. Allowed types are: BWBarrier, BWPrimeBarrier and ExpBarrier.

Definition at line 477 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setBelleNRAlpha ( Double_t  alpha)
inline

Set the alpha parameter for new Belle non-resonant components.

Parameters
[in]alphathe alpha parameter value

Definition at line 303 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setDataEventNo ( UInt_t  iEvt)
virtual

Load the data for a given event.

Parameters
[in]iEvtthe number of the event

Reimplemented from LauAbsDPDynamics.

Definition at line 1706 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::setFFTerm ( UInt_t  index,
Double_t  realPart,
Double_t  imagPart 
)
protectedvirtual

Set the dynamic part of the amplitude for a given amplitude component at the current point in the Dalitz plot.

Parameters
[in]indexthe index of the amplitude component
[in]realPartthe real part of the amplitude
[in]imagPartthe imaginary part of the amplitude

Definition at line 1411 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::setFlatteg1 ( Double_t  g1)
inline

Set the g1 parameter for new Flatte components.

Parameters
[in]g1the coupling constant for channel 1

Definition at line 427 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setFlatteg2 ( Double_t  g2)
inline

Set the g2 parameter for new Flatte components.

Parameters
[in]g2the coupling constant for channel 2

Definition at line 433 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setFlatteParameters ( Double_t  g1,
Double_t  g2 
)
inline

Set the parameters for new Flatte components.

Parameters
[in]g1the coupling constant for channel 1
[in]g2the coupling constant for channel 2

Definition at line 421 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setIntegralBinWidths ( Double_t  m13BinWidth,
Double_t  m23BinWidth 
)
virtual

Set the widths of the bins to use when integrating across the Dalitz plot.

Parameters
[in]m13BinWidththe bin width to use when integrating over m13
[in]m23BinWidththe bin width to use when integrating over m23

Definition at line 1154 of file LauIsobarDynamics.cc.

void LauIsobarDynamics::setIntFileName ( const TString &  fileName)
inline

Set the name of the file to which to save the results of the integrals.

Parameters
[in]fileNamethe name of the file

Definition at line 67 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSBackgroundMag ( Double_t  B)
inline

Set the background magnitude for new LASS components.

Parameters
[in]Bthe background magnitude

Definition at line 360 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSBackgroundPhase ( Double_t  phiB)
inline

Set the background phase for new LASS components.

Parameters
[in]phiBthe background phase

Definition at line 366 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSCutOff ( Double_t  cutOff)
inline

Set the cutoff value for new LASS components.

Parameters
[in]cutOffthe cutoff value

Definition at line 372 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSEffectiveRange ( Double_t  r)
inline

Set the effective range for new LASS components.

Parameters
[in]rthe effective range

Definition at line 342 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSParameters ( Double_t  a,
Double_t  r,
Double_t  R = 1.0,
Double_t  phiR = 0.0,
Double_t  B = 1.0,
Double_t  phiB = 0.0,
Double_t  cutOff = 1.8 
)
inline

Set the parameters for new LASS components.

Parameters
[in]athe scattering length
[in]rthe effective range
[in]Rthe resonance magnitude
[in]phiRthe resonance phase
[in]Bthe background magnitude
[in]phiBthe background phase
[in]cutOffthe cutoff value

Definition at line 321 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSResonanceMag ( Double_t  R)
inline

Set the resonance magnitude for new LASS components.

Parameters
[in]Rthe resonance magnitude

Definition at line 348 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSResonancePhase ( Double_t  phiR)
inline

Set the resonance phase for new LASS components.

Parameters
[in]phiRthe resonance phase

Definition at line 354 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setLASSScatteringLength ( Double_t  a)
inline

Set the scattering length for new LASS components.

Parameters
[in]athe scattering length

Definition at line 336 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setParBarrierRadius ( Double_t  radius)
inline

Set the radius of the barrier factor due to the parent to use for new amplitude components.

Parameters
[in]radiusthe barrier radius

Definition at line 471 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::setResBarrierRadius ( Double_t  radius)
inline

Set the radius of the barrier factor due to the resonance to use for new amplitude components.

Parameters
[in]radiusthe barrier radius

Definition at line 465 of file LauIsobarDynamics.hh.

void LauIsobarDynamics::writeIntegralsFile ( )
protectedvirtual

Write the results of the integrals (and related information) to a file.

Definition at line 330 of file LauIsobarDynamics.cc.

Member Data Documentation

Double_t LauIsobarDynamics::ASq_
private

The value of A squared for the current event.

Definition at line 649 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::aSqMaxSet_
private

The maximum allowed value of A squared.

Definition at line 685 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::aSqMaxVar_
private

The maximum value of A squared that has been seen so far while generating.

Definition at line 688 of file LauIsobarDynamics.hh.

LauAbsResonance::BarrierType LauIsobarDynamics::barrierType_
private

The type of the barrier factor for new amplitude components.

Definition at line 754 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::BelleNRAlpha_
private

The alpha parameter for new Belle non-resonant components.

Definition at line 691 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeFlatteParameterg1_
private

Whether the default value of the Flatte parameter g1 has been changed.

Definition at line 742 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeFlatteParameterg2_
private

Whether the default value of the Flatte parameter g2 has been changed.

Definition at line 745 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSBackgroundMag_
private

Whether the default value of the LASS background magnitude has been changed.

Definition at line 727 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSBackgroundPhase_
private

Whether the default value of the LASS background phase has been changed.

Definition at line 730 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSCutOff_
private

Whether the default value of the LASS cutoff has been changed.

Definition at line 733 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSEffectiveRange_
private

Whether the default value of the LASS effective range has been changed.

Definition at line 718 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSResonanceMag_
private

Whether the default value of the LASS resonance magnitude has been changed.

Definition at line 721 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSResonancePhase_
private

Whether the default value of the LASS resonance phase has been changed.

Definition at line 724 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::changeLASSScatteringLength_
private

Whether the default value of the LASS scattering length has been changed.

Definition at line 715 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::eff_
private

The efficiency at the current point in the Dalitz plot.

Definition at line 640 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::evtLike_
private

The normalised likelihood for the current event.

Definition at line 652 of file LauIsobarDynamics.hh.

std::vector<LauComplex> LauIsobarDynamics::ff_
private

The dynamic part of the amplitude for each amplitude component at the current point in the Dalitz plot.

Definition at line 670 of file LauIsobarDynamics.hh.

std::vector< std::vector<LauComplex> > LauIsobarDynamics::fifjEffSum_
private

The event-by-event running total of efficiency corrected amplitude cross terms for each pair of amplitude components.

Calculated as the sum of ff_[i]*ff_[j]*efficiency for all events

Definition at line 661 of file LauIsobarDynamics.hh.

std::vector< std::vector<LauComplex> > LauIsobarDynamics::fifjSum_
private

The event-by-event running total of the amplitude cross terms for each pair of amplitude components.

Calculated as the sum of ff_[i]*ff_[j] for all events

Definition at line 667 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::FlatteParameterg1_
private

The constant parameter g1 for new Flatte components.

Definition at line 736 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::FlatteParameterg2_
private

The constant parameter g2 for new Flatte components.

Definition at line 739 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::flipHelicity_
private

The helicity flip flag for new amplitude components.

Definition at line 757 of file LauIsobarDynamics.hh.

std::vector<Double_t> LauIsobarDynamics::fNorm_
private

The normalisation factors for the dynamic parts of the amplitude for each amplitude component.

Definition at line 676 of file LauIsobarDynamics.hh.

std::vector<Double_t> LauIsobarDynamics::fSqSum_
private

The event-by-event running total of the dynamical amplitude squared for each amplitude component.

Definition at line 673 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::integralsDone_
private

Whether the integrals have been performed.

Definition at line 616 of file LauIsobarDynamics.hh.

TString LauIsobarDynamics::intFileName_
private

The name of the file to save integrals to.

Definition at line 619 of file LauIsobarDynamics.hh.

Int_t LauIsobarDynamics::iterationsMax_
private

The maximum allowed number of attempts when generating an event.

Definition at line 679 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::jacobian_
private

The Jacobian, for the transformation into square DP coordinates at the current point in the Dalitz plot.

Definition at line 646 of file LauIsobarDynamics.hh.

KMPropMap LauIsobarDynamics::kMatrixPropagators_
private

The K-matrix propagators.

Definition at line 595 of file LauIsobarDynamics.hh.

KMStringMap LauIsobarDynamics::kMatrixPropSet_
private

The names of the M-matrix components in the model mapped to their propagators.

Definition at line 598 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSBackgroundMag_
private

The background magnitude for new LASS components.

Definition at line 706 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSBackgroundPhase_
private

The background phase for new LASS components.

Definition at line 709 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSCutOff_
private

The cutoff value for new LASS components.

Definition at line 712 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSEffectiveRange_
private

The effective range for new LASS components.

Definition at line 697 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSResonanceMag_
private

The resonance magnitude for new LASS components.

Definition at line 700 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSResonancePhase_
private

The resonance phase for new LASS components.

Definition at line 703 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::LASSScatteringLength_
private

The scattering length for new LASS components.

Definition at line 694 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::m13BinWidth_
private

The bin width to use when integrating over m13.

Definition at line 622 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::m13Sq_
private

The invariant mass squared of the first and third daughters.

Definition at line 628 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::m23BinWidth_
private

The bin width to use when integrating over m23.

Definition at line 625 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::m23Sq_
private

The invariant mass squared of the second and third daughters.

Definition at line 631 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::mPrime_
private

The square Dalitz plot coordinate, m'.

Definition at line 634 of file LauIsobarDynamics.hh.

Int_t LauIsobarDynamics::nSigGenLoop_
private

The number of unsucessful attempts to generate an event so far.

Definition at line 682 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::parBarrierRadius_
private

The radius of the parent barrier factor for new amplitude components.

Definition at line 751 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::resBarrierRadius_
private

The radius of the resonance barrier factor for new amplitude components.

Definition at line 748 of file LauIsobarDynamics.hh.

std::vector<Int_t> LauIsobarDynamics::resIntAmp_
private

The index within the resonance maker for each amplitude component.

Definition at line 604 of file LauIsobarDynamics.hh.

std::vector<Int_t> LauIsobarDynamics::resPairAmp_
private

The index of the daughter not produced by the resonance for each amplitude component.

Definition at line 607 of file LauIsobarDynamics.hh.

std::vector<TString> LauIsobarDynamics::resTypAmp_
private

The resonance types of all of the amplitude components.

Definition at line 601 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::scfFraction_
private

The fraction of events that are poorly reconstructed (the self cross feed fraction) at the current point in the Dalitz plot.

Definition at line 643 of file LauIsobarDynamics.hh.

std::vector<LauAbsResonance*> LauIsobarDynamics::sigResonances_
private

The resonances in the model.

Definition at line 592 of file LauIsobarDynamics.hh.

Bool_t LauIsobarDynamics::symmetricalDP_
private

Whether the Dalitz plot is symmetrical.

Definition at line 613 of file LauIsobarDynamics.hh.

Double_t LauIsobarDynamics::thPrime_
private

The square Dalitz plot coordinate theta'.

Definition at line 637 of file LauIsobarDynamics.hh.

LauComplex LauIsobarDynamics::totAmp_
private

The total amplitude for the current event.

Definition at line 655 of file LauIsobarDynamics.hh.

std::vector<Int_t> LauIsobarDynamics::typDaug_
private

The PDG codes of the daughters.

Definition at line 610 of file LauIsobarDynamics.hh.


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