laura is hosted by Hepforge, IPPP Durham
Laura++  v3r1
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 - 2014.
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 #include <set>
19 #include <vector>
20 
21 #include "TFile.h"
22 #include "TRandom.h"
23 #include "TSystem.h"
24 
25 #include "LauAbsEffModel.hh"
26 #include "LauAbsResonance.hh"
27 #include "LauAbsIncohRes.hh"
28 #include "LauBelleNR.hh"
29 #include "LauBelleSymNR.hh"
30 #include "LauCacheData.hh"
31 #include "LauConstants.hh"
32 #include "LauDaughters.hh"
34 #include "LauFitDataTree.hh"
35 #include "LauIsobarDynamics.hh"
36 #include "LauKinematics.hh"
37 #include "LauKMatrixProdPole.hh"
38 #include "LauKMatrixProdSVP.hh"
39 #include "LauKMatrixPropagator.hh"
40 #include "LauKMatrixPropFactory.hh"
41 #include "LauNRAmplitude.hh"
42 #include "LauPrint.hh"
43 #include "LauRandom.hh"
44 #include "LauResonanceMaker.hh"
45 
47 
48 // for Kpipi: only one scfFraction 2D histogram is needed
50  daughters_(daughters),
51  kinematics_(daughters_ ? daughters_->getKinematics() : 0),
52  effModel_(effModel),
53  nAmp_(0),
54  nIncohAmp_(0),
55  DPNorm_(0.0),
56  DPRate_("DPRate", 0.0, 0.0, 1000.0),
57  meanDPEff_("meanDPEff", 0.0, 0.0, 1.0),
58  currentEvent_(0),
59  symmetricalDP_(kFALSE),
60  fullySymmetricDP_(kFALSE),
61  integralsDone_(kFALSE),
62  normalizationSchemeDone_(kFALSE),
63  intFileName_("integ.dat"),
64  m13BinWidth_(0.005),
65  m23BinWidth_(0.005),
66  mPrimeBinWidth_(0.001),
67  thPrimeBinWidth_(0.001),
68  narrowWidth_(0.020),
69  binningFactor_(100.0),
70  m13Sq_(0.0),
71  m23Sq_(0.0),
72  mPrime_(0.0),
73  thPrime_(0.0),
74  tagCat_(-1),
75  eff_(1.0),
76  scfFraction_(0.0),
77  jacobian_(0.0),
78  ASq_(0.0),
79  evtLike_(0.0),
80  iterationsMax_(100000),
81  nSigGenLoop_(0),
82  aSqMaxSet_(1.25),
83  aSqMaxVar_(0.0),
84  flipHelicity_(kTRUE),
85  recalcNormalisation_(kFALSE)
86 {
87  if (daughters != 0) {
88  symmetricalDP_ = daughters->gotSymmetricalDP();
89  fullySymmetricDP_ = daughters->gotFullySymmetricDP();
90  typDaug_.push_back(daughters->getTypeDaug1());
91  typDaug_.push_back(daughters->getTypeDaug2());
92  typDaug_.push_back(daughters->getTypeDaug3());
93  }
94 
95  if (scfFractionModel != 0) {
96  scfFractionModel_[0] = scfFractionModel;
97  }
98 
99  sigResonances_.clear();
100  sigIncohResonances_.clear();
101  kMatrixPropagators_.clear();
102  kMatrixPropSet_.clear();
103  extraParameters_.clear();
104 }
105 
106 // for Kspipi, we need a scfFraction 2D histogram for each tagging category. They are provided by the map.
107 // Also, we need to know the place that the tagging category of the current event occupies in the data structure inputFitTree
109  daughters_(daughters),
110  kinematics_(daughters_ ? daughters_->getKinematics() : 0),
111  effModel_(effModel),
112  scfFractionModel_(scfFractionModel),
113  nAmp_(0),
114  nIncohAmp_(0),
115  DPNorm_(0.0),
116  DPRate_("DPRate", 0.0, 0.0, 1000.0),
117  meanDPEff_("meanDPEff", 0.0, 0.0, 1.0),
118  currentEvent_(0),
119  symmetricalDP_(kFALSE),
120  fullySymmetricDP_(kFALSE),
121  integralsDone_(kFALSE),
122  normalizationSchemeDone_(kFALSE),
123  intFileName_("integ.dat"),
124  m13BinWidth_(0.005),
125  m23BinWidth_(0.005),
126  mPrimeBinWidth_(0.001),
127  thPrimeBinWidth_(0.001),
128  narrowWidth_(0.020),
129  binningFactor_(100.0),
130  m13Sq_(0.0),
131  m23Sq_(0.0),
132  mPrime_(0.0),
133  thPrime_(0.0),
134  tagCat_(-1),
135  eff_(1.0),
136  scfFraction_(0.0),
137  jacobian_(0.0),
138  ASq_(0.0),
139  evtLike_(0.0),
140  iterationsMax_(100000),
141  nSigGenLoop_(0),
142  aSqMaxSet_(1.25),
143  aSqMaxVar_(0.0),
144  flipHelicity_(kTRUE),
145  recalcNormalisation_(kFALSE)
146 {
147  // Constructor for the isobar signal model
148  if (daughters != 0) {
149  symmetricalDP_ = daughters->gotSymmetricalDP();
150  fullySymmetricDP_ = daughters->gotFullySymmetricDP();
151  typDaug_.push_back(daughters->getTypeDaug1());
152  typDaug_.push_back(daughters->getTypeDaug2());
153  typDaug_.push_back(daughters->getTypeDaug3());
154  }
155 
156  sigResonances_.clear();
157  sigIncohResonances_.clear();
158  kMatrixPropagators_.clear();
159  kMatrixPropSet_.clear();
160  extraParameters_.clear();
161 }
162 
164 {
165  extraParameters_.clear();
166 
167  for ( std::vector<LauCacheData*>::iterator iter = data_.begin(); iter != data_.end(); ++iter ) {
168  delete (*iter);
169  }
170  data_.clear();
171 
172  for (std::vector<LauDPPartialIntegralInfo*>::iterator it = dpPartialIntegralInfo_.begin(); it != dpPartialIntegralInfo_.end(); ++it)
173  {
174  delete (*it);
175  }
176  dpPartialIntegralInfo_.clear();
177 }
178 
180 {
181  for (UInt_t i = 0; i < nAmp_; i++) {
182 
183  fSqSum_[i] = 0.0;
184  fSqEffSum_[i] = 0.0;
185  fNorm_[i] = 0.0;
186  ff_[i].zero();
187 
188  for (UInt_t j = 0; j < nAmp_; j++) {
189 
190  fifjEffSum_[i][j].zero();
191  fifjSum_[i][j].zero();
192 
193  }
194  }
195 
196  for (UInt_t i = 0; i < nIncohAmp_; i++) {
197  fSqSum_[i+nAmp_] = 0.0;
198  fSqEffSum_[i+nAmp_] = 0.0;
199  fNorm_[i+nAmp_] = 0.0;
200  incohInten_[i] = 0.0;
201  }
202 }
203 
205 {
206  if ( recalcNormalisation_ == kFALSE ) {
207  return;
208  }
209 
210  // We need to calculate the normalisation constants for the
211  // Dalitz plot generation/fitting.
212 
213  integralsDone_ = kFALSE;
214 
215  this->resetNormVectors();
217  this->calcDPNormalisation();
218 
219  integralsDone_ = kTRUE;
220 }
221 
223 {
224  // Loop through the resonance parameters and see which ones have changed
225  // For those that have changed mark the corresponding resonance(s) as needing to be re-evaluated
226 
227  integralsToBeCalculated_.clear();
228 
229  const UInt_t nResPars = resonancePars_.size();
230  for ( UInt_t iPar(0); iPar < nResPars; ++iPar ) {
231  const Double_t newValue = resonancePars_[iPar]->value();
232  if ( newValue != resonanceParValues_[iPar] ) {
233  resonanceParValues_[iPar] = newValue;
234 
235  const std::vector<UInt_t>& indices = resonanceParResIndex_[iPar];
236  std::vector<UInt_t>::const_iterator indexIter = indices.begin();
237  const std::vector<UInt_t>::const_iterator indexEnd = indices.end();
238  for( ; indexIter != indexEnd; ++indexIter) {
239  integralsToBeCalculated_.insert(*indexIter);
240  }
241  }
242  }
243 }
244 
245 void LauIsobarDynamics::initialise(const std::vector<LauComplex>& coeffs)
246 {
247  // Check whether we have a valid set of integration constants for
248  // the normalisation of the signal likelihood function.
249  this->initialiseVectors();
250 
251  // Mark the DP integrals as undetermined
252  integralsDone_ = kFALSE;
253 
254  // Initialise all resonance models
255  resonancePars_.clear();
256  resonanceParValues_.clear();
257  resonanceParResIndex_.clear();
258 
259  std::set<LauParameter*> uniqueResPars;
260  UInt_t resIndex(0);
261  for ( std::vector<LauAbsResonance*>::iterator resIter = sigResonances_.begin(); resIter != sigResonances_.end(); ++resIter ) {
262  (*resIter)->initialise();
263 
264  // Check if this resonance has floating parameters
265  // Append all unique parameters to our list
266  const std::vector<LauParameter*>& resPars = (*resIter)->getFloatingParameters();
267 
268  for ( std::vector<LauParameter*>::const_iterator parIter = resPars.begin(); parIter != resPars.end(); ++parIter ) {
269  if ( uniqueResPars.insert( *parIter ).second ) {
270  // This parameter has not already been added to
271  // the list of unique ones. Add it, its value
272  // and its associated resonance ID to the
273  // appropriate lists.
274  resonancePars_.push_back( *parIter );
275  resonanceParValues_.push_back( (*parIter)->value() );
276  std::vector<UInt_t> resIndices( 1, resIndex );
277  resonanceParResIndex_.push_back( resIndices );
278  } else {
279  // This parameter has already been added to the
280  // list of unique ones. However, we still need
281  // to indicate that this resonance should be
282  // associated with it.
283  std::vector<LauParameter*>::const_iterator uniqueParIter = resonancePars_.begin();
284  std::vector<std::vector<UInt_t> >::iterator indicesIter = resonanceParResIndex_.begin();
285  while( (*uniqueParIter) != (*parIter) ) {
286  ++uniqueParIter;
287  ++indicesIter;
288  }
289  ( *indicesIter ).push_back( resIndex );
290  }
291  }
292 
293  ++resIndex;
294  }
295 
296  for ( std::vector<LauAbsIncohRes*>::iterator resIter = sigIncohResonances_.begin(); resIter != sigIncohResonances_.end(); ++resIter ) {
297  (*resIter)->initialise();
298 
299  // Check if this resonance has floating parameters
300  // Append all unique parameters to our list
301  const std::vector<LauParameter*>& resPars = (*resIter)->getFloatingParameters();
302 
303  for ( std::vector<LauParameter*>::const_iterator parIter = resPars.begin(); parIter != resPars.end(); ++parIter ) {
304  if ( uniqueResPars.insert( *parIter ).second ) {
305  // This parameter has not already been added to
306  // the list of unique ones. Add it, its value
307  // and its associated resonance ID to the
308  // appropriate lists.
309  resonancePars_.push_back( *parIter );
310  resonanceParValues_.push_back( (*parIter)->value() );
311  std::vector<UInt_t> resIndices( 1, resIndex );
312  resonanceParResIndex_.push_back( resIndices );
313  } else {
314  // This parameter has already been added to the
315  // list of unique ones. However, we still need
316  // to indicate that this resonance should be
317  // associated with it.
318  std::vector<LauParameter*>::const_iterator uniqueParIter = resonancePars_.begin();
319  std::vector<std::vector<UInt_t> >::iterator indicesIter = resonanceParResIndex_.begin();
320  while( (*uniqueParIter) != (*parIter) ) {
321  ++uniqueParIter;
322  ++indicesIter;
323  }
324  ( *indicesIter ).push_back( resIndex );
325  }
326  }
327 
328  ++resIndex;
329  }
330 
331  if ( resonancePars_.empty() ) {
332  recalcNormalisation_ = kFALSE;
333  } else {
334  recalcNormalisation_ = kTRUE;
335  }
336 
337  // Print summary of what we have so far to screen
338  this->initSummary();
339 
340  if ( nAmp_+nIncohAmp_ == 0 ) {
341  std::cout << "INFO in LauIsobarDynamics::initialise : No contributions to DP model, not performing normalisation integrals." << std::endl;
342  } else {
343 
344  // We need to calculate the normalisation constants for the Dalitz plot generation/fitting.
345 
346  std::cout<<"INFO in LauIsobarDynamics::initialise : Starting special run to generate the integrals for normalising the PDF..."<<std::endl;
347 
348  // Since this is the initialisation, we need to calculate everything for every resonance
349  integralsToBeCalculated_.clear();
350  for ( UInt_t i(0); i < nAmp_+nIncohAmp_; ++i ) {
351  integralsToBeCalculated_.insert(i);
352  }
353 
354  // Calculate and cache the normalisations of each resonance _dynamic_ amplitude
355  // (e.g. Breit-Wigner contribution, not from the complex coefficients).
356  // These are stored in fNorm_[i].
357  // fSqSum[i] is the total of the dynamical amplitude squared for a given resonance, i.
358  // We require that:
359  // |fNorm_[i]|^2 * |fSqSum[i]|^2 = 1,
360  // i.e. fNorm_[i] normalises each resonance contribution to give the same number of
361  // events in the DP, accounting for the total DP area and the dynamics of the resonance.
362  this->calcDPNormalisation();
363 
364  // Write the integrals to a file (mainly for debugging purposes)
365  this->writeIntegralsFile();
366  }
367 
368 
369  integralsDone_ = kTRUE;
370 
371  std::cout << std::setprecision(10);
372 
373  std::cout<<"INFO in LauIsobarDynamics::initialise : Summary of the integrals:"<<std::endl;
374 
375  for (UInt_t i = 0; i < nAmp_+nIncohAmp_; i++) {
376  std::cout<<" fNorm["<<i<<"] = "<<fNorm_[i]<<std::endl;
377  std::cout<<" fSqSum["<<i<<"] = "<<fSqSum_[i]<<std::endl;
378  std::cout<<" fSqEffSum["<<i<<"] = "<<fSqEffSum_[i]<<std::endl;
379  }
380 
381  for (UInt_t i = 0; i < nAmp_; i++) {
382  for (UInt_t j = 0; j < nAmp_; j++) {
383  std::cout<<" fifjEffSum["<<i<<"]["<<j<<"] = "<<fifjEffSum_[i][j];
384  }
385  std::cout<<std::endl;
386  }
387 
388  for (UInt_t i = 0; i < nAmp_; i++) {
389  for (UInt_t j = 0; j < nAmp_; j++) {
390  std::cout<<" fifjSum["<<i<<"]["<<j<<"] = "<<fifjSum_[i][j];
391  }
392  std::cout<<std::endl;
393  }
394 
395  // Calculate the initial fit fractions (for later comparison with Toy MC, if required)
396  this->updateCoeffs(coeffs);
397  this->calcExtraInfo(kTRUE);
398 
399  for (UInt_t i = 0; i < nAmp_; i++) {
400  for (UInt_t j = i; j < nAmp_; j++) {
401  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial fit fraction for amplitude ("<<i<<","<<j<<") = "<<fitFrac_[i][j].genValue()<<std::endl;
402  }
403  }
404 
405  for (UInt_t i = 0; i < nIncohAmp_; i++) {
406  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial fit fraction for incoherent amplitude ("<<i<<","<<i<<") = "<<fitFrac_[i+nAmp_][i+nAmp_].genValue()<<std::endl;
407  }
408 
409  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial efficiency = "<<meanDPEff_.initValue()<<std::endl;
410 
411  std::cout<<"INFO in LauIsobarDynamics::initialise : Initial DPRate = "<<DPRate_.initValue()<<std::endl;
412 }
413 
415 {
416 
417  UInt_t i(0);
418  TString nameP = daughters_->getNameParent();
419  TString name1 = daughters_->getNameDaug1();
420  TString name2 = daughters_->getNameDaug2();
421  TString name3 = daughters_->getNameDaug3();
422 
423  std::cout<<"INFO in LauIsobarDynamics::initSummary : We are going to do a DP with "<<nameP<<" going to "<<name1<<" "<<name2<<" "<<name3<<std::endl;
424 
425  std::cout<<" : For the following resonance combinations:"<<std::endl;
426 
427  std::cout<<" : In tracks 2 and 3:"<<std::endl;
428  for (i = 0; i < nAmp_; i++) {
429  if (resPairAmp_[i] == 1) {
430  std::cout<<" : A"<<i<<": "<<resTypAmp_[i]<<" to "<<name2<<", "<< name3<<std::endl;
431  }
432  }
433 
434  for (i = 0; i < nIncohAmp_; i++) {
435  if (incohResPairAmp_[i] == 1) {
436  std::cout<<" : A"<<nAmp_+i<<": "<<incohResTypAmp_[i]<<" (incoherent) to "<<name2<<", "<< name3<<std::endl;
437  }
438  }
439 
440  std::cout<<" : In tracks 1 and 3:"<<std::endl;
441  for (i = 0; i < nAmp_; i++) {
442  if (resPairAmp_[i] == 2) {
443  std::cout<<" : A"<<i<<": "<<resTypAmp_[i]<<" to "<<name1<<", "<< name3<<std::endl;
444  }
445  }
446 
447  for (i = 0; i < nIncohAmp_; i++) {
448  if (incohResPairAmp_[i] == 2) {
449  std::cout<<" : A"<<nAmp_+i<<": "<<incohResTypAmp_[i]<<" (incoherent) to "<<name1<<", "<< name3<<std::endl;
450  }
451  }
452 
453  std::cout<<" : In tracks 1 and 2:"<<std::endl;
454  for (i = 0; i < nAmp_; i++) {
455  if (resPairAmp_[i] == 3) {
456  std::cout<<" : A"<<i<<": "<<resTypAmp_[i]<<" to "<<name1<<", "<< name2<<std::endl;
457  }
458  }
459 
460  for (i = 0; i < nIncohAmp_; i++) {
461  if (incohResPairAmp_[i] == 3) {
462  std::cout<<" : A"<<nAmp_+i<<": "<<incohResTypAmp_[i]<<" (incoherent) to "<<name1<<", "<< name2<<std::endl;
463  }
464  }
465 
466  for (i = 0; i < nAmp_; i++) {
467  if (resPairAmp_[i] == 0) {
468  std::cout<<" : and a non-resonant amplitude:"<<std::endl;
469  std::cout<<" : A"<<i<<std::endl;
470  }
471  }
472 
473  std::cout<<std::endl;
474 }
475 
477 {
478  fSqSum_.clear(); fSqSum_.resize(nAmp_+nIncohAmp_);
479  fSqEffSum_.clear(); fSqEffSum_.resize(nAmp_+nIncohAmp_);
480  fifjEffSum_.clear(); fifjEffSum_.resize(nAmp_);
481  fifjSum_.clear(); fifjSum_.resize(nAmp_);
482  ff_.clear(); ff_.resize(nAmp_);
483  incohInten_.clear(); incohInten_.resize(nIncohAmp_);
484  fNorm_.clear(); fNorm_.resize(nAmp_+nIncohAmp_);
485  fitFrac_.clear(); fitFrac_.resize(nAmp_+nIncohAmp_);
487 
488  LauComplex null(0.0, 0.0);
489 
490  for (UInt_t i = 0; i < nAmp_; i++) {
491 
492  fSqSum_[i] = 0.0; fNorm_[i] = 0.0;
493  ff_[i] = null;
494  fifjEffSum_[i].resize(nAmp_);
495  fifjSum_[i].resize(nAmp_);
496  fitFrac_[i].resize(nAmp_+nIncohAmp_);
497  fitFracEffUnCorr_[i].resize(nAmp_+nIncohAmp_);
498 
499  for (UInt_t j = 0; j < nAmp_; j++) {
500  fifjEffSum_[i][j] = null;
501  fifjSum_[i][j] = null;
502  fitFrac_[i][j].valueAndRange(0.0, -100.0, 100.0);
503  fitFracEffUnCorr_[i][j].valueAndRange(0.0, -100.0, 100.0);
504  }
505  for (UInt_t j = nAmp_; j < nAmp_+nIncohAmp_; j++) {
506  fitFrac_[i][j].valueAndRange(0.0, -100.0, 100.0);
507  fitFracEffUnCorr_[i][j].valueAndRange(0.0, -100.0, 100.0);
508  }
509  }
510  for (UInt_t i = nAmp_; i < nAmp_+nIncohAmp_; i++) {
511  fSqSum_[i] = 0.0; fNorm_[i] = 0.0;
512  incohInten_[i-nAmp_] = 0.0;
513  fitFrac_[i].resize(nAmp_+nIncohAmp_);
514  fitFracEffUnCorr_[i].resize(nAmp_+nIncohAmp_);
515 
516  for (UInt_t j = 0; j < nAmp_+nIncohAmp_; j++) {
517  fitFrac_[i][j].valueAndRange(0.0, -100.0, 100.0);
518  fitFracEffUnCorr_[i][j].valueAndRange(0.0, -100.0, 100.0);
519  }
520  }
521 
522  UInt_t nKMatrixProps = kMatrixPropagators_.size();
523  extraParameters_.clear();
524  extraParameters_.resize( nKMatrixProps );
525 
526  for ( UInt_t i = 0; i < nKMatrixProps; ++i ) {
527  extraParameters_[i].valueAndRange(0.0, -100.0, 100.0);
528  }
529 }
530 
531 
533 {
534  // Routine to end integration calculation for loglike normalisation.
535  // This writes out the normalisation integral output into the file named
536  // outputFileName.
537  std::cout<<"INFO in LauIsobarDynamics::writeIntegralsFile : Writing integral output to integrals file "<<intFileName_.Data()<<std::endl;
538 
539  UInt_t i(0), j(0);
540  std::ofstream getChar(intFileName_.Data());
541 
542  getChar << std::setprecision(10);
543 
544  // Write out daughter types (pi, pi0, K, K0s?)
545  for (i = 0; i < 3; i++) {
546  getChar << typDaug_[i] << " ";
547  }
548 
549  // Write out number of resonances in the Dalitz plot model
550  getChar << nAmp_ << std::endl;
551 
552  // Write out the resonances
553  for (i = 0; i < nAmp_; i++) {
554  getChar << resTypAmp_[i] << " ";
555  }
556 
557  getChar << std::endl;
558 
559  // Write out the resonance model types (BW, RelBW etc...)
560  for (i = 0; i < nAmp_; i++) {
561  LauAbsResonance* theResonance = sigResonances_[i];
562  Int_t resModelInt = theResonance->getResonanceModel();
563  getChar << resModelInt << " ";
564  }
565 
566  getChar << std::endl;
567 
568  // Write out the track pairings for each resonance. This is specified
569  // by the resPairAmpInt integer in the addResonance function.
570  for (i = 0; i < nAmp_; i++) {
571  getChar << resPairAmp_[i] << " ";
572  }
573 
574  getChar << std::endl;
575 
576  // Write out the fSqSum = |ff|^2, where ff = resAmp()
577  for (i = 0; i < nAmp_; i++) {
578  getChar << fSqSum_[i] << " ";
579  }
580 
581  getChar << std::endl;
582 
583  // Similar to fSqSum, but with the efficiency term included.
584  for (i = 0; i < nAmp_; i++) {
585  getChar << fSqEffSum_[i] << " ";
586  }
587 
588  getChar << std::endl;
589 
590  // Write out the f_i*f_j_conj*eff values = resAmp_i*resAmp_j_conj*eff.
591  // Note that only the top half of the i*j "matrix" is required, as it
592  // is symmetric w.r.t i, j.
593  for (i = 0; i < nAmp_; i++) {
594  for (j = i; j < nAmp_; j++) {
595  getChar << fifjEffSum_[i][j] << " ";
596  }
597  }
598 
599  getChar << std::endl;
600 
601  // Similar to fifjEffSum, but without the efficiency term included.
602  for (i = 0; i < nAmp_; i++) {
603  for (j = i; j < nAmp_; j++) {
604  getChar << fifjSum_[i][j] << " ";
605  }
606  }
607 
608  getChar << std::endl;
609 
610  // Write out number of incoherent resonances in the Dalitz plot model
611  getChar << nIncohAmp_ << std::endl;
612 
613  // Write out the incoherent resonances
614  for (i = 0; i < nIncohAmp_; i++) {
615  getChar << incohResTypAmp_[i] << " ";
616  }
617 
618  getChar << std::endl;
619 
620  // Write out the incoherent resonance model types (BW, RelBW etc...)
621  for (i = 0; i < nIncohAmp_; i++) {
622  LauAbsResonance* theResonance = sigIncohResonances_[i];
623  Int_t resModelInt = theResonance->getResonanceModel();
624  getChar << resModelInt << " ";
625  }
626 
627  getChar << std::endl;
628 
629  // Write out the track pairings for each incoherent resonance. This is specified
630  // by the resPairAmpInt integer in the addIncohResonance function.
631  for (i = 0; i < nIncohAmp_; i++) {
632  getChar << incohResPairAmp_[i] << " ";
633  }
634 
635  getChar << std::endl;
636 
637  // Write out the fSqSum = |ff|^2, where |ff|^2 = incohResAmp()
638  for (i = nAmp_; i < nAmp_+nIncohAmp_; i++) {
639  getChar << fSqSum_[i] << " ";
640  }
641 
642  getChar << std::endl;
643 
644  // Similar to fSqSum, but with the efficiency term included.
645  for (i = nAmp_; i < nAmp_+nIncohAmp_; i++) {
646  getChar << fSqEffSum_[i] << " ";
647  }
648 
649  getChar << std::endl;
650 
651 }
652 
654 {
655  // Function to add a resonance in a Dalitz plot.
656  // No check is made w.r.t flavour and charge conservation rules, and so
657  // the user is responsible for checking the internal consistency of
658  // their function statements with these laws. For example, the program
659  // will not prevent the user from asking for a rho resonance in a K-pi
660  // pair or a K* resonance in a pi-pi pair.
661  // However, to assist the user, a summary of the resonant structure requested
662  // by the user is printed before the program runs. It is important to check this
663  // information when you first define your Dalitz plot model before doing
664  // any fitting/generating.
665  // Arguments are: resonance name, integer to specify the resonance track pairing
666  // (1 => m_23, 2 => m_13, 3 => m_12), i.e. the bachelor track number.
667  // The third argument resType specifies whether the resonance is a Breit-Wigner (BW)
668  // Relativistic Breit-Wigner (RelBW) or Flatte distribution (Flatte), for example.
669 
670  if( LauAbsResonance::isIncoherentModel(resType) == true ) {
671  std::cerr<<"ERROR in LauIsobarDynamics::addResonance : Resonance type \""<<resType<<"\" not allowed for a coherent resonance"<<std::endl;
672  return 0;
673  }
674 
675  LauResonanceMaker& resonanceMaker = LauResonanceMaker::get();
676  LauAbsResonance *theResonance = resonanceMaker.getResonance(daughters_, resName, resPairAmpInt, resType, bwCategory, bwType);
677 
678  if (theResonance == 0) {
679  std::cerr<<"ERROR in LauIsobarDynamics::addResonance : Couldn't create the resonance \""<<resName<<"\""<<std::endl;
680  return 0;
681  }
682 
683  // implement the helicity flip here
684  if (flipHelicity_ && daughters_->getCharge(resPairAmpInt) == 0 && daughters_->getChargeParent() == 0 && daughters_->getTypeParent() > 0 ) {
685  if ( ( resPairAmpInt == 1 && TMath::Abs(daughters_->getTypeDaug2()) == TMath::Abs(daughters_->getTypeDaug3()) ) ||
686  ( resPairAmpInt == 2 && TMath::Abs(daughters_->getTypeDaug1()) == TMath::Abs(daughters_->getTypeDaug3()) ) ||
687  ( resPairAmpInt == 3 && TMath::Abs(daughters_->getTypeDaug1()) == TMath::Abs(daughters_->getTypeDaug2()) ) ) {
688  theResonance->flipHelicity(kTRUE);
689  }
690  }
691 
692  // Set the resonance name and what track is the bachelor
693  TString resonanceName = theResonance->getResonanceName();
694  resTypAmp_.push_back(resonanceName);
695 
696  // Always force the non-resonant amplitude pair to have resPairAmp = 0
697  // in case the user chooses the wrong number.
698  if ( resType == LauAbsResonance::FlatNR ||
699  resType == LauAbsResonance::NRModel ) {
700  std::cout<<"INFO in LauIsobarDynamics::addResonance : Setting resPairAmp to 0 for "<<resonanceName<<" contribution."<<std::endl;
701  resPairAmp_.push_back(0);
702  } else {
703  resPairAmp_.push_back(resPairAmpInt);
704  }
705 
706  // Increment the number of resonance amplitudes we have so far
707  ++nAmp_;
708 
709  // Finally, add the resonance object to the internal array
710  sigResonances_.push_back(theResonance);
711 
712  std::cout<<"INFO in LauIsobarDynamics::addResonance : Successfully added resonance. Total number of resonances so far = "<<nAmp_<<std::endl;
713 
714  return theResonance;
715 }
716 
717 LauAbsResonance* LauIsobarDynamics::addIncoherentResonance(const TString& resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType)
718 {
719  // Function to add an incoherent resonance in a Dalitz plot.
720  // No check is made w.r.t flavour and charge conservation rules, and so
721  // the user is responsible for checking the internal consistency of
722  // their function statements with these laws. For example, the program
723  // will not prevent the user from asking for a rho resonance in a K-pi
724  // pair or a K* resonance in a pi-pi pair.
725  // However, to assist the user, a summary of the resonant structure requested
726  // by the user is printed before the program runs. It is important to check this
727  // information when you first define your Dalitz plot model before doing
728  // any fitting/generating.
729  // Arguments are: resonance name, integer to specify the resonance track pairing
730  // (1 => m_23, 2 => m_13, 3 => m_12), i.e. the bachelor track number.
731  // The third argument resType specifies the shape of the resonance
732  // Gaussian (GaussIncoh), for example.
733 
734  if( LauAbsResonance::isIncoherentModel(resType) == false ) {
735  std::cerr<<"ERROR in LauIsobarDynamics::addIncohResonance : Resonance type \""<<resType<<"\" not allowed for an incoherent resonance"<<std::endl;
736  return 0;
737  }
738 
739  LauResonanceMaker& resonanceMaker = LauResonanceMaker::get();
740  LauAbsIncohRes *theResonance = dynamic_cast<LauAbsIncohRes*>( resonanceMaker.getResonance(daughters_, resName, resPairAmpInt, resType) );
741 
742  if (theResonance == 0) {
743  std::cerr<<"ERROR in LauIsobarDynamics::addIncohResonance : Couldn't create the resonance \""<<resName<<"\""<<std::endl;
744  return 0;
745  }
746 
747  // Set the resonance name and what track is the bachelor
748  TString resonanceName = theResonance->getResonanceName();
749  incohResTypAmp_.push_back(resonanceName);
750  incohResPairAmp_.push_back(resPairAmpInt);
751 
752  // Increment the number of resonance amplitudes we have so far
753  ++nIncohAmp_;
754 
755  // Finally, add the resonance object to the internal array
756  sigIncohResonances_.push_back(theResonance);
757 
758  std::cout<<"INFO in LauIsobarDynamics::addIncohResonance : Successfully added incoherent resonance. Total number of incoherent resonances so far = "<<nIncohAmp_<<std::endl;
759 
760  return theResonance;
761 }
762 
763 void LauIsobarDynamics::defineKMatrixPropagator(const TString& propName, const TString& paramFileName, Int_t resPairAmpInt,
764  Int_t nChannels, Int_t nPoles, Int_t rowIndex)
765 {
766  // Define the K-matrix propagator. The resPairAmpInt integer specifies which mass combination should be used
767  // for the invariant mass-squared variable "s". The pole masses and coupling constants are defined in the
768  // paramFileName parameter file.
769  // The number of channels and poles are defined by the nChannels and nPoles integers, respectively.
770  // The integer rowIndex specifies which row of the propagator should be used when
771  // summing over all amplitude channels: S-wave will be the first row, so rowIndex = 1.
772 
773  // Check that the rowIndex is valid
774  if (rowIndex < 1 || rowIndex > nChannels) {
775  std::cerr << "ERROR in LauIsobarDynamics::defineKMatrixPropagator. The rowIndex, which is set to "
776  << rowIndex << ", must be between 1 and the number of channels "
777  << nChannels << std::endl;
778  gSystem->Exit(EXIT_FAILURE);
779  }
780 
781  TString propagatorName(propName), parameterFile(paramFileName);
782 
783  LauKMatrixPropagator* thePropagator = LauKMatrixPropFactory::getInstance()->getPropagator(propagatorName, parameterFile,
784  resPairAmpInt, nChannels,
785  nPoles, rowIndex);
786  kMatrixPropagators_[propagatorName] = thePropagator;
787 
788 }
789 
790 void LauIsobarDynamics::addKMatrixProdPole(const TString& poleName, const TString& propName, Int_t poleIndex)
791 {
792 
793  // Add a K-matrix production pole term, using the K-matrix propagator given by the propName.
794  // Here, poleIndex is the integer specifying the pole number.
795 
796  // First, find the K-matrix propagator.
797  KMPropMap::iterator mapIter = kMatrixPropagators_.find(propName);
798  if (mapIter != kMatrixPropagators_.end()) {
799 
800  LauKMatrixPropagator* thePropagator = mapIter->second;
801 
802  // Make sure the pole index is valid
803  Int_t nPoles = thePropagator->getNPoles();
804  if (poleIndex < 1 || poleIndex > nPoles) {
805  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdPole : The pole index "<<poleIndex
806  <<" is not between 1 and "<<nPoles<<". Not adding production pole "<<poleName
807  <<" for K-matrix propagator "<<propName<<std::endl;
808  return;
809  }
810 
811  // Now add the K-matrix production pole amplitude to the vector of LauAbsResonance pointers.
812  Int_t resPairAmpInt = thePropagator->getResPairAmpInt();
813  LauAbsResonance* prodPole = new LauKMatrixProdPole(poleName, poleIndex, resPairAmpInt, thePropagator, daughters_);
814 
815  resTypAmp_.push_back(poleName);
816  resPairAmp_.push_back(resPairAmpInt);
817 
818  ++nAmp_;
819  sigResonances_.push_back(prodPole);
820 
821  // Also store the propName-poleName pair for calculating total fit fractions later on
822  // (avoiding the need to use dynamic casts to check which resonances are of the K-matrix type)
823  kMatrixPropSet_[poleName] = propName;
824 
825  std::cout<<"INFO in LauIsobarDynamics::addKMatrixProdPole : Successfully added K-matrix production pole term. Total number of resonances so far = "<<nAmp_<<std::endl;
826 
827  } else {
828 
829  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdPole : The propagator of the name "<<propName
830  <<" could not be found for the production pole "<<poleName<<std::endl;
831 
832  }
833 
834 }
835 
836 
837 void LauIsobarDynamics::addKMatrixProdSVP(const TString& SVPName, const TString& propName, Int_t channelIndex)
838 {
839 
840  // Add a K-matrix production "slowly-varying part" (SVP) term, using the K-matrix propagator
841  // given by the propName. Here, channelIndex is the integer specifying the channel number.
842 
843  // First, find the K-matrix propagator.
844  KMPropMap::iterator mapIter = kMatrixPropagators_.find(propName);
845  if (mapIter != kMatrixPropagators_.end()) {
846 
847  LauKMatrixPropagator* thePropagator = mapIter->second;
848 
849  // Make sure the channel index is valid
850  Int_t nChannels = thePropagator->getNChannels();
851  if (channelIndex < 1 || channelIndex > nChannels) {
852  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdSVP : The channel index "<<channelIndex
853  <<" is not between 1 and "<<nChannels<<". Not adding production slowly-varying part "<<SVPName
854  <<" for K-matrix propagator "<<propName<<std::endl;
855  return;
856  }
857 
858  // Now add the K-matrix production SVP amplitude to the vector of LauAbsResonance pointers.
859  Int_t resPairAmpInt = thePropagator->getResPairAmpInt();
860  LauAbsResonance* prodSVP = new LauKMatrixProdSVP(SVPName, channelIndex, resPairAmpInt, thePropagator, daughters_);
861 
862  resTypAmp_.push_back(SVPName);
863  resPairAmp_.push_back(resPairAmpInt);
864 
865  ++nAmp_;
866  sigResonances_.push_back(prodSVP);
867 
868  // Also store the SVPName-propName pair for calculating total fit fractions later on
869  // (avoiding the need to use dynamic casts to check which resonances are of the K-matrix type)
870  kMatrixPropSet_[SVPName] = propName;
871 
872  std::cout<<"INFO in LauIsobarDynamics::addKMatrixProdSVP : Successfully added K-matrix production slowly-varying (SVP) term. Total number of resonances so far = "<<nAmp_<<std::endl;
873 
874  } else {
875 
876  std::cerr<<"ERROR in LauIsobarDynamics::addKMatrixProdSVP : The propagator of the name "<<propName
877  <<" could not be found for the production slowly-varying part "<<SVPName<<std::endl;
878 
879  }
880 }
881 
882 Int_t LauIsobarDynamics::resonanceIndex(const TString& resName) const
883 {
884  Int_t index(0);
885  const LauAbsResonance* theResonance(0);
886 
887  for (std::vector<LauAbsResonance*>::const_iterator iter=sigResonances_.begin(); iter!=sigResonances_.end(); ++iter) {
888  theResonance = (*iter);
889  if (theResonance != 0) {
890  const TString& resString = theResonance->getResonanceName();
891  if (resString == resName) {
892  return index;
893  }
894  }
895  ++index;
896  }
897 
898  for (std::vector<LauAbsIncohRes*>::const_iterator iter=sigIncohResonances_.begin(); iter!=sigIncohResonances_.end(); ++iter) {
899  theResonance = (*iter);
900  if (theResonance != 0) {
901  const TString& resString = theResonance->getResonanceName();
902  if (resString == resName) {
903  return index;
904  }
905  }
906  ++index;
907  }
908 
909  return -1;
910 }
911 
912 Bool_t LauIsobarDynamics::hasResonance(const TString& resName) const
913 {
914  const Int_t index = this->resonanceIndex(resName);
915  if (index < 0) {
916  return kFALSE;
917  } else {
918  return kTRUE;
919  }
920 }
921 
922 const LauAbsResonance* LauIsobarDynamics::getResonance(const UInt_t resIndex) const
923 {
924  if ( resIndex < this->getnCohAmp() ) {
925  return sigResonances_[resIndex];
926  } else if ( resIndex < this->getnTotAmp() ) {
927  return sigIncohResonances_[ resIndex - nAmp_ ];
928  } else {
929  std::cerr<<"ERROR in LauIsobarDynamics::getResonance : Couldn't find resonance with index \""<<resIndex<<"\" in the model."<<std::endl;
930  return 0;
931  }
932 }
933 
935 {
936  if ( resIndex < this->getnCohAmp() ) {
937  return sigResonances_[resIndex];
938  } else if ( resIndex < this->getnTotAmp() ) {
939  return sigIncohResonances_[ resIndex - nAmp_ ];
940  } else {
941  std::cerr<<"ERROR in LauIsobarDynamics::getResonance : Couldn't find resonance with index \""<<resIndex<<"\" in the model."<<std::endl;
942  return 0;
943  }
944 }
945 
947 {
948  const Int_t index = this->resonanceIndex( resName );
949  if ( index < 0 ) {
950  std::cerr<<"ERROR in LauIsobarDynamics::findResonance : Couldn't find resonance with name \""<<resName<<"\" in the model."<<std::endl;
951  return 0;
952  } else {
953  return this->getResonance( index );
954  }
955 }
956 
957 const LauAbsResonance* LauIsobarDynamics::findResonance(const TString& resName) const
958 {
959  const Int_t index = this->resonanceIndex( resName );
960  if ( index < 0 ) {
961  std::cerr<<"ERROR in LauIsobarDynamics::findResonance : Couldn't find resonance with name \""<<resName<<"\" in the model."<<std::endl;
962  return 0;
963  } else {
964  return this->getResonance( index );
965  }
966 }
967 
968 void LauIsobarDynamics::removeCharge(TString& string) const
969 {
970  Ssiz_t index = string.Index("+");
971  if (index != -1) {
972  string.Remove(index,1);
973  }
974  index = string.Index("-");
975  if (index != -1) {
976  string.Remove(index,1);
977  }
978 }
979 
981 {
984  }
985 
986  for (std::vector<LauDPPartialIntegralInfo*>::iterator it = dpPartialIntegralInfo_.begin(); it != dpPartialIntegralInfo_.end(); ++it)
987  {
988  this->calcDPPartialIntegral( *it );
989  }
990 
991  for (UInt_t i = 0; i < nAmp_+nIncohAmp_; ++i) {
992  fNorm_[i] = 0.0;
993  if (fSqSum_[i] > 0.0) {fNorm_[i] = TMath::Sqrt(1.0/(fSqSum_[i]));}
994  }
995 }
996 
998 {
999  if ( ! dpPartialIntegralInfo_.empty() ) {
1000  std::cerr<<"ERROR in LauIsobarDynamics::calcDPNormalisationScheme : Scheme already stored!"<<std::endl;
1001  return;
1002  }
1003 
1004  // The precision for the Gauss-Legendre weights
1005  const Double_t precision(1e-6);
1006 
1007  // Get the rectangle that encloses the DP
1008  Double_t minm13 = kinematics_->getm13Min();
1009  Double_t maxm13 = kinematics_->getm13Max();
1010  Double_t minm23 = kinematics_->getm23Min();
1011  Double_t maxm23 = kinematics_->getm23Max();
1012  Double_t minm12 = kinematics_->getm12Min();
1013  Double_t maxm12 = kinematics_->getm12Max();
1014 
1015  // Find out whether we have narrow resonances in the DP (defined here as width < 20 MeV).
1016  std::multimap<Double_t,Double_t> m13NarrowRes;
1017  std::multimap<Double_t,Double_t> m23NarrowRes;
1018  std::multimap<Double_t,Double_t> m12NarrowRes;
1019  for ( std::vector<LauAbsResonance*>::const_iterator iter = sigResonances_.begin(); iter != sigResonances_.end(); ++iter ) {
1020  Double_t width = (*iter)->getWidth();
1021  if ( width > narrowWidth_ || width == 0.0 ) { continue; }
1022  Double_t mass = (*iter)->getMass();
1023  Int_t pair = (*iter)->getPairInt();
1024  TString name = (*iter)->getResonanceName();
1025  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Found narrow resonance: "<<name<<", mass = "<<mass<<", width = "<<width<<", pair int = "<<pair<<std::endl;
1026  if ( pair == 1 ) {
1027  if ( mass < minm23 || mass > maxm23 ){ continue; }
1028  m23NarrowRes.insert( std::make_pair(width,mass) );
1029  } else if ( pair == 2 ) {
1030  if ( mass < minm13 || mass > maxm13 ){ continue; }
1031  m13NarrowRes.insert( std::make_pair(width,mass) );
1032  } else if ( pair == 3 ) {
1033  if ( mass < minm12 || mass > maxm12 ){ continue; }
1034  m12NarrowRes.insert( std::make_pair(width,mass) );
1035  } else {
1036  std::cerr<<"WARNING in LauIsobarDynamics::calcDPNormalisationScheme : strange pair integer, "<<pair<<", for resonance \""<<(*iter)->getResonanceName()<<std::endl;
1037  }
1038  }
1039 
1040  for ( std::vector<LauAbsIncohRes*>::const_iterator iter = sigIncohResonances_.begin(); iter != sigIncohResonances_.end(); ++iter ) {
1041  Double_t width = (*iter)->getWidth();
1042  if ( width > narrowWidth_ || width == 0.0 ) { continue; }
1043  Double_t mass = (*iter)->getMass();
1044  Int_t pair = (*iter)->getPairInt();
1045  TString name = (*iter)->getResonanceName();
1046  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Found narrow resonance: "<<name<<", mass = "<<mass<<", width = "<<width<<", pair int = "<<pair<<std::endl;
1047  if ( pair == 1 ) {
1048  if ( mass < minm23 || mass > maxm23 ){
1049  std::cout<<" : But its pole is outside the kinematically allowed range, so will not consider it narrow for the purposes of integration"<<std::endl;
1050  } else {
1051  m23NarrowRes.insert( std::make_pair(width,mass) );
1052  }
1053  } else if ( pair == 2 ) {
1054  if ( mass < minm13 || mass > maxm13 ){
1055  std::cout<<" : But its pole is outside the kinematically allowed range, so will not consider it narrow for the purposes of integration"<<std::endl;
1056  } else {
1057  m13NarrowRes.insert( std::make_pair(width,mass) );
1058  }
1059  } else if ( pair == 3 ) {
1060  if ( mass < minm12 || mass > maxm12 ){
1061  std::cout<<" : But its pole is outside the kinematically allowed range, so will not consider it narrow for the purposes of integration"<<std::endl;
1062  } else {
1063  m12NarrowRes.insert( std::make_pair(width,mass) );
1064  }
1065  } else {
1066  std::cerr<<"WARNING in LauIsobarDynamics::calcDPNormalisationScheme : strange pair integer, "<<pair<<", for resonance \""<<(*iter)->getResonanceName()<<std::endl;
1067  }
1068  }
1069 
1070  // Find the narrowest resonance in each mass pairing
1071  const Bool_t e12 = m12NarrowRes.empty();
1072  const Bool_t e13 = m13NarrowRes.empty();
1073  const Bool_t e23 = m23NarrowRes.empty();
1074  //const UInt_t s12 = e12 ? 0 : m12NarrowRes.size();
1075  const UInt_t s13 = e13 ? 0 : m13NarrowRes.size();
1076  const UInt_t s23 = e23 ? 0 : m23NarrowRes.size();
1077  //const Double_t w12 = e12 ? DBL_MAX : m12NarrowRes.begin()->first / binningFactor_;
1078  const Double_t w13 = e13 ? DBL_MAX : m13NarrowRes.begin()->first / binningFactor_;
1079  const Double_t w23 = e23 ? DBL_MAX : m23NarrowRes.begin()->first / binningFactor_;
1080 
1081  // Start off with default bin widths (5 MeV for m13 vs m23 and 0.001 for m' vs theta')
1082  Double_t m13BinWidth = m13BinWidth_;
1083  Double_t m23BinWidth = m23BinWidth_;
1084  Double_t mPrimeBinWidth = mPrimeBinWidth_;
1085  Double_t thPrimeBinWidth = thPrimeBinWidth_;
1086 
1087  // Depending on how many narrow resonances we have and where they
1088  // are we adopt different approaches
1089  if ( e12 && e13 && e23 ) {
1090  // If we have no narrow resonances just integrate the whole
1091  // DP with the standard bin widths
1092  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : No narrow resonances found, integrating over whole Dalitz plot..."<<std::endl;
1093  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1094  } else if ( ! e12 ) {
1095  // We have at least one narrow resonance in m12
1096  // Switch to using the square DP for the integration
1097  // TODO - for the time being just use a single, reasonably fine by default and tunable, grid
1098  // - can later consider whether there's a need to split up the mPrime axis into regions around particularly narrow resonances in m12
1099  // - but it seems that this isn't really needed since even the default tune gives a good resolution for most narrow resonances such as phi / chi_c0
1100  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One or more narrow resonances found in m12, integrating over whole square Dalitz plot with bin widths of "<<mPrimeBinWidth<<" in mPrime and "<<thPrimeBinWidth<<" in thetaPrime..."<<std::endl;
1101  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(0.0, 1.0, 0.0, 1.0, mPrimeBinWidth, thPrimeBinWidth, precision, nAmp_, nIncohAmp_, kTRUE, kinematics_));
1102  } else if ( s13==1 && e23 ) {
1103  // We have a single narrow resonance in m13
1104  // Divide the plot into 3 regions: the resonance band and
1105  // the two areas either side.
1106  Double_t mass = m13NarrowRes.begin()->second;
1107  Double_t width = m13NarrowRes.begin()->first;
1108  Double_t resMin = mass - 5.0*width;
1109  Double_t resMax = mass + 5.0*width;
1110  // if the resonance is close to threshold just go from
1111  // threshold to resMax, otherwise treat threshold to resMin
1112  // as a separate region
1113  if ( resMin < (minm13+50.0*m13BinWidth_) ) {
1114  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One narrow resonance found in m13, close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
1115  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1116  m13BinWidth = w13;
1117  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1118  } else {
1119  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One narrow resonance found in m13, dividing Dalitz plot into three regions..."<<std::endl;
1120  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1121  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1122  m13BinWidth = w13;
1123  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin, resMax, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1124  }
1125  } else if ( s13==2 && e23 ) {
1126  // We have a two narrow resonances in m13
1127  // Divide the plot into 5 regions: the resonance bands,
1128  // the two areas either side and the area in between.
1129  std::multimap<Double_t,Double_t> massordered;
1130  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1131  massordered.insert( std::make_pair( iter->second, iter->first ) );
1132  }
1133  Double_t res1Mass = massordered.begin()->first;
1134  Double_t res1Width = massordered.begin()->second;
1135  Double_t res1Min = res1Mass - 5.0*res1Width;
1136  Double_t res1Max = res1Mass + 5.0*res1Width;
1137  Double_t res2Mass = (++(massordered.begin()))->first;
1138  Double_t res2Width = (++(massordered.begin()))->second;
1139  Double_t res2Min = res2Mass - 5.0*res2Width;
1140  Double_t res2Max = res2Mass + 5.0*res2Width;
1141  // if the resonance is close to threshold just go from
1142  // threshold to resMax, otherwise treat threshold to resMin
1143  // as a separate region
1144  if ( res1Min < (minm13+50.0*m13BinWidth_) ) {
1145  if ( res1Max > res2Min ) {
1146  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found in m13, both close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
1147  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1148  m13BinWidth = TMath::Min(res1Width,res2Width)/100.0;
1149  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1150  } else {
1151  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found in m13, one close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
1152  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res1Max, res2Min, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1153  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1154  m13BinWidth = res1Width/100.0;
1155  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, res1Max, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1156  m13BinWidth = res2Width/100.0;
1157  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res2Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1158  }
1159  } else {
1160  if ( res1Max > res2Min ) {
1161  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found close together in m13, dividing Dalitz plot into three regions..."<<std::endl;
1162  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, res1Min, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1163  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1164  m13BinWidth = TMath::Min(res1Width,res2Width)/100.0;
1165  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res1Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1166  } else {
1167  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found in m13, dividing Dalitz plot into five regions..."<<std::endl;
1168  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, res1Min, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1169  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res1Max, res2Min, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1170  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res2Max, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1171  m13BinWidth = res1Width/100.0;
1172  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res1Min, res1Max, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1173  m13BinWidth = res2Width/100.0;
1174  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(res2Min, res2Max, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1175  }
1176  }
1177  } else if ( s23==1 && e13 ) {
1178  // We have a single narrow resonance in m23
1179  // Divide the plot into 3 regions: the resonance band and
1180  // the two areas either side.
1181  Double_t mass = m23NarrowRes.begin()->second;
1182  Double_t width = m23NarrowRes.begin()->first;
1183  Double_t resMin = mass - 5.0*width;
1184  Double_t resMax = mass + 5.0*width;
1185  // if the resonance is close to threshold just go from
1186  // threshold to resMax, otherwise treat threshold to resMin
1187  // as a separate region
1188  if ( resMin < (minm23+50.0*m23BinWidth_) ) {
1189  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One narrow resonance found in m23, close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
1190  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1191  m23BinWidth = w23;
1192  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, resMax, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1193  } else {
1194  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One narrow resonance found in m23, dividing Dalitz plot into three regions..."<<std::endl;
1195  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, resMin, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1196  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1197  m23BinWidth = w23;
1198  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, resMin, resMax, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1199  }
1200  } else if ( s23==2 && e13 ) {
1201  // We have a two narrow resonances in m23
1202  // Divide the plot into 5 regions: the resonance bands,
1203  // the two areas either side and the area in between.
1204  std::multimap<Double_t,Double_t> massordered;
1205  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1206  massordered.insert( std::make_pair( iter->second, iter->first ) );
1207  }
1208  Double_t res1Mass = massordered.begin()->first;
1209  Double_t res1Width = massordered.begin()->second;
1210  Double_t res1Min = res1Mass - 5.0*res1Width;
1211  Double_t res1Max = res1Mass + 5.0*res1Width;
1212  Double_t res2Mass = (++(massordered.begin()))->first;
1213  Double_t res2Width = (++(massordered.begin()))->second;
1214  Double_t res2Min = res2Mass - 5.0*res2Width;
1215  Double_t res2Max = res2Mass + 5.0*res2Width;
1216  // if the resonance is close to threshold just go from
1217  // threshold to resMax, otherwise treat threshold to resMin
1218  // as a separate region
1219  if ( res1Min < (minm23+50.0*m23BinWidth_) ) {
1220  if ( res1Max > res2Min ) {
1221  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found in m23, both close to threshold, dividing Dalitz plot into two regions..."<<std::endl;
1222  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1223  m23BinWidth = TMath::Min(res1Width,res2Width)/100.0;
1224  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, res2Max, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1225  } else {
1226  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found in m23, one close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
1227  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res1Max, res2Min, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1228  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1229  m23BinWidth = res1Width/100.0;
1230  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, res1Max, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1231  m23BinWidth = res2Width/100.0;
1232  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res2Min, res2Max, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1233  }
1234  } else {
1235  if ( res1Max > res2Min ) {
1236  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found close together in m23, dividing Dalitz plot into three regions..."<<std::endl;
1237  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, res1Min, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1238  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1239  m23BinWidth = TMath::Min(res1Width,res2Width)/100.0;
1240  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res1Min, res2Max, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1241  } else {
1242  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Two narrow resonances found in m23, dividing Dalitz plot into five regions..."<<std::endl;
1243  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, res1Min, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1244  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res1Max, res2Min, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1245  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res2Max, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1246  m23BinWidth = res1Width/100.0;
1247  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res1Min, res1Max, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1248  m23BinWidth = res2Width/100.0;
1249  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, res2Min, res2Max, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1250  }
1251  }
1252  } else if ( s13==1 && s23==1 ) {
1253  // We have a single narrow resonance in both m13 and m23
1254  // Divide the plot into 9 regions: the point where the
1255  // resonance bands cross, the four other parts of the bands
1256  // and the four remaining areas of the DP.
1257  Double_t mass23 = m23NarrowRes.begin()->second;
1258  Double_t width23 = m23NarrowRes.begin()->first;
1259  Double_t resMin23 = mass23 - 5.0*width23;
1260  Double_t resMax23 = mass23 + 5.0*width23;
1261  Double_t mass13 = m13NarrowRes.begin()->second;
1262  Double_t width13 = m13NarrowRes.begin()->first;
1263  Double_t resMin13 = mass13 - 5.0*width13;
1264  Double_t resMax13 = mass13 + 5.0*width13;
1265  // if either resonance is close to threshold just go from
1266  // threshold to resMax, otherwise treat threshold to resMin
1267  // as a separate region
1268  if ( resMin13 < (minm13+50.0*m13BinWidth_) && resMin23 < (minm23+50.0*m23BinWidth_) ) {
1269  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One narrow resonance found in m13 and one in m23, both close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
1270  m13BinWidth = m13BinWidth_;
1271  m23BinWidth = m23BinWidth_;
1272  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1273  m13BinWidth = m13BinWidth_;
1274  m23BinWidth = w23;
1275  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1276  m13BinWidth = w13;
1277  m23BinWidth = m23BinWidth_;
1278  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1279  m13BinWidth = w13;
1280  m23BinWidth = w23;
1281  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1282  } else if ( resMin13 < (minm13+50.0*m13BinWidth_) ) {
1283  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : 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;
1284  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1285  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1286  m13BinWidth = m13BinWidth_;
1287  m23BinWidth = w23;
1288  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1289  m13BinWidth = w13;
1290  m23BinWidth = m23BinWidth_;
1291  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1292  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1293  m13BinWidth = w13;
1294  m23BinWidth = w23;
1295  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1296  } else if ( resMin23 < (minm23+50.0*m23BinWidth_) ) {
1297  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : 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;
1298  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1299  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1300  m13BinWidth = m13BinWidth_;
1301  m23BinWidth = w23;
1302  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1303  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1304  m13BinWidth = w13;
1305  m23BinWidth = m23BinWidth_;
1306  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1307  m13BinWidth = w13;
1308  m23BinWidth = w23;
1309  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1310  } else {
1311  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : One narrow resonance found in both m13 and m23, neither close to threshold, dividing Dalitz plot into nine regions..."<<std::endl;
1312  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1313  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1314  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1315  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1316  m13BinWidth = m13BinWidth_;
1317  m23BinWidth = w23;
1318  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1319  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1320  m13BinWidth = w13;
1321  m23BinWidth = m23BinWidth_;
1322  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1323  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1324  m13BinWidth = w13;
1325  m23BinWidth = w23;
1326  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1327  }
1328  } else if ( e23 && s13>1 ) {
1329  // We have multiple narrow resonances in m13 only.
1330  // Divide the plot into 2 regions: threshold to the most
1331  // massive of the narrow resonances, and the rest
1332  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in m13, dividing Dalitz plot into two regions..."<<std::endl;
1333  Double_t mass = 0.0;
1334  Double_t width = 0.0;
1335  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1336  if ( mass < iter->second ) {
1337  mass = iter->second;
1338  width = iter->first;
1339  }
1340  }
1341  Double_t resMax = mass + 5.0*width;
1342  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax, maxm13, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1343  m13BinWidth = w13;
1344  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax, minm23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1345  } else if ( e13 && s23>1 ) {
1346  // We have multiple narrow resonances in m23 only.
1347  // Divide the plot into 2 regions: threshold to the most
1348  // massive of the narrow resonances, and the rest
1349  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in m23, dividing Dalitz plot into two regions..."<<std::endl;
1350  Double_t mass = 0.0;
1351  Double_t width = 0.0;
1352  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1353  if ( mass < iter->second ) {
1354  mass = iter->second;
1355  width = iter->first;
1356  }
1357  }
1358  Double_t resMax = mass + 5.0*width;
1359  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, resMax, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1360  m23BinWidth = w23;
1361  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, maxm13, minm23, resMax, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1362  } else if ( s13==1 && s23>1 ) {
1363  // We've got a single narrow resonance in m13 and multiple
1364  // narrow resonances in m23. Divide the plot into 6 regions.
1365  Double_t mass23 = 0.0;
1366  Double_t width23 = 0.0;
1367  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1368  if ( mass23 < iter->second ) {
1369  mass23 = iter->second;
1370  width23 = iter->first;
1371  }
1372  }
1373  Double_t resMax23 = mass23 + 5.0*width23;
1374  Double_t mass13 = m13NarrowRes.begin()->second;
1375  Double_t width13 = m13NarrowRes.begin()->first;
1376  Double_t resMin13 = mass13 - 5.0*width13;
1377  Double_t resMax13 = mass13 + 5.0*width13;
1378  // if the m13 resonance is close to threshold just go from
1379  // threshold to resMax, otherwise treat threshold to resMin
1380  // as a separate region
1381  if ( resMin13 < (minm13+50.0*m13BinWidth_) ) {
1382  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in m23 and one in m13, close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
1383  m13BinWidth = m13BinWidth_;
1384  m23BinWidth = m23BinWidth_;
1385  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1386  m13BinWidth = m13BinWidth_;
1387  m23BinWidth = w23;
1388  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1389  m13BinWidth = w13;
1390  m23BinWidth = m23BinWidth_;
1391  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1392  m13BinWidth = w13;
1393  m23BinWidth = w23;
1394  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1395  } else {
1396  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in m23 and one in m13, not close to threshold, dividing Dalitz plot into six regions..."<<std::endl;
1397  m13BinWidth = m13BinWidth_;
1398  m23BinWidth = m23BinWidth_;
1399  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1400  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1401  m13BinWidth = m13BinWidth_;
1402  m23BinWidth = w23;
1403  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMin13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1404  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1405  m13BinWidth = w13;
1406  m23BinWidth = m23BinWidth_;
1407  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1408  m13BinWidth = w13;
1409  m23BinWidth = w23;
1410  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMin13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1411  }
1412  } else if ( s13>1 && s23==1 ) {
1413  // We've got a single narrow resonance in m23 and multiple
1414  // narrow resonances in m13. Divide the plot into 6 regions.
1415  Double_t mass13 = 0.0;
1416  Double_t width13 = 0.0;
1417  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1418  if ( mass13 < iter->second ) {
1419  mass13 = iter->second;
1420  width13 = iter->first;
1421  }
1422  }
1423  Double_t resMax13 = mass13 + 5.0*width13;
1424  Double_t mass23 = m23NarrowRes.begin()->second;
1425  Double_t width23 = m23NarrowRes.begin()->first;
1426  Double_t resMin23 = mass23 - 5.0*width23;
1427  Double_t resMax23 = mass23 + 5.0*width23;
1428  // if the m23 resonance is close to threshold just go from
1429  // threshold to resMax, otherwise treat threshold to resMin
1430  // as a separate region
1431  if ( resMin23 < (minm23+50.0*m23BinWidth_) ) {
1432  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in m13 and one in m23, close to threshold, dividing Dalitz plot into four regions..."<<std::endl;
1433  m13BinWidth = m13BinWidth_;
1434  m23BinWidth = m23BinWidth_;
1435  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1436  m13BinWidth = m13BinWidth_;
1437  m23BinWidth = w23;
1438  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1439  m13BinWidth = w13;
1440  m23BinWidth = m23BinWidth_;
1441  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1442  m13BinWidth = w13;
1443  m23BinWidth = w23;
1444  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1445  } else {
1446  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in m13 and one in m23, not close to threshold, dividing Dalitz plot into six regions..."<<std::endl;
1447  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1448  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1449  m13BinWidth = m13BinWidth_;
1450  m23BinWidth = w23;
1451  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1452  m13BinWidth = w13;
1453  m23BinWidth = m23BinWidth_;
1454  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, minm23, resMin23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1455  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1456  m13BinWidth = w13;
1457  m23BinWidth = w23;
1458  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMin23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1459  }
1460  } else {
1461  // We've got multiple narrow resonances in both m13 and m23.
1462  // Divide the plot into 4 regions.
1463  std::cout<<"INFO in LauIsobarDynamics::calcDPNormalisationScheme : Multiple narrow resonances found in both m13 and m23, dividing Dalitz plot into four regions..."<<std::endl;
1464  Double_t mass13 = 0.0;
1465  Double_t width13 = 0.0;
1466  for ( std::map<Double_t,Double_t>::const_iterator iter = m13NarrowRes.begin(); iter != m13NarrowRes.end(); ++iter ) {
1467  if ( mass13 < iter->second ) {
1468  mass13 = iter->second;
1469  width13 = iter->first;
1470  }
1471  }
1472  Double_t resMax13 = mass13 + 5.0*width13;
1473 
1474  Double_t mass23 = 0.0;
1475  Double_t width23 = 0.0;
1476  for ( std::map<Double_t,Double_t>::const_iterator iter = m23NarrowRes.begin(); iter != m23NarrowRes.end(); ++iter ) {
1477  if ( mass23 < iter->second ) {
1478  mass23 = iter->second;
1479  width23 = iter->first;
1480  }
1481  }
1482  Double_t resMax23 = mass23 + 5.0*width23;
1483 
1484  m13BinWidth = m13BinWidth_;
1485  m23BinWidth = m23BinWidth_;
1486  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1487  m13BinWidth = m13BinWidth_;
1488  m23BinWidth = w23;
1489  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(resMax13, maxm13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1490  m13BinWidth = w13;
1491  m23BinWidth = m23BinWidth_;
1492  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, resMax23, maxm23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1493  m13BinWidth = w13;
1494  m23BinWidth = w23;
1495  dpPartialIntegralInfo_.push_back(new LauDPPartialIntegralInfo(minm13, resMax13, minm23, resMax23, m13BinWidth, m23BinWidth, precision, nAmp_, nIncohAmp_));
1496  }
1497 
1498  normalizationSchemeDone_ = kTRUE;
1499 }
1500 
1501 void LauIsobarDynamics::setIntegralBinWidths(const Double_t m13BinWidth, const Double_t m23BinWidth,
1502  const Double_t mPrimeBinWidth, const Double_t thPrimeBinWidth)
1503 {
1504  // Set the bin widths for the m13 vs m23 integration grid
1505  m13BinWidth_ = m13BinWidth;
1506  m23BinWidth_ = m23BinWidth;
1507 
1508  // Set the bin widths for the m' vs theta' integration grid
1509  mPrimeBinWidth_ = mPrimeBinWidth;
1510  thPrimeBinWidth_ = thPrimeBinWidth;
1511 }
1512 
1514 {
1515  // Calculate the integrals for all parts of the amplitude in the given region of the DP
1516 
1517  const Bool_t squareDP = intInfo->getSquareDP();
1518  const UInt_t nm13Points = intInfo->getnm13Points();
1519  const UInt_t nm23Points = intInfo->getnm23Points();
1520 
1521  //Double_t dpArea(0.0);
1522  for (UInt_t i = 0; i < nm13Points; ++i) {
1523 
1524  const Double_t m13 = intInfo->getM13Value(i);
1525  const Double_t m13Sq = m13*m13;
1526 
1527  for (UInt_t j = 0; j < nm23Points; ++j) {
1528 
1529  const Double_t m23 = intInfo->getM23Value(j);
1530  const Double_t m23Sq = m23*m23;
1531  const Double_t weight = intInfo->getWeight(i,j);
1532 
1533  // Calculate the integral contributions for each resonance.
1534  // Only points within the DP area contribute.
1535  // This also calculates the total DP area as a check.
1536  // NB if squareDP is true, m13 and m23 are actually mPrime and thetaPrime
1537  Bool_t withinDP = squareDP ? kinematics_->withinSqDPLimits(m13, m23) : kinematics_->withinDPLimits(m13Sq, m23Sq);
1538  if (withinDP == kTRUE) {
1539 
1540  if ( squareDP ) {
1541  // NB m13 and m23 are actually mPrime and thetaPrime
1542  kinematics_->updateSqDPKinematics(m13, m23);
1543  } else {
1544  kinematics_->updateKinematics(m13Sq, m23Sq);
1545  }
1546 
1547  this->calculateAmplitudes(intInfo, i, j);
1548 
1549  this->addGridPointToIntegrals(weight);
1550 
1551  // Increment total DP area
1552  //dpArea += weight;
1553  }
1554 
1555  } // j weights loop
1556  } // i weights loop
1557 
1558  // Print out DP area to check whether we have a sensible output
1559  //std::cout<<" : dpArea = "<<dpArea<<std::endl;
1560 }
1561 
1562 void LauIsobarDynamics::calculateAmplitudes( LauDPPartialIntegralInfo* intInfo, const UInt_t m13Point, const UInt_t m23Point )
1563 {
1564  const std::set<UInt_t>::const_iterator intEnd = integralsToBeCalculated_.end();
1565 
1566  for (UInt_t iAmp = 0; iAmp < nAmp_; ++iAmp) {
1567 
1568  if ( integralsToBeCalculated_.find(iAmp) != intEnd ) {
1569  // Calculate the dynamics for this resonance
1570  this->resAmp(iAmp);
1571 
1572  // Store the new value in the integration info object
1573  intInfo->storeAmplitude( m13Point, m23Point, iAmp, ff_[iAmp] );
1574  } else {
1575  // Retrieve the cached value of the amplitude
1576  ff_[iAmp] = intInfo->getAmplitude( m13Point, m23Point, iAmp );
1577  }
1578  }
1579 
1580  for (UInt_t iAmp = 0; iAmp < nIncohAmp_; ++iAmp) {
1581 
1582  if ( integralsToBeCalculated_.find(iAmp+nAmp_) != intEnd ) {
1583  // Calculate the dynamics for this resonance
1584  this->incohResAmp(iAmp);
1585 
1586  // Store the new value in the integration info object
1587  intInfo->storeIntensity( m13Point, m23Point, iAmp, incohInten_[iAmp] );
1588  } else {
1589  // Retrieve the cached value of the amplitude
1590  incohInten_[iAmp] = intInfo->getIntensity( m13Point, m23Point, iAmp );
1591  }
1592  }
1593 
1594  // If we haven't cached the data, then we need to find out the efficiency.
1595  eff_ = this->retrieveEfficiency();
1596 
1597  intInfo->storeEfficiency( m13Point, m23Point, eff_ );
1598 }
1599 
1601 {
1602  std::set<UInt_t>::const_iterator iter = integralsToBeCalculated_.begin();
1603  const std::set<UInt_t>::const_iterator intEnd = integralsToBeCalculated_.end();
1604 
1605  for ( ; iter != intEnd; ++iter) {
1606 
1607  // Calculate the dynamics for this resonance
1608  if(*iter < nAmp_) {
1609  this->resAmp(*iter);
1610  } else {
1611  this->incohResAmp(*iter-nAmp_);
1612  }
1613  }
1614 
1615  // If we haven't cached the data, then we need to find out the efficiency.
1616  eff_ = this->retrieveEfficiency();
1617 }
1618 
1619 void LauIsobarDynamics::calcTotalAmp(const Bool_t useEff)
1620 {
1621  // Reset the total amplitude to zero
1622  totAmp_.zero();
1623 
1624  // Loop over all signal amplitudes
1625  LauComplex ATerm;
1626  for (UInt_t i = 0; i < nAmp_; ++i) {
1627 
1628  // Get the partial complex amplitude - (mag, phase)*(resonance dynamics)
1629  ATerm = Amp_[i]*ff_[i];
1630 
1631  // Scale this contribution by its relative normalisation w.r.t. the whole dynamics
1632  ATerm.rescale(fNorm_[i]);
1633 
1634  // Add this partial amplitude to the sum
1635  totAmp_ += ATerm;
1636 
1637  } // Loop over amplitudes
1638 
1639  // |Sum of partial amplitudes|^2
1640  ASq_ = totAmp_.abs2();
1641 
1642  for (UInt_t i = 0; i < nIncohAmp_; ++i) {
1643 
1644  // Get the partial complex amplitude - (mag, phase)
1645  ATerm = Amp_[i+nAmp_];
1646 
1647  // Scale this contribution by its relative normalisation w.r.t. the whole dynamics
1648  ATerm.rescale(fNorm_[i+nAmp_]);
1649 
1650  // Add this partial amplitude to the sum
1651  ASq_ += ATerm.abs2()*incohInten_[i];
1652 
1653  }
1654 
1655  // Apply the efficiency correction for this event.
1656  // Multiply the amplitude squared sum by the DP efficiency
1657  if ( useEff ) {
1658  ASq_ *= eff_;
1659  }
1660 }
1661 
1663 {
1664  // Combine the Gauss-Legendre weight with the efficiency
1665  const Double_t effWeight = eff_*weight;
1666 
1667  LauComplex fifjEffSumTerm;
1668  LauComplex fifjSumTerm;
1669 
1670  // Calculates the half-matrix of amplitude-squared and interference
1671  // terms (dynamical part only)
1672  // Add the values at this point on the integration grid to the sums
1673  // (one weighted only by the integration weights, one also weighted by
1674  // the efficiency)
1675  for (UInt_t i = 0; i < nAmp_; ++i) {
1676 
1677  // Add the dynamical amplitude squared for this resonance.
1678  Double_t fSqVal = ff_[i].abs2();
1679  fSqSum_[i] += fSqVal*weight;
1680  fSqEffSum_[i] += fSqVal*effWeight;
1681 
1682  for (UInt_t j = i; j < nAmp_; ++j) {
1683 
1684  fifjEffSumTerm = fifjSumTerm = ff_[i]*ff_[j].conj();
1685 
1686  fifjEffSumTerm.rescale(effWeight);
1687  fifjEffSum_[i][j] += fifjEffSumTerm;
1688 
1689  fifjSumTerm.rescale(weight);
1690  fifjSum_[i][j] += fifjSumTerm;
1691  }
1692  }
1693  for (UInt_t i = 0; i < nIncohAmp_; ++i) {
1694 
1695  // Add the dynamical amplitude squared for this resonance.
1696  Double_t fSqVal = incohInten_[i];
1697  fSqSum_[i+nAmp_] += fSqVal*weight;
1698  fSqEffSum_[i+nAmp_] += fSqVal*effWeight;
1699  }
1700 }
1701 
1702 void LauIsobarDynamics::resAmp(const UInt_t index)
1703 {
1704  // Routine to calculate the resonance dynamics (amplitude)
1705  // using the appropriate Breit-Wigner/Form Factors.
1706 
1707  if ( index >= nAmp_ ) {
1708  std::cerr<<"ERROR in LauIsobarDynamics::resAmp : index = "<<index<<" is not within the range 0 to "<<nAmp_-1<<std::endl;
1709  return;
1710  }
1711 
1712  // Get the signal resonance from the stored vector
1713  LauAbsResonance* sigResonance = sigResonances_[index];
1714 
1715  if (sigResonance == 0) {
1716  std::cerr<<"ERROR in LauIsobarDynamics::resAmp : Couldn't retrieve valid resonance with index = "<<index<<std::endl;
1717  ff_[index] = LauComplex(0.0, 0.0);
1718  return;
1719  }
1720 
1721  ff_[index] = sigResonance->amplitude(kinematics_);
1722 
1723  // If we have a symmetrical Dalitz plot, flip the m_13^2 and m_23^2
1724  // variables, recalculate the dynamics, and combine both contributions.
1725  if (symmetricalDP_ == kTRUE && sigResonance->preSymmetrised() == kFALSE ) {
1727 
1728  ff_[index] += sigResonance->amplitude(kinematics_);
1729 
1730  // Flip the m_13^2 and m_23^2 variables back to their original values
1732  }
1733 
1734  // If we have a fully symmetric DP we need to calculate the amplitude at 6 points
1735  else if (fullySymmetricDP_ == kTRUE) {
1736  // rotate and evaluate
1738  ff_[index] += sigResonance->amplitude(kinematics_);
1739 
1740  // rotate and evaluate
1742  ff_[index] += sigResonance->amplitude(kinematics_);
1743 
1744  // rotate, flip and evaluate
1747  ff_[index] += sigResonance->amplitude(kinematics_);
1748 
1749  // rotate and evaluate
1751  ff_[index] += sigResonance->amplitude(kinematics_);
1752 
1753  // rotate and evaluate
1755  ff_[index] += sigResonance->amplitude(kinematics_);
1756 
1757  // rotate and flip to get us back to where we started
1760  }
1761 }
1762 
1763 void LauIsobarDynamics::incohResAmp(const UInt_t index)
1764 {
1765  // Routine to calculate the resonance dynamics (amplitude)
1766  // using the appropriate Breit-Wigner/Form Factors.
1767 
1768  if ( index >= nIncohAmp_ ) {
1769  std::cerr<<"ERROR in LauIsobarDynamics::incohResAmp : index = "<<index<<" is not within the range 0 to "<<nIncohAmp_-1<<std::endl;
1770  return;
1771  }
1772 
1773  // Get the signal resonance from the stored vector
1774  LauAbsIncohRes* sigResonance = sigIncohResonances_[index];
1775 
1776  if (sigResonance == 0) {
1777  std::cerr<<"ERROR in LauIsobarDynamics::incohResAmp : Couldn't retrieve valid incoherent resonance with index = "<<index<<std::endl;
1778  incohInten_[index] = 0.0;
1779  return;
1780  }
1781 
1782  LauComplex ff = sigResonance->amplitude(kinematics_);
1783  incohInten_[index] = sigResonance->intensityFactor(kinematics_)*ff.abs2();
1784 
1785  // If we have a symmetrical Dalitz plot, flip the m_13^2 and m_23^2
1786  // variables, recalculate the dynamics, and combine both contributions.
1787  if (symmetricalDP_ == kTRUE) {
1789 
1790  ff = sigResonance->amplitude(kinematics_);
1791  incohInten_[index] += sigResonance->intensityFactor(kinematics_)*ff.abs2();
1792 
1793  // Flip the m_13^2 and m_23^2 variables back to their original values
1795  }
1796 
1797  // If we have a fully symmetric DP we need to calculate the amplitude at 6 points
1798  else if (fullySymmetricDP_ == kTRUE) {
1799  // rotate and evaluate
1801  ff = sigResonance->amplitude(kinematics_);
1802  incohInten_[index] += sigResonance->intensityFactor(kinematics_)*ff.abs2();
1803 
1804  // rotate and evaluate
1806  ff = sigResonance->amplitude(kinematics_);
1807  incohInten_[index] += sigResonance->intensityFactor(kinematics_)*ff.abs2();
1808 
1809  // rotate, flip and evaluate
1812  ff = sigResonance->amplitude(kinematics_);
1813  incohInten_[index] += sigResonance->intensityFactor(kinematics_)*ff.abs2();
1814 
1815  // rotate and evaluate
1817  ff = sigResonance->amplitude(kinematics_);
1818  incohInten_[index] += sigResonance->intensityFactor(kinematics_)*ff.abs2();
1819 
1820  // rotate and evaluate
1822  ff = sigResonance->amplitude(kinematics_);
1823  incohInten_[index] += sigResonance->intensityFactor(kinematics_)*ff.abs2();
1824 
1825  // rotate and flip to get us back to where we started
1828  }
1829 
1830 }
1831 
1832 void LauIsobarDynamics::setFFTerm(const UInt_t index, const Double_t realPart, const Double_t imagPart)
1833 {
1834  // Function to set the internal ff term (normally calculated using resAmp(index).
1835  if ( index >= nAmp_ ) {
1836  std::cerr<<"ERROR in LauIsobarDynamics::setFFTerm : index = "<<index<<" is not within the range 0 to "<<nAmp_-1<<std::endl;
1837  return;
1838  }
1839 
1840  ff_[index].setRealImagPart( realPart, imagPart );
1841 }
1842 
1843 void LauIsobarDynamics::setIncohIntenTerm(const UInt_t index, const Double_t value)
1844 {
1845  // Function to set the internal incoherent intensity term (normally calculated using incohResAmp(index).
1846  if ( index >= nIncohAmp_ ) {
1847  std::cerr<<"ERROR in LauIsobarDynamics::setFFTerm : index = "<<index<<" is not within the range 0 to "<<nIncohAmp_-1<<std::endl;
1848  return;
1849  }
1850 
1851  incohInten_[index] = value;
1852 }
1853 
1854 void LauIsobarDynamics::calcExtraInfo(const Bool_t init)
1855 {
1856  // This method calculates the fit fractions, mean efficiency and total DP rate
1857 
1858  Double_t fifjEffTot(0.0), fifjTot(0.0);
1859  UInt_t i, j;
1860  for (i = 0; i < nAmp_; i++) {
1861 
1862  // Calculate the diagonal terms
1863  TString name = "A"; name += i; name += "Sq_FitFrac";
1864  fitFrac_[i][i].name(name);
1865 
1866  name += "EffUnCorr";
1867  fitFracEffUnCorr_[i][i].name(name);
1868 
1869  Double_t fifjSumReal = fifjSum_[i][i].re();
1870  Double_t sumTerm = Amp_[i].abs2()*fifjSumReal*fNorm_[i]*fNorm_[i];
1871  fifjTot += sumTerm;
1872 
1873  Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
1874  Double_t sumEffTerm = Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
1875  fifjEffTot += sumEffTerm;
1876 
1877  fitFrac_[i][i].value(sumTerm);
1878  fitFracEffUnCorr_[i][i].value(sumEffTerm);
1879  }
1880 
1881  for (i = 0; i < nAmp_; i++) {
1882  for (j = i+1; j < nAmp_+nIncohAmp_; j++) {
1883  // Calculate the cross-terms
1884  TString name = "A"; name += i; name += "A"; name += j; name += "_FitFrac";
1885  fitFrac_[i][j].name(name);
1886 
1887  name += "EffUnCorr";
1888  fitFracEffUnCorr_[i][j].name(name);
1889 
1890  if(j >= nAmp_) {
1891  //Set off-diagonal incoherent terms to zero
1892  fitFrac_[i][j].value(0.);
1893  fitFracEffUnCorr_[i][j].value(0.);
1894  continue;
1895  }
1896 
1897  LauComplex AmpjConj = Amp_[j].conj();
1898  LauComplex AmpTerm = Amp_[i]*AmpjConj;
1899 
1900  Double_t crossTerm = 2.0*(AmpTerm*fifjSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1901  fifjTot += crossTerm;
1902 
1903  Double_t crossEffTerm = 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
1904  fifjEffTot += crossEffTerm;
1905 
1906  fitFrac_[i][j].value(crossTerm);
1907  fitFracEffUnCorr_[i][j].value(crossEffTerm);
1908  }
1909  }
1910 
1911  for (i = nAmp_; i < nAmp_+nIncohAmp_; i++) {
1912 
1913  // Calculate the incoherent terms
1914  TString name = "A"; name += i; name += "Sq_FitFrac";
1915  fitFrac_[i][i].name(name);
1916 
1917  name += "EffUnCorr";
1918  fitFracEffUnCorr_[i][i].name(name);
1919 
1920  Double_t sumTerm = Amp_[i].abs2()*fSqSum_[i]*fNorm_[i]*fNorm_[i];
1921  fifjTot += sumTerm;
1922 
1923  Double_t sumEffTerm = Amp_[i].abs2()*fSqEffSum_[i]*fNorm_[i]*fNorm_[i];
1924  fifjEffTot += sumEffTerm;
1925 
1926  fitFrac_[i][i].value(sumTerm);
1927  fitFracEffUnCorr_[i][i].value(sumEffTerm);
1928  }
1929 
1930  for (i = nAmp_; i < nAmp_+nIncohAmp_; i++) {
1931  for (j = i+1; j < nAmp_+nIncohAmp_; j++) {
1932  //Set off-diagonal incoherent terms to zero
1933  TString name = "A"; name += i; name += "A"; name += j; name += "_FitFrac";
1934  fitFrac_[i][j].name(name);
1935 
1936  name += "EffUnCorr";
1937  fitFracEffUnCorr_[i][j].name(name);
1938 
1939  fitFrac_[i][j].value(0.);
1940  fitFracEffUnCorr_[i][j].value(0.);
1941  }
1942  }
1943 
1944  if (TMath::Abs(fifjTot) > 1e-10) {
1945  meanDPEff_.value(fifjEffTot/fifjTot);
1946  if (init) {
1949  }
1950  }
1951  DPRate_.value(fifjTot);
1952  if (init) {
1955  }
1956 
1957  // Now divide the fitFraction sums by the overall integral
1958  for (i = 0; i < nAmp_+nIncohAmp_; i++) {
1959  for (j = i; j < nAmp_+nIncohAmp_; j++) {
1960  // Get the actual fractions by dividing by the total DP rate
1961  Double_t fitFracVal = fitFrac_[i][j].value();
1962  fitFracVal /= fifjTot;
1963  fitFrac_[i][j].value( fitFracVal );
1964 
1965  Double_t fitFracEffUnCorrVal = fitFracEffUnCorr_[i][j].value();
1966  fitFracEffUnCorrVal /= fifjEffTot;
1967  fitFracEffUnCorr_[i][j].value( fitFracEffUnCorrVal );
1968 
1969  if (init) {
1970  fitFrac_[i][j].genValue( fitFrac_[i][j].value() );
1971  fitFrac_[i][j].initValue( fitFrac_[i][j].value() );
1972  fitFracEffUnCorr_[i][j].genValue( fitFracEffUnCorr_[i][j].value() );
1973  fitFracEffUnCorr_[i][j].initValue( fitFracEffUnCorr_[i][j].value() );
1974  }
1975  }
1976  }
1977 
1978  // Work out total fit fraction over all K-matrix components (for each propagator)
1979  KMPropMap::iterator mapIter;
1980  Int_t propInt(0);
1981 
1982  for (mapIter = kMatrixPropagators_.begin(); mapIter != kMatrixPropagators_.end(); ++mapIter) {
1983 
1984  LauKMatrixPropagator* thePropagator = mapIter->second;
1985 
1986  TString propName = thePropagator->getName();
1987 
1988  // Now loop over all resonances and find those which are K-matrix components for this propagator
1989  Double_t kMatrixTotFitFrac(0.0);
1990 
1991  for (i = 0; i < nAmp_; i++) {
1992 
1993  Bool_t gotKMRes1 = this->gotKMatrixMatch(i, propName);
1994  if (gotKMRes1 == kFALSE) {continue;}
1995 
1996  Double_t fifjSumReal = fifjSum_[i][i].re();
1997  Double_t sumTerm = Amp_[i].abs2()*fifjSumReal*fNorm_[i]*fNorm_[i];
1998 
1999  //Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
2000  //Double_t sumEffTerm = Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
2001 
2002  kMatrixTotFitFrac += sumTerm;
2003 
2004  for (j = i+1; j < nAmp_; j++) {
2005 
2006  Bool_t gotKMRes2 = this->gotKMatrixMatch(j, propName);
2007  if (gotKMRes2 == kFALSE) {continue;}
2008 
2009  LauComplex AmpjConj = Amp_[j].conj();
2010  LauComplex AmpTerm = Amp_[i]*AmpjConj;
2011 
2012  Double_t crossTerm = 2.0*(AmpTerm*fifjSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
2013  //Double_t crossEffTerm = 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
2014 
2015  kMatrixTotFitFrac += crossTerm;
2016 
2017  }
2018 
2019  }
2020 
2021  kMatrixTotFitFrac /= fifjTot;
2022 
2023  TString parName("KMatrixTotFF_"); parName += propInt;
2024  extraParameters_[propInt].name( parName );
2025  extraParameters_[propInt].value(kMatrixTotFitFrac);
2026  if (init) {
2027  extraParameters_[propInt].genValue(kMatrixTotFitFrac);
2028  extraParameters_[propInt].initValue(kMatrixTotFitFrac);
2029  }
2030 
2031  std::cout<<"INFO in LauIsobarDynamics::calcExtraInfo : Total K-matrix fit fraction for propagator "<<propName<<" is "<<kMatrixTotFitFrac<<std::endl;
2032 
2033  ++propInt;
2034 
2035  }
2036 
2037 
2038 }
2039 
2040 Bool_t LauIsobarDynamics::gotKMatrixMatch(UInt_t resAmpInt, const TString& propName) const
2041 {
2042 
2043  Bool_t gotMatch(kFALSE);
2044 
2045  if (resAmpInt >= nAmp_) {return kFALSE;}
2046 
2047  const LauAbsResonance* theResonance = sigResonances_[resAmpInt];
2048 
2049  if (theResonance == 0) {return kFALSE;}
2050 
2051  Int_t resModelInt = theResonance->getResonanceModel();
2052 
2053  if (resModelInt == LauAbsResonance::KMatrix) {
2054 
2055  TString resName = theResonance->getResonanceName();
2056 
2057  KMStringMap::const_iterator kMPropSetIter = kMatrixPropSet_.find(resName);
2058 
2059  if (kMPropSetIter != kMatrixPropSet_.end()) {
2060  TString kmPropString = kMPropSetIter->second;
2061  if (kmPropString == propName) {gotMatch = kTRUE;}
2062  }
2063 
2064  }
2065 
2066  return gotMatch;
2067 
2068 }
2069 
2071 {
2072  // Calculate the normalisation for the log-likelihood function.
2073  DPNorm_ = 0.0;
2074 
2075  for (UInt_t i = 0; i < nAmp_; i++) {
2076  // fifjEffSum is the contribution from the term involving the resonance
2077  // dynamics (f_i for resonance i) and the efficiency term.
2078  Double_t fifjEffSumReal = fifjEffSum_[i][i].re();
2079  // We need to normalise this contribution w.r.t. the complete dynamics in the DP.
2080  // Hence we scale by the fNorm_i factor (squared), which is calculated by the
2081  // initialise() function, when the normalisation integrals are calculated and cached.
2082  // We also include the complex amplitude squared to get the total normalisation
2083  // contribution from this resonance.
2084  DPNorm_ += Amp_[i].abs2()*fifjEffSumReal*fNorm_[i]*fNorm_[i];
2085  }
2086 
2087  // We now come to the cross-terms (between resonances i and j) in the normalisation.
2088  for (UInt_t i = 0; i < nAmp_; i++) {
2089  for (UInt_t j = i+1; j < nAmp_; j++) {
2090  LauComplex AmpjConj = Amp_[j].conj();
2091  LauComplex AmpTerm = Amp_[i]*AmpjConj;
2092  // Again, fifjEffSum is the contribution from the term involving the resonance
2093  // dynamics (f_i*f_j_conjugate) and the efficiency cross term.
2094  // Also include the relative normalisation between these two resonances w.r.t. the
2095  // total DP dynamical structure (fNorm_i and fNorm_j) and the complex
2096  // amplitude squared (mag,phase) part.
2097  DPNorm_ += 2.0*(AmpTerm*fifjEffSum_[i][j]).re()*fNorm_[i]*fNorm_[j];
2098  }
2099  }
2100 
2101  for (UInt_t i = 0; i < nIncohAmp_; i++) {
2102  DPNorm_ += Amp_[i+nAmp_].abs2()*fSqEffSum_[i+nAmp_]*fNorm_[i+nAmp_]*fNorm_[i+nAmp_];
2103  }
2104 
2105  return DPNorm_;
2106 }
2107 
2109 {
2110  // Routine to generate a signal event according to the Dalitz plot
2111  // model we have defined.
2112 
2113  // We need to make sure to calculate everything for every resonance
2114  integralsToBeCalculated_.clear();
2115  for ( UInt_t i(0); i < nAmp_+nIncohAmp_; ++i ) {
2116  integralsToBeCalculated_.insert(i);
2117  }
2118 
2119  nSigGenLoop_ = 0;
2120  Bool_t generatedSig(kFALSE);
2121 
2122  while (generatedSig == kFALSE && nSigGenLoop_ < iterationsMax_) {
2123 
2124  // Generates uniform DP phase-space distribution
2125  Double_t m13Sq(0.0), m23Sq(0.0);
2126  kinematics_->genFlatPhaseSpace(m13Sq, m23Sq);
2127 
2128  // If we're in a symmetrical DP then we should only generate events in one half
2129  // TODO - what do we do for fully symmetric?
2130  if ( symmetricalDP_ && !fullySymmetricDP_ && m13Sq > m23Sq ) {
2131  Double_t tmpSq = m13Sq;
2132  m13Sq = m23Sq;
2133  m23Sq = tmpSq;
2134  }
2135 
2136  // Calculate the amplitudes and total amplitude for the given DP point
2137  this->calcLikelihoodInfo(m13Sq, m23Sq);
2138 
2139 
2140  // Throw the random number and check it against the ratio of ASq and the accept/reject ceiling
2141  const Double_t randNo = LauRandom::randomFun()->Rndm();
2142  if (randNo > ASq_/aSqMaxSet_) {
2143  ++nSigGenLoop_;
2144  } else {
2145  generatedSig = kTRUE;
2146  nSigGenLoop_ = 0;
2147 
2148  // Keep a note of the maximum ASq that we've found
2149  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
2150  }
2151 
2152  } // while loop
2153 
2154  // Check that all is well with the generation
2155  Bool_t sigGenOK(kTRUE);
2156  if (GenOK != this->checkToyMC(kTRUE,kFALSE)) {
2157  sigGenOK = kFALSE;
2158  }
2159 
2160  return sigGenOK;
2161 }
2162 
2163 LauIsobarDynamics::ToyMCStatus LauIsobarDynamics::checkToyMC(Bool_t printErrorMessages, Bool_t printInfoMessages)
2164 {
2165  // Check whether we have generated the toy MC OK.
2166  ToyMCStatus ok(GenOK);
2167 
2168  if (nSigGenLoop_ >= iterationsMax_) {
2169  // Exceeded maximum allowed iterations - the generation is too inefficient
2170  if (printErrorMessages) {
2171  std::cerr<<"WARNING in LauIsobarDynamics::checkToyMC : More than "<<iterationsMax_<<" iterations performed and no event accepted."<<std::endl;
2172  }
2173 
2174  if ( aSqMaxSet_ > 1.01 * aSqMaxVar_ ) {
2175  if (printErrorMessages) {
2176  std::cerr<<" : |A|^2 maximum was set to "<<aSqMaxSet_<<" but this appears to be too high."<<std::endl;
2177  std::cerr<<" : Maximum value of |A|^2 found so far = "<<aSqMaxVar_<<std::endl;
2178  std::cerr<<" : The value of |A|^2 maximum will be decreased and the generation restarted."<<std::endl;
2179  }
2180  aSqMaxSet_ = 1.01 * aSqMaxVar_;
2181  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : |A|^2 max reset to "<<aSqMaxSet_<<std::endl;
2182  ok = MaxIterError;
2183  } else {
2184  if (printErrorMessages) {
2185  std::cerr<<" : |A|^2 maximum was set to "<<aSqMaxSet_<<", which seems to be correct for the given model."<<std::endl;
2186  std::cerr<<" : However, the generation is very inefficient - please check your model."<<std::endl;
2187  std::cerr<<" : The maximum number of iterations will be increased and the generation restarted."<<std::endl;
2188  }
2189  iterationsMax_ *= 2;
2190  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : max number of iterations reset to "<<iterationsMax_<<std::endl;
2191  ok = MaxIterError;
2192  }
2193  } else if (aSqMaxVar_ > aSqMaxSet_) {
2194  // Found a value of ASq higher than the accept/reject ceiling - the generation is biased
2195  if (printErrorMessages) {
2196  std::cerr<<"WARNING in LauIsobarDynamics::checkToyMC : |A|^2 maximum was set to "<<aSqMaxSet_<<" but a value exceeding this was found: "<<aSqMaxVar_<<std::endl;
2197  std::cerr<<" : Run was invalid, as any generated MC will be biased, according to the accept/reject method!"<<std::endl;
2198  std::cerr<<" : The value of |A|^2 maximum be reset to be > "<<aSqMaxVar_<<" and the generation restarted."<<std::endl;
2199  }
2200  aSqMaxSet_ = 1.01 * aSqMaxVar_;
2201  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : |A|^2 max reset to "<<aSqMaxSet_<<std::endl;
2202  ok = ASqMaxError;
2203  } else if (printInfoMessages) {
2204  std::cout<<"INFO in LauIsobarDynamics::checkToyMC : aSqMaxSet = "<<aSqMaxSet_<<" and aSqMaxVar = "<<aSqMaxVar_<<std::endl;
2205  }
2206 
2207  return ok;
2208 }
2209 
2211 {
2212  // Retrieve the data for event iEvt
2213  if (data_.size() > iEvt) {
2214  currentEvent_ = data_[iEvt];
2215  } else {
2216  std::cerr<<"ERROR in LauIsobarDynamics::setDataEventNo : Event index too large: "<<iEvt<<" >= "<<data_.size()<<"."<<std::endl;
2217  }
2218 
2225  scfFraction_ = currentEvent_->retrieveScfFraction(); // These two are necessary, even though the dynamics don't actually use scfFraction_ or jacobian_,
2226  jacobian_ = currentEvent_->retrieveJacobian(); // since this is at the heart of the caching mechanism.
2227 }
2228 
2230 {
2231  // Calculate the likelihood and associated info
2232  // for the given event using cached information
2233  evtLike_ = 0.0;
2234 
2235  // retrieve the cached dynamics from the tree:
2236  // realAmp, imagAmp for each resonance plus efficiency, scf fraction and jacobian
2237  this->setDataEventNo(iEvt);
2238 
2239  // use realAmp and imagAmp to create the resonance amplitudes
2240  const std::vector<Double_t>& realAmp = currentEvent_->retrieveRealAmp();
2241  const std::vector<Double_t>& imagAmp = currentEvent_->retrieveImagAmp();
2242  const std::vector<Double_t>& incohInten = currentEvent_->retrieveIncohIntensities();
2243  for (UInt_t i = 0; i < nAmp_; i++) {
2244  this->setFFTerm(i, realAmp[i], imagAmp[i]);
2245  }
2246  for (UInt_t i = 0; i < nIncohAmp_; i++) {
2247  this->setIncohIntenTerm(i, incohInten[i]);
2248  }
2249 
2250  // Update the dynamics - calculates totAmp_ and then ASq_ = totAmp_.abs2() * eff_
2251  // All calculated using cached information on the individual amplitudes and efficiency.
2252  this->calcTotalAmp(kTRUE);
2253 
2254  // Calculate the normalised matrix element squared value
2255  if (DPNorm_ > 1e-10) {
2256  evtLike_ = ASq_/DPNorm_;
2257  }
2258 }
2259 
2260 void LauIsobarDynamics::calcLikelihoodInfo(const Double_t m13Sq, const Double_t m23Sq)
2261 {
2262  this->calcLikelihoodInfo(m13Sq, m23Sq, -1);
2263 }
2264 
2265 void LauIsobarDynamics::calcLikelihoodInfo(const Double_t m13Sq, const Double_t m23Sq, const Int_t tagCat)
2266 {
2267  // Calculate the likelihood and associated info
2268  // for the given point in the Dalitz plot
2269  // Also retrieves the SCF fraction in the bin where the event lies (done
2270  // here to cache it along with the the rest of the DP quantities, like eff)
2271  // The jacobian for the square DP is calculated here for the same reason.
2272  evtLike_ = 0.0;
2273 
2274  // update the kinematics for the specified DP point
2275  kinematics_->updateKinematics(m13Sq, m23Sq);
2276 
2277  // calculate the jacobian and the scfFraction to cache them later
2278  scfFraction_ = this->retrieveScfFraction(tagCat);
2279  if (kinematics_->squareDP() == kTRUE) {
2281  }
2282 
2283  // calculate the ff_ terms and retrieves eff_ from the efficiency model
2284  this->calculateAmplitudes();
2285  // then calculate totAmp_ and finally ASq_ = totAmp_.abs2() * eff_
2286  this->calcTotalAmp(kTRUE);
2287 
2288  // Calculate the normalised matrix element squared value
2289  if (DPNorm_ > 1e-10) {
2290  evtLike_ = ASq_/DPNorm_;
2291  }
2292 }
2293 
2295 {
2296  if ( recalcNormalisation_ == kFALSE ) {
2297  return;
2298  }
2299 
2300  const UInt_t nEvents = data_.size();
2301 
2302  std::set<UInt_t>::const_iterator iter = integralsToBeCalculated_.begin();
2303  const std::set<UInt_t>::const_iterator intEnd = integralsToBeCalculated_.end();
2304 
2305  for (UInt_t iEvt = 0; iEvt < nEvents; ++iEvt) {
2306 
2307  currentEvent_ = data_[iEvt];
2308 
2309  std::vector<Double_t>& realAmp = currentEvent_->retrieveRealAmp();
2310  std::vector<Double_t>& imagAmp = currentEvent_->retrieveImagAmp();
2311  std::vector<Double_t>& incohInten = currentEvent_->retrieveIncohIntensities();
2312 
2313  const Double_t m13Sq = currentEvent_->retrievem13Sq();
2314  const Double_t m23Sq = currentEvent_->retrievem23Sq();
2315  const Int_t tagCat = currentEvent_->retrieveTagCat();
2316 
2317  this->calcLikelihoodInfo(m13Sq, m23Sq, tagCat);
2318 
2319  for ( iter = integralsToBeCalculated_.begin(); iter != intEnd; ++iter) {
2320  const UInt_t i = *iter;
2321  if(*iter < nAmp_) {
2322  realAmp[i] = ff_[i].re();
2323  imagAmp[i] = ff_[i].im();
2324  } else {
2325  incohInten[i-nAmp_] = incohInten_[i-nAmp_];
2326  }
2327  }
2328  }
2329 }
2330 
2332 {
2333  // In LauFitDataTree, the first two variables should always be m13^2 and m23^2.
2334  // Other variables follow thus: charge/flavour tag prob, etc.
2335 
2336  // Since this is the first caching, we need to make sure to calculate everything for every resonance
2337  integralsToBeCalculated_.clear();
2338  for ( UInt_t i(0); i < nAmp_+nIncohAmp_; ++i ) {
2339  integralsToBeCalculated_.insert(i);
2340  }
2341 
2342  UInt_t nBranches = inputFitTree.nBranches();
2343 
2344  if (nBranches < 2) {
2345  std::cerr<<"ERROR in LauIsobarDynamics::fillDataTree : Expecting at least 2 variables " <<"in input data tree, but there are "<<nBranches<<"!\n";
2346  std::cerr<<" : Make sure you have the right number of variables in your input data file!"<<std::endl;
2347  gSystem->Exit(EXIT_FAILURE);
2348  }
2349 
2350  // Data structure that will cache the variables required to
2351  // calculate the signal likelihood for this experiment
2352  for ( std::vector<LauCacheData*>::iterator iter = data_.begin(); iter != data_.end(); ++iter ) {
2353  delete (*iter);
2354  }
2355  data_.clear();
2356 
2357  Double_t m13Sq(0.0), m23Sq(0.0);
2358  Double_t mPrime(0.0), thPrime(0.0);
2359  Int_t tagCat(-1);
2360  std::vector<Double_t> realAmp(nAmp_), imagAmp(nAmp_);
2361  Double_t eff(0.0), scfFraction(0.0), jacobian(0.0);
2362 
2363  UInt_t nEvents = inputFitTree.nEvents() + inputFitTree.nFakeEvents();
2364 
2365  data_.reserve(nEvents);
2366 
2367  for (UInt_t iEvt = 0; iEvt < nEvents; ++iEvt) {
2368 
2369  const LauFitData& dataValues = inputFitTree.getData(iEvt);
2370  LauFitData::const_iterator iter = dataValues.find("m13Sq");
2371  m13Sq = iter->second;
2372  iter = dataValues.find("m23Sq");
2373  m23Sq = iter->second;
2374 
2375  // is there more than one tagging category?
2376  // if so then we need to know the category from the data
2377  if (scfFractionModel_.size()>1) {
2378  iter = dataValues.find("tagCat");
2379  tagCat = static_cast<Int_t>(iter->second);
2380  }
2381 
2382  // calculates the amplitudes and total amplitude for the given DP point
2383  // tagging category not needed by dynamics, but to find out the scfFraction
2384  this->calcLikelihoodInfo(m13Sq, m23Sq, tagCat);
2385 
2386  // extract the real and imaginary parts of the ff_ terms for storage
2387  for (UInt_t i = 0; i < nAmp_; i++) {
2388  realAmp[i] = ff_[i].re();
2389  imagAmp[i] = ff_[i].im();
2390  }
2391 
2392  if ( kinematics_->squareDP() ) {
2393  mPrime = kinematics_->getmPrime();
2394  thPrime = kinematics_->getThetaPrime();
2395  }
2396 
2397  eff = this->getEvtEff();
2398  scfFraction = this->getEvtScfFraction();
2399  jacobian = this->getEvtJacobian();
2400 
2401  // store the data for each event in the list
2402  data_.push_back( new LauCacheData() );
2403  data_[iEvt]->storem13Sq(m13Sq);
2404  data_[iEvt]->storem23Sq(m23Sq);
2405  data_[iEvt]->storemPrime(mPrime);
2406  data_[iEvt]->storethPrime(thPrime);
2407  data_[iEvt]->storeTagCat(tagCat);
2408  data_[iEvt]->storeEff(eff);
2409  data_[iEvt]->storeScfFraction(scfFraction);
2410  data_[iEvt]->storeJacobian(jacobian);
2411  data_[iEvt]->storeRealAmp(realAmp);
2412  data_[iEvt]->storeImagAmp(imagAmp);
2413  data_[iEvt]->storeIncohIntensities(incohInten_);
2414  }
2415 }
2416 
2418 {
2419  // Select the event (kinematics_) using an accept/reject method based on the
2420  // ratio of the current value of ASq to the maximal value.
2421  Bool_t accepted(kFALSE);
2422 
2423  // calculate the ff_ terms and retrieves eff_ from the efficiency model
2424  this->calculateAmplitudes();
2425  // then calculate totAmp_ and finally ASq_ = totAmp_.abs2() (without the efficiency correction!)
2426  this->calcTotalAmp(kFALSE);
2427 
2428  // Compare the ASq value with the maximal value (set by the user)
2429  if (LauRandom::randomFun()->Rndm() < ASq_/aSqMaxSet_) {
2430  accepted = kTRUE;
2431  }
2432 
2433  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
2434 
2435  return accepted;
2436 }
2437 
2439 {
2440  // calculate the ff_ terms and retrieves eff_ from the efficiency model
2441  this->calculateAmplitudes();
2442  // then calculate totAmp_ and finally ASq_ = totAmp_.abs2() (without the efficiency correction!)
2443  this->calcTotalAmp(kFALSE);
2444 
2445  if (ASq_ > aSqMaxVar_) {aSqMaxVar_ = ASq_;}
2446 
2447  // return the event weight = the value of the squared amplitude divided
2448  // by the user-defined ceiling
2449  return ASq_ / aSqMaxSet_;
2450 }
2451 
2452 void LauIsobarDynamics::updateCoeffs(const std::vector<LauComplex>& coeffs)
2453 {
2454  // Check that the number of coeffs is correct
2455  if (coeffs.size() != this->getnTotAmp()) {
2456  std::cerr << "ERROR in LauIsobarDynamics::updateCoeffs : Expected " << this->getnTotAmp() << " but got " << coeffs.size() << std::endl;
2457  gSystem->Exit(EXIT_FAILURE);
2458  }
2459 
2460  // Now check if the coeffs have changed
2461  Bool_t changed = (Amp_ != coeffs);
2462  if (changed) {
2463  // Copy the coeffs
2464  Amp_ = coeffs;
2465  }
2466 
2467  // TODO should perhaps keep track of whether the resonance parameters have changed here and if none of those and none of the coeffs have changed then we don't need to update the norm
2468 
2469  // Update the total normalisation for the signal likelihood
2470  this->calcSigDPNorm();
2471 }
2472 
2473 TString LauIsobarDynamics::getConjResName(const TString& resName) const
2474 {
2475  // Get the name of the charge conjugate resonance
2476  TString conjName(resName);
2477 
2478  Ssiz_t index1 = resName.Index("+");
2479  Ssiz_t index2 = resName.Index("-");
2480  if (index1 != -1) {
2481  conjName.Replace(index1, 1, "-");
2482  } else if (index2 != -1) {
2483  conjName.Replace(index2, 1, "+");
2484  }
2485 
2486  return conjName;
2487 
2488 }
2489 
2491 {
2492  Double_t eff(1.0);
2493  if (effModel_ != 0) {
2495  }
2496  return eff;
2497 }
2498 
2500 {
2501  Double_t scfFraction(0.0);
2502 
2503  // scf model and eff model are exactly the same, functionally
2504  // so we use an instance of LauAbsEffModel, and the method
2505  // calcEfficiency actually calculates the scf fraction
2506  if (tagCat == -1) {
2507  if (!scfFractionModel_.empty()) {
2508  scfFraction = scfFractionModel_[0]->calcEfficiency(kinematics_);
2509  }
2510  } else {
2511  scfFraction = scfFractionModel_[tagCat]->calcEfficiency(kinematics_);
2512  }
2513  return scfFraction;
2514 }
2515 
KMStringMap kMatrixPropSet_
The names of the M-matrix components in the model mapped to their propagators.
std::vector< LauComplex > Amp_
The complex coefficients for the amplitude components.
void modifyDataTree()
Recache the amplitude values for those that have changed.
Int_t nSigGenLoop_
The number of unsucessful attempts to generate an event so far.
TString getNameDaug1() const
Get name of the first daughter particle.
Double_t getm23Max() const
Get the m23 maximum defined as (mParent - m1)
Double_t retrievemPrime() const
Retrieve the square Dalitz plot coordinate, m&#39;.
Bool_t squareDP() const
Are the square Dalitz plot co-ordinates being calculated?
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
File containing declaration of LauKMatrixPropFactory class.
File containing declaration of LauFitDataTree class.
File containing declaration of LauNRAmplitude class.
void setIntegralBinWidths(const Double_t m13BinWidth, const Double_t m23BinWidth, const Double_t mPrimeBinWidth=0.001, const Double_t thPrimeBinWidth=0.001)
Set the widths of the bins to use when integrating across the Dalitz plot or square Dalitz plot...
Bool_t recalcNormalisation_
Flag to recalculate the normalisation.
void storeIntensity(const UInt_t m13Point, const UInt_t m23Point, const UInt_t iAmp, const Double_t intensity)
Store the intensity for the given grid point and intensity index.
Double_t retrievem23Sq() const
Retrieve the invariant mass squared of the second and third daugthers.
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 updateCoeffs(const std::vector< LauComplex > &coeffs)
Update the complex coefficients for the resonances.
void updateSqDPKinematics(Double_t mPrime, Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m&#39; and theta&#39;.
std::set< UInt_t > integralsToBeCalculated_
Resonance indices for which the amplitudes and integrals should be recalculated.
Double_t m13BinWidth_
The bin width to use when integrating over m13.
Bool_t flipHelicity_
The helicity flip flag for new amplitude components.
File containing declaration of LauAbsEffModel class.
ClassImp(LauAbsCoeffSet)
std::vector< Int_t > resPairAmp_
The index of the daughter not produced by the resonance for each amplitude component.
Bool_t hasResonance(const TString &resName) const
Check whether this model includes a named resonance.
const LauAbsResonance * findResonance(const TString &resName) const
Retrieve the named resonance.
std::vector< LauParameter > extraParameters_
any extra parameters/quantities (e.g. K-matrix total fit fractions)
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.
Double_t getEvtJacobian() const
Retrieve the Jacobian, for the transformation into square DP coordinates, for the current event...
File containing declaration of LauBelleNR class.
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.
std::vector< Double_t > fSqEffSum_
The event-by-event running total of the dynamical amplitude squared for each amplitude component...
File containing declaration of LauKMatrixProdSVP class.
void calcDPNormalisationScheme()
Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
LauAbsEffModel * effModel_
The efficiency model across the Dalitz plot.
void calcExtraInfo(const Bool_t init=kFALSE)
Calculate the fit fractions, mean efficiency and total DP rate.
std::vector< Int_t > incohResPairAmp_
The index of the daughter not produced by the resonance for each incoherent amplitude component...
Double_t DPNorm_
The normalisation factor for the log-likelihood function.
LauTagCatScfFractionModelMap scfFractionModel_
The self cross feed fraction models across the Dalitz plot.
LauParArray fitFrac_
The fit fractions for the amplitude components.
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...
Abstract class for defining incoherent resonant amplitude models.
Double_t m13Sq_
The invariant mass squared of the first and third daughters.
LauAbsResonance * getResonance(const LauDaughters *daughters, const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory=LauBlattWeisskopfFactor::Default, const LauBlattWeisskopfFactor::BarrierType bwType=LauBlattWeisskopfFactor::BWPrimeBarrier)
Create a resonance.
Int_t iterationsMax_
The maximum allowed number of attempts when generating an event.
LauComplex totAmp_
The total amplitude for the current event.
TString getConjResName(const TString &resName) const
Retrieve the name of the charge conjugate of a named resonance.
Double_t aSqMaxSet_
The maximum allowed value of A squared.
std::vector< Double_t > resonanceParValues_
List of floating resonance parameter values from previous calculation.
File containing declaration of LauBelleSymNR class.
Bool_t symmetricalDP_
Whether the Dalitz plot is symmetrical.
File containing declaration of LauPrint class.
Pure abstract base class for defining the efficiency description across the signal Dalitz plot...
void setIncohIntenTerm(const UInt_t index, const Double_t value)
Set the dynamic part of the intensity for a given incoherent amplitude component at the current point...
Int_t getNChannels() const
Get the number of channels.
Double_t getm13Min() const
Get the m13 minimum defined as (m1 + m3)
Class for defining (a section of) the Dalitz plot integration binning scheme.
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
std::map< TString, Double_t > LauFitData
Type for holding event data.
void calcDPNormalisation()
Calculate the Dalitz plot normalisation integrals across the whole Dalitz plot.
TString getNameDaug2() const
Get name of the second daughter particle.
Double_t getM23Value(const UInt_t m23Point) const
Retrieve the m23 value at the given grid point.
void recalculateNormalisation()
recalculate Normalization
std::vector< TString > incohResTypAmp_
The resonance types of all of the incoherent amplitude components.
std::vector< LauAbsIncohRes * > sigIncohResonances_
The incoherent resonances in the model.
UInt_t nBranches() const
Obtain the number of branches in the tree.
Double_t getEvtScfFraction() const
Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the ...
void resAmp(const UInt_t index)
Calculate the dynamic part of the amplitude for a given component at the current point in the Dalitz ...
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.
LauParameter DPRate_
The overall Dalitz plot rate.
void calcLikelihoodInfo(const UInt_t iEvt)
Calculate the likelihood (and all associated information) for the given event number.
ToyMCStatus
The possible statuses for toy MC generation.
KMPropMap kMatrixPropagators_
The K-matrix propagators.
std::vector< Double_t > fSqSum_
The event-by-event running total of the dynamical amplitude squared for each amplitude component...
Bool_t flipHelicity() const
Get the helicity flip flag.
File containing declaration of LauKinematics class.
UInt_t getnm23Points() const
Retrieve the number of bins in m23.
File containing declaration of LauKMatrixPropagator class.
Double_t getm13Max() const
Get the m13 maximum defined as (mParent - m2)
UInt_t getnm13Points() const
Retrieve the number of bins in m13.
Class for defining a K-matrix production pole amplitude term.
Double_t getm12Min() const
Get the m12 minimum defined as (m1 + m2)
const LauAbsResonance * getResonance(const UInt_t resIndex) const
Retrieve a resonance by its index.
UInt_t nAmp_
The number of amplitude components.
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...
Singleton factory class for creating resonances.
void updateKinematics(Double_t m13Sq, Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
void initSummary()
Print a summary of the model to be used.
Bool_t fullySymmetricDP_
Whether the Dalitz plot is fully symmetric.
void fillDataTree(const LauFitDataTree &fitDataTree)
Fill the internal data structure that caches the resonance dynamics.
void resetNormVectors()
Zero the various values used to store integrals info.
Double_t ASq_
The value of A squared for the current event.
Double_t calcSqDPJacobian(const Double_t mPrime, const Double_t thPrime) const
Calculate the Jacobian for the transformation m23^2, m13^2 -&gt; m&#39;, theta&#39; (square DP) at the given poin...
Double_t abs2() const
Obtain the square of the absolute value of the complex number.
Definition: LauComplex.hh:232
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:66
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.
TString intFileName_
The name of the file to save integrals to.
Int_t getResPairAmpInt() const
Get the DP axis identifier.
File containing declaration of LauCacheData class.
Bool_t getSquareDP() const
Retrieve the square DP flag.
virtual LauResonanceModel getResonanceModel() const =0
Get the resonance model type.
Double_t getm12Max() const
Get the m12 maximum defined as (mParent - m3)
void addKMatrixProdPole(const TString &poleName, const TString &propName, Int_t poleIndex)
Add a K-matrix production pole term to the model.
void calcTotalAmp(const Bool_t useEff)
Calculate the total Dalitz plot amplitude at the current point in the Dalitz plot.
Double_t retrieveScfFraction(Int_t tagCat)
Obtain the self cross feed fraction of the current event from the model.
std::vector< TString > resTypAmp_
The resonance types of all of the amplitude components.
const std::vector< Double_t > & retrieveIncohIntensities() const
Retrieve the incoherent intensities.
const std::vector< Double_t > & retrieveRealAmp() const
Retrieve the real parts of the amplitudes.
Int_t retrieveTagCat() const
Retrieve the tagging category.
LauIsobarDynamics(LauDaughters *daughters, LauAbsEffModel *effModel, LauAbsEffModel *scfFractionModel=0)
Constructor.
Int_t tagCat_
The tagging category.
LauCacheData * currentEvent_
The cached data for the current event.
void initialiseVectors()
Initialise the internal storage for this model.
Double_t jacobian_
The Jacobian, for the transformation into square DP coordinates at the current point in the Dalitz pl...
LauDaughters * daughters_
The daughters of the decay.
Bool_t withinSqDPLimits(Double_t mPrime, Double_t thetaPrime) const
Check whether a given (m&#39;,theta&#39;) point is within the kinematic limits of the Dalitz plot...
void flipAndUpdateKinematics()
Flips the DP variables m13^2 &lt;-&gt; m23^2 and recalculates all kinematic quantities. ...
std::map< Int_t, LauAbsEffModel * > LauTagCatScfFractionModelMap
The type used for containing multiple self cross feed fraction models for different categories (e...
Double_t narrowWidth_
The value below which a resonance width is considered to be narrow.
Double_t getM13Value(const UInt_t m13Point) const
Retrieve the m13 value at the given grid point.
void incohResAmp(const UInt_t index)
Calculate the dynamic part of the intensity for a given incoherent component at the current point in ...
File containing declaration of LauDPPartialIntegralInfo class.
std::vector< std::vector< LauComplex > > fifjSum_
The event-by-event running total of the amplitude cross terms for each pair of amplitude components...
Double_t getm23Min() const
Get the m23 minimum defined as (m2 + m3)
Double_t retrieveScfFraction() const
Retrieve the fraction of poorly reconstructed events (the so-called self cross feed fraction) ...
Double_t retrievethPrime() const
Retrieve the square Dalitz plot coordinate, theta&#39;.
UInt_t getnTotAmp() const
Retrieve the total number of amplitude components.
virtual Double_t calcEfficiency(const LauKinematics *kinematics) const =0
Determine the efficiency for a given point in the Dalitz plot.
Double_t calcSigDPNorm()
Calculate the normalisation factor for the log-likelihood function.
Bool_t withinDPLimits(Double_t m13Sq, Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot...
Int_t getTypeDaug3() const
Get PDG code of the third daughter particle.
void addGridPointToIntegrals(const Double_t weight)
Add the amplitude values (with the appropriate weight) at the current grid point to the running integ...
Double_t thPrimeBinWidth_
The bin width to use when integrating over thetaPrime.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
File containing LauRandom namespace.
Bool_t generate()
Generate a toy MC signal event.
LauAbsResonance * addResonance(const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory=LauBlattWeisskopfFactor::Default, const LauBlattWeisskopfFactor::BarrierType bwType=LauBlattWeisskopfFactor::BWPrimeBarrier)
Add a resonance to the Dalitz plot.
void rotateAndUpdateKinematics()
Cyclically rotates the DP variables (m12 -&gt; m23, m23 -&gt; m13, m13 -&gt; m12) and recalculates all kinemat...
void removeCharge(TString &string) const
Remove the charge from the given particle name.
std::vector< Int_t > typDaug_
The PDG codes of the daughters.
LauResonanceModel
Define the allowed resonance types.
Double_t eff_
The efficiency at the current point in the Dalitz plot.
void calculateAmplitudes()
Calculate the amplitudes for all resonances for the current kinematics.
Int_t getTypeDaug2() const
Get PDG code of the second daughter particle.
BarrierType
Define the allowed types of barrier factors.
Double_t getWeight(const UInt_t m13Point, const UInt_t m23Point) const
Retrieve the weight for the given grid point.
UInt_t nFakeEvents() const
Retrieve the number of fake events.
Int_t resonanceIndex(const TString &resName) const
Retrieve the index for the given resonance.
Bool_t gotReweightedEvent()
Calculates whether an event with the current kinematics should be accepted in order to produce a dist...
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
TString getNameParent() const
Get name of the parent particle.
Double_t binningFactor_
The factor relating the width of the narrowest resonance and the binning size.
void addKMatrixProdSVP(const TString &SVPName, const TString &propName, Int_t channelIndex)
Add a K-matrix slowly-varying part (SVP) term to the model.
std::vector< LauCacheData * > data_
The cached data for all events.
Bool_t gotFullySymmetricDP() const
Is Dalitz plot fully symmetric, i.e. 3 identical particles.
Definition: LauDaughters.hh:72
void rescale(Double_t scaleVal)
Scale this by a factor.
Definition: LauComplex.hh:285
Int_t getChargeParent() const
Get charge of the parent particle.
void findIntegralsToBeRecalculated()
Determine which amplitudes and integrals need to be recalculated.
void zero()
Set both real and imaginary part to zero.
Definition: LauComplex.hh:321
Double_t initValue() const
The initial value of the parameter.
Class for defining signal dynamics using the isobar model.
Double_t getEventWeight()
Calculate the acceptance rate, for events with the current kinematics, when generating events accordi...
LauAbsResonance * addIncoherentResonance(const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType)
Add an incoherent resonance to the Dalitz plot.
File containing LauConstants namespace.
void writeIntegralsFile()
Write the results of the integrals (and related information) to a file.
File containing declaration of LauAbsResonance class.
void calcDPPartialIntegral(LauDPPartialIntegralInfo *intInfo)
Calculate the Dalitz plot normalisation integrals over a given range.
Double_t mPrimeBinWidth_
The bin width to use when integrating over mPrime.
ToyMCStatus checkToyMC(Bool_t printErrorMessages=kTRUE, Bool_t printInfoMessages=kFALSE)
Check the status of the toy MC generation.
LauKinematics * kinematics_
The kinematics of the decay.
UInt_t getnCohAmp() const
Retrieve the number of coherent amplitude components.
Double_t getThetaPrime() const
Get theta&#39; value.
Class for defining a complex number.
Definition: LauComplex.hh:47
Double_t getEvtEff() const
Retrieve the efficiency for the current event.
UInt_t nIncohAmp_
The number of incoherent amplitude components.
Bool_t gotKMatrixMatch(UInt_t resAmpInt, const TString &propName) const
Check whether a resonance is a K-matrix component of a given propagator.
std::vector< LauDPPartialIntegralInfo * > dpPartialIntegralInfo_
The storage of the integration scheme.
Double_t getIntensity(const UInt_t m13Point, const UInt_t m23Point, const UInt_t iAmp) const
Retrieve the intensity for the given grid point and intensity index.
void genFlatPhaseSpace(Double_t &m13Sq, Double_t &m23Sq) const
Routine to generate events flat in phase-space.
TString getName() const
Get the propagator name.
void setDataEventNo(UInt_t iEvt)
Load the data for a given event.
Double_t value() const
The value of the parameter.
Bool_t normalizationSchemeDone_
Whether the scheme for the integration has been determined.
UInt_t nEvents() const
Retrieve the number of events.
std::vector< LauParameter * > resonancePars_
List of floating resonance parameters.
std::vector< std::vector< UInt_t > > resonanceParResIndex_
Indices in sigResonances_ to point to the corresponding signal resonance(s) for each floating paramet...
BlattWeisskopfCategory
Define resonance categories that will share common barrier factor radii.
void storeEfficiency(const UInt_t m13Point, const UInt_t m23Point, const Double_t efficiency)
Store the efficiency for the given grid point.
const LauComplex & getAmplitude(const UInt_t m13Point, const UInt_t m23Point, const UInt_t iAmp) const
Retrieve the amplitude for the given grid point and amplitude index.
LauParameter meanDPEff_
The mean efficiency across the Dalitz plot.
virtual Bool_t preSymmetrised() const
Is the amplitude pre-symmetrised?
static LauResonanceMaker & get()
Get the factory instance.
Bool_t integralsDone_
Whether the integrals have been performed.
std::vector< Double_t > incohInten_
The dynamic part of the intensity for each incoherent amplitude component at the current point in the...
virtual LauComplex amplitude(const LauKinematics *kinematics)
Calculate the complex amplitude.
void storeAmplitude(const UInt_t m13Point, const UInt_t m23Point, const UInt_t iAmp, const LauComplex &amplitude)
Store the amplitude for the given grid point and amplitude index.
std::vector< LauComplex > ff_
The dynamic part of the amplitude for each amplitude component at the current point in the Dalitz plo...
static bool isIncoherentModel(LauResonanceModel model)
Is the resonance model incoherent?
void initialise(const std::vector< LauComplex > &coeffs)
Initialise the Dalitz plot dynamics.
TString getNameDaug3() const
Get name of the third daughter particle.
Double_t genValue() const
The value generated for the parameter.
Double_t mPrime_
The square Dalitz plot coordinate, m&#39;.
std::vector< Double_t > fNorm_
The normalisation factors for the dynamic parts of the amplitude for each amplitude component...
Class for defining a K-matrix production &quot;slowly-varying part&quot; (SVP) amplitude.
Class to store the input fit variables.
Double_t retrieveEff() const
Retrieve the efficiency.
File containing declaration of LauAbsIncohRes class.
std::vector< LauAbsResonance * > sigResonances_
The resonances in the model.
virtual Double_t intensityFactor(const LauKinematics *kinematics)=0
Get intensity factor.
void setFFTerm(const UInt_t index, const Double_t realPart, const Double_t imagPart)
Set the dynamic part of the amplitude for a given amplitude component at the current point in the Dal...
Double_t retrieveEfficiency()
Obtain the efficiency of the current event from the model.
LauParArray fitFracEffUnCorr_
The efficiency-uncorrected fit fractions for the amplitude components.
Class for defining a K-matrix propagator.