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