laura is hosted by Hepforge, IPPP Durham
Laura++  v2r1p1
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  if ( res1Max > res2Min ) {
831  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m13, both close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
832  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
833  m13BinWidth = TMath::Min(res1Width,res2Width)/100.0;
834  this->calcDPPartialIntegral(minm13, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
835  } else {
836  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m13, one close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
837  this->calcDPPartialIntegral(res1Max, res2Min, minm23, maxm23, m13BinWidth, m23BinWidth);
838  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
839  m13BinWidth = res1Width/100.0;
840  this->calcDPPartialIntegral(minm13, res1Max, minm23, maxm23, m13BinWidth, m23BinWidth);
841  m13BinWidth = res2Width/100.0;
842  this->calcDPPartialIntegral(res2Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
843  }
844  } else {
845  if ( res1Max > res2Min ) {
846  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found close together in m13, dividing Dalitz plot into three regions..."<<std::endl;
847  this->calcDPPartialIntegral(minm13, res1Min, minm23, maxm23, m13BinWidth, m23BinWidth);
848  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
849  m13BinWidth = TMath::Min(res1Width,res2Width)/100.0;
850  this->calcDPPartialIntegral(res1Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
851  } else {
852  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m13, dividing Dalitz plot into five regions..."<<std::endl;
853  this->calcDPPartialIntegral(minm13, res1Min, minm23, maxm23, m13BinWidth, m23BinWidth);
854  this->calcDPPartialIntegral(res1Max, res2Min, minm23, maxm23, m13BinWidth, m23BinWidth);
855  this->calcDPPartialIntegral(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
856  m13BinWidth = res1Width/100.0;
857  this->calcDPPartialIntegral(res1Min, res1Max, minm23, maxm23, m13BinWidth, m23BinWidth);
858  m13BinWidth = res2Width/100.0;
859  this->calcDPPartialIntegral(res2Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth);
860  }
861  }
862  } else if ( s23==1 && e13 ) {
863  // We have a single narrow resonance in m23
864  // Divide the plot into 3 regions: the resonance band and
865  // the two areas either side.
866  Double_t mass = m23NarrowRes.begin()->second;
867  Double_t width = m23NarrowRes.begin()->first;
868  Double_t resMin = mass - 5.0*width;
869  Double_t resMax = mass + 5.0*width;
870  // if the resonance is close to threshold just go from
871  // threshold to resMax, otherwise treat threshold to resMin
872  // as a separate region
873  if ( resMin < (minm23+50.0*m23BinWidth_) ) {
874  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m23, close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
875  this->calcDPPartialIntegral(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth);
876  m23BinWidth = w23;
877  this->calcDPPartialIntegral(minm13, maxm13, minm23, resMax, m13BinWidth, m23BinWidth);
878  } else {
879  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : One narrow resonance found in m23, dividing Dalitz plot into three regions..."<<std::endl;
880  this->calcDPPartialIntegral(minm13, maxm13, minm23, resMin, m13BinWidth, m23BinWidth);
881  this->calcDPPartialIntegral(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth);
882  m23BinWidth = w23;
883  this->calcDPPartialIntegral(minm13, maxm13, resMin, resMax, m13BinWidth, m23BinWidth);
884  }
885  } else if ( s23==2 && e13 ) {
886  // We have a two narrow resonances in m23
887  // Divide the plot into 5 regions: the resonance bands,
888  // the two areas either side and the area in between.
889  std::multimap<Double_t,Double_t> massordered;
890  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
891  massordered.insert( std::make_pair( iter->second, iter->first ) );
892  }
893  Double_t res1Mass = massordered.begin()->first;
894  Double_t res1Width = massordered.begin()->second;
895  Double_t res1Min = res1Mass - 5.0*res1Width;
896  Double_t res1Max = res1Mass + 5.0*res1Width;
897  Double_t res2Mass = (++(massordered.begin()))->first;
898  Double_t res2Width = (++(massordered.begin()))->second;
899  Double_t res2Min = res2Mass - 5.0*res2Width;
900  Double_t res2Max = res2Mass + 5.0*res2Width;
901  // if the resonance is close to threshold just go from
902  // threshold to resMax, otherwise treat threshold to resMin
903  // as a separate region
904  if ( res1Min < (minm23+50.0*m23BinWidth_) ) {
905  if ( res1Max > res2Min ) {
906  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m23, both close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
907  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
908  m23BinWidth = TMath::Min(res1Width,res2Width)/100.0;
909  this->calcDPPartialIntegral(minm13, maxm13, minm23, res2Max, m13BinWidth, m23BinWidth);
910  } else {
911  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m23, one close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
912  this->calcDPPartialIntegral(minm13, maxm13, res1Max, res2Min, m13BinWidth, m23BinWidth);
913  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
914  m23BinWidth = res1Width/100.0;
915  this->calcDPPartialIntegral(minm13, maxm13, minm23, res1Max, m13BinWidth, m23BinWidth);
916  m23BinWidth = res2Width/100.0;
917  this->calcDPPartialIntegral(minm13, maxm13, res2Min, res2Max, m13BinWidth, m23BinWidth);
918  }
919  } else {
920  if ( res1Max > res2Min ) {
921  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found close together in m23, dividing Dalitz plot into three regions..."<<std::endl;
922  this->calcDPPartialIntegral(minm13, maxm13, minm23, res1Min, m13BinWidth, m23BinWidth);
923  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
924  m23BinWidth = TMath::Min(res1Width,res2Width)/100.0;
925  this->calcDPPartialIntegral(minm13, maxm13, res1Min, res2Max, m13BinWidth, m23BinWidth);
926  } else {
927  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Two narrow resonances found in m23, dividing Dalitz plot into five regions..."<<std::endl;
928  this->calcDPPartialIntegral(minm13, maxm13, minm23, res1Min, m13BinWidth, m23BinWidth);
929  this->calcDPPartialIntegral(minm13, maxm13, res1Max, res2Min, m13BinWidth, m23BinWidth);
930  this->calcDPPartialIntegral(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth);
931  m23BinWidth = res1Width/100.0;
932  this->calcDPPartialIntegral(minm13, maxm13, res1Min, res1Max, m13BinWidth, m23BinWidth);
933  m23BinWidth = res2Width/100.0;
934  this->calcDPPartialIntegral(minm13, maxm13, res2Min, res2Max, m13BinWidth, m23BinWidth);
935  }
936  }
937  } else if ( s13==1 && s23==1 ) {
938  // We have a single narrow resonance in both m13 and m23
939  // Divide the plot into 9 regions: the point where the
940  // resonance bands cross, the four other parts of the bands
941  // and the four remaining areas of the DP.
942  Double_t mass23 = m23NarrowRes.begin()->second;
943  Double_t width23 = m23NarrowRes.begin()->first;
944  Double_t resMin23 = mass23 - 5.0*width23;
945  Double_t resMax23 = mass23 + 5.0*width23;
946  Double_t mass13 = m13NarrowRes.begin()->second;
947  Double_t width13 = m13NarrowRes.begin()->first;
948  Double_t resMin13 = mass13 - 5.0*width13;
949  Double_t resMax13 = mass13 + 5.0*width13;
950  // if either resonance is close to threshold just go from
951  // threshold to resMax, otherwise treat threshold to resMin
952  // as a separate region
953  if ( resMin13 < (minm13+50.0*m13BinWidth_) && resMin23 < (minm23+50.0*m23BinWidth_) ) {
954  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;
955  m13BinWidth = m13BinWidth_;
956  m23BinWidth = m23BinWidth_;
957  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
958  m13BinWidth = m13BinWidth_;
959  m23BinWidth = w23;
960  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
961  m13BinWidth = w13;
962  m23BinWidth = m23BinWidth_;
963  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
964  m13BinWidth = w13;
965  m23BinWidth = w23;
966  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
967  } else if ( resMin13 < (minm13+50.0*m13BinWidth_) ) {
968  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;
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(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth);
974  m13BinWidth = w13;
975  m23BinWidth = m23BinWidth_;
976  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth);
977  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
978  m13BinWidth = w13;
979  m23BinWidth = w23;
980  this->calcDPPartialIntegral(minm13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth);
981  } else if ( resMin23 < (minm23+50.0*m23BinWidth_) ) {
982  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;
983  this->calcDPPartialIntegral(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth);
984  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
985  m13BinWidth = m13BinWidth_;
986  m23BinWidth = w23;
987  this->calcDPPartialIntegral(minm13, resMin13, minm23, resMax23, m13BinWidth, m23BinWidth);
988  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
989  m13BinWidth = w13;
990  m23BinWidth = m23BinWidth_;
991  this->calcDPPartialIntegral(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
992  m13BinWidth = w13;
993  m23BinWidth = w23;
994  this->calcDPPartialIntegral(resMin13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
995  } else {
996  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;
997  this->calcDPPartialIntegral(minm13, resMin13, minm23, resMin23, m13BinWidth, m23BinWidth);
998  this->calcDPPartialIntegral(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth);
999  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth);
1000  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1001  m13BinWidth = m13BinWidth_;
1002  m23BinWidth = w23;
1003  this->calcDPPartialIntegral(minm13, resMin13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1004  this->calcDPPartialIntegral(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1005  m13BinWidth = w13;
1006  m23BinWidth = m23BinWidth_;
1007  this->calcDPPartialIntegral(resMin13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth);
1008  this->calcDPPartialIntegral(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1009  m13BinWidth = w13;
1010  m23BinWidth = w23;
1011  this->calcDPPartialIntegral(resMin13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1012  }
1013  } else if ( e23 && s13>1 ) {
1014  // We have multiple narrow resonances in m13 only.
1015  // Divide the plot into 2 regions: threshold to the most
1016  // massive of the narrow resonances, and the rest
1017  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m13, dividing Dalitz plot into two regions..."<<std::endl;
1018  Double_t mass = 0.0;
1019  Double_t width = 0.0;
1020  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1021  if ( mass < iter->second ) {
1022  mass = iter->second;
1023  width = iter->first;
1024  }
1025  }
1026  Double_t resMax = mass + 5.0*width;
1027  this->calcDPPartialIntegral(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth);
1028  m13BinWidth = w13;
1029  this->calcDPPartialIntegral(minm13, resMax, minm23, maxm23, m13BinWidth, m23BinWidth);
1030  } else if ( e13 && s23>1 ) {
1031  // We have multiple narrow resonances in m23 only.
1032  // Divide the plot into 2 regions: threshold to the most
1033  // massive of the narrow resonances, and the rest
1034  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in m23, dividing Dalitz plot into two regions..."<<std::endl;
1035  Double_t mass = 0.0;
1036  Double_t width = 0.0;
1037  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1038  if ( mass < iter->second ) {
1039  mass = iter->second;
1040  width = iter->first;
1041  }
1042  }
1043  Double_t resMax = mass + 5.0*width;
1044  this->calcDPPartialIntegral(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth);
1045  m23BinWidth = w23;
1046  this->calcDPPartialIntegral(minm13, maxm13, minm23, resMax, m13BinWidth, m23BinWidth);
1047  } else if ( s13==1 && s23>1 ) {
1048  // We've got a single narrow resonance in m13 and multiple
1049  // narrow resonances in m23. Divide the plot into 6 regions.
1050  Double_t mass23 = 0.0;
1051  Double_t width23 = 0.0;
1052  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1053  if ( mass23 < iter->second ) {
1054  mass23 = iter->second;
1055  width23 = iter->first;
1056  }
1057  }
1058  Double_t resMax23 = mass23 + 5.0*width23;
1059  Double_t mass13 = m13NarrowRes.begin()->second;
1060  Double_t width13 = m13NarrowRes.begin()->first;
1061  Double_t resMin13 = mass13 - 5.0*width13;
1062  Double_t resMax13 = mass13 + 5.0*width13;
1063  // if the m13 resonance is close to threshold just go from
1064  // threshold to resMax, otherwise treat threshold to resMin
1065  // as a separate region
1066  if ( resMin13 < (minm13+50.0*m13BinWidth_) ) {
1067  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;
1068  m13BinWidth = m13BinWidth_;
1069  m23BinWidth = m23BinWidth_;
1070  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1071  m13BinWidth = m13BinWidth_;
1072  m23BinWidth = w23;
1073  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
1074  m13BinWidth = w13;
1075  m23BinWidth = m23BinWidth_;
1076  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1077  m13BinWidth = w13;
1078  m23BinWidth = w23;
1079  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1080  } else {
1081  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;
1082  m13BinWidth = m13BinWidth_;
1083  m23BinWidth = m23BinWidth_;
1084  this->calcDPPartialIntegral(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1085  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1086  m13BinWidth = m13BinWidth_;
1087  m23BinWidth = w23;
1088  this->calcDPPartialIntegral(minm13, resMin13, minm23, resMax23, m13BinWidth, m23BinWidth);
1089  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
1090  m13BinWidth = w13;
1091  m23BinWidth = m23BinWidth_;
1092  this->calcDPPartialIntegral(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1093  m13BinWidth = w13;
1094  m23BinWidth = w23;
1095  this->calcDPPartialIntegral(resMin13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1096  }
1097  } else if ( s13>1 && s23==1 ) {
1098  // We've got a single narrow resonance in m23 and multiple
1099  // narrow resonances in m13. Divide the plot into 6 regions.
1100  Double_t mass13 = 0.0;
1101  Double_t width13 = 0.0;
1102  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1103  if ( mass13 < iter->second ) {
1104  mass13 = iter->second;
1105  width13 = iter->first;
1106  }
1107  }
1108  Double_t resMax13 = mass13 + 5.0*width13;
1109  Double_t mass23 = m23NarrowRes.begin()->second;
1110  Double_t width23 = m23NarrowRes.begin()->first;
1111  Double_t resMin23 = mass23 - 5.0*width23;
1112  Double_t resMax23 = mass23 + 5.0*width23;
1113  // if the m23 resonance is close to threshold just go from
1114  // threshold to resMax, otherwise treat threshold to resMin
1115  // as a separate region
1116  if ( resMin23 < (minm23+50.0*m23BinWidth_) ) {
1117  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;
1118  m13BinWidth = m13BinWidth_;
1119  m23BinWidth = m23BinWidth_;
1120  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1121  m13BinWidth = m13BinWidth_;
1122  m23BinWidth = w23;
1123  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
1124  m13BinWidth = w13;
1125  m23BinWidth = m23BinWidth_;
1126  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1127  m13BinWidth = w13;
1128  m23BinWidth = w23;
1129  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1130  } else {
1131  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;
1132  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth);
1133  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1134  m13BinWidth = m13BinWidth_;
1135  m23BinWidth = w23;
1136  this->calcDPPartialIntegral(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1137  m13BinWidth = w13;
1138  m23BinWidth = m23BinWidth_;
1139  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth);
1140  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1141  m13BinWidth = w13;
1142  m23BinWidth = w23;
1143  this->calcDPPartialIntegral(minm13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth);
1144  }
1145  } else {
1146  // We've got multiple narrow resonances in both m13 and m23.
1147  // Divide the plot into 4 regions.
1148  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisation : Multiple narrow resonances found in both m13 and m23, dividing Dalitz plot into four regions..."<<std::endl;
1149  Double_t mass13 = 0.0;
1150  Double_t width13 = 0.0;
1151  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1152  if ( mass13 < iter->second ) {
1153  mass13 = iter->second;
1154  width13 = iter->first;
1155  }
1156  }
1157  Double_t resMax13 = mass13 + 5.0*width13;
1158 
1159  Double_t mass23 = 0.0;
1160  Double_t width23 = 0.0;
1161  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1162  if ( mass23 < iter->second ) {
1163  mass23 = iter->second;
1164  width23 = iter->first;
1165  }
1166  }
1167  Double_t resMax23 = mass23 + 5.0*width23;
1168 
1169  m13BinWidth = m13BinWidth_;
1170  m23BinWidth = m23BinWidth_;
1171  this->calcDPPartialIntegral(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1172  m13BinWidth = m13BinWidth_;
1173  m23BinWidth = w23;
1174  this->calcDPPartialIntegral(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth);
1175  m13BinWidth = w13;
1176  m23BinWidth = m23BinWidth_;
1177  this->calcDPPartialIntegral(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth);
1178  m13BinWidth = w13;
1179  m23BinWidth = w23;
1180  this->calcDPPartialIntegral(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth);
1181  }
1182 }
1183 
1184 void LauIsobarDynamics::setIntegralBinWidths(Double_t m13BinWidth, Double_t m23BinWidth)
1185 {
1186  // Specify whether we're going to use Gauss-Legendre integration to calculate the normalisation
1187  // integrals, and the bin widths we require for the m13 and m23 axes. Note that the integration
1188  // is done over m13, m23 space, with the appropriate Jacobian applied, and not m13^2, m23^2 space.
1189  // The default bin widths in m13 and m23 space are 5 MeV.
1190 
1191  m13BinWidth_ = m13BinWidth;
1192  m23BinWidth_ = m23BinWidth;
1193 }
1194 
1195 void LauIsobarDynamics::calcDPPartialIntegral(Double_t minm13, Double_t maxm13, Double_t minm23, Double_t maxm23,
1196  Double_t m13BinWidth, Double_t m23BinWidth)
1197 {
1198  // Calculate the total DP area, as well as finding the normalisation terms for
1199  // the signal resonances
1200 
1201  Int_t i(0), j(0);
1202  Double_t precision(1e-6);
1203 
1204  Double_t meanm13 = 0.5*(minm13 + maxm13);
1205  Double_t rangem13 = maxm13 - minm13;
1206  Double_t halfRangem13 = 0.5*rangem13;
1207 
1208  Double_t meanm23 = 0.5*(minm23 + maxm23);
1209  Double_t rangem23 = maxm23 - minm23;
1210  Double_t halfRangem23 = 0.5*rangem23;
1211 
1212  Double_t intFactor = halfRangem13*halfRangem23;
1213 
1214  // Choose smallest of mass ranges to set number of bins, given specified bin width
1215  Int_t nm13Points = static_cast<Int_t>((rangem13/m13BinWidth));
1216  Int_t nm23Points = static_cast<Int_t>((rangem23/m23BinWidth));
1217 
1218  // Avoid integral if we have no points in either x or y space
1219  if (nm13Points == 0 || nm23Points == 0) {return;}
1220 
1221  std::cout<<"INFO in LauIsobarDynamics::calcDPPartialIntegral : nm13Points = "<<nm13Points<<", nm23Points = "<<nm23Points<<std::endl;
1222  std::cout<<" : m13BinWidth = "<<m13BinWidth<<", m23BinWidth = "<<m23BinWidth<<std::endl;
1223  std::cout<<" : Integrating over m13 = "<<minm13<<" to "<<maxm13<<", m23 = "<<minm23<<" to "<<maxm23<<std::endl;
1224 
1225  LauIntegrals DPIntegrals(precision);
1226  std::vector<Double_t> m13Weights, m23Weights;
1227  std::vector<Double_t> m13Abscissas, m23Abscissas;
1228 
1229  DPIntegrals.calcGaussLegendreWeights(nm13Points, m13Abscissas, m13Weights);
1230  DPIntegrals.calcGaussLegendreWeights(nm23Points, m23Abscissas, m23Weights);
1231 
1232  Int_t nm13Weights = static_cast<Int_t>(m13Weights.size());
1233  Int_t nm23Weights = static_cast<Int_t>(m23Weights.size());
1234 
1235  //std::cout<<" : nm13Weights = "<<nm13Weights<<", nm23Weights = "<<nm23Weights<<std::endl;
1236  // Print out abscissas and weights for the integration
1237  Double_t totm13Weight(0.0), totm23Weight(0.0);
1238  for (i = 0; i < nm13Weights; i++) {
1239  totm13Weight += m13Weights[i];
1240  }
1241  for (i = 0; i < nm23Weights; i++) {
1242  totm23Weight += m23Weights[i];
1243  }
1244  std::cout<<" : totm13Weight = "<<totm13Weight<<", totm23Weight = "<<totm23Weight<<std::endl;
1245 
1246  std::vector<Double_t> m13(nm13Weights), m23(nm23Weights);
1247  std::vector<Double_t> m13Sq(nm13Weights), m23Sq(nm23Weights);
1248 
1249  // Use same number of abscissas for x and y co-ordinates
1250  Int_t m = (nm13Weights + 1)/2;
1251  for (i = 0; i < m; i++) {
1252 
1253  Int_t ii = nm13Weights - 1 - i; // symmetric i index
1254 
1255  Double_t dm13 = halfRangem13*m13Abscissas[i];
1256  Double_t m13Val = meanm13 - dm13;
1257  m13[i] = m13Val;
1258  m13Sq[i] = m13Val*m13Val;
1259 
1260  m13Val = meanm13 + dm13;
1261  m13[ii] = m13Val;
1262  m13Sq[ii] = m13Val*m13Val;
1263 
1264  }
1265 
1266  m = (nm23Weights +1)/2;
1267  for (i = 0; i < m; i++) {
1268 
1269  Int_t ii = nm23Weights - 1 - i; // symmetric i index
1270 
1271  Double_t dm23 = halfRangem23*m23Abscissas[i];
1272  Double_t m23Val = meanm23 - dm23;
1273  m23[i] = m23Val;
1274  m23Sq[i] = m23Val*m23Val;
1275  m23Val = meanm23 + dm23;
1276  m23[ii] = m23Val;
1277  m23Sq[ii] = m23Val*m23Val;
1278  }
1279 
1280  // Now compute the integral
1281  Double_t dpArea(0.0);
1282  for (i = 0; i < nm13Weights; i++) {
1283 
1284  for (j = 0; j < nm23Weights; j++) {
1285 
1286  Double_t weight = m13Weights[i]*m23Weights[j];
1287  Double_t Jacobian = 4.0*m13[i]*m23[j];
1288  weight *= (Jacobian*intFactor);
1289 
1290  // Calculate the integral contributions for each resonance.
1291  // Only resonances within the DP area contribute.
1292  // This also calculates the total DP area as a check.
1293  Bool_t withinDP = kinematics_->withinDPLimits(m13Sq[i], m23Sq[j]);
1294  if (withinDP == kTRUE) {
1295 
1296  kinematics_->updateKinematics(m13Sq[i], m23Sq[j]);
1297  this->dynamics(kFALSE, weight);
1298  // Increment total DP area
1299  dpArea += weight;
1300 
1301  }
1302 
1303  } // j weights loop
1304  } // i weights loop
1305 
1306  // Print out DP area to check whether we have a sensible output
1307  std::cout<<" : dpArea = "<<dpArea<<std::endl;
1308 
1309 }
1310 
1311 void LauIsobarDynamics::dynamics(Bool_t cacheResData, Double_t weight, Bool_t useEff)
1312 {
1313  // Routine that calculates the Dalitz plot amplitude, incorporating
1314  // resonance dynamics (using resAmp()) and any interference between them.
1315  // Used by the fit() and sigGen() functions.
1316 
1317  UInt_t i(0), j(0);
1318 
1319  // Reset the total amplitude to zero
1320  totAmp_.zero();
1321 
1322  // Loop over the number of resonance amplitudes defined in the model
1323  // Have we already calculated this for this event (during fit?)
1324  // Or do we have a resonance that has varying pole mass/width/other factors?
1325  if (cacheResData == kFALSE) {
1326 
1327  for (i = 0; i < nAmp_; i++) {
1328 
1329  // Calculate the dynamics for this resonance, using the resAmp function.
1330  ff_[i] = resAmp(i);
1331 
1332  //std::cout<<"ff_["<<i<<"] = "<<ff_[i]<<std::endl;
1333 
1334  // If we have a symmetrical Dalitz plot, flip the m_13^2 and m_23^2
1335  // variables, recalculate the dynamics, and average both contributions.
1336  // Although, the factor of 0.5 cancels out with the fact that
1337  // the resonance appears on both sides of the Dalitz plot.
1338  if (symmetricalDP_ == kTRUE) { // was tieSg_ == 12
1340  ff_[i] += resAmp(i);
1341  // Flip the m_13^2 and m_23^2 variables back to their original values
1343  }
1344 
1345  } // Loop over amplitudes
1346 
1347  // If we haven't cached the data, then we need to find out the efficiency.
1348  eff_ = this->retrieveEfficiency();
1349 
1350  } // Already cached data?
1351 
1352  if (integralsDone_ == kTRUE) {
1353 
1354  // Loop over all signal amplitudes
1355  LauComplex ATerm;
1356  for (i = 0; i < nAmp_; i++) {
1357 
1358  // Get the partial complex amplitude - (mag, phase)*(resonance dynamics)
1359  ATerm = Amp_[i]*ff_[i];
1360  // Scale this contribution by its relative normalisation w.r.t. the whole dynamics
1361  ATerm.rescale(fNorm_[i]);
1362 
1363  // Add this partial amplitude to the sum
1364  //std::cout<<"For i = "<<i<<", ATerm = "<<ATerm<<", Amp = "<<Amp_[i]<<", ff = "<<ff_[i]<<std::endl;
1365  totAmp_ += ATerm;
1366 
1367  } // Loop over amplitudes
1368 
1369  // |Sum of partial amplitudes|^2
1370  ASq_ = totAmp_.abs2();
1371 
1372  // Apply the efficiency correction for this event.
1373  // Multiply the amplitude squared sum by the DP efficiency
1374  if ( useEff ) {
1375  ASq_ *= eff_;
1376  }
1377 
1378  } else { // integrals not done
1379 
1380  // Find the efficiency correction to be applied for this event.
1381  eff_ = this->retrieveEfficiency();
1382 
1383  Double_t effWeight = eff_*weight;
1384 
1385  // Need this for integrals for normalisation of likelihood function.
1386  LauComplex fifjEffSumTerm;
1387  LauComplex fifjSumTerm;
1388  for (i = 0; i < nAmp_; i++) {
1389 
1390  // Add the dynamical amplitude squared for this resonance.
1391  Double_t fSqVal = ff_[i].abs2();
1392  fSqSum_[i] += fSqVal*weight;
1393 
1394  for (j = i; j < nAmp_; j++) {
1395 
1396  fifjEffSumTerm = fifjSumTerm = ff_[i]*ff_[j].conj();
1397 
1398  fifjEffSumTerm.rescale(effWeight);
1399  fifjEffSum_[i][j] += fifjEffSumTerm;
1400 
1401  fifjSumTerm.rescale(weight);
1402  fifjSum_[i][j] += fifjSumTerm;
1403  }
1404  }
1405  }
1406 }
1407 
1409 {
1410  // Routine to calculate the resonance dynamics (amplitude)
1411  // using the appropriate Breit-Wigner/Form Factors.
1412  // Called by the dynamics() function.
1413 
1414  // Get the signal resonance from the stored vector
1415  LauAbsResonance* sigResonance = sigResonances_[index];
1416 
1417  if (sigResonance == 0) {
1418  std::cout<<"ERROR in LauIsobarDynamics::resAmp : Couldn't retrieve resonance with index = "<<index<<std::endl;
1419  return LauComplex(0.0, 0.0);
1420  }
1421 
1422  // Get the integer index of the resonance.
1423  Int_t resInt = resIntAmp_[index];
1424 
1425  LauComplex resAmplitude(0.0, 0.0);
1426 
1427  if (resInt < 0 || resInt >= static_cast<Int_t>(this->getnDefinedResonances())) {
1428 
1429  std::cout<<"ERROR in LauIsobarDynamics::resAmp : Probably bad resonance name."<<std::endl;
1430  resAmplitude = LauComplex(0.0, 0.0);
1431 
1432  } else {
1433 
1434  resAmplitude = sigResonance->amplitude(kinematics_);
1435 
1436  }
1437 
1438  return resAmplitude;
1439 }
1440 
1441 void LauIsobarDynamics::setFFTerm(UInt_t index, Double_t realPart, Double_t imagPart)
1442 {
1443  // Function to set the internal ff = resAmp() term.
1444  if ( index >= nAmp_ ) {
1445  std::cerr<<"ERROR in LauIsobarDynamics::setFFTerm : index = "<<index<<" is not within the range 0 and "<<nAmp_-1<<". Bailing out."<<std::endl;
1446  return;
1447  }
1448 
1449  ff_[index].setRealImagPart( realPart, imagPart );
1450 }
1451 
1453 {
1454  // This method calculates the fit fractions, mean efficiency and total DP rate
1455 
1456  Double_t fifjEffTot(0.0), fifjTot(0.0);
1457  UInt_t i, j;
1458  for (i = 0; i < nAmp_; i++) {
1459 
1460  // Calculate the diagonal terms
1461  TString name = "A"; name += i; name += "Sq_FitFrac";
1462  fitFrac_[i][i].name(name);
1463 
1464  Double_t fifjSumReal = fifjSum_[i][i].re();
1465  Double_t sumTerm = Amp_[i].abs2()*fifjSumReal*fNorm_[i]*fNorm_[i];
1466  fifjTot += sumTerm;
1467 
1468  Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1469  Double_t sumEffTerm = Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1470  fifjEffTot += sumEffTerm;
1471 
1472  fitFrac_[i][i] = sumTerm;
1473  }
1474 
1475  for (i = 0; i < nAmp_; i++) {
1476  for (j = i+1; j < nAmp_; j++) {
1477  // Calculate the cross-terms
1478  TString name = "A"; name += i; name += "A"; name += j; name += "_FitFrac";
1479  fitFrac_[i][j].name(name);
1480 
1481  LauComplex AmpjConj = Amp_[j].conj();
1482  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1483 
1484  Double_t crossTerm = 2.0*(AmpTerm*fifjSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1485  fifjTot += crossTerm;
1486 
1487  Double_t crossEffTerm = 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1488  fifjEffTot += crossEffTerm;
1489 
1490  fitFrac_[i][j] = crossTerm;
1491  }
1492  }
1493 
1494  if (TMath::Abs(fifjTot) > 1e-10) {
1495  meanDPEff_ = fifjEffTot/fifjTot;
1496  if (init) {
1499  }
1500  }
1501  DPRate_ = fifjTot;
1502  if (init) {
1505  }
1506 
1507  // Now divide the fitFraction sums by the overall integral
1508  for (i = 0; i < nAmp_; i++) {
1509  for (j = i; j < nAmp_; j++) {
1510  // Get the actual fractions by dividing by the total DP rate
1511  fitFrac_[i][j] /= fifjTot;
1512  if (init) {
1513  fitFrac_[i][j].genValue( fitFrac_[i][j].value() );
1514  fitFrac_[i][j].initValue( fitFrac_[i][j].value() );
1515  }
1516  }
1517  }
1518 
1519  // Work out total fit fraction over all K-matrix components (for each propagator)
1520  KMPropMap::iterator mapIter;
1521  Int_t propInt(0);
1522 
1523  for (mapIter = kMatrixPropagators_.begin(); mapIter != kMatrixPropagators_.end(); ++mapIter) {
1524 
1525  LauKMatrixPropagator* thePropagator = mapIter->second;
1526 
1527  TString propName = thePropagator->getName();
1528 
1529  // Now loop over all resonances and find those which are K-matrix components for this propagator
1530  Double_t kMatrixTotFitFrac(0.0);
1531 
1532  for (i = 0; i < nAmp_; i++) {
1533 
1534  Bool_t gotKMRes1 = this->gotKMatrixMatch(i, propName);
1535  if (gotKMRes1 == kFALSE) {continue;}
1536 
1537  Double_t fifjSumReal = fifjSum_[i][i].re();
1538  Double_t sumTerm = Amp_[i].abs2()*fifjSumReal*fNorm_[i]*fNorm_[i];
1539 
1540  //Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1541  //Double_t sumEffTerm = Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1542 
1543  kMatrixTotFitFrac += sumTerm;
1544 
1545  for (j = i+1; j < nAmp_; j++) {
1546 
1547  Bool_t gotKMRes2 = this->gotKMatrixMatch(j, propName);
1548  if (gotKMRes2 == kFALSE) {continue;}
1549 
1550  LauComplex AmpjConj = Amp_[j].conj();
1551  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1552 
1553  Double_t crossTerm = 2.0*(AmpTerm*fifjSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1554  //Double_t crossEffTerm = 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1555 
1556  kMatrixTotFitFrac += crossTerm;
1557 
1558  }
1559 
1560  }
1561 
1562  kMatrixTotFitFrac /= fifjTot;
1563 
1564  TString parName("KMatrixTotFF_"); parName += propInt;
1565  extraParameters_[propInt].name( parName );
1566  extraParameters_[propInt] = kMatrixTotFitFrac;
1567  if (init) {
1568  extraParameters_[propInt].genValue(kMatrixTotFitFrac);
1569  extraParameters_[propInt].initValue(kMatrixTotFitFrac);
1570  }
1571 
1572  std::cout<<"INFO in LauIsobarDynamics::calcExtraInfo : Total K-matrix fit fraction for propagator "<<propName<<" is "<<kMatrixTotFitFrac<<std::endl;
1573 
1574  ++propInt;
1575 
1576  }
1577 
1578 
1579 }
1580 
1581 Bool_t LauIsobarDynamics::gotKMatrixMatch(UInt_t resAmpInt, const TString& propName) const
1582 {
1583 
1584  Bool_t gotMatch(kFALSE);
1585 
1586  if (resAmpInt >= nAmp_) {return kFALSE;}
1587 
1588  const LauAbsResonance* theResonance = sigResonances_[resAmpInt];
1589 
1590  if (theResonance == 0) {return kFALSE;}
1591 
1592  Int_t resModelInt = theResonance->getResonanceModel();
1593 
1594  if (resModelInt == LauAbsResonance::KMatrix) {
1595 
1596  TString resName = theResonance->getResonanceName();
1597 
1598  KMStringMap::const_iterator kMPropSetIter = kMatrixPropSet_.find(resName);
1599 
1600  if (kMPropSetIter != kMatrixPropSet_.end()) {
1601  TString kmPropString = kMPropSetIter->second;
1602  if (kmPropString == propName) {gotMatch = kTRUE;}
1603  }
1604 
1605  }
1606 
1607  return gotMatch;
1608 
1609 }
1610 
1612 {
1613  // Calculate the normalisation for the log-likelihood function.
1614  DPNorm_ = 0.0;
1615 
1616  for (UInt_t i = 0; i < nAmp_; i++) {
1617  // fifjEffSum is the contribution from the term involving the resonance
1618  // dynamics (f_i for resonance i) and the efficiency term.
1619  Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1620  // We need to normalise this contribution w.r.t. the complete dynamics in the DP.
1621  // Hence we scale by the fNorm_i factor (squared), which is calculated by the
1622  // initialise() function, when the normalisation integrals are calculated and cached.
1623  // We also include the complex amplitude squared to get the total normalisation
1624  // contribution from this resonance.
1625  DPNorm_ += Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1626  }
1627 
1628  // We now come to the cross-terms (between resonances i and j) in the normalisation.
1629  for (UInt_t i = 0; i < nAmp_; i++) {
1630  for (UInt_t j = i+1; j < nAmp_; j++) {
1631  LauComplex AmpjConj = Amp_[j].conj();
1632  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1633  // Again, fifjEffSum is the contribution from the term involving the resonance
1634  // dynamics (f_i*f_j_conjugate) and the efficiency cross term.
1635  // Also include the relative normalisation between these two resonances w.r.t. the
1636  // total DP dynamical structure (fNorm_i and fNorm_j) and the complex
1637  // amplitude squared (mag,phase) part.
1638  DPNorm_ += 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1639  }
1640  }
1641 
1642  return DPNorm_;
1643 }
1644 
1646 {
1647  // Routine to generate a signal event according to the Dalitz plot
1648  // model we have defined.
1649 
1650  nSigGenLoop_ = 0;
1651  Bool_t generatedSig(kFALSE);
1652 
1653  while (generatedSig == kFALSE && nSigGenLoop_ < iterationsMax_) {
1654 
1655  // Generates uniform DP phase-space distribution
1656  Double_t m13Sq(0.0), m23Sq(0.0);
1657  kinematics_->genFlatPhaseSpace(m13Sq, m23Sq);
1658 
1659  // If we're in a symmetrical DP then we should only generate events in one half
1660  if ( symmetricalDP_ && m13Sq > m23Sq ) {
1661  Double_t tmpSq = m13Sq;
1662  m13Sq = m23Sq;
1663  m23Sq = tmpSq;
1664  }
1665 
1666  // calculates the amplitudes and total amplitude for the given DP point
1667  this->calcLikelihoodInfo(m13Sq, m23Sq);
1668 
1669  if (integralsDone_ == kTRUE) {
1670 
1671  // No need for efficiency correction here. It is already done in
1672  // the dynamics() function (unlike the Fortran version).
1673 
1674  // Very important line to avoid bias in MC generation for the accept/reject method.
1675  // Make sure that the total amplitude squared is below some number, given
1676  // by aSqMaxSet_. If it is, then the event is valid.
1677  // Otherwise, go through another toy MC loop until we can generate the event
1678  // OK, or until we reach the maximum iteration limit.
1679 
1680  Double_t randNo = LauRandom::randomFun()->Rndm();
1681 
1682  if (randNo > ASq_/aSqMaxSet_) {
1683  nSigGenLoop_++;
1684  } else {
1685  generatedSig = kTRUE;
1686  nSigGenLoop_ = 0;
1687  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
1688  }
1689 
1690  } else {
1691  // For toy MC numerical integration only
1692  generatedSig = kTRUE;
1693  }
1694  } // while loop
1695 
1696  Bool_t sigGenOK(kTRUE);
1697  if (GenOK != this->checkToyMC(kFALSE,kFALSE)) {
1698  sigGenOK = kFALSE;
1699  }
1700 
1701  return sigGenOK;
1702 }
1703 
1704 LauAbsDPDynamics::ToyMCStatus LauIsobarDynamics::checkToyMC(Bool_t printErrorMessages, Bool_t printInfoMessages)
1705 {
1706  // Check whether we have generated the toy MC OK.
1707  ToyMCStatus ok(GenOK);
1708 
1709  if (nSigGenLoop_ >= iterationsMax_) {
1710  if (printErrorMessages) {
1711  std::cerr<<"WARNING in LauIsobarDynamics::checkToyMC : More than "<<iterationsMax_<<" iterations required."<<std::endl;
1712  std::cerr<<" : Try to decrease the maximum allowed value of the total amplitude squared using the "
1713  <<"LauIsobarDynamics::setASqMaxValue(Double_t) function and re-run."<<std::endl;
1714  std::cerr<<" : This needs to be, perhaps significantly, less than "<<aSqMaxSet_<<std::endl;
1715  std::cerr<<" : Maximum value of ASq so far = "<<aSqMaxVar_<<std::endl;
1716  }
1717  aSqMaxSet_ = 1.01 * aSqMaxVar_;
1718  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : |A|^2 max reset to "<<aSqMaxSet_<<std::endl;
1719  ok = MaxIterError;
1720  } else if (aSqMaxVar_ > aSqMaxSet_) {
1721  if (printErrorMessages) {
1722  std::cerr<<"WARNING in LauIsobarDynamics::checkToyMC : aSqMaxSet_ was set to "<<aSqMaxSet_<<" but actual aSqMax was "<<aSqMaxVar_<<std::endl;
1723  std::cerr<<" : Run was invalid, as any generated MC will be biased, according to the accept/reject method!"<<std::endl;
1724  std::cerr<<" : Please set aSqMaxSet >= "<<aSqMaxVar_<<" using the LauIsobarDynamics::setASqMaxValue(Double_t) function and re-run."<<std::endl;
1725  }
1726  aSqMaxSet_ = 1.01 * aSqMaxVar_;
1727  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : |A|^2 max reset to "<<aSqMaxSet_<<std::endl;
1728  ok = ASqMaxError;
1729  } else if (printInfoMessages) {
1730  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : aSqMaxSet = "<<aSqMaxSet_<<" and aSqMaxVar = "<<aSqMaxVar_<<std::endl;
1731  }
1732 
1733  return ok;
1734 }
1735 
1737 {
1744  scfFraction_ = currentEvent_->retrieveScfFraction(); // These two are necessary, even though the dynamics don't actually use scfFraction_ or jacobian_,
1745  jacobian_ = currentEvent_->retrieveJacobian(); // since this is at the heart of the caching mechanism.
1746 }
1747 
1749 {
1750  // Calculate the likelihood and associated info
1751  // for the given event using cached information
1752  evtLike_ = 0.0;
1753 
1754  // retrieve the cached dynamics from the tree:
1755  // realAmp, imagAmp for each resonance plus efficiency, scf fraction and jacobian
1756  this->setDataEventNo(iEvt);
1757 
1758  // use realAmp and imagAmp to create the resonance amplitudes
1759  for (UInt_t i = 0; i < nAmp_; i++) {
1761  }
1762 
1763  // Update the dynamics - calculates totAmp_ and then ASq_ = totAmp_.abs2() * eff_
1764  // All calculated using cached information on the individual amplitudes and efficiency.
1765  this->dynamics(kTRUE, 1.0);
1766 
1767  // Calculate the normalised matrix element squared value
1768  if (DPNorm_ > 1e-10) {
1769  evtLike_ = ASq_/DPNorm_;
1770  }
1771 }
1772 
1773 void LauIsobarDynamics::calcLikelihoodInfo(Double_t m13Sq, Double_t m23Sq)
1774 {
1775  this->calcLikelihoodInfo(m13Sq, m23Sq, -1);
1776 }
1777 
1778 void LauIsobarDynamics::calcLikelihoodInfo(Double_t m13Sq, Double_t m23Sq, Int_t tagCat)
1779 {
1780  // Calculate the likelihood and associated info
1781  // for the given point in the Dalitz plot
1782  // Also retrieves the SCF fraction in the bin where the event lies (done
1783  // here to cache it along with the the rest of the DP quantities, like eff)
1784  // The jacobian for the square DP is calculated here for the same reason.
1785  evtLike_ = 0.0;
1786 
1787  // update the kinematics for the specified DP point
1788  kinematics_->updateKinematics(m13Sq, m23Sq);
1789 
1790  // calculate the jacobian and the scfFraction to cache them later
1791  scfFraction_ = this->retrieveScfFraction(tagCat);
1792  if (kinematics_->squareDP() == kTRUE) {
1793  // If cacheResData == kFALSE, updateKinematics has been called before dynamics(). Then get Jacobian.
1795  }
1796 
1797  // calculates the ff_ terms and retrieves eff_ from the efficiency model
1798  // then calculates totAmp_ and finally ASq_ = totAmp_.abs2() * eff_
1799  this->dynamics(kFALSE, 1.0);
1800 
1801  // Calculate the normalised matrix element squared value
1802  if (DPNorm_ > 1e-10) {
1803  evtLike_ = ASq_/DPNorm_;
1804  }
1805 }
1806 
1808 {
1809  // In LauFitDataTree, the first two variables should always be m13^2 and m23^2.
1810  // Other variables follow thus: charge/flavour tag prob, etc.
1811 
1812  UInt_t nBranches = inputFitTree.nBranches();
1813 
1814  if (nBranches < 2) {
1815  std::cout<<"ERROR in LauIsobarDynamics::fillDataTree : Expecting at least 2 variables "
1816  <<"in input data tree, but there are "<<nBranches<<"! Make sure you have "
1817  <<"the right number of variables in your input data file!"<<std::endl;
1818  gSystem->Exit(-1);
1819  }
1820 
1821  // Data structure that will cache the variables required to
1822  // calculate the signal likelihood for this experiment
1823  for ( std::vector<LauCacheData*>::iterator iter = data_.begin(); iter != data_.end(); ++iter ) {
1824  delete (*iter);
1825  }
1826  data_.clear();
1827 
1828  Double_t m13Sq(0.0), m23Sq(0.0);
1829  Double_t mPrime(0.0), thPrime(0.0);
1830  Int_t tagCat(-1);
1831  std::vector<Double_t> realAmp(nAmp_), imagAmp(nAmp_);
1832  Double_t eff(0.0), scfFraction(0.0), jacobian(0.0);
1833 
1834  UInt_t nEvents = inputFitTree.nEvents() + inputFitTree.nFakeEvents();
1835 
1836  data_.reserve(nEvents);
1837 
1838  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1839 
1840  const LauFitData& dataValues = inputFitTree.getData(iEvt);
1841  LauFitData::const_iterator iter = dataValues.find("m13Sq");
1842  m13Sq = iter->second;
1843  iter = dataValues.find("m23Sq");
1844  m23Sq = iter->second;
1845 
1846  // is there more than one tagging category?
1847  // if so then we need to know the category from the data
1848  if (scfFractionModel_.size()>1) {
1849  iter = dataValues.find("tagCat");
1850  tagCat = static_cast<Int_t>(iter->second);
1851  }
1852 
1853  // calculates the amplitudes and total amplitude for the given DP point
1854  // tagging category not needed by dynamics, but to find out the scfFraction
1855  this->calcLikelihoodInfo(m13Sq, m23Sq, tagCat);
1856 
1857  // extract the real and imaginary parts of the ff_ terms for storage
1858  for (UInt_t i = 0; i < nAmp_; i++) {
1859  realAmp[i] = ff_[i].re();
1860  imagAmp[i] = ff_[i].im();
1861  }
1862 
1863  if ( kinematics_->squareDP() ) {
1864  mPrime = kinematics_->getmPrime();
1865  thPrime = kinematics_->getThetaPrime();
1866  }
1867 
1868  eff = this->getEvtEff();
1869  scfFraction = this->getEvtScfFraction();
1870  jacobian = this->getEvtJacobian();
1871 
1872  // store the data for each event in the list
1873  data_.push_back( new LauCacheData() );
1874  data_[iEvt]->storem13Sq(m13Sq);
1875  data_[iEvt]->storem23Sq(m23Sq);
1876  data_[iEvt]->storemPrime(mPrime);
1877  data_[iEvt]->storethPrime(thPrime);
1878  data_[iEvt]->storeEff(eff);
1879  data_[iEvt]->storeScfFraction(scfFraction);
1880  data_[iEvt]->storeJacobian(jacobian);
1881  data_[iEvt]->storeRealAmp(realAmp);
1882  data_[iEvt]->storeImagAmp(imagAmp);
1883  }
1884 }
1885 
1887 {
1888  // Select the event (kinematics_) using an accept/reject method based on the
1889  // ratio of the current value of ASq to the maximal value.
1890  Bool_t accepted(kFALSE);
1891 
1892  this->dynamics(kFALSE, 1.0, kFALSE);
1893 
1894  // Compare the ASq value with the maximal value (set by the user)
1895  if (LauRandom::randomFun()->Rndm() < ASq_/aSqMaxSet_) {
1896  accepted = kTRUE;
1897  }
1898 
1899  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
1900 
1901  return accepted;
1902 }
1903 
1905 {
1906  // calculate the dynamics from the current kinematics
1907  this->dynamics(kFALSE, 1.0, kFALSE);
1908 
1909  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
1910 
1911  // return the event weight = the value of the squared amplitude divided
1912  // by the user-defined ceiling
1913  return ASq_ / aSqMaxSet_;
1914 }
1915 
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.