laura is hosted by Hepforge, IPPP Durham
Laura++  v2r2p1
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, const Bool_t upperHalf, const Bool_t fluctuateBins, 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, fluctuateBins, 0.0, 0.0, upperHalf, daughters->squareDP() );
159 
160  scfMap_ = scfMap;
161 
162  useSCF_ = kTRUE;
163  useSCFHist_ = kTRUE;
164 }
165 
166 void LauSimpleFitModel::splitSignalComponent( const Double_t scfFrac, const 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  // Add the fit fraction that has not been corrected for the efficiency for each signal component
510  if (fitFracEffUnCorr_.size() != nSigComp_) {
511  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << fitFracEffUnCorr_.size() << std::endl;
512  gSystem->Exit(EXIT_FAILURE);
513  }
514  for (UInt_t i(0); i<nSigComp_; ++i) {
515  if (fitFracEffUnCorr_[i].size() != nSigComp_) {
516  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << fitFracEffUnCorr_[i].size() << std::endl;
517  gSystem->Exit(EXIT_FAILURE);
518  }
519  }
520 
521  for (UInt_t i = 0; i < nSigComp_; i++) {
522  for (UInt_t j = i; j < nSigComp_; j++) {
523  extraVars.push_back(fitFrac_[i][j]);
524  extraVars.push_back(fitFracEffUnCorr_[i][j]);
525  }
526  }
527 
528  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
529  std::vector<LauParameter> extraParams = sigDPModel_->getExtraParameters();
530  std::vector<LauParameter>::iterator extraIter;
531  for (extraIter = extraParams.begin(); extraIter != extraParams.end(); ++extraIter) {
532  LauParameter extraParameter = (*extraIter);
533  extraVars.push_back(extraParameter);
534  }
535 
536  // Now add in the DP efficiency value
537  Double_t initMeanEff = sigDPModel_->getMeanEff().initValue();
538  meanEff_.value(initMeanEff); meanEff_.initValue(initMeanEff); meanEff_.genValue(initMeanEff);
539  extraVars.push_back(meanEff_);
540 
541  // Also add in the DP rate
542  Double_t initDPRate = sigDPModel_->getDPRate().initValue();
543  dpRate_.value(initDPRate); dpRate_.initValue(initDPRate); dpRate_.genValue(initDPRate);
544  extraVars.push_back(dpRate_);
545 }
546 
547 void LauSimpleFitModel::finaliseFitResults(const TString& tablePrefixName)
548 {
549  // Retrieve parameters from the fit results for calculations and toy generation
550  // and eventually store these in output root ntuples/text files
551 
552  // Now take the fit parameters and update them as necessary
553  // e.g. to make mag > 0.0 and phase in the right range.
554  // This function will also calculate any other values, such as the
555  // fit fractions, using any errors provided by fitParErrors as appropriate.
556  // Also obtain the pull values: (measured - generated)/(average error)
557 
558  if (this->useDP() == kTRUE) {
559  for (UInt_t i = 0; i < nSigComp_; i++) {
560  // Check whether we have mag > 0.0, and phase in the right range
561  coeffPars_[i]->finaliseValues();
562  }
563  }
564 
565  // update the pulls on the events
566  if (this->doEMLFit()) {
568  }
569  if (useSCF_ && !useSCFHist_) {
571  }
572  if (usingBkgnd_ == kTRUE) {
573  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
574  (*iter)->updatePull();
575  }
576  }
577 
578  // Update the pulls on all the extra PDFs' parameters
580  if (useSCF_ == kTRUE) {
582  }
583  if (usingBkgnd_ == kTRUE) {
584  for (LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
585  this->updateFitParameters(*iter);
586  }
587  }
588 
589  // Fill the fit results to the ntuple for current experiment
590 
591  // update the coefficients and then calculate the fit fractions
592  if (this->useDP() == kTRUE) {
593  this->updateCoeffs();
596 
598  if (fitFrac.size() != nSigComp_) {
599  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFrac.size() << std::endl;
600  gSystem->Exit(EXIT_FAILURE);
601  }
602  for (UInt_t i(0); i<nSigComp_; ++i) {
603  if (fitFrac[i].size() != nSigComp_) {
604  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFrac[i].size() << std::endl;
605  gSystem->Exit(EXIT_FAILURE);
606  }
607  }
608 
610  if (fitFracEffUnCorr.size() != nSigComp_) {
611  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFracEffUnCorr.size() << std::endl;
612  gSystem->Exit(EXIT_FAILURE);
613  }
614  for (UInt_t i(0); i<nSigComp_; ++i) {
615  if (fitFracEffUnCorr[i].size() != nSigComp_) {
616  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << fitFracEffUnCorr[i].size() << std::endl;
617  gSystem->Exit(EXIT_FAILURE);
618  }
619  }
620 
621  for (UInt_t i(0); i<nSigComp_; ++i) {
622  for (UInt_t j(i); j<nSigComp_; ++j) {
623  fitFrac_[i][j].value(fitFrac[i][j].value());
624  fitFracEffUnCorr_[i][j].value(fitFracEffUnCorr[i][j].value());
625  }
626  }
627 
630 
631  this->clearExtraVarVectors();
632  LauParameterList& extraVars = this->extraPars();
633 
634  // Then store the final fit parameters, and any extra parameters for
635  // the signal model (e.g. fit fractions)
636  for (UInt_t i(0); i<nSigComp_; ++i) {
637  for (UInt_t j(i); j<nSigComp_; ++j) {
638  extraVars.push_back(fitFrac_[i][j]);
639  extraVars.push_back(fitFracEffUnCorr_[i][j]);
640  }
641  }
642 
643  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
644  std::vector<LauParameter> extraParams = sigDPModel_->getExtraParameters();
645  std::vector<LauParameter>::iterator extraIter;
646  for (extraIter = extraParams.begin(); extraIter != extraParams.end(); ++extraIter) {
647  LauParameter extraParameter = (*extraIter);
648  extraVars.push_back(extraParameter);
649  }
650 
651  // Now add in the DP efficiency value
652  extraVars.push_back(meanEff_);
653 
654  // Also add in the DP rate
655  extraVars.push_back(dpRate_);
656 
657  this->printFitFractions(std::cout);
658  }
659 
660  const LauParameterPList& fitVars = this->fitPars();
661  const LauParameterList& extraVars = this->extraPars();
662  LauFitNtuple* ntuple = this->fitNtuple();
663  ntuple->storeParsAndErrors(fitVars, extraVars);
664 
665  // find out the correlation matrix for the parameters
666  ntuple->storeCorrMatrix(this->iExpt(), this->nll(), this->fitStatus(), this->covarianceMatrix());
667 
668  // Fill the data into ntuple
669  ntuple->updateFitNtuple();
670 
671  // Print out the partial fit fractions, phases and the
672  // averaged efficiency, reweighted by the dynamics (and anything else)
673  if (this->writeLatexTable()) {
674  TString sigOutFileName(tablePrefixName);
675  sigOutFileName += "_"; sigOutFileName += this->iExpt(); sigOutFileName += "Expt.tex";
676  this->writeOutTable(sigOutFileName);
677  }
678 }
679 
680 void LauSimpleFitModel::printFitFractions(std::ostream& output)
681 {
682  // Print out Fit Fractions, total DP rate and mean efficiency
683  for (UInt_t i = 0; i < nSigComp_; i++) {
684  output << "FitFraction for component " << i << " (" << coeffPars_[i]->name() << ") = " << fitFrac_[i][i] << std::endl;
685  }
686  output << "Overall DP rate (integral of matrix element squared) = " << dpRate_ << std::endl;
687  output << "Average efficiency weighted by whole DP dynamics = " << meanEff_ << std::endl;
688 }
689 
690 void LauSimpleFitModel::writeOutTable(const TString& outputFile)
691 {
692  // Write out the results of the fit to a tex-readable table
693  // TODO - need to include the yields in this table
694  std::ofstream fout(outputFile);
695  LauPrint print;
696 
697  std::cout << "INFO in LauSimpleFitModel::writeOutTable : Writing out results of the fit to the tex file " << outputFile << std::endl;
698 
699  if (this->useDP() == kTRUE) {
700  // print the fit coefficients in one table
701  coeffPars_.front()->printTableHeading(fout);
702  for (UInt_t i = 0; i < nSigComp_; i++) {
703  coeffPars_[i]->printTableRow(fout);
704  }
705  fout << "\\hline" << std::endl;
706  fout << "\\end{tabular}" << std::endl << std::endl;
707 
708  // print the fit fractions in another
709  fout << "\\begin{tabular}{|l|c|}" << std::endl;
710  fout << "\\hline" << std::endl;
711  fout << "Component & FitFraction \\\\" << std::endl;
712  fout << "\\hline" << std::endl;
713  Double_t fitFracSum(0.0);
714  for (UInt_t i = 0; i < nSigComp_; i++) {
715  TString resName = coeffPars_[i]->name();
716  resName = resName.ReplaceAll("_", "\\_");
717  fout << resName << " & $";
718  Double_t fitFrac = fitFrac_[i][i].value();
719  fitFracSum += fitFrac;
720  print.printFormat(fout, fitFrac);
721  fout << "$ \\\\" << std::endl;
722  }
723  fout << "\\hline" << std::endl;
724 
725  // Also print out sum of fit fractions
726  fout << "Fit Fraction Sum & $";
727  print.printFormat(fout, fitFracSum);
728  fout << "$ \\\\" << std::endl;
729  fout << "\\hline" << std::endl;
730 
731  fout << "DP rate & $";
732  print.printFormat(fout, dpRate_.value());
733  fout << "$ \\\\" << std::endl;
734  fout << "\\hline" << std::endl;
735 
736  fout << "$< \\varepsilon >$ & $";
737  print.printFormat(fout, meanEff_.value());
738  fout << "$ \\\\" << std::endl;
739  fout << "\\hline" << std::endl;
740  fout << "\\end{tabular}" << std::endl << std::endl;
741  }
742 
743  if (!signalPdfs_.empty()) {
744  fout << "\\begin{tabular}{|l|c|}" << std::endl;
745  fout << "\\hline" << std::endl;
746  if (useSCF_ == kTRUE) {
747  fout << "\\Extra TM Signal PDFs' Parameters: & \\\\" << std::endl;
748  } else {
749  fout << "\\Extra Signal PDFs' Parameters: & \\\\" << std::endl;
750  }
751  this->printFitParameters(signalPdfs_, fout);
752  if (useSCF_ == kTRUE && !scfPdfs_.empty()) {
753  fout << "\\hline" << std::endl;
754  fout << "\\Extra SCF Signal PDFs' Parameters: & \\\\" << std::endl;
755  this->printFitParameters(scfPdfs_, fout);
756  }
757  if (usingBkgnd_ == kTRUE && !bkgndPdfs_.empty()) {
758  fout << "\\hline" << std::endl;
759  fout << "\\Extra Background PDFs' Parameters: & \\\\" << std::endl;
760  for (LauBkgndPdfsList::const_iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
761  this->printFitParameters(*iter, fout);
762  }
763  }
764  fout << "\\hline \n\\end{tabular}" << std::endl << std::endl;
765  }
766 }
767 
769 {
770  // Update the number of signal events to be total-sum(background events)
771  this->updateSigEvents();
772 
773  // Check whether we want to have randomised initial fit parameters for the signal model
774  if (this->useRandomInitFitPars() == kTRUE) {
775  std::cout << "INFO in LauSimpleFitModel::checkInitFitParams : Setting random parameters for the signal DP model" << std::endl;
776  this->randomiseInitFitPars();
777  }
778 }
779 
781 {
782  // Only randomise those parameters that are not fixed!
783  std::cout << "INFO in LauSimpleFitModel::randomiseInitFitPars : Randomising the initial fit magnitudes and phases of the resonances..." << std::endl;
784 
785  for (UInt_t i = 0; i < nSigComp_; i++) {
786  coeffPars_[i]->randomiseInitValues();
787  }
788 }
789 
791 {
792  // Determine the number of events to generate for each hypothesis
793  // If we're smearing then smear each one individually
794 
795  LauGenInfo nEvtsGen;
796 
797  // Signal
798  Double_t evtWeight(1.0);
799  Int_t nEvts = TMath::FloorNint(signalEvents_->genValue());
800  if ( nEvts < 0 ) {
801  evtWeight = -1.0;
802  nEvts = TMath::Abs( nEvts );
803  }
804  if (this->doPoissonSmearing()) {
805  nEvts = LauRandom::randomFun()->Poisson(nEvts);
806  }
807  nEvtsGen["signal"] = std::make_pair( nEvts, evtWeight );
808 
809  // Backgrounds
810  const UInt_t nBkgnds = this->nBkgndClasses();
811  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
812  const TString& bkgndClass = this->bkgndClassName(bkgndID);
813  const LauParameter* evtsPar = bkgndEvents_[bkgndID];
814  evtWeight = 1.0;
815  nEvts = TMath::FloorNint( evtsPar->genValue() );
816  if ( nEvts < 0 ) {
817  evtWeight = -1.0;
818  nEvts = TMath::Abs( nEvts );
819  }
820  if (this->doPoissonSmearing()) {
821  nEvts = LauRandom::randomFun()->Poisson(nEvts);
822  }
823  nEvtsGen[bkgndClass] = std::make_pair( nEvts, evtWeight );
824  }
825 
826  return nEvtsGen;
827 }
828 
830 {
831  // Routine to generate toy Monte Carlo events according to the various models we have defined.
832 
833  // Determine the number of events to generate for each hypothesis
834  LauGenInfo nEvts = this->eventsToGenerate();
835 
836  Bool_t genOK(kTRUE);
837  Int_t evtNum(0);
838 
839  const UInt_t nBkgnds = this->nBkgndClasses();
840  std::vector<TString> bkgndClassNames(nBkgnds);
841  std::vector<TString> bkgndClassNamesGen(nBkgnds);
842  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
843  TString name( this->bkgndClassName(iBkgnd) );
844  bkgndClassNames[iBkgnd] = name;
845  bkgndClassNamesGen[iBkgnd] = "gen"+name;
846  }
847 
848  const Bool_t storeSCFTruthInfo = ( useSCF_ || ( this->enableEmbedding() && signalTree_ != 0 && signalTree_->haveBranch("mcMatch") ) );
849 
850  // Loop over the hypotheses and generate the requested number of events for each one
851  for (LauGenInfo::const_iterator iter = nEvts.begin(); iter != nEvts.end(); ++iter) {
852 
853  const TString& type(iter->first);
854  Int_t nEvtsGen( iter->second.first );
855  Double_t evtWeight( iter->second.second );
856 
857  for (Int_t iEvt(0); iEvt<nEvtsGen; ++iEvt) {
858 
859  this->setGenNtupleDoubleBranchValue( "evtWeight", evtWeight );
860  // Add efficiency information
861  this->setGenNtupleDoubleBranchValue( "efficiency", 1 );
862 
863  if (type == "signal") {
864  this->setGenNtupleIntegerBranchValue("genSig",1);
865  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
866  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], 0 );
867  }
868  genOK = this->generateSignalEvent();
869  this->setGenNtupleDoubleBranchValue( "efficiency", sigDPModel_->getEvtEff() );
870  } else {
871  this->setGenNtupleIntegerBranchValue("genSig",0);
872  if ( storeSCFTruthInfo ) {
873  this->setGenNtupleIntegerBranchValue("genTMSig",0);
874  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
875  }
876  UInt_t bkgndID(0);
877  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
878  Int_t gen(0);
879  if ( bkgndClassNames[iBkgnd] == type ) {
880  gen = 1;
881  bkgndID = iBkgnd;
882  }
883  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], gen );
884  }
885  genOK = this->generateBkgndEvent(bkgndID);
886  }
887 
888  if (!genOK) {
889  // If there was a problem with the generation then break out and return.
890  // The problem model will have adjusted itself so that all should be OK next time.
891  break;
892  }
893 
894  if (this->useDP() == kTRUE) {
895  this->setDPBranchValues();
896  }
897 
898  // Store the event number (within this experiment)
899  this->setGenNtupleIntegerBranchValue("iEvtWithinExpt",evtNum);
900  // and then increment it
901  ++evtNum;
902 
903  this->fillGenNtupleBranches();
904  if (iEvt%500 == 0) {
905  std::cout << "INFO in LauSimpleFitModel::genExpt : Generated event number " << iEvt << " out of " << nEvtsGen << " " << type << " events." << std::endl;
906  }
907  }
908 
909  if (!genOK) {
910  break;
911  }
912  }
913 
914  if (this->useDP() && genOK) {
915 
916  sigDPModel_->checkToyMC(kTRUE,kTRUE);
917 
918  // Get the fit fractions if they're to be written into the latex table
919  if (this->writeLatexTable()) {
921  if (fitFrac.size() != nSigComp_) {
922  std::cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: " << fitFrac.size() << std::endl;
923  gSystem->Exit(EXIT_FAILURE);
924  }
925  for (UInt_t i(0); i<nSigComp_; ++i) {
926  if (fitFrac[i].size() != nSigComp_) {
927  std::cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: " << fitFrac[i].size() << std::endl;
928  gSystem->Exit(EXIT_FAILURE);
929  }
930  }
931 
932  for (UInt_t i(0); i<nSigComp_; ++i) {
933  for (UInt_t j = i; j < nSigComp_; j++) {
934  fitFrac_[i][j].value(fitFrac[i][j].value());
935  }
936  }
939  }
940  }
941 
942  // If we're reusing embedded events or if the generation is being
943  // reset then clear the lists of used events
944  if (reuseSignal_ || !genOK) {
945  if (signalTree_) {
947  }
948  }
949  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
950  LauEmbeddedData* data = bkgndTree_[bkgndID];
951  if (reuseBkgnd_[bkgndID] || !genOK) {
952  if (data) {
953  data->clearUsedList();
954  }
955  }
956  }
957 
958  return genOK;
959 }
960 
962 {
963  // Generate signal event
964  Bool_t genOK(kTRUE);
965  Bool_t genSCF(kFALSE);
966 
967  if (this->useDP()) {
968  if (this->enableEmbedding() && signalTree_) {
969  if (useReweighting_) {
970  // Select a (random) event from the generated data. Then store the
971  // reconstructed DP co-ords, together with other pdf information,
972  // as the embedded data.
974  } else {
975 
977  }
978  if (signalTree_->haveBranch("mcMatch")) {
979  Int_t match = static_cast<Int_t>(signalTree_->getValue("mcMatch"));
980  if (match) {
981  this->setGenNtupleIntegerBranchValue("genTMSig",1);
982  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
983  genSCF = kFALSE;
984  } else {
985  this->setGenNtupleIntegerBranchValue("genTMSig",0);
986  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
987  genSCF = kTRUE;
988  }
989  }
990  } else {
991  genOK = sigDPModel_->generate();
992  if ( genOK && useSCF_ ) {
993  Double_t frac(0.0);
994  if ( useSCFHist_ ) {
996  } else {
997  frac = scfFrac_.genValue();
998  }
999  if ( frac < LauRandom::randomFun()->Rndm() ) {
1000  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1001  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1002  genSCF = kFALSE;
1003  } else {
1004  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1005  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1006  genSCF = kTRUE;
1007 
1008  // Optionally smear the DP position
1009  // of the SCF event
1010  if ( scfMap_ != 0 ) {
1011  Double_t xCoord(0.0), yCoord(0.0);
1012  if ( kinematics_->squareDP() ) {
1013  xCoord = kinematics_->getmPrime();
1014  yCoord = kinematics_->getThetaPrime();
1015  } else {
1016  xCoord = kinematics_->getm13Sq();
1017  yCoord = kinematics_->getm23Sq();
1018  }
1019 
1020  // Find the bin number where this event is generated
1021  Int_t binNo = scfMap_->binNumber( xCoord, yCoord );
1022 
1023  // Retrieve the migration histogram
1024  TH2* histo = scfMap_->trueHist( binNo );
1025 
1026  LauAbsEffModel * effModel = sigDPModel_->getEffModel();
1027  do {
1028  // Get a random point from the histogram
1029  histo->GetRandom2( xCoord, yCoord );
1030 
1031  // Update the kinematics
1032  if ( kinematics_->squareDP() ) {
1033  kinematics_->updateSqDPKinematics( xCoord, yCoord );
1034  } else {
1035  kinematics_->updateKinematics( xCoord, yCoord );
1036  }
1037  } while ( ! effModel->passVeto( kinematics_ ) );
1038  }
1039  }
1040  }
1041  }
1042  } else {
1043  if (this->enableEmbedding() && signalTree_) {
1045  if (signalTree_->haveBranch("mcMatch")) {
1046  Int_t match = static_cast<Int_t>(signalTree_->getValue("mcMatch"));
1047  if (match) {
1048  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1049  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1050  genSCF = kFALSE;
1051  } else {
1052  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1053  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1054  genSCF = kTRUE;
1055  }
1056  }
1057  } else if ( useSCF_ ) {
1058  Double_t frac = scfFrac_.genValue();
1059  if ( frac < LauRandom::randomFun()->Rndm() ) {
1060  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1061  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1062  genSCF = kFALSE;
1063  } else {
1064  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1065  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1066  genSCF = kTRUE;
1067  }
1068  }
1069  }
1070  if (genOK) {
1071  if ( useSCF_ ) {
1072  if ( genSCF ) {
1074  } else {
1076  }
1077  } else {
1079  }
1080  }
1081  // Check for problems with the embedding
1082  if (this->enableEmbedding() && signalTree_ && (signalTree_->nEvents() == signalTree_->nUsedEvents())) {
1083  std::cerr << "WARNING in LauSimpleFitModel::generateSignalEvent : Source of embedded signal events used up, clearing the list of used events." << std::endl;
1085  }
1086  return genOK;
1087 }
1088 
1090 {
1091  // Generate background event
1092  Bool_t genOK(kTRUE);
1093 
1094  LauAbsBkgndDPModel* model = bkgndDPModels_[bkgndID];
1095  LauEmbeddedData* embeddedData(0);
1096  if (this->enableEmbedding()) {
1097  embeddedData = bkgndTree_[bkgndID];
1098  }
1099  LauPdfList* extraPdfs = &bkgndPdfs_[bkgndID];
1100  if (this->useDP()) {
1101  if (embeddedData) {
1102  embeddedData->getEmbeddedEvent(kinematics_);
1103  } else {
1104  if (model == 0) {
1105  std::cerr << "ERROR in LauSimpleFitModel::generateBkgndEvent : Can't find the DP model for background class \"" << this->bkgndClassName(bkgndID) << "\"." << std::endl;
1106  gSystem->Exit(EXIT_FAILURE);
1107  }
1108  genOK = model->generate();
1109  }
1110  } else {
1111  if (embeddedData) {
1112  embeddedData->getEmbeddedEvent(0);
1113  }
1114  }
1115  if (genOK) {
1116  this->generateExtraPdfValues(extraPdfs, embeddedData);
1117  }
1118  // Check for problems with the embedding
1119  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1120  std::cerr << "WARNING in LauSimpleFitModel::generateBkgndEvent : Source of embedded " << this->bkgndClassName(bkgndID) << " events used up, clearing the list of used events." << std::endl;
1121  embeddedData->clearUsedList();
1122  }
1123  return genOK;
1124 }
1125 
1127 {
1128  // Setup the required ntuple branches
1129  this->addGenNtupleDoubleBranch("evtWeight");
1130  this->addGenNtupleIntegerBranch("genSig");
1131  this->addGenNtupleDoubleBranch("efficiency");
1132  if ( useSCF_ || ( this->enableEmbedding() && signalTree_ != 0 && signalTree_->haveBranch("mcMatch") ) ) {
1133  this->addGenNtupleIntegerBranch("genTMSig");
1134  this->addGenNtupleIntegerBranch("genSCFSig");
1135  }
1136  const UInt_t nBkgnds = this->nBkgndClasses();
1137  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1138  TString name( this->bkgndClassName(iBkgnd) );
1139  name.Prepend("gen");
1140  this->addGenNtupleIntegerBranch(name);
1141  }
1142  if (this->useDP() == kTRUE) {
1143  this->addGenNtupleDoubleBranch("m12");
1144  this->addGenNtupleDoubleBranch("m23");
1145  this->addGenNtupleDoubleBranch("m13");
1146  this->addGenNtupleDoubleBranch("m12Sq");
1147  this->addGenNtupleDoubleBranch("m23Sq");
1148  this->addGenNtupleDoubleBranch("m13Sq");
1149  this->addGenNtupleDoubleBranch("cosHel12");
1150  this->addGenNtupleDoubleBranch("cosHel23");
1151  this->addGenNtupleDoubleBranch("cosHel13");
1152  if (kinematics_->squareDP()) {
1153  this->addGenNtupleDoubleBranch("mPrime");
1154  this->addGenNtupleDoubleBranch("thPrime");
1155  }
1156  }
1157  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1158  std::vector<TString> varNames = (*pdf_iter)->varNames();
1159  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1160  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1161  this->addGenNtupleDoubleBranch( (*var_iter) );
1162  }
1163  }
1164  }
1165 }
1166 
1168 {
1169  // Store all the DP information
1176  this->setGenNtupleDoubleBranchValue("cosHel12", kinematics_->getc12());
1177  this->setGenNtupleDoubleBranchValue("cosHel23", kinematics_->getc23());
1178  this->setGenNtupleDoubleBranchValue("cosHel13", kinematics_->getc13());
1179  if (kinematics_->squareDP()) {
1182  }
1183 }
1184 
1186 {
1187  // Generate from the extra PDFs
1188  if (extraPdfs) {
1189  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1190  LauFitData genValues;
1191  if (embeddedData) {
1192  genValues = embeddedData->getValues( (*pdf_iter)->varNames() );
1193  } else {
1194  genValues = (*pdf_iter)->generate(kinematics_);
1195  }
1196  for ( LauFitData::const_iterator var_iter = genValues.begin(); var_iter != genValues.end(); ++var_iter ) {
1197  TString varName = var_iter->first;
1198  if ( varName != "m13Sq" && varName != "m23Sq" ) {
1199  Double_t value = var_iter->second;
1200  this->setGenNtupleDoubleBranchValue(varName,value);
1201  }
1202  }
1203  }
1204  }
1205 }
1206 
1208 {
1209  // Update the total normalisation for the signal likelihood
1210  if (this->useDP() == kTRUE) {
1211  this->updateCoeffs();
1213  }
1214 
1215  // Update the signal events from the background events if not doing an extended fit
1216  if ( !this->doEMLFit() && !signalEvents_->fixed() ) {
1217  this->updateSigEvents();
1218  }
1219 }
1220 
1222 {
1223  // The background parameters will have been set from Minuit.
1224  // We need to update the signal events using these.
1225  Double_t nTotEvts = this->eventsPerExpt();
1226 
1227  signalEvents_->range(-2.0*nTotEvts,2.0*nTotEvts);
1228  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1229  (*iter)->range(-2.0*nTotEvts,2.0*nTotEvts);
1230  }
1231 
1232  if ( signalEvents_->fixed() ) {
1233  return;
1234  }
1235 
1236  // Subtract background events (if any) from signal.
1237  Double_t signalEvents = nTotEvts;
1238  if ( usingBkgnd_ == kTRUE ) {
1239  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1240  signalEvents -= (*iter)->value();
1241  }
1242  }
1243  signalEvents_->value(signalEvents);
1244 }
1245 
1247 {
1248  // Fill the internal data trees of the signal and backgrond models.
1249 
1250  LauFitDataTree* inputFitData = this->fitData();
1251 
1252  // First the Dalitz plot variables (m_ij^2)
1253  if (this->useDP() == kTRUE) {
1254 
1255  // need to append SCF smearing bins before caching DP amplitudes
1256  if ( scfMap_ != 0 ) {
1257  this->appendBinCentres( inputFitData );
1258  }
1259  sigDPModel_->fillDataTree(*inputFitData);
1260 
1261  if (usingBkgnd_ == kTRUE) {
1262  for (LauBkgndDPModelList::iterator iter = bkgndDPModels_.begin(); iter != bkgndDPModels_.end(); ++iter) {
1263  (*iter)->fillDataTree(*inputFitData);
1264  }
1265  }
1266  }
1267 
1268  // ...and then the extra PDFs
1269  this->cacheInfo(signalPdfs_, *inputFitData);
1270  this->cacheInfo(scfPdfs_, *inputFitData);
1271  for (LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter) {
1272  this->cacheInfo((*iter), *inputFitData);
1273  }
1274 
1275  // and finally the SCF fractions and jacobians
1276  if ( useSCF_ && useSCFHist_ ) {
1277  if ( !inputFitData->haveBranch( "m13Sq" ) || !inputFitData->haveBranch( "m23Sq" ) ) {
1278  std::cerr << "ERROR in LauSimpleFitModel::cacheInputFitVars : Input data does not contain DP branches and so can't cache the SCF fraction." << std::endl;
1279  gSystem->Exit(EXIT_FAILURE);
1280  }
1281 
1282  UInt_t nEvents = inputFitData->nEvents();
1283  recoSCFFracs_.clear();
1284  recoSCFFracs_.reserve( nEvents );
1285  if ( kinematics_->squareDP() ) {
1286  recoJacobians_.clear();
1287  recoJacobians_.reserve( nEvents );
1288  }
1289  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1290  const LauFitData& dataValues = inputFitData->getData(iEvt);
1291  LauFitData::const_iterator m13_iter = dataValues.find("m13Sq");
1292  LauFitData::const_iterator m23_iter = dataValues.find("m23Sq");
1293  kinematics_->updateKinematics( m13_iter->second, m23_iter->second );
1294  Double_t scfFrac = scfFracHist_->calcEfficiency( kinematics_ );
1295  recoSCFFracs_.push_back( scfFrac );
1296  if ( kinematics_->squareDP() ) {
1298  }
1299  }
1300  }
1301 }
1302 
1304 {
1305  // We'll be caching the DP amplitudes and efficiencies of the centres of the true bins.
1306  // To do so, we attach some fake points at the end of inputData, the number of the entry
1307  // minus the total number of events corresponding to the number of the histogram for that
1308  // given true bin in the LauScfMap object. (What this means is that when Laura is provided with
1309  // the LauScfMap object by the user, it's the latter who has to make sure that it contains the
1310  // right number of histograms and in exactly the right order!)
1311 
1312  // Get the x and y co-ordinates of the bin centres
1313  std::vector<Double_t> binCentresXCoords;
1314  std::vector<Double_t> binCentresYCoords;
1315  scfMap_->listBinCentres(binCentresXCoords, binCentresYCoords);
1316 
1317  // The SCF histograms could be in square Dalitz plot histograms.
1318  // The dynamics takes normal Dalitz plot coords, so we might have to convert them back.
1319  Bool_t sqDP = kinematics_->squareDP();
1320  UInt_t nBins = binCentresXCoords.size();
1321  fakeSCFFracs_.clear();
1322  fakeSCFFracs_.reserve( nBins );
1323  if ( sqDP ) {
1324  fakeJacobians_.clear();
1325  fakeJacobians_.reserve( nBins );
1326  }
1327 
1328  for (UInt_t iBin = 0; iBin < nBins; ++iBin) {
1329 
1330  if ( sqDP ) {
1331  kinematics_->updateSqDPKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1332  binCentresXCoords[iBin] = kinematics_->getm13Sq();
1333  binCentresYCoords[iBin] = kinematics_->getm23Sq();
1335  } else {
1336  kinematics_->updateKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1337  }
1338 
1340  }
1341 
1342  // Set up inputFitVars_ object to hold the fake events
1343  inputData->appendFakePoints(binCentresXCoords,binCentresYCoords);
1344 }
1345 
1347 {
1348  // Get the DP likelihood for signal and backgrounds
1349  this->getEvtDPLikelihood(iEvt);
1350 
1351  // Get the combined extra PDFs likelihood for signal and backgrounds
1352  this->getEvtExtraLikelihoods(iEvt);
1353 
1354  // If appropriate, combine the TM and SCF likelihoods
1355  Double_t sigLike = sigDPLike_ * sigExtraLike_;
1356  if ( useSCF_ ) {
1357  Double_t scfFrac(0.0);
1358  if (useSCFHist_) {
1359  scfFrac = recoSCFFracs_[iEvt];
1360  } else {
1361  scfFrac = scfFrac_.value();
1362  }
1363  sigLike *= (1.0 - scfFrac);
1364  if ( (scfMap_ != 0) && (this->useDP() == kTRUE) ) {
1365  // if we're smearing the SCF DP PDF then the SCF frac
1366  // is already included in the SCF DP likelihood
1367  sigLike += (scfDPLike_ * scfExtraLike_);
1368  } else {
1369  sigLike += (scfFrac * scfDPLike_ * scfExtraLike_);
1370  }
1371  }
1372 
1373  // Construct the total event likelihood
1374  Double_t likelihood = signalEvents_->value() * sigLike;
1375  const UInt_t nBkgnds = this->nBkgndClasses();
1376  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1377  likelihood += (bkgndEvents_[bkgndID]->value() * bkgndDPLike_[bkgndID] * bkgndExtraLike_[bkgndID]);
1378  }
1379 
1380  return likelihood;
1381 }
1382 
1384 {
1385  Double_t eventSum(0.0);
1386  eventSum += signalEvents_->value();
1387  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1388  eventSum += (*iter)->value();
1389  }
1390  return eventSum;
1391 }
1392 
1394 {
1395  // Function to return the signal and background likelihoods for the
1396  // Dalitz plot for the given event evtNo.
1397 
1398  if (this->useDP() == kTRUE) {
1401 
1402  if ( useSCF_ == kTRUE ) {
1403  if ( scfMap_ == 0 ) {
1404  // we're not smearing the SCF DP position
1405  // so the likelihood is the same as the TM
1407  } else {
1408  // calculate the smeared SCF DP likelihood
1409  scfDPLike_ = this->getEvtSCFDPLikelihood(iEvt);
1410  }
1411  }
1412 
1413  const UInt_t nBkgnds = this->nBkgndClasses();
1414  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1415  if (usingBkgnd_ == kTRUE) {
1416  bkgndDPLike_[bkgndID] = bkgndDPModels_[bkgndID]->getLikelihood(iEvt);
1417  } else {
1418  bkgndDPLike_[bkgndID] = 0.0;
1419  }
1420  }
1421  } else {
1422  // There's always going to be a term in the likelihood for the
1423  // signal, so we'd better not zero it.
1424  sigDPLike_ = 1.0;
1425  scfDPLike_ = 1.0;
1426 
1427  const UInt_t nBkgnds = this->nBkgndClasses();
1428  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1429  if (usingBkgnd_ == kTRUE) {
1430  bkgndDPLike_[bkgndID] = 1.0;
1431  } else {
1432  bkgndDPLike_[bkgndID] = 0.0;
1433  }
1434  }
1435  }
1436 }
1437 
1439 {
1440  Double_t scfDPLike(0.0);
1441 
1442  Double_t recoJacobian(1.0);
1443  Double_t xCoord(0.0);
1444  Double_t yCoord(0.0);
1445  Bool_t squareDP = kinematics_->squareDP();
1446  if ( squareDP ) {
1447  xCoord = sigDPModel_->getEvtmPrime();
1448  yCoord = sigDPModel_->getEvtthPrime();
1449  recoJacobian = recoJacobians_[iEvt];
1450  } else {
1451  xCoord = sigDPModel_->getEvtm13Sq();
1452  yCoord = sigDPModel_->getEvtm23Sq();
1453  }
1454 
1455  // Find the bin that our reco event falls in
1456  Int_t recoBin = scfMap_->binNumber( xCoord, yCoord );
1457 
1458  // Find out which true Bins contribute to the given reco bin
1459  const std::vector<Int_t>* trueBins = scfMap_->trueBins(recoBin);
1460 
1461  Int_t nDataEvents = this->eventsPerExpt();
1462 
1463  // Loop over the true bins
1464  for (std::vector<Int_t>::const_iterator iter = trueBins->begin(); iter != trueBins->end(); ++iter)
1465  {
1466  Int_t trueBin = (*iter);
1467 
1468  // prob of a true event in the given true bin migrating to the reco bin
1469  Double_t pRecoGivenTrue = scfMap_->prob( recoBin, trueBin );
1470 
1471  // We've cached the DP amplitudes and the efficiency for the
1472  // true bin centres, just after the data points
1473  sigDPModel_->calcLikelihoodInfo( nDataEvents + trueBin );
1474  Double_t pTrue = sigDPModel_->getEvtLikelihood();
1475 
1476  // Get the cached SCF fraction (and jacobian if we're using the square DP)
1477  Double_t scfFraction = fakeSCFFracs_[ trueBin ];
1478  Double_t jacobian(1.0);
1479  if ( squareDP ) {
1480  jacobian = fakeJacobians_[ trueBin ];
1481  }
1482 
1483  scfDPLike += pTrue * jacobian * scfFraction * pRecoGivenTrue;
1484 
1485  }
1486 
1487  // Divide by the reco jacobian
1488  scfDPLike /= recoJacobian;
1489 
1490  return scfDPLike;
1491 }
1492 
1494 {
1495  // Function to return the signal and background likelihoods for the
1496  // extra variables for the given event evtNo.
1497 
1498  sigExtraLike_ = 1.0;
1499 
1500  for (LauPdfList::iterator iter = signalPdfs_.begin(); iter != signalPdfs_.end(); ++iter) {
1501  (*iter)->calcLikelihoodInfo(iEvt);
1502  sigExtraLike_ *= (*iter)->getLikelihood();
1503  }
1504  if (useSCF_) {
1505  scfExtraLike_ = 1.0;
1506  for (LauPdfList::iterator iter = scfPdfs_.begin(); iter != scfPdfs_.end(); ++iter) {
1507  (*iter)->calcLikelihoodInfo(iEvt);
1508  scfExtraLike_ *= (*iter)->getLikelihood();
1509  }
1510  }
1511 
1512  const UInt_t nBkgnds = this->nBkgndClasses();
1513  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1514  if (usingBkgnd_) {
1515  bkgndExtraLike_[bkgndID] = 1.0;
1516  LauPdfList& pdfList = bkgndPdfs_[bkgndID];
1517  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1518  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1519  bkgndExtraLike_[bkgndID] *= (*pdf_iter)->getLikelihood();
1520  }
1521  } else {
1522  bkgndExtraLike_[bkgndID] = 0.0;
1523  }
1524  }
1525 }
1526 
1528 {
1529  coeffs_.clear();
1530  coeffs_.reserve(nSigComp_);
1531  for (UInt_t i = 0; i < nSigComp_; i++) {
1532  coeffs_.push_back(coeffPars_[i]->particleCoeff());
1533  }
1534 }
1535 
1537 {
1538  // add branches for storing the experiment number and the number of
1539  // the event within the current experiment
1540  this->addSPlotNtupleIntegerBranch("iExpt");
1541  this->addSPlotNtupleIntegerBranch("iEvtWithinExpt");
1542 
1543  // Store the efficiency of the event (for inclusive BF calculations).
1544  if (this->storeDPEff()) {
1545  this->addSPlotNtupleDoubleBranch("efficiency");
1546  if ( sigDPModel_->usingScfModel() ) {
1547  this->addSPlotNtupleDoubleBranch("scffraction");
1548  }
1549  }
1550 
1551  // Store the total event likelihood for each species.
1552  if (useSCF_) {
1553  this->addSPlotNtupleDoubleBranch("sigTMTotalLike");
1554  this->addSPlotNtupleDoubleBranch("sigSCFTotalLike");
1555  this->addSPlotNtupleDoubleBranch("sigSCFFrac");
1556  } else {
1557  this->addSPlotNtupleDoubleBranch("sigTotalLike");
1558  }
1559  if (usingBkgnd_) {
1560  const UInt_t nBkgnds = this->nBkgndClasses();
1561  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1562  TString name( this->bkgndClassName(iBkgnd) );
1563  name += "TotalLike";
1564  this->addSPlotNtupleDoubleBranch(name);
1565  }
1566  }
1567 
1568  // Store the DP likelihoods
1569  if (this->useDP()) {
1570  if (useSCF_) {
1571  this->addSPlotNtupleDoubleBranch("sigTMDPLike");
1572  this->addSPlotNtupleDoubleBranch("sigSCFDPLike");
1573  } else {
1574  this->addSPlotNtupleDoubleBranch("sigDPLike");
1575  }
1576  if (usingBkgnd_) {
1577  const UInt_t nBkgnds = this->nBkgndClasses();
1578  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1579  TString name( this->bkgndClassName(iBkgnd) );
1580  name += "DPLike";
1581  this->addSPlotNtupleDoubleBranch(name);
1582  }
1583  }
1584  }
1585 
1586  // Store the likelihoods for each extra PDF
1587  if (useSCF_) {
1588  this->addSPlotNtupleBranches(&signalPdfs_, "sigTM");
1589  this->addSPlotNtupleBranches(&scfPdfs_, "sigSCF");
1590  } else {
1591  this->addSPlotNtupleBranches(&signalPdfs_, "sig");
1592  }
1593  if (usingBkgnd_) {
1594  const UInt_t nBkgnds = this->nBkgndClasses();
1595  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1596  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1597  const LauPdfList* pdfList = &(bkgndPdfs_[iBkgnd]);
1598  this->addSPlotNtupleBranches(pdfList, bkgndClass);
1599  }
1600  }
1601 }
1602 
1603 void LauSimpleFitModel::addSPlotNtupleBranches(const LauPdfList* extraPdfs, const TString& prefix)
1604 {
1605  if (extraPdfs) {
1606  // Loop through each of the PDFs
1607  for (LauPdfList::const_iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1608 
1609  // Count the number of input variables that are not
1610  // DP variables (used in the case where there is DP
1611  // dependence for e.g. MVA)
1612  UInt_t nVars(0);
1613  std::vector<TString> varNames = (*pdf_iter)->varNames();
1614  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1615  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1616  ++nVars;
1617  }
1618  }
1619 
1620  if ( nVars == 1 ) {
1621  // If the PDF only has one variable then
1622  // simply add one branch for that variable
1623  TString varName = (*pdf_iter)->varName();
1624  TString name(prefix);
1625  name += varName;
1626  name += "Like";
1627  this->addSPlotNtupleDoubleBranch(name);
1628  } else if ( nVars == 2 ) {
1629  // If the PDF has two variables then we
1630  // need a branch for them both together and
1631  // branches for each
1632  TString allVars("");
1633  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1634  allVars += (*var_iter);
1635  TString name(prefix);
1636  name += (*var_iter);
1637  name += "Like";
1638  this->addSPlotNtupleDoubleBranch(name);
1639  }
1640  TString name(prefix);
1641  name += allVars;
1642  name += "Like";
1643  this->addSPlotNtupleDoubleBranch(name);
1644  } else {
1645  std::cerr << "WARNING in LauSimpleFitModel::addSPlotNtupleBranches : Can't yet deal with 3D PDFs." << std::endl;
1646  }
1647  }
1648  }
1649 }
1650 
1651 Double_t LauSimpleFitModel::setSPlotNtupleBranchValues(LauPdfList* extraPdfs, const TString& prefix, UInt_t iEvt)
1652 {
1653  // Store the PDF value for each variable in the list
1654  Double_t totalLike(1.0);
1655  Double_t extraLike(0.0);
1656  if (extraPdfs) {
1657  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1658 
1659  // calculate the likelihood for this event
1660  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1661  extraLike = (*pdf_iter)->getLikelihood();
1662  totalLike *= extraLike;
1663 
1664  // Count the number of input variables that are not
1665  // DP variables (used in the case where there is DP
1666  // dependence for e.g. MVA)
1667  UInt_t nVars(0);
1668  std::vector<TString> varNames = (*pdf_iter)->varNames();
1669  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1670  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1671  ++nVars;
1672  }
1673  }
1674 
1675  if ( nVars == 1 ) {
1676  // If the PDF only has one variable then
1677  // simply store the value for that variable
1678  TString varName = (*pdf_iter)->varName();
1679  TString name(prefix);
1680  name += varName;
1681  name += "Like";
1682  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
1683  } else if ( nVars == 2 ) {
1684  // If the PDF has two variables then we
1685  // store the value for them both together
1686  // and for each on their own
1687  TString allVars("");
1688  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1689  allVars += (*var_iter);
1690  TString name(prefix);
1691  name += (*var_iter);
1692  name += "Like";
1693  Double_t indivLike = (*pdf_iter)->getLikelihood( (*var_iter) );
1694  this->setSPlotNtupleDoubleBranchValue(name, indivLike);
1695  }
1696  TString name(prefix);
1697  name += allVars;
1698  name += "Like";
1699  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
1700  } else {
1701  std::cerr << "WARNING in LauSimpleFitModel::setSPlotNtupleBranchValues : Can't yet deal with 3D PDFs." << std::endl;
1702  }
1703  }
1704  }
1705  return totalLike;
1706 }
1707 
1709 {
1710  LauSPlot::NameSet nameSet;
1711  if (this->useDP()) {
1712  nameSet.insert("DP");
1713  }
1714  // Loop through all the signal PDFs
1715  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1716  // Loop over the variables involved in each PDF
1717  std::vector<TString> varNames = (*pdf_iter)->varNames();
1718  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1719  // If they are not DP coordinates then add them
1720  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1721  nameSet.insert( (*var_iter) );
1722  }
1723  }
1724  }
1725  return nameSet;
1726 }
1727 
1729 {
1730  LauSPlot::NumbMap numbMap;
1731  if (!signalEvents_->fixed() && this->doEMLFit()) {
1732  numbMap["sig"] = signalEvents_->genValue();
1733  }
1734  if ( usingBkgnd_ ) {
1735  const UInt_t nBkgnds = this->nBkgndClasses();
1736  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1737  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1738  const LauParameter* par = bkgndEvents_[iBkgnd];
1739  if (!par->fixed()) {
1740  numbMap[bkgndClass] = par->genValue();
1741  }
1742  }
1743  }
1744  return numbMap;
1745 }
1746 
1748 {
1749  LauSPlot::NumbMap numbMap;
1750  if (signalEvents_->fixed() && this->doEMLFit()) {
1751  numbMap["sig"] = signalEvents_->genValue();
1752  }
1753  if ( usingBkgnd_ ) {
1754  const UInt_t nBkgnds = this->nBkgndClasses();
1755  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1756  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1757  const LauParameter* par = bkgndEvents_[iBkgnd];
1758  if (par->fixed()) {
1759  numbMap[bkgndClass] = par->genValue();
1760  }
1761  }
1762  }
1763  return numbMap;
1764 }
1765 
1767 {
1768  LauSPlot::TwoDMap twodimMap;
1769 
1770  for (LauPdfList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end(); ++pdf_iter) {
1771  // Count the number of input variables that are not DP variables
1772  UInt_t nVars(0);
1773  std::vector<TString> varNames = (*pdf_iter)->varNames();
1774  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1775  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1776  ++nVars;
1777  }
1778  }
1779  if ( nVars == 2 ) {
1780  if (useSCF_) {
1781  twodimMap.insert( std::make_pair( "sigTM", std::make_pair( varNames[0], varNames[1] ) ) );
1782  } else {
1783  twodimMap.insert( std::make_pair( "sig", std::make_pair( varNames[0], varNames[1] ) ) );
1784  }
1785  }
1786  }
1787 
1788  if ( useSCF_ ) {
1789  for (LauPdfList::const_iterator pdf_iter = scfPdfs_.begin(); pdf_iter != scfPdfs_.end(); ++pdf_iter) {
1790  // Count the number of input variables that are not DP variables
1791  UInt_t nVars(0);
1792  std::vector<TString> varNames = (*pdf_iter)->varNames();
1793  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1794  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1795  ++nVars;
1796  }
1797  }
1798  if ( nVars == 2 ) {
1799  twodimMap.insert( std::make_pair( "sigSCF", std::make_pair( varNames[0], varNames[1] ) ) );
1800  }
1801  }
1802  }
1803 
1804  if (usingBkgnd_) {
1805  const UInt_t nBkgnds = this->nBkgndClasses();
1806  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1807  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1808  const LauPdfList& pdfList = bkgndPdfs_[iBkgnd];
1809  for (LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1810  // Count the number of input variables that are not DP variables
1811  UInt_t nVars(0);
1812  std::vector<TString> varNames = (*pdf_iter)->varNames();
1813  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1814  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1815  ++nVars;
1816  }
1817  }
1818  if ( nVars == 2 ) {
1819  twodimMap.insert( std::make_pair( bkgndClass, std::make_pair( varNames[0], varNames[1] ) ) );
1820  }
1821  }
1822  }
1823  }
1824 
1825  return twodimMap;
1826 }
1827 
1829 {
1830  std::cout << "INFO in LauSimpleFitModel::storePerEvtLlhds : Storing per-event likelihood values..." << std::endl;
1831 
1832  // if we've not been using the DP model then we need to cache all
1833  // the info here so that we can get the efficiency from it
1834 
1835  LauFitDataTree* inputFitData = this->fitData();
1836 
1837  if (!this->useDP() && this->storeDPEff()) {
1839  sigDPModel_->fillDataTree(*inputFitData);
1840  }
1841 
1842  UInt_t evtsPerExpt(this->eventsPerExpt());
1843 
1844  for (UInt_t iEvt = 0; iEvt < evtsPerExpt; ++iEvt) {
1845 
1846  this->setSPlotNtupleIntegerBranchValue("iExpt",this->iExpt());
1847  this->setSPlotNtupleIntegerBranchValue("iEvtWithinExpt",iEvt);
1848 
1849  // the DP information
1850  this->getEvtDPLikelihood(iEvt);
1851  if (this->storeDPEff()) {
1852  if (!this->useDP()) {
1854  }
1855  this->setSPlotNtupleDoubleBranchValue("efficiency",sigDPModel_->getEvtEff());
1856  if ( sigDPModel_->usingScfModel() ) {
1858  }
1859  }
1860  if (this->useDP()) {
1862  if (useSCF_) {
1864  this->setSPlotNtupleDoubleBranchValue("sigTMDPLike",sigDPLike_);
1865  this->setSPlotNtupleDoubleBranchValue("sigSCFDPLike",scfDPLike_);
1866  } else {
1867  this->setSPlotNtupleDoubleBranchValue("sigDPLike",sigDPLike_);
1868  }
1869  if (usingBkgnd_) {
1870  const UInt_t nBkgnds = this->nBkgndClasses();
1871  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1872  TString name = this->bkgndClassName(iBkgnd);
1873  name += "DPLike";
1874  this->setSPlotNtupleDoubleBranchValue(name,bkgndDPLike_[iBkgnd]);
1875  }
1876  }
1877  } else {
1878  sigTotalLike_ = 1.0;
1879  if (useSCF_) {
1880  scfTotalLike_ = 1.0;
1881  }
1882  if (usingBkgnd_) {
1883  const UInt_t nBkgnds = this->nBkgndClasses();
1884  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1885  bkgndTotalLike_[iBkgnd] = 1.0;
1886  }
1887  }
1888  }
1889 
1890  // the signal PDF values
1891  if ( useSCF_ ) {
1892  sigTotalLike_ *= this->setSPlotNtupleBranchValues(&signalPdfs_, "sigTM", iEvt);
1893  scfTotalLike_ *= this->setSPlotNtupleBranchValues(&scfPdfs_, "sigSCF", iEvt);
1894  } else {
1895  sigTotalLike_ *= this->setSPlotNtupleBranchValues(&signalPdfs_, "sig", iEvt);
1896  }
1897 
1898  // the background PDF values
1899  if (usingBkgnd_) {
1900  const UInt_t nBkgnds = this->nBkgndClasses();
1901  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1902  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
1903  LauPdfList& pdfs = bkgndPdfs_[iBkgnd];
1904  bkgndTotalLike_[iBkgnd] *= this->setSPlotNtupleBranchValues(&(pdfs), bkgndClass, iEvt);
1905  }
1906  }
1907 
1908  // the total likelihoods
1909  if (useSCF_) {
1910  Double_t scfFrac(0.0);
1911  if ( useSCFHist_ ) {
1912  scfFrac = recoSCFFracs_[iEvt];
1913  } else {
1914  scfFrac = scfFrac_.value();
1915  }
1916  this->setSPlotNtupleDoubleBranchValue("sigSCFFrac",scfFrac);
1917  sigTotalLike_ *= ( 1.0 - scfFrac );
1918  if ( scfMap_ == 0 ) {
1919  scfTotalLike_ *= scfFrac;
1920  }
1921  this->setSPlotNtupleDoubleBranchValue("sigTMTotalLike",sigTotalLike_);
1922  this->setSPlotNtupleDoubleBranchValue("sigSCFTotalLike",scfTotalLike_);
1923  } else {
1924  this->setSPlotNtupleDoubleBranchValue("sigTotalLike",sigTotalLike_);
1925  }
1926  if (usingBkgnd_) {
1927  const UInt_t nBkgnds = this->nBkgndClasses();
1928  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1929  TString name = this->bkgndClassName(iBkgnd);
1930  name += "TotalLike";
1932  }
1933  }
1934  // fill the tree
1935  this->fillSPlotNtupleBranches();
1936  }
1937  std::cout << "INFO in LauSimpleFitModel::storePerEvtLlhds : Finished storing per-event likelihood values." << std::endl;
1938 }
1939 
1940 void LauSimpleFitModel::embedSignal(const TString& fileName, const TString& treeName,
1941  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
1942  Bool_t useReweighting)
1943 {
1944  if (signalTree_) {
1945  std::cerr << "ERROR in LauSimpleFitModel::embedSignal : Already embedding signal from a file." << std::endl;
1946  return;
1947  }
1948 
1949  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
1950  std::cerr << "WARNING in LauSimpleFitModel::embedSignal : Conflicting options provided, will not reuse events at all." << std::endl;
1951  reuseEventsWithinExperiment = kFALSE;
1952  }
1953 
1954  signalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
1955 
1956  Bool_t dataOK = signalTree_->findBranches();
1957  if (!dataOK) {
1958  delete signalTree_; signalTree_ = 0;
1959  std::cerr << "ERROR in LauSimpleFitModel::embedSignal : Problem creating data tree for embedding." << std::endl;
1960  return;
1961  }
1962 
1963  reuseSignal_ = reuseEventsWithinEnsemble;
1964  useReweighting_ = useReweighting;
1965 
1966  if (this->enableEmbedding() == kFALSE) {
1967  this->enableEmbedding(kTRUE);
1968  }
1969 }
1970 
1971 void LauSimpleFitModel::embedBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
1972  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
1973 {
1974  if ( ! this->validBkgndClass( bkgndClass ) ) {
1975  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
1976  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << std::endl;
1977  return;
1978  }
1979 
1980  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
1981 
1982  if (bkgndTree_[bkgndID]) {
1983  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Already embedding background from a file." << std::endl;
1984  return;
1985  }
1986 
1987  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
1988  std::cerr << "WARNING in LauSimpleFitModel::embedBkgnd : Conflicting options provided, will not reuse events at all." << std::endl;
1989  reuseEventsWithinExperiment = kFALSE;
1990  }
1991 
1992  bkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
1993 
1994  Bool_t dataOK = bkgndTree_[bkgndID]->findBranches();
1995  if (!dataOK) {
1996  delete bkgndTree_[bkgndID]; bkgndTree_[bkgndID] = 0;
1997  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Problem creating data tree for embedding." << std::endl;
1998  return;
1999  }
2000 
2001  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
2002 
2003  if (this->enableEmbedding() == kFALSE) {
2004  this->enableEmbedding(kTRUE);
2005  }
2006 }
2007 
2008 void LauSimpleFitModel::weightEvents( const TString& dataFileName, const TString& dataTreeName )
2009 {
2010  // Routine to provide weights for events that are uniformly distributed
2011  // in the DP (or square DP) so as to reproduce the given DP model
2012 
2013  if ( kinematics_->squareDP() ) {
2014  std::cout << "INFO in LauSimpleFitModel::weightEvents : will create weights assuming events were generated flat in the square DP" << std::endl;
2015  } else {
2016  std::cout << "INFO in LauSimpleFitModel::weightEvents : will create weights assuming events were generated flat in phase space" << std::endl;
2017  }
2018 
2019  // This reads in the given dataFile and creates an input
2020  // fit data tree that stores them for all events and experiments.
2021  Bool_t dataOK = this->cacheFitData(dataFileName,dataTreeName);
2022  if (!dataOK) {
2023  std::cerr << "ERROR in LauSimpleFitModel::weightEvents : Problem caching the data." << std::endl;
2024  return;
2025  }
2026 
2027  LauFitDataTree* inputFitData = this->fitData();
2028 
2029  if ( ! inputFitData->haveBranch( "m13Sq_MC" ) || ! inputFitData->haveBranch( "m23Sq_MC" ) ) {
2030  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Cannot find MC truth DP coordinate branches in supplied data, aborting." << std::endl;
2031  return;
2032  }
2033 
2034  // Create the ntuple to hold the DP weights
2035  TString weightsFileName( dataFileName );
2036  Ssiz_t index = weightsFileName.Last('.');
2037  weightsFileName.Insert( index, "_DPweights" );
2038  LauGenNtuple * weightsTuple = new LauGenNtuple( weightsFileName, dataTreeName );
2039  weightsTuple->addIntegerBranch("iExpt");
2040  weightsTuple->addIntegerBranch("iEvtWithinExpt");
2041  weightsTuple->addDoubleBranch("dpModelWeight");
2042 
2043  UInt_t iExpmt = this->iExpt();
2044  UInt_t nExpmt = this->nExpt();
2045  UInt_t firstExpmt = this->firstExpt();
2046  for (iExpmt = firstExpmt; iExpmt < (firstExpmt+nExpmt); ++iExpmt) {
2047 
2048  inputFitData->readExperimentData(iExpmt);
2049  UInt_t nEvents = inputFitData->nEvents();
2050 
2051  if (nEvents < 1) {
2052  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Zero events in experiment " << iExpmt << ", skipping..." << std::endl;
2053  continue;
2054  }
2055 
2056  weightsTuple->setIntegerBranchValue( "iExpt", iExpmt );
2057 
2058  // Calculate and store the weights for the events in this experiment
2059  for ( UInt_t iEvent(0); iEvent < nEvents; ++iEvent ) {
2060 
2061  weightsTuple->setIntegerBranchValue( "iEvtWithinExpt", iEvent );
2062 
2063  const LauFitData& evtData = inputFitData->getData( iEvent );
2064 
2065  Double_t m13Sq_MC = evtData.find("m13Sq_MC")->second;
2066  Double_t m23Sq_MC = evtData.find("m23Sq_MC")->second;
2067 
2068  Double_t dpModelWeight(0.0);
2069 
2070  if ( kinematics_->withinDPLimits( m13Sq_MC, m23Sq_MC ) ) {
2071 
2072  kinematics_->updateKinematics( m13Sq_MC, m23Sq_MC );
2073  dpModelWeight = sigDPModel_->getEventWeight();
2074 
2075  if ( kinematics_->squareDP() ) {
2076  dpModelWeight *= kinematics_->calcSqDPJacobian();
2077  }
2078 
2079  if (LauAbsDPDynamics::GenOK != sigDPModel_->checkToyMC(kTRUE,kFALSE)) {
2080  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Problem in calculating the weight, aborting." << std::endl;
2081  delete weightsTuple;
2082  return;
2083  }
2084  }
2085 
2086  weightsTuple->setDoubleBranchValue( "dpModelWeight", dpModelWeight );
2087  weightsTuple->fillBranches();
2088  }
2089 
2090  }
2091 
2092  weightsTuple->buildIndex( "iExpt", "iEvtWithinExpt" );
2093  weightsTuple->addFriendTree(dataFileName, dataTreeName);
2094  weightsTuple->writeOutGenResults();
2095 
2096  delete weightsTuple;
2097 }
2098 
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:55
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.
void splitSignalComponent(const TH2 *dpHisto, const Bool_t upperHalf=kFALSE, const Bool_t fluctuateBins=kFALSE, LauScfMap *scfMap=0)
Split the signal component into well-reconstructed and mis-reconstructed parts.
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()
Calculate things that depend on the fit parameters after they have been updated by Minuit...
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...
Pure abstract base class for defining the efficiency description across the signal Dalitz plot...
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 normalisation parameters (i.e. yields)
Double_t calcEfficiency(const LauKinematics *kinematics) const
Determine the efficiency for a given point in the Dalitz plot.
Definition: LauEffModel.cc:261
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)
Append fake data points to the inputData for each bin in the SCF smearing matrix. ...
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.
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.
LauAbsEffModel * getEffModel()
Retrieve the model for the efficiency across the Dalitz plot.
virtual void addSPlotNtupleDoubleBranch(const TString &name)
Add a branch to the sPlot tree for storing a double.
LauParArray fitFracEffUnCorr_
Fit fractions (uncorrected for the efficiency)
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.
const LauParArray & getFitFractionsEfficiencyUncorrected() const
Retrieve the fit fractions for the amplitude components.
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
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.
virtual Bool_t passVeto(const LauKinematics *kinematics) const =0
Determine whether the given DP position is outside the vetoes.
Bool_t generateSignalEvent()
Generate signal event.