laura is hosted by Hepforge, IPPP Durham
Laura++  v2r2p1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauIsobarDynamics.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2005 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 #include <iomanip>
17 #include <fstream>
18 
19 #include "TFile.h"
20 #include "TRandom.h"
21 #include "TSystem.h"
22 
23 #include "LauAbsResonance.hh"
24 #include "LauBelleNR.hh"
25 #include "LauBelleSymNR.hh"
26 #include "LauCacheData.hh"
27 #include "LauConstants.hh"
28 #include "LauDaughters.hh"
29 #include "LauAbsEffModel.hh"
30 #include "LauFitDataTree.hh"
31 #include "LauGounarisSakuraiRes.hh"
32 #include "LauIntegrals.hh"
33 #include "LauIsobarDynamics.hh"
34 #include "LauKinematics.hh"
35 #include "LauKMatrixPropagator.hh"
36 #include "LauKMatrixPropFactory.hh"
37 #include "LauKMatrixProdPole.hh"
38 #include "LauKMatrixProdSVP.hh"
39 #include "LauNRAmplitude.hh"
40 #include "LauPolNR.hh"
41 #include "LauPrint.hh"
42 #include "LauRandom.hh"
43 #include "LauRelBreitWignerRes.hh"
44 #include "LauResonanceMaker.hh"
45 
47 
48 
49 // for Kpipi: only one scfFraction 2D histogram is needed
51  LauAbsDPDynamics(daughters, effModel, scfFractionModel),
52  symmetricalDP_(kFALSE),
53  integralsDone_(kFALSE),
54  intFileName_("integ.dat"),
55  m13BinWidth_(0.005),
56  m23BinWidth_(0.005),
57  m13Sq_(0.0),
58  m23Sq_(0.0),
59  mPrime_(0.0),
60  thPrime_(0.0),
61  eff_(1.0),
62  scfFraction_(0.0),
63  jacobian_(0.0),
64  ASq_(0.0),
65  evtLike_(0.0),
66  iterationsMax_(100000),
67  nSigGenLoop_(0),
68  aSqMaxSet_(1.25),
69  aSqMaxVar_(0.0),
70  setBarrierRadius_(kFALSE),
71  resBarrierRadius_(4.0),
72  parBarrierRadius_(4.0),
73  barrierType_( LauAbsResonance::BWPrimeBarrier ),
74  flipHelicity_(kTRUE)
75 {
76  // Constructor for the isobar signal model
77  if (daughters != 0) {
78  symmetricalDP_ = daughters->gotSymmetricalDP();
79  typDaug_.push_back(daughters->getTypeDaug1());
80  typDaug_.push_back(daughters->getTypeDaug2());
81  typDaug_.push_back(daughters->getTypeDaug3());
82  }
83  sigResonances_.clear();
84  kMatrixPropagators_.clear();
85  kMatrixPropSet_.clear();
86 }
87 
88 // for Kspipi, we need a scfFraction 2D histogram for each tagging category. They are provided by the map.
89 // Also, we need to know the place that the tagging category of the current event occupies in the data structure inputFitTree
91  LauAbsDPDynamics(daughters, effModel, scfFractionModel),
92  symmetricalDP_(kFALSE),
93  integralsDone_(kFALSE),
94  intFileName_("integ.dat"),
95  m13BinWidth_(0.005),
96  m23BinWidth_(0.005),
97  m13Sq_(0.0),
98  m23Sq_(0.0),
99  mPrime_(0.0),
100  thPrime_(0.0),
101  eff_(0.0),
102  scfFraction_(0.0),
103  jacobian_(0.0),
104  ASq_(0.0),
105  evtLike_(0.0),
106  iterationsMax_(100000),
107  nSigGenLoop_(0),
108  aSqMaxSet_(1.25),
109  aSqMaxVar_(0.0),
110  setBarrierRadius_(kFALSE),
111  resBarrierRadius_(4.0),
112  parBarrierRadius_(4.0),
113  barrierType_( LauAbsResonance::BWPrimeBarrier ),
114  flipHelicity_(kTRUE)
115 {
116  // Constructor for the isobar signal model
117  if (daughters != 0) {
118  symmetricalDP_ = daughters->gotSymmetricalDP();
119  typDaug_.push_back(daughters->getTypeDaug1());
120  typDaug_.push_back(daughters->getTypeDaug2());
121  typDaug_.push_back(daughters->getTypeDaug3());
122  }
123  sigResonances_.clear();
124  kMatrixPropagators_.clear();
125  kMatrixPropSet_.clear();
126 }
127 
129 {
130 }
131 
132 void LauIsobarDynamics::initialise(const std::vector<LauComplex>& coeffs)
133 {
134  // Check whether we have a valid set of integration constants for
135  // the normalisation of the signal likelihood function.
136  this->initialiseVectors();
137 
138  // Mark the DP integrals as undetermined
139  integralsDone_ = kFALSE;
140 
141  // Initialise all resonance models
142  for ( std::vector<LauAbsResonance*>::iterator iter = sigResonances_.begin(); iter != sigResonances_.end(); ++iter ) {
143  (*iter)->initialise();
144  }
145 
146  // Print summary of what we have so far to screen
147  this->initSummary();
148 
149  if ( nAmp_ == 0 ) {
150  std::cout << "INFO in LauIsobarDynamics::initialise : No contributions to DP model, not performing normalisation integrals." << std::endl;
151  } else {
152 
153  // We need to calculate the normalisation constants for the
154  // Dalitz plot generation/fitting.
155  std::cout<<"INFO in LauIsobarDynamics::initialise : Starting special run to generate the integrals for normalising the PDF..."<<std::endl;
156  this->calcDPNormalisation();
157 
158  // Write the integrals to a file (mainly for debugging purposes)
159  this->writeIntegralsFile();
160  }
161 
162 
163  integralsDone_ = kTRUE;
164 
165  std::cout << std::setprecision(10);
166 
167  // Calculate and cache the relative normalisations of each resonance _dynamic_ amplitude
168  // (e.g. Breit-Wigner contribution, not from the complex amplitude/phase) w.r.t. the
169  // total DP amplitude. These are stored in fNorm_[i].
170  // The normalisation uses fSqSum[i], which is calculated within the dynamics() function,
171  // which has already been called by the calcDPNomalisation() function above to create the
172  // normalisation integrals.
173  // fSqSum[i] is the event-by-event running total of the dynamical amplitude
174  // squared for a given resonance, i. We require that:
175  // |fNorm_[i]|^2 * |fSqSum[i]|^2 = 1,
176  // i.e. fNorm_[i] normalises each resonance contribution to give the same number of
177  // events in the DP, accounting for the total DP area and the
178  // dynamics of the resonance.
179 
180  std::cout<<"INFO in LauIsobarDynamics::initialise : Summary of the integrals:"<<std::endl;
181  for (UInt_t i = 0; i < nAmp_; i++) {
182  fNorm_[i] = 0.0;
183  if (fSqSum_[i] > 0.0) {fNorm_[i] = TMath::Sqrt(1.0/(fSqSum_[i]));}
184  std::cout<<" fNorm["<<i<<"] = "<<fNorm_[i]<<std::endl;
185  std::cout<<" fSqSum["<<i<<"] = "<<fSqSum_[i]<<std::endl;
186  }
187 
188  for (UInt_t i = 0; i < nAmp_; i++) {
189  for (UInt_t j = 0; j < nAmp_; j++) {
190  std::cout<<" fifjEffSum["<<i<<"]["<<j<<"] = "<<fifjEffSum_[i][j];
191  }
192  std::cout<<std::endl;
193  }
194 
195  for (UInt_t i = 0; i < nAmp_; i++) {
196  for (UInt_t j = 0; j < nAmp_; j++) {
197  std::cout<<" fifjSum["<<i<<"]["<<j<<"] = "<<fifjSum_[i][j];
198  }
199  std::cout<<std::endl;
200  }
201 
202  // Calculate the initial fit fractions (for later comparison with Toy MC, if required)
203  this->updateCoeffs(coeffs);
204  this->calcExtraInfo(kTRUE);
205  for (UInt_t i = 0; i < nAmp_; i++) {
206  for (UInt_t j = i; j < nAmp_; j++) {
207  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial fit fraction for amplitude ("<<i<<","<<j<<") = "<<fitFrac_[i][j].genValue()<<std::endl;
208  }
209  }
210 
211  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial efficiency = "<<meanDPEff_.initValue()<<std::endl;
212 
213  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial DPRate = "<<DPRate_.initValue()<<std::endl;
214 }
215 
217 {
218 
219  UInt_t i(0);
220  TString nameP = daughters_->getNameParent();
221  TString name1 = daughters_->getNameDaug1();
222  TString name2 = daughters_->getNameDaug2();
223  TString name3 = daughters_->getNameDaug3();
224 
225  std::cout<<"INFO in LauIsobarDynamics::initSummary : We are going to do a DP with "<<nameP<<" going to "<<name1<<" "<<name2<<" "<<name3<<std::endl;
226 
227  std::cout<<" : For the following resonance combinations:"<<std::endl;
228 
229  std::cout<<" : In tracks 2 and 3:"<<std::endl;
230  for (i = 0; i < nAmp_; i++) {
231  if (resPairAmp_[i] == 1) {
232  std::cout<<" : "<<resTypAmp_[i]<<" to "<<name2<<", "<< name3<<std::endl;
233  }
234  }
235 
236  std::cout<<std::endl;
237 
238  std::cout<<" : In tracks 1 and 3:"<<std::endl;
239  for (i = 0; i < nAmp_; i++) {
240  if (resPairAmp_[i] == 2) {
241  std::cout<<" : "<<resTypAmp_[i]<<" to "<<name1<<", "<< name3<<std::endl;
242  }
243  }
244 
245  std::cout<<std::endl;
246 
247  std::cout<<" : In tracks 1 and 2:"<<std::endl;
248  for (i = 0; i < nAmp_; i++) {
249  if (resPairAmp_[i] == 3) {
250  std::cout<<" : "<<resTypAmp_[i]<<" to "<<name1<<", "<< name2<<std::endl;
251  }
252  }
253 
254  std::cout<<std::endl;
255 
256  for (i = 0; i < nAmp_; i++) {
257  if (resPairAmp_[i] == 0) {
258  std::cout<<" : and a non-resonant amplitude."<<std::endl;
259  }
260  }
261 
262  std::cout<<std::endl;
263 }
264 
266 {
267  fSqSum_.clear(); fSqSum_.resize(nAmp_);
268  fifjEffSum_.clear(); fifjEffSum_.resize(nAmp_);
269  fifjSum_.clear(); fifjSum_.resize(nAmp_);
270  ff_.clear(); ff_.resize(nAmp_);
271  fNorm_.clear(); fNorm_.resize(nAmp_);
272  fitFrac_.clear(); fitFrac_.resize(nAmp_);
273  fitFracEffUnCorr_.clear(); fitFracEffUnCorr_.resize(nAmp_);
274 
275  LauComplex null(0.0, 0.0);
276 
277  for (UInt_t i = 0; i < nAmp_; i++) {
278 
279  fSqSum_[i] = 0.0; fNorm_[i] = 0.0;
280  ff_[i] = null;
281  fifjEffSum_[i].resize(nAmp_);
282  fifjSum_[i].resize(nAmp_);
283  fitFrac_[i].resize(nAmp_);
284  fitFracEffUnCorr_[i].resize(nAmp_);
285 
286  for (UInt_t j = 0; j < nAmp_; j++) {
287  fifjEffSum_[i][j] = null;
288  fifjSum_[i][j] = null;
289  fitFrac_[i][j].valueAndRange(0.0, -100.0, 100.0);
290  fitFracEffUnCorr_[i][j].valueAndRange(0.0, -100.0, 100.0);
291  }
292  }
293 
294  UInt_t nKMatrixProps = kMatrixPropagators_.size();
295  extraParameters_.clear();
296  extraParameters_.resize( nKMatrixProps );
297 
298  for ( UInt_t i = 0; i < nKMatrixProps; ++i ) {
299  extraParameters_[i].valueAndRange(0.0, -100.0, 100.0);
300  }
301 }
302 
303 
305 {
306  // Routine to end integration calculation for loglike normalisation.
307  // This writes out the normalisation integral output into the file named
308  // outputFileName.
309  std::cout<<"INFO in LauIsobarDynamics::writeIntegralsFile : Writing integral output to integrals file "<<intFileName_.Data()<<std::endl;
310 
311  UInt_t i(0), j(0);
312  std::ofstream getChar(intFileName_.Data());
313 
314  getChar << std::setprecision(10);
315 
316  // Write out daughter types (pi, pi0, K, K0s?)
317  for (i = 0; i < 3; i++) {
318  getChar << typDaug_[i] << " ";
319  }
320 
321  // Write out number of resonances in the Dalitz plot model
322  getChar << nAmp_ << std::endl;
323 
324  // Write out the resonances
325  for (i = 0; i < nAmp_; i++) {
326  getChar << resTypAmp_[i] << " ";
327  }
328 
329  getChar << std::endl;
330 
331  // Write out the resonance model types (BW, RelBW etc...)
332  for (i = 0; i < nAmp_; i++) {
333  LauAbsResonance* theResonance = sigResonances_[i];
334  Int_t resModelInt = theResonance->getResonanceModel();
335  getChar << resModelInt << " ";
336  }
337 
338  getChar << std::endl;
339 
340  // Write out the track pairings for each resonance. This is specified
341  // by the resPairAmpInt integer in the addResonance function.
342  for (i = 0; i < nAmp_; i++) {
343  getChar << resPairAmp_[i] << " ";
344  }
345 
346  getChar << std::endl;
347 
348  // Write out the fSqSum = |ff|^2, where ff = resAmp()
349  for (i = 0; i < nAmp_; i++) {
350  getChar << fSqSum_[i] << " ";
351  }
352 
353  getChar << std::endl;
354 
355  // Write out the f_i*f_j_conj*eff values = resAmp_i*resAmp_j_conj*eff.
356  // Note that only the top half of the i*j "matrix" is required, as it
357  // is symmetric w.r.t i, j.
358  for (i = 0; i < nAmp_; i++) {
359  for (j = i; j < nAmp_; j++) {
360  getChar << fifjEffSum_[i][j] << " ";
361  }
362  }
363 
364  getChar << std::endl;
365 
366  // Similar to fifjEffSum, but without the efficiency term included.
367  for (i = 0; i < nAmp_; i++) {
368  for (j = i; j < nAmp_; j++) {
369  getChar << fifjSum_[i][j] << " ";
370  }
371  }
372 
373  getChar << std::endl;
374 
375 }
376 
377 LauAbsResonance* LauIsobarDynamics::addResonance(const TString& resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType)
378 {
379  // Function to add a resonance in a Dalitz plot.
380  // No check is made w.r.t flavour and charge conservation rules, and so
381  // the user is responsible for checking the internal consistency of
382  // their function statements with these laws. For example, the program
383  // will not prevent the user from asking for a rho resonance in a K-pi
384  // pair or a K* resonance in a pi-pi pair.
385  // However, to assist the user, a summary of the resonant structure requested
386  // by the user is printed before the program runs. It is important to check this
387  // information when you first define your Dalitz plot model before doing
388  // any fitting/generating.
389  // Arguments are: resonance name, integer to specify the resonance track pairing
390  // (1 => m_23, 2 => m_13, 3 => m_12), i.e. the bachelor track number.
391  // The third argument resType specifies whether the resonance is a Breit-Wigner (BW)
392  // Relativistic Breit-Wigner (RelBW) or Flatte distribution (Flatte), for example.
393 
394  LauAbsResonance *theResonance =
395  resonanceMaker_->getResonance(resName, resPairAmpInt, resType);
396 
397  if (theResonance == 0) {
398  std::cerr<<"ERROR in LauIsobarDynamics::addResonance : Couldn't create the resonance \""<<resName<<"\""<<std::endl;
399  return 0;
400  }
401 
402  // implement the helicity flip here
403  if (flipHelicity_ && daughters_->getCharge(resPairAmpInt) == 0) {
404  if ( daughters_->getChargeParent() == 0 && daughters_->getTypeParent() > 0 ) {
405  theResonance->flipHelicity(kTRUE);
406  }
407  }
408 
409  // Set the Blatt-Weisskopf barrier factors as appropriate
410  if (setBarrierRadius_) {
412  }
413 
414  // Set the resonance name and what track is the bachelor
415  TString resonanceName = theResonance->getResonanceName();
416  resTypAmp_.push_back(resonanceName);
417  resIntAmp_.push_back(resonanceMaker_->resTypeInt(resonanceName));
418 
419  // Always force the non-resonant amplitude pair to have resPairAmp = 0
420  // in case the user chooses the wrong number.
421  if ( resType == LauAbsResonance::FlatNR ||
422  resType == LauAbsResonance::BelleSymNR ||
423  resType == LauAbsResonance::NRModel ) {
424  std::cout<<"INFO in LauIsobarDynamics::addResonance : Setting resPairAmp to 0 for "<<resonanceName<<" contribution."<<std::endl;
425  resPairAmp_.push_back(0);
426  } else {
427  resPairAmp_.push_back(resPairAmpInt);
428  }
429 
430  // Increment the number of resonance amplitudes we have so far
431  ++nAmp_;
432 
433  // Finally, add the resonance object to the internal array
434  sigResonances_.push_back(theResonance);
435 
436  std::cout<<"INFO in LauIsobarDynamics::addResonance : Successfully added resonance. Total number of resonances so far = "<<nAmp_<<std::endl;
437 
438  return theResonance;
439 }
440 
441 void LauIsobarDynamics::defineKMatrixPropagator(const TString& propName, const TString& paramFileName, Int_t resPairAmpInt,
442  Int_t nChannels, Int_t nPoles, Int_t rowIndex)
443 {
444  // Define the K-matrix propagator. The resPairAmpInt integer specifies which mass combination should be used
445  // for the invariant mass-squared variable "s". The pole masses and coupling constants are defined in the
446  // paramFileName parameter file. The number of channels and poles are defined by the nChannels and nPoles integers, respectively.
447  // The integer rowIndex specifies which row of the propagator should be used when
448  // summing over all amplitude channels: S-wave will be the first row, so rowIndex = 1.
449 
450  if (rowIndex < 1) {
451  std::cerr<<"Error in defineKMatrixPropagator: rowIndex must be > 0 but is equal to "<<rowIndex<<std::endl;
452  return;
453  }
454 
455  TString propagatorName(propName), parameterFile(paramFileName);
456 
457  LauKMatrixPropagator* thePropagator = LauKMatrixPropFactory::getInstance()->getPropagator(propagatorName, parameterFile,
458  resPairAmpInt, nChannels,
459  nPoles, rowIndex);
460  kMatrixPropagators_[propagatorName] = thePropagator;
461 
462 }
463 
464 void LauIsobarDynamics::addKMatrixProdPole(const TString& poleName, const TString& propName, Int_t poleIndex)
465 {
466 
467  // Add a K-matrix production pole term, using the K-matrix propagator given by the propName.
468  // Here, poleIndex is the integer specifying the pole number.
469 
470  // First, find the K-matrix propagator.
471  KMPropMap::iterator mapIter = kMatrixPropagators_.find(propName);
472  if (mapIter != kMatrixPropagators_.end()) {
473 
474  LauKMatrixPropagator* thePropagator = mapIter->second;
475 
476  // Make sure the pole index is valid
477  Int_t nPoles = thePropagator->getNPoles();
478  if (poleIndex < 1 || poleIndex > nPoles) {
479  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdPole : The pole index "<<poleIndex
480  <<" is not between 1 and "<<nPoles<<". Not adding production pole "<<poleName
481  <<" for K-matrix propagator "<<propName<<std::endl;
482  return;
483  }
484 
485  // Now add the K-matrix production pole amplitude to the vector of LauAbsResonance pointers.
486  Int_t resPairAmpInt = thePropagator->getResPairAmpInt();
487  LauAbsResonance* prodPole = new LauKMatrixProdPole(poleName, poleIndex, resPairAmpInt, thePropagator, daughters_);
488 
489  resTypAmp_.push_back(poleName);
490  resIntAmp_.push_back(0);
491  resPairAmp_.push_back(resPairAmpInt);
492 
493  nAmp_++;
494  sigResonances_.push_back(prodPole);
495 
496  // Also store the propName-poleName pair for calculating total fit fractions later on
497  // (avoiding the need to use dynamic casts to check which resonances are of the K-matrix type)
498  kMatrixPropSet_[poleName] = propName;
499 
500  std::cout<<"INFO in LauIsobarDynamics::addKMatrixProdPole : Successfully added K-matrix production pole term. Total number of resonances so far = "<<nAmp_<<std::endl;
501 
502  } else {
503 
504  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdPole : The propagator of the name "<<propName
505  <<" could not be found for the production pole "<<poleName<<std::endl;
506 
507  }
508 
509 }
510 
511 
512 void LauIsobarDynamics::addKMatrixProdSVP(const TString& SVPName, const TString& propName, Int_t channelIndex)
513 {
514 
515  // Add a K-matrix production "slowly-varying part" (SVP) term, using the K-matrix propagator
516  // given by the propName. Here, channelIndex is the integer specifying the channel number.
517 
518  // First, find the K-matrix propagator.
519  KMPropMap::iterator mapIter = kMatrixPropagators_.find(propName);
520  if (mapIter != kMatrixPropagators_.end()) {
521 
522  LauKMatrixPropagator* thePropagator = mapIter->second;
523 
524  // Make sure the channel index is valid
525  Int_t nChannels = thePropagator->getNChannels();
526  if (channelIndex < 1 || channelIndex > nChannels) {
527  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdSVP : The channel index "<<channelIndex
528  <<" is not between 1 and "<<nChannels<<". Not adding production slowly-varying part "<<SVPName
529  <<" for K-matrix propagator "<<propName<<std::endl;
530  return;
531  }
532 
533  // Now add the K-matrix production SVP amplitude to the vector of LauAbsResonance pointers.
534  Int_t resPairAmpInt = thePropagator->getResPairAmpInt();
535  LauAbsResonance* prodSVP = new LauKMatrixProdSVP(SVPName, channelIndex, resPairAmpInt, thePropagator, daughters_);
536 
537  resTypAmp_.push_back(SVPName);
538  resIntAmp_.push_back(0);
539  resPairAmp_.push_back(resPairAmpInt);
540 
541  nAmp_++;
542  sigResonances_.push_back(prodSVP);
543 
544  // Also store the SVPName-propName pair for calculating total fit fractions later on
545  // (avoiding the need to use dynamic casts to check which resonances are of the K-matrix type)
546  kMatrixPropSet_[SVPName] = propName;
547 
548  std::cout<<"INFO in LauIsobarDynamics::addKMatrixProdSVP : Successfully added K-matrix production slowly-varying (SVP) term. Total number of resonances so far = "<<nAmp_<<std::endl;
549 
550  } else {
551 
552  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdSVP : The propagator of the name "<<propName
553  <<" could not be found for the production slowly-varying part "<<SVPName<<std::endl;
554 
555  }
556 }
557 
559 {
560  TString testName(name);
561  testName.ToLower();
562  LauAbsResonance* theResonance(0);
563 
564  for (std::vector<LauAbsResonance*>::iterator iter=sigResonances_.begin(); iter!=sigResonances_.end(); ++iter) {
565  theResonance = (*iter);
566  if (theResonance != 0) {
567  TString resString = theResonance->getResonanceName();
568  resString.ToLower();
569  if (resString.BeginsWith(testName, TString::kExact)) {
570  return theResonance;
571  }
572  }
573  }
574 
575  std::cerr<<"ERROR in LauIsobarDynamics::findResonance : Couldn't find resonance \""<<name<<"\" in the model."<<std::endl;
576  return 0;
577 }
578 
580 {
581  TString testName(name);
582  testName.ToLower();
583  const LauAbsResonance* theResonance(0);
584 
585  for (std::vector<LauAbsResonance*>::const_iterator iter=sigResonances_.begin(); iter!=sigResonances_.end(); ++iter) {
586  theResonance = (*iter);
587  if (theResonance != 0) {
588  TString resString = theResonance->getResonanceName();
589  resString.ToLower();
590  if (resString.BeginsWith(testName, TString::kExact)) {
591  return theResonance;
592  }
593  }
594  }
595 
596  std::cerr<<"ERROR in LauIsobarDynamics::findResonance : Couldn't find resonance \""<<name<<"\" in the model."<<std::endl;
597  return 0;
598 }
599 
600 void LauIsobarDynamics::removeCharge(TString& string) const
601 {
602  Ssiz_t index = string.Index("+");
603  if (index != -1) {
604  string.Remove(index,1);
605  }
606  index = string.Index("-");
607  if (index != -1) {
608  string.Remove(index,1);
609  }
610 }
611 
612 void LauIsobarDynamics::changeResonance(const TString& resName, Double_t newMass, Double_t newWidth, Int_t newSpin)
613 {
614  // Change the mass, width or spin of a resonance.
615 
616  if (newMass > 0.0 || newWidth > 0.0 || newSpin > -1) {
617  std::cerr<<"ERROR in LauIsobarDynamics::changeResonance : mass, width and spin parameters all out of range."<<std::endl;
618  return;
619  }
620 
621  LauAbsResonance* theRes = this->findResonance(resName);
622  if (theRes != 0) {
623  theRes->changeResonance(newMass, newWidth, newSpin);
624  }
625 
626 }
627 
629 {
630  // Use Gauss-Legendre quadrature integration
631 
632  // Get the rectangle that encloses the DP
633  Double_t minm13 = kinematics_->getm13Min();
634  Double_t maxm13 = kinematics_->getm13Max();
635  Double_t minm23 = kinematics_->getm23Min();
636  Double_t maxm23 = kinematics_->getm23Max();
637  Double_t minm12 = kinematics_->getm12Min();
638  Double_t maxm12 = kinematics_->getm12Max();
639 
640  // Find out whether we have narrow resonances in the DP (defined here as width < 20 MeV).
641  std::multimap<Double_t,Double_t> m13NarrowRes;
642  std::multimap<Double_t,Double_t> m23NarrowRes;
643  std::multimap<Double_t,Double_t> m12NarrowRes;
644  for ( std::vector<LauAbsResonance*>::const_iterator iter = sigResonances_.begin(); iter != sigResonances_.end(); ++iter ) {
645  Double_t width = (*iter)->getWidth();
646  if ( width > 0.020 || width == 0.0 ) { continue; }
647  Double_t mass = (*iter)->getMass();
648  Int_t pair = (*iter)->getPairInt();
649  TString name = (*iter)->getResonanceName();
650  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Found narrow resonance: "<<name<<", mass = "<<mass<<", width = "<<width<<", pair int = "<<pair<<std::endl;
651  if ( pair == 1 ) {
652  if ( mass < minm23 || mass > maxm23 ){ continue; }
653  m23NarrowRes.insert( std::make_pair(width,mass) );
654  } else if ( pair == 2 ) {
655  if ( mass < minm13 || mass > maxm13 ){ continue; }
656  m13NarrowRes.insert( std::make_pair(width,mass) );
657  } else if ( pair == 3 ) {
658  if ( mass < minm12 || mass > maxm12 ){ continue; }
659  m12NarrowRes.insert( std::make_pair(width,mass) );
660  } else {
661  std::cerr<<"WARNING in LauIsobarDynamics::calcDPNormalisation : strange pair integer, "<<pair<<", for resonance \""<<(*iter)->getResonanceName()<<std::endl;
662  }
663  }
664 
665  // Find the narrowest resonance in each mass pairing
666  const Bool_t e12 = m12NarrowRes.empty();
667  const Bool_t e13 = m13NarrowRes.empty();
668  const Bool_t e23 = m23NarrowRes.empty();
669  //const UInt_t s12 = e12 ? 0 : m12NarrowRes.size();
670  const UInt_t s13 = e13 ? 0 : m13NarrowRes.size();
671  const UInt_t s23 = e23 ? 0 : m23NarrowRes.size();
672  const Double_t w12 = e12 ? DBL_MAX : m12NarrowRes.begin()->first / 100.0;
673  const Double_t w13 = e13 ? DBL_MAX : m13NarrowRes.begin()->first / 100.0;
674  const Double_t w23 = e23 ? DBL_MAX : m23NarrowRes.begin()->first / 100.0;
675 
676  // Start off with default bin width (5 MeV)
677  Double_t m13BinWidth = m13BinWidth_;
678  Double_t m23BinWidth = m23BinWidth_;
679 
680  // Depending on how many narrow resonances we have and where they
681  // are we adopt different approaches
682  if ( e12 && e13 && e23 ) {
683  // If we have no narrow resonances just integrate the whole
684  // DP with the standard bin widths
685  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : No narrow resonances found, integrating over whole Dalitz plot..."<<std::endl;
686  this->calcDPPartialIntegral(minm13, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
687  } else if ( ! e12 ) {
688  // If we have a narrow resonance on the diagonal then we'll have to
689  // just use a narrow bin width over the whole DP (1/10 of the width
690  // of the narrowest resonance in any mass pair)
691  m13BinWidth = m23BinWidth = 10.0*TMath::Min( w12, TMath::Min( w13, w23 ) );
692  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One or more narrow resonances found in m12, integrating over whole Dalitz plot with bin width of "<<m13BinWidth<<" GeV/c2..."<<std::endl;
693  this->calcDPPartialIntegral(minm13, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
694  } else if ( s13==1 && e23 ) {
695  // We have a single narrow resonance in m13
696  // Divide the plot into 3 regions: the resonance band and
697  // the two areas either side.
698  Double_t mass = m13NarrowRes.begin()->second;
699  Double_t width = m13NarrowRes.begin()->first;
700  Double_t resMin = mass - 5.0*width;
701  Double_t resMax = mass + 5.0*width;
702  // if the resonance is close to threshold just go from
703  // threshold to resMax, otherwise treat threshold to resMin
704  // as a separate region
705  if ( resMin < (minm13+50.0*m13BinWidth_) ) {
706  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m13, close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
707  this->calcDPPartialIntegral(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
708  m13BinWidth = w13;
709  this->calcDPPartialIntegral(minm13, resMax, minm23, maxm23, m13BinWidth, m23BinWidth);
710  } else {
711  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m13, dividing Dalitz plot into three regions..."<<std::endl;
712  this->calcDPPartialIntegral(minm13, resMin, minm23, maxm23, m13BinWidth, m23BinWidth);
713  this->calcDPPartialIntegral(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
714  m13BinWidth = w13;
715  this->calcDPPartialIntegral(resMin, resMax, minm23, maxm23, m13BinWidth, m23BinWidth);
716  }
717  } else if ( s13==2 && e23 ) {
718  // We have a two narrow resonances in m13
719  // Divide the plot into 5 regions: the resonance bands,
720  // the two areas either side and the area in between.
721  std::multimap<Double_t,Double_t> massordered;
722  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
723  massordered.insert( std::make_pair( iter->second, iter->first ) );
724  }
725  Double_t res1Mass = massordered.begin()->first;
726  Double_t res1Width = massordered.begin()->second;
727  Double_t res1Min = res1Mass - 5.0*res1Width;
728  Double_t res1Max = res1Mass + 5.0*res1Width;
729  Double_t res2Mass = (++(massordered.begin()))->first;
730  Double_t res2Width = (++(massordered.begin()))->second;
731  Double_t res2Min = res2Mass - 5.0*res2Width;
732  Double_t res2Max = res2Mass + 5.0*res2Width;
733  // if the resonance is close to threshold just go from
734  // threshold to resMax, otherwise treat threshold to resMin
735  // as a separate region
736  if ( res1Min < (minm13+50.0*m13BinWidth_) ) {
737  if ( res1Max > res2Min ) {
738  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m13, both close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
739  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
740  m13BinWidth = TMath::Min(res1Width,res2Width)/100.0;
741  this->calcDPPartialIntegral(minm13, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
742  } else {
743  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m13, one close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
744  this->calcDPPartialIntegral(res1Max, res2Min, minm23, maxm23, m13BinWidth, m23BinWidth);
745  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
746  m13BinWidth = res1Width/100.0;
747  this->calcDPPartialIntegral(minm13, res1Max, minm23, maxm23, m13BinWidth, m23BinWidth);
748  m13BinWidth = res2Width/100.0;
749  this->calcDPPartialIntegral(res2Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
750  }
751  } else {
752  if ( res1Max > res2Min ) {
753  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found close together in m13, dividing Dalitz plot into three regions..."<<std::endl;
754  this->calcDPPartialIntegral(minm13, res1Min, minm23, maxm23, m13BinWidth, m23BinWidth);
755  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
756  m13BinWidth = TMath::Min(res1Width,res2Width)/100.0;
757  this->calcDPPartialIntegral(res1Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
758  } else {
759  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m13, dividing Dalitz plot into five regions..."<<std::endl;
760  this->calcDPPartialIntegral(minm13, res1Min, minm23, maxm23, m13BinWidth, m23BinWidth);
761  this->calcDPPartialIntegral(res1Max, res2Min, minm23, maxm23, m13BinWidth, m23BinWidth);
762  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
763  m13BinWidth = res1Width/100.0;
764  this->calcDPPartialIntegral(res1Min, res1Max, minm23, maxm23, m13BinWidth, m23BinWidth);
765  m13BinWidth = res2Width/100.0;
766  this->calcDPPartialIntegral(res2Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
767  }
768  }
769  } else if ( s23==1 && e13 ) {
770  // We have a single narrow resonance in m23
771  // Divide the plot into 3 regions: the resonance band and
772  // the two areas either side.
773  Double_t mass = m23NarrowRes.begin()->second;
774  Double_t width = m23NarrowRes.begin()->first;
775  Double_t resMin = mass - 5.0*width;
776  Double_t resMax = mass + 5.0*width;
777  // if the resonance is close to threshold just go from
778  // threshold to resMax, otherwise treat threshold to resMin
779  // as a separate region
780  if ( resMin < (minm23+50.0*m23BinWidth_) ) {
781  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m23, close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
782  this->calcDPPartialIntegral(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth);
783  m23BinWidth = w23;
784  this->calcDPPartialIntegral(minm13, maxm13, minm23, resMax, m13BinWidth, m23BinWidth);
785  } else {
786  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m23, dividing Dalitz plot into three regions..."<<std::endl;
787  this->calcDPPartialIntegral(minm13, maxm13, minm23, resMin, m13BinWidth, m23BinWidth);
788  this->calcDPPartialIntegral(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth);
789  m23BinWidth = w23;
790  this->calcDPPartialIntegral(minm13, maxm13, resMin, resMax, m13BinWidth, m23BinWidth);
791  }
792  } else if ( s23==2 && e13 ) {
793  // We have a two narrow resonances in m23
794  // Divide the plot into 5 regions: the resonance bands,
795  // the two areas either side and the area in between.
796  std::multimap<Double_t,Double_t> massordered;
797  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
798  massordered.insert( std::make_pair( iter->second, iter->first ) );
799  }
800  Double_t res1Mass = massordered.begin()->first;
801  Double_t res1Width = massordered.begin()->second;
802  Double_t res1Min = res1Mass - 5.0*res1Width;
803  Double_t res1Max = res1Mass + 5.0*res1Width;
804  Double_t res2Mass = (++(massordered.begin()))->first;
805  Double_t res2Width = (++(massordered.begin()))->second;
806  Double_t res2Min = res2Mass - 5.0*res2Width;
807  Double_t res2Max = res2Mass + 5.0*res2Width;
808  // if the resonance is close to threshold just go from
809  // threshold to resMax, otherwise treat threshold to resMin
810  // as a separate region
811  if ( res1Min < (minm23+50.0*m23BinWidth_) ) {
812  if ( res1Max > res2Min ) {
813  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m23, both close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
814  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
815  m23BinWidth = TMath::Min(res1Width,res2Width)/100.0;
816  this->calcDPPartialIntegral(minm13, maxm13, minm23, res2Max, m13BinWidth, m23BinWidth);
817  } else {
818  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m23, one close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
819  this->calcDPPartialIntegral(minm13, maxm13, res1Max, res2Min, m13BinWidth, m23BinWidth);
820  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
821  m23BinWidth = res1Width/100.0;
822  this->calcDPPartialIntegral(minm13, maxm13, minm23, res1Max, m13BinWidth, m23BinWidth);
823  m23BinWidth = res2Width/100.0;
824  this->calcDPPartialIntegral(minm13, maxm13, res2Min, res2Max, m13BinWidth, m23BinWidth);
825  }
826  } else {
827  if ( res1Max > res2Min ) {
828  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found close together in m23, dividing Dalitz plot into three regions..."<<std::endl;
829  this->calcDPPartialIntegral(minm13, maxm13, minm23, res1Min, m13BinWidth, m23BinWidth);
830  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
831  m23BinWidth = TMath::Min(res1Width,res2Width)/100.0;
832  this->calcDPPartialIntegral(minm13, maxm13, res1Min, res2Max, m13BinWidth, m23BinWidth);
833  } else {
834  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m23, dividing Dalitz plot into five regions..."<<std::endl;
835  this->calcDPPartialIntegral(minm13, maxm13, minm23, res1Min, m13BinWidth, m23BinWidth);
836  this->calcDPPartialIntegral(minm13, maxm13, res1Max, res2Min, m13BinWidth, m23BinWidth);
837  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
838  m23BinWidth = res1Width/100.0;
839  this->calcDPPartialIntegral(minm13, maxm13, res1Min, res1Max, m13BinWidth, m23BinWidth);
840  m23BinWidth = res2Width/100.0;
841  this->calcDPPartialIntegral(minm13, maxm13, res2Min, res2Max, m13BinWidth, m23BinWidth);
842  }
843  }
844  } else if ( s13==1 && s23==1 ) {
845  // We have a single narrow resonance in both m13 and m23
846  // Divide the plot into 9 regions: the point where the
847  // resonance bands cross, the four other parts of the bands
848  // and the four remaining areas of the DP.
849  Double_t mass23 = m23NarrowRes.begin()->second;
850  Double_t width23 = m23NarrowRes.begin()->first;
851  Double_t resMin23 = mass23 - 5.0*width23;
852  Double_t resMax23 = mass23 + 5.0*width23;
853  Double_t mass13 = m13NarrowRes.begin()->second;
854  Double_t width13 = m13NarrowRes.begin()->first;
855  Double_t resMin13 = mass13 - 5.0*width13;
856  Double_t resMax13 = mass13 + 5.0*width13;
857  // if either resonance is close to threshold just go from
858  // threshold to resMax, otherwise treat threshold to resMin
859  // as a separate region
860  if ( resMin13 < (minm13+50.0*m13BinWidth_) && resMin23 < (minm23+50.0*m23BinWidth_) ) {
861  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m13 and one in m23, both close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
862  m13BinWidth = m13BinWidth_;
863  m23BinWidth = m23BinWidth_;
864  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
865  m13BinWidth = m13BinWidth_;
866  m23BinWidth = w23;
867  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
868  m13BinWidth = w13;
869  m23BinWidth = m23BinWidth_;
870  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
871  m13BinWidth = w13;
872  m23BinWidth = w23;
873  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
874  } else if ( resMin13 < (minm13+50.0*m13BinWidth_) ) {
875  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m13, close to threshold, and one in m23, not close to threshold, dividing Dalitz plot into six regions..."<<std::endl;
876  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth);
877  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
878  m13BinWidth = m13BinWidth_;
879  m23BinWidth = w23;
880  this->calcDPPartialIntegral(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth);
881  m13BinWidth = w13;
882  m23BinWidth = m23BinWidth_;
883  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth);
884  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
885  m13BinWidth = w13;
886  m23BinWidth = w23;
887  this->calcDPPartialIntegral(minm13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth);
888  } else if ( resMin23 < (minm23+50.0*m23BinWidth_) ) {
889  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m23, close to threshold, and one in m13, not close to threshold, dividing Dalitz plot into six regions..."<<std::endl;
890  this->calcDPPartialIntegral(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth);
891  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
892  m13BinWidth = m13BinWidth_;
893  m23BinWidth = w23;
894  this->calcDPPartialIntegral(minm13, resMin13, minm23, resMax23, m13BinWidth, m23BinWidth);
895  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
896  m13BinWidth = w13;
897  m23BinWidth = m23BinWidth_;
898  this->calcDPPartialIntegral(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
899  m13BinWidth = w13;
900  m23BinWidth = w23;
901  this->calcDPPartialIntegral(resMin13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
902  } else {
903  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in both m13 and m23, neither close to threshold, dividing Dalitz plot into nine regions..."<<std::endl;
904  this->calcDPPartialIntegral(minm13, resMin13, minm23, resMin23, m13BinWidth, m23BinWidth);
905  this->calcDPPartialIntegral(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth);
906  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth);
907  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
908  m13BinWidth = m13BinWidth_;
909  m23BinWidth = w23;
910  this->calcDPPartialIntegral(minm13, resMin13, resMin23, resMax23, m13BinWidth, m23BinWidth);
911  this->calcDPPartialIntegral(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth);
912  m13BinWidth = w13;
913  m23BinWidth = m23BinWidth_;
914  this->calcDPPartialIntegral(resMin13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth);
915  this->calcDPPartialIntegral(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
916  m13BinWidth = w13;
917  m23BinWidth = w23;
918  this->calcDPPartialIntegral(resMin13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth);
919  }
920  } else if ( e23 && s13>1 ) {
921  // We have multiple narrow resonances in m13 only.
922  // Divide the plot into 2 regions: threshold to the most
923  // massive of the narrow resonances, and the rest
924  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m13, dividing Dalitz plot into two regions..."<<std::endl;
925  Double_t mass = 0.0;
926  Double_t width = 0.0;
927  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
928  if ( mass < iter->second ) {
929  mass = iter->second;
930  width = iter->first;
931  }
932  }
933  Double_t resMax = mass + 5.0*width;
934  this->calcDPPartialIntegral(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
935  m13BinWidth = w13;
936  this->calcDPPartialIntegral(minm13, resMax, minm23, maxm23, m13BinWidth, m23BinWidth);
937  } else if ( e13 && s23>1 ) {
938  // We have multiple narrow resonances in m23 only.
939  // Divide the plot into 2 regions: threshold to the most
940  // massive of the narrow resonances, and the rest
941  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m23, dividing Dalitz plot into two regions..."<<std::endl;
942  Double_t mass = 0.0;
943  Double_t width = 0.0;
944  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
945  if ( mass < iter->second ) {
946  mass = iter->second;
947  width = iter->first;
948  }
949  }
950  Double_t resMax = mass + 5.0*width;
951  this->calcDPPartialIntegral(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth);
952  m23BinWidth = w23;
953  this->calcDPPartialIntegral(minm13, maxm13, minm23, resMax, m13BinWidth, m23BinWidth);
954  } else if ( s13==1 && s23>1 ) {
955  // We've got a single narrow resonance in m13 and multiple
956  // narrow resonances in m23. Divide the plot into 6 regions.
957  Double_t mass23 = 0.0;
958  Double_t width23 = 0.0;
959  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
960  if ( mass23 < iter->second ) {
961  mass23 = iter->second;
962  width23 = iter->first;
963  }
964  }
965  Double_t resMax23 = mass23 + 5.0*width23;
966  Double_t mass13 = m13NarrowRes.begin()->second;
967  Double_t width13 = m13NarrowRes.begin()->first;
968  Double_t resMin13 = mass13 - 5.0*width13;
969  Double_t resMax13 = mass13 + 5.0*width13;
970  // if the m13 resonance is close to threshold just go from
971  // threshold to resMax, otherwise treat threshold to resMin
972  // as a separate region
973  if ( resMin13 < (minm13+50.0*m13BinWidth_) ) {
974  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m23 and one in m13, close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
975  m13BinWidth = m13BinWidth_;
976  m23BinWidth = m23BinWidth_;
977  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
978  m13BinWidth = m13BinWidth_;
979  m23BinWidth = w23;
980  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
981  m13BinWidth = w13;
982  m23BinWidth = m23BinWidth_;
983  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
984  m13BinWidth = w13;
985  m23BinWidth = w23;
986  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
987  } else {
988  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m23 and one in m13, not close to threshold, dividing Dalitz plot into six regions..."<<std::endl;
989  m13BinWidth = m13BinWidth_;
990  m23BinWidth = m23BinWidth_;
991  this->calcDPPartialIntegral(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth);
992  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
993  m13BinWidth = m13BinWidth_;
994  m23BinWidth = w23;
995  this->calcDPPartialIntegral(minm13, resMin13, minm23, resMax23, m13BinWidth, m23BinWidth);
996  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
997  m13BinWidth = w13;
998  m23BinWidth = m23BinWidth_;
999  this->calcDPPartialIntegral(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1000  m13BinWidth = w13;
1001  m23BinWidth = w23;
1002  this->calcDPPartialIntegral(resMin13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1003  }
1004  } else if ( s13>1 && s23==1 ) {
1005  // We've got a single narrow resonance in m23 and multiple
1006  // narrow resonances in m13. Divide the plot into 6 regions.
1007  Double_t mass13 = 0.0;
1008  Double_t width13 = 0.0;
1009  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1010  if ( mass13 < iter->second ) {
1011  mass13 = iter->second;
1012  width13 = iter->first;
1013  }
1014  }
1015  Double_t resMax13 = mass13 + 5.0*width13;
1016  Double_t mass23 = m23NarrowRes.begin()->second;
1017  Double_t width23 = m23NarrowRes.begin()->first;
1018  Double_t resMin23 = mass23 - 5.0*width23;
1019  Double_t resMax23 = mass23 + 5.0*width23;
1020  // if the m23 resonance is close to threshold just go from
1021  // threshold to resMax, otherwise treat threshold to resMin
1022  // as a separate region
1023  if ( resMin23 < (minm23+50.0*m23BinWidth_) ) {
1024  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m13 and one in m23, close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
1025  m13BinWidth = m13BinWidth_;
1026  m23BinWidth = m23BinWidth_;
1027  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1028  m13BinWidth = m13BinWidth_;
1029  m23BinWidth = w23;
1030  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
1031  m13BinWidth = w13;
1032  m23BinWidth = m23BinWidth_;
1033  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1034  m13BinWidth = w13;
1035  m23BinWidth = w23;
1036  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1037  } else {
1038  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m13 and one in m23, not close to threshold, dividing Dalitz plot into six regions..."<<std::endl;
1039  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth);
1040  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1041  m13BinWidth = m13BinWidth_;
1042  m23BinWidth = w23;
1043  this->calcDPPartialIntegral(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1044  m13BinWidth = w13;
1045  m23BinWidth = m23BinWidth_;
1046  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth);
1047  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1048  m13BinWidth = w13;
1049  m23BinWidth = w23;
1050  this->calcDPPartialIntegral(minm13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1051  }
1052  } else {
1053  // We've got multiple narrow resonances in both m13 and m23.
1054  // Divide the plot into 4 regions.
1055  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in both m13 and m23, dividing Dalitz plot into four regions..."<<std::endl;
1056  Double_t mass13 = 0.0;
1057  Double_t width13 = 0.0;
1058  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1059  if ( mass13 < iter->second ) {
1060  mass13 = iter->second;
1061  width13 = iter->first;
1062  }
1063  }
1064  Double_t resMax13 = mass13 + 5.0*width13;
1065 
1066  Double_t mass23 = 0.0;
1067  Double_t width23 = 0.0;
1068  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1069  if ( mass23 < iter->second ) {
1070  mass23 = iter->second;
1071  width23 = iter->first;
1072  }
1073  }
1074  Double_t resMax23 = mass23 + 5.0*width23;
1075 
1076  m13BinWidth = m13BinWidth_;
1077  m23BinWidth = m23BinWidth_;
1078  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1079  m13BinWidth = m13BinWidth_;
1080  m23BinWidth = w23;
1081  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
1082  m13BinWidth = w13;
1083  m23BinWidth = m23BinWidth_;
1084  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1085  m13BinWidth = w13;
1086  m23BinWidth = w23;
1087  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1088  }
1089 }
1090 
1091 void LauIsobarDynamics::setIntegralBinWidths(Double_t m13BinWidth, Double_t m23BinWidth)
1092 {
1093  // Specify whether we're going to use Gauss-Legendre integration to calculate the normalisation
1094  // integrals, and the bin widths we require for the m13 and m23 axes. Note that the integration
1095  // is done over m13, m23 space, with the appropriate Jacobian applied, and not m13^2, m23^2 space.
1096  // The default bin widths in m13 and m23 space are 5 MeV.
1097 
1098  m13BinWidth_ = m13BinWidth;
1099  m23BinWidth_ = m23BinWidth;
1100 }
1101 
1102 void LauIsobarDynamics::calcDPPartialIntegral(Double_t minm13, Double_t maxm13, Double_t minm23, Double_t maxm23,
1103  Double_t m13BinWidth, Double_t m23BinWidth)
1104 {
1105  // Calculate the total DP area, as well as finding the normalisation terms for
1106  // the signal resonances
1107 
1108  Int_t i(0), j(0);
1109  Double_t precision(1e-6);
1110 
1111  Double_t meanm13 = 0.5*(minm13 + maxm13);
1112  Double_t rangem13 = maxm13 - minm13;
1113  Double_t halfRangem13 = 0.5*rangem13;
1114 
1115  Double_t meanm23 = 0.5*(minm23 + maxm23);
1116  Double_t rangem23 = maxm23 - minm23;
1117  Double_t halfRangem23 = 0.5*rangem23;
1118 
1119  Double_t intFactor = halfRangem13*halfRangem23;
1120 
1121  // Choose smallest of mass ranges to set number of bins, given specified bin width
1122  Int_t nm13Points = static_cast<Int_t>((rangem13/m13BinWidth));
1123  Int_t nm23Points = static_cast<Int_t>((rangem23/m23BinWidth));
1124 
1125  // Avoid integral if we have no points in either x or y space
1126  if (nm13Points == 0 || nm23Points == 0) {return;}
1127 
1128  std::cout<<"INFO in LauIsobarDynamics::calcDPPartialIntegral : nm13Points = "<<nm13Points<<", nm23Points = "<<nm23Points<<std::endl;
1129  std::cout<<" : m13BinWidth = "<<m13BinWidth<<", m23BinWidth = "<<m23BinWidth<<std::endl;
1130  std::cout<<" : Integrating over m13 = "<<minm13<<" to "<<maxm13<<", m23 = "<<minm23<<" to "<<maxm23<<std::endl;
1131 
1132  LauIntegrals DPIntegrals(precision);
1133  std::vector<Double_t> m13Weights, m23Weights;
1134  std::vector<Double_t> m13Abscissas, m23Abscissas;
1135 
1136  DPIntegrals.calcGaussLegendreWeights(nm13Points, m13Abscissas, m13Weights);
1137  DPIntegrals.calcGaussLegendreWeights(nm23Points, m23Abscissas, m23Weights);
1138 
1139  Int_t nm13Weights = static_cast<Int_t>(m13Weights.size());
1140  Int_t nm23Weights = static_cast<Int_t>(m23Weights.size());
1141 
1142  //std::cout<<" : nm13Weights = "<<nm13Weights<<", nm23Weights = "<<nm23Weights<<std::endl;
1143  // Print out abscissas and weights for the integration
1144  Double_t totm13Weight(0.0), totm23Weight(0.0);
1145  for (i = 0; i < nm13Weights; i++) {
1146  totm13Weight += m13Weights[i];
1147  }
1148  for (i = 0; i < nm23Weights; i++) {
1149  totm23Weight += m23Weights[i];
1150  }
1151  std::cout<<" : totm13Weight = "<<totm13Weight<<", totm23Weight = "<<totm23Weight<<std::endl;
1152 
1153  std::vector<Double_t> m13(nm13Weights), m23(nm23Weights);
1154  std::vector<Double_t> m13Sq(nm13Weights), m23Sq(nm23Weights);
1155 
1156  // Use same number of abscissas for x and y co-ordinates
1157  Int_t m = (nm13Weights + 1)/2;
1158  for (i = 0; i < m; i++) {
1159 
1160  Int_t ii = nm13Weights - 1 - i; // symmetric i index
1161 
1162  Double_t dm13 = halfRangem13*m13Abscissas[i];
1163  Double_t m13Val = meanm13 - dm13;
1164  m13[i] = m13Val;
1165  m13Sq[i] = m13Val*m13Val;
1166 
1167  m13Val = meanm13 + dm13;
1168  m13[ii] = m13Val;
1169  m13Sq[ii] = m13Val*m13Val;
1170 
1171  }
1172 
1173  m = (nm23Weights +1)/2;
1174  for (i = 0; i < m; i++) {
1175 
1176  Int_t ii = nm23Weights - 1 - i; // symmetric i index
1177 
1178  Double_t dm23 = halfRangem23*m23Abscissas[i];
1179  Double_t m23Val = meanm23 - dm23;
1180  m23[i] = m23Val;
1181  m23Sq[i] = m23Val*m23Val;
1182  m23Val = meanm23 + dm23;
1183  m23[ii] = m23Val;
1184  m23Sq[ii] = m23Val*m23Val;
1185  }
1186 
1187  // Now compute the integral
1188  Double_t dpArea(0.0);
1189  for (i = 0; i < nm13Weights; i++) {
1190 
1191  for (j = 0; j < nm23Weights; j++) {
1192 
1193  Double_t weight = m13Weights[i]*m23Weights[j];
1194  Double_t Jacobian = 4.0*m13[i]*m23[j];
1195  weight *= (Jacobian*intFactor);
1196 
1197  // Calculate the integral contributions for each resonance.
1198  // Only resonances within the DP area contribute.
1199  // This also calculates the total DP area as a check.
1200  Bool_t withinDP = kinematics_->withinDPLimits(m13Sq[i], m23Sq[j]);
1201  if (withinDP == kTRUE) {
1202 
1203  kinematics_->updateKinematics(m13Sq[i], m23Sq[j]);
1204  this->dynamics(kFALSE, weight);
1205  // Increment total DP area
1206  dpArea += weight;
1207 
1208  }
1209 
1210  } // j weights loop
1211  } // i weights loop
1212 
1213  // Print out DP area to check whether we have a sensible output
1214  std::cout<<" : dpArea = "<<dpArea<<std::endl;
1215 
1216 }
1217 
1218 void LauIsobarDynamics::dynamics(Bool_t cacheResData, Double_t weight, Bool_t useEff)
1219 {
1220  // Routine that calculates the Dalitz plot amplitude, incorporating
1221  // resonance dynamics (using resAmp()) and any interference between them.
1222  // Used by the fit() and sigGen() functions.
1223 
1224  UInt_t i(0), j(0);
1225 
1226  // Reset the total amplitude to zero
1227  totAmp_.zero();
1228 
1229  // Loop over the number of resonance amplitudes defined in the model
1230  // Have we already calculated this for this event (during fit?)
1231  // Or do we have a resonance that has varying pole mass/width/other factors?
1232  if (cacheResData == kFALSE) {
1233 
1234  for (i = 0; i < nAmp_; i++) {
1235 
1236  // Calculate the dynamics for this resonance, using the resAmp function.
1237  ff_[i] = resAmp(i);
1238 
1239  //std::cout<<"ff_["<<i<<"] = "<<ff_[i]<<std::endl;
1240 
1241  // If we have a symmetrical Dalitz plot, flip the m_13^2 and m_23^2
1242  // variables, recalculate the dynamics, and average both contributions.
1243  // Although, the factor of 0.5 cancels out with the fact that
1244  // the resonance appears on both sides of the Dalitz plot.
1245  if (symmetricalDP_ == kTRUE) { // was tieSg_ == 12
1247  ff_[i] += resAmp(i);
1248  // Flip the m_13^2 and m_23^2 variables back to their original values
1250  }
1251 
1252  } // Loop over amplitudes
1253 
1254  // If we haven't cached the data, then we need to find out the efficiency.
1255  eff_ = this->retrieveEfficiency();
1256 
1257  } // Already cached data?
1258 
1259  if (integralsDone_ == kTRUE) {
1260 
1261  // Loop over all signal amplitudes
1262  LauComplex ATerm;
1263  for (i = 0; i < nAmp_; i++) {
1264 
1265  // Get the partial complex amplitude - (mag, phase)*(resonance dynamics)
1266  ATerm = Amp_[i]*ff_[i];
1267  // Scale this contribution by its relative normalisation w.r.t. the whole dynamics
1268  ATerm.rescale(fNorm_[i]);
1269 
1270  // Add this partial amplitude to the sum
1271  //std::cout<<"For i = "<<i<<", ATerm = "<<ATerm<<", Amp = "<<Amp_[i]<<", ff = "<<ff_[i]<<std::endl;
1272  totAmp_ += ATerm;
1273 
1274  } // Loop over amplitudes
1275 
1276  // |Sum of partial amplitudes|^2
1277  ASq_ = totAmp_.abs2();
1278 
1279  // Apply the efficiency correction for this event.
1280  // Multiply the amplitude squared sum by the DP efficiency
1281  if ( useEff ) {
1282  ASq_ *= eff_;
1283  }
1284 
1285  } else { // integrals not done
1286 
1287  // Find the efficiency correction to be applied for this event.
1288  eff_ = this->retrieveEfficiency();
1289 
1290  Double_t effWeight = eff_*weight;
1291 
1292  // Need this for integrals for normalisation of likelihood function.
1293  LauComplex fifjEffSumTerm;
1294  LauComplex fifjSumTerm;
1295  for (i = 0; i < nAmp_; i++) {
1296 
1297  // Add the dynamical amplitude squared for this resonance.
1298  Double_t fSqVal = ff_[i].abs2();
1299  fSqSum_[i] += fSqVal*weight;
1300 
1301  for (j = i; j < nAmp_; j++) {
1302 
1303  fifjEffSumTerm = fifjSumTerm = ff_[i]*ff_[j].conj();
1304 
1305  fifjEffSumTerm.rescale(effWeight);
1306  fifjEffSum_[i][j] += fifjEffSumTerm;
1307 
1308  fifjSumTerm.rescale(weight);
1309  fifjSum_[i][j] += fifjSumTerm;
1310  }
1311  }
1312  }
1313 }
1314 
1316 {
1317  // Routine to calculate the resonance dynamics (amplitude)
1318  // using the appropriate Breit-Wigner/Form Factors.
1319  // Called by the dynamics() function.
1320 
1321  // Get the signal resonance from the stored vector
1322  LauAbsResonance* sigResonance = sigResonances_[index];
1323 
1324  if (sigResonance == 0) {
1325  std::cout<<"ERROR in LauIsobarDynamics::resAmp : Couldn't retrieve resonance with index = "<<index<<std::endl;
1326  return LauComplex(0.0, 0.0);
1327  }
1328 
1329  // Get the integer index of the resonance.
1330  Int_t resInt = resIntAmp_[index];
1331 
1332  LauComplex resAmplitude(0.0, 0.0);
1333 
1334  if (resInt < 0 || resInt >= static_cast<Int_t>(this->getnDefinedResonances())) {
1335 
1336  std::cout<<"ERROR in LauIsobarDynamics::resAmp : Probably bad resonance name."<<std::endl;
1337  resAmplitude = LauComplex(0.0, 0.0);
1338 
1339  } else {
1340 
1341  resAmplitude = sigResonance->amplitude(kinematics_);
1342 
1343  }
1344 
1345  return resAmplitude;
1346 }
1347 
1348 void LauIsobarDynamics::setFFTerm(UInt_t index, Double_t realPart, Double_t imagPart)
1349 {
1350  // Function to set the internal ff = resAmp() term.
1351  if ( index >= nAmp_ ) {
1352  std::cerr<<"ERROR in LauIsobarDynamics::setFFTerm : index = "<<index<<" is not within the range 0 and "<<nAmp_-1<<". Bailing out."<<std::endl;
1353  return;
1354  }
1355 
1356  ff_[index].setRealImagPart( realPart, imagPart );
1357 }
1358 
1360 {
1361  // This method calculates the fit fractions, mean efficiency and total DP rate
1362 
1363  Double_t fifjEffTot(0.0), fifjTot(0.0);
1364  UInt_t i, j;
1365  for (i = 0; i < nAmp_; i++) {
1366 
1367  // Calculate the diagonal terms
1368  TString name = "A"; name += i; name += "Sq_FitFrac";
1369  fitFrac_[i][i].name(name);
1370 
1371  name += "EffUnCorr";
1372  fitFracEffUnCorr_[i][i].name(name);
1373 
1374  Double_t fifjSumReal = fifjSum_[i][i].re();
1375  Double_t sumTerm = Amp_[i].abs2()*fifjSumReal*fNorm_[i]*fNorm_[i];
1376  fifjTot += sumTerm;
1377 
1378  Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1379  Double_t sumEffTerm = Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1380  fifjEffTot += sumEffTerm;
1381 
1382  fitFrac_[i][i] = sumTerm;
1383  fitFracEffUnCorr_[i][i] = sumEffTerm;
1384  }
1385 
1386  for (i = 0; i < nAmp_; i++) {
1387  for (j = i+1; j < nAmp_; j++) {
1388  // Calculate the cross-terms
1389  TString name = "A"; name += i; name += "A"; name += j; name += "_FitFrac";
1390  fitFrac_[i][j].name(name);
1391 
1392  name += "EffUnCorr";
1393  fitFracEffUnCorr_[i][j].name(name);
1394 
1395  LauComplex AmpjConj = Amp_[j].conj();
1396  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1397 
1398  Double_t crossTerm = 2.0*(AmpTerm*fifjSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1399  fifjTot += crossTerm;
1400 
1401  Double_t crossEffTerm = 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1402  fifjEffTot += crossEffTerm;
1403 
1404  fitFrac_[i][j] = crossTerm;
1405  fitFracEffUnCorr_[i][j] = crossEffTerm;
1406  }
1407  }
1408 
1409  if (TMath::Abs(fifjTot) > 1e-10) {
1410  meanDPEff_ = fifjEffTot/fifjTot;
1411  if (init) {
1414  }
1415  }
1416  DPRate_ = fifjTot;
1417  if (init) {
1420  }
1421 
1422  // Now divide the fitFraction sums by the overall integral
1423  for (i = 0; i < nAmp_; i++) {
1424  for (j = i; j < nAmp_; j++) {
1425  // Get the actual fractions by dividing by the total DP rate
1426  fitFrac_[i][j] /= fifjTot;
1427  fitFracEffUnCorr_[i][j] /= fifjEffTot;
1428  if (init) {
1429  fitFrac_[i][j].genValue( fitFrac_[i][j].value() );
1430  fitFrac_[i][j].initValue( fitFrac_[i][j].value() );
1431  fitFracEffUnCorr_[i][j].genValue( fitFracEffUnCorr_[i][j].value() );
1432  fitFracEffUnCorr_[i][j].initValue( fitFracEffUnCorr_[i][j].value() );
1433  }
1434  }
1435  }
1436 
1437  // Work out total fit fraction over all K-matrix components (for each propagator)
1438  KMPropMap::iterator mapIter;
1439  Int_t propInt(0);
1440 
1441  for (mapIter = kMatrixPropagators_.begin(); mapIter != kMatrixPropagators_.end(); ++mapIter) {
1442 
1443  LauKMatrixPropagator* thePropagator = mapIter->second;
1444 
1445  TString propName = thePropagator->getName();
1446 
1447  // Now loop over all resonances and find those which are K-matrix components for this propagator
1448  Double_t kMatrixTotFitFrac(0.0);
1449 
1450  for (i = 0; i < nAmp_; i++) {
1451 
1452  Bool_t gotKMRes1 = this->gotKMatrixMatch(i, propName);
1453  if (gotKMRes1 == kFALSE) {continue;}
1454 
1455  Double_t fifjSumReal = fifjSum_[i][i].re();
1456  Double_t sumTerm = Amp_[i].abs2()*fifjSumReal*fNorm_[i]*fNorm_[i];
1457 
1458  //Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1459  //Double_t sumEffTerm = Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1460 
1461  kMatrixTotFitFrac += sumTerm;
1462 
1463  for (j = i+1; j < nAmp_; j++) {
1464 
1465  Bool_t gotKMRes2 = this->gotKMatrixMatch(j, propName);
1466  if (gotKMRes2 == kFALSE) {continue;}
1467 
1468  LauComplex AmpjConj = Amp_[j].conj();
1469  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1470 
1471  Double_t crossTerm = 2.0*(AmpTerm*fifjSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1472  //Double_t crossEffTerm = 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1473 
1474  kMatrixTotFitFrac += crossTerm;
1475 
1476  }
1477 
1478  }
1479 
1480  kMatrixTotFitFrac /= fifjTot;
1481 
1482  TString parName("KMatrixTotFF_"); parName += propInt;
1483  extraParameters_[propInt].name( parName );
1484  extraParameters_[propInt] = kMatrixTotFitFrac;
1485  if (init) {
1486  extraParameters_[propInt].genValue(kMatrixTotFitFrac);
1487  extraParameters_[propInt].initValue(kMatrixTotFitFrac);
1488  }
1489 
1490  std::cout<<"INFO in LauIsobarDynamics::calcExtraInfo : Total K-matrix fit fraction for propagator "<<propName<<" is "<<kMatrixTotFitFrac<<std::endl;
1491 
1492  ++propInt;
1493 
1494  }
1495 
1496 
1497 }
1498 
1499 Bool_t LauIsobarDynamics::gotKMatrixMatch(UInt_t resAmpInt, const TString& propName) const
1500 {
1501 
1502  Bool_t gotMatch(kFALSE);
1503 
1504  if (resAmpInt >= nAmp_) {return kFALSE;}
1505 
1506  const LauAbsResonance* theResonance = sigResonances_[resAmpInt];
1507 
1508  if (theResonance == 0) {return kFALSE;}
1509 
1510  Int_t resModelInt = theResonance->getResonanceModel();
1511 
1512  if (resModelInt == LauAbsResonance::KMatrix) {
1513 
1514  TString resName = theResonance->getResonanceName();
1515 
1516  KMStringMap::const_iterator kMPropSetIter = kMatrixPropSet_.find(resName);
1517 
1518  if (kMPropSetIter != kMatrixPropSet_.end()) {
1519  TString kmPropString = kMPropSetIter->second;
1520  if (kmPropString == propName) {gotMatch = kTRUE;}
1521  }
1522 
1523  }
1524 
1525  return gotMatch;
1526 
1527 }
1528 
1530 {
1531  // Calculate the normalisation for the log-likelihood function.
1532  DPNorm_ = 0.0;
1533 
1534  for (UInt_t i = 0; i < nAmp_; i++) {
1535  // fifjEffSum is the contribution from the term involving the resonance
1536  // dynamics (f_i for resonance i) and the efficiency term.
1537  Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1538  // We need to normalise this contribution w.r.t. the complete dynamics in the DP.
1539  // Hence we scale by the fNorm_i factor (squared), which is calculated by the
1540  // initialise() function, when the normalisation integrals are calculated and cached.
1541  // We also include the complex amplitude squared to get the total normalisation
1542  // contribution from this resonance.
1543  DPNorm_ += Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1544  }
1545 
1546  // We now come to the cross-terms (between resonances i and j) in the normalisation.
1547  for (UInt_t i = 0; i < nAmp_; i++) {
1548  for (UInt_t j = i+1; j < nAmp_; j++) {
1549  LauComplex AmpjConj = Amp_[j].conj();
1550  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1551  // Again, fifjEffSum is the contribution from the term involving the resonance
1552  // dynamics (f_i*f_j_conjugate) and the efficiency cross term.
1553  // Also include the relative normalisation between these two resonances w.r.t. the
1554  // total DP dynamical structure (fNorm_i and fNorm_j) and the complex
1555  // amplitude squared (mag,phase) part.
1556  DPNorm_ += 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1557  }
1558  }
1559 
1560  return DPNorm_;
1561 }
1562 
1564 {
1565  // Routine to generate a signal event according to the Dalitz plot
1566  // model we have defined.
1567 
1568  nSigGenLoop_ = 0;
1569  Bool_t generatedSig(kFALSE);
1570 
1571  while (generatedSig == kFALSE && nSigGenLoop_ < iterationsMax_) {
1572 
1573  // Generates uniform DP phase-space distribution
1574  Double_t m13Sq(0.0), m23Sq(0.0);
1575  kinematics_->genFlatPhaseSpace(m13Sq, m23Sq);
1576 
1577  // If we're in a symmetrical DP then we should only generate events in one half
1578  if ( symmetricalDP_ && m13Sq > m23Sq ) {
1579  Double_t tmpSq = m13Sq;
1580  m13Sq = m23Sq;
1581  m23Sq = tmpSq;
1582  }
1583 
1584  // calculates the amplitudes and total amplitude for the given DP point
1585  this->calcLikelihoodInfo(m13Sq, m23Sq);
1586 
1587  if (integralsDone_ == kTRUE) {
1588 
1589  // No need for efficiency correction here. It is already done in
1590  // the dynamics() function (unlike the Fortran version).
1591 
1592  // Very important line to avoid bias in MC generation for the accept/reject method.
1593  // Make sure that the total amplitude squared is below some number, given
1594  // by aSqMaxSet_. If it is, then the event is valid.
1595  // Otherwise, go through another toy MC loop until we can generate the event
1596  // OK, or until we reach the maximum iteration limit.
1597 
1598  Double_t randNo = LauRandom::randomFun()->Rndm();
1599 
1600  if (randNo > ASq_/aSqMaxSet_) {
1601  nSigGenLoop_++;
1602  } else {
1603  generatedSig = kTRUE;
1604  nSigGenLoop_ = 0;
1605  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
1606  }
1607 
1608  } else {
1609  // For toy MC numerical integration only
1610  generatedSig = kTRUE;
1611  }
1612  } // while loop
1613 
1614  Bool_t sigGenOK(kTRUE);
1615  if (GenOK != this->checkToyMC(kTRUE,kFALSE)) {
1616  sigGenOK = kFALSE;
1617  }
1618 
1619  return sigGenOK;
1620 }
1621 
1622 LauAbsDPDynamics::ToyMCStatus LauIsobarDynamics::checkToyMC(Bool_t printErrorMessages, Bool_t printInfoMessages)
1623 {
1624  // Check whether we have generated the toy MC OK.
1625  ToyMCStatus ok(GenOK);
1626 
1627  if (nSigGenLoop_ >= iterationsMax_) {
1628  if (printErrorMessages) {
1629  std::cerr<<"WARNING in LauIsobarDynamics::checkToyMC : More than "<<iterationsMax_<<" iterations performed and no event accepted."<<std::endl;
1630  }
1631 
1632  if ( aSqMaxSet_ > 1.01 * aSqMaxVar_ ) {
1633  if (printErrorMessages) {
1634  std::cerr<<" : |A|^2 maximum was set to "<<aSqMaxSet_<<" but this appears to be too high."<<std::endl;
1635  std::cerr<<" : Maximum value of |A|^2 found so far = "<<aSqMaxVar_<<std::endl;
1636  std::cerr<<" : The value of |A|^2 maximum will be decreased and the generation restarted."<<std::endl;
1637  }
1638  aSqMaxSet_ = 1.01 * aSqMaxVar_;
1639  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : |A|^2 max reset to "<<aSqMaxSet_<<std::endl;
1640  ok = MaxIterError;
1641  } else {
1642  if (printErrorMessages) {
1643  std::cerr<<" : |A|^2 maximum was set to "<<aSqMaxSet_<<", which seems to be correct for the given model."<<std::endl;
1644  std::cerr<<" : However, the generation is very inefficient - please check your model."<<std::endl;
1645  std::cerr<<" : The maximum number of iterations will be increased and the generation restarted."<<std::endl;
1646  }
1647  iterationsMax_ *= 2;
1648  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : max number of iterations reset to "<<iterationsMax_<<std::endl;
1649  ok = MaxIterError;
1650  }
1651  } else if (aSqMaxVar_ > aSqMaxSet_) {
1652  if (printErrorMessages) {
1653  std::cerr<<"WARNING in LauIsobarDynamics::checkToyMC : |A|^2 maximum was set to "<<aSqMaxSet_<<" but a value exceeding this was found: "<<aSqMaxVar_<<std::endl;
1654  std::cerr<<" : Run was invalid, as any generated MC will be biased, according to the accept/reject method!"<<std::endl;
1655  std::cerr<<" : The value of |A|^2 maximum be reset to be > "<<aSqMaxVar_<<" and the generation restarted."<<std::endl;
1656  }
1657  aSqMaxSet_ = 1.01 * aSqMaxVar_;
1658  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : |A|^2 max reset to "<<aSqMaxSet_<<std::endl;
1659  ok = ASqMaxError;
1660  } else if (printInfoMessages) {
1661  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : aSqMaxSet = "<<aSqMaxSet_<<" and aSqMaxVar = "<<aSqMaxVar_<<std::endl;
1662  }
1663 
1664  return ok;
1665 }
1666 
1668 {
1675  scfFraction_ = currentEvent_->retrieveScfFraction(); // These two are necessary, even though the dynamics don't actually use scfFraction_ or jacobian_,
1676  jacobian_ = currentEvent_->retrieveJacobian(); // since this is at the heart of the caching mechanism.
1677 }
1678 
1680 {
1681  // Calculate the likelihood and associated info
1682  // for the given event using cached information
1683  evtLike_ = 0.0;
1684 
1685  // retrieve the cached dynamics from the tree:
1686  // realAmp, imagAmp for each resonance plus efficiency, scf fraction and jacobian
1687  this->setDataEventNo(iEvt);
1688 
1689  // use realAmp and imagAmp to create the resonance amplitudes
1690  for (UInt_t i = 0; i < nAmp_; i++) {
1692  }
1693 
1694  // Update the dynamics - calculates totAmp_ and then ASq_ = totAmp_.abs2() * eff_
1695  // All calculated using cached information on the individual amplitudes and efficiency.
1696  this->dynamics(kTRUE, 1.0);
1697 
1698  // Calculate the normalised matrix element squared value
1699  if (DPNorm_ > 1e-10) {
1700  evtLike_ = ASq_/DPNorm_;
1701  }
1702 }
1703 
1704 void LauIsobarDynamics::calcLikelihoodInfo(Double_t m13Sq, Double_t m23Sq)
1705 {
1706  this->calcLikelihoodInfo(m13Sq, m23Sq, -1);
1707 }
1708 
1709 void LauIsobarDynamics::calcLikelihoodInfo(Double_t m13Sq, Double_t m23Sq, Int_t tagCat)
1710 {
1711  // Calculate the likelihood and associated info
1712  // for the given point in the Dalitz plot
1713  // Also retrieves the SCF fraction in the bin where the event lies (done
1714  // here to cache it along with the the rest of the DP quantities, like eff)
1715  // The jacobian for the square DP is calculated here for the same reason.
1716  evtLike_ = 0.0;
1717 
1718  // update the kinematics for the specified DP point
1719  kinematics_->updateKinematics(m13Sq, m23Sq);
1720 
1721  // calculate the jacobian and the scfFraction to cache them later
1722  scfFraction_ = this->retrieveScfFraction(tagCat);
1723  if (kinematics_->squareDP() == kTRUE) {
1724  // If cacheResData == kFALSE, updateKinematics has been called before dynamics(). Then get Jacobian.
1726  }
1727 
1728  // calculates the ff_ terms and retrieves eff_ from the efficiency model
1729  // then calculates totAmp_ and finally ASq_ = totAmp_.abs2() * eff_
1730  this->dynamics(kFALSE, 1.0);
1731 
1732  // Calculate the normalised matrix element squared value
1733  if (DPNorm_ > 1e-10) {
1734  evtLike_ = ASq_/DPNorm_;
1735  }
1736 }
1737 
1739 {
1740  // In LauFitDataTree, the first two variables should always be m13^2 and m23^2.
1741  // Other variables follow thus: charge/flavour tag prob, etc.
1742 
1743  UInt_t nBranches = inputFitTree.nBranches();
1744 
1745  if (nBranches < 2) {
1746  std::cout<<"ERROR in LauIsobarDynamics::fillDataTree : Expecting at least 2 variables "
1747  <<"in input data tree, but there are "<<nBranches<<"! Make sure you have "
1748  <<"the right number of variables in your input data file!"<<std::endl;
1749  gSystem->Exit(-1);
1750  }
1751 
1752  // Data structure that will cache the variables required to
1753  // calculate the signal likelihood for this experiment
1754  for ( std::vector<LauCacheData*>::iterator iter = data_.begin(); iter != data_.end(); ++iter ) {
1755  delete (*iter);
1756  }
1757  data_.clear();
1758 
1759  Double_t m13Sq(0.0), m23Sq(0.0);
1760  Double_t mPrime(0.0), thPrime(0.0);
1761  Int_t tagCat(-1);
1762  std::vector<Double_t> realAmp(nAmp_), imagAmp(nAmp_);
1763  Double_t eff(0.0), scfFraction(0.0), jacobian(0.0);
1764 
1765  UInt_t nEvents = inputFitTree.nEvents() + inputFitTree.nFakeEvents();
1766 
1767  data_.reserve(nEvents);
1768 
1769  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1770 
1771  const LauFitData& dataValues = inputFitTree.getData(iEvt);
1772  LauFitData::const_iterator iter = dataValues.find("m13Sq");
1773  m13Sq = iter->second;
1774  iter = dataValues.find("m23Sq");
1775  m23Sq = iter->second;
1776 
1777  // is there more than one tagging category?
1778  // if so then we need to know the category from the data
1779  if (scfFractionModel_.size()>1) {
1780  iter = dataValues.find("tagCat");
1781  tagCat = static_cast<Int_t>(iter->second);
1782  }
1783 
1784  // calculates the amplitudes and total amplitude for the given DP point
1785  // tagging category not needed by dynamics, but to find out the scfFraction
1786  this->calcLikelihoodInfo(m13Sq, m23Sq, tagCat);
1787 
1788  // extract the real and imaginary parts of the ff_ terms for storage
1789  for (UInt_t i = 0; i < nAmp_; i++) {
1790  realAmp[i] = ff_[i].re();
1791  imagAmp[i] = ff_[i].im();
1792  }
1793 
1794  if ( kinematics_->squareDP() ) {
1795  mPrime = kinematics_->getmPrime();
1796  thPrime = kinematics_->getThetaPrime();
1797  }
1798 
1799  eff = this->getEvtEff();
1800  scfFraction = this->getEvtScfFraction();
1801  jacobian = this->getEvtJacobian();
1802 
1803  // store the data for each event in the list
1804  data_.push_back( new LauCacheData() );
1805  data_[iEvt]->storem13Sq(m13Sq);
1806  data_[iEvt]->storem23Sq(m23Sq);
1807  data_[iEvt]->storemPrime(mPrime);
1808  data_[iEvt]->storethPrime(thPrime);
1809  data_[iEvt]->storeEff(eff);
1810  data_[iEvt]->storeScfFraction(scfFraction);
1811  data_[iEvt]->storeJacobian(jacobian);
1812  data_[iEvt]->storeRealAmp(realAmp);
1813  data_[iEvt]->storeImagAmp(imagAmp);
1814  }
1815 }
1816 
1818 {
1819  // Select the event (kinematics_) using an accept/reject method based on the
1820  // ratio of the current value of ASq to the maximal value.
1821  Bool_t accepted(kFALSE);
1822 
1823  this->dynamics(kFALSE, 1.0, kFALSE);
1824 
1825  // Compare the ASq value with the maximal value (set by the user)
1826  if (LauRandom::randomFun()->Rndm() < ASq_/aSqMaxSet_) {
1827  accepted = kTRUE;
1828  }
1829 
1830  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
1831 
1832  return accepted;
1833 }
1834 
1836 {
1837  // calculate the dynamics from the current kinematics
1838  this->dynamics(kFALSE, 1.0, kFALSE);
1839 
1840  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
1841 
1842  // return the event weight = the value of the squared amplitude divided
1843  // by the user-defined ceiling
1844  return ASq_ / aSqMaxSet_;
1845 }
1846 
KMStringMap kMatrixPropSet_
The names of the M-matrix components in the model mapped to their propagators.
virtual Double_t getEvtJacobian() const
Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event...
Int_t nSigGenLoop_
The number of unsucessful attempts to generate an event so far.
TString getNameDaug1() const
Get name of the first daughter particle.
Double_t getm23Max() const
Get the m23 maximum defined as (mParent - m1)
Class for defining the abstract interface for signal Dalitz plot dynamics.
Double_t retrievemPrime() const
Retrieve the square Dalitz plot coordinate, m&#39;.
Bool_t squareDP() const
Are the square Dalitz plot co-ordinates being calculated?
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
File containing declaration of LauKMatrixPropFactory class.
File containing declaration of LauFitDataTree class.
File containing declaration of LauNRAmplitude class.
Double_t retrievem23Sq() const
Retrieve the invariant mass squared of the second and third daugthers.
Double_t resBarrierRadius_
The radius of the resonance barrier factor for new amplitude components.
virtual LauComplex resAmp(Int_t index)
Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz ...
virtual ~LauIsobarDynamics()
Destructor.
Int_t getTypeDaug1() const
Get PDG code of the first daughter particle.
const TString & getResonanceName() const
Get the name of the resonance.
void changeResonance(const Double_t newMass, const Double_t newWidth, const Int_t newSpin)
Allow the mass, width and spin of the resonance to be changed.
Double_t m13BinWidth_
The bin width to use when integrating over m13.
Bool_t flipHelicity_
The helicity flip flag for new amplitude components.
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.
File containing declaration of LauAbsEffModel class.
ClassImp(LauAbsCoeffSet)
LauParArray fitFracEffUnCorr_
The efficiency-uncorrected fit fractions for the amplitude components.
std::vector< Int_t > resPairAmp_
The index of the daughter not produced by the resonance for each amplitude component.
LauParameter DPRate_
The overall Dalitz plot rate.
File containing declaration of LauKMatrixProdPole class.
const TString & name() const
The parameter name.
Double_t getmPrime() const
Get m&#39; value.
Int_t getNPoles() const
Get the number of poles.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:33
const std::vector< Double_t > & retrieveImagAmp() const
Retrieve the imaginary parts of the amplitudes.
Int_t getTypeParent() const
Get PDG code of the parent particle.
Int_t getCharge(Int_t resPairAmpInt) const
Get charge of a particular two-daughter combination.
File containing declaration of LauBelleNR class.
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.
File containing declaration of LauKMatrixProdSVP class.
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.
File containing declaration of LauDaughters class.
Double_t aSqMaxVar_
The maximum value of A squared that has been seen so far while generating.
Double_t scfFraction_
The fraction of events that are poorly reconstructed (the self cross feed fraction) at the current po...
LauAbsResonance * getResonance(const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType)
Create a resonance.
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.
File containing declaration of LauGounarisSakuraiRes class.
LauAbsResonance::BarrierType barrierType_
The type of the barrier factor for new amplitude components.
Double_t aSqMaxSet_
The maximum allowed value of A squared.
File containing declaration of LauBelleSymNR class.
virtual Double_t retrieveEfficiency()
Obtain the efficiency of the current event from the model.
Bool_t symmetricalDP_
Whether the Dalitz plot is symmetrical.
File containing declaration of LauPrint class.
Pure abstract base class for defining the efficiency description across the signal Dalitz plot...
Int_t getNChannels() const
Get the number of channels.
Double_t getm13Min() const
Get the m13 minimum defined as (m1 + m3)
Double_t parBarrierRadius_
The radius of the parent barrier factor for new amplitude components.
Double_t retrieveJacobian() const
Retrieve the Jacobian for the transformation into square-Dalitz-plot coordinates. ...
Class to contain cached data relating to an event.
Definition: LauCacheData.hh:30
Class for performing numerical integration routines.
Definition: LauIntegrals.hh:30
std::map< TString, Double_t > LauFitData
Type for holding event data.
virtual void calcDPNormalisation()
Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
TString getNameDaug2() const
Get name of the second daughter particle.
UInt_t nBranches() const
Obtain the number of branches in the tree.
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 Dal...
std::vector< LauComplex > Amp_
The complex coefficients for the amplitude components.
LauKMatrixPropagator * getPropagator(const TString &name, const TString &paramFileName, Int_t resPairAmpInt, Int_t nChannels, Int_t nPoles, Int_t rowIndex)
Retrieve the propagator if it already exists, otherwise create one.
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...
Int_t resTypeInt(const TString &name) const
Retrieve the integer index for the specified resonance.
Bool_t flipHelicity() const
Get the helicity flip flag.
File containing declaration of LauKinematics class.
File containing declaration of LauRelBreitWignerRes class.
virtual Double_t getEvtEff() const
Retrieve the efficiency for the current event.
File containing declaration of LauKMatrixPropagator class.
Double_t getm13Max() const
Get the m13 maximum defined as (mParent - m2)
virtual LauAbsResonance * addResonance(const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType)
Add a resonance to the Dalitz plot.
std::vector< LauCacheData * > data_
The cached data for all events.
Class for defining a K-matrix production pole amplitude term.
Double_t getm12Min() const
Get the m12 minimum defined as (m1 + m2)
Double_t retrievem13Sq() const
Retrieve the invariant mass squared of the first and third daugthers.
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.
File containing declaration of LauIsobarDynamics class.
std::vector< std::vector< LauComplex > > fifjEffSum_
The event-by-event running total of efficiency corrected amplitude cross terms for each pair of ampli...
void updateKinematics(Double_t m13Sq, Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
virtual void initSummary()
Print a summary of the model to be used.
virtual void fillDataTree(const LauFitDataTree &fitDataTree)
Fill the internal data structure that caches the resonance dynamics.
virtual LauAbsResonance * findResonance(const TString &name)
Retrieve the named resonance.
std::map< Int_t, LauAbsEffModel * > LauTagCatScfFractionModelMap
The type used for containing multiple self cross feed fraction models for different categories (e...
Double_t ASq_
The value of A squared for the current event.
Double_t abs2() const
Obtain the square of the absolute value of the complex number.
Definition: LauComplex.hh:229
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric.
Definition: LauDaughters.hh:66
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.
std::vector< LauParameter > extraParameters_
any extra parameters/quantities (e.g. K-matrix total fit fractions)
static LauKMatrixPropFactory * getInstance()
Get a static instance of this factory class. Only one is created per application. ...
File containing declaration of LauResonanceMaker class.
Double_t m23BinWidth_
The bin width to use when integrating over m23.
void calcGaussLegendreWeights(const Int_t numPoints, std::vector< Double_t > &abscissas, std::vector< Double_t > &weights)
Calculate the Gauss-Legendre weights.
Definition: LauIntegrals.cc:37
TString intFileName_
The name of the file to save integrals to.
Int_t getResPairAmpInt() const
Get the DP axis identifier.
File containing declaration of LauCacheData class.
UInt_t nAmp_
The number of amplitude components.
virtual LauResonanceModel getResonanceModel() const =0
Get the resonance model type.
Double_t getm12Max() const
Get the m12 maximum defined as (mParent - m3)
virtual void addKMatrixProdPole(const TString &poleName, const TString &propName, Int_t poleIndex)
Add a K-matrix production pole term to the model.
std::vector< TString > resTypAmp_
The resonance types of all of the amplitude components.
virtual void calcLikelihoodInfo(UInt_t iEvt)
Calculate the likelihood (and all associated information) for the given event number.
const std::vector< Double_t > & retrieveRealAmp() const
Retrieve the real parts of the amplitudes.
std::vector< Int_t > resIntAmp_
The index within the resonance maker for each amplitude component.
LauIsobarDynamics(LauDaughters *daughters, LauAbsEffModel *effModel, LauAbsEffModel *scfFractionModel=0)
Constructor.
LauParArray fitFrac_
The fit fractions for the amplitude components.
virtual 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...
void flipAndUpdateKinematics()
Flips the DP variables m13^2 &lt;-&gt; m23^2 and recalculates all kinematic quantities. ...
std::vector< std::vector< LauComplex > > fifjSum_
The event-by-event running total of the amplitude cross terms for each pair of amplitude components...
Double_t getm23Min() const
Get the m23 minimum defined as (m2 + m3)
UInt_t getnDefinedResonances() const
Retrieve the number of defined resonances in the resonance maker.
Double_t retrieveScfFraction() const
Retrieve the fraction of poorly reconstructed events (the so-called self cross feed fraction) ...
Double_t retrievethPrime() const
Retrieve the square Dalitz plot coordinate, theta&#39;.
virtual Double_t retrieveScfFraction(Int_t tagCat)
Obtain the self cross feed fraction of the current event from the model.
Bool_t setBarrierRadius_
Should the radii of the resonance barrier factors be adjusted for new amplitude components.
virtual Double_t calcSigDPNorm()
Calculate the normalisation factor for the log-likelihood function.
virtual Double_t getEvtScfFraction() const
Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the ...
Bool_t withinDPLimits(Double_t m13Sq, Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot...
Int_t getTypeDaug3() const
Get PDG code of the third daughter particle.
LauDaughters * daughters_
The daughters of the decay.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
File containing LauRandom namespace.
virtual Bool_t generate()
Generate a toy MC signal event.
virtual 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.
Double_t eff_
The efficiency at the current point in the Dalitz plot.
Int_t getTypeDaug2() const
Get PDG code of the second daughter particle.
UInt_t nFakeEvents() const
Retrieve the number of fake events.
virtual 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...
TString getNameParent() const
Get name of the parent particle.
virtual void updateCoeffs(const std::vector< LauComplex > &coeffs)
Update the complex coefficients for the resonances.
virtual void addKMatrixProdSVP(const TString &SVPName, const TString &propName, Int_t channelIndex)
Add a K-matrix slowly-varying part (SVP) term to the model.
void rescale(Double_t scaleVal)
Scale this by a factor.
Definition: LauComplex.hh:282
Int_t getChargeParent() const
Get charge of the parent particle.
void zero()
Set both real and imaginary part to zero.
Definition: LauComplex.hh:318
Double_t DPNorm_
The normalisation factor for the log-likelihood function.
Double_t initValue() const
The initial value of the parameter.
Class for defining signal dynamics using the isobar model.
virtual void setDataEventNo(UInt_t iEvt)
Load the data for a given event.
virtual Double_t getEventWeight()
Calculate the acceptance rate, for events with the current kinematics, when generating events accordi...
File containing LauConstants namespace.
ToyMCStatus
The possible statuses for toy MC generation.
virtual void writeIntegralsFile()
Write the results of the integrals (and related information) to a file.
virtual ToyMCStatus checkToyMC(Bool_t printErrorMessages=kTRUE, Bool_t printInfoMessages=kFALSE)
Check the status of the toy MC generation.
File containing declaration of LauAbsResonance class.
virtual void setIntegralBinWidths(Double_t m13BinWidth, Double_t m23BinWidth)
Set the widths of the bins to use when integrating across the Dalitz plot.
LauResonanceMaker * resonanceMaker_
Object to create resonances.
LauKinematics * kinematics_
The kinematics of the decay.
Double_t getThetaPrime() const
Get theta&#39; value.
Class for defining a complex number.
Definition: LauComplex.hh:47
Bool_t gotKMatrixMatch(UInt_t resAmpInt, const TString &propName) const
Check whether a resonance is a K-matrix component of a given propagator.
void genFlatPhaseSpace(Double_t &m13Sq, Double_t &m23Sq) const
Routine to generate events flat in phase-space.
File containing declaration of LauIntegrals class.
TString getName() const
Get the propagator name.
virtual void setDataEventNo(UInt_t iEvt)
Load the data for a given event.
LauCacheData * currentEvent_
The cached data for the current event.
Double_t value() const
The value of the parameter.
File containing declaration of LauPolNR class.
UInt_t nEvents() const
Retrieve the number of events.
LauParameter meanDPEff_
The mean efficiency across the Dalitz plot.
Double_t calcSqDPJacobian()
Calculate the Jacobian for the transformation m23^2, m13^2 -&gt; m&#39;, theta&#39; (square DP) ...
Bool_t integralsDone_
Whether the integrals have been performed.
virtual LauComplex amplitude(const LauKinematics *kinematics)
Calculate the complex amplitude.
std::vector< LauComplex > ff_
The dynamic part of the amplitude for each amplitude component at the current point in the Dalitz plo...
virtual void initialise(const std::vector< LauComplex > &coeffs)
Initialise the Dalitz plot dynamics.
TString getNameDaug3() const
Get name of the third daughter particle.
Double_t genValue() const
The value generated for the parameter.
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...
Class for defining a K-matrix production &quot;slowly-varying part&quot; (SVP) amplitude.
Class to store the input fit variables.
Double_t retrieveEff() const
Retrieve the efficiency.
std::vector< LauAbsResonance * > sigResonances_
The resonances in the model.
virtual void calcExtraInfo(Bool_t init=kFALSE)
Calculate the fit fractions, mean efficiency and total DP rate.
LauTagCatScfFractionModelMap scfFractionModel_
The self cross feed fraction models across the Dalitz plot.
Class for defining a K-matrix propagator.
virtual void setBarrierRadii(const Double_t resRadius, const Double_t parRadius, const BarrierType type)
Set the form factor model and parameters.