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