laura is hosted by Hepforge, IPPP Durham
Laura++  v2r1
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 
19 #include "TFile.h"
20 #include "TH2.h"
21 #include "TMinuit.h"
22 #include "TRandom.h"
23 #include "TSystem.h"
24 #include "TVirtualFitter.h"
25 
26 #include "LauAbsBkgndDPModel.hh"
27 #include "LauAbsCoeffSet.hh"
28 #include "LauAbsDPDynamics.hh"
29 #include "LauAbsPdf.hh"
30 #include "LauComplex.hh"
31 #include "LauConstants.hh"
32 #include "LauDaughters.hh"
33 #include "LauEmbeddedData.hh"
34 #include "LauEffModel.hh"
35 #include "LauFitNtuple.hh"
36 #include "LauGenNtuple.hh"
37 #include "LauKinematics.hh"
38 #include "LauPrint.hh"
39 #include "LauRandom.hh"
40 #include "LauScfMap.hh"
41 #include "LauSimpleFitModel.hh"
42 
44 
45 
47  sigDPModel_(sigModel),
48  kinematics_(sigModel ? sigModel->getKinematics() : 0),
49  usingBkgnd_(kFALSE),
50  nSigComp_(0),
51  nSigDPPar_(0),
52  nExtraPdfPar_(0),
53  nNormPar_(0),
54  meanEff_("meanEff",0.0,0.0,1.0),
55  dpRate_("dpRate",0.0,0.0,100.0),
56  //signalEvents_("signalEvents",1.0,0.0,1.0),
57  signalEvents_(0),
58  useSCF_(kFALSE),
59  useSCFHist_(kFALSE),
60  scfFrac_("scfFrac",0.0,0.0,1.0),
61  scfFracHist_(0),
62  scfMap_(0),
63  compareFitData_(kFALSE),
64  signalTree_(0),
65  reuseSignal_(kFALSE),
66  useReweighting_(kFALSE),
67  sigDPLike_(0.0),
68  scfDPLike_(0.0),
69  sigExtraLike_(0.0),
70  scfExtraLike_(0.0),
71  sigTotalLike_(0.0),
72  scfTotalLike_(0.0)
73 {
74 }
75 
77 {
78  delete signalTree_;
79  for (LauBkgndEmbDataList::iterator iter = bkgndTree_.begin(); iter != bkgndTree_.end(); ++iter) {
80  delete (*iter);
81  }
82  delete scfFracHist_;
83 }
84 
86 {
87  UInt_t nBkgnds = this->nBkgndClasses();
88  bkgndDPModels_.resize( nBkgnds );
89  bkgndPdfs_.resize( nBkgnds );
90  bkgndEvents_.resize( nBkgnds );
91  bkgndTree_.resize( nBkgnds );
92  reuseBkgnd_.resize( nBkgnds );
93  bkgndDPLike_.resize( nBkgnds );
94  bkgndExtraLike_.resize( nBkgnds );
95  bkgndTotalLike_.resize( nBkgnds );
96 }
97 
99 {
100  if ( nSigEvents == 0 ) {
101  std::cerr << "ERROR in LauSimpleFitModel::setNSigEvents : The signal yield LauParameter pointer is null." << std::endl;
102  gSystem->Exit(EXIT_FAILURE);
103  }
104  if ( signalEvents_ != 0 ) {
105  std::cerr << "ERROR in LauSimpleFitModel::setNSigEvents : You are trying to overwrite the signal yield." << std::endl;
106  return;
107  }
108 
109  signalEvents_ = nSigEvents;
110  TString name = signalEvents_->name();
111  if ( ! name.Contains("signalEvents") && !( name.BeginsWith("signal") && name.EndsWith("Events") ) ) {
112  signalEvents_->name("signalEvents");
113  }
114  Double_t value = nSigEvents->value();
115  signalEvents_->range(-2.0*(TMath::Abs(value)+1.0), 2.0*(TMath::Abs(value)+1.0));
116 }
117 
119 {
120  if ( nBkgndEvents == 0 ) {
121  std::cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : The background yield LauParameter pointer is null." << std::endl;
122  gSystem->Exit(EXIT_FAILURE);
123  }
124 
125  if ( ! this->validBkgndClass( nBkgndEvents->name() ) ) {
126  std::cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : Invalid background class \"" << nBkgndEvents->name() << "\"." << std::endl;
127  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << std::endl;
128  gSystem->Exit(EXIT_FAILURE);
129  }
130 
131  UInt_t bkgndID = this->bkgndClassID( nBkgndEvents->name() );
132 
133  if ( bkgndEvents_[bkgndID] != 0 ) {
134  std::cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : You are trying to overwrite the background yield." << std::endl;
135  return;
136  }
137 
138  bkgndEvents_[bkgndID] = nBkgndEvents;
139  bkgndEvents_[bkgndID]->name( nBkgndEvents->name()+"Events" );
140  Double_t value = nBkgndEvents->value();
141  bkgndEvents_[bkgndID]->range(-2.0*(TMath::Abs(value)+1.0), 2.0*(TMath::Abs(value)+1.0));
142 }
143 
144 void LauSimpleFitModel::splitSignalComponent( const TH2* dpHisto, Bool_t upperHalf, LauScfMap* scfMap )
145 {
146  if ( useSCF_ == kTRUE ) {
147  std::cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : Have already setup SCF." << std::endl;
148  return;
149  }
150 
151  if ( dpHisto == 0 ) {
152  std::cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : The histogram pointer is null." << std::endl;
153  return;
154  }
155 
156  LauDaughters* daughters = sigDPModel_->getDaughters();
157  scfFracHist_ = new LauEffModel( daughters, 0 );
158  scfFracHist_->setEffHisto( dpHisto, kTRUE, kFALSE, 0.0, 0.0, upperHalf, daughters->squareDP() );
159 
160  scfMap_ = scfMap;
161 
162  useSCF_ = kTRUE;
163  useSCFHist_ = kTRUE;
164 }
165 
166 void LauSimpleFitModel::splitSignalComponent( Double_t scfFrac, Bool_t fixed )
167 {
168  if ( useSCF_ == kTRUE ) {
169  std::cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : Have already setup SCF." << std::endl;
170  return;
171  }
172 
173  scfFrac_.range( 0.0, 1.0 );
174  scfFrac_.value( scfFrac ); scfFrac_.initValue( scfFrac ); scfFrac_.genValue( scfFrac );
175  scfFrac_.fixed( fixed );
176 
177  useSCF_ = kTRUE;
178  useSCFHist_ = kFALSE;
179 }
180 
181 void LauSimpleFitModel::setBkgndDPModel(const TString& bkgndClass, LauAbsBkgndDPModel* bkgndDPModel)
182 {
183  if (bkgndDPModel == 0) {
184  std::cerr << "ERROR in LauSimpleFitModel::setBkgndDPModel : The model pointer is null." << std::endl;
185  return;
186  }
187 
188  // check that this background name is valid
189  if ( ! this->validBkgndClass( bkgndClass ) ) {
190  std::cerr << "ERROR in LauSimpleFitModel::setBkgndDPModel : Invalid background class \"" << bkgndClass << "\"." << std::endl;
191  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << std::endl;
192  return;
193  }
194 
195  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
196  bkgndDPModels_[bkgndID] = bkgndDPModel;
197 
198  usingBkgnd_ = kTRUE;
199 }
200 
202 {
203  if (pdf==0) {
204  std::cerr << "ERROR in LauSimpleFitModel::setSignalPdf : The PDF pointer is null." << std::endl;
205  return;
206  }
207  signalPdfs_.push_back(pdf);
208 }
209 
211 {
212  if (pdf==0) {
213  std::cerr << "ERROR in LauSimpleFitModel::setSCFPdf : The PDF pointer is null." << std::endl;
214  return;
215  }
216  scfPdfs_.push_back(pdf);
217 }
218 
219 void LauSimpleFitModel::setBkgndPdf(const TString& bkgndClass, LauAbsPdf* pdf)
220 {
221  if (pdf == 0) {
222  std::cerr << "ERROR in LauSimpleFitModel::setBkgndPdf : The PDF pointer is null." << std::endl;
223  return;
224  }
225 
226  // check that this background name is valid
227  if ( ! this->validBkgndClass( bkgndClass ) ) {
228  std::cerr << "ERROR in LauSimpleFitModel::setBkgndPdf : Invalid background class \"" << bkgndClass << "\"." << std::endl;
229  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << std::endl;
230  return;
231  }
232 
233  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
234  bkgndPdfs_[bkgndID].push_back(pdf);
235 
236  usingBkgnd_ = kTRUE;
237 }
238 
240 {
241  // Is there a component called compName in the signal model?
242  TString compName(coeffSet->name());
243  Bool_t ok = sigDPModel_->hasResonance(compName);
244  if (!ok) {
245  std::cerr << "ERROR in LauSimpleFitModel::setAmpCoeffSet : Signal DP model doesn't contain component \"" << compName << "\"." << std::endl;
246  return;
247  }
248 
249  // Do we already have it in our list of names?
250  for (std::vector<LauAbsCoeffSet*>::const_iterator iter = coeffPars_.begin(); iter != coeffPars_.end(); ++iter) {
251  if ((*iter)->name() == compName) {
252  std::cerr << "ERROR in LauSimpleFitModel::setAmpCoeffSet : Have already set coefficients for \"" << compName << "\"." << std::endl;
253  return;
254  }
255  }
256 
257  coeffSet->index(nSigComp_);
258  coeffPars_.push_back(coeffSet);
259 
260  ++nSigComp_;
261 
262  std::cout << "INFO in LauSimpleFitModel::setAmpCoeffSet : Added coefficients for component \"" << compName << "\" to the fit model." << std::endl;
263  coeffSet->printParValues();
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  std::cerr << "ERROR in LauSimpleFitModel::initialise : The pointer to the signal DP model is null.\n";
274  std::cerr << " : Removing the Dalitz Plot from the model." << std::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  std::cerr << "ERROR in LauSimpleFitModel::initialise : The pointer to one of the background DP models is null.\n";
281  std::cerr << " : Removing the Dalitz Plot from the model." << std::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  std::cerr << "ERROR in LauSimpleFitModel::initialise : There are " << nsigpdfvars << " TM signal PDF variables but " << nscfpdfvars << " SCF signal PDF variables." << std::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  std::cerr << "ERROR in LauSimpleFitModel::initialise : There are " << nsigpdfvars << " signal PDF variables but " << nbkgndpdfvars << " bkgnd PDF variables." << std::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  std::cerr << "ERROR in LauSimpleFitModel::initialise : Number of fit parameters not of expected size." << std::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  std::cerr << "ERROR in LauSimpleFitModel::initialise : Signal model doesn't exist for any variable." << std::endl;
366  gSystem->Exit(EXIT_FAILURE);
367  }
368 
369  this->setExtraNtupleVars();
370 }
371 
373 {
374  std::cout << "INFO in LauSimpleFitModel::initialiseDPModels : Initialising signal DP model" << std::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  std::cout << "INFO in LauSimpleFitModel::setSignalDPParameters : Setting the initial fit parameters for the signal DP model." << std::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  std::cerr << "ERROR in LauSimpleFitModel::setSignalDPParameters : Number of signal DP components with magnitude and phase set not right." << std::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  std::cerr << "ERROR in LauSimpleFitModel::setFitNEvents : Signal yield not defined." << std::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  std::cerr << "ERROR in LauSimpleFitModel::setFitNEvents : Background yield not defined." << std::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  std::cout << "INFO in LauSimpleFitModel::setFitNEvents : Initialising number of events for signal and background components..." << std::endl;
462  // add the signal events to the list of fit parameters
463  fitVars.push_back(signalEvents_);
464  ++nNormPar_;
465  } else {
466  std::cout << "INFO in LauSimpleFitModel::setFitNEvents : Initialising number of events for background components (and hence signal)..." << std::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  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << fitFrac_.size() << std::endl;
499  gSystem->Exit(EXIT_FAILURE);
500  }
501  for (UInt_t i(0); i<nSigComp_; ++i) {
502  if (fitFrac_[i].size() != nSigComp_) {
503  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << fitFrac_[i].size() << std::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  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFrac.size() << std::endl;
586  gSystem->Exit(EXIT_FAILURE);
587  }
588  for (UInt_t i(0); i<nSigComp_; ++i) {
589  if (fitFrac[i].size() != nSigComp_) {
590  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFrac[i].size() << std::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(std::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(), this->covarianceMatrix());
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] << std::endl;
657  }
658  output << "Overall DP rate (integral of matrix element squared) = " << dpRate_ << std::endl;
659  output << "Average efficiency weighted by whole DP dynamics = " << meanEff_ << std::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  std::cout << "INFO in LauSimpleFitModel::writeOutTable : Writing out results of the fit to the tex file " << outputFile << std::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" << std::endl;
678  fout << "\\end{tabular}" << std::endl << std::endl;
679 
680  // print the fit fractions in another
681  fout << "\\begin{tabular}{|l|c|}" << std::endl;
682  fout << "\\hline" << std::endl;
683  fout << "Component & FitFraction \\\\" << std::endl;
684  fout << "\\hline" << std::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 << "$ \\\\" << std::endl;
694  }
695  fout << "\\hline" << std::endl;
696 
697  // Also print out sum of fit fractions
698  fout << "Fit Fraction Sum & $";
699  print.printFormat(fout, fitFracSum);
700  fout << "$ \\\\" << std::endl;
701  fout << "\\hline" << std::endl;
702 
703  fout << "DP rate & $";
704  print.printFormat(fout, dpRate_.value());
705  fout << "$ \\\\" << std::endl;
706  fout << "\\hline" << std::endl;
707 
708  fout << "$< \\varepsilon >$ & $";
709  print.printFormat(fout, meanEff_.value());
710  fout << "$ \\\\" << std::endl;
711  fout << "\\hline" << std::endl;
712  fout << "\\end{tabular}" << std::endl << std::endl;
713  }
714 
715  if (!signalPdfs_.empty()) {
716  fout << "\\begin{tabular}{|l|c|}" << std::endl;
717  fout << "\\hline" << std::endl;
718  if (useSCF_ == kTRUE) {
719  fout << "\\Extra TM Signal PDFs' Parameters: & \\\\" << std::endl;
720  } else {
721  fout << "\\Extra Signal PDFs' Parameters: & \\\\" << std::endl;
722  }
723  this->printFitParameters(signalPdfs_, fout);
724  if (useSCF_ == kTRUE && !scfPdfs_.empty()) {
725  fout << "\\hline" << std::endl;
726  fout << "\\Extra SCF Signal PDFs' Parameters: & \\\\" << std::endl;
727  this->printFitParameters(scfPdfs_, fout);
728  }
729  if (usingBkgnd_ == kTRUE && !bkgndPdfs_.empty()) {
730  fout << "\\hline" << std::endl;
731  fout << "\\Extra Background PDFs' Parameters: & \\\\" << std::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}" << std::endl << std::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  std::cout << "INFO in LauSimpleFitModel::checkInitFitParams : Setting random parameters for the signal DP model" << std::endl;
748  this->randomiseInitFitPars();
749  }
750 }
751 
753 {
754  // Only randomise those parameters that are not fixed!
755  std::cout << "INFO in LauSimpleFitModel::randomiseInitFitPars : Randomising the initial fit magnitudes and phases of the resonances..." << std::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  // Add efficiency information
833  this->setGenNtupleDoubleBranchValue( "efficiency", 1 );
834 
835  if (type == "signal") {
836  this->setGenNtupleIntegerBranchValue("genSig",1);
837  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
838  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], 0 );
839  }
840  genOK = this->generateSignalEvent();
841  this->setGenNtupleDoubleBranchValue( "efficiency", sigDPModel_->getEvtEff() );
842  } else {
843  this->setGenNtupleIntegerBranchValue("genSig",0);
844  if ( storeSCFTruthInfo ) {
845  this->setGenNtupleIntegerBranchValue("genTMSig",0);
846  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
847  }
848  UInt_t bkgndID(0);
849  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
850  Int_t gen(0);
851  if ( bkgndClassNames[iBkgnd] == type ) {
852  gen = 1;
853  bkgndID = iBkgnd;
854  }
855  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], gen );
856  }
857  genOK = this->generateBkgndEvent(bkgndID);
858  }
859 
860  if (!genOK) {
861  // If there was a problem with the generation then break out and return.
862  // The problem model will have adjusted itself so that all should be OK next time.
863  break;
864  }
865 
866  if (this->useDP() == kTRUE) {
867  this->setDPBranchValues();
868  }
869 
870  // Store the event number (within this experiment)
871  this->setGenNtupleIntegerBranchValue("iEvtWithinExpt",evtNum);
872  // and then increment it
873  ++evtNum;
874 
875  this->fillGenNtupleBranches();
876  if (iEvt%500 == 0) {
877  std::cout << "INFO in LauSimpleFitModel::genExpt : Generated event number " << iEvt << " out of " << nEvtsGen << " " << type << " events." << std::endl;
878  }
879  }
880 
881  if (!genOK) {
882  break;
883  }
884  }
885 
886  if (this->useDP() && genOK) {
887 
888  sigDPModel_->checkToyMC(kTRUE,kTRUE);
889 
890  // Get the fit fractions if they're to be written into the latex table
891  if (this->writeLatexTable()) {
893  if (fitFrac.size() != nSigComp_) {
894  std::cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: " << fitFrac.size() << std::endl;
895  gSystem->Exit(EXIT_FAILURE);
896  }
897  for (UInt_t i(0); i<nSigComp_; ++i) {
898  if (fitFrac[i].size() != nSigComp_) {
899  std::cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: " << fitFrac[i].size() << std::endl;
900  gSystem->Exit(EXIT_FAILURE);
901  }
902  }
903 
904  for (UInt_t i(0); i<nSigComp_; ++i) {
905  for (UInt_t j = i; j < nSigComp_; j++) {
906  fitFrac_[i][j].value(fitFrac[i][j].value());
907  }
908  }
911  }
912  }
913 
914  // If we're reusing embedded events or if the generation is being
915  // reset then clear the lists of used events
916  if (reuseSignal_ || !genOK) {
917  if (signalTree_) {
919  }
920  }
921  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
922  LauEmbeddedData* data = bkgndTree_[bkgndID];
923  if (reuseBkgnd_[bkgndID] || !genOK) {
924  if (data) {
925  data->clearUsedList();
926  }
927  }
928  }
929 
930  return genOK;
931 }
932 
934 {
935  // Generate signal event
936  Bool_t genOK(kTRUE);
937  Bool_t genSCF(kFALSE);
938 
939  if (this->useDP()) {
940  if (this->enableEmbedding() && signalTree_) {
941  if (useReweighting_) {
942  // Select a (random) event from the generated data. Then store the
943  // reconstructed DP co-ords, together with other pdf information,
944  // as the embedded data.
946  } else {
947 
949  }
950  if (signalTree_->haveBranch("mcMatch")) {
951  Int_t match = static_cast<Int_t>(signalTree_->getValue("mcMatch"));
952  if (match) {
953  this->setGenNtupleIntegerBranchValue("genTMSig",1);
954  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
955  genSCF = kFALSE;
956  } else {
957  this->setGenNtupleIntegerBranchValue("genTMSig",0);
958  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
959  genSCF = kTRUE;
960  }
961  }
962  } else {
963  genOK = sigDPModel_->generate();
964  if ( genOK && useSCF_ ) {
965  Double_t frac(0.0);
966  if ( useSCFHist_ ) {
968  } else {
969  frac = scfFrac_.genValue();
970  }
971  if ( frac < LauRandom::randomFun()->Rndm() ) {
972  this->setGenNtupleIntegerBranchValue("genTMSig",1);
973  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
974  genSCF = kFALSE;
975  } else {
976  this->setGenNtupleIntegerBranchValue("genTMSig",0);
977  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
978  genSCF = kTRUE;
979 
980  // Optionally smear the DP position
981  // of the SCF event
982  if ( scfMap_ != 0 ) {
983  Double_t xCoord(0.0), yCoord(0.0);
984  if ( kinematics_->squareDP() ) {
985  xCoord = kinematics_->getmPrime();
986  yCoord = kinematics_->getThetaPrime();
987  } else {
988  xCoord = kinematics_->getm13Sq();
989  yCoord = kinematics_->getm23Sq();
990  }
991 
992  // Find the bin number where this event is generated
993  Int_t binNo = scfMap_->binNumber( xCoord, yCoord );
994 
995  // Retrieve the migration histogram
996  TH2* histo = scfMap_->trueHist( binNo );
997 
998  LauEffModel * effModel = sigDPModel_->getEffModel();
999  do {
1000  // Get a random point from the histogram
1001  histo->GetRandom2( xCoord, yCoord );
1002 
1003  // Update the kinematics
1004  if ( kinematics_->squareDP() ) {
1005  kinematics_->updateSqDPKinematics( xCoord, yCoord );
1006  } else {
1007  kinematics_->updateKinematics( xCoord, yCoord );
1008  }
1009  } while ( ! effModel->passVeto( kinematics_ ) );
1010  }
1011  }
1012  }
1013  }
1014  } else {
1015  if (this->enableEmbedding() && signalTree_) {
1017  if (signalTree_->haveBranch("mcMatch")) {
1018  Int_t match = static_cast<Int_t>(signalTree_->getValue("mcMatch"));
1019  if (match) {
1020  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1021  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1022  genSCF = kFALSE;
1023  } else {
1024  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1025  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1026  genSCF = kTRUE;
1027  }
1028  }
1029  } else if ( useSCF_ ) {
1030  Double_t frac = scfFrac_.genValue();
1031  if ( frac < LauRandom::randomFun()->Rndm() ) {
1032  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1033  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1034  genSCF = kFALSE;
1035  } else {
1036  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1037  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1038  genSCF = kTRUE;
1039  }
1040  }
1041  }
1042  if (genOK) {
1043  if ( useSCF_ ) {
1044  if ( genSCF ) {
1046  } else {
1048  }
1049  } else {
1051  }
1052  }
1053  // Check for problems with the embedding
1054  if (this->enableEmbedding() && signalTree_ && (signalTree_->nEvents() == signalTree_->nUsedEvents())) {
1055  std::cerr << "WARNING in LauSimpleFitModel::generateSignalEvent : Source of embedded signal events used up, clearing the list of used events." << std::endl;
1057  }
1058  return genOK;
1059 }
1060 
1062 {
1063  // Generate background event
1064  Bool_t genOK(kTRUE);
1065 
1066  LauAbsBkgndDPModel* model = bkgndDPModels_[bkgndID];
1067  LauEmbeddedData* embeddedData(0);
1068  if (this->enableEmbedding()) {
1069  embeddedData = bkgndTree_[bkgndID];
1070  }
1071  LauPdfList* extraPdfs = &bkgndPdfs_[bkgndID];
1072  if (this->useDP()) {
1073  if (embeddedData) {
1074  embeddedData->getEmbeddedEvent(kinematics_);
1075  } else {
1076  if (model == 0) {
1077  std::cerr << "ERROR in LauSimpleFitModel::generateBkgndEvent : Can't find the DP model for background class \"" << this->bkgndClassName(bkgndID) << "\"." << std::endl;
1078  gSystem->Exit(EXIT_FAILURE);
1079  }
1080  genOK = model->generate();
1081  }
1082  } else {
1083  if (embeddedData) {
1084  embeddedData->getEmbeddedEvent(0);
1085  }
1086  }
1087  if (genOK) {
1088  this->generateExtraPdfValues(extraPdfs, embeddedData);
1089  }
1090  // Check for problems with the embedding
1091  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1092  std::cerr << "WARNING in LauSimpleFitModel::generateBkgndEvent : Source of embedded " << this->bkgndClassName(bkgndID) << " events used up, clearing the list of used events." << std::endl;
1093  embeddedData->clearUsedList();
1094  }
1095  return genOK;
1096 }
1097 
1099 {
1100  // Setup the required ntuple branches
1101  this->addGenNtupleDoubleBranch("evtWeight");
1102  this->addGenNtupleIntegerBranch("genSig");
1103  this->addGenNtupleDoubleBranch("efficiency");
1104  if ( useSCF_ || ( this->enableEmbedding() && signalTree_ != 0 && signalTree_->haveBranch("mcMatch") ) ) {
1105  this->addGenNtupleIntegerBranch("genTMSig");
1106  this->addGenNtupleIntegerBranch("genSCFSig");
1107  }
1108  const UInt_t nBkgnds = this->nBkgndClasses();
1109  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1110  TString name( this->bkgndClassName(iBkgnd) );
1111  name.Prepend("gen");
1112  this->addGenNtupleIntegerBranch(name);
1113  }
1114  if (this->useDP() == kTRUE) {
1115  this->addGenNtupleDoubleBranch("m12");
1116  this->addGenNtupleDoubleBranch("m23");
1117  this->addGenNtupleDoubleBranch("m13");
1118  this->addGenNtupleDoubleBranch("m12Sq");
1119  this->addGenNtupleDoubleBranch("m23Sq");
1120  this->addGenNtupleDoubleBranch("m13Sq");
1121  this->addGenNtupleDoubleBranch("cosHel12");
1122  this->addGenNtupleDoubleBranch("cosHel23");
1123  this->addGenNtupleDoubleBranch("cosHel13");
1124  if (kinematics_->squareDP()) {
1125  this->addGenNtupleDoubleBranch("mPrime");
1126  this->addGenNtupleDoubleBranch("thPrime");
1127  }
1128  }
1129  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1130  std::vector<TString> varNames = (*pdf_iter)->varNames();
1131  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1132  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1133  this->addGenNtupleDoubleBranch( (*var_iter) );
1134  }
1135  }
1136  }
1137 }
1138 
1140 {
1141  // Store all the DP information
1148  this->setGenNtupleDoubleBranchValue("cosHel12", kinematics_->getc12());
1149  this->setGenNtupleDoubleBranchValue("cosHel23", kinematics_->getc23());
1150  this->setGenNtupleDoubleBranchValue("cosHel13", kinematics_->getc13());
1151  if (kinematics_->squareDP()) {
1154  }
1155 }
1156 
1158 {
1159  // Generate from the extra PDFs
1160  if (extraPdfs) {
1161  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1162  LauFitData genValues;
1163  if (embeddedData) {
1164  genValues = embeddedData->getValues( (*pdf_iter)->varNames() );
1165  } else {
1166  genValues = (*pdf_iter)->generate(kinematics_);
1167  }
1168  for ( LauFitData::const_iterator var_iter = genValues.begin(); var_iter != genValues.end(); ++var_iter ) {
1169  TString varName = var_iter->first;
1170  if ( varName != "m13Sq" && varName != "m23Sq" ) {
1171  Double_t value = var_iter->second;
1172  this->setGenNtupleDoubleBranchValue(varName,value);
1173  }
1174  }
1175  }
1176  }
1177 }
1178 
1180 {
1181  // Update the total normalisation for the signal likelihood
1182  if (this->useDP() == kTRUE) {
1183  this->updateCoeffs();
1185  }
1186 
1187  // Update the signal events from the background events if not doing an extended fit
1188  if ( !this->doEMLFit() && !signalEvents_->fixed() ) {
1189  this->updateSigEvents();
1190  }
1191 }
1192 
1194 {
1195  // The background parameters will have been set from Minuit.
1196  // We need to update the signal events using these.
1197  Double_t nTotEvts = this->eventsPerExpt();
1198 
1199  signalEvents_->range(-2.0*nTotEvts,2.0*nTotEvts);
1200  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1201  (*iter)->range(-2.0*nTotEvts,2.0*nTotEvts);
1202  }
1203 
1204  if ( signalEvents_->fixed() ) {
1205  return;
1206  }
1207 
1208  // Subtract background events (if any) from signal.
1209  Double_t signalEvents = nTotEvts;
1210  if ( usingBkgnd_ == kTRUE ) {
1211  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1212  signalEvents -= (*iter)->value();
1213  }
1214  }
1215  signalEvents_->value(signalEvents);
1216 }
1217 
1219 {
1220  // Fill the internal data trees of the signal and backgrond models.
1221 
1222  LauFitDataTree* inputFitData = this->fitData();
1223 
1224  // First the Dalitz plot variables (m_ij^2)
1225  if (this->useDP() == kTRUE) {
1226 
1227  // need to append SCF smearing bins before caching DP amplitudes
1228  if ( scfMap_ != 0 ) {
1229  this->appendBinCentres( inputFitData );
1230  }
1231  sigDPModel_->fillDataTree(*inputFitData);
1232 
1233  if (usingBkgnd_ == kTRUE) {
1234  for (LauBkgndDPModelList::iterator iter = bkgndDPModels_.begin(); iter != bkgndDPModels_.end(); ++iter) {
1235  (*iter)->fillDataTree(*inputFitData);
1236  }
1237  }
1238  }
1239 
1240  // ...and then the extra PDFs
1241  this->cacheInfo(signalPdfs_, *inputFitData);
1242  this->cacheInfo(scfPdfs_, *inputFitData);
1243  for (LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
1244  this->cacheInfo((*iter), *inputFitData);
1245  }
1246 
1247  // and finally the SCF fractions and jacobians
1248  if ( useSCF_ && useSCFHist_ ) {
1249  if ( !inputFitData->haveBranch( "m13Sq" ) || !inputFitData->haveBranch( "m23Sq" ) ) {
1250  std::cerr << "ERROR in LauSimpleFitModel::cacheInputFitVars : Input data does not contain DP branches and so can't cache the SCF fraction." << std::endl;
1251  gSystem->Exit(EXIT_FAILURE);
1252  }
1253 
1254  UInt_t nEvents = inputFitData->nEvents();
1255  recoSCFFracs_.clear();
1256  recoSCFFracs_.reserve( nEvents );
1257  if ( kinematics_->squareDP() ) {
1258  recoJacobians_.clear();
1259  recoJacobians_.reserve( nEvents );
1260  }
1261  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1262  const LauFitData& dataValues = inputFitData->getData(iEvt);
1263  LauFitData::const_iterator m13_iter = dataValues.find("m13Sq");
1264  LauFitData::const_iterator m23_iter = dataValues.find("m23Sq");
1265  kinematics_->updateKinematics( m13_iter->second, m23_iter->second );
1266  Double_t scfFrac = scfFracHist_->calcEfficiency( kinematics_ );
1267  recoSCFFracs_.push_back( scfFrac );
1268  if ( kinematics_->squareDP() ) {
1270  }
1271  }
1272  }
1273 }
1274 
1276 {
1277  // We'll be caching the DP amplitudes and efficiencies of the centres of the true bins.
1278  // To do so, we attach some fake points at the end of inputData, the number of the entry
1279  // minus the total number of events corresponding to the number of the histogram for that
1280  // given true bin in the LauScfMap object. (What this means is that when Laura is provided with
1281  // the LauScfMap object by the user, it's the latter who has to make sure that it contains the
1282  // right number of histograms and in exactly the right order!)
1283 
1284  // Get the x and y co-ordinates of the bin centres
1285  std::vector<Double_t> binCentresXCoords;
1286  std::vector<Double_t> binCentresYCoords;
1287  scfMap_->listBinCentres(binCentresXCoords, binCentresYCoords);
1288 
1289  // The SCF histograms could be in square Dalitz plot histograms.
1290  // The dynamics takes normal Dalitz plot coords, so we might have to convert them back.
1291  Bool_t sqDP = kinematics_->squareDP();
1292  UInt_t nBins = binCentresXCoords.size();
1293  fakeSCFFracs_.clear();
1294  fakeSCFFracs_.reserve( nBins );
1295  if ( sqDP ) {
1296  fakeJacobians_.clear();
1297  fakeJacobians_.reserve( nBins );
1298  }
1299 
1300  for (UInt_t iBin = 0; iBin < nBins; ++iBin) {
1301 
1302  if ( sqDP ) {
1303  kinematics_->updateSqDPKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1304  binCentresXCoords[iBin] = kinematics_->getm13Sq();
1305  binCentresYCoords[iBin] = kinematics_->getm23Sq();
1307  } else {
1308  kinematics_->updateKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1309  }
1310 
1312  }
1313 
1314  // Set up inputFitVars_ object to hold the fake events
1315  inputData->appendFakePoints(binCentresXCoords,binCentresYCoords);
1316 }
1317 
1319 {
1320  // Get the DP likelihood for signal and backgrounds
1321  this->getEvtDPLikelihood(iEvt);
1322 
1323  // Get the combined extra PDFs likelihood for signal and backgrounds
1324  this->getEvtExtraLikelihoods(iEvt);
1325 
1326  // If appropriate, combine the TM and SCF likelihoods
1327  Double_t sigLike = sigDPLike_ * sigExtraLike_;
1328  if ( useSCF_ ) {
1329  Double_t scfFrac(0.0);
1330  if (useSCFHist_) {
1331  scfFrac = recoSCFFracs_[iEvt];
1332  } else {
1333  scfFrac = scfFrac_.value();
1334  }
1335  sigLike *= (1.0 - scfFrac);
1336  if ( (scfMap_ != 0) && (this->useDP() == kTRUE) ) {
1337  // if we're smearing the SCF DP PDF then the SCF frac
1338  // is already included in the SCF DP likelihood
1339  sigLike += (scfDPLike_ * scfExtraLike_);
1340  } else {
1341  sigLike += (scfFrac * scfDPLike_ * scfExtraLike_);
1342  }
1343  }
1344 
1345  // Construct the total event likelihood
1346  Double_t likelihood = signalEvents_->value() * sigLike;
1347  const UInt_t nBkgnds = this->nBkgndClasses();
1348  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1349  likelihood += (bkgndEvents_[bkgndID]->value() * bkgndDPLike_[bkgndID] * bkgndExtraLike_[bkgndID]);
1350  }
1351 
1352  return likelihood;
1353 }
1354 
1356 {
1357  Double_t eventSum(0.0);
1358  eventSum += signalEvents_->value();
1359  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1360  eventSum += (*iter)->value();
1361  }
1362  return eventSum;
1363 }
1364 
1366 {
1367  // Function to return the signal and background likelihoods for the
1368  // Dalitz plot for the given event evtNo.
1369 
1370  if (this->useDP() == kTRUE) {
1373 
1374  if ( useSCF_ == kTRUE ) {
1375  if ( scfMap_ == 0 ) {
1376  // we're not smearing the SCF DP position
1377  // so the likelihood is the same as the TM
1379  } else {
1380  // calculate the smeared SCF DP likelihood
1381  scfDPLike_ = this->getEvtSCFDPLikelihood(iEvt);
1382  }
1383  }
1384 
1385  const UInt_t nBkgnds = this->nBkgndClasses();
1386  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1387  if (usingBkgnd_ == kTRUE) {
1388  bkgndDPLike_[bkgndID] = bkgndDPModels_[bkgndID]->getLikelihood(iEvt);
1389  } else {
1390  bkgndDPLike_[bkgndID] = 0.0;
1391  }
1392  }
1393  } else {
1394  // There's always going to be a term in the likelihood for the
1395  // signal, so we'd better not zero it.
1396  sigDPLike_ = 1.0;
1397  scfDPLike_ = 1.0;
1398 
1399  const UInt_t nBkgnds = this->nBkgndClasses();
1400  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1401  if (usingBkgnd_ == kTRUE) {
1402  bkgndDPLike_[bkgndID] = 1.0;
1403  } else {
1404  bkgndDPLike_[bkgndID] = 0.0;
1405  }
1406  }
1407  }
1408 }
1409 
1411 {
1412  Double_t scfDPLike(0.0);
1413 
1414  Double_t recoJacobian(1.0);
1415  Double_t xCoord(0.0);
1416  Double_t yCoord(0.0);
1417  Bool_t squareDP = kinematics_->squareDP();
1418  if ( squareDP ) {
1419  xCoord = sigDPModel_->getEvtmPrime();
1420  yCoord = sigDPModel_->getEvtthPrime();
1421  recoJacobian = recoJacobians_[iEvt];
1422  } else {
1423  xCoord = sigDPModel_->getEvtm13Sq();
1424  yCoord = sigDPModel_->getEvtm23Sq();
1425  }
1426 
1427  // Find the bin that our reco event falls in
1428  Int_t recoBin = scfMap_->binNumber( xCoord, yCoord );
1429 
1430  // Find out which true Bins contribute to the given reco bin
1431  const std::vector<Int_t>* trueBins = scfMap_->trueBins(recoBin);
1432 
1433  Int_t nDataEvents = this->eventsPerExpt();
1434 
1435  // Loop over the true bins
1436  for (std::vector<Int_t>::const_iterator iter = trueBins->begin(); iter != trueBins->end(); ++iter)
1437  {
1438  Int_t trueBin = (*iter);
1439 
1440  // prob of a true event in the given true bin migrating to the reco bin
1441  Double_t pRecoGivenTrue = scfMap_->prob( recoBin, trueBin );
1442 
1443  // We've cached the DP amplitudes and the efficiency for the
1444  // true bin centres, just after the data points
1445  sigDPModel_->calcLikelihoodInfo( nDataEvents + trueBin );
1446  Double_t pTrue = sigDPModel_->getEvtLikelihood();
1447 
1448  // Get the cached SCF fraction (and jacobian if we're using the square DP)
1449  Double_t scfFraction = fakeSCFFracs_[ trueBin ];
1450  Double_t jacobian(1.0);
1451  if ( squareDP ) {
1452  jacobian = fakeJacobians_[ trueBin ];
1453  }
1454 
1455  scfDPLike += pTrue * jacobian * scfFraction * pRecoGivenTrue;
1456 
1457  }
1458 
1459  // Divide by the reco jacobian
1460  scfDPLike /= recoJacobian;
1461 
1462  return scfDPLike;
1463 }
1464 
1466 {
1467  // Function to return the signal and background likelihoods for the
1468  // extra variables for the given event evtNo.
1469 
1470  sigExtraLike_ = 1.0;
1471 
1472  for (LauPdfList::iterator iter = signalPdfs_.begin(); iter != signalPdfs_.end(); ++iter) {
1473  (*iter)->calcLikelihoodInfo(iEvt);
1474  sigExtraLike_ *= (*iter)->getLikelihood();
1475  }
1476  if (useSCF_) {
1477  scfExtraLike_ = 1.0;
1478  for (LauPdfList::iterator iter = scfPdfs_.begin(); iter != scfPdfs_.end(); ++iter) {
1479  (*iter)->calcLikelihoodInfo(iEvt);
1480  scfExtraLike_ *= (*iter)->getLikelihood();
1481  }
1482  }
1483 
1484  const UInt_t nBkgnds = this->nBkgndClasses();
1485  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1486  if (usingBkgnd_) {
1487  bkgndExtraLike_[bkgndID] = 1.0;
1488  LauPdfList& pdfList = bkgndPdfs_[bkgndID];
1489  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1490  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1491  bkgndExtraLike_[bkgndID] *= (*pdf_iter)->getLikelihood();
1492  }
1493  } else {
1494  bkgndExtraLike_[bkgndID] = 0.0;
1495  }
1496  }
1497 }
1498 
1500 {
1501  coeffs_.clear();
1502  coeffs_.reserve(nSigComp_);
1503  for (UInt_t i = 0; i < nSigComp_; i++) {
1504  coeffs_.push_back(coeffPars_[i]->particleCoeff());
1505  }
1506 }
1507 
1509 {
1510  // add branches for storing the experiment number and the number of
1511  // the event within the current experiment
1512  this->addSPlotNtupleIntegerBranch("iExpt");
1513  this->addSPlotNtupleIntegerBranch("iEvtWithinExpt");
1514 
1515  // Store the efficiency of the event (for inclusive BF calculations).
1516  if (this->storeDPEff()) {
1517  this->addSPlotNtupleDoubleBranch("efficiency");
1518  if ( sigDPModel_->usingScfModel() ) {
1519  this->addSPlotNtupleDoubleBranch("scffraction");
1520  }
1521  }
1522 
1523  // Store the total event likelihood for each species.
1524  if (useSCF_) {
1525  this->addSPlotNtupleDoubleBranch("sigTMTotalLike");
1526  this->addSPlotNtupleDoubleBranch("sigSCFTotalLike");
1527  this->addSPlotNtupleDoubleBranch("sigSCFFrac");
1528  } else {
1529  this->addSPlotNtupleDoubleBranch("sigTotalLike");
1530  }
1531  if (usingBkgnd_) {
1532  const UInt_t nBkgnds = this->nBkgndClasses();
1533  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1534  TString name( this->bkgndClassName(iBkgnd) );
1535  name += "TotalLike";
1536  this->addSPlotNtupleDoubleBranch(name);
1537  }
1538  }
1539 
1540  // Store the DP likelihoods
1541  if (this->useDP()) {
1542  if (useSCF_) {
1543  this->addSPlotNtupleDoubleBranch("sigTMDPLike");
1544  this->addSPlotNtupleDoubleBranch("sigSCFDPLike");
1545  } else {
1546  this->addSPlotNtupleDoubleBranch("sigDPLike");
1547  }
1548  if (usingBkgnd_) {
1549  const UInt_t nBkgnds = this->nBkgndClasses();
1550  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1551  TString name( this->bkgndClassName(iBkgnd) );
1552  name += "DPLike";
1553  this->addSPlotNtupleDoubleBranch(name);
1554  }
1555  }
1556  }
1557 
1558  // Store the likelihoods for each extra PDF
1559  if (useSCF_) {
1560  this->addSPlotNtupleBranches(&signalPdfs_, "sigTM");
1561  this->addSPlotNtupleBranches(&scfPdfs_, "sigSCF");
1562  } else {
1563  this->addSPlotNtupleBranches(&signalPdfs_, "sig");
1564  }
1565  if (usingBkgnd_) {
1566  const UInt_t nBkgnds = this->nBkgndClasses();
1567  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1568  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1569  const LauPdfList* pdfList = &(bkgndPdfs_[iBkgnd]);
1570  this->addSPlotNtupleBranches(pdfList, bkgndClass);
1571  }
1572  }
1573 }
1574 
1575 void LauSimpleFitModel::addSPlotNtupleBranches(const LauPdfList* extraPdfs, const TString& prefix)
1576 {
1577  if (extraPdfs) {
1578  // Loop through each of the PDFs
1579  for (LauPdfList::const_iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1580 
1581  // Count the number of input variables that are not
1582  // DP variables (used in the case where there is DP
1583  // dependence for e.g. MVA)
1584  UInt_t nVars(0);
1585  std::vector<TString> varNames = (*pdf_iter)->varNames();
1586  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1587  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1588  ++nVars;
1589  }
1590  }
1591 
1592  if ( nVars == 1 ) {
1593  // If the PDF only has one variable then
1594  // simply add one branch for that variable
1595  TString varName = (*pdf_iter)->varName();
1596  TString name(prefix);
1597  name += varName;
1598  name += "Like";
1599  this->addSPlotNtupleDoubleBranch(name);
1600  } else if ( nVars == 2 ) {
1601  // If the PDF has two variables then we
1602  // need a branch for them both together and
1603  // branches for each
1604  TString allVars("");
1605  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1606  allVars += (*var_iter);
1607  TString name(prefix);
1608  name += (*var_iter);
1609  name += "Like";
1610  this->addSPlotNtupleDoubleBranch(name);
1611  }
1612  TString name(prefix);
1613  name += allVars;
1614  name += "Like";
1615  this->addSPlotNtupleDoubleBranch(name);
1616  } else {
1617  std::cerr << "WARNING in LauSimpleFitModel::addSPlotNtupleBranches : Can't yet deal with 3D PDFs." << std::endl;
1618  }
1619  }
1620  }
1621 }
1622 
1623 Double_t LauSimpleFitModel::setSPlotNtupleBranchValues(LauPdfList* extraPdfs, const TString& prefix, UInt_t iEvt)
1624 {
1625  // Store the PDF value for each variable in the list
1626  Double_t totalLike(1.0);
1627  Double_t extraLike(0.0);
1628  if (extraPdfs) {
1629  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1630 
1631  // calculate the likelihood for this event
1632  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1633  extraLike = (*pdf_iter)->getLikelihood();
1634  totalLike *= extraLike;
1635 
1636  // Count the number of input variables that are not
1637  // DP variables (used in the case where there is DP
1638  // dependence for e.g. MVA)
1639  UInt_t nVars(0);
1640  std::vector<TString> varNames = (*pdf_iter)->varNames();
1641  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1642  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1643  ++nVars;
1644  }
1645  }
1646 
1647  if ( nVars == 1 ) {
1648  // If the PDF only has one variable then
1649  // simply store the value for that variable
1650  TString varName = (*pdf_iter)->varName();
1651  TString name(prefix);
1652  name += varName;
1653  name += "Like";
1654  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
1655  } else if ( nVars == 2 ) {
1656  // If the PDF has two variables then we
1657  // store the value for them both together
1658  // and for each on their own
1659  TString allVars("");
1660  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1661  allVars += (*var_iter);
1662  TString name(prefix);
1663  name += (*var_iter);
1664  name += "Like";
1665  Double_t indivLike = (*pdf_iter)->getLikelihood( (*var_iter) );
1666  this->setSPlotNtupleDoubleBranchValue(name, indivLike);
1667  }
1668  TString name(prefix);
1669  name += allVars;
1670  name += "Like";
1671  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
1672  } else {
1673  std::cerr << "WARNING in LauSimpleFitModel::setSPlotNtupleBranchValues : Can't yet deal with 3D PDFs." << std::endl;
1674  }
1675  }
1676  }
1677  return totalLike;
1678 }
1679 
1681 {
1682  LauSPlot::NameSet nameSet;
1683  if (this->useDP()) {
1684  nameSet.insert("DP");
1685  }
1686  // Loop through all the signal PDFs
1687  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1688  // Loop over the variables involved in each PDF
1689  std::vector<TString> varNames = (*pdf_iter)->varNames();
1690  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1691  // If they are not DP coordinates then add them
1692  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1693  nameSet.insert( (*var_iter) );
1694  }
1695  }
1696  }
1697  return nameSet;
1698 }
1699 
1701 {
1702  LauSPlot::NumbMap numbMap;
1703  if (!signalEvents_->fixed() && this->doEMLFit()) {
1704  numbMap["sig"] = signalEvents_->genValue();
1705  }
1706  if ( usingBkgnd_ ) {
1707  const UInt_t nBkgnds = this->nBkgndClasses();
1708  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1709  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1710  const LauParameter* par = bkgndEvents_[iBkgnd];
1711  if (!par->fixed()) {
1712  numbMap[bkgndClass] = par->genValue();
1713  }
1714  }
1715  }
1716  return numbMap;
1717 }
1718 
1720 {
1721  LauSPlot::NumbMap numbMap;
1722  if (signalEvents_->fixed() && this->doEMLFit()) {
1723  numbMap["sig"] = signalEvents_->genValue();
1724  }
1725  if ( usingBkgnd_ ) {
1726  const UInt_t nBkgnds = this->nBkgndClasses();
1727  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1728  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1729  const LauParameter* par = bkgndEvents_[iBkgnd];
1730  if (par->fixed()) {
1731  numbMap[bkgndClass] = par->genValue();
1732  }
1733  }
1734  }
1735  return numbMap;
1736 }
1737 
1739 {
1740  LauSPlot::TwoDMap twodimMap;
1741 
1742  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1743  // Count the number of input variables that are not DP variables
1744  UInt_t nVars(0);
1745  std::vector<TString> varNames = (*pdf_iter)->varNames();
1746  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1747  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1748  ++nVars;
1749  }
1750  }
1751  if ( nVars == 2 ) {
1752  if (useSCF_) {
1753  twodimMap.insert( std::make_pair( "sigTM", std::make_pair( varNames[0], varNames[1] ) ) );
1754  } else {
1755  twodimMap.insert( std::make_pair( "sig", std::make_pair( varNames[0], varNames[1] ) ) );
1756  }
1757  }
1758  }
1759 
1760  if ( useSCF_ ) {
1761  for (LauPdfList::const_iterator pdf_iter = scfPdfs_.begin(); pdf_iter != scfPdfs_.end(); ++pdf_iter) {
1762  // Count the number of input variables that are not DP variables
1763  UInt_t nVars(0);
1764  std::vector<TString> varNames = (*pdf_iter)->varNames();
1765  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1766  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1767  ++nVars;
1768  }
1769  }
1770  if ( nVars == 2 ) {
1771  twodimMap.insert( std::make_pair( "sigSCF", std::make_pair( varNames[0], varNames[1] ) ) );
1772  }
1773  }
1774  }
1775 
1776  if (usingBkgnd_) {
1777  const UInt_t nBkgnds = this->nBkgndClasses();
1778  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1779  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1780  const LauPdfList& pdfList = bkgndPdfs_[iBkgnd];
1781  for (LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1782  // Count the number of input variables that are not DP variables
1783  UInt_t nVars(0);
1784  std::vector<TString> varNames = (*pdf_iter)->varNames();
1785  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1786  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1787  ++nVars;
1788  }
1789  }
1790  if ( nVars == 2 ) {
1791  twodimMap.insert( std::make_pair( bkgndClass, std::make_pair( varNames[0], varNames[1] ) ) );
1792  }
1793  }
1794  }
1795  }
1796 
1797  return twodimMap;
1798 }
1799 
1801 {
1802  std::cout << "INFO in LauSimpleFitModel::storePerEvtLlhds : Storing per-event likelihood values..." << std::endl;
1803 
1804  // if we've not been using the DP model then we need to cache all
1805  // the info here so that we can get the efficiency from it
1806 
1807  LauFitDataTree* inputFitData = this->fitData();
1808 
1809  if (!this->useDP() && this->storeDPEff()) {
1811  sigDPModel_->fillDataTree(*inputFitData);
1812  }
1813 
1814  UInt_t evtsPerExpt(this->eventsPerExpt());
1815 
1816  for (UInt_t iEvt = 0; iEvt < evtsPerExpt; ++iEvt) {
1817 
1818  this->setSPlotNtupleIntegerBranchValue("iExpt",this->iExpt());
1819  this->setSPlotNtupleIntegerBranchValue("iEvtWithinExpt",iEvt);
1820 
1821  // the DP information
1822  this->getEvtDPLikelihood(iEvt);
1823  if (this->storeDPEff()) {
1824  if (!this->useDP()) {
1826  }
1827  this->setSPlotNtupleDoubleBranchValue("efficiency",sigDPModel_->getEvtEff());
1828  if ( sigDPModel_->usingScfModel() ) {
1830  }
1831  }
1832  if (this->useDP()) {
1834  if (useSCF_) {
1836  this->setSPlotNtupleDoubleBranchValue("sigTMDPLike",sigDPLike_);
1837  this->setSPlotNtupleDoubleBranchValue("sigSCFDPLike",scfDPLike_);
1838  } else {
1839  this->setSPlotNtupleDoubleBranchValue("sigDPLike",sigDPLike_);
1840  }
1841  if (usingBkgnd_) {
1842  const UInt_t nBkgnds = this->nBkgndClasses();
1843  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1844  TString name = this->bkgndClassName(iBkgnd);
1845  name += "DPLike";
1846  this->setSPlotNtupleDoubleBranchValue(name,bkgndDPLike_[iBkgnd]);
1847  }
1848  }
1849  } else {
1850  sigTotalLike_ = 1.0;
1851  if (useSCF_) {
1852  scfTotalLike_ = 1.0;
1853  }
1854  if (usingBkgnd_) {
1855  const UInt_t nBkgnds = this->nBkgndClasses();
1856  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1857  bkgndTotalLike_[iBkgnd] = 1.0;
1858  }
1859  }
1860  }
1861 
1862  // the signal PDF values
1863  if ( useSCF_ ) {
1864  sigTotalLike_ *= this->setSPlotNtupleBranchValues(&signalPdfs_, "sigTM", iEvt);
1865  scfTotalLike_ *= this->setSPlotNtupleBranchValues(&scfPdfs_, "sigSCF", iEvt);
1866  } else {
1867  sigTotalLike_ *= this->setSPlotNtupleBranchValues(&signalPdfs_, "sig", iEvt);
1868  }
1869 
1870  // the background PDF values
1871  if (usingBkgnd_) {
1872  const UInt_t nBkgnds = this->nBkgndClasses();
1873  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1874  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1875  LauPdfList& pdfs = bkgndPdfs_[iBkgnd];
1876  bkgndTotalLike_[iBkgnd] *= this->setSPlotNtupleBranchValues(&(pdfs), bkgndClass, iEvt);
1877  }
1878  }
1879 
1880  // the total likelihoods
1881  if (useSCF_) {
1882  Double_t scfFrac(0.0);
1883  if ( useSCFHist_ ) {
1884  scfFrac = recoSCFFracs_[iEvt];
1885  } else {
1886  scfFrac = scfFrac_.value();
1887  }
1888  this->setSPlotNtupleDoubleBranchValue("sigSCFFrac",scfFrac);
1889  sigTotalLike_ *= ( 1.0 - scfFrac );
1890  if ( scfMap_ == 0 ) {
1891  scfTotalLike_ *= scfFrac;
1892  }
1893  this->setSPlotNtupleDoubleBranchValue("sigTMTotalLike",sigTotalLike_);
1894  this->setSPlotNtupleDoubleBranchValue("sigSCFTotalLike",scfTotalLike_);
1895  } else {
1896  this->setSPlotNtupleDoubleBranchValue("sigTotalLike",sigTotalLike_);
1897  }
1898  if (usingBkgnd_) {
1899  const UInt_t nBkgnds = this->nBkgndClasses();
1900  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1901  TString name = this->bkgndClassName(iBkgnd);
1902  name += "TotalLike";
1904  }
1905  }
1906  // fill the tree
1907  this->fillSPlotNtupleBranches();
1908  }
1909  std::cout << "INFO in LauSimpleFitModel::storePerEvtLlhds : Finished storing per-event likelihood values." << std::endl;
1910 }
1911 
1912 void LauSimpleFitModel::embedSignal(const TString& fileName, const TString& treeName,
1913  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
1914  Bool_t useReweighting)
1915 {
1916  if (signalTree_) {
1917  std::cerr << "ERROR in LauSimpleFitModel::embedSignal : Already embedding signal from a file." << std::endl;
1918  return;
1919  }
1920 
1921  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
1922  std::cerr << "WARNING in LauSimpleFitModel::embedSignal : Conflicting options provided, will not reuse events at all." << std::endl;
1923  reuseEventsWithinExperiment = kFALSE;
1924  }
1925 
1926  signalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
1927 
1928  Bool_t dataOK = signalTree_->findBranches();
1929  if (!dataOK) {
1930  delete signalTree_; signalTree_ = 0;
1931  std::cerr << "ERROR in LauSimpleFitModel::embedSignal : Problem creating data tree for embedding." << std::endl;
1932  return;
1933  }
1934 
1935  reuseSignal_ = reuseEventsWithinEnsemble;
1936  useReweighting_ = useReweighting;
1937 
1938  if (this->enableEmbedding() == kFALSE) {
1939  this->enableEmbedding(kTRUE);
1940  }
1941 }
1942 
1943 void LauSimpleFitModel::embedBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
1944  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
1945 {
1946  if ( ! this->validBkgndClass( bkgndClass ) ) {
1947  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
1948  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << std::endl;
1949  return;
1950  }
1951 
1952  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
1953 
1954  if (bkgndTree_[bkgndID]) {
1955  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Already embedding background from a file." << std::endl;
1956  return;
1957  }
1958 
1959  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
1960  std::cerr << "WARNING in LauSimpleFitModel::embedBkgnd : Conflicting options provided, will not reuse events at all." << std::endl;
1961  reuseEventsWithinExperiment = kFALSE;
1962  }
1963 
1964  bkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
1965 
1966  Bool_t dataOK = bkgndTree_[bkgndID]->findBranches();
1967  if (!dataOK) {
1968  delete bkgndTree_[bkgndID]; bkgndTree_[bkgndID] = 0;
1969  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Problem creating data tree for embedding." << std::endl;
1970  return;
1971  }
1972 
1973  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
1974 
1975  if (this->enableEmbedding() == kFALSE) {
1976  this->enableEmbedding(kTRUE);
1977  }
1978 }
1979 
1980 void LauSimpleFitModel::weightEvents( const TString& dataFileName, const TString& dataTreeName )
1981 {
1982  // Routine to provide weights for events that are uniformly distributed
1983  // in the DP (or square DP) so as to reproduce the given DP model
1984 
1985  if ( kinematics_->squareDP() ) {
1986  std::cout << "INFO in LauSimpleFitModel::weightEvents : will create weights assuming events were generated flat in the square DP" << std::endl;
1987  } else {
1988  std::cout << "INFO in LauSimpleFitModel::weightEvents : will create weights assuming events were generated flat in phase space" << std::endl;
1989  }
1990 
1991  // This reads in the given dataFile and creates an input
1992  // fit data tree that stores them for all events and experiments.
1993  Bool_t dataOK = this->cacheFitData(dataFileName,dataTreeName);
1994  if (!dataOK) {
1995  std::cerr << "ERROR in LauSimpleFitModel::weightEvents : Problem caching the data." << std::endl;
1996  return;
1997  }
1998 
1999  LauFitDataTree* inputFitData = this->fitData();
2000 
2001  if ( ! inputFitData->haveBranch( "m13Sq_MC" ) || ! inputFitData->haveBranch( "m23Sq_MC" ) ) {
2002  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Cannot find MC truth DP coordinate branches in supplied data, aborting." << std::endl;
2003  return;
2004  }
2005 
2006  // Create the ntuple to hold the DP weights
2007  TString weightsFileName( dataFileName );
2008  Ssiz_t index = weightsFileName.Last('.');
2009  weightsFileName.Insert( index, "_DPweights" );
2010  LauGenNtuple * weightsTuple = new LauGenNtuple( weightsFileName, dataTreeName );
2011  weightsTuple->addIntegerBranch("iExpt");
2012  weightsTuple->addIntegerBranch("iEvtWithinExpt");
2013  weightsTuple->addDoubleBranch("dpModelWeight");
2014 
2015  UInt_t iExpmt = this->iExpt();
2016  UInt_t nExpmt = this->nExpt();
2017  UInt_t firstExpmt = this->firstExpt();
2018  for (iExpmt = firstExpmt; iExpmt < (firstExpmt+nExpmt); ++iExpmt) {
2019 
2020  inputFitData->readExperimentData(iExpmt);
2021  UInt_t nEvents = inputFitData->nEvents();
2022 
2023  if (nEvents < 1) {
2024  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Zero events in experiment " << iExpmt << ", skipping..." << std::endl;
2025  continue;
2026  }
2027 
2028  weightsTuple->setIntegerBranchValue( "iExpt", iExpmt );
2029 
2030  // Calculate and store the weights for the events in this experiment
2031  for ( UInt_t iEvent(0); iEvent < nEvents; ++iEvent ) {
2032 
2033  weightsTuple->setIntegerBranchValue( "iEvtWithinExpt", iEvent );
2034 
2035  const LauFitData& evtData = inputFitData->getData( iEvent );
2036 
2037  Double_t m13Sq_MC = evtData.find("m13Sq_MC")->second;
2038  Double_t m23Sq_MC = evtData.find("m23Sq_MC")->second;
2039 
2040  Double_t dpModelWeight(0.0);
2041 
2042  if ( kinematics_->withinDPLimits( m13Sq_MC, m23Sq_MC ) ) {
2043 
2044  kinematics_->updateKinematics( m13Sq_MC, m23Sq_MC );
2045  dpModelWeight = sigDPModel_->getEventWeight();
2046 
2047  if ( kinematics_->squareDP() ) {
2048  dpModelWeight *= kinematics_->calcSqDPJacobian();
2049  }
2050 
2051  if (LauAbsDPDynamics::GenOK != sigDPModel_->checkToyMC(kTRUE,kFALSE)) {
2052  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Problem in calculating the weight, aborting." << std::endl;
2053  delete weightsTuple;
2054  return;
2055  }
2056  }
2057 
2058  weightsTuple->setDoubleBranchValue( "dpModelWeight", dpModelWeight );
2059  weightsTuple->fillBranches();
2060  }
2061 
2062  }
2063 
2064  weightsTuple->buildIndex( "iExpt", "iEvtWithinExpt" );
2065  weightsTuple->addFriendTree(dataFileName, dataTreeName);
2066  weightsTuple->writeOutGenResults();
2067 
2068  delete weightsTuple;
2069 }
2070 
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.
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.
ClassImp(LauAbsCoeffSet)
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.
virtual void printParValues() const =0
Print the current values of the parameters.
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 storeCorrMatrix(UInt_t iExpt, Double_t NLL, Int_t fitStatus, const TMatrixD &covMatrix)
Store the correlation matrix and other fit information.
Definition: LauFitNtuple.cc:64
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:33
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:42
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:41
const TMatrixD & covarianceMatrix() const
Access the fit covariance matrix.
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.