laura is hosted by Hepforge, IPPP Durham
Laura++  v1r0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauSimpleFitModel.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 #include <iomanip>
17 #include <fstream>
18 using std::cout;
19 using std::cerr;
20 using std::endl;
21 using std::setprecision;
22 
23 #include "TFile.h"
24 #include "TH2.h"
25 #include "TMinuit.h"
26 #include "TRandom.h"
27 #include "TSystem.h"
28 #include "TVirtualFitter.h"
29 
30 #include "LauAbsBkgndDPModel.hh"
31 #include "LauAbsCoeffSet.hh"
32 #include "LauAbsDPDynamics.hh"
33 #include "LauAbsPdf.hh"
34 #include "LauComplex.hh"
35 #include "LauConstants.hh"
36 #include "LauDaughters.hh"
37 #include "LauEmbeddedData.hh"
38 #include "LauEffModel.hh"
39 #include "LauFitNtuple.hh"
40 #include "LauGenNtuple.hh"
41 #include "LauKinematics.hh"
42 #include "LauPrint.hh"
43 #include "LauRandom.hh"
44 #include "LauScfMap.hh"
45 #include "LauSimpleFitModel.hh"
46 
47 ClassImp(LauSimpleFitModel)
48 
49 
51  sigDPModel_(sigModel),
52  kinematics_(sigModel ? sigModel->getKinematics() : 0),
53  usingBkgnd_(kFALSE),
54  nSigComp_(0),
55  nSigDPPar_(0),
56  nExtraPdfPar_(0),
57  nNormPar_(0),
58  meanEff_("meanEff",0.0,0.0,1.0),
59  dpRate_("dpRate",0.0,0.0,100.0),
60  signalEvents_("signalEvents",1.0,0.0,1.0),
61  useSCF_(kFALSE),
62  useSCFHist_(kFALSE),
63  scfFrac_("scfFrac",0.0,0.0,1.0),
64  scfFracHist_(0),
65  scfMap_(0),
66  compareFitData_(kFALSE),
67  signalTree_(0),
68  reuseSignal_(kFALSE),
69  useReweighting_(kFALSE),
70  sigDPLike_(0.0),
71  scfDPLike_(0.0),
72  sigExtraLike_(0.0),
73  scfExtraLike_(0.0),
74  sigTotalLike_(0.0),
75  scfTotalLike_(0.0)
76 {
77 }
78 
80 {
81  delete signalTree_;
82  for (LauBkgndEmbDataList::iterator iter = bkgndTree_.begin(); iter != bkgndTree_.end(); ++iter) {
83  delete (*iter);
84  }
85  delete scfFracHist_;
86 }
87 
89 {
90  UInt_t nBkgnds = this->nBkgndClasses();
91  bkgndDPModels_.resize( nBkgnds );
92  bkgndPdfs_.resize( nBkgnds );
93  bkgndEvents_.resize( nBkgnds );
94  bkgndTree_.resize( nBkgnds );
95  reuseBkgnd_.resize( nBkgnds );
96  bkgndDPLike_.resize( nBkgnds );
97  bkgndExtraLike_.resize( nBkgnds );
98  bkgndTotalLike_.resize( nBkgnds );
99 }
100 
101 void LauSimpleFitModel::setNSigEvents(Double_t nSigEvents, Bool_t fixSigEvents)
102 {
103  signalEvents_.range(-2.0*(TMath::Abs(nSigEvents)+1.0), 2.0*(TMath::Abs(nSigEvents)+1.0));
104  signalEvents_.value(nSigEvents);
105  signalEvents_.initValue(nSigEvents);
106  signalEvents_.genValue(nSigEvents);
107  signalEvents_.fixed(fixSigEvents);
108 }
109 
110 void LauSimpleFitModel::setNBkgndEvents(const TString& bkgndClass, Double_t nBkgndEvents, Bool_t fixBkgndEvents)
111 {
112  if ( ! this->validBkgndClass( bkgndClass ) ) {
113  cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : Invalid background class \"" << bkgndClass << "\"." << std::endl;
114  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
115  return;
116  }
117 
118  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
119 
120  bkgndEvents_[bkgndID].name( bkgndClass+"Events" );
121  bkgndEvents_[bkgndID].range(-2.0*(TMath::Abs(nBkgndEvents)+1.0), 2.0*(TMath::Abs(nBkgndEvents)+1.0));
122  bkgndEvents_[bkgndID].value(nBkgndEvents);
123  bkgndEvents_[bkgndID].initValue(nBkgndEvents);
124  bkgndEvents_[bkgndID].genValue(nBkgndEvents);
125  bkgndEvents_[bkgndID].fixed(fixBkgndEvents);
126 }
127 
128 void LauSimpleFitModel::splitSignalComponent( const TH2* dpHisto, Bool_t upperHalf, LauScfMap* scfMap )
129 {
130  if ( useSCF_ == kTRUE ) {
131  cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : Have already setup SCF." << endl;
132  return;
133  }
134 
135  if ( dpHisto == 0 ) {
136  cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : The histogram pointer is null." << endl;
137  return;
138  }
139 
140  LauDaughters* daughters = sigDPModel_->getDaughters();
141  scfFracHist_ = new LauEffModel( daughters, 0 );
142  scfFracHist_->setEffHisto( dpHisto, kTRUE, kFALSE, 0.0, 0.0, upperHalf, daughters->squareDP() );
143 
144  scfMap_ = scfMap;
145 
146  useSCF_ = kTRUE;
147  useSCFHist_ = kTRUE;
148 }
149 
150 void LauSimpleFitModel::splitSignalComponent( Double_t scfFrac, Bool_t fixed )
151 {
152  if ( useSCF_ == kTRUE ) {
153  cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : Have already setup SCF." << endl;
154  return;
155  }
156 
157  scfFrac_.range( 0.0, 1.0 );
158  scfFrac_.value( scfFrac ); scfFrac_.initValue( scfFrac ); scfFrac_.genValue( scfFrac );
159  scfFrac_.fixed( fixed );
160 
161  useSCF_ = kTRUE;
162  useSCFHist_ = kFALSE;
163 }
164 
165 void LauSimpleFitModel::setBkgndDPModel(const TString& bkgndClass, LauAbsBkgndDPModel* bkgndDPModel)
166 {
167  if (bkgndDPModel == 0) {
168  cerr << "ERROR in LauSimpleFitModel::setBkgndDPModel : The model pointer is null." << endl;
169  return;
170  }
171 
172  // check that this background name is valid
173  if ( ! this->validBkgndClass( bkgndClass ) ) {
174  cerr << "ERROR in LauSimpleFitModel::setBkgndDPModel : Invalid background class \"" << bkgndClass << "\"." << std::endl;
175  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
176  return;
177  }
178 
179  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
180  bkgndDPModels_[bkgndID] = bkgndDPModel;
181 
182  usingBkgnd_ = kTRUE;
183 }
184 
186 {
187  if (pdf==0) {
188  cerr << "ERROR in LauSimpleFitModel::setSignalPdf : The PDF pointer is null." << endl;
189  return;
190  }
191  signalPdfs_.push_back(pdf);
192 }
193 
195 {
196  if (pdf==0) {
197  cerr << "ERROR in LauSimpleFitModel::setSCFPdf : The PDF pointer is null." << endl;
198  return;
199  }
200  scfPdfs_.push_back(pdf);
201 }
202 
203 void LauSimpleFitModel::setBkgndPdf(const TString& bkgndClass, LauAbsPdf* pdf)
204 {
205  if (pdf == 0) {
206  cerr << "ERROR in LauSimpleFitModel::setBkgndPdf : The PDF pointer is null." << endl;
207  return;
208  }
209 
210  // check that this background name is valid
211  if ( ! this->validBkgndClass( bkgndClass ) ) {
212  cerr << "ERROR in LauSimpleFitModel::setBkgndPdf : Invalid background class \"" << bkgndClass << "\"." << std::endl;
213  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
214  return;
215  }
216 
217  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
218  bkgndPdfs_[bkgndID].push_back(pdf);
219 
220  usingBkgnd_ = kTRUE;
221 }
222 
224 {
225  // Is there a component called compName in the signal model?
226  TString compName(coeffSet->name());
227  Bool_t ok = sigDPModel_->hasResonance(compName);
228  if (!ok) {
229  cerr << "ERROR in LauSimpleFitModel::setAmpCoeffSet : Signal DP model doesn't contain component \"" << compName << "\"." << endl;
230  return;
231  }
232 
233  // Do we already have it in our list of names?
234  for (std::vector<LauAbsCoeffSet*>::const_iterator iter = coeffPars_.begin(); iter != coeffPars_.end(); ++iter) {
235  if ((*iter)->name() == compName) {
236  cerr << "ERROR in LauSimpleFitModel::setAmpCoeffSet : Have already set coefficients for \"" << compName << "\"." << endl;
237  return;
238  }
239  }
240 
241  coeffSet->index(nSigComp_);
242  coeffPars_.push_back(coeffSet);
243 
244  ++nSigComp_;
245 
246  cout << "Added coefficients for component \"" << compName << "\" to the fit model." << endl;
247 }
248 
249 
251 {
252  // First of all check that, we have all the Dalitz-plot models
253 
254  if (this->useDP()) {
255  if ( sigDPModel_ == 0 ) {
256  cerr << "ERROR in LauSimpleFitModel::initialise : The pointer to the signal DP model is null.\n";
257  cerr << " : Removing the Dalitz Plot from the model." << endl;
258  this->useDP(kFALSE);
259  }
260  if ( usingBkgnd_ ) {
261  for (LauBkgndDPModelList::const_iterator dpmodel_iter = bkgndDPModels_.begin(); dpmodel_iter != bkgndDPModels_.end(); ++dpmodel_iter ) {
262  if ( (*dpmodel_iter) == 0 ) {
263  cerr << "ERROR in LauSimpleFitModel::initialise : The pointer to one of the background DP models is null.\n";
264  cerr << " : Removing the Dalitz Plot from the model." << endl;
265  this->useDP(kFALSE);
266  break;
267  }
268  }
269  }
270  }
271 
272  // Next check that, if a given component is being used we've got the
273  // right number of PDFs for all the variables involved
274  // TODO - should probably check variable names and so on as well
275 
276  UInt_t nsigpdfvars(0);
277  for ( LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter ) {
278  std::vector<TString> varNames = (*pdf_iter)->varNames();
279  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
280  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
281  ++nsigpdfvars;
282  }
283  }
284  }
285  if (useSCF_) {
286  UInt_t nscfpdfvars(0);
287  for ( LauPdfList::const_iterator pdf_iter = scfPdfs_.begin(); pdf_iter != scfPdfs_.end(); ++pdf_iter ) {
288  std::vector<TString> varNames = (*pdf_iter)->varNames();
289  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
290  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
291  ++nscfpdfvars;
292  }
293  }
294  }
295  if (nscfpdfvars != nsigpdfvars) {
296  cerr << "ERROR in LauSimpleFitModel::initialise : There are " << nsigpdfvars << " TM signal PDF variables but " << nscfpdfvars << " SCF signal PDF variables." << endl;
297  gSystem->Exit(EXIT_FAILURE);
298  }
299  }
300  if (usingBkgnd_) {
301  for (LauBkgndPdfsList::const_iterator bgclass_iter = bkgndPdfs_.begin(); bgclass_iter != bkgndPdfs_.end(); ++bgclass_iter) {
302  UInt_t nbkgndpdfvars(0);
303  const LauPdfList& pdfList = (*bgclass_iter);
304  for ( LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter ) {
305  std::vector<TString> varNames = (*pdf_iter)->varNames();
306  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
307  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
308  ++nbkgndpdfvars;
309  }
310  }
311  }
312  if (nbkgndpdfvars != nsigpdfvars) {
313  cerr << "ERROR in LauSimpleFitModel::initialise : There are " << nsigpdfvars << " signal PDF variables but " << nbkgndpdfvars << " bkgnd PDF variables." << endl;
314  gSystem->Exit(EXIT_FAILURE);
315  }
316  }
317  }
318 
319  // Clear the vectors of parameter information so we can start from scratch
320  this->clearFitParVectors();
321 
322  // Set the fit parameters for signal and background models
323  this->setSignalDPParameters();
324 
325  // Set the fit parameters for the various extra PDFs
326  this->setExtraPdfParameters();
327 
328  // Set the initial bg and signal events
329  this->setFitNEvents();
330 
331  // Check that we have the expected number of fit variables
332  const LauParameterPList& fitVars = this->fitPars();
333  if (fitVars.size() != (nSigDPPar_ + nExtraPdfPar_ + nNormPar_)) {
334  cerr << "ERROR in LauSimpleFitModel::initialise : Number of fit parameters not of expected size." << endl;
335  gSystem->Exit(EXIT_FAILURE);
336  }
337 
338  // From the initial parameter values calculate the coefficients
339  // so they can be passed to the signal model
340  this->updateCoeffs();
341 
342  // Initialisation
343  if (this->useDP() == kTRUE) {
344  this->initialiseDPModels();
345  }
346 
347  if (!this->useDP() && signalPdfs_.empty()) {
348  cerr << "ERROR in LauSimpleFitModel::initialise : Signal model doesn't exist for any variable." << endl;
349  gSystem->Exit(EXIT_FAILURE);
350  }
351 
352  this->setExtraNtupleVars();
353 }
354 
356 {
357  cout << "Initialising signal DP model" << endl;
359 
360  if (usingBkgnd_ == kTRUE) {
361  for (LauBkgndDPModelList::iterator iter = bkgndDPModels_.begin(); iter != bkgndDPModels_.end(); ++iter) {
362  (*iter)->initialise();
363  }
364  }
365 }
366 
368 {
369  // Set the fit parameters for the signal model.
370 
371  nSigDPPar_ = 0;
372  if ( ! this->useDP() ) {
373  return;
374  }
375 
376  cout << "INFO in LauSimpleFitModel::setSignalDPParameters : Setting the initial fit parameters for the signal DP model." << endl;
377 
378  // Need to check that the number of components we have and that the dynamics has matches up
379  UInt_t nAmp = sigDPModel_->getnAmp();
380  if (nAmp != nSigComp_) {
381  cerr << "ERROR in LauSimpleFitModel::setSignalDPParameters : Number of signal DP components with magnitude and phase set not right." << endl;
382  gSystem->Exit(EXIT_FAILURE);
383  }
384 
385 
386  // Place signal model parameters in vector of fit variables
387  LauParameterPList& fitVars = this->fitPars();
388  for (UInt_t i = 0; i < nSigComp_; i++) {
389  LauParameterPList pars = coeffPars_[i]->getParameters();
390  for (LauParameterPList::iterator iter = pars.begin(); iter != pars.end(); ++iter) {
391  if ( !(*iter)->clone() ) {
392  fitVars.push_back(*iter);
393  ++nSigDPPar_;
394  }
395  }
396  }
397 }
398 
400 {
401  // Include all the parameters of the various PDFs in the fit
402  // NB all of them are passed to the fit, even though some have been fixed through parameter.fixed(kTRUE)
403  // With the new "cloned parameter" scheme only "original" parameters are passed to the fit.
404  // Their clones are updated automatically when the originals are updated.
405 
406  nExtraPdfPar_ = 0;
407 
409 
410  if (useSCF_ == kTRUE) {
412  }
413 
414  if (usingBkgnd_ == kTRUE) {
415  for (LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
416  nExtraPdfPar_ += this->addFitParameters(*iter);
417  }
418  }
419 }
420 
422 {
423  nNormPar_ = 0;
424 
425  // initialise the total number of events to be the sum of all the hypotheses
426  Double_t nTotEvts = signalEvents_.value();
427  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
428  nTotEvts += iter->value();
429  }
430  this->eventsPerExpt(TMath::FloorNint(nTotEvts));
431 
432  LauParameterPList& fitVars = this->fitPars();
433 
434  // if doing an extended ML fit add the number of signal events into the fit parameters
435  if (this->doEMLFit()) {
436  cout << "INFO in LauSimpleFitModel::setFitNEvents : Initialising number of events for signal and background components..." << endl;
437  // add the signal events to the list of fit parameters
438  fitVars.push_back(&signalEvents_);
439  ++nNormPar_;
440  } else {
441  cout << "INFO in LauSimpleFitModel::setFitNEvents : Initialising number of events for background components (and hence signal)..." << endl;
442  }
443 
444  if (useSCF_ && !useSCFHist_) {
445  fitVars.push_back(&scfFrac_);
446  ++nNormPar_;
447  }
448 
449  if (usingBkgnd_ == kTRUE) {
450  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
451  LauParameter& parameter = (*iter);
452  fitVars.push_back(&parameter);
453  ++nNormPar_;
454  }
455  }
456 }
457 
459 {
460  // Set-up other parameters derived from the fit results, e.g. fit fractions.
461  if (this->useDP() != kTRUE) {
462  return;
463  }
464 
465  // First clear the vectors so we start from scratch
466  this->clearExtraVarVectors();
467 
468  LauParameterList& extraVars = this->extraPars();
469 
470  // Add a fit fraction for each signal component
472  if (fitFrac_.size() != nSigComp_) {
473  cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << fitFrac_.size() << endl;
474  gSystem->Exit(EXIT_FAILURE);
475  }
476  for (UInt_t i(0); i<nSigComp_; ++i) {
477  if (fitFrac_[i].size() != nSigComp_) {
478  cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << fitFrac_[i].size() << endl;
479  gSystem->Exit(EXIT_FAILURE);
480  }
481  }
482 
483  for (UInt_t i = 0; i < nSigComp_; i++) {
484  for (UInt_t j = i; j < nSigComp_; j++) {
485  extraVars.push_back(fitFrac_[i][j]);
486  }
487  }
488 
489  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
490  std::vector<LauParameter> extraParams = sigDPModel_->getExtraParameters();
491  std::vector<LauParameter>::iterator extraIter;
492  for (extraIter = extraParams.begin(); extraIter != extraParams.end(); ++extraIter) {
493  LauParameter extraParameter = (*extraIter);
494  extraVars.push_back(extraParameter);
495  }
496 
497  // Now add in the DP efficiency value
498  Double_t initMeanEff = sigDPModel_->getMeanEff().initValue();
499  meanEff_.value(initMeanEff); meanEff_.initValue(initMeanEff); meanEff_.genValue(initMeanEff);
500  extraVars.push_back(meanEff_);
501 
502  // Also add in the DP rate
503  Double_t initDPRate = sigDPModel_->getDPRate().initValue();
504  dpRate_.value(initDPRate); dpRate_.initValue(initDPRate); dpRate_.genValue(initDPRate);
505  extraVars.push_back(dpRate_);
506 }
507 
508 void LauSimpleFitModel::finaliseFitResults(const TString& tablePrefixName)
509 {
510  // Retrieve parameters from the fit results for calculations and toy generation
511  // and eventually store these in output root ntuples/text files
512 
513  // Now take the fit parameters and update them as necessary
514  // e.g. to make mag > 0.0 and phase in the right range.
515  // This function will also calculate any other values, such as the
516  // fit fractions, using any errors provided by fitParErrors as appropriate.
517  // Also obtain the pull values: (measured - generated)/(average error)
518 
519  if (this->useDP() == kTRUE) {
520  for (UInt_t i = 0; i < nSigComp_; i++) {
521  // Check whether we have mag > 0.0, and phase in the right range
522  coeffPars_[i]->finaliseValues();
523  }
524  }
525 
526  // update the pulls on the events
527  if (this->doEMLFit()) {
529  }
530  if (useSCF_ && !useSCFHist_) {
532  }
533  if (usingBkgnd_ == kTRUE) {
534  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
535  iter->updatePull();
536  }
537  }
538 
539  // Update the pulls on all the extra PDFs' parameters
541  if (useSCF_ == kTRUE) {
543  }
544  if (usingBkgnd_ == kTRUE) {
545  for (LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
546  this->updateFitParameters(*iter);
547  }
548  }
549 
550  // Fill the fit results to the ntuple for current experiment
551 
552  // update the coefficients and then calculate the fit fractions
553  if (this->useDP() == kTRUE) {
554  this->updateCoeffs();
557 
559  if (fitFrac.size() != nSigComp_) {
560  cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFrac.size() << endl;
561  gSystem->Exit(EXIT_FAILURE);
562  }
563  for (UInt_t i(0); i<nSigComp_; ++i) {
564  if (fitFrac[i].size() != nSigComp_) {
565  cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFrac[i].size() << endl;
566  gSystem->Exit(EXIT_FAILURE);
567  }
568  }
569 
570  for (UInt_t i(0); i<nSigComp_; ++i) {
571  for (UInt_t j(i); j<nSigComp_; ++j) {
572  fitFrac_[i][j].value(fitFrac[i][j].value());
573  }
574  }
575 
578 
579  this->clearExtraVarVectors();
580  LauParameterList& extraVars = this->extraPars();
581 
582  // Then store the final fit parameters, and any extra parameters for
583  // the signal model (e.g. fit fractions)
584  for (UInt_t i(0); i<nSigComp_; ++i) {
585  for (UInt_t j(i); j<nSigComp_; ++j) {
586  extraVars.push_back(fitFrac_[i][j]);
587  }
588  }
589 
590  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
591  std::vector<LauParameter> extraParams = sigDPModel_->getExtraParameters();
592  std::vector<LauParameter>::iterator extraIter;
593  for (extraIter = extraParams.begin(); extraIter != extraParams.end(); ++extraIter) {
594  LauParameter extraParameter = (*extraIter);
595  extraVars.push_back(extraParameter);
596  }
597 
598  // Now add in the DP efficiency value
599  extraVars.push_back(meanEff_);
600 
601  // Also add in the DP rate
602  extraVars.push_back(dpRate_);
603 
604  this->printFitFractions(cout);
605  }
606 
607  const LauParameterPList& fitVars = this->fitPars();
608  const LauParameterList& extraVars = this->extraPars();
609  LauFitNtuple* ntuple = this->fitNtuple();
610  ntuple->storeParsAndErrors(fitVars, extraVars);
611 
612  // find out the correlation matrix for the parameters
613  ntuple->storeCorrMatrix(this->iExpt(), this->nll(), this->fitStatus());
614 
615  // Fill the data into ntuple
616  ntuple->updateFitNtuple();
617 
618  // Print out the partial fit fractions, phases and the
619  // averaged efficiency, reweighted by the dynamics (and anything else)
620  if (this->writeLatexTable()) {
621  TString sigOutFileName(tablePrefixName);
622  sigOutFileName += "_"; sigOutFileName += this->iExpt(); sigOutFileName += "Expt.tex";
623  this->writeOutTable(sigOutFileName);
624  }
625 }
626 
627 void LauSimpleFitModel::printFitFractions(std::ostream& output)
628 {
629  // Print out Fit Fractions, total DP rate and mean efficiency
630  for (UInt_t i = 0; i < nSigComp_; i++) {
631  output << "FitFraction for component " << i << " (" << coeffPars_[i]->name() << ") = " << fitFrac_[i][i] << endl;
632  }
633  output << "Overall DP rate (integral of matrix element squared) = " << dpRate_ << endl;
634  output << "Average efficiency weighted by whole DP dynamics = " << meanEff_ << endl;
635 }
636 
637 void LauSimpleFitModel::writeOutTable(const TString& outputFile)
638 {
639  // Write out the results of the fit to a tex-readable table
640  // TODO - need to include the yields in this table
641  std::ofstream fout(outputFile);
642  LauPrint print;
643 
644  cout << "Writing out results of the fit to the tex file " << outputFile << endl;
645 
646  if (this->useDP() == kTRUE) {
647  // print the fit coefficients in one table
648  coeffPars_.front()->printTableHeading(fout);
649  for (UInt_t i = 0; i < nSigComp_; i++) {
650  coeffPars_[i]->printTableRow(fout);
651  }
652  fout << "\\hline" << endl;
653  fout << "\\end{tabular}" << endl << endl;
654 
655  // print the fit fractions in another
656  fout << "\\begin{tabular}{|l|c|}" << endl;
657  fout << "\\hline" << endl;
658  fout << "Component & FitFraction \\\\" << endl;
659  fout << "\\hline" << endl;
660  Double_t fitFracSum(0.0);
661  for (UInt_t i = 0; i < nSigComp_; i++) {
662  TString resName = coeffPars_[i]->name();
663  resName = resName.ReplaceAll("_", "\\_");
664  fout << resName << " & $";
665  Double_t fitFrac = fitFrac_[i][i].value();
666  fitFracSum += fitFrac;
667  print.printFormat(fout, fitFrac);
668  fout << "$ \\\\" << endl;
669  }
670  fout << "\\hline" << endl;
671 
672  // Also print out sum of fit fractions
673  fout << "Fit Fraction Sum & $";
674  print.printFormat(fout, fitFracSum);
675  fout << "$ \\\\" << endl;
676  fout << "\\hline" << endl;
677 
678  fout << "DP rate & $";
679  print.printFormat(fout, dpRate_.value());
680  fout << "$ \\\\" << endl;
681  fout << "\\hline" << endl;
682 
683  fout << "$< \\varepsilon >$ & $";
684  print.printFormat(fout, meanEff_.value());
685  fout << "$ \\\\" << endl;
686  fout << "\\hline" << endl;
687  fout << "\\end{tabular}" << endl << endl;
688  }
689 
690  if (!signalPdfs_.empty()) {
691  fout << "\\begin{tabular}{|l|c|}" << endl;
692  fout << "\\hline" << endl;
693  if (useSCF_ == kTRUE) {
694  fout << "\\Extra TM Signal PDFs' Parameters: & \\\\" << endl;
695  } else {
696  fout << "\\Extra Signal PDFs' Parameters: & \\\\" << endl;
697  }
698  this->printFitParameters(signalPdfs_, fout);
699  if (useSCF_ == kTRUE && !scfPdfs_.empty()) {
700  fout << "\\hline" << endl;
701  fout << "\\Extra SCF Signal PDFs' Parameters: & \\\\" << endl;
702  this->printFitParameters(scfPdfs_, fout);
703  }
704  if (usingBkgnd_ == kTRUE && !bkgndPdfs_.empty()) {
705  fout << "\\hline" << endl;
706  fout << "\\Extra Background PDFs' Parameters: & \\\\" << endl;
707  for (LauBkgndPdfsList::const_iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
708  this->printFitParameters(*iter, fout);
709  }
710  }
711  fout << "\\hline \n\\end{tabular}" << endl << endl;
712  }
713 }
714 
716 {
717  // Update the number of signal events to be total-sum(background events)
718  this->updateSigEvents();
719 
720  // Check whether we want to have randomised initial fit parameters for the signal model
721  if (this->useRandomInitFitPars() == kTRUE) {
722  cout << "Setting random parameters for the signal DP model" << endl;
723  this->randomiseInitFitPars();
724  }
725 }
726 
728 {
729  // Only randomise those parameters that are not fixed!
730  cout << "Randomising the initial fit magnitudes and phases of the resonances..." << endl;
731 
732  for (UInt_t i = 0; i < nSigComp_; i++) {
733  coeffPars_[i]->randomiseInitValues();
734  }
735 }
736 
738 {
739  // Determine the number of events to generate for each hypothesis
740  // If we're smearing then smear each one individually
741 
742  LauGenInfo nEvtsGen;
743 
744  // Signal
745  Double_t evtWeight(1.0);
746  Int_t nEvts = TMath::FloorNint(signalEvents_.genValue());
747  if ( nEvts < 0 ) {
748  evtWeight = -1.0;
749  nEvts = TMath::Abs( nEvts );
750  }
751  if (this->doPoissonSmearing()) {
752  nEvts = LauRandom::randomFun()->Poisson(nEvts);
753  }
754  nEvtsGen["signal"] = std::make_pair( nEvts, evtWeight );
755 
756  // Backgrounds
757  const UInt_t nBkgnds = this->nBkgndClasses();
758  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
759  const TString& bkgndClass = this->bkgndClassName(bkgndID);
760  const LauParameter& evtsPar = bkgndEvents_[bkgndID];
761  evtWeight = 1.0;
762  nEvts = TMath::FloorNint( evtsPar.genValue() );
763  if ( nEvts < 0 ) {
764  evtWeight = -1.0;
765  nEvts = TMath::Abs( nEvts );
766  }
767  if (this->doPoissonSmearing()) {
768  nEvts = LauRandom::randomFun()->Poisson(nEvts);
769  }
770  nEvtsGen[bkgndClass] = std::make_pair( nEvts, evtWeight );
771  }
772 
773  return nEvtsGen;
774 }
775 
777 {
778  // Routine to generate toy Monte Carlo events according to the various models we have defined.
779 
780  // Determine the number of events to generate for each hypothesis
781  LauGenInfo nEvts = this->eventsToGenerate();
782 
783  Bool_t genOK(kTRUE);
784  Int_t evtNum(0);
785 
786  const UInt_t nBkgnds = this->nBkgndClasses();
787  std::vector<TString> bkgndClassNames(nBkgnds);
788  std::vector<TString> bkgndClassNamesGen(nBkgnds);
789  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
790  TString name( this->bkgndClassName(iBkgnd) );
791  bkgndClassNames[iBkgnd] = name;
792  bkgndClassNamesGen[iBkgnd] = "gen"+name;
793  }
794 
795  const Bool_t storeSCFTruthInfo = ( useSCF_ || ( this->enableEmbedding() && signalTree_ != 0 && signalTree_->haveBranch("mcMatch") ) );
796 
797  // Loop over the hypotheses and generate the requested number of events for each one
798  for (LauGenInfo::const_iterator iter = nEvts.begin(); iter != nEvts.end(); ++iter) {
799 
800  const TString& type(iter->first);
801  Int_t nEvtsGen( iter->second.first );
802  Double_t evtWeight( iter->second.second );
803 
804  for (Int_t iEvt(0); iEvt<nEvtsGen; ++iEvt) {
805 
806  this->setGenNtupleDoubleBranchValue( "evtWeight", evtWeight );
807 
808  if (type == "signal") {
809  this->setGenNtupleIntegerBranchValue("genSig",1);
810  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
811  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], 0 );
812  }
813  genOK = this->generateSignalEvent();
814  } else {
815  this->setGenNtupleIntegerBranchValue("genSig",0);
816  if ( storeSCFTruthInfo ) {
817  this->setGenNtupleIntegerBranchValue("genTMSig",0);
818  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
819  }
820  UInt_t bkgndID(0);
821  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
822  Int_t gen(0);
823  if ( bkgndClassNames[iBkgnd] == type ) {
824  gen = 1;
825  bkgndID = iBkgnd;
826  }
827  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], gen );
828  }
829  genOK = this->generateBkgndEvent(bkgndID);
830  }
831 
832  if (!genOK) {
833  // If there was a problem with the generation then break out and return.
834  // The problem model will have adjusted itself so that all should be OK next time.
835  break;
836  }
837 
838  if (this->useDP() == kTRUE) {
839  this->setDPBranchValues();
840  }
841 
842  // Store the event number (within this experiment)
843  this->setGenNtupleIntegerBranchValue("iEvtWithinExpt",evtNum);
844  // and then increment it
845  ++evtNum;
846 
847  this->fillGenNtupleBranches();
848  if (iEvt%500 == 0) {
849  cout << "INFO in LauSimpleFitModel::genExpt : Generated event number " << iEvt << " out of " << nEvtsGen << " " << type << " events." << endl;
850  }
851  }
852 
853  if (!genOK) {
854  break;
855  }
856  }
857 
858  if (this->useDP() && genOK) {
859 
860  sigDPModel_->checkToyMC(kTRUE,kTRUE);
861 
862  // Get the fit fractions if they're to be written into the latex table
863  if (this->writeLatexTable()) {
865  if (fitFrac.size() != nSigComp_) {
866  cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: " << fitFrac.size() << endl;
867  gSystem->Exit(EXIT_FAILURE);
868  }
869  for (UInt_t i(0); i<nSigComp_; ++i) {
870  if (fitFrac[i].size() != nSigComp_) {
871  cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: " << fitFrac[i].size() << endl;
872  gSystem->Exit(EXIT_FAILURE);
873  }
874  }
875 
876  for (UInt_t i(0); i<nSigComp_; ++i) {
877  for (UInt_t j = i; j < nSigComp_; j++) {
878  fitFrac_[i][j].value(fitFrac[i][j].value());
879  }
880  }
883  }
884  }
885 
886  // If we're reusing embedded events or if the generation is being
887  // reset then clear the lists of used events
888  if (reuseSignal_ || !genOK) {
889  if (signalTree_) {
891  }
892  }
893  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
894  LauEmbeddedData* data = bkgndTree_[bkgndID];
895  if (reuseBkgnd_[bkgndID] || !genOK) {
896  if (data) {
897  data->clearUsedList();
898  }
899  }
900  }
901 
902  return genOK;
903 }
904 
906 {
907  // Generate signal event
908  Bool_t genOK(kTRUE);
909  Bool_t genSCF(kFALSE);
910 
911  if (this->useDP()) {
912  if (this->enableEmbedding() && signalTree_) {
913  if (useReweighting_) {
914  // Select a (random) event from the generated data. Then store the
915  // reconstructed DP co-ords, together with other pdf information,
916  // as the embedded data.
918  } else {
919 
921  }
922  if (signalTree_->haveBranch("mcMatch")) {
923  Int_t match = static_cast<Int_t>(signalTree_->getValue("mcMatch"));
924  if (match) {
925  this->setGenNtupleIntegerBranchValue("genTMSig",1);
926  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
927  genSCF = kFALSE;
928  } else {
929  this->setGenNtupleIntegerBranchValue("genTMSig",0);
930  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
931  genSCF = kTRUE;
932  }
933  }
934  } else {
935  genOK = sigDPModel_->generate();
936  if ( genOK && useSCF_ ) {
937  Double_t frac(0.0);
938  if ( useSCFHist_ ) {
940  } else {
941  frac = scfFrac_.genValue();
942  }
943  if ( frac < LauRandom::randomFun()->Rndm() ) {
944  this->setGenNtupleIntegerBranchValue("genTMSig",1);
945  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
946  genSCF = kFALSE;
947  } else {
948  this->setGenNtupleIntegerBranchValue("genTMSig",0);
949  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
950  genSCF = kTRUE;
951 
952  // Optionally smear the DP position
953  // of the SCF event
954  if ( scfMap_ != 0 ) {
955  Double_t xCoord(0.0), yCoord(0.0);
956  if ( kinematics_->squareDP() ) {
957  xCoord = kinematics_->getmPrime();
958  yCoord = kinematics_->getThetaPrime();
959  } else {
960  xCoord = kinematics_->getm13Sq();
961  yCoord = kinematics_->getm23Sq();
962  }
963 
964  // Find the bin number where this event is generated
965  Int_t binNo = scfMap_->binNumber( xCoord, yCoord );
966 
967  // Retrieve the migration histogram
968  TH2* histo = scfMap_->trueHist( binNo );
969 
970  LauEffModel * effModel = sigDPModel_->getEffModel();
971  do {
972  // Get a random point from the histogram
973  histo->GetRandom2( xCoord, yCoord );
974 
975  // Update the kinematics
976  if ( kinematics_->squareDP() ) {
977  kinematics_->updateSqDPKinematics( xCoord, yCoord );
978  } else {
979  kinematics_->updateKinematics( xCoord, yCoord );
980  }
981  } while ( ! effModel->passVeto( kinematics_ ) );
982  }
983  }
984  }
985  }
986  } else {
987  if (this->enableEmbedding() && signalTree_) {
989  if (signalTree_->haveBranch("mcMatch")) {
990  Int_t match = static_cast<Int_t>(signalTree_->getValue("mcMatch"));
991  if (match) {
992  this->setGenNtupleIntegerBranchValue("genTMSig",1);
993  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
994  genSCF = kFALSE;
995  } else {
996  this->setGenNtupleIntegerBranchValue("genTMSig",0);
997  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
998  genSCF = kTRUE;
999  }
1000  }
1001  } else if ( useSCF_ ) {
1002  Double_t frac = scfFrac_.genValue();
1003  if ( frac < LauRandom::randomFun()->Rndm() ) {
1004  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1005  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1006  genSCF = kFALSE;
1007  } else {
1008  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1009  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1010  genSCF = kTRUE;
1011  }
1012  }
1013  }
1014  if (genOK) {
1015  if ( useSCF_ ) {
1016  if ( genSCF ) {
1018  } else {
1020  }
1021  } else {
1023  }
1024  }
1025  // Check for problems with the embedding
1026  if (this->enableEmbedding() && signalTree_ && (signalTree_->nEvents() == signalTree_->nUsedEvents())) {
1027  cerr << "WARNING in LauSimpleFitModel::generateSignalEvent : Source of embedded signal events used up, clearing the list of used events." << endl;
1029  }
1030  return genOK;
1031 }
1032 
1034 {
1035  // Generate background event
1036  Bool_t genOK(kTRUE);
1037 
1038  LauAbsBkgndDPModel* model = bkgndDPModels_[bkgndID];
1039  LauEmbeddedData* embeddedData(0);
1040  if (this->enableEmbedding()) {
1041  embeddedData = bkgndTree_[bkgndID];
1042  }
1043  LauPdfList* extraPdfs = &bkgndPdfs_[bkgndID];
1044  if (this->useDP()) {
1045  if (embeddedData) {
1046  embeddedData->getEmbeddedEvent(kinematics_);
1047  } else {
1048  if (model == 0) {
1049  cerr << "ERROR in LauSimpleFitModel::generateBkgndEvent : Can't find the DP model for background class \"" << this->bkgndClassName(bkgndID) << "\"." << endl;
1050  gSystem->Exit(EXIT_FAILURE);
1051  }
1052  genOK = model->generate();
1053  }
1054  } else {
1055  if (embeddedData) {
1056  embeddedData->getEmbeddedEvent(0);
1057  }
1058  }
1059  if (genOK) {
1060  this->generateExtraPdfValues(extraPdfs, embeddedData);
1061  }
1062  // Check for problems with the embedding
1063  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1064  cerr << "WARNING in LauSimpleFitModel::generateBkgndEvent : Source of embedded " << this->bkgndClassName(bkgndID) << " events used up, clearing the list of used events." << endl;
1065  embeddedData->clearUsedList();
1066  }
1067  return genOK;
1068 }
1069 
1071 {
1072  // Setup the required ntuple branches
1073  this->addGenNtupleDoubleBranch("evtWeight");
1074  this->addGenNtupleIntegerBranch("genSig");
1075  if ( useSCF_ || ( this->enableEmbedding() && signalTree_ != 0 && signalTree_->haveBranch("mcMatch") ) ) {
1076  this->addGenNtupleIntegerBranch("genTMSig");
1077  this->addGenNtupleIntegerBranch("genSCFSig");
1078  }
1079  const UInt_t nBkgnds = this->nBkgndClasses();
1080  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1081  TString name( this->bkgndClassName(iBkgnd) );
1082  name.Prepend("gen");
1083  this->addGenNtupleIntegerBranch(name);
1084  }
1085  if (this->useDP() == kTRUE) {
1086  this->addGenNtupleDoubleBranch("m12");
1087  this->addGenNtupleDoubleBranch("m23");
1088  this->addGenNtupleDoubleBranch("m13");
1089  this->addGenNtupleDoubleBranch("m12Sq");
1090  this->addGenNtupleDoubleBranch("m23Sq");
1091  this->addGenNtupleDoubleBranch("m13Sq");
1092  this->addGenNtupleDoubleBranch("cosHel12");
1093  this->addGenNtupleDoubleBranch("cosHel23");
1094  this->addGenNtupleDoubleBranch("cosHel13");
1095  if (kinematics_->squareDP()) {
1096  this->addGenNtupleDoubleBranch("mPrime");
1097  this->addGenNtupleDoubleBranch("thPrime");
1098  }
1099  }
1100  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1101  std::vector<TString> varNames = (*pdf_iter)->varNames();
1102  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1103  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1104  this->addGenNtupleDoubleBranch( (*var_iter) );
1105  }
1106  }
1107  }
1108 }
1109 
1111 {
1112  // Store all the DP information
1119  this->setGenNtupleDoubleBranchValue("cosHel12", kinematics_->getc12());
1120  this->setGenNtupleDoubleBranchValue("cosHel23", kinematics_->getc23());
1121  this->setGenNtupleDoubleBranchValue("cosHel13", kinematics_->getc13());
1122  if (kinematics_->squareDP()) {
1125  }
1126 }
1127 
1129 {
1130  // Generate from the extra PDFs
1131  if (extraPdfs) {
1132  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1133  LauFitData genValues;
1134  if (embeddedData) {
1135  genValues = embeddedData->getValues( (*pdf_iter)->varNames() );
1136  } else {
1137  genValues = (*pdf_iter)->generate(kinematics_);
1138  }
1139  for ( LauFitData::const_iterator var_iter = genValues.begin(); var_iter != genValues.end(); ++var_iter ) {
1140  TString varName = var_iter->first;
1141  if ( varName != "m13Sq" && varName != "m23Sq" ) {
1142  Double_t value = var_iter->second;
1143  this->setGenNtupleDoubleBranchValue(varName,value);
1144  }
1145  }
1146  }
1147  }
1148 }
1149 
1151 {
1152  // Update the total normalisation for the signal likelihood
1153  if (this->useDP() == kTRUE) {
1154  this->updateCoeffs();
1156  }
1157 
1158  // Update the signal events from the background events if not doing an extended fit
1159  if ( !this->doEMLFit() && !signalEvents_.fixed() ) {
1160  this->updateSigEvents();
1161  }
1162 }
1163 
1165 {
1166  // The background parameters will have been set from Minuit.
1167  // We need to update the signal events using these.
1168  Double_t nTotEvts = this->eventsPerExpt();
1169 
1170  signalEvents_.range(-2.0*nTotEvts,2.0*nTotEvts);
1171  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1172  (*iter).range(-2.0*nTotEvts,2.0*nTotEvts);
1173  }
1174 
1175  if ( signalEvents_.fixed() ) {
1176  return;
1177  }
1178 
1179  // Subtract background events (if any) from signal.
1180  Double_t signalEvents = nTotEvts;
1181  if ( usingBkgnd_ == kTRUE ) {
1182  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1183  signalEvents -= (*iter);
1184  }
1185  }
1186  signalEvents_.value(signalEvents);
1187 }
1188 
1190 {
1191  // Fill the internal data trees of the signal and backgrond models.
1192 
1193  LauFitDataTree* inputFitData = this->fitData();
1194 
1195  // First the Dalitz plot variables (m_ij^2)
1196  if (this->useDP() == kTRUE) {
1197 
1198  // need to append SCF smearing bins before caching DP amplitudes
1199  if ( scfMap_ != 0 ) {
1200  this->appendBinCentres( inputFitData );
1201  }
1202  sigDPModel_->fillDataTree(*inputFitData);
1203 
1204  if (usingBkgnd_ == kTRUE) {
1205  for (LauBkgndDPModelList::iterator iter = bkgndDPModels_.begin(); iter != bkgndDPModels_.end(); ++iter) {
1206  (*iter)->fillDataTree(*inputFitData);
1207  }
1208  }
1209  }
1210 
1211  // ...and then the extra PDFs
1212  this->cacheInfo(signalPdfs_, *inputFitData);
1213  this->cacheInfo(scfPdfs_, *inputFitData);
1214  for (LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
1215  this->cacheInfo((*iter), *inputFitData);
1216  }
1217 
1218  // and finally the SCF fractions and jacobians
1219  if ( useSCF_ && useSCFHist_ ) {
1220  if ( !inputFitData->haveBranch( "m13Sq" ) || !inputFitData->haveBranch( "m23Sq" ) ) {
1221  cerr << "ERROR in LauSimpleFitModel::cacheInputFitVars : Input data does not contain DP branches and so can't cache the SCF fraction." << endl;
1222  gSystem->Exit(EXIT_FAILURE);
1223  }
1224 
1225  UInt_t nEvents = inputFitData->nEvents();
1226  recoSCFFracs_.clear();
1227  recoSCFFracs_.reserve( nEvents );
1228  if ( kinematics_->squareDP() ) {
1229  recoJacobians_.clear();
1230  recoJacobians_.reserve( nEvents );
1231  }
1232  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1233  const LauFitData& dataValues = inputFitData->getData(iEvt);
1234  LauFitData::const_iterator m13_iter = dataValues.find("m13Sq");
1235  LauFitData::const_iterator m23_iter = dataValues.find("m23Sq");
1236  kinematics_->updateKinematics( m13_iter->second, m23_iter->second );
1237  Double_t scfFrac = scfFracHist_->calcEfficiency( kinematics_ );
1238  recoSCFFracs_.push_back( scfFrac );
1239  if ( kinematics_->squareDP() ) {
1241  }
1242  }
1243  }
1244 }
1245 
1247 {
1248  // We'll be caching the DP amplitudes and efficiencies of the centres of the true bins.
1249  // To do so, we attach some fake points at the end of inputData, the number of the entry
1250  // minus the total number of events corresponding to the number of the histogram for that
1251  // given true bin in the LauScfMap object. (What this means is that when Laura is provided with
1252  // the LauScfMap object by the user, it's the latter who has to make sure that it contains the
1253  // right number of histograms and in exactly the right order!)
1254 
1255  // Get the x and y co-ordinates of the bin centres
1256  std::vector<Double_t> binCentresXCoords;
1257  std::vector<Double_t> binCentresYCoords;
1258  scfMap_->listBinCentres(binCentresXCoords, binCentresYCoords);
1259 
1260  // The SCF histograms could be in square Dalitz plot histograms.
1261  // The dynamics takes normal Dalitz plot coords, so we might have to convert them back.
1262  Bool_t sqDP = kinematics_->squareDP();
1263  UInt_t nBins = binCentresXCoords.size();
1264  fakeSCFFracs_.clear();
1265  fakeSCFFracs_.reserve( nBins );
1266  if ( sqDP ) {
1267  fakeJacobians_.clear();
1268  fakeJacobians_.reserve( nBins );
1269  }
1270 
1271  for (UInt_t iBin = 0; iBin < nBins; ++iBin) {
1272 
1273  if ( sqDP ) {
1274  kinematics_->updateSqDPKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1275  binCentresXCoords[iBin] = kinematics_->getm13Sq();
1276  binCentresYCoords[iBin] = kinematics_->getm23Sq();
1278  } else {
1279  kinematics_->updateKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1280  }
1281 
1283  }
1284 
1285  // Set up inputFitVars_ object to hold the fake events
1286  inputData->appendFakePoints(binCentresXCoords,binCentresYCoords);
1287 }
1288 
1290 {
1291  // Get the DP likelihood for signal and backgrounds
1292  this->getEvtDPLikelihood(iEvt);
1293 
1294  // Get the combined extra PDFs likelihood for signal and backgrounds
1295  this->getEvtExtraLikelihoods(iEvt);
1296 
1297  // If appropriate, combine the TM and SCF likelihoods
1298  Double_t sigLike = sigDPLike_ * sigExtraLike_;
1299  if ( useSCF_ ) {
1300  Double_t scfFrac(0.0);
1301  if (useSCFHist_) {
1302  scfFrac = recoSCFFracs_[iEvt];
1303  } else {
1304  scfFrac = scfFrac_.value();
1305  }
1306  sigLike *= (1.0 - scfFrac);
1307  if ( (scfMap_ != 0) && (this->useDP() == kTRUE) ) {
1308  // if we're smearing the SCF DP PDF then the SCF frac
1309  // is already included in the SCF DP likelihood
1310  sigLike += (scfDPLike_ * scfExtraLike_);
1311  } else {
1312  sigLike += (scfFrac * scfDPLike_ * scfExtraLike_);
1313  }
1314  }
1315 
1316  // Construct the total event likelihood
1317  Double_t likelihood = signalEvents_ * sigLike;
1318  const UInt_t nBkgnds = this->nBkgndClasses();
1319  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1320  likelihood += (bkgndEvents_[bkgndID] * bkgndDPLike_[bkgndID] * bkgndExtraLike_[bkgndID]);
1321  }
1322 
1323  return likelihood;
1324 }
1325 
1327 {
1328  Double_t eventSum(0.0);
1329  eventSum += signalEvents_;
1330  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1331  eventSum += (*iter);
1332  }
1333  return eventSum;
1334 }
1335 
1337 {
1338  // Function to return the signal and background likelihoods for the
1339  // Dalitz plot for the given event evtNo.
1340 
1341  if (this->useDP() == kTRUE) {
1344 
1345  if ( useSCF_ == kTRUE ) {
1346  if ( scfMap_ == 0 ) {
1347  // we're not smearing the SCF DP position
1348  // so the likelihood is the same as the TM
1350  } else {
1351  // calculate the smeared SCF DP likelihood
1352  scfDPLike_ = this->getEvtSCFDPLikelihood(iEvt);
1353  }
1354  }
1355 
1356  const UInt_t nBkgnds = this->nBkgndClasses();
1357  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1358  if (usingBkgnd_ == kTRUE) {
1359  bkgndDPLike_[bkgndID] = bkgndDPModels_[bkgndID]->getLikelihood(iEvt);
1360  } else {
1361  bkgndDPLike_[bkgndID] = 0.0;
1362  }
1363  }
1364  } else {
1365  // There's always going to be a term in the likelihood for the
1366  // signal, so we'd better not zero it.
1367  sigDPLike_ = 1.0;
1368  scfDPLike_ = 1.0;
1369 
1370  const UInt_t nBkgnds = this->nBkgndClasses();
1371  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1372  if (usingBkgnd_ == kTRUE) {
1373  bkgndDPLike_[bkgndID] = 1.0;
1374  } else {
1375  bkgndDPLike_[bkgndID] = 0.0;
1376  }
1377  }
1378  }
1379 }
1380 
1382 {
1383  Double_t scfDPLike(0.0);
1384 
1385  Double_t recoJacobian(1.0);
1386  Double_t xCoord(0.0);
1387  Double_t yCoord(0.0);
1388  Bool_t squareDP = kinematics_->squareDP();
1389  if ( squareDP ) {
1390  xCoord = sigDPModel_->getEvtmPrime();
1391  yCoord = sigDPModel_->getEvtthPrime();
1392  recoJacobian = recoJacobians_[iEvt];
1393  } else {
1394  xCoord = sigDPModel_->getEvtm13Sq();
1395  yCoord = sigDPModel_->getEvtm23Sq();
1396  }
1397 
1398  // Find the bin that our reco event falls in
1399  Int_t recoBin = scfMap_->binNumber( xCoord, yCoord );
1400 
1401  // Find out which true Bins contribute to the given reco bin
1402  const std::vector<Int_t>* trueBins = scfMap_->trueBins(recoBin);
1403 
1404  Int_t nDataEvents = this->eventsPerExpt();
1405 
1406  // Loop over the true bins
1407  for (std::vector<Int_t>::const_iterator iter = trueBins->begin(); iter != trueBins->end(); ++iter)
1408  {
1409  Int_t trueBin = (*iter);
1410 
1411  // prob of a true event in the given true bin migrating to the reco bin
1412  Double_t pRecoGivenTrue = scfMap_->prob( recoBin, trueBin );
1413 
1414  // We've cached the DP amplitudes and the efficiency for the
1415  // true bin centres, just after the data points
1416  sigDPModel_->calcLikelihoodInfo( nDataEvents + trueBin );
1417  Double_t pTrue = sigDPModel_->getEvtLikelihood();
1418 
1419  // Get the cached SCF fraction (and jacobian if we're using the square DP)
1420  Double_t scfFraction = fakeSCFFracs_[ trueBin ];
1421  Double_t jacobian(1.0);
1422  if ( squareDP ) {
1423  jacobian = fakeJacobians_[ trueBin ];
1424  }
1425 
1426  scfDPLike += pTrue * jacobian * scfFraction * pRecoGivenTrue;
1427 
1428  }
1429 
1430  // Divide by the reco jacobian
1431  scfDPLike /= recoJacobian;
1432 
1433  return scfDPLike;
1434 }
1435 
1437 {
1438  // Function to return the signal and background likelihoods for the
1439  // extra variables for the given event evtNo.
1440 
1441  sigExtraLike_ = 1.0;
1442 
1443  for (LauPdfList::iterator iter = signalPdfs_.begin(); iter != signalPdfs_.end(); ++iter) {
1444  (*iter)->calcLikelihoodInfo(iEvt);
1445  sigExtraLike_ *= (*iter)->getLikelihood();
1446  }
1447  if (useSCF_) {
1448  scfExtraLike_ = 1.0;
1449  for (LauPdfList::iterator iter = scfPdfs_.begin(); iter != scfPdfs_.end(); ++iter) {
1450  (*iter)->calcLikelihoodInfo(iEvt);
1451  scfExtraLike_ *= (*iter)->getLikelihood();
1452  }
1453  }
1454 
1455  const UInt_t nBkgnds = this->nBkgndClasses();
1456  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1457  if (usingBkgnd_) {
1458  bkgndExtraLike_[bkgndID] = 1.0;
1459  LauPdfList& pdfList = bkgndPdfs_[bkgndID];
1460  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1461  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1462  bkgndExtraLike_[bkgndID] *= (*pdf_iter)->getLikelihood();
1463  }
1464  } else {
1465  bkgndExtraLike_[bkgndID] = 0.0;
1466  }
1467  }
1468 }
1469 
1471 {
1472  coeffs_.clear();
1473  coeffs_.reserve(nSigComp_);
1474  for (UInt_t i = 0; i < nSigComp_; i++) {
1475  coeffs_.push_back(coeffPars_[i]->particleCoeff());
1476  }
1477 }
1478 
1480 {
1481  // add branches for storing the experiment number and the number of
1482  // the event within the current experiment
1483  this->addSPlotNtupleIntegerBranch("iExpt");
1484  this->addSPlotNtupleIntegerBranch("iEvtWithinExpt");
1485 
1486  // Store the efficiency of the event (for inclusive BF calculations).
1487  if (this->storeDPEff()) {
1488  this->addSPlotNtupleDoubleBranch("efficiency");
1489  if ( sigDPModel_->usingScfModel() ) {
1490  this->addSPlotNtupleDoubleBranch("scffraction");
1491  }
1492  }
1493 
1494  // Store the total event likelihood for each species.
1495  if (useSCF_) {
1496  this->addSPlotNtupleDoubleBranch("sigTMTotalLike");
1497  this->addSPlotNtupleDoubleBranch("sigSCFTotalLike");
1498  this->addSPlotNtupleDoubleBranch("sigSCFFrac");
1499  } else {
1500  this->addSPlotNtupleDoubleBranch("sigTotalLike");
1501  }
1502  if (usingBkgnd_) {
1503  const UInt_t nBkgnds = this->nBkgndClasses();
1504  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1505  TString name( this->bkgndClassName(iBkgnd) );
1506  name += "TotalLike";
1507  this->addSPlotNtupleDoubleBranch(name);
1508  }
1509  }
1510 
1511  // Store the DP likelihoods
1512  if (this->useDP()) {
1513  if (useSCF_) {
1514  this->addSPlotNtupleDoubleBranch("sigTMDPLike");
1515  this->addSPlotNtupleDoubleBranch("sigSCFDPLike");
1516  } else {
1517  this->addSPlotNtupleDoubleBranch("sigDPLike");
1518  }
1519  if (usingBkgnd_) {
1520  const UInt_t nBkgnds = this->nBkgndClasses();
1521  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1522  TString name( this->bkgndClassName(iBkgnd) );
1523  name += "DPLike";
1524  this->addSPlotNtupleDoubleBranch(name);
1525  }
1526  }
1527  }
1528 
1529  // Store the likelihoods for each extra PDF
1530  if (useSCF_) {
1531  this->addSPlotNtupleBranches(&signalPdfs_, "sigTM");
1532  this->addSPlotNtupleBranches(&scfPdfs_, "sigSCF");
1533  } else {
1534  this->addSPlotNtupleBranches(&signalPdfs_, "sig");
1535  }
1536  if (usingBkgnd_) {
1537  const UInt_t nBkgnds = this->nBkgndClasses();
1538  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1539  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1540  const LauPdfList* pdfList = &(bkgndPdfs_[iBkgnd]);
1541  this->addSPlotNtupleBranches(pdfList, bkgndClass);
1542  }
1543  }
1544 }
1545 
1546 void LauSimpleFitModel::addSPlotNtupleBranches(const LauPdfList* extraPdfs, const TString& prefix)
1547 {
1548  if (extraPdfs) {
1549  // Loop through each of the PDFs
1550  for (LauPdfList::const_iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1551 
1552  // Count the number of input variables that are not
1553  // DP variables (used in the case where there is DP
1554  // dependence for e.g. MVA)
1555  UInt_t nVars(0);
1556  std::vector<TString> varNames = (*pdf_iter)->varNames();
1557  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1558  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1559  ++nVars;
1560  }
1561  }
1562 
1563  if ( nVars == 1 ) {
1564  // If the PDF only has one variable then
1565  // simply add one branch for that variable
1566  TString varName = (*pdf_iter)->varName();
1567  TString name(prefix);
1568  name += varName;
1569  name += "Like";
1570  this->addSPlotNtupleDoubleBranch(name);
1571  } else if ( nVars == 2 ) {
1572  // If the PDF has two variables then we
1573  // need a branch for them both together and
1574  // branches for each
1575  TString allVars("");
1576  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1577  allVars += (*var_iter);
1578  TString name(prefix);
1579  name += (*var_iter);
1580  name += "Like";
1581  this->addSPlotNtupleDoubleBranch(name);
1582  }
1583  TString name(prefix);
1584  name += allVars;
1585  name += "Like";
1586  this->addSPlotNtupleDoubleBranch(name);
1587  } else {
1588  cerr << "WARNING in LauSimpleFitModel::addSPlotNtupleBranches : Can't yet deal with 3D PDFs." << endl;
1589  }
1590  }
1591  }
1592 }
1593 
1594 Double_t LauSimpleFitModel::setSPlotNtupleBranchValues(LauPdfList* extraPdfs, const TString& prefix, UInt_t iEvt)
1595 {
1596  // Store the PDF value for each variable in the list
1597  Double_t totalLike(1.0);
1598  Double_t extraLike(0.0);
1599  if (extraPdfs) {
1600  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1601 
1602  // calculate the likelihood for this event
1603  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1604  extraLike = (*pdf_iter)->getLikelihood();
1605  totalLike *= extraLike;
1606 
1607  // Count the number of input variables that are not
1608  // DP variables (used in the case where there is DP
1609  // dependence for e.g. MVA)
1610  UInt_t nVars(0);
1611  std::vector<TString> varNames = (*pdf_iter)->varNames();
1612  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1613  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1614  ++nVars;
1615  }
1616  }
1617 
1618  if ( nVars == 1 ) {
1619  // If the PDF only has one variable then
1620  // simply store the value for that variable
1621  TString varName = (*pdf_iter)->varName();
1622  TString name(prefix);
1623  name += varName;
1624  name += "Like";
1625  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
1626  } else if ( nVars == 2 ) {
1627  // If the PDF has two variables then we
1628  // store the value for them both together
1629  // and for each on their own
1630  TString allVars("");
1631  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1632  allVars += (*var_iter);
1633  TString name(prefix);
1634  name += (*var_iter);
1635  name += "Like";
1636  Double_t indivLike = (*pdf_iter)->getLikelihood( (*var_iter) );
1637  this->setSPlotNtupleDoubleBranchValue(name, indivLike);
1638  }
1639  TString name(prefix);
1640  name += allVars;
1641  name += "Like";
1642  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
1643  } else {
1644  cerr << "WARNING in LauSimpleFitModel::setSPlotNtupleBranchValues : Can't yet deal with 3D PDFs." << endl;
1645  }
1646  }
1647  }
1648  return totalLike;
1649 }
1650 
1652 {
1653  LauSPlot::NameSet nameSet;
1654  if (this->useDP()) {
1655  nameSet.insert("DP");
1656  }
1657  // Loop through all the signal PDFs
1658  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1659  // Loop over the variables involved in each PDF
1660  std::vector<TString> varNames = (*pdf_iter)->varNames();
1661  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1662  // If they are not DP coordinates then add them
1663  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1664  nameSet.insert( (*var_iter) );
1665  }
1666  }
1667  }
1668  return nameSet;
1669 }
1670 
1672 {
1673  LauSPlot::NumbMap numbMap;
1674  if (!signalEvents_.fixed() && this->doEMLFit()) {
1675  numbMap["sig"] = signalEvents_.genValue();
1676  }
1677  if ( usingBkgnd_ ) {
1678  const UInt_t nBkgnds = this->nBkgndClasses();
1679  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1680  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1681  const LauParameter& par = bkgndEvents_[iBkgnd];
1682  if (!par.fixed()) {
1683  numbMap[bkgndClass] = par.genValue();
1684  }
1685  }
1686  }
1687  return numbMap;
1688 }
1689 
1691 {
1692  LauSPlot::NumbMap numbMap;
1693  if (signalEvents_.fixed() && this->doEMLFit()) {
1694  numbMap["sig"] = signalEvents_.genValue();
1695  }
1696  if ( usingBkgnd_ ) {
1697  const UInt_t nBkgnds = this->nBkgndClasses();
1698  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1699  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1700  const LauParameter& par = bkgndEvents_[iBkgnd];
1701  if (par.fixed()) {
1702  numbMap[bkgndClass] = par.genValue();
1703  }
1704  }
1705  }
1706  return numbMap;
1707 }
1708 
1710 {
1711  LauSPlot::TwoDMap twodimMap;
1712 
1713  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1714  // Count the number of input variables that are not DP variables
1715  UInt_t nVars(0);
1716  std::vector<TString> varNames = (*pdf_iter)->varNames();
1717  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1718  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1719  ++nVars;
1720  }
1721  }
1722  if ( nVars == 2 ) {
1723  if (useSCF_) {
1724  twodimMap.insert( std::make_pair( "sigTM", std::make_pair( varNames[0], varNames[1] ) ) );
1725  } else {
1726  twodimMap.insert( std::make_pair( "sig", std::make_pair( varNames[0], varNames[1] ) ) );
1727  }
1728  }
1729  }
1730 
1731  if ( useSCF_ ) {
1732  for (LauPdfList::const_iterator pdf_iter = scfPdfs_.begin(); pdf_iter != scfPdfs_.end(); ++pdf_iter) {
1733  // Count the number of input variables that are not DP variables
1734  UInt_t nVars(0);
1735  std::vector<TString> varNames = (*pdf_iter)->varNames();
1736  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1737  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1738  ++nVars;
1739  }
1740  }
1741  if ( nVars == 2 ) {
1742  twodimMap.insert( std::make_pair( "sigSCF", std::make_pair( varNames[0], varNames[1] ) ) );
1743  }
1744  }
1745  }
1746 
1747  if (usingBkgnd_) {
1748  const UInt_t nBkgnds = this->nBkgndClasses();
1749  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1750  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1751  const LauPdfList& pdfList = bkgndPdfs_[iBkgnd];
1752  for (LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1753  // Count the number of input variables that are not DP variables
1754  UInt_t nVars(0);
1755  std::vector<TString> varNames = (*pdf_iter)->varNames();
1756  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1757  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1758  ++nVars;
1759  }
1760  }
1761  if ( nVars == 2 ) {
1762  twodimMap.insert( std::make_pair( bkgndClass, std::make_pair( varNames[0], varNames[1] ) ) );
1763  }
1764  }
1765  }
1766  }
1767 
1768  return twodimMap;
1769 }
1770 
1772 {
1773  cout << "Storing per-event likelihood values..." << endl;
1774 
1775  // if we've not been using the DP model then we need to cache all
1776  // the info here so that we can get the efficiency from it
1777 
1778  LauFitDataTree* inputFitData = this->fitData();
1779 
1780  if (!this->useDP() && this->storeDPEff()) {
1782  sigDPModel_->fillDataTree(*inputFitData);
1783  }
1784 
1785  UInt_t evtsPerExpt(this->eventsPerExpt());
1786 
1787  for (UInt_t iEvt = 0; iEvt < evtsPerExpt; ++iEvt) {
1788 
1789  this->setSPlotNtupleIntegerBranchValue("iExpt",this->iExpt());
1790  this->setSPlotNtupleIntegerBranchValue("iEvtWithinExpt",iEvt);
1791 
1792  // the DP information
1793  this->getEvtDPLikelihood(iEvt);
1794  if (this->storeDPEff()) {
1795  if (!this->useDP()) {
1797  }
1798  this->setSPlotNtupleDoubleBranchValue("efficiency",sigDPModel_->getEvtEff());
1799  if ( sigDPModel_->usingScfModel() ) {
1801  }
1802  }
1803  if (this->useDP()) {
1805  if (useSCF_) {
1807  this->setSPlotNtupleDoubleBranchValue("sigTMDPLike",sigDPLike_);
1808  this->setSPlotNtupleDoubleBranchValue("sigSCFDPLike",scfDPLike_);
1809  } else {
1810  this->setSPlotNtupleDoubleBranchValue("sigDPLike",sigDPLike_);
1811  }
1812  if (usingBkgnd_) {
1813  const UInt_t nBkgnds = this->nBkgndClasses();
1814  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1815  TString name = this->bkgndClassName(iBkgnd);
1816  name += "DPLike";
1817  this->setSPlotNtupleDoubleBranchValue(name,bkgndDPLike_[iBkgnd]);
1818  }
1819  }
1820  } else {
1821  sigTotalLike_ = 1.0;
1822  if (useSCF_) {
1823  scfTotalLike_ = 1.0;
1824  }
1825  if (usingBkgnd_) {
1826  const UInt_t nBkgnds = this->nBkgndClasses();
1827  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1828  bkgndTotalLike_[iBkgnd] = 1.0;
1829  }
1830  }
1831  }
1832 
1833  // the signal PDF values
1834  if ( useSCF_ ) {
1835  sigTotalLike_ *= this->setSPlotNtupleBranchValues(&signalPdfs_, "sigTM", iEvt);
1836  scfTotalLike_ *= this->setSPlotNtupleBranchValues(&scfPdfs_, "sigSCF", iEvt);
1837  } else {
1838  sigTotalLike_ *= this->setSPlotNtupleBranchValues(&signalPdfs_, "sig", iEvt);
1839  }
1840 
1841  // the background PDF values
1842  if (usingBkgnd_) {
1843  const UInt_t nBkgnds = this->nBkgndClasses();
1844  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1845  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1846  LauPdfList& pdfs = bkgndPdfs_[iBkgnd];
1847  bkgndTotalLike_[iBkgnd] *= this->setSPlotNtupleBranchValues(&(pdfs), bkgndClass, iEvt);
1848  }
1849  }
1850 
1851  // the total likelihoods
1852  if (useSCF_) {
1853  Double_t scfFrac(0.0);
1854  if ( useSCFHist_ ) {
1855  scfFrac = recoSCFFracs_[iEvt];
1856  } else {
1857  scfFrac = scfFrac_.value();
1858  }
1859  this->setSPlotNtupleDoubleBranchValue("sigSCFFrac",scfFrac);
1860  sigTotalLike_ *= ( 1.0 - scfFrac );
1861  if ( scfMap_ == 0 ) {
1862  scfTotalLike_ *= scfFrac;
1863  }
1864  this->setSPlotNtupleDoubleBranchValue("sigTMTotalLike",sigTotalLike_);
1865  this->setSPlotNtupleDoubleBranchValue("sigSCFTotalLike",scfTotalLike_);
1866  } else {
1867  this->setSPlotNtupleDoubleBranchValue("sigTotalLike",sigTotalLike_);
1868  }
1869  if (usingBkgnd_) {
1870  const UInt_t nBkgnds = this->nBkgndClasses();
1871  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1872  TString name = this->bkgndClassName(iBkgnd);
1873  name += "TotalLike";
1875  }
1876  }
1877  // fill the tree
1878  this->fillSPlotNtupleBranches();
1879  }
1880  cout << "Finished storing per-event likelihood values." << endl;
1881 }
1882 
1883 void LauSimpleFitModel::embedSignal(const TString& fileName, const TString& treeName,
1884  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
1885  Bool_t useReweighting)
1886 {
1887  if (signalTree_) {
1888  cerr << "ERROR in LauSimpleFitModel::embedSignal : Already embedding signal from a file." << endl;
1889  return;
1890  }
1891 
1892  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
1893  cerr << "WARNING in LauSimpleFitModel::embedSignal : Conflicting options provided, will not reuse events at all." << endl;
1894  reuseEventsWithinExperiment = kFALSE;
1895  }
1896 
1897  signalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
1898 
1899  Bool_t dataOK = signalTree_->findBranches();
1900  if (!dataOK) {
1901  delete signalTree_; signalTree_ = 0;
1902  cerr << "ERROR in LauSimpleFitModel::embedSignal : Problem creating data tree for embedding." << endl;
1903  return;
1904  }
1905 
1906  reuseSignal_ = reuseEventsWithinEnsemble;
1907  useReweighting_ = useReweighting;
1908 
1909  if (this->enableEmbedding() == kFALSE) {
1910  this->enableEmbedding(kTRUE);
1911  }
1912 }
1913 
1914 void LauSimpleFitModel::embedBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
1915  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
1916 {
1917  if ( ! this->validBkgndClass( bkgndClass ) ) {
1918  cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
1919  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
1920  return;
1921  }
1922 
1923  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
1924 
1925  if (bkgndTree_[bkgndID]) {
1926  cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Already embedding background from a file." << endl;
1927  return;
1928  }
1929 
1930  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
1931  cerr << "WARNING in LauSimpleFitModel::embedBkgnd : Conflicting options provided, will not reuse events at all." << endl;
1932  reuseEventsWithinExperiment = kFALSE;
1933  }
1934 
1935  bkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
1936 
1937  Bool_t dataOK = bkgndTree_[bkgndID]->findBranches();
1938  if (!dataOK) {
1939  delete bkgndTree_[bkgndID]; bkgndTree_[bkgndID] = 0;
1940  cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Problem creating data tree for embedding." << endl;
1941  return;
1942  }
1943 
1944  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
1945 
1946  if (this->enableEmbedding() == kFALSE) {
1947  this->enableEmbedding(kTRUE);
1948  }
1949 }
1950 
1951 void LauSimpleFitModel::weightEvents( const TString& dataFileName, const TString& dataTreeName )
1952 {
1953  // Routine to provide weights for events that are uniformly distributed
1954  // in the DP (or square DP) so as to reproduce the given DP model
1955 
1956  if ( kinematics_->squareDP() ) {
1957  std::cout << "INFO in LauSimpleFitModel::weightEvents : will create weights assuming events were generated flat in the square DP" << std::endl;
1958  } else {
1959  std::cout << "INFO in LauSimpleFitModel::weightEvents : will create weights assuming events were generated flat in phase space" << std::endl;
1960  }
1961 
1962  // This reads in the given dataFile and creates an input
1963  // fit data tree that stores them for all events and experiments.
1964  Bool_t dataOK = this->cacheFitData(dataFileName,dataTreeName);
1965  if (!dataOK) {
1966  std::cerr << "ERROR in LauSimpleFitModel::weightEvents : Problem caching the data." << std::endl;
1967  return;
1968  }
1969 
1970  LauFitDataTree* inputFitData = this->fitData();
1971 
1972  if ( ! inputFitData->haveBranch( "m13Sq_MC" ) || ! inputFitData->haveBranch( "m23Sq_MC" ) ) {
1973  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Cannot find MC truth DP coordinate branches in supplied data, aborting." << std::endl;
1974  return;
1975  }
1976 
1977  // Create the ntuple to hold the DP weights
1978  TString weightsFileName( dataFileName );
1979  Ssiz_t index = weightsFileName.Last('.');
1980  weightsFileName.Insert( index, "_DPweights" );
1981  LauGenNtuple * weightsTuple = new LauGenNtuple( weightsFileName, dataTreeName );
1982  weightsTuple->addIntegerBranch("iExpt");
1983  weightsTuple->addIntegerBranch("iEvtWithinExpt");
1984  weightsTuple->addDoubleBranch("dpModelWeight");
1985 
1986  UInt_t iExpmt = this->iExpt();
1987  UInt_t nExpmt = this->nExpt();
1988  UInt_t firstExpmt = this->firstExpt();
1989  for (iExpmt = firstExpmt; iExpmt < (firstExpmt+nExpmt); ++iExpmt) {
1990 
1991  inputFitData->readExperimentData(iExpmt);
1992  UInt_t nEvents = inputFitData->nEvents();
1993 
1994  if (nEvents < 1) {
1995  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Zero events in experiment " << iExpmt << ", skipping..." << std::endl;
1996  continue;
1997  }
1998 
1999  weightsTuple->setIntegerBranchValue( "iExpt", iExpmt );
2000 
2001  // Calculate and store the weights for the events in this experiment
2002  for ( UInt_t iEvent(0); iEvent < nEvents; ++iEvent ) {
2003 
2004  weightsTuple->setIntegerBranchValue( "iEvtWithinExpt", iEvent );
2005 
2006  const LauFitData& evtData = inputFitData->getData( iEvent );
2007 
2008  Double_t m13Sq_MC = evtData.find("m13Sq_MC")->second;
2009  Double_t m23Sq_MC = evtData.find("m23Sq_MC")->second;
2010 
2011  Double_t dpModelWeight(0.0);
2012 
2013  if ( kinematics_->withinDPLimits( m13Sq_MC, m23Sq_MC ) ) {
2014 
2015  kinematics_->updateKinematics( m13Sq_MC, m23Sq_MC );
2016  dpModelWeight = sigDPModel_->getEventWeight();
2017 
2018  if ( kinematics_->squareDP() ) {
2019  dpModelWeight *= kinematics_->calcSqDPJacobian();
2020  }
2021 
2022  if (LauAbsDPDynamics::GenOK != sigDPModel_->checkToyMC(kTRUE,kFALSE)) {
2023  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Problem in calculating the weight, aborting." << std::endl;
2024  delete weightsTuple;
2025  return;
2026  }
2027  }
2028 
2029  weightsTuple->setDoubleBranchValue( "dpModelWeight", dpModelWeight );
2030  weightsTuple->fillBranches();
2031  }
2032 
2033  }
2034 
2035  weightsTuple->buildIndex( "iExpt", "iEvtWithinExpt" );
2036  weightsTuple->addFriendTree(dataFileName, dataTreeName);
2037  weightsTuple->writeOutGenResults();
2038 
2039  delete weightsTuple;
2040 }
2041 
void setSignalPdf(LauAbsPdf *pdf)
Set the signal PDF for a given variable.
Double_t range() const
The range allowed for the parameter.
LauAbsDPDynamics * sigDPModel_
The signal Dalitz plot model.
void setBkgndDPModel(const TString &bkgndClass, LauAbsBkgndDPModel *bkgndModel)
Set the background DP models.
Bool_t reuseSignal_
Boolean to reuse signal events.
void setDPBranchValues()
Store all of the DP information.
virtual void addSPlotNtupleIntegerBranch(const TString &name)
Add a branch to the sPlot tree for storing an integer.
Class for defining the abstract interface for signal Dalitz plot dynamics.
The abstract interface for a background Dalitz plot model.
void storeCorrMatrix(UInt_t iExpt, Double_t NLL, Int_t fitStatus)
Store the correlation matrix and other fit information.
Definition: LauFitNtuple.cc:61
Bool_t squareDP() const
Are the square Dalitz plot co-ordinates being calculated?
UInt_t firstExpt() const
Obtain the number of the first experiment.
Double_t getc23() const
Get the cosine of the helicity angle theta23.
LauParameter meanEff_
The mean efficiency.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
Double_t sigExtraLike_
Signal likelihood from extra PDFs.
std::vector< LauParameter > LauParameterList
List of parameters.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
UInt_t nEvents() const
Retrieve the number of events.
UInt_t eventsPerExpt() const
Obtain the total number of events in the current experiment.
virtual Double_t getEvtthPrime() const =0
Retrieve the square Dalitz plot coordinate, theta&#39;, for the current event.
Bool_t writeLatexTable() const
Determine whether writing out of the latex table is enabled.
Double_t sigDPLike_
Signal likelihood value.
void randomiseInitFitPars()
Randomise the initial fit parameters.
LauEffModel * scfFracHist_
The histogram giving the DP-dependence of the SCF fraction.
virtual void setNBkgndEvents(const TString &bkgndClass, Double_t nBkgndEvents, Bool_t fixBkgndEvents=kFALSE)
Set the background event yield(s)
void updateSqDPKinematics(Double_t mPrime, Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m&#39; and theta&#39;.
virtual Double_t getTotEvtLikelihood(UInt_t iEvt)
Get the total likelihood for each event.
void setEffHisto(const TH2 *effHisto, Bool_t useInterpolation=kTRUE, Bool_t fluctuateBins=kFALSE, Double_t avEff=-1.0, Double_t absError=-1.0, Bool_t useUpperHalfOnly=kFALSE, Bool_t squareDP=kFALSE)
Set the efficiency variation across the phase space using a predetermined 2D histogram.
Definition: LauEffModel.cc:60
void cacheInfo(LauPdfList &pdfList, const LauFitDataTree &theData)
Have all PDFs in the list cache the data.
virtual void setNSigEvents(Double_t nSigEvents, Bool_t fixSigEvents=kFALSE)
Set the signal event yield.
virtual void setGenNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the gen tree.
std::vector< Double_t > recoJacobians_
The cached values of the sqDP jacobians for each event.
void generateExtraPdfValues(LauPdfList *extraPdfs, LauEmbeddedData *embeddedData)
Generate from the extra PDFs.
Double_t getm23() const
Get the m23 invariant mass.
virtual void initialise()
Initialise the fit.
const TString & name() const
The parameter name.
Double_t getmPrime() const
Get m&#39; value.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:33
UInt_t addFitParameters(LauPdfList &pdfList)
Add parameters of the PDFs in the list to the list of all fit parameters.
std::map< TString, std::pair< Int_t, Double_t > > LauGenInfo
Define a map to be used to store a category name and numbers.
Bool_t useDP() const
Is the Dalitz plot term in the likelihood.
LauDaughters * getDaughters()
Retrieve the daughters.
Double_t getc13() const
Get the cosine of the helicity angle theta13.
File containing declaration of LauAbsCoeffSet class.
void clearUsedList()
Clear the list of used events.
virtual void weightEvents(const TString &dataFileName, const TString &dataTreeName)
Weight events based on the DP model.
TH2 * trueHist(Int_t trueBin)
Retrieve the migration histogram for trueBin.
Definition: LauScfMap.cc:176
virtual void setAmpCoeffSet(LauAbsCoeffSet *coeffSet)
Set the DP amplitude coefficients.
std::multimap< TString, std::pair< TString, TString > > TwoDMap
Type to associate the name of the species that have 2D PDFs with the names of the two variables invol...
Definition: LauSPlot.hh:68
UInt_t nExtraPdfPar_
Number of extra PDF parameters.
Int_t buildIndex(const TString &majorName, const TString &minorName="0")
Create an index table using leaves of the tree.
Double_t scfDPLike_
SCF likelihood value.
void listBinCentres(std::vector< Double_t > &xCoords, std::vector< Double_t > &yCoords) const
Create lists of the co-ordinates of the centres of true bins.
Definition: LauScfMap.cc:114
File containing declaration of LauDaughters class.
const std::vector< Int_t > * trueBins(Int_t recoBin) const
Find which true bins contribute to the given reco bin.
Definition: LauScfMap.cc:155
virtual void setGenNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the gen tree.
virtual void finaliseFitResults(const TString &tablePrefixName)
Get the fit results and store them.
virtual void propagateParUpdates()
Bool_t useSCFHist_
Is the SCF fraction DP-dependent.
void addFriendTree(const TString &rootFileName, const TString &rootTreeName)
Add a friend tree.
File containing declaration of LauScfMap class.
LauBkgndReuseEventsList reuseBkgnd_
Vector of booleans to reuse background events.
std::vector< std::vector< LauParameter > > LauParArray
Type to define an array of parameters.
File containing declaration of LauSimpleFitModel class.
LauParameter getDPRate() const
Retrieve the overall Dalitz plot rate.
Int_t fitStatus() const
Access the fit status information.
File containing declaration of LauPrint class.
virtual void calcLikelihoodInfo(Double_t m13Sq, Double_t m23Sq)=0
Calculate the likelihood (and all associated information) given values of the Dalitz plot coordinates...
void setBkgndPdf(const TString &bkgndClass, LauAbsPdf *pdf)
Set the background PDF.
std::vector< LauAbsPdf * > LauPdfList
List of Pdfs.
virtual void updateCoeffs()
Update the coefficients.
std::map< TString, Double_t > LauFitData
Type for holding event data.
void fillBranches()
Fill branches in the ntuple.
virtual void getEvtDPLikelihood(UInt_t iEvt)
Calculate the signal and background likelihoods for the DP for a given event.
LauEmbeddedData * signalTree_
The signal event tree.
void setSCFPdf(LauAbsPdf *pdf)
Set the SCF PDF for a given variable.
Double_t getm23Sq() const
Get the m23 invariant mass square.
Bool_t useReweighting_
Boolean to use reweighting.
virtual void printFitFractions(std::ostream &output)
Print the fit fractions, total DP rate and mean efficiency.
Bool_t squareDP() const
Determine to use or not the square Dalitz plot.
Definition: LauDaughters.hh:72
virtual Bool_t hasResonance(const TString &resName) const
Check whether this model includes a named resonance.
void addSPlotNtupleBranches(const LauPdfList *extraPdfs, const TString &prefix)
Add sPlot branches for the extra PDFs.
virtual LauSPlot::TwoDMap twodimPDFs() const
Returns the species and variables for all 2D PDFs in the fit.
std::vector< LauParameter > getExtraParameters()
Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions)
Class to define various output print commands.
Definition: LauPrint.hh:29
Double_t getValue(const TString &name) const
Get the value of a specified branch.
std::vector< LauComplex > coeffs_
The complex coefficients.
virtual void fillGenNtupleBranches()
Fill the gen tuple branches.
std::vector< Double_t > bkgndExtraLike_
Background likelihood value(s) from extra PDFs.
void updateFitNtuple()
Update the fit ntuple.
File containing declaration of LauKinematics class.
virtual void setSPlotNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the sPlot tree.
LauParameter getMeanEff() const
Retrieve the mean efficiency across the Dalitz plot.
File containing declaration of LauEmbeddedData class.
LauBkgndYieldList bkgndEvents_
Background yield(s)
void printFitParameters(const LauPdfList &pdfList, std::ostream &fout) const
Print the fit parameters for all PDFs in the list.
LauPdfList signalPdfs_
The signal PDFs.
Bool_t haveBranch(const TString &name) const
Boolean to determine whether branch exists.
void setExtraPdfParameters()
Set the fit parameters for the extra PDFs.
Bool_t generateBkgndEvent(UInt_t bkgndID)
Generate background event.
Double_t scfTotalLike_
Total SCF likelihood.
void readExperimentData(UInt_t iExpt)
Read events only for the given experiment.
void updateSigEvents()
Update the signal events after Minuit sets background parameters.
virtual UInt_t index() const
Retrieve the index number of the coefficient set.
virtual LauSPlot::NumbMap fixdSpeciesNames() const
Returns the names and yields of species that are fixed in the fit.
void updateKinematics(Double_t m13Sq, Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
LauParameter dpRate_
The average DP rate.
const LauParameterList & extraPars() const
Access the extra variables.
std::vector< Double_t > bkgndDPLike_
Background likelihood value(s)
std::vector< LauParameter * > LauParameterPList
List of parameter pointers.
virtual void fillDataTree(const LauFitDataTree &inputFitTree)=0
Obtain data from a fit tree.
UInt_t nNormPar_
Number of parameters.
Bool_t passVeto(const LauKinematics *kinematics) const
Determine whether the given DP position is outside the vetoes.
Definition: LauEffModel.cc:134
Double_t calcEfficiency(const LauKinematics *kinematics) const
Determine the efficiency for a given point in the Dalitz plot.
Definition: LauEffModel.cc:95
void setupGenNtupleBranches()
Setup the required ntuple branches.
void addIntegerBranch(const TString &name)
Add integer branch to tree.
Definition: LauGenNtuple.cc:73
void appendBinCentres(LauFitDataTree *inputData)
virtual void calcExtraInfo(Bool_t init=kFALSE)=0
Calculate the fit fractions, mean efficiency and total DP rate.
virtual void cacheInputFitVars()
Read in the input fit data variables, e.g. m13Sq and m23Sq.
Abstract interface to the fitting and toy MC model.
void setSignalDPParameters()
Set the fit parameters for the DP model.
UInt_t iExpt() const
Obtain the number of the current experiment.
virtual Double_t getEvtScfFraction() const =0
Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) in the D...
virtual void addGenNtupleIntegerBranch(const TString &name)
Add a branch to the gen tree for storing an integer.
const LauFitNtuple * fitNtuple() const
Access the fit ntuple.
std::vector< LauAbsCoeffSet * > coeffPars_
Magnitudes and Phases.
Int_t binNumber(Double_t xCoord, Double_t yCoord) const
Find the global bin number for the given co-ordinates.
Definition: LauScfMap.cc:144
virtual LauSPlot::NumbMap freeSpeciesNames() const
Returns the names and yields of species that are free in the fit.
UInt_t nExpt() const
Obtain the number of experiments.
virtual LauSPlot::NameSet variableNames() const
Returns the names of all variables in the fit.
Bool_t doPoissonSmearing() const
Determine whether Poisson smearing is enabled for the toy MC generation.
virtual void getEvtExtraLikelihoods(UInt_t iEvt)
Determine the signal and background likelihood for the extra variables for a given event...
LauBkgndPdfsList bkgndPdfs_
The background PDFs.
Double_t getc12() const
Get the cosine of the helicity angle theta12.
Double_t getm13() const
Get the m13 invariant mass.
std::map< TString, Double_t > NumbMap
Type to associate a category name with a double precision number, e.g. a yield or PDF value for a giv...
Definition: LauSPlot.hh:62
Bool_t doEMLFit() const
Determine whether an extended maximum likelihood fit it being performed.
Bool_t useSCF_
Is the signal split into TM and SCF.
void splitSignalComponent(const TH2 *dpHisto, Bool_t upperHalf=kFALSE, LauScfMap *scfMap=0)
Split the signal component into well-reconstructed and mis-reconstructed parts.
File containing declaration of LauComplex class.
LauFitData getValues(const std::vector< TString > &names) const
Get values of specified branches.
virtual void setupBkgndVectors()
Define the length of the background vectors.
LauGenInfo eventsToGenerate()
Determine the number of events to generate for each hypothesis.
Class to store the results from the toy MC generation into an ntuple.
Definition: LauGenNtuple.hh:32
Class for defining the abstract interface for complex coefficient classes.
Double_t nll() const
Access the current NLL value.
virtual Double_t getEvtLikelihood() const =0
Retrieve the likelihood for the current event.
std::vector< Double_t > fakeJacobians_
The cached values of the sqDP jacobians for each true bin.
void embedSignal(const TString &fileName, const TString &treeName, Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment=kFALSE, Bool_t useReweighting=kFALSE)
Embed full simulation events for the signal, rather than generating toy from the PDFs.
std::vector< Double_t > bkgndTotalLike_
Total background likelihood(s)
void clearFitParVectors()
Clear the vectors containing fit parameters.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:31
LauBkgndDPModelList bkgndDPModels_
The background Dalitz Plot model.
Class to store the data for embedding in toy experiments.
virtual Double_t getEvtmPrime() const =0
Retrieve the square Dalitz plot coordinate, m&#39;, for the current event.
virtual void addSPlotNtupleDoubleBranch(const TString &name)
Add a branch to the sPlot tree for storing a double.
void clearExtraVarVectors()
Clear the vectors containing extra ntuple variables.
void getEmbeddedEvent(LauKinematics *kinematics)
Retrieve an event from the data sample.
Bool_t withinDPLimits(Double_t m13Sq, Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot...
LauScfMap * scfMap_
The smearing matrix for the SCF DP PDF.
UInt_t nSigComp_
Number of signal components.
void appendFakePoints(const std::vector< Double_t > &xCoords, const std::vector< Double_t > &yCoords)
Add fake events to the data.
virtual ToyMCStatus checkToyMC(Bool_t printErrorMessages=kTRUE, Bool_t printInfoMessages=kFALSE)=0
Check the status of the toy MC generation.
File containing declaration of LauAbsPdf class.
Bool_t storeDPEff() const
Determine whether the efficiency information should be stored in the sPlot ntuple.
virtual Double_t getEvtSCFDPLikelihood(UInt_t iEvt)
Calculate the SCF likelihood for the DP for a given event.
UInt_t bkgndClassID(const TString &className) const
The number assigned to a background class.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
File containing LauRandom namespace.
UInt_t getnAmp() const
Retrieve the number of amplitude components.
Bool_t useRandomInitFitPars() const
Determine whether the initial values of the fit parameters, in particular the isobar coefficient para...
File containing declaration of LauEffModel class.
const LauFitDataTree * fitData() const
Access the data store.
File containing declaration of LauAbsDPDynamics class.
void setIntegerBranchValue(const TString &name, Int_t value)
Set value of an integer branch.
Definition: LauGenNtuple.cc:91
Double_t setSPlotNtupleBranchValues(LauPdfList *extraPdfs, const TString &prefix, UInt_t iEvt)
Set the branches for the sPlot ntuple with extra PDFs.
Bool_t validBkgndClass(const TString &className) const
Check if the given background class is in the list.
virtual Double_t getEvtm23Sq() const =0
Retrieve the invariant mass squared of the second and third daughters in the current event...
virtual void writeOutTable(const TString &outputFile)
Write the fit results in latex table format.
Double_t getm12Sq() const
Get the m12 invariant mass square.
Define a Dalitz plot according to the isobar model.
LauPdfList scfPdfs_
The SCF PDFs.
Class that implements the efficiency description across the signal Dalitz plot.
Definition: LauEffModel.hh:37
Bool_t cacheFitData(const TString &dataFileName, const TString &dataTreeName)
Store variables from the input file into the internal data storage.
UInt_t nUsedEvents() const
Retrieve the number of events that have already been sampled.
Bool_t usingBkgnd_
Background boolean.
virtual void updateCoeffs(const std::vector< LauComplex > &coeffs)
Update the complex coefficients for the resonances.
Double_t getm13Sq() const
Get the m13 invariant mass square.
void storeParsAndErrors(const std::vector< LauParameter * > &fitVars, const std::vector< LauParameter > &extraVars)
Store parameters and their errors.
std::vector< Double_t > recoSCFFracs_
The cached values of the SCF fraction for each event.
virtual void setSPlotNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the sPlot tree.
void embedBkgnd(const TString &bkgndClass, const TString &fileName, const TString &treeName, Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment=kFALSE)
Embed full simulation events for the given background class, rather than generating toy from the PDFs...
Double_t initValue() const
The initial value of the parameter.
Double_t scfExtraLike_
SCF likelihood from extra PDFs.
virtual Double_t getEvtm13Sq() const =0
Retrieve the invariant mass squared of the first and third daughters in the current event...
Bool_t getReweightedEvent(LauAbsDPDynamics *dynamics)
Retrieve an event from the data sample, applying an accept/reject based on the given DP model...
Double_t getm12() const
Get the m12 invariant mass.
UInt_t nSigDPPar_
Number of signal DP parameters.
virtual void fillSPlotNtupleBranches()
Fill the sPlot tuple.
virtual Double_t getEventWeight()=0
Calculate the acceptance rate, for events with the current kinematics, when generating events accordi...
File containing LauConstants namespace.
LauParameter signalEvents_
Signal yield.
std::vector< Double_t > fakeSCFFracs_
The cached values of the SCF fraction for each bin centre.
LauParArray fitFrac_
Fit fractions.
void writeOutGenResults()
Write out the results from the generation.
void printFormat(std::ostream &stream, Double_t value) const
Method to choose the printing format to a specified level of precision.
Definition: LauPrint.cc:32
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:41
Double_t sigTotalLike_
Total signal likelihood.
Bool_t haveBranch(const TString &name) const
Check if the named branch is stored.
Class for representing the 4D smearing matrix for mis-reconstructed signal (self cross feed) ...
Definition: LauScfMap.hh:29
LauParameter scfFrac_
The (global) SCF fraction parameter.
virtual ~LauSimpleFitModel()
Destructor.
Double_t getThetaPrime() const
Get theta&#39; value.
void updatePull()
Call to update the bias and pull values.
UInt_t nBkgndClasses() const
Returns the number of background classes.
const TString & bkgndClassName(UInt_t classID) const
Get the name of a background class from the number.
Bool_t enableEmbedding() const
Determine whether embedding of events is enabled in the generation.
virtual TString name() const
Retrieve the name of the coefficient set.
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:40
const LauParameterPList & fitPars() const
Access the fit variables.
Double_t value() const
The value of the parameter.
const LauParArray & getFitFractions() const
Retrieve the fit fractions for the amplitude components.
UInt_t nEvents() const
Retrieve the number of events.
virtual void addGenNtupleDoubleBranch(const TString &name)
Add a branch to the gen tree for storing a double.
Double_t calcSqDPJacobian()
Calculate the Jacobian for the transformation m23^2, m13^2 -&gt; m&#39;, theta&#39; (square DP) ...
virtual void setupSPlotNtupleBranches()
Add branches to store experiment number and the event number within the experiment.
void updateFitParameters(LauPdfList &pdfList)
Update the fit parameters for all PDFs in the list.
virtual void checkInitFitParams()
Check the initial fit parameters.
virtual Double_t getEvtEff() const =0
Retrieve the efficiency for the current event.
void setExtraNtupleVars()
Set-up other parameters that are derived from the fit results, e.g. fit fractions.
void addDoubleBranch(const TString &name)
Add double branch to tree.
Definition: LauGenNtuple.cc:82
virtual void storePerEvtLlhds()
Store the per event likelihood values.
virtual Bool_t genExpt()
Toy MC generation and fitting overloaded functions.
Bool_t usingScfModel()
Check whether a self cross feed fraction model is being used.
Bool_t findBranches()
Find and read the branches in data tree.
File containing declaration of LauAbsBkgndDPModel class.
void setFitNEvents()
Set the initial yields.
Double_t prob(Int_t recoBin, Int_t trueBin) const
Probability of a true event in the given true bin migrating to the reco bin.
Definition: LauScfMap.cc:165
LauEffModel * getEffModel()
Retrieve the model for the efficiency across the Dalitz plot.
Double_t genValue() const
The value generated for the parameter.
virtual Double_t getEventSum() const
Get the total number of events.
LauKinematics * kinematics_
The Dalitz plot kinematics object.
virtual void initialise(const std::vector< LauComplex > &coeffs)=0
Initialise the Dalitz plot dynamics.
Class to store the input fit variables.
std::set< TString > NameSet
Type to store names, e.g. of the discriminating/control variables.
Definition: LauSPlot.hh:59
File containing declaration of LauFitNtuple class.
virtual void initialiseDPModels()
Initialise the signal DP model.
File containing declaration of LauGenNtuple class.
void setDoubleBranchValue(const TString &name, Double_t value)
Set value of a double branch.
Definition: LauGenNtuple.cc:96
LauBkgndEmbDataList bkgndTree_
The background event tree.
virtual Bool_t generate()=0
Generate a toy MC event from the model.
virtual Bool_t generate()=0
Generate a toy MC signal event.
Bool_t generateSignalEvent()
Generate signal event.