|
Laura++
3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
|
Class for defining signal dynamics using the isobar model.
More...
#include <LauIsobarDynamics.hh>
|
| LauIsobarDynamics (LauDaughters *daughters, LauAbsEffModel *effModel, LauAbsEffModel *scfFractionModel=0) |
| Constructor. More...
|
|
| LauIsobarDynamics (LauDaughters *daughters, LauAbsEffModel *effModel, LauTagCatScfFractionModelMap scfFractionModel) |
| Constructor. More...
|
|
virtual | ~LauIsobarDynamics () |
| Destructor.
|
|
void | initialise (const std::vector< LauComplex > &coeffs) |
| Initialise the Dalitz plot dynamics. More...
|
|
void | recalculateNormalisation () |
| recalculate Normalization
|
|
void | setIntFileName (const TString &fileName) |
| Set the name of the file to which to save the results of the integrals. More...
|
|
void | setIntegralBinWidths (const Double_t m13BinWidth, const Double_t m23BinWidth, const Double_t mPrimeBinWidth=0.001, const Double_t thPrimeBinWidth=0.001) |
| Set the widths of the bins to use when integrating across the Dalitz plot or square Dalitz plot. More...
|
|
void | setNarrowResonanceThreshold (const Double_t narrowWidth) |
| Set the value below which a resonance width is considered to be narrow. More...
|
|
void | setIntegralBinningFactor (const Double_t binningFactor) |
| Set the factor relating the width of a narrow resonance and the binning size in its integration region. More...
|
|
void | forceSymmetriseIntegration (const Bool_t force) |
| Force the symmetrisation of the integration in m13 <-> m23 for non-symmetric but flavour-conjugate final states. More...
|
|
LauAbsResonance * | addResonance (const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory=LauBlattWeisskopfFactor::Default) |
| Add a resonance to the Dalitz plot. More...
|
|
LauAbsResonance * | addIncoherentResonance (const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType) |
| Add an incoherent resonance to the Dalitz plot. More...
|
|
LauKMatrixPropagator * | defineKMatrixPropagator (const TString &propName, const TString ¶mFileName, Int_t resPairAmpInt, Int_t nChannels, Int_t nPoles, Int_t rowIndex=1) |
| Define a new K-matrix Propagator. More...
|
|
void | addKMatrixProdPole (const TString &poleName, const TString &propName, Int_t poleIndex, Bool_t useProdAdler=kFALSE) |
| Add a K-matrix production pole term to the model. More...
|
|
void | addKMatrixProdSVP (const TString &SVPName, const TString &propName, Int_t channelIndex, Bool_t useProdAdler=kFALSE) |
| Add a K-matrix slowly-varying part (SVP) term 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...
|
|
Bool_t | generate () |
| Generate a toy MC signal event. More...
|
|
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...
|
|
void | calcLikelihoodInfo (const UInt_t iEvt) |
| Calculate the likelihood (and all associated information) for the given event number. More...
|
|
void | calcLikelihoodInfo (const Double_t m13Sq, const Double_t m23Sq) |
| Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates. More...
|
|
void | calcLikelihoodInfo (const Double_t m13Sq, const Double_t m23Sq, const Int_t tagCat) |
| Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates and the tagging category. More...
|
|
void | calcExtraInfo (const Bool_t init=kFALSE) |
| Calculate the fit fractions, mean efficiency and total DP rate. More...
|
|
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...
|
|
Double_t | getEventWeight () |
| Calculate the acceptance rate, for events with the current kinematics, when generating events according to the model. More...
|
|
const LauComplex & | getEvtDPAmp () const |
| Retrieve the total amplitude for the current event. More...
|
|
Double_t | getEvtm13Sq () const |
| Retrieve the invariant mass squared of the first and third daughters in the current event. More...
|
|
Double_t | getEvtm23Sq () const |
| Retrieve the invariant mass squared of the second and third daughters in the current event. More...
|
|
Double_t | getEvtmPrime () const |
| Retrieve the square Dalitz plot coordinate, m', for the current event. More...
|
|
Double_t | getEvtthPrime () const |
| Retrieve the square Dalitz plot coordinate, theta', for the current event. More...
|
|
Double_t | getEvtEff () const |
| Retrieve the efficiency for the current event. More...
|
|
Double_t | getEvtScfFraction () const |
| Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the current event. More...
|
|
Double_t | getEvtJacobian () const |
| Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event. More...
|
|
Double_t | getEvtIntensity () const |
| Retrieve the total intensity multiplied by the efficiency for the current event. More...
|
|
Double_t | getEvtLikelihood () const |
| Retrieve the likelihood for the current event. More...
|
|
LauComplex | getDynamicAmp (const 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...
|
|
LauComplex | getFullAmplitude (const Int_t resID) const |
| Retrieve the Amplitude of resonance resID. 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...
|
|
void | fillDataTree (const LauFitDataTree &fitDataTree) |
| Fill the internal data structure that caches the resonance dynamics. More...
|
|
void | modifyDataTree () |
| Recache the amplitude values for those that have changed.
|
|
Bool_t | hasResonance (const TString &resName) const |
| Check whether this model includes a named resonance. More...
|
|
Int_t | resonanceIndex (const TString &resName) const |
| Retrieve the index for the given resonance. More...
|
|
TString | getConjResName (const TString &resName) const |
| Retrieve the name of the charge conjugate of a named resonance. More...
|
|
const LauAbsResonance * | findResonance (const TString &resName) const |
| Retrieve the named resonance. More...
|
|
const LauAbsResonance * | getResonance (const UInt_t resIndex) const |
| Retrieve a resonance by its index. More...
|
|
void | updateCoeffs (const std::vector< LauComplex > &coeffs) |
| Update the complex coefficients for the resonances. More...
|
|
void | collateResonanceParameters () |
| Collate the resonance parameters to initialise (or re-initialise) the model. More...
|
|
void | flipHelicityForCPEigenstates (const Bool_t boolean) |
| Set the helicity flip flag for new amplitude components. More...
|
|
const LauParameter & | getMeanEff () const |
| Retrieve the mean efficiency across the Dalitz plot. More...
|
|
const LauParameter & | getDPRate () const |
| Retrieve the overall Dalitz plot rate. More...
|
|
const LauParArray & | getFitFractions () const |
| Retrieve the fit fractions for the amplitude components. More...
|
|
const LauParArray & | getFitFractionsEfficiencyUncorrected () const |
| Retrieve the fit fractions for the amplitude components. More...
|
|
UInt_t | getnTotAmp () const |
| Retrieve the total number of amplitude components. More...
|
|
UInt_t | getnCohAmp () const |
| Retrieve the number of coherent amplitude components. More...
|
|
UInt_t | getnIncohAmp () const |
| Retrieve the number of incoherent amplitude components. More...
|
|
Double_t | getDPNorm () const |
| Retrieve the normalisation factor for the log-likelihood function. More...
|
|
const LauDaughters * | getDaughters () const |
| Retrieve the daughters. More...
|
|
const LauKinematics * | getKinematics () const |
| Retrieve the Dalitz plot kinematics. More...
|
|
LauKinematics * | getKinematics () |
| Retrieve the Dalitz plot kinematics. More...
|
|
const LauAbsEffModel * | getEffModel () const |
| Retrieve the model for the efficiency across the Dalitz plot. More...
|
|
Bool_t | usingScfModel () const |
| Check whether a self cross feed fraction model is being used. More...
|
|
const std::vector< LauParameter > & | getExtraParameters () const |
| Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions) More...
|
|
std::vector< LauParameter * > & | getFloatingParameters () |
| Retrieve the floating parameters of the resonance models. More...
|
|
void | calculateRhoOmegaFitFractions (const Bool_t calcFF) |
| Whether to calculate separate rho and omega fit-fractions from LauRhoOmegaMix.
|
|
|
void | initSummary () |
| Print a summary of the model to be used.
|
|
void | initialiseVectors () |
| Initialise the internal storage for this model.
|
|
void | resetNormVectors () |
| Zero the various values used to store integrals info.
|
|
void | calcDPNormalisation () |
| Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
|
|
std::vector< std::pair< Double_t, Double_t > > | formGapsFromRegions (const std::vector< std::pair< Double_t, Double_t >> ®ions, const Double_t min, const Double_t max) const |
| Form the regions that are produced by the spaces between narrow resonances. More...
|
|
void | cullNullRegions (std::vector< LauDPPartialIntegralInfo * > ®ions) const |
| Removes entries in the vector of LauDPPartialIntegralInfo* that are null. More...
|
|
LauDPPartialIntegralInfo * | newDPIntegrationRegion (const Double_t minm13, const Double_t maxm13, const Double_t minm23, const Double_t maxm23, const Double_t m13BinWidth, const Double_t m23BinWidth, const Double_t precision, const UInt_t nAmp, const UInt_t nIncohAmp) const |
| Wrapper for LauDPPartialIntegralInfo constructor. More...
|
|
void | correctDPOverlap (std::vector< std::pair< Double_t, Double_t >> ®ions, const std::vector< Double_t > &binnings) const |
| Correct regions to ensure that the finest integration grid takes precedence. More...
|
|
std::vector< LauDPPartialIntegralInfo * > | m23IntegrationRegions (const std::vector< std::pair< Double_t, Double_t >> &m13Regions, const std::vector< std::pair< Double_t, Double_t >> &m23Regions, const std::vector< Double_t > &m13Binnings, const std::vector< Double_t > &m23Binnings, const Double_t precision, const Double_t defaultBinning) const |
| Create the integration grid objects for the m23 narrow resonance regions, including the overlap regions with the m13 narrow resonances. More...
|
|
std::vector< LauDPPartialIntegralInfo * > | m13IntegrationRegions (const std::vector< std::pair< Double_t, Double_t >> &m13Regions, const std::vector< std::pair< Double_t, Double_t >> &m23Regions, const std::vector< Double_t > &m13Binnings, const Double_t precision, const Double_t defaultBinning) const |
| Create the integration grid objects for the m13 narrow resonance regions, excluding the overlap regions with the m23 narrow resonances. More...
|
|
void | calcDPNormalisationScheme () |
| Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
|
|
void | findIntegralsToBeRecalculated () |
| Determine which amplitudes and integrals need to be recalculated.
|
|
void | calcDPPartialIntegral (LauDPPartialIntegralInfo *intInfo) |
| Calculate the Dalitz plot normalisation integrals over a given range. More...
|
|
void | writeIntegralsFile () |
| Write the results of the integrals (and related information) to a file.
|
|
void | setFFTerm (const UInt_t index, const Double_t realPart, const Double_t imagPart) |
| Set the dynamic part of the amplitude for a given amplitude component at the current point in the Dalitz plot. More...
|
|
void | setIncohIntenTerm (const UInt_t index, const Double_t value) |
| Set the dynamic part of the intensity for a given incoherent amplitude component at the current point in the Dalitz plot. More...
|
|
void | calculateAmplitudes () |
| Calculate the amplitudes for all resonances for the current kinematics.
|
|
void | calculateAmplitudes (LauDPPartialIntegralInfo *intInfo, const UInt_t m13Point, const UInt_t m23Point) |
| Calculate or retrieve the cached value of the amplitudes for all resonances at the specified integration grid point. More...
|
|
void | addGridPointToIntegrals (const Double_t weight) |
| Add the amplitude values (with the appropriate weight) at the current grid point to the running integral values. More...
|
|
void | calcTotalAmp (const Bool_t useEff) |
| Calculate the total Dalitz plot amplitude at the current point in the Dalitz plot. More...
|
|
Double_t | retrieveEfficiency () |
| Obtain the efficiency of the current event from the model. More...
|
|
Double_t | retrieveScfFraction (Int_t tagCat) |
| Obtain the self cross feed fraction of the current event from the model. More...
|
|
void | setASqMaxVarValue (Double_t value) |
| Set the maximum of A squared that has been found. More...
|
|
Double_t | calcSigDPNorm () |
| Calculate the normalisation factor for the log-likelihood function. More...
|
|
LauComplex | resAmp (const UInt_t index) |
| Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz plot. More...
|
|
Double_t | incohResAmp (const UInt_t index) |
| Calculate the dynamic part of the intensity for a given incoherent component at the current point in the Dalitz plot. More...
|
|
void | setDataEventNo (UInt_t iEvt) |
| Load the data for a given event. More...
|
|
LauAbsResonance * | findResonance (const TString &resName) |
| Retrieve the named resonance. More...
|
|
LauAbsResonance * | getResonance (const UInt_t resIndex) |
| Retrieve a resonance by its index. More...
|
|
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...
|
|
|
typedef std::map< TString, LauKMatrixPropagator * > | KMPropMap |
| The type used for containing the K-matrix propagators.
|
|
typedef std::map< TString, TString > | KMStringMap |
| The type used for mapping K-matrix components to their propagators.
|
|
|
LauDaughters * | daughters_ |
| The daughters of the decay.
|
|
LauKinematics * | kinematics_ |
| The kinematics of the decay.
|
|
LauAbsEffModel * | effModel_ |
| The efficiency model across the Dalitz plot.
|
|
LauTagCatScfFractionModelMap | scfFractionModel_ |
| The self cross feed fraction models across the Dalitz plot. More...
|
|
UInt_t | nAmp_ |
| The number of amplitude components.
|
|
UInt_t | nIncohAmp_ |
| The number of incoherent amplitude components.
|
|
std::vector< LauComplex > | Amp_ |
| The complex coefficients for the amplitude components.
|
|
Double_t | DPNorm_ |
| The normalisation factor for the log-likelihood function.
|
|
LauParArray | fitFrac_ |
| The fit fractions for the amplitude components.
|
|
LauParArray | fitFracEffUnCorr_ |
| The efficiency-uncorrected fit fractions for the amplitude components.
|
|
LauParameter | DPRate_ |
| The overall Dalitz plot rate.
|
|
LauParameter | meanDPEff_ |
| The mean efficiency across the Dalitz plot.
|
|
std::vector< LauCacheData * > | data_ |
| The cached data for all events.
|
|
LauCacheData * | currentEvent_ |
| The cached data for the current event.
|
|
std::vector< LauParameter > | extraParameters_ |
| any extra parameters/quantities (e.g. K-matrix total fit fractions)
|
|
std::vector< LauAbsResonance * > | sigResonances_ |
| The resonances in the model.
|
|
std::vector< LauAbsIncohRes * > | sigIncohResonances_ |
| The incoherent resonances in the model.
|
|
KMPropMap | kMatrixPropagators_ |
| The K-matrix propagators.
|
|
KMStringMap | kMatrixPropSet_ |
| The names of the M-matrix components in the model mapped to their propagators.
|
|
std::vector< TString > | resTypAmp_ |
| The resonance types of all of the amplitude components.
|
|
std::vector< Int_t > | resPairAmp_ |
| The index of the daughter not produced by the resonance for each amplitude component.
|
|
std::vector< TString > | incohResTypAmp_ |
| The resonance types of all of the incoherent amplitude components.
|
|
std::vector< Int_t > | incohResPairAmp_ |
| The index of the daughter not produced by the resonance for each incoherent amplitude component.
|
|
std::vector< Int_t > | typDaug_ |
| The PDG codes of the daughters.
|
|
Bool_t | symmetricalDP_ |
| Whether the Dalitz plot is symmetrical.
|
|
Bool_t | fullySymmetricDP_ |
| Whether the Dalitz plot is fully symmetric.
|
|
Bool_t | flavConjDP_ |
| Whether the Dalitz plot is a flavour-conjugate final state.
|
|
Bool_t | integralsDone_ |
| Whether the integrals have been performed.
|
|
Bool_t | normalizationSchemeDone_ |
| Whether the scheme for the integration has been determined.
|
|
Bool_t | forceSymmetriseIntegration_ |
| Force the symmetrisation of the integration in m13 <-> m23 for non-symmetric but flavour-conjugate final states.
|
|
std::vector< LauDPPartialIntegralInfo * > | dpPartialIntegralInfo_ |
| The storage of the integration scheme.
|
|
TString | intFileName_ |
| The name of the file to save integrals to.
|
|
Double_t | m13BinWidth_ |
| The bin width to use when integrating over m13.
|
|
Double_t | m23BinWidth_ |
| The bin width to use when integrating over m23.
|
|
Double_t | mPrimeBinWidth_ |
| The bin width to use when integrating over mPrime.
|
|
Double_t | thPrimeBinWidth_ |
| The bin width to use when integrating over thetaPrime.
|
|
Double_t | narrowWidth_ |
| The value below which a resonance width is considered to be narrow.
|
|
Double_t | binningFactor_ |
| The factor relating the width of the narrowest resonance and the binning size.
|
|
Double_t | m13Sq_ |
| The invariant mass squared of the first and third daughters.
|
|
Double_t | m23Sq_ |
| The invariant mass squared of the second and third daughters.
|
|
Double_t | mPrime_ |
| The square Dalitz plot coordinate, m'.
|
|
Double_t | thPrime_ |
| The square Dalitz plot coordinate theta'.
|
|
Int_t | tagCat_ |
| The tagging category.
|
|
Double_t | eff_ |
| The efficiency at the current point in the Dalitz plot.
|
|
Double_t | scfFraction_ |
| The fraction of events that are poorly reconstructed (the self cross feed fraction) at the current point in the Dalitz plot.
|
|
Double_t | jacobian_ |
| The Jacobian, for the transformation into square DP coordinates at the current point in the Dalitz plot.
|
|
Double_t | ASq_ |
| The value of A squared for the current event.
|
|
Double_t | evtLike_ |
| The normalised likelihood for the current event.
|
|
LauComplex | totAmp_ |
| The total amplitude for the current event.
|
|
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< LauComplex > | ff_ |
| The dynamic part of the amplitude for each amplitude component at the current point in the Dalitz plot.
|
|
std::vector< Double_t > | incohInten_ |
| The dynamic part of the intensity for each incoherent amplitude component at the current point in the Dalitz plot.
|
|
std::vector< Double_t > | fSqSum_ |
| The event-by-event running total of the dynamical amplitude squared for each amplitude component.
|
|
std::vector< Double_t > | fSqEffSum_ |
| The event-by-event running total of the dynamical amplitude squared for each amplitude component.
|
|
std::vector< Double_t > | fNorm_ |
| The normalisation factors for the dynamic parts of the amplitude for each amplitude component.
|
|
Int_t | iterationsMax_ |
| The maximum allowed number of attempts when generating an event.
|
|
Int_t | nSigGenLoop_ |
| The number of unsucessful attempts to generate an event so far.
|
|
Double_t | aSqMaxSet_ |
| The maximum allowed value of A squared.
|
|
Double_t | aSqMaxVar_ |
| The maximum value of A squared that has been seen so far while generating.
|
|
Bool_t | aSqMaxAuto_ { kTRUE } |
| Flag to generate aSqMaxSet_ once generate is called.
|
|
Bool_t | flipHelicity_ |
| The helicity flip flag for new amplitude components.
|
|
Bool_t | recalcNormalisation_ |
| Flag to recalculate the normalisation.
|
|
std::vector< LauParameter * > | resonancePars_ |
| List of floating resonance parameters.
|
|
std::vector< Double_t > | resonanceParValues_ |
| List of floating resonance parameter values from previous calculation.
|
|
std::vector< std::vector< UInt_t > > | resonanceParResIndex_ |
| Indices in sigResonances_ to point to the corresponding signal resonance(s) for each floating parameter.
|
|
std::set< UInt_t > | integralsToBeCalculated_ |
| Resonance indices for which the amplitudes and integrals should be recalculated.
|
|
Bool_t | calculateRhoOmegaFitFractions_ |
| Whether to calculate separate rho and omega fit fractions from the LauRhoOmegaMix model.
|
|
Class for defining signal dynamics using the isobar model.
Definition at line 53 of file LauIsobarDynamics.hh.
◆ ToyMCStatus
The possible statuses for toy MC generation.
Enumerator |
---|
GenOK | Generation completed OK
|
MaxIterError | Maximum allowed number of iterations completed without success (ASqMax is too high)
|
ASqMaxError | An amplitude squared value was returned that was larger than the maximum expected (ASqMax is too low)
|
Definition at line 60 of file LauIsobarDynamics.hh.
◆ LauIsobarDynamics() [1/2]
Constructor.
- Parameters
-
[in] | daughters | the daughters of the decay |
[in] | effModel | the model to describe the efficiency across the Dalitz plot |
[in] | scfFractionModel | the model to describe the fraction of poorly constructed events (the self cross feed fraction) across the Dalitz plot |
Definition at line 65 of file LauIsobarDynamics.cc.
◆ LauIsobarDynamics() [2/2]
Constructor.
- Parameters
-
[in] | daughters | the daughters of the decay |
[in] | effModel | the model to describe the efficiency across the Dalitz plot |
[in] | scfFractionModel | the 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 130 of file LauIsobarDynamics.cc.
◆ addGridPointToIntegrals()
void LauIsobarDynamics::addGridPointToIntegrals |
( |
const Double_t |
weight | ) |
|
|
protected |
Add the amplitude values (with the appropriate weight) at the current grid point to the running integral values.
- Parameters
-
[in] | weight | the weight to apply |
Definition at line 2080 of file LauIsobarDynamics.cc.
◆ addIncoherentResonance()
Add an incoherent resonance to the Dalitz plot.
NB the stored order of resonances is:
- Parameters
-
[in] | resName | the name of the resonant particle |
[in] | resPairAmpInt | the index of the daughter not produced by the resonance |
[in] | resType | the model for the resonance dynamics |
- Returns
- the newly created resonance
Definition at line 798 of file LauIsobarDynamics.cc.
◆ addKMatrixProdPole()
void LauIsobarDynamics::addKMatrixProdPole |
( |
const TString & |
poleName, |
|
|
const TString & |
propName, |
|
|
Int_t |
poleIndex, |
|
|
Bool_t |
useProdAdler = kFALSE |
|
) |
| |
Add a K-matrix production pole term to the model.
NB the stored order of resonances is:
- Parameters
-
[in] | poleName | the name of the pole |
[in] | propName | the name of the propagator to use |
[in] | poleIndex | the index of the pole within the propagator |
[in] | useProdAdler | boolean to turn on/off the production Adler zero factor (default = off) |
Definition at line 887 of file LauIsobarDynamics.cc.
◆ addKMatrixProdSVP()
void LauIsobarDynamics::addKMatrixProdSVP |
( |
const TString & |
SVPName, |
|
|
const TString & |
propName, |
|
|
Int_t |
channelIndex, |
|
|
Bool_t |
useProdAdler = kFALSE |
|
) |
| |
Add a K-matrix slowly-varying part (SVP) term to the model.
NB the stored order of resonances is:
- Parameters
-
[in] | SVPName | the name of the term |
[in] | propName | the name of the propagator to use |
[in] | channelIndex | the index of the channel within the propagator |
[in] | useProdAdler | boolean to turn on/off the production Adler zero factor (default = off) |
Definition at line 943 of file LauIsobarDynamics.cc.
◆ addResonance()
Add a resonance to the Dalitz plot.
NB the stored order of resonances is:
- Parameters
-
[in] | resName | the name of the resonant particle |
[in] | resPairAmpInt | the index of the daughter not produced by the resonance |
[in] | resType | the model for the resonance dynamics |
[in] | bwCategory | the Blatt-Weisskopf barrier factor category |
- Returns
- the newly created resonance
Definition at line 722 of file LauIsobarDynamics.cc.
◆ calcDPPartialIntegral()
Calculate the Dalitz plot normalisation integrals over a given range.
- Parameters
-
[in] | intInfo | the integration information object |
Definition at line 1711 of file LauIsobarDynamics.cc.
◆ calcExtraInfo()
void LauIsobarDynamics::calcExtraInfo |
( |
const Bool_t |
init = kFALSE | ) |
|
Calculate the fit fractions, mean efficiency and total DP rate.
- Parameters
-
[in] | init | whether the calculated values should be stored as the initial/generated values or the fitted values |
Definition at line 2199 of file LauIsobarDynamics.cc.
◆ calcLikelihoodInfo() [1/3]
void LauIsobarDynamics::calcLikelihoodInfo |
( |
const Double_t |
m13Sq, |
|
|
const Double_t |
m23Sq |
|
) |
| |
Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates.
- Parameters
-
[in] | m13Sq | the invariant mass squared of the first and third daughters |
[in] | m23Sq | the invariant mass squared of the second and third daughters |
Definition at line 2756 of file LauIsobarDynamics.cc.
◆ calcLikelihoodInfo() [2/3]
void LauIsobarDynamics::calcLikelihoodInfo |
( |
const Double_t |
m13Sq, |
|
|
const Double_t |
m23Sq, |
|
|
const Int_t |
tagCat |
|
) |
| |
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] | m13Sq | the invariant mass squared of the first and third daughters |
[in] | m23Sq | the invariant mass squared of the second and third daughters |
[in] | tagCat | the tagging category |
Definition at line 2761 of file LauIsobarDynamics.cc.
◆ calcLikelihoodInfo() [3/3]
void LauIsobarDynamics::calcLikelihoodInfo |
( |
const UInt_t |
iEvt | ) |
|
Calculate the likelihood (and all associated information) for the given event number.
- Parameters
-
Definition at line 2725 of file LauIsobarDynamics.cc.
◆ calcSigDPNorm()
Double_t LauIsobarDynamics::calcSigDPNorm |
( |
| ) |
|
|
protected |
Calculate the normalisation factor for the log-likelihood function.
- Returns
- the normalisation factor
Definition at line 2527 of file LauIsobarDynamics.cc.
◆ calcTotalAmp()
void LauIsobarDynamics::calcTotalAmp |
( |
const Bool_t |
useEff | ) |
|
|
protected |
Calculate the total Dalitz plot amplitude at the current point in the Dalitz plot.
- Parameters
-
[in] | useEff | whether to apply efficiency corrections |
Definition at line 2038 of file LauIsobarDynamics.cc.
◆ calculateAmplitudes()
void LauIsobarDynamics::calculateAmplitudes |
( |
LauDPPartialIntegralInfo * |
intInfo, |
|
|
const UInt_t |
m13Point, |
|
|
const UInt_t |
m23Point |
|
) |
| |
|
protected |
Calculate or retrieve the cached value of the amplitudes for all resonances at the specified integration grid point.
- Parameters
-
[in,out] | intInfo | the integration information object |
[in] | m13Point | the grid index in m13 |
[in] | m23Point | the grid index in m23 |
Definition at line 1761 of file LauIsobarDynamics.cc.
◆ checkToyMC()
Check the status of the toy MC generation.
- Parameters
-
[in] | printErrorMessages | whether error messages should be printed |
[in] | printInfoMessages | whether info messages should be printed |
- Returns
- the status of the toy MC generation
Definition at line 2637 of file LauIsobarDynamics.cc.
◆ collateResonanceParameters()
void LauIsobarDynamics::collateResonanceParameters |
( |
| ) |
|
Collate the resonance parameters to initialise (or re-initialise) the model.
NB: This has been factored out of the initialise() method to allow for use in the importation of parameters in LauAbsFitModel
Definition at line 273 of file LauIsobarDynamics.cc.
◆ correctDPOverlap()
void LauIsobarDynamics::correctDPOverlap |
( |
std::vector< std::pair< Double_t, Double_t >> & |
regions, |
|
|
const std::vector< Double_t > & |
binnings |
|
) |
| const |
|
protected |
Correct regions to ensure that the finest integration grid takes precedence.
- Parameters
-
[in] | regions | the windows in invariant mass |
[in] | binnings | the corresponding binnings for each window |
Definition at line 1152 of file LauIsobarDynamics.cc.
◆ cullNullRegions()
Removes entries in the vector of LauDPPartialIntegralInfo* that are null.
- Parameters
-
[in] | regions | the list of region pointers |
Definition at line 1146 of file LauIsobarDynamics.cc.
◆ defineKMatrixPropagator()
LauKMatrixPropagator * LauIsobarDynamics::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.
- Parameters
-
[in] | propName | the name of the propagator |
[in] | paramFileName | the file that defines the propagator |
[in] | resPairAmpInt | the index of the bachelor |
[in] | nChannels | the number of channels |
[in] | nPoles | the number of poles |
[in] | rowIndex | the index of the row to be used when summing over all amplitude channels: S-wave corresponds to rowIndex = 1. |
Definition at line 851 of file LauIsobarDynamics.cc.
◆ fillDataTree()
void LauIsobarDynamics::fillDataTree |
( |
const LauFitDataTree & |
fitDataTree | ) |
|
Fill the internal data structure that caches the resonance dynamics.
- Parameters
-
[in] | fitDataTree | the data source |
Definition at line 2829 of file LauIsobarDynamics.cc.
◆ findResonance() [1/2]
LauAbsResonance * LauIsobarDynamics::findResonance |
( |
const TString & |
resName | ) |
|
|
protected |
Retrieve the named resonance.
- Parameters
-
[in] | resName | the name of the resonance to retrieve |
- Returns
- the requested resonance
Definition at line 1069 of file LauIsobarDynamics.cc.
◆ findResonance() [2/2]
const LauAbsResonance * LauIsobarDynamics::findResonance |
( |
const TString & |
resName | ) |
const |
Retrieve the named resonance.
- Parameters
-
[in] | resName | the name of the resonance to retrieve |
- Returns
- the requested resonance
Definition at line 1081 of file LauIsobarDynamics.cc.
◆ flipHelicityForCPEigenstates()
void LauIsobarDynamics::flipHelicityForCPEigenstates |
( |
const Bool_t |
boolean | ) |
|
|
inline |
Set the helicity flip flag for new amplitude components.
- Parameters
-
[in] | boolean | the helicity flip flag |
Definition at line 483 of file LauIsobarDynamics.hh.
◆ forceSymmetriseIntegration()
void LauIsobarDynamics::forceSymmetriseIntegration |
( |
const Bool_t |
force | ) |
|
|
inline |
Force the symmetrisation of the integration in m13 <-> m23 for non-symmetric but flavour-conjugate final states.
This can be necessary for time-dependent fits (where interference terms between A and Abar need to be integrated)
- Parameters
-
[in] | force | toggle forcing symmetrisation of the integration for apparently flavour-conjugate final states |
Definition at line 151 of file LauIsobarDynamics.hh.
◆ formGapsFromRegions()
std::vector< std::pair< Double_t, Double_t > > LauIsobarDynamics::formGapsFromRegions |
( |
const std::vector< std::pair< Double_t, Double_t >> & |
regions, |
|
|
const Double_t |
min, |
|
|
const Double_t |
max |
|
) |
| const |
|
protected |
Form the regions that are produced by the spaces between narrow resonances.
- Parameters
-
[in] | regions | the regions defined around narrow resonances |
[in] | min | the minimum value of the invariant mass |
[in] | max | the maximum value of the invariant mass |
- Returns
- vector of pointers to LauDPPartialIntegralInfo objects that contain the individual coarse regions
Definition at line 1125 of file LauIsobarDynamics.cc.
◆ generate()
Bool_t LauIsobarDynamics::generate |
( |
| ) |
|
Generate a toy MC signal event.
- Returns
- kTRUE if the event is successfully generated, kFALSE otherwise
Definition at line 2566 of file LauIsobarDynamics.cc.
◆ getASqMaxSetValue()
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 250 of file LauIsobarDynamics.hh.
◆ getASqMaxVarValue()
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 256 of file LauIsobarDynamics.hh.
◆ getConjResName()
TString LauIsobarDynamics::getConjResName |
( |
const TString & |
resName | ) |
const |
Retrieve the name of the charge conjugate of a named resonance.
- Parameters
-
- Returns
- the name of the charge conjugate
Definition at line 2974 of file LauIsobarDynamics.cc.
◆ getDaughters()
const LauDaughters* LauIsobarDynamics::getDaughters |
( |
| ) |
const |
|
inline |
◆ getDPNorm()
Double_t LauIsobarDynamics::getDPNorm |
( |
| ) |
const |
|
inline |
Retrieve the normalisation factor for the log-likelihood function.
- Returns
- the normalisation factor
Definition at line 534 of file LauIsobarDynamics.hh.
◆ getDPRate()
Retrieve the overall Dalitz plot rate.
- Returns
- the overall Dalitz plot rate
Definition at line 495 of file LauIsobarDynamics.hh.
◆ getDynamicAmp()
LauComplex LauIsobarDynamics::getDynamicAmp |
( |
const Int_t |
resID | ) |
const |
|
inline |
Retrieve the normalised dynamic part of the amplitude of the given amplitude component at the current point in the Dalitz plot.
- Parameters
-
[in] | resID | the index of the component within the model |
- Returns
- the amplitude of the given component
Definition at line 389 of file LauIsobarDynamics.hh.
◆ getEffModel()
Retrieve the model for the efficiency across the Dalitz plot.
- Returns
- the efficiency model
Definition at line 558 of file LauIsobarDynamics.hh.
◆ getEventWeight()
Double_t LauIsobarDynamics::getEventWeight |
( |
| ) |
|
Calculate the acceptance rate, for events with the current kinematics, when generating events according to the model.
- Returns
- the weight for the current kinematics
Definition at line 2940 of file LauIsobarDynamics.cc.
◆ getEvtDPAmp()
const LauComplex& LauIsobarDynamics::getEvtDPAmp |
( |
| ) |
const |
|
inline |
Retrieve the total amplitude for the current event.
- Returns
- the total amplitude
Definition at line 325 of file LauIsobarDynamics.hh.
◆ getEvtEff()
Double_t LauIsobarDynamics::getEvtEff |
( |
| ) |
const |
|
inline |
Retrieve the efficiency for the current event.
- Returns
- the efficiency for the current event
Definition at line 355 of file LauIsobarDynamics.hh.
◆ getEvtIntensity()
Double_t LauIsobarDynamics::getEvtIntensity |
( |
| ) |
const |
|
inline |
Retrieve the total intensity multiplied by the efficiency for the current event.
- Returns
- the total intensity multiplied by the efficiency for the current event
Definition at line 373 of file LauIsobarDynamics.hh.
◆ getEvtJacobian()
Double_t LauIsobarDynamics::getEvtJacobian |
( |
| ) |
const |
|
inline |
Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event.
- Returns
- the Jacobian for the current event
Definition at line 367 of file LauIsobarDynamics.hh.
◆ getEvtLikelihood()
Double_t LauIsobarDynamics::getEvtLikelihood |
( |
| ) |
const |
|
inline |
Retrieve the likelihood for the current event.
The likelihood is the normalised total intensity: evtLike_ = ASq_/DPNorm_
- Returns
- the likelihood for the current event
Definition at line 382 of file LauIsobarDynamics.hh.
◆ getEvtm13Sq()
Double_t LauIsobarDynamics::getEvtm13Sq |
( |
| ) |
const |
|
inline |
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
Definition at line 331 of file LauIsobarDynamics.hh.
◆ getEvtm23Sq()
Double_t LauIsobarDynamics::getEvtm23Sq |
( |
| ) |
const |
|
inline |
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
Definition at line 337 of file LauIsobarDynamics.hh.
◆ getEvtmPrime()
Double_t LauIsobarDynamics::getEvtmPrime |
( |
| ) |
const |
|
inline |
Retrieve the square Dalitz plot coordinate, m', for the current event.
- Returns
- the square Dalitz plot coordinate, m', for the current event
Definition at line 343 of file LauIsobarDynamics.hh.
◆ getEvtScfFraction()
Double_t LauIsobarDynamics::getEvtScfFraction |
( |
| ) |
const |
|
inline |
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
Definition at line 361 of file LauIsobarDynamics.hh.
◆ getEvtthPrime()
Double_t LauIsobarDynamics::getEvtthPrime |
( |
| ) |
const |
|
inline |
Retrieve the square Dalitz plot coordinate, theta', for the current event.
- Returns
- the square Dalitz plot coordinate, theta', for the current event
Definition at line 349 of file LauIsobarDynamics.hh.
◆ getExtraParameters()
const std::vector<LauParameter>& LauIsobarDynamics::getExtraParameters |
( |
| ) |
const |
|
inline |
Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions)
- Returns
- any extra parameters
Definition at line 570 of file LauIsobarDynamics.hh.
◆ getFiFjEffSum()
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 414 of file LauIsobarDynamics.hh.
◆ getFiFjSum()
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 408 of file LauIsobarDynamics.hh.
◆ getFitFractions()
const LauParArray& LauIsobarDynamics::getFitFractions |
( |
| ) |
const |
|
inline |
Retrieve the fit fractions for the amplitude components.
- Returns
- the fit fractions
Definition at line 501 of file LauIsobarDynamics.hh.
◆ getFitFractionsEfficiencyUncorrected()
const LauParArray& LauIsobarDynamics::getFitFractionsEfficiencyUncorrected |
( |
| ) |
const |
|
inline |
Retrieve the fit fractions for the amplitude components.
- Returns
- the fit fractions
Definition at line 507 of file LauIsobarDynamics.hh.
◆ getFloatingParameters()
std::vector<LauParameter*>& LauIsobarDynamics::getFloatingParameters |
( |
| ) |
|
|
inline |
Retrieve the floating parameters of the resonance models.
- Returns
- the list of floating parameters
Definition at line 576 of file LauIsobarDynamics.hh.
◆ getFNorm()
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 420 of file LauIsobarDynamics.hh.
◆ getFullAmplitude()
LauComplex LauIsobarDynamics::getFullAmplitude |
( |
const Int_t |
resID | ) |
const |
|
inline |
Retrieve the Amplitude of resonance resID.
- Parameters
-
[in] | resID | the index of the component within the model |
- Returns
- the amplitude of the given component
Definition at line 399 of file LauIsobarDynamics.hh.
◆ getKinematics() [1/2]
Retrieve the Dalitz plot kinematics.
- Returns
- the Dalitz plot kinematics
Definition at line 552 of file LauIsobarDynamics.hh.
◆ getKinematics() [2/2]
Retrieve the Dalitz plot kinematics.
- Returns
- the Dalitz plot kinematics
Definition at line 546 of file LauIsobarDynamics.hh.
◆ getMeanEff()
Retrieve the mean efficiency across the Dalitz plot.
- Returns
- the mean efficiency across the Dalitz plot
Definition at line 489 of file LauIsobarDynamics.hh.
◆ getnCohAmp()
UInt_t LauIsobarDynamics::getnCohAmp |
( |
| ) |
const |
|
inline |
Retrieve the number of coherent amplitude components.
- Returns
- the number of coherent amplitude components
Definition at line 522 of file LauIsobarDynamics.hh.
◆ getnIncohAmp()
UInt_t LauIsobarDynamics::getnIncohAmp |
( |
| ) |
const |
|
inline |
Retrieve the number of incoherent amplitude components.
- Returns
- the number of incoherent amplitude components
Definition at line 528 of file LauIsobarDynamics.hh.
◆ getnTotAmp()
UInt_t LauIsobarDynamics::getnTotAmp |
( |
| ) |
const |
|
inline |
Retrieve the total number of amplitude components.
- Returns
- the total number of amplitude components
Definition at line 516 of file LauIsobarDynamics.hh.
◆ getResonance() [1/2]
Retrieve a resonance by its index.
- Parameters
-
[in] | resIndex | the index of the resonance to retrieve |
- Returns
- the requested resonance
Definition at line 1056 of file LauIsobarDynamics.cc.
◆ getResonance() [2/2]
const LauAbsResonance * LauIsobarDynamics::getResonance |
( |
const UInt_t |
resIndex | ) |
const |
Retrieve a resonance by its index.
- Parameters
-
[in] | resIndex | the index of the resonance to retrieve |
- Returns
- the requested resonance
Definition at line 1043 of file LauIsobarDynamics.cc.
◆ gotKMatrixMatch()
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] | resAmpInt | the index of the resonance within the model |
[in] | propName | the name of the K-matrix propagator |
- Returns
- true if the resonance is a component of the given propagator, otherwise return false
Definition at line 2493 of file LauIsobarDynamics.cc.
◆ gotReweightedEvent()
Bool_t LauIsobarDynamics::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.
Uses the accept/reject method.
- Returns
- kTRUE if the event has been accepted, kFALSE otherwise
Definition at line 2917 of file LauIsobarDynamics.cc.
◆ hasResonance()
Bool_t LauIsobarDynamics::hasResonance |
( |
const TString & |
resName | ) |
const |
Check whether this model includes a named resonance.
- Parameters
-
- Returns
- true if the resonance is present, false otherwise
Definition at line 1033 of file LauIsobarDynamics.cc.
◆ incohResAmp()
Double_t LauIsobarDynamics::incohResAmp |
( |
const UInt_t |
index | ) |
|
|
protected |
Calculate the dynamic part of the intensity for a given incoherent component at the current point in the Dalitz plot.
- Parameters
-
[in] | index | the index of the incoherent component within the model |
Definition at line 2147 of file LauIsobarDynamics.cc.
◆ initialise()
void LauIsobarDynamics::initialise |
( |
const std::vector< LauComplex > & |
coeffs | ) |
|
Initialise the Dalitz plot dynamics.
- Parameters
-
[in] | coeffs | the complex coefficients for the resonances |
Definition at line 361 of file LauIsobarDynamics.cc.
◆ m13IntegrationRegions()
std::vector< LauDPPartialIntegralInfo * > LauIsobarDynamics::m13IntegrationRegions |
( |
const std::vector< std::pair< Double_t, Double_t >> & |
m13Regions, |
|
|
const std::vector< std::pair< Double_t, Double_t >> & |
m23Regions, |
|
|
const std::vector< Double_t > & |
m13Binnings, |
|
|
const Double_t |
precision, |
|
|
const Double_t |
defaultBinning |
|
) |
| const |
|
protected |
Create the integration grid objects for the m13 narrow resonance regions, excluding the overlap regions with the m23 narrow resonances.
The regions will have a m13Binnings x defaultBinning grid. The overlap regions are created by the m23IntegrationRegions function.
- Parameters
-
[in] | m13Regions | the limits of each narrow-resonance region in m13 |
[in] | m23Regions | the limits of each narrow-resonance region in m23 |
[in] | m13Binnings | the binning of each narrow-resonance region in m13 |
[in] | precision | the precision required for the Gauss-Legendre weights |
[in] | defaultBinning | the binning used in the bulk of the phase space |
- Returns
- vector of pointers to LauDPPartialIntegralInfo objects that contain the individual regions
Definition at line 1172 of file LauIsobarDynamics.cc.
◆ m23IntegrationRegions()
std::vector< LauDPPartialIntegralInfo * > LauIsobarDynamics::m23IntegrationRegions |
( |
const std::vector< std::pair< Double_t, Double_t >> & |
m13Regions, |
|
|
const std::vector< std::pair< Double_t, Double_t >> & |
m23Regions, |
|
|
const std::vector< Double_t > & |
m13Binnings, |
|
|
const std::vector< Double_t > & |
m23Binnings, |
|
|
const Double_t |
precision, |
|
|
const Double_t |
defaultBinning |
|
) |
| const |
|
protected |
Create the integration grid objects for the m23 narrow resonance regions, including the overlap regions with the m13 narrow resonances.
The overlap regions will have an m13Binnings x m23Binnings grid. The other regions will have a defaultBinning x m23Binnings grid.
- Parameters
-
[in] | m13Regions | the limits of each narrow-resonance region in m13 |
[in] | m23Regions | the limits of each narrow-resonance region in m23 |
[in] | m13Binnings | the binning of each narrow-resonance region in m13 |
[in] | m23Binnings | the binning of each narrow-resonance region in m23 |
[in] | precision | the precision required for the Gauss-Legendre weights |
[in] | defaultBinning | the binning used in the bulk of the phase space |
- Returns
- vector of pointers to LauDPPartialIntegralInfo objects that contain the individual regions
Definition at line 1247 of file LauIsobarDynamics.cc.
◆ maxGenIterations()
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 276 of file LauIsobarDynamics.hh.
◆ newDPIntegrationRegion()
LauDPPartialIntegralInfo * LauIsobarDynamics::newDPIntegrationRegion |
( |
const Double_t |
minm13, |
|
|
const Double_t |
maxm13, |
|
|
const Double_t |
minm23, |
|
|
const Double_t |
maxm23, |
|
|
const Double_t |
m13BinWidth, |
|
|
const Double_t |
m23BinWidth, |
|
|
const Double_t |
precision, |
|
|
const UInt_t |
nAmp, |
|
|
const UInt_t |
nIncohAmp |
|
) |
| const |
|
protected |
Wrapper for LauDPPartialIntegralInfo constructor.
- Parameters
-
[in] | minm13 | the minimum of the m13 range |
[in] | maxm13 | the maximum of the m13 range |
[in] | minm23 | the minimum of the m23 range |
[in] | maxm23 | the maximum of the m23 range |
[in] | m13BinWidth | the m13 bin width |
[in] | m23BinWidth | the m23 bin width |
[in] | precision | the precision required for the Gauss-Legendre weights |
[in] | nAmp | the number of coherent amplitude components |
[in] | nIncohAmp | the number of incoherent amplitude components |
- Returns
- 0 if the integration region has no internal points, otherwise returns a pointer to the newly constructed LauDPPartialIntegralInfo object
Definition at line 1336 of file LauIsobarDynamics.cc.
◆ removeCharge()
void LauIsobarDynamics::removeCharge |
( |
TString & |
string | ) |
const |
|
protected |
Remove the charge from the given particle name.
- Parameters
-
[in,out] | string | the particle name |
Definition at line 1093 of file LauIsobarDynamics.cc.
◆ resAmp()
LauComplex LauIsobarDynamics::resAmp |
( |
const UInt_t |
index | ) |
|
|
protected |
Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz plot.
- Parameters
-
[in] | index | the index of the amplitude component within the model |
Definition at line 2120 of file LauIsobarDynamics.cc.
◆ resonanceIndex()
Int_t LauIsobarDynamics::resonanceIndex |
( |
const TString & |
resName | ) |
const |
Retrieve the index for the given resonance.
- Parameters
-
- Returns
- the index of the resonance if it is present, -1 otherwise
Definition at line 999 of file LauIsobarDynamics.cc.
◆ retrieveEfficiency()
Double_t LauIsobarDynamics::retrieveEfficiency |
( |
| ) |
|
|
protected |
Obtain the efficiency of the current event from the model.
- Returns
- the efficiency
Definition at line 2990 of file LauIsobarDynamics.cc.
◆ retrieveScfFraction()
Double_t LauIsobarDynamics::retrieveScfFraction |
( |
Int_t |
tagCat | ) |
|
|
protected |
Obtain the self cross feed fraction of the current event from the model.
- Parameters
-
[in] | tagCat | the tagging category of the current event |
- Returns
- the self cross feed fraction
Definition at line 2999 of file LauIsobarDynamics.cc.
◆ setASqMaxValue()
void LauIsobarDynamics::setASqMaxValue |
( |
Double_t |
value | ) |
|
|
inline |
Set the maximum value of A squared to be used in the accept/reject.
Disables the automatic determination of ASqMax
- Parameters
-
Definition at line 240 of file LauIsobarDynamics.hh.
◆ setASqMaxVarValue()
void LauIsobarDynamics::setASqMaxVarValue |
( |
Double_t |
value | ) |
|
|
inlineprotected |
◆ setDataEventNo()
void LauIsobarDynamics::setDataEventNo |
( |
UInt_t |
iEvt | ) |
|
|
protected |
Load the data for a given event.
- Parameters
-
[in] | iEvt | the number of the event |
Definition at line 2704 of file LauIsobarDynamics.cc.
◆ setFFTerm()
void LauIsobarDynamics::setFFTerm |
( |
const UInt_t |
index, |
|
|
const Double_t |
realPart, |
|
|
const Double_t |
imagPart |
|
) |
| |
|
protected |
Set the dynamic part of the amplitude for a given amplitude component at the current point in the Dalitz plot.
- Parameters
-
[in] | index | the index of the amplitude component |
[in] | realPart | the real part of the amplitude |
[in] | imagPart | the imaginary part of the amplitude |
Definition at line 2175 of file LauIsobarDynamics.cc.
◆ setIncohIntenTerm()
void LauIsobarDynamics::setIncohIntenTerm |
( |
const UInt_t |
index, |
|
|
const Double_t |
value |
|
) |
| |
|
protected |
Set the dynamic part of the intensity for a given incoherent amplitude component at the current point in the Dalitz plot.
- Parameters
-
[in] | index | the index of the incoherent amplitude component |
[in] | value | the intensity |
Definition at line 2187 of file LauIsobarDynamics.cc.
◆ setIntegralBinningFactor()
void LauIsobarDynamics::setIntegralBinningFactor |
( |
const Double_t |
binningFactor | ) |
|
|
inline |
Set the factor relating the width of a narrow resonance and the binning size in its integration region.
Narrow resonances trigger different integration behaviour - dividing the DP into regions where a finer binning is used. This can cause high memory usage, so use this method and LauIsobarDynamics::setNarrowResonanceThreshold to tune this behaviour, if needed.
- Parameters
-
[in] | binningFactor | the factor by which the resonance width is divided to obtain the bin size (defaults to 100) |
Definition at line 140 of file LauIsobarDynamics.hh.
◆ setIntegralBinWidths()
void LauIsobarDynamics::setIntegralBinWidths |
( |
const Double_t |
m13BinWidth, |
|
|
const Double_t |
m23BinWidth, |
|
|
const Double_t |
mPrimeBinWidth = 0.001 , |
|
|
const Double_t |
thPrimeBinWidth = 0.001 |
|
) |
| |
Set the widths of the bins to use when integrating across the Dalitz plot or square Dalitz plot.
Specify the bin widths required when performing the DP integration. Note that the integration is not performed in m13^2 vs m23^2 space but in either m13 vs m23 space or mPrime vs thetaPrime space, with the appropriate Jacobian applied. The default bin widths in m13 vs m23 space are 0.005 GeV. The default bin widths in mPrime vs thetaPrime space are 0.001.
- Parameters
-
[in] | m13BinWidth | the bin width to use when integrating over m13 |
[in] | m23BinWidth | the bin width to use when integrating over m23 |
[in] | mPrimeBinWidth | the bin width to use when integrating over mPrime |
[in] | thPrimeBinWidth | the bin width to use when integrating over thetaPrime |
Definition at line 1697 of file LauIsobarDynamics.cc.
◆ setIntFileName()
void LauIsobarDynamics::setIntFileName |
( |
const TString & |
fileName | ) |
|
|
inline |
Set the name of the file to which to save the results of the integrals.
- Parameters
-
[in] | fileName | the name of the file |
Definition at line 102 of file LauIsobarDynamics.hh.
◆ setNarrowResonanceThreshold()
void LauIsobarDynamics::setNarrowResonanceThreshold |
( |
const Double_t |
narrowWidth | ) |
|
|
inline |
Set the value below which a resonance width is considered to be narrow.
Narrow resonances trigger different integration behaviour - dividing the DP into regions where a finer binning is used. This can cause high memory usage, so use this method and LauIsobarDynamics::setIntegralBinningFactor to tune this behaviour, if needed.
- Parameters
-
[in] | narrowWidth | the value below which a resonance is considered to be narrow (defaults to 0.02 GeV/c2) |
Definition at line 131 of file LauIsobarDynamics.hh.
◆ updateCoeffs()
void LauIsobarDynamics::updateCoeffs |
( |
const std::vector< LauComplex > & |
coeffs | ) |
|
Update the complex coefficients for the resonances.
- Parameters
-
[in] | coeffs | the new set of coefficients |
Definition at line 2952 of file LauIsobarDynamics.cc.
◆ usingScfModel()
Bool_t LauIsobarDynamics::usingScfModel |
( |
| ) |
const |
|
inline |
Check whether a self cross feed fraction model is being used.
- Returns
- true if a self cross feed fraction model is being used, false otherwise
Definition at line 564 of file LauIsobarDynamics.hh.
◆ fifjEffSum_
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 980 of file LauIsobarDynamics.hh.
◆ fifjSum_
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 986 of file LauIsobarDynamics.hh.
◆ scfFractionModel_
The self cross feed fraction models across the Dalitz plot.
These model the fraction of signal events that are poorly reconstructed (the self cross feed fraction) as a function of Dalitz plot position. If the self cross feed is depependent on the tagging category then seperate models can be defined.
Definition at line 839 of file LauIsobarDynamics.hh.
The documentation for this class was generated from the following files:
|