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