laura is hosted by Hepforge, IPPP Durham
Laura++  v3r4
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauIsobarDynamics.hh
Go to the documentation of this file.
1 
2 /*
3 Copyright 2005 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
33 #ifndef LAU_ISOBAR_DYNAMICS
34 #define LAU_ISOBAR_DYNAMICS
35 
36 #include <set>
37 #include <vector>
38 
39 #include "TString.h"
40 
41 #include "LauAbsResonance.hh"
42 #include "LauComplex.hh"
43 
44 class LauCacheData;
45 class LauDaughters;
46 class LauAbsEffModel;
47 class LauAbsIncohRes;
48 class LauFitDataTree;
51 class LauKinematics;
52 
54 
55  public:
57  typedef std::map<Int_t,LauAbsEffModel*> LauTagCatScfFractionModelMap;
58 
60  enum ToyMCStatus {
64  };
65 
67 
72  LauIsobarDynamics(LauDaughters* daughters, LauAbsEffModel* effModel, LauAbsEffModel* scfFractionModel = 0);
73 
75 
80  LauIsobarDynamics(LauDaughters* daughters, LauAbsEffModel* effModel, LauTagCatScfFractionModelMap scfFractionModel);
81 
83  virtual ~LauIsobarDynamics();
84 
86 
89  void initialise(const std::vector<LauComplex>& coeffs);
90 
93 
95 
98  inline void setIntFileName(const TString& fileName) {intFileName_ = fileName;}
99 
100  // Integration
102 
115  void setIntegralBinWidths(const Double_t m13BinWidth, const Double_t m23BinWidth,
116  const Double_t mPrimeBinWidth = 0.001, const Double_t thPrimeBinWidth = 0.001);
117 
119 
125  void setNarrowResonanceThreshold(const Double_t narrowWidth) { narrowWidth_ = narrowWidth; }
126 
128 
134  void setIntegralBinningFactor(const Double_t binningFactor) { binningFactor_ = binningFactor; }
135 
137 
142  void forceSymmetriseIntegration(const Bool_t force) { forceSymmetriseIntegration_ = force; }
143 
145 
156  LauAbsResonance* addResonance(const TString& resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory = LauBlattWeisskopfFactor::Default);
157 
159 
169  LauAbsResonance* addIncoherentResonance(const TString& resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType);
170 
172 
180  void defineKMatrixPropagator(const TString& propName, const TString& paramFileName,
181  Int_t resPairAmpInt, Int_t nChannels, Int_t nPoles, Int_t rowIndex = 1);
182 
184 
194  void addKMatrixProdPole(const TString& poleName, const TString& propName, Int_t poleIndex, Bool_t useProdAdler = kFALSE);
195 
197 
207  void addKMatrixProdSVP(const TString& SVPName, const TString& propName, Int_t channelIndex, Bool_t useProdAdler = kFALSE);
208 
210 
213  inline void setASqMaxValue(Double_t value) {aSqMaxSet_ = value;}
214 
216 
219  inline Double_t getASqMaxSetValue() const { return aSqMaxSet_; }
220 
222 
225  inline Double_t getASqMaxVarValue() const { return aSqMaxVar_; }
226 
228 
231  Bool_t generate();
232 
234 
239  ToyMCStatus checkToyMC(Bool_t printErrorMessages = kTRUE, Bool_t printInfoMessages = kFALSE);
240 
242 
245  inline Int_t maxGenIterations() const {return iterationsMax_;}
246 
248 
251  void calcLikelihoodInfo(const UInt_t iEvt);
252 
254 
258  void calcLikelihoodInfo(const Double_t m13Sq, const Double_t m23Sq);
259 
261 
267  void calcLikelihoodInfo(const Double_t m13Sq, const Double_t m23Sq, const Int_t tagCat);
268 
270 
273  void calcExtraInfo(const Bool_t init = kFALSE);
274 
276 
280  Bool_t gotReweightedEvent();
281 
283 
286  Double_t getEventWeight();
287 
289 
292  inline const LauComplex& getEvtDPAmp() const {return totAmp_;}
293 
295 
298  inline Double_t getEvtm13Sq() const {return m13Sq_;}
299 
301 
304  inline Double_t getEvtm23Sq() const {return m23Sq_;}
305 
307 
310  inline Double_t getEvtmPrime() const {return mPrime_;}
311 
313 
316  inline Double_t getEvtthPrime() const {return thPrime_;}
317 
319 
322  inline Double_t getEvtEff() const {return eff_;}
323 
325 
328  inline Double_t getEvtScfFraction() const {return scfFraction_;}
329 
331 
334  inline Double_t getEvtJacobian() const {return jacobian_;}
335 
337 
340  inline Double_t getEvtIntensity() const {return ASq_;}
341 
343 
349  inline Double_t getEvtLikelihood() const {return evtLike_;}
350 
352 
356  inline LauComplex getDynamicAmp(const Int_t resID) const {return ff_[resID].scale(fNorm_[resID]);}
357 
359 
363  inline LauComplex getFullAmplitude(const Int_t resID) const {return Amp_[resID] * this->getDynamicAmp(resID);}
364 
366 
369  inline const std::vector< std::vector<LauComplex> >& getFiFjSum() const {return fifjSum_;}
370 
372 
375  inline const std::vector< std::vector<LauComplex> >& getFiFjEffSum() const {return fifjEffSum_;}
376 
378 
381  inline const std::vector<Double_t>& getFNorm() const {return fNorm_;}
382 
384 
387  void fillDataTree(const LauFitDataTree& fitDataTree);
388 
390  void modifyDataTree();
391 
393 
397  Bool_t hasResonance(const TString& resName) const;
398 
400 
404  Int_t resonanceIndex(const TString& resName) const;
405 
407 
411  TString getConjResName(const TString& resName) const;
412 
414 
418  const LauAbsResonance* findResonance(const TString& resName) const;
419 
421 
425  const LauAbsResonance* getResonance(const UInt_t resIndex) const;
426 
428 
431  void updateCoeffs(const std::vector<LauComplex>& coeffs);
432 
434 
437  inline void flipHelicityForCPEigenstates(Bool_t boolean) {flipHelicity_ = boolean;}
438 
440 
443  inline const LauParameter& getMeanEff() const {return meanDPEff_;}
444 
446 
449  inline const LauParameter& getDPRate() const {return DPRate_;}
450 
452 
455  inline const LauParArray& getFitFractions() const {return fitFrac_;}
456 
458 
462 
464 
467  inline UInt_t getnTotAmp() const {return nAmp_+nIncohAmp_;}
468 
470 
473  inline UInt_t getnCohAmp() const {return nAmp_;}
474 
476 
479  inline UInt_t getnIncohAmp() const {return nIncohAmp_;}
480 
482 
485  inline Double_t getDPNorm() const {return DPNorm_;}
486 
488 
491  inline const LauDaughters* getDaughters() const {return daughters_;}
492 
494 
497  inline const LauKinematics* getKinematics() const {return kinematics_;}
498 
500 
504 
506 
509  inline const LauAbsEffModel* getEffModel() const {return effModel_;}
510 
512 
515  inline Bool_t usingScfModel() const { return ! scfFractionModel_.empty(); }
516 
518 
521  inline const std::vector<LauParameter>& getExtraParameters() const {return extraParameters_;}
522 
524 
527  inline std::vector<LauParameter*>& getFloatingParameters() {return resonancePars_;}
528 
529  protected:
531  void initSummary();
532 
534  void initialiseVectors();
535 
537  void resetNormVectors();
538 
540  void calcDPNormalisation();
541 
543 
549  std::vector< std::pair<Double_t,Double_t> > formGapsFromRegions(const std::vector< std::pair<Double_t,Double_t> >& regions, const Double_t min, const Double_t max) const;
550 
552 
555  void cullNullRegions(std::vector<LauDPPartialIntegralInfo*>& regions) const;
556 
558 
570  LauDPPartialIntegralInfo* newDPIntegrationRegion(const Double_t minm13, const Double_t maxm13,
571  const Double_t minm23, const Double_t maxm23,
572  const Double_t m13BinWidth, const Double_t m23BinWidth,
573  const Double_t precision,
574  const UInt_t nAmp,
575  const UInt_t nIncohAmp) const;
576 
578 
582  void correctDPOverlap(std::vector< std::pair<Double_t,Double_t> >& regions, const std::vector<Double_t>& binnings) const;
583 
585 
597  std::vector<LauDPPartialIntegralInfo*> m23IntegrationRegions(const std::vector< std::pair<Double_t,Double_t> >& m13Regions,
598  const std::vector< std::pair<Double_t,Double_t> >& m23Regions,
599  const std::vector<Double_t>& m13Binnings,
600  const std::vector<Double_t>& m23Binnings,
601  const Double_t precision,
602  const Double_t defaultBinning) const;
603 
605 
616  std::vector<LauDPPartialIntegralInfo*> m13IntegrationRegions(const std::vector< std::pair<Double_t,Double_t> >& m13Regions,
617  const std::vector< std::pair<Double_t,Double_t> >& m23Regions,
618  const std::vector<Double_t>& m13Binnings,
619  const Double_t precision,
620  const Double_t defaultBinning) const;
621 
624 
627 
629 
633 
635  void writeIntegralsFile();
636 
638 
643  void setFFTerm(const UInt_t index, const Double_t realPart, const Double_t imagPart);
644 
646 
650  void setIncohIntenTerm(const UInt_t index, const Double_t value);
651 
653  void calculateAmplitudes();
654 
656 
661  void calculateAmplitudes( LauDPPartialIntegralInfo* intInfo, const UInt_t m13Point, const UInt_t m23Point );
662 
664 
667  void addGridPointToIntegrals(const Double_t weight);
668 
670 
673  void calcTotalAmp(const Bool_t useEff);
674 
676 
679  Double_t retrieveEfficiency();
680 
682 
686  Double_t retrieveScfFraction(Int_t tagCat);
687 
689 
692  inline void setASqMaxVarValue(Double_t value) {aSqMaxVar_ = value;}
693 
695 
698  Double_t calcSigDPNorm();
699 
701 
704  LauComplex resAmp(const UInt_t index);
705 
707 
710  Double_t incohResAmp(const UInt_t index);
711 
713 
716  void setDataEventNo(UInt_t iEvt);
717 
719 
723  LauAbsResonance* findResonance(const TString& resName);
724 
726 
730  LauAbsResonance* getResonance(const UInt_t resIndex);
731 
733 
736  void removeCharge(TString& string) const;
737 
739 
744  Bool_t gotKMatrixMatch(UInt_t resAmpInt, const TString& propName) const;
745 
746  private:
749 
752 
754  typedef std::map<TString, LauKMatrixPropagator*> KMPropMap;
755 
757  typedef std::map<TString, TString> KMStringMap;
758 
761 
764 
767 
769 
774 
776  UInt_t nAmp_;
777 
779  UInt_t nIncohAmp_;
780 
782  std::vector<LauComplex> Amp_;
783 
785  Double_t DPNorm_;
786 
789 
792 
795 
798 
800  std::vector<LauCacheData*> data_;
801 
804 
806  std::vector<LauParameter> extraParameters_;
807 
809  std::vector<LauAbsResonance*> sigResonances_;
810 
812  std::vector<LauAbsIncohRes*> sigIncohResonances_;
813 
816 
819 
821  std::vector<TString> resTypAmp_;
822 
824  std::vector<Int_t> resPairAmp_;
825 
827  std::vector<TString> incohResTypAmp_;
828 
830  std::vector<Int_t> incohResPairAmp_;
831 
833  std::vector<Int_t> typDaug_;
834 
837 
840 
842  Bool_t flavConjDP_;
843 
846 
849 
852 
854  std::vector<LauDPPartialIntegralInfo*> dpPartialIntegralInfo_;
855 
857  TString intFileName_;
858 
860  Double_t m13BinWidth_;
861 
863  Double_t m23BinWidth_;
864 
866  Double_t mPrimeBinWidth_;
867 
870 
872  Double_t narrowWidth_;
873 
875  Double_t binningFactor_;
876 
878  Double_t m13Sq_;
879 
881  Double_t m23Sq_;
882 
884  Double_t mPrime_;
885 
887  Double_t thPrime_;
888 
890  Int_t tagCat_;
891 
893  Double_t eff_;
894 
896  Double_t scfFraction_;
897 
899  Double_t jacobian_;
900 
902  Double_t ASq_;
903 
905  Double_t evtLike_;
906 
909 
911 
914  std::vector< std::vector<LauComplex> > fifjEffSum_;
915 
917 
920  std::vector< std::vector<LauComplex> > fifjSum_;
921 
923  std::vector<LauComplex> ff_;
924 
926  std::vector<Double_t> incohInten_;
927 
929  std::vector<Double_t> fSqSum_;
930 
932  std::vector<Double_t> fSqEffSum_;
933 
935  std::vector<Double_t> fNorm_;
936 
939 
942 
944  Double_t aSqMaxSet_;
945 
947  Double_t aSqMaxVar_;
948 
951 
954 
956  std::vector<LauParameter*> resonancePars_;
957 
959  std::vector<Double_t> resonanceParValues_;
960 
962  std::vector< std::vector<UInt_t> > resonanceParResIndex_;
963 
965  std::set<UInt_t> integralsToBeCalculated_;
966 
967  ClassDef(LauIsobarDynamics,0)
968 };
969 
970 #endif
KMStringMap kMatrixPropSet_
The names of the M-matrix components in the model mapped to their propagators.
std::vector< LauComplex > Amp_
The complex coefficients for the amplitude components.
void modifyDataTree()
Recache the amplitude values for those that have changed.
Int_t nSigGenLoop_
The number of unsucessful attempts to generate an event so far.
Double_t getASqMaxVarValue() const
Retrieve the maximum of A squared that has been found while generating.
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 regio...
const LauParameter & getDPRate() const
Retrieve the overall Dalitz plot rate.
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 regio...
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...
Bool_t recalcNormalisation_
Flag to recalculate the normalisation.
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.
virtual ~LauIsobarDynamics()
Destructor.
const std::vector< Double_t > & getFNorm() const
Retrieve the normalisation factors for the dynamic parts of the amplitudes for all of the amplitude c...
void updateCoeffs(const std::vector< LauComplex > &coeffs)
Update the complex coefficients for the resonances.
std::set< UInt_t > integralsToBeCalculated_
Resonance indices for which the amplitudes and integrals should be recalculated.
Double_t m13BinWidth_
The bin width to use when integrating over m13.
Bool_t flipHelicity_
The helicity flip flag for new amplitude components.
Double_t incohResAmp(const UInt_t index)
Calculate the dynamic part of the intensity for a given incoherent component at the current point in ...
std::vector< std::pair< Double_t, Double_t > > formGapsFromRegions(const std::vector< std::pair< Double_t, Double_t > > &regions, const Double_t min, const Double_t max) const
Form the regions that are produced by the spaces between narrow resonances.
std::vector< Int_t > resPairAmp_
The index of the daughter not produced by the resonance for each amplitude component.
Bool_t hasResonance(const TString &resName) const
Check whether this model includes a named resonance.
const LauParameter & getMeanEff() const
Retrieve the mean efficiency across the Dalitz plot.
const LauAbsResonance * findResonance(const TString &resName) const
Retrieve the named resonance.
std::vector< LauParameter > extraParameters_
any extra parameters/quantities (e.g. K-matrix total fit fractions)
const LauKinematics * getKinematics() const
Retrieve the Dalitz plot kinematics.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getEvtJacobian() const
Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event...
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.
Double_t getDPNorm() const
Retrieve the normalisation factor for the log-likelihood function.
std::vector< Double_t > fSqEffSum_
The event-by-event running total of the dynamical amplitude squared for each amplitude component...
void calcDPNormalisationScheme()
Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
LauAbsEffModel * effModel_
The efficiency model across the Dalitz plot.
void calcExtraInfo(const Bool_t init=kFALSE)
Calculate the fit fractions, mean efficiency and total DP rate.
std::vector< Int_t > incohResPairAmp_
The index of the daughter not produced by the resonance for each incoherent amplitude component...
void setASqMaxValue(Double_t value)
Set the maximum value of A squared to be used in the accept/reject.
Double_t DPNorm_
The normalisation factor for the log-likelihood function.
LauTagCatScfFractionModelMap scfFractionModel_
The self cross feed fraction models across the Dalitz plot.
LauParArray fitFrac_
The fit fractions for the amplitude components.
Double_t aSqMaxVar_
The maximum value of A squared that has been seen so far while generating.
void forceSymmetriseIntegration(const Bool_t force)
Force the symmetrisation of the integration in m13 &lt;-&gt; m23 for non-symmetric but flavour-conjugate fi...
Double_t scfFraction_
The fraction of events that are poorly reconstructed (the self cross feed fraction) at the current po...
Abstract class for defining incoherent resonant amplitude models.
Double_t m13Sq_
The invariant mass squared of the first and third daughters.
Int_t iterationsMax_
The maximum allowed number of attempts when generating an event.
LauComplex totAmp_
The total amplitude for the current event.
Double_t getEvtthPrime() const
Retrieve the square Dalitz plot coordinate, theta&#39;, for the current event.
Double_t getASqMaxSetValue() const
Retrieve the maximum value of A squared to be used in the accept/reject.
TString getConjResName(const TString &resName) const
Retrieve the name of the charge conjugate of a named resonance.
std::vector< std::vector< LauParameter > > LauParArray
Type to define an array of parameters.
Double_t aSqMaxSet_
The maximum allowed value of A squared.
std::vector< Double_t > resonanceParValues_
List of floating resonance parameter values from previous calculation.
Bool_t symmetricalDP_
Whether the Dalitz plot is symmetrical.
UInt_t getnIncohAmp() const
Retrieve the number of incoherent amplitude components.
Pure abstract base class for defining the efficiency description across the signal Dalitz plot...
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...
Class for defining (a section of) the Dalitz plot integration binning scheme.
Class to contain cached data relating to an event.
Definition: LauCacheData.hh:44
void calcDPNormalisation()
Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
LauComplex getFullAmplitude(const Int_t resID) const
Retrieve the Amplitude of resonance resID.
const LauParArray & getFitFractions() const
Retrieve the fit fractions for the amplitude components.
void recalculateNormalisation()
recalculate Normalization
std::vector< TString > incohResTypAmp_
The resonance types of all of the incoherent amplitude components.
std::vector< LauAbsIncohRes * > sigIncohResonances_
The incoherent resonances in the model.
void setASqMaxVarValue(Double_t value)
Set the maximum of A squared that has been found.
Double_t getEvtScfFraction() const
Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the ...
LauParameter DPRate_
The overall Dalitz plot rate.
void calcLikelihoodInfo(const UInt_t iEvt)
Calculate the likelihood (and all associated information) for the given event number.
ToyMCStatus
The possible statuses for toy MC generation.
KMPropMap kMatrixPropagators_
The K-matrix propagators.
std::vector< Double_t > fSqSum_
The event-by-event running total of the dynamical amplitude squared for each amplitude component...
LauComplex resAmp(const UInt_t index)
Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz ...
Bool_t usingScfModel() const
Check whether a self cross feed fraction model is being used.
const LauAbsResonance * getResonance(const UInt_t resIndex) const
Retrieve a resonance by its index.
Double_t getEvtmPrime() const
Retrieve the square Dalitz plot coordinate, m&#39;, for the current event.
UInt_t nAmp_
The number of amplitude components.
Double_t m23Sq_
The invariant mass squared of the second and third daughters.
Double_t thPrime_
The square Dalitz plot coordinate theta&#39;.
Double_t evtLike_
The normalised likelihood for the current event.
void setIntegralBinningFactor(const Double_t binningFactor)
Set the factor relating the width of a narrow resonance and the binning size in its integration regio...
std::vector< std::vector< LauComplex > > fifjEffSum_
The event-by-event running total of efficiency corrected amplitude cross terms for each pair of ampli...
void initSummary()
Print a summary of the model to be used.
Bool_t fullySymmetricDP_
Whether the Dalitz plot is fully symmetric.
void fillDataTree(const LauFitDataTree &fitDataTree)
Fill the internal data structure that caches the resonance dynamics.
void resetNormVectors()
Zero the various values used to store integrals info.
Double_t ASq_
The value of A squared for the current event.
Bool_t flavConjDP_
Whether the Dalitz plot is a flavour-conjugate final state.
const std::vector< std::vector< LauComplex > > & getFiFjSum() const
Retrieve the event-by-event running totals of amplitude cross terms for all pairs of amplitude compon...
Double_t getEvtm23Sq() const
Retrieve the invariant mass squared of the second and third daughters in the current event...
Double_t m23BinWidth_
The bin width to use when integrating over m23.
TString intFileName_
The name of the file to save integrals to.
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.
std::vector< LauParameter * > & getFloatingParameters()
Retrieve the floating parameters of the resonance models.
void setNarrowResonanceThreshold(const Double_t narrowWidth)
Set the value below which a resonance width is considered to be narrow.
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.
void calcTotalAmp(const Bool_t useEff)
Calculate the total Dalitz plot amplitude at the current point in the Dalitz plot.
const LauComplex & getEvtDPAmp() const
Retrieve the total amplitude for the current event.
Double_t retrieveScfFraction(Int_t tagCat)
Obtain the self cross feed fraction of the current event from the model.
std::vector< TString > resTypAmp_
The resonance types of all of the amplitude components.
LauKinematics * getKinematics()
Retrieve the Dalitz plot kinematics.
LauIsobarDynamics(LauDaughters *daughters, LauAbsEffModel *effModel, LauAbsEffModel *scfFractionModel=0)
Constructor.
Int_t tagCat_
The tagging category.
LauCacheData * currentEvent_
The cached data for the current event.
void initialiseVectors()
Initialise the internal storage for this model.
Double_t jacobian_
The Jacobian, for the transformation into square DP coordinates at the current point in the Dalitz pl...
LauDaughters * daughters_
The daughters of the decay.
std::map< Int_t, LauAbsEffModel * > LauTagCatScfFractionModelMap
The type used for containing multiple self cross feed fraction models for different categories (e...
const LauDaughters * getDaughters() const
Retrieve the daughters.
Double_t narrowWidth_
The value below which a resonance width is considered to be narrow.
File containing declaration of LauComplex class.
const LauParArray & getFitFractionsEfficiencyUncorrected() const
Retrieve the fit fractions for the amplitude components.
std::vector< std::vector< LauComplex > > fifjSum_
The event-by-event running total of the amplitude cross terms for each pair of amplitude components...
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
UInt_t getnTotAmp() const
Retrieve the total number of amplitude components.
Double_t calcSigDPNorm()
Calculate the normalisation factor for the log-likelihood function.
const LauAbsEffModel * getEffModel() const
Retrieve the model for the efficiency across the Dalitz plot.
void addGridPointToIntegrals(const Double_t weight)
Add the amplitude values (with the appropriate weight) at the current grid point to the running integ...
Double_t getEvtIntensity() const
Retrieve the total intensity multiplied by the efficiency for the current event.
Double_t thPrimeBinWidth_
The bin width to use when integrating over thetaPrime.
Bool_t generate()
Generate a toy MC signal event.
void removeCharge(TString &string) const
Remove the charge from the given particle name.
std::vector< Int_t > typDaug_
The PDG codes of the daughters.
LauResonanceModel
Define the allowed resonance types.
void flipHelicityForCPEigenstates(Bool_t boolean)
Set the helicity flip flag for new amplitude components.
Double_t eff_
The efficiency at the current point in the Dalitz plot.
void calculateAmplitudes()
Calculate the amplitudes for all resonances for the current kinematics.
void correctDPOverlap(std::vector< std::pair< Double_t, Double_t > > &regions, const std::vector< Double_t > &binnings) const
Correct regions to ensure that the finest integration grid takes precedence.
Int_t maxGenIterations() const
Retrieve the maximum number of iterations allowed when generating an event.
Int_t resonanceIndex(const TString &resName) const
Retrieve the index for the given resonance.
Bool_t gotReweightedEvent()
Calculates whether an event with the current kinematics should be accepted in order to produce a dist...
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
Double_t binningFactor_
The factor relating the width of the narrowest resonance and the binning size.
std::vector< LauCacheData * > data_
The cached data for all events.
void findIntegralsToBeRecalculated()
Determine which amplitudes and integrals need to be recalculated.
Class for defining signal dynamics using the isobar model.
Double_t getEventWeight()
Calculate the acceptance rate, for events with the current kinematics, when generating events accordi...
void setIntFileName(const TString &fileName)
Set the name of the file to which to save the results of the integrals.
LauAbsResonance * addIncoherentResonance(const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType)
Add an incoherent resonance to the Dalitz plot.
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.
void writeIntegralsFile()
Write the results of the integrals (and related information) to a file.
File containing declaration of LauAbsResonance class.
void calcDPPartialIntegral(LauDPPartialIntegralInfo *intInfo)
Calculate the Dalitz plot normalisation integrals over a given range.
Double_t mPrimeBinWidth_
The bin width to use when integrating over mPrime.
ToyMCStatus checkToyMC(Bool_t printErrorMessages=kTRUE, Bool_t printInfoMessages=kFALSE)
Check the status of the toy MC generation.
LauKinematics * kinematics_
The kinematics of the decay.
UInt_t getnCohAmp() const
Retrieve the number of coherent amplitude components.
LauIsobarDynamics & operator=(const LauIsobarDynamics &rhs)
Copy assignment operator (not implemented)
Class for defining a complex number.
Definition: LauComplex.hh:61
Double_t getEvtEff() const
Retrieve the efficiency for the current event.
UInt_t nIncohAmp_
The number of incoherent amplitude components.
Bool_t gotKMatrixMatch(UInt_t resAmpInt, const TString &propName) const
Check whether a resonance is a K-matrix component of a given propagator.
Double_t getEvtLikelihood() const
Retrieve the likelihood for the current event.
std::vector< LauDPPartialIntegralInfo * > dpPartialIntegralInfo_
The storage of the integration scheme.
std::map< TString, TString > KMStringMap
The type used for mapping K-matrix components to their propagators.
void cullNullRegions(std::vector< LauDPPartialIntegralInfo * > &regions) const
Removes entries in the vector of LauDPPartialIntegralInfo* that are null.
void setDataEventNo(UInt_t iEvt)
Load the data for a given event.
std::map< TString, LauKMatrixPropagator * > KMPropMap
The type used for containing the K-matrix propagators.
Class for calculating 3-body kinematic quantities.
Bool_t normalizationSchemeDone_
Whether the scheme for the integration has been determined.
std::vector< LauParameter * > resonancePars_
List of floating resonance parameters.
std::vector< std::vector< UInt_t > > resonanceParResIndex_
Indices in sigResonances_ to point to the corresponding signal resonance(s) for each floating paramet...
BlattWeisskopfCategory
Define resonance categories that will share common barrier factor radii.
LauParameter meanDPEff_
The mean efficiency across the Dalitz plot.
Bool_t integralsDone_
Whether the integrals have been performed.
std::vector< Double_t > incohInten_
The dynamic part of the intensity for each incoherent amplitude component at the current point in the...
std::vector< LauComplex > ff_
The dynamic part of the amplitude for each amplitude component at the current point in the Dalitz plo...
const std::vector< LauParameter > & getExtraParameters() const
Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions)
Bool_t forceSymmetriseIntegration_
Force the symmetrisation of the integration in m13 &lt;-&gt; m23 for non-symmetric but flavour-conjugate fi...
void initialise(const std::vector< LauComplex > &coeffs)
Initialise the Dalitz plot dynamics.
const std::vector< std::vector< LauComplex > > & getFiFjEffSum() const
Retrieve the event-by-event running totals of efficiency corrected amplitude cross terms for all pair...
Double_t getEvtm13Sq() const
Retrieve the invariant mass squared of the first and third daughters in the current event...
Double_t mPrime_
The square Dalitz plot coordinate, m&#39;.
std::vector< Double_t > fNorm_
The normalisation factors for the dynamic parts of the amplitude for each amplitude component...
LauComplex getDynamicAmp(const Int_t resID) const
Retrieve the normalised dynamic part of the amplitude of the given amplitude component at the current...
Class to store the input fit variables.
std::vector< LauAbsResonance * > sigResonances_
The resonances in the model.
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 Dal...
Double_t retrieveEfficiency()
Obtain the efficiency of the current event from the model.
LauParArray fitFracEffUnCorr_
The efficiency-uncorrected fit fractions for the amplitude components.
Class for defining a K-matrix propagator.