laura is hosted by Hepforge, IPPP Durham
Laura++  v1r2
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauCPFitModel.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 #include <vector>
19 using std::cout;
20 using std::cerr;
21 using std::endl;
22 using std::setprecision;
23 using std::vector;
24 
25 #include "TVirtualFitter.h"
26 #include "TSystem.h"
27 #include "TMinuit.h"
28 #include "TRandom.h"
29 #include "TFile.h"
30 #include "TMath.h"
31 #include "TH2.h"
32 
33 #include "LauAbsBkgndDPModel.hh"
34 #include "LauAbsCoeffSet.hh"
35 #include "LauAbsDPDynamics.hh"
36 #include "LauAbsPdf.hh"
37 #include "LauAsymmCalc.hh"
38 #include "LauComplex.hh"
39 #include "LauConstants.hh"
40 #include "LauCPFitModel.hh"
41 #include "LauDaughters.hh"
42 #include "LauEffModel.hh"
43 #include "LauEmbeddedData.hh"
44 #include "LauFitNtuple.hh"
45 #include "LauKinematics.hh"
46 #include "LauPrint.hh"
47 #include "LauRandom.hh"
48 #include "LauScfMap.hh"
49 
50 ClassImp(LauCPFitModel)
51 
52 
53 LauCPFitModel::LauCPFitModel(LauAbsDPDynamics* negModel, LauAbsDPDynamics* posModel, Bool_t tagged, const TString& tagVarName) : LauAbsFitModel(),
54  negSigModel_(negModel), posSigModel_(posModel),
55  negKinematics_(negModel ? negModel->getKinematics() : 0),
56  posKinematics_(posModel ? posModel->getKinematics() : 0),
57  usingBkgnd_(kFALSE),
58  nSigComp_(0),
59  nSigDPPar_(0),
60  nExtraPdfPar_(0),
61  nNormPar_(0),
62  negMeanEff_("negMeanEff",0.0,0.0,1.0), posMeanEff_("posMeanEff",0.0,0.0,1.0),
63  negDPRate_("negDPRate",0.0,0.0,100.0), posDPRate_("posDPRate",0.0,0.0,100.0),
64  signalEvents_(0),
65  signalAsym_(0),
66  forceAsym_(kFALSE),
67  tagged_(tagged),
68  tagVarName_(tagVarName),
69  curEvtCharge_(0),
70  useSCF_(kFALSE),
71  useSCFHist_(kFALSE),
72  scfFrac_("scfFrac",0.0,0.0,1.0),
73  scfFracHist_(0),
74  scfMap_(0),
75  compareFitData_(kFALSE),
76  negParent_("B-"), posParent_("B+"),
77  negSignalTree_(0), posSignalTree_(0),
78  reuseSignal_(kFALSE),
79  useNegReweighting_(kFALSE), usePosReweighting_(kFALSE),
80  sigDPLike_(0.0),
81  scfDPLike_(0.0),
82  sigExtraLike_(0.0),
83  scfExtraLike_(0.0),
84  sigTotalLike_(0.0),
85  scfTotalLike_(0.0)
86 {
87  LauDaughters* negDaug = negSigModel_->getDaughters();
88  if (negDaug != 0) {negParent_ = negDaug->getNameParent();}
89  LauDaughters* posDaug = posSigModel_->getDaughters();
90  if (posDaug != 0) {posParent_ = posDaug->getNameParent();}
91 }
92 
94 {
95  delete negSignalTree_;
96  delete posSignalTree_;
97  for (LauBkgndEmbDataList::iterator iter = negBkgndTree_.begin(); iter != negBkgndTree_.end(); ++iter) {
98  delete (*iter);
99  }
100  for (LauBkgndEmbDataList::iterator iter = posBkgndTree_.begin(); iter != posBkgndTree_.end(); ++iter) {
101  delete (*iter);
102  }
103  delete scfFracHist_;
104 }
105 
107 {
108  UInt_t nBkgnds = this->nBkgndClasses();
109  negBkgndDPModels_.resize( nBkgnds );
110  posBkgndDPModels_.resize( nBkgnds );
111  negBkgndPdfs_.resize( nBkgnds );
112  posBkgndPdfs_.resize( nBkgnds );
113  bkgndEvents_.resize( nBkgnds );
114  bkgndAsym_.resize( nBkgnds );
115  negBkgndTree_.resize( nBkgnds );
116  posBkgndTree_.resize( nBkgnds );
117  reuseBkgnd_.resize( nBkgnds );
118  bkgndDPLike_.resize( nBkgnds );
119  bkgndExtraLike_.resize( nBkgnds );
120  bkgndTotalLike_.resize( nBkgnds );
121 }
122 
124 {
125  if ( nSigEvents == 0 ) {
126  cerr << "ERROR in LauCPFitModel::setNSigEvents : The LauParameter pointer is null." << endl;
127  gSystem->Exit(EXIT_FAILURE);
128  }
129  if ( signalEvents_ != 0 ) {
130  cerr << "ERROR in LauCPFitModel::setNSigEvents : You are trying to overwrite the signal yield." << endl;
131  return;
132  }
133  if ( signalAsym_ != 0 ) {
134  cerr << "ERROR in LauCPFitModel::setNSigEvents : You are trying to overwrite the signal asymmetry." << endl;
135  return;
136  }
137 
138  signalEvents_ = nSigEvents;
139  signalEvents_->name("signalEvents");
140  Double_t value = nSigEvents->value();
141  signalEvents_->range(-2.0*(TMath::Abs(value)+1.0), 2.0*(TMath::Abs(value)+1.0));
142 
143  signalAsym_ = new LauParameter("signalAsym",0.0,-1.0,1.0,kTRUE);
144 }
145 
146 void LauCPFitModel::setNSigEvents( LauParameter* nSigEvents, LauParameter* sigAsym, Bool_t forceAsym )
147 {
148  if ( nSigEvents == 0 ) {
149  cerr << "ERROR in LauCPFitModel::setNSigEvents : The event LauParameter pointer is null." << endl;
150  gSystem->Exit(EXIT_FAILURE);
151  }
152  if ( sigAsym == 0 ) {
153  cerr << "ERROR in LauCPFitModel::setNSigEvents : The asym LauParameter pointer is null." << endl;
154  gSystem->Exit(EXIT_FAILURE);
155  }
156 
157  if ( signalEvents_ != 0 ) {
158  cerr << "ERROR in LauCPFitModel::setNSigEvents : You are trying to overwrite the signal yield." << endl;
159  return;
160  }
161  if ( signalAsym_ != 0 ) {
162  cerr << "ERROR in LauCPFitModel::setNSigEvents : You are trying to overwrite the signal asymmetry." << endl;
163  return;
164  }
165 
166  signalEvents_ = nSigEvents;
167  signalEvents_->name("signalEvents");
168  Double_t value = nSigEvents->value();
169  signalEvents_->range(-2.0*(TMath::Abs(value)+1.0), 2.0*(TMath::Abs(value)+1.0));
170 
171  signalAsym_ = sigAsym;
172  signalAsym_->name("signalAsym");
173  signalAsym_->range(-1.0,1.0);
174 
175  forceAsym_ = forceAsym;
176 }
177 
179 {
180  if ( nBkgndEvents == 0 ) {
181  cerr << "ERROR in LauCPFitModel::setNBgkndEvents : The background yield LauParameter pointer is null." << endl;
182  gSystem->Exit(EXIT_FAILURE);
183  }
184 
185  if ( ! this->validBkgndClass( nBkgndEvents->name() ) ) {
186  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : Invalid background class \"" << nBkgndEvents->name() << "\"." << std::endl;
187  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
188  gSystem->Exit(EXIT_FAILURE);
189  }
190 
191  UInt_t bkgndID = this->bkgndClassID( nBkgndEvents->name() );
192 
193  if ( bkgndEvents_[bkgndID] != 0 ) {
194  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : You are trying to overwrite the background yield." << endl;
195  return;
196  }
197 
198  if ( bkgndAsym_[bkgndID] != 0 ) {
199  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : You are trying to overwrite the background asymmetry." << endl;
200  return;
201  }
202 
203  bkgndEvents_[bkgndID] = nBkgndEvents;
204  bkgndEvents_[bkgndID]->name( nBkgndEvents->name()+"Events" );
205  Double_t value = nBkgndEvents->value();
206  bkgndEvents_[bkgndID]->range(-2.0*(TMath::Abs(value)+1.0), 2.0*(TMath::Abs(value)+1.0));
207 
208  bkgndAsym_[bkgndID] = new LauParameter(nBkgndEvents->name()+"Asym",0.0,-1.0,1.0,kTRUE);
209 }
210 
212 {
213  if ( nBkgndEvents == 0 ) {
214  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : The background yield LauParameter pointer is null." << endl;
215  gSystem->Exit(EXIT_FAILURE);
216  }
217 
218  if ( bkgndAsym == 0 ) {
219  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : The background asym LauParameter pointer is null." << endl;
220  gSystem->Exit(EXIT_FAILURE);
221  }
222 
223  if ( ! this->validBkgndClass( nBkgndEvents->name() ) ) {
224  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : Invalid background class \"" << nBkgndEvents->name() << "\"." << std::endl;
225  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
226  gSystem->Exit(EXIT_FAILURE);
227  }
228 
229  UInt_t bkgndID = this->bkgndClassID( nBkgndEvents->name() );
230 
231  if ( bkgndEvents_[bkgndID] != 0 ) {
232  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : You are trying to overwrite the background yield." << endl;
233  return;
234  }
235 
236  if ( bkgndAsym_[bkgndID] != 0 ) {
237  cerr << "ERROR in LauCPFitModel::setNBkgndEvents : You are trying to overwrite the background asymmetry." << endl;
238  return;
239  }
240 
241  bkgndEvents_[bkgndID] = nBkgndEvents;
242  bkgndEvents_[bkgndID]->name( nBkgndEvents->name()+"Events" );
243  Double_t value = nBkgndEvents->value();
244  bkgndEvents_[bkgndID]->range(-2.0*(TMath::Abs(value)+1.0), 2.0*(TMath::Abs(value)+1.0));
245 
246  bkgndAsym_[bkgndID] = bkgndAsym;
247  bkgndAsym_[bkgndID]->name( nBkgndEvents->name()+"Asym" );
248  bkgndAsym_[bkgndID]->range(-1.0,1.0);
249 }
250 
251 void LauCPFitModel::splitSignalComponent( const TH2* dpHisto, Bool_t upperHalf, LauScfMap* scfMap )
252 {
253  if ( useSCF_ == kTRUE ) {
254  cerr << "ERROR in LauCPFitModel::splitSignalComponent : Have already setup SCF." << endl;
255  return;
256  }
257 
258  if ( dpHisto == 0 ) {
259  cerr << "ERROR in LauCPFitModel::splitSignalComponent : The histogram pointer is null." << endl;
260  return;
261  }
262 
263  LauDaughters* daughters = negSigModel_->getDaughters();
264  scfFracHist_ = new LauEffModel( daughters, 0 );
265  scfFracHist_->setEffHisto( dpHisto, kTRUE, kFALSE, 0.0, 0.0, upperHalf, daughters->squareDP() );
266 
267  scfMap_ = scfMap;
268 
269  useSCF_ = kTRUE;
270  useSCFHist_ = kTRUE;
271 }
272 
273 void LauCPFitModel::splitSignalComponent( Double_t scfFrac, Bool_t fixed )
274 {
275  if ( useSCF_ == kTRUE ) {
276  cerr << "ERROR in LauCPFitModel::splitSignalComponent : Have already setup SCF." << endl;
277  return;
278  }
279 
280  scfFrac_.range( 0.0, 1.0 );
281  scfFrac_.value( scfFrac ); scfFrac_.initValue( scfFrac ); scfFrac_.genValue( scfFrac );
282  scfFrac_.fixed( fixed );
283 
284  useSCF_ = kTRUE;
285  useSCFHist_ = kFALSE;
286 }
287 
288 void LauCPFitModel::setBkgndDPModels(const TString& bkgndClass, LauAbsBkgndDPModel* negModel, LauAbsBkgndDPModel* posModel)
289 {
290  if ((negModel==0) || (posModel==0)) {
291  cerr << "ERROR in LauCPFitModel::setBkgndDPModels : One or both of the model pointers is null." << endl;
292  return;
293  }
294 
295  // check that this background name is valid
296  if ( ! this->validBkgndClass( bkgndClass) ) {
297  cerr << "ERROR in LauCPFitModel::setBkgndDPModel : Invalid background class \"" << bkgndClass << "\"." << std::endl;
298  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
299  return;
300  }
301 
302  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
303  negBkgndDPModels_[bkgndID] = negModel;
304  posBkgndDPModels_[bkgndID] = posModel;
305 
306  usingBkgnd_ = kTRUE;
307 }
308 
310 {
311  if ( tagged_ ) {
312  if (negPdf==0 || posPdf==0) {
313  cerr << "ERROR in LauCPFitModel::setSignalPdfs : One or both of the PDF pointers is null." << endl;
314  return;
315  }
316  } else {
317  // if we're doing an untagged analysis we will only use the negative PDFs
318  if ( negPdf==0 ) {
319  cerr << "ERROR in LauCPFitModel::setSignalPdfs : The negative PDF pointer is null." << endl;
320  return;
321  }
322  if ( posPdf!=0 ) {
323  cerr << "WARNING in LauCPFitModel::setSignalPdfs : Doing an untagged fit so will not use the positive PDF." << endl;
324  }
325  }
326  negSignalPdfs_.push_back(negPdf);
327  posSignalPdfs_.push_back(posPdf);
328 }
329 
331 {
332  if ( tagged_ ) {
333  if (negPdf==0 || posPdf==0) {
334  cerr << "ERROR in LauCPFitModel::setSCFPdfs : One or both of the PDF pointers is null." << endl;
335  return;
336  }
337  } else {
338  // if we're doing an untagged analysis we will only use the negative PDFs
339  if ( negPdf==0 ) {
340  cerr << "ERROR in LauCPFitModel::setSCFPdfs : The negative PDF pointer is null." << endl;
341  return;
342  }
343  if ( posPdf!=0 ) {
344  cerr << "WARNING in LauCPFitModel::setSCFPdfs : Doing an untagged fit so will not use the positive PDF." << endl;
345  }
346  }
347  negScfPdfs_.push_back(negPdf);
348  posScfPdfs_.push_back(posPdf);
349 }
350 
351 void LauCPFitModel::setBkgndPdfs(const TString& bkgndClass, LauAbsPdf* negPdf, LauAbsPdf* posPdf)
352 {
353  if ( tagged_ ) {
354  if (negPdf==0 || posPdf==0) {
355  cerr << "ERROR in LauCPFitModel::setBkgndPdfs : One or both of the PDF pointers is null." << endl;
356  return;
357  }
358  } else {
359  // if we're doing an untagged analysis we will only use the negative PDFs
360  if ( negPdf==0 ) {
361  cerr << "ERROR in LauCPFitModel::setBkgndPdfs : The negative PDF pointer is null." << endl;
362  return;
363  }
364  if ( posPdf!=0 ) {
365  cerr << "WARNING in LauCPFitModel::setBkgndPdfs : Doing an untagged fit so will not use the positive PDF." << endl;
366  }
367  }
368 
369  // check that this background name is valid
370  if ( ! this->validBkgndClass( bkgndClass ) ) {
371  cerr << "ERROR in LauCPFitModel::setBkgndPdfs : Invalid background class \"" << bkgndClass << "\"." << std::endl;
372  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
373  return;
374  }
375 
376  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
377  negBkgndPdfs_[bkgndID].push_back(negPdf);
378  posBkgndPdfs_[bkgndID].push_back(posPdf);
379 
380  usingBkgnd_ = kTRUE;
381 }
382 
384 {
385  // Is there a component called compName in the signal model?
386  TString compName(coeffSet->name());
387  Bool_t negOK = negSigModel_->hasResonance(compName);
388  TString conjName = negSigModel_->getConjResName(compName);
389  Bool_t posOK = posSigModel_->hasResonance(conjName);
390  if (!negOK) {
391  cout << "ERROR in LauCPFitModel::setMagPhase : " << negParent_ << " signal DP model doesn't contain component \"" << compName << "\"." << endl;
392  return;
393  }
394  if (!posOK) {
395  cout << "ERROR in LauCPFitModel::setMagPhase : " << posParent_ << " signal DP model doesn't contain component \"" << conjName << "\"." << endl;
396  return;
397  }
398 
399  // Do we already have it in our list of names?
400  for (std::vector<LauAbsCoeffSet*>::const_iterator iter=coeffPars_.begin(); iter!=coeffPars_.end(); ++iter) {
401  if ((*iter)->name() == compName) {
402  cerr << "ERROR in LauCPFitModel::setAmpCoeffSet : Have already set coefficients for \"" << compName << "\"." << endl;
403  return;
404  }
405  }
406 
407  coeffSet->index(nSigComp_);
408  coeffPars_.push_back(coeffSet);
409 
410  TString parName = coeffSet->baseName(); parName += "FitFracAsym";
411  fitFracAsymm_.push_back(LauParameter(parName, 0.0, -1.0, 1.0));
412 
413  acp_.push_back(coeffSet->acp());
414 
415  ++nSigComp_;
416 
417  cout << "Added coefficients for component \"" << compName << "\" to the fit model." << endl;
418 }
419 
421 {
422  // First of all check that, if a given component is being used,
423  // we've got the PDFs for all the variables involved
424  if ( this->useDP() ) {
425  if ((negSigModel_ == 0) || (posSigModel_ == 0)) {
426  cerr << "ERROR in LauCPFitModel::initialise : the pointer to one (neg or pos) of the signal DP models is null.\n";
427  cerr << " : Removing the Dalitz Plot from the model." << endl;
428  this->useDP(kFALSE);
429  }
430  if ( usingBkgnd_ ) {
431  if ( negBkgndDPModels_.empty() || posBkgndDPModels_.empty() ) {
432  cerr << "ERROR in LauCPFitModel::initialise : No background DP models found.\n";
433  cerr << " : Removing the Dalitz plot from the model." << endl;
434  this->useDP(kFALSE);
435  }
436  for (LauBkgndDPModelList::const_iterator dpmodel_iter = negBkgndDPModels_.begin(); dpmodel_iter != negBkgndDPModels_.end(); ++dpmodel_iter ) {
437  if ( (*dpmodel_iter) == 0 ) {
438  cerr << "ERROR in LauCPFitModel::initialise : The pointer to one of the background DP models is null.\n";
439  cerr << " : Removing the Dalitz Plot from the model." << endl;
440  this->useDP(kFALSE);
441  break;
442  }
443  }
444  for (LauBkgndDPModelList::const_iterator dpmodel_iter = posBkgndDPModels_.begin(); dpmodel_iter != posBkgndDPModels_.end(); ++dpmodel_iter ) {
445  if ( (*dpmodel_iter) == 0 ) {
446  cerr << "ERROR in LauCPFitModel::initialise : The pointer to one of the background DP models is null.\n";
447  cerr << " : Removing the Dalitz Plot from the model." << endl;
448  this->useDP(kFALSE);
449  break;
450  }
451  }
452  }
453  }
454 
455  // Next check that, if a given component is being used we've got the
456  // right number of PDFs for all the variables involved
457  // TODO - should probably check variable names and so on as well
458 
459  UInt_t nsigpdfvars(0);
460  for ( LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter ) {
461  std::vector<TString> varNames = (*pdf_iter)->varNames();
462  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
463  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
464  ++nsigpdfvars;
465  }
466  }
467  }
468  if (useSCF_) {
469  UInt_t nscfpdfvars(0);
470  for ( LauPdfList::const_iterator pdf_iter = negScfPdfs_.begin(); pdf_iter != negScfPdfs_.end(); ++pdf_iter ) {
471  std::vector<TString> varNames = (*pdf_iter)->varNames();
472  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
473  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
474  ++nscfpdfvars;
475  }
476  }
477  }
478  if (nscfpdfvars != nsigpdfvars) {
479  cerr << "ERROR in LauCPFitModel::initialise : There are " << nsigpdfvars << " TM signal PDF variables but " << nscfpdfvars << " SCF signal PDF variables." << endl;
480  gSystem->Exit(EXIT_FAILURE);
481  }
482  }
483  if (usingBkgnd_) {
484  for (LauBkgndPdfsList::const_iterator bgclass_iter = negBkgndPdfs_.begin(); bgclass_iter != negBkgndPdfs_.end(); ++bgclass_iter) {
485  UInt_t nbkgndpdfvars(0);
486  const LauPdfList& pdfList = (*bgclass_iter);
487  for ( LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter ) {
488  std::vector<TString> varNames = (*pdf_iter)->varNames();
489  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
490  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
491  ++nbkgndpdfvars;
492  }
493  }
494  }
495  if (nbkgndpdfvars != nsigpdfvars) {
496  cerr << "ERROR in LauCPFitModel::initialise : There are " << nsigpdfvars << " signal PDF variables but " << nbkgndpdfvars << " bkgnd PDF variables." << endl;
497  gSystem->Exit(EXIT_FAILURE);
498  }
499  }
500  }
501 
502  // Clear the vectors of parameter information so we can start from scratch
503  this->clearFitParVectors();
504 
505  // Set the fit parameters for signal and background models
506  this->setSignalDPParameters();
507 
508  // Set the fit parameters for the various extra PDFs
509  this->setExtraPdfParameters();
510 
511  // Set the initial bg and signal events
512  this->setFitNEvents();
513 
514  // Check that we have the expected number of fit variables
515  const LauParameterPList& fitVars = this->fitPars();
516  if (fitVars.size() != (nSigDPPar_ + nExtraPdfPar_ + nNormPar_)) {
517  cerr << "ERROR in LauCPFitModel::initialise : Number of fit parameters not of expected size. Exiting" << endl;
518  gSystem->Exit(EXIT_FAILURE);
519  }
520 
521  // From the initial parameter values calculate the coefficients
522  // so they can be passed to the signal model
523  this->updateCoeffs();
524 
525  // Initialisation
526  if (this->useDP() == kTRUE) {
527  this->initialiseDPModels();
528  }
529 
530  if (!this->useDP() && negSignalPdfs_.empty()) {
531  cerr << "ERROR in LauCPFitModel::initialise : Signal model doesn't exist for any variable." << endl;
532  gSystem->Exit(EXIT_FAILURE);
533  }
534 
535  this->setExtraNtupleVars();
536 }
537 
539 {
540  cout << "INFO in LauCPFitModel::initialiseDPModels : Initialising signal DP model" << endl;
543 
544  if (usingBkgnd_ == kTRUE) {
545  for (LauBkgndDPModelList::iterator iter = negBkgndDPModels_.begin(); iter != negBkgndDPModels_.end(); ++iter) {
546  (*iter)->initialise();
547  }
548  for (LauBkgndDPModelList::iterator iter = posBkgndDPModels_.begin(); iter != posBkgndDPModels_.end(); ++iter) {
549  (*iter)->initialise();
550  }
551  }
552 }
553 
555 {
556  // Set the fit parameters for the signal model.
557 
558  nSigDPPar_ = 0;
559  if ( ! this->useDP() ) {
560  return;
561  }
562 
563  cout << "INFO in LauCPFitModel::setSignalDPParameters : Setting the initial fit parameters for the signal DP model." << endl;
564 
565  // Need to check that the number of components we have and that the dynamics has matches up
566  UInt_t nNegAmp = negSigModel_->getnAmp();
567  UInt_t nPosAmp = posSigModel_->getnAmp();
568  if ( nNegAmp != nPosAmp ) {
569  cerr << "ERROR in LauCPFitModel::setSignalDPParameters : Unequal number of signal DP components in the negative and positive models: " << nNegAmp << " != " << nPosAmp << endl;
570  gSystem->Exit(EXIT_FAILURE);
571  }
572  if ( nNegAmp != nSigComp_ ) {
573  cerr << "ERROR in LauCPFitModel::setSignalDPParameters : Number of signal DP components in the model (" << nNegAmp << ") not equal to number of coefficients supplied (" << nSigComp_ << ")." << endl;
574  gSystem->Exit(EXIT_FAILURE);
575  }
576 
577 
578  // Place signal model parameters in vector of fit variables
579  LauParameterPList& fitVars = this->fitPars();
580  for (UInt_t i = 0; i < nSigComp_; i++) {
581  LauParameterPList pars = coeffPars_[i]->getParameters();
582  for (LauParameterPList::iterator iter = pars.begin(); iter != pars.end(); ++iter) {
583  if ( !(*iter)->clone() ) {
584  fitVars.push_back(*iter);
585  ++nSigDPPar_;
586  }
587  }
588  }
589 }
590 
592 {
593  // Include all the parameters of the PDF in the fit
594  // NB all of them are passed to the fit, even though some have been fixed through parameter.fixed(kTRUE)
595  // With the new "cloned parameter" scheme only "original" parameters are passed to the fit.
596  // Their clones are updated automatically when the originals are updated.
597 
598  nExtraPdfPar_ = 0;
599 
601  if ( tagged_ ) {
603  }
604 
605  if (useSCF_ == kTRUE) {
607  if ( tagged_ ) {
609  }
610  }
611 
612  if (usingBkgnd_ == kTRUE) {
613  for (LauBkgndPdfsList::iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
614  nExtraPdfPar_ += this->addFitParameters(*iter);
615  }
616  if ( tagged_ ) {
617  for (LauBkgndPdfsList::iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
618  nExtraPdfPar_ += this->addFitParameters(*iter);
619  }
620  }
621  }
622 }
623 
625 {
626  if ( signalEvents_ == 0 ) {
627  cerr << "ERROR in LauCPFitModel::setFitNEvents : Signal yield not defined." << endl;
628  return;
629  }
630  nNormPar_ = 0;
631 
632  // initialise the total number of events to be the sum of all the hypotheses
633  Double_t nTotEvts = signalEvents_->value();
634  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
635  nTotEvts += (*iter)->value();
636  if ( (*iter) == 0 ) {
637  cerr << "ERROR in LauCPFitModel::setFitNEvents : Background yield not defined." << endl;
638  return;
639  }
640  }
641  this->eventsPerExpt(TMath::FloorNint(nTotEvts));
642 
643  LauParameterPList& fitVars = this->fitPars();
644 
645  // if doing an extended ML fit add the number of signal events into the fit parameters
646  if (this->doEMLFit()) {
647  cout << "INFO in LauCPFitModel::setFitNEvents : Initialising number of events for signal and background components..." << endl;
648  // add the signal fraction to the list of fit parameters
649  fitVars.push_back(signalEvents_);
650  ++nNormPar_;
651  } else {
652  cout << "INFO in LauCPFitModel::setFitNEvents : Initialising number of events for background components (and hence signal)..." << endl;
653  }
654 
655  // if not using the DP in the model we need an explicit signal asymmetry parameter
656  if (this->useDP() == kFALSE) {
657  fitVars.push_back(signalAsym_);
658  ++nNormPar_;
659  }
660 
661  if (useSCF_ && !useSCFHist_) {
662  fitVars.push_back(&scfFrac_);
663  ++nNormPar_;
664  }
665 
666  if (usingBkgnd_ == kTRUE) {
667  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
668  LauParameter* parameter = (*iter);
669  fitVars.push_back(parameter);
670  ++nNormPar_;
671  }
672  for (LauBkgndYieldList::iterator iter = bkgndAsym_.begin(); iter != bkgndAsym_.end(); ++iter) {
673  LauParameter* parameter = (*iter);
674  fitVars.push_back(parameter);
675  ++nNormPar_;
676  }
677  }
678 }
679 
681 {
682  // Set-up other parameters derived from the fit results, e.g. fit fractions.
683 
684  if (this->useDP() != kTRUE) {
685  return;
686  }
687 
688  // First clear the vectors so we start from scratch
689  this->clearExtraVarVectors();
690 
691  // Add the positive and negative fit fractions for each signal component
693  if (negFitFrac_.size() != nSigComp_) {
694  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << negFitFrac_.size() << endl;
695  gSystem->Exit(EXIT_FAILURE);
696  }
697  for (UInt_t i(0); i<nSigComp_; ++i) {
698  if (negFitFrac_[i].size() != nSigComp_) {
699  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << negFitFrac_[i].size() << endl;
700  gSystem->Exit(EXIT_FAILURE);
701  }
702  }
703 
704  LauParameterList& extraVars = this->extraPars();
705 
706  for (UInt_t i(0); i<nSigComp_; ++i) {
707  for (UInt_t j(i); j<nSigComp_; ++j) {
708  TString name = negFitFrac_[i][j].name();
709  name.Insert( name.Index("FitFrac"), "Neg" );
710  negFitFrac_[i][j].name(name);
711  extraVars.push_back(negFitFrac_[i][j]);
712  }
713  }
714 
716  if (posFitFrac_.size() != nSigComp_) {
717  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << posFitFrac_.size() << endl;
718  gSystem->Exit(EXIT_FAILURE);
719  }
720  for (UInt_t i(0); i<nSigComp_; ++i) {
721  if (posFitFrac_[i].size() != nSigComp_) {
722  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << posFitFrac_[i].size() << endl;
723  gSystem->Exit(EXIT_FAILURE);
724  }
725  }
726  for (UInt_t i(0); i<nSigComp_; ++i) {
727  for (UInt_t j(i); j<nSigComp_; ++j) {
728  TString name = posFitFrac_[i][j].name();
729  name.Insert( name.Index("FitFrac"), "Pos" );
730  posFitFrac_[i][j].name(name);
731  extraVars.push_back(posFitFrac_[i][j]);
732  }
733  }
734 
735  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
736  std::vector<LauParameter> negExtraPars = negSigModel_->getExtraParameters();
737  std::vector<LauParameter>::iterator negExtraIter;
738  for (negExtraIter = negExtraPars.begin(); negExtraIter != negExtraPars.end(); ++negExtraIter) {
739  LauParameter negExtraParameter = (*negExtraIter);
740  extraVars.push_back(negExtraParameter);
741  }
742 
743  std::vector<LauParameter> posExtraPars = posSigModel_->getExtraParameters();
744  std::vector<LauParameter>::iterator posExtraIter;
745  for (posExtraIter = posExtraPars.begin(); posExtraIter != posExtraPars.end(); ++posExtraIter) {
746  LauParameter posExtraParameter = (*posExtraIter);
747  extraVars.push_back(posExtraParameter);
748  }
749 
750  // Now add in the DP efficiency value
751  Double_t initMeanEff = negSigModel_->getMeanEff().initValue();
752  negMeanEff_.value(initMeanEff);
753  negMeanEff_.genValue(initMeanEff);
754  negMeanEff_.initValue(initMeanEff);
755  extraVars.push_back(negMeanEff_);
756 
757  initMeanEff = posSigModel_->getMeanEff().initValue();
758  posMeanEff_.value(initMeanEff);
759  posMeanEff_.genValue(initMeanEff);
760  posMeanEff_.initValue(initMeanEff);
761  extraVars.push_back(posMeanEff_);
762 
763  // Also add in the DP rates
764  Double_t initDPRate = negSigModel_->getDPRate().initValue();
765  negDPRate_.value(initDPRate);
766  negDPRate_.genValue(initDPRate);
767  negDPRate_.initValue(initDPRate);
768  extraVars.push_back(negDPRate_);
769 
770  initDPRate = posSigModel_->getDPRate().initValue();
771  posDPRate_.value(initDPRate);
772  posDPRate_.genValue(initDPRate);
773  posDPRate_.initValue(initDPRate);
774  extraVars.push_back(posDPRate_);
775 
776  // Calculate the CPC and CPV Fit Fractions, ACPs and FitFrac asymmetries
777  this->calcExtraFractions(kTRUE);
778  this->calcAsymmetries(kTRUE);
779 
780  // Add the CP violating and CP conserving fit fractions for each signal component
781  for (UInt_t i = 0; i < nSigComp_; i++) {
782  for (UInt_t j = i; j < nSigComp_; j++) {
783  extraVars.push_back(CPVFitFrac_[i][j]);
784  }
785  }
786  for (UInt_t i = 0; i < nSigComp_; i++) {
787  for (UInt_t j = i; j < nSigComp_; j++) {
788  extraVars.push_back(CPCFitFrac_[i][j]);
789  }
790  }
791 
792  // Add the Fit Fraction asymmetry for each signal component
793  for (UInt_t i = 0; i < nSigComp_; i++) {
794  extraVars.push_back(fitFracAsymm_[i]);
795  }
796 
797  // Add the calculated CP asymmetry for each signal component
798  for (UInt_t i = 0; i < nSigComp_; i++) {
799  extraVars.push_back(acp_[i]);
800  }
801 }
802 
803 void LauCPFitModel::calcExtraFractions(Bool_t initValues)
804 {
805  // Calculate the CP-conserving and CP-violating fit fractions
806 
807  if (initValues) {
808  // create the structure
809  CPCFitFrac_.clear();
810  CPVFitFrac_.clear();
811  CPCFitFrac_.resize(nSigComp_);
812  CPVFitFrac_.resize(nSigComp_);
813  for (UInt_t i(0); i<nSigComp_; ++i) {
814  CPCFitFrac_[i].resize(nSigComp_);
815  CPVFitFrac_[i].resize(nSigComp_);
816 
817  for (UInt_t j(i); j<nSigComp_; ++j) {
818  TString name = negFitFrac_[i][j].name();
819  name.Replace( name.Index("Neg"), 3, "CPC" );
820  CPCFitFrac_[i][j].name( name );
821  CPCFitFrac_[i][j].valueAndRange( 0.0, -100.0, 100.0 );
822 
823  name = negFitFrac_[i][j].name();
824  name.Replace( name.Index("Neg"), 3, "CPV" );
825  CPVFitFrac_[i][j].name( name );
826  CPVFitFrac_[i][j].valueAndRange( 0.0, -100.0, 100.0 );
827  }
828  }
829  }
830 
831  Double_t denom = negDPRate_ + posDPRate_;
832 
833  for (UInt_t i(0); i<nSigComp_; ++i) {
834  for (UInt_t j(i); j<nSigComp_; ++j) {
835 
836  Double_t negTerm = negFitFrac_[i][j]*negDPRate_;
837  Double_t posTerm = posFitFrac_[i][j]*posDPRate_;
838 
839  Double_t cpcFitFrac = (negTerm + posTerm)/denom;
840  Double_t cpvFitFrac = (negTerm - posTerm)/denom;
841 
842  CPCFitFrac_[i][j] = cpcFitFrac;
843  CPVFitFrac_[i][j] = cpvFitFrac;
844 
845  if (initValues) {
846  CPCFitFrac_[i][j].genValue(cpcFitFrac);
847  CPCFitFrac_[i][j].initValue(cpcFitFrac);
848 
849  CPVFitFrac_[i][j].genValue(cpvFitFrac);
850  CPVFitFrac_[i][j].initValue(cpvFitFrac);
851  }
852  }
853  }
854 }
855 
856 void LauCPFitModel::calcAsymmetries(Bool_t initValues)
857 {
858  // Calculate the CP asymmetries
859  // Also calculate the fit fraction asymmetries
860 
861  for (UInt_t i = 0; i < nSigComp_; i++) {
862 
863  acp_[i] = coeffPars_[i]->acp();
864 
865  LauAsymmCalc asymmCalc(negFitFrac_[i][i].value(), posFitFrac_[i][i].value());
866  Double_t asym = asymmCalc.getAsymmetry();
867  fitFracAsymm_[i] = asym;
868  if (initValues) {
869  fitFracAsymm_[i].genValue(asym);
870  fitFracAsymm_[i].initValue(asym);
871  }
872  }
873 }
874 
875 void LauCPFitModel::finaliseFitResults(const TString& tablePrefixName)
876 {
877  // Retrieve parameters from the fit results for calculations and toy generation
878  // and eventually store these in output root ntuples/text files
879 
880  // Now take the fit parameters and update them as necessary
881  // i.e. to make mag > 0.0, phase in the right range.
882  // This function will also calculate any other values, such as the
883  // fit fractions, using any errors provided by fitParErrors as appropriate.
884  // Also obtain the pull values: (measured - generated)/(average error)
885 
886  if (this->useDP() == kTRUE) {
887  for (UInt_t i = 0; i < nSigComp_; ++i) {
888  // Check whether we have "a/b > 0.0", and phases in the right range
889  coeffPars_[i]->finaliseValues();
890  }
891  }
892 
893  // update the pulls on the event fractions and asymmetries
894  if (this->doEMLFit()) {
896  }
897  if (this->useDP() == kFALSE) {
899  }
900  if (useSCF_ && !useSCFHist_) {
902  }
903  if (usingBkgnd_ == kTRUE) {
904  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
905  (*iter)->updatePull();
906  }
907  for (LauBkgndYieldList::iterator iter = bkgndAsym_.begin(); iter != bkgndAsym_.end(); ++iter) {
908  (*iter)->updatePull();
909  }
910  }
911 
912  // Update the pulls on all the extra PDFs' parameters
915  if (useSCF_ == kTRUE) {
918  }
919  if (usingBkgnd_ == kTRUE) {
920  for (LauBkgndPdfsList::iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
921  this->updateFitParameters(*iter);
922  }
923  for (LauBkgndPdfsList::iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
924  this->updateFitParameters(*iter);
925  }
926  }
927 
928  // Fill the fit results to the ntuple
929 
930  // update the coefficients and then calculate the fit fractions and ACP's
931  if (this->useDP() == kTRUE) {
932  this->updateCoeffs();
935 
936  LauParArray negFitFrac = negSigModel_->getFitFractions();
937  if (negFitFrac.size() != nSigComp_) {
938  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << negFitFrac.size() << endl;
939  gSystem->Exit(EXIT_FAILURE);
940  }
941  for (UInt_t i(0); i<nSigComp_; ++i) {
942  if (negFitFrac[i].size() != nSigComp_) {
943  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << negFitFrac[i].size() << endl;
944  gSystem->Exit(EXIT_FAILURE);
945  }
946  }
947 
948  LauParArray posFitFrac = posSigModel_->getFitFractions();
949  if (posFitFrac.size() != nSigComp_) {
950  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << posFitFrac.size() << endl;
951  gSystem->Exit(EXIT_FAILURE);
952  }
953  for (UInt_t i(0); i<nSigComp_; ++i) {
954  if (posFitFrac[i].size() != nSigComp_) {
955  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << posFitFrac[i].size() << endl;
956  gSystem->Exit(EXIT_FAILURE);
957  }
958  }
959 
960  for (UInt_t i(0); i<nSigComp_; ++i) {
961  for (UInt_t j(i); j<nSigComp_; ++j) {
962  negFitFrac_[i][j].value(negFitFrac[i][j].value());
963  posFitFrac_[i][j].value(posFitFrac[i][j].value());
964  }
965  }
966 
971 
972  this->calcExtraFractions();
973  this->calcAsymmetries();
974 
975  // Then store the final fit parameters, and any extra parameters for
976  // the signal model (e.g. fit fractions, FF asymmetries, ACPs, mean efficiency and DP rate)
977 
978  this->clearExtraVarVectors();
979  LauParameterList& extraVars = this->extraPars();
980 
981  // Add the positive and negative fit fractions for each signal component
982  for (UInt_t i(0); i<nSigComp_; ++i) {
983  for (UInt_t j(i); j<nSigComp_; ++j) {
984  extraVars.push_back(negFitFrac_[i][j]);
985  }
986  }
987  for (UInt_t i(0); i<nSigComp_; ++i) {
988  for (UInt_t j(i); j<nSigComp_; ++j) {
989  extraVars.push_back(posFitFrac_[i][j]);
990  }
991  }
992 
993  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
994  std::vector<LauParameter> negExtraPars = negSigModel_->getExtraParameters();
995  std::vector<LauParameter>::iterator negExtraIter;
996  for (negExtraIter = negExtraPars.begin(); negExtraIter != negExtraPars.end(); ++negExtraIter) {
997  LauParameter negExtraParameter = (*negExtraIter);
998  extraVars.push_back(negExtraParameter);
999  }
1000 
1001  std::vector<LauParameter> posExtraPars = posSigModel_->getExtraParameters();
1002  std::vector<LauParameter>::iterator posExtraIter;
1003  for (posExtraIter = posExtraPars.begin(); posExtraIter != posExtraPars.end(); ++posExtraIter) {
1004  LauParameter posExtraParameter = (*posExtraIter);
1005  extraVars.push_back(posExtraParameter);
1006  }
1007 
1008  extraVars.push_back(negMeanEff_);
1009  extraVars.push_back(posMeanEff_);
1010  extraVars.push_back(negDPRate_);
1011  extraVars.push_back(posDPRate_);
1012 
1013  for (UInt_t i = 0; i < nSigComp_; i++) {
1014  for (UInt_t j(i); j<nSigComp_; ++j) {
1015  extraVars.push_back(CPVFitFrac_[i][j]);
1016  }
1017  }
1018  for (UInt_t i = 0; i < nSigComp_; i++) {
1019  for (UInt_t j(i); j<nSigComp_; ++j) {
1020  extraVars.push_back(CPCFitFrac_[i][j]);
1021  }
1022  }
1023  for (UInt_t i(0); i<nSigComp_; ++i) {
1024  extraVars.push_back(fitFracAsymm_[i]);
1025  }
1026  for (UInt_t i(0); i<nSigComp_; ++i) {
1027  extraVars.push_back(acp_[i]);
1028  }
1029 
1030  this->printFitFractions(cout);
1031  this->printAsymmetries(cout);
1032  }
1033 
1034  const LauParameterPList& fitVars = this->fitPars();
1035  const LauParameterList& extraVars = this->extraPars();
1036  LauFitNtuple* ntuple = this->fitNtuple();
1037  ntuple->storeParsAndErrors(fitVars, extraVars);
1038 
1039  // find out the correlation matrix for the parameters
1040  ntuple->storeCorrMatrix(this->iExpt(), this->nll(), this->fitStatus());
1041 
1042  // Fill the data into ntuple
1043  ntuple->updateFitNtuple();
1044 
1045  // Print out the partial fit fractions, phases and the
1046  // averaged efficiency, reweighted by the dynamics (and anything else)
1047  if (this->writeLatexTable()) {
1048  TString sigOutFileName(tablePrefixName);
1049  sigOutFileName += "_"; sigOutFileName += this->iExpt(); sigOutFileName += "Expt.tex";
1050  this->writeOutTable(sigOutFileName);
1051  }
1052 }
1053 
1054 void LauCPFitModel::printFitFractions(std::ostream& output)
1055 {
1056  // Print out Fit Fractions, total DP rate and mean efficiency
1057  // First for the B- events
1058  for (UInt_t i = 0; i < nSigComp_; i++) {
1059  const TString compName(coeffPars_[i]->name());
1060  output << negParent_ << " FitFraction for component " << i << " (" << compName << ") = " << negFitFrac_[i][i] << endl;
1061  }
1062  output << negParent_ << " overall DP rate (integral of matrix element squared) = " << negDPRate_ << endl;
1063  output << negParent_ << " average efficiency weighted by whole DP dynamics = " << negMeanEff_ << endl;
1064 
1065  // Then for the positive sample
1066  for (UInt_t i = 0; i < nSigComp_; i++) {
1067  const TString compName(coeffPars_[i]->name());
1068  const TString conjName(negSigModel_->getConjResName(compName));
1069  output << posParent_ << " FitFraction for component " << i << " (" << conjName << ") = " << posFitFrac_[i][i] << endl;
1070  }
1071  output << posParent_ << " overall DP rate (integral of matrix element squared) = " << posDPRate_ << endl;
1072  output << posParent_ << " average efficiency weighted by whole DP dynamics = " << posMeanEff_ << endl;
1073 }
1074 
1075 void LauCPFitModel::printAsymmetries(std::ostream& output)
1076 {
1077  for (UInt_t i = 0; i < nSigComp_; i++) {
1078  const TString compName(coeffPars_[i]->name());
1079  output << "Fit Fraction asymmetry for component " << i << " (" << compName << ") = " << fitFracAsymm_[i] << endl;
1080  }
1081  for (UInt_t i = 0; i < nSigComp_; i++) {
1082  const TString compName(coeffPars_[i]->name());
1083  output << "ACP for component " << i << " (" << compName << ") = " << acp_[i] << endl;
1084  }
1085 }
1086 
1087 void LauCPFitModel::writeOutTable(const TString& outputFile)
1088 {
1089  // Write out the results of the fit to a tex-readable table
1090  // TODO - need to include the yields in this table
1091  std::ofstream fout(outputFile);
1092  LauPrint print;
1093 
1094  cout << "Writing out results of the fit to the tex file " << outputFile << endl;
1095 
1096  if (this->useDP() == kTRUE) {
1097  // print the fit coefficients in one table
1098  coeffPars_.front()->printTableHeading(fout);
1099  for (UInt_t i = 0; i < nSigComp_; i++) {
1100  coeffPars_[i]->printTableRow(fout);
1101  }
1102  fout << "\\hline" << endl;
1103  fout << "\\end{tabular}" << endl << endl;
1104 
1105  // print the fit fractions and asymmetries in another
1106  fout << "\\begin{tabular}{|l|c|c|c|c|}" << endl;
1107  fout << "\\hline" << endl;
1108  fout << "Component & " << negParent_ << " Fit Fraction & " << posParent_ << " Fit Fraction & Fit Fraction Asymmetry & ACP \\\\" << endl;
1109  fout << "\\hline" << endl;
1110 
1111  Double_t negFitFracSum(0.0);
1112  Double_t posFitFracSum(0.0);
1113 
1114  for (UInt_t i = 0; i < nSigComp_; i++) {
1115  TString resName = coeffPars_[i]->name();
1116  resName = resName.ReplaceAll("_", "\\_");
1117 
1118  Double_t negFitFrac = negFitFrac_[i][i].value();
1119  Double_t posFitFrac = posFitFrac_[i][i].value();
1120  negFitFracSum += negFitFrac;
1121  posFitFracSum += posFitFrac;
1122 
1123  Double_t fitFracAsymm = fitFracAsymm_[i].value();
1124 
1125  Double_t acp = acp_[i].value();
1126  Double_t acpErr = acp_[i].error();
1127 
1128  fout << resName << " & $";
1129 
1130  print.printFormat(fout, negFitFrac);
1131  fout << "$ & $";
1132  print.printFormat(fout, posFitFrac);
1133  fout << "$ & $";
1134 
1135  print.printFormat(fout, fitFracAsymm);
1136  fout << "$ & $";
1137 
1138  print.printFormat(fout, acp);
1139  fout << " \\pm ";
1140  print.printFormat(fout, acpErr);
1141  fout << "$ \\\\" << endl;
1142  }
1143  fout << "\\hline" << endl;
1144 
1145  // Also print out sum of fit fractions
1146  fout << "Fit Fraction Sum & $";
1147  print.printFormat(fout, negFitFracSum);
1148  fout << "$ & $";
1149  print.printFormat(fout, posFitFracSum);
1150  fout << "$ & & \\\\" << endl;
1151  fout << "\\hline" << endl;
1152 
1153  fout << "DP rate & $";
1154  print.printFormat(fout, negDPRate_.value());
1155  fout << "$ & $";
1156  print.printFormat(fout, posDPRate_.value());
1157  fout << "$ & & \\\\" << endl;
1158 
1159  fout << "$< \\varepsilon > $ & $";
1160  print.printFormat(fout, negMeanEff_.value());
1161  fout << "$ & $";
1162  print.printFormat(fout, posMeanEff_.value());
1163  fout << "$ & & \\\\" << endl;
1164  fout << "\\hline" << endl;
1165  fout << "\\end{tabular}" << endl << endl;
1166  }
1167 
1168  if (!negSignalPdfs_.empty()) {
1169  fout << "\\begin{tabular}{|l|c|}" << endl;
1170  fout << "\\hline" << endl;
1171  if (useSCF_ == kTRUE) {
1172  fout << "\\Extra TM Signal PDFs' Parameters: & \\\\" << endl;
1173  } else {
1174  fout << "\\Extra Signal PDFs' Parameters: & \\\\" << endl;
1175  }
1176  this->printFitParameters(negSignalPdfs_, fout);
1177  if ( tagged_ ) {
1178  this->printFitParameters(posSignalPdfs_, fout);
1179  }
1180  if (useSCF_ == kTRUE && !negScfPdfs_.empty()) {
1181  fout << "\\hline" << endl;
1182  fout << "\\Extra SCF Signal PDFs' Parameters: & \\\\" << endl;
1183  this->printFitParameters(negScfPdfs_, fout);
1184  if ( tagged_ ) {
1185  this->printFitParameters(posScfPdfs_, fout);
1186  }
1187  }
1188  if (usingBkgnd_ == kTRUE && !negBkgndPdfs_.empty()) {
1189  fout << "\\hline" << endl;
1190  fout << "\\Extra Background PDFs' Parameters: & \\\\" << endl;
1191  for (LauBkgndPdfsList::const_iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
1192  this->printFitParameters(*iter, fout);
1193  }
1194  if ( tagged_ ) {
1195  for (LauBkgndPdfsList::const_iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
1196  this->printFitParameters(*iter, fout);
1197  }
1198  }
1199  }
1200  fout << "\\hline \n\\end{tabular}" << endl << endl;
1201  }
1202 }
1203 
1205 {
1206  // Update the number of signal events to be total-sum(background events)
1207  this->updateSigEvents();
1208 
1209  // Check whether we want to have randomised initial fit parameters for the signal model
1210  if (this->useRandomInitFitPars() == kTRUE) {
1211  cout << "Setting random parameters for the signal model" << endl;
1212  this->randomiseInitFitPars();
1213  }
1214 }
1215 
1217 {
1218  // Only randomise those parameters that are not fixed!
1219  cout << "Randomising the initial fit magnitudes and phases of the components..." << endl;
1220 
1221  for (UInt_t i = 0; i < nSigComp_; i++) {
1222  coeffPars_[i]->randomiseInitValues();
1223  }
1224 }
1225 
1227 {
1228  // Determine the number of events to generate for each hypothesis
1229  // If we're smearing then smear each one individually
1230 
1231  LauGenInfo nEvtsGen;
1232 
1233  // Signal
1234  Double_t evtWeight(1.0);
1235  Double_t nEvts = signalEvents_->genValue();
1236  if ( nEvts < 0.0 ) {
1237  evtWeight = -1.0;
1238  nEvts = TMath::Abs( nEvts );
1239  }
1240  Double_t asym(0.0);
1241  Double_t sigAsym(0.0);
1242  // need to include this as an alternative in case the DP isn't in the model
1243  if ( !this->useDP() || forceAsym_ ) {
1244  sigAsym = signalAsym_->genValue();
1245  } else {
1246  Double_t negRate = negSigModel_->getDPNorm();
1247  Double_t posRate = posSigModel_->getDPNorm();
1248  if (negRate+posRate>1e-30) {
1249  sigAsym = (negRate-posRate)/(negRate+posRate);
1250  }
1251  }
1252  asym = sigAsym;
1253  Int_t nPosEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 - asym)) + 0.5);
1254  Int_t nNegEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 + asym)) + 0.5);
1255  if (this->doPoissonSmearing()) {
1256  nNegEvts = LauRandom::randomFun()->Poisson(nNegEvts);
1257  nPosEvts = LauRandom::randomFun()->Poisson(nPosEvts);
1258  }
1259  nEvtsGen[std::make_pair("signal",-1)] = std::make_pair(nNegEvts,evtWeight);
1260  nEvtsGen[std::make_pair("signal",+1)] = std::make_pair(nPosEvts,evtWeight);
1261 
1262  // backgrounds
1263  const UInt_t nBkgnds = this->nBkgndClasses();
1264  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1265  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1266  const LauParameter* evtsPar = bkgndEvents_[bkgndID];
1267  const LauParameter* asymPar = bkgndAsym_[bkgndID];
1268  evtWeight = 1.0;
1269  nEvts = TMath::FloorNint( evtsPar->genValue() );
1270  if ( nEvts < 0 ) {
1271  evtWeight = -1.0;
1272  nEvts = TMath::Abs( nEvts );
1273  }
1274  asym = asymPar->genValue();
1275  nPosEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 - asym)) + 0.5);
1276  nNegEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 + asym)) + 0.5);
1277  if (this->doPoissonSmearing()) {
1278  nNegEvts = LauRandom::randomFun()->Poisson(nNegEvts);
1279  nPosEvts = LauRandom::randomFun()->Poisson(nPosEvts);
1280  }
1281  nEvtsGen[std::make_pair(bkgndClass,-1)] = std::make_pair(nNegEvts,evtWeight);
1282  nEvtsGen[std::make_pair(bkgndClass,+1)] = std::make_pair(nPosEvts,evtWeight);
1283  }
1284 
1285  cout << "Generating toy MC with:" << endl;
1286  cout << "Signal asymmetry = " << sigAsym << " and number of signal events = " << signalEvents_->genValue() << endl;
1287  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1288  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1289  const LauParameter* evtsPar = bkgndEvents_[bkgndID];
1290  const LauParameter* asymPar = bkgndAsym_[bkgndID];
1291  cout << bkgndClass << " asymmetry = " << asymPar->genValue() << " and number of " << bkgndClass << " events = " << evtsPar->genValue() << endl;
1292  }
1293 
1294  return nEvtsGen;
1295 }
1296 
1298 {
1299  // Routine to generate toy Monte Carlo events according to the various models we have defined.
1300 
1301  // Determine the number of events to generate for each hypothesis
1302  LauGenInfo nEvts = this->eventsToGenerate();
1303 
1304  Bool_t genOK(kTRUE);
1305  Int_t evtNum(0);
1306 
1307  const UInt_t nBkgnds = this->nBkgndClasses();
1308  std::vector<TString> bkgndClassNames(nBkgnds);
1309  std::vector<TString> bkgndClassNamesGen(nBkgnds);
1310  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1311  TString name( this->bkgndClassName(iBkgnd) );
1312  bkgndClassNames[iBkgnd] = name;
1313  bkgndClassNamesGen[iBkgnd] = "gen"+name;
1314  }
1315 
1316  const Bool_t storeSCFTruthInfo = ( useSCF_ || ( this->enableEmbedding() &&
1317  negSignalTree_ != 0 && negSignalTree_->haveBranch("mcMatch") &&
1318  posSignalTree_ != 0 && posSignalTree_->haveBranch("mcMatch") ) );
1319 
1320  // Loop over the hypotheses and generate the requested number of events for each one
1321  for (LauGenInfo::const_iterator iter = nEvts.begin(); iter != nEvts.end(); ++iter) {
1322 
1323  const TString& type(iter->first.first);
1324  curEvtCharge_ = iter->first.second;
1325  Double_t evtWeight( iter->second.second );
1326  Int_t nEvtsGen( iter->second.first );
1327 
1328  for (Int_t iEvt(0); iEvt<nEvtsGen; ++iEvt) {
1329 
1330  this->setGenNtupleDoubleBranchValue( "evtWeight", evtWeight );
1331 
1332  if (type == "signal") {
1333  this->setGenNtupleIntegerBranchValue("genSig",1);
1334  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1335  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], 0 );
1336  }
1337  genOK = this->generateSignalEvent();
1338  } else {
1339  this->setGenNtupleIntegerBranchValue("genSig",0);
1340  if ( storeSCFTruthInfo ) {
1341  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1342  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1343  }
1344  UInt_t bkgndID(0);
1345  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1346  Int_t gen(0);
1347  if ( bkgndClassNames[iBkgnd] == type ) {
1348  gen = 1;
1349  bkgndID = iBkgnd;
1350  }
1351  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], gen );
1352  }
1353  genOK = this->generateBkgndEvent(bkgndID);
1354  }
1355  if (!genOK) {
1356  // If there was a problem with the generation then break out and return.
1357  // The problem model will have adjusted itself so that all should be OK next time.
1358  break;
1359  }
1360  if (this->useDP() == kTRUE) {
1361  this->setDPBranchValues();
1362  }
1363 
1364  // Store the event charge
1366 
1367  // Store the event number (within this experiment)
1368  // and then increment it
1369  this->setGenNtupleIntegerBranchValue("iEvtWithinExpt",evtNum);
1370  ++evtNum;
1371 
1372  this->fillGenNtupleBranches();
1373  if (iEvt%500 == 0) {cout << "Generated event number " << iEvt << " out of " << nEvtsGen << " " << type << " events." << endl;}
1374  }
1375 
1376  if (!genOK) {
1377  break;
1378  }
1379  }
1380 
1381  if (this->useDP() && genOK) {
1382 
1383  negSigModel_->checkToyMC(kTRUE,kTRUE);
1384  posSigModel_->checkToyMC(kTRUE,kTRUE);
1385 
1386  // Get the fit fractions if they're to be written into the latex table
1387  if (this->writeLatexTable()) {
1388  LauParArray negFitFrac = negSigModel_->getFitFractions();
1389  if (negFitFrac.size() != nSigComp_) {
1390  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << negFitFrac.size() << endl;
1391  gSystem->Exit(EXIT_FAILURE);
1392  }
1393  for (UInt_t i(0); i<nSigComp_; ++i) {
1394  if (negFitFrac[i].size() != nSigComp_) {
1395  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << negFitFrac[i].size() << endl;
1396  gSystem->Exit(EXIT_FAILURE);
1397  }
1398  }
1399  LauParArray posFitFrac = posSigModel_->getFitFractions();
1400  if (posFitFrac.size() != nSigComp_) {
1401  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << posFitFrac.size() << endl;
1402  gSystem->Exit(EXIT_FAILURE);
1403  }
1404  for (UInt_t i(0); i<nSigComp_; ++i) {
1405  if (posFitFrac[i].size() != nSigComp_) {
1406  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << posFitFrac[i].size() << endl;
1407  gSystem->Exit(EXIT_FAILURE);
1408  }
1409  }
1410 
1411  for (UInt_t i(0); i<nSigComp_; ++i) {
1412  for (UInt_t j(i); j<nSigComp_; ++j) {
1413  negFitFrac_[i][j].value(negFitFrac[i][j].value());
1414  posFitFrac_[i][j].value(posFitFrac[i][j].value());
1415  }
1416  }
1421  }
1422  }
1423 
1424  // If we're reusing embedded events or if the generation is being
1425  // reset then clear the lists of used events
1426  if (reuseSignal_ || !genOK) {
1427  if (negSignalTree_) {
1429  }
1430  if (posSignalTree_) {
1432  }
1433  }
1434  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1435  LauEmbeddedData* data = negBkgndTree_[bkgndID];
1436  if (reuseBkgnd_[bkgndID] || !genOK) {
1437  if (data) {
1438  data->clearUsedList();
1439  }
1440  }
1441  }
1442  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1443  LauEmbeddedData* data = posBkgndTree_[bkgndID];
1444  if (reuseBkgnd_[bkgndID] || !genOK) {
1445  if (data) {
1446  data->clearUsedList();
1447  }
1448  }
1449  }
1450 
1451  return genOK;
1452 }
1453 
1455 {
1456  // Generate signal event
1457  Bool_t genOK(kTRUE);
1458  Bool_t genSCF(kFALSE);
1459 
1460  LauAbsDPDynamics* model(0);
1461  LauKinematics* kinematics(0);
1462  LauEmbeddedData* embeddedData(0);
1463  LauPdfList* sigPdfs(0);
1464  LauPdfList* scfPdfs(0);
1465 
1466  Bool_t doReweighting(kFALSE);
1467 
1468  if (curEvtCharge_<0) {
1469  model = negSigModel_;
1470  kinematics = negKinematics_;
1471  sigPdfs = &negSignalPdfs_;
1472  scfPdfs = &negScfPdfs_;
1473  if (this->enableEmbedding()) {
1474  embeddedData = negSignalTree_;
1475  doReweighting = useNegReweighting_;
1476  }
1477  } else {
1478  model = posSigModel_;
1479  kinematics = posKinematics_;
1480  if ( tagged_ ) {
1481  sigPdfs = &posSignalPdfs_;
1482  scfPdfs = &posScfPdfs_;
1483  } else {
1484  sigPdfs = &negSignalPdfs_;
1485  scfPdfs = &negScfPdfs_;
1486  }
1487  if (this->enableEmbedding()) {
1488  embeddedData = posSignalTree_;
1489  doReweighting = usePosReweighting_;
1490  }
1491  }
1492 
1493  if (this->useDP()) {
1494  if (embeddedData) {
1495  if (doReweighting) {
1496  // Select a (random) event from the generated data. Then store the
1497  // reconstructed DP co-ords, together with other pdf information,
1498  // as the embedded data.
1499  genOK = embeddedData->getReweightedEvent(model);
1500  } else {
1501  // Just get the information of a (randomly) selected event in the
1502  // embedded data
1503  embeddedData->getEmbeddedEvent(kinematics);
1504  }
1505  genSCF = this->storeSignalMCMatch( embeddedData );
1506  } else {
1507  genOK = model->generate();
1508  if ( genOK && useSCF_ ) {
1509  Double_t frac(0.0);
1510  if ( useSCFHist_ ) {
1511  frac = scfFracHist_->calcEfficiency( kinematics );
1512  } else {
1513  frac = scfFrac_.genValue();
1514  }
1515  if ( frac < LauRandom::randomFun()->Rndm() ) {
1516  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1517  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1518  genSCF = kFALSE;
1519  } else {
1520  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1521  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1522  genSCF = kTRUE;
1523 
1524  // Optionally smear the DP position
1525  // of the SCF event
1526  if ( scfMap_ != 0 ) {
1527  Double_t xCoord(0.0), yCoord(0.0);
1528  if ( kinematics->squareDP() ) {
1529  xCoord = kinematics->getmPrime();
1530  yCoord = kinematics->getThetaPrime();
1531  } else {
1532  xCoord = kinematics->getm13Sq();
1533  yCoord = kinematics->getm23Sq();
1534  }
1535 
1536  // Find the bin number where this event is generated
1537  Int_t binNo = scfMap_->binNumber( xCoord, yCoord );
1538 
1539  // Retrieve the migration histogram
1540  TH2* histo = scfMap_->trueHist( binNo );
1541 
1542  LauEffModel * effModel = model->getEffModel();
1543  do {
1544  // Get a random point from the histogram
1545  histo->GetRandom2( xCoord, yCoord );
1546 
1547  // Update the kinematics
1548  if ( kinematics->squareDP() ) {
1549  kinematics->updateSqDPKinematics( xCoord, yCoord );
1550  } else {
1551  kinematics->updateKinematics( xCoord, yCoord );
1552  }
1553  } while ( ! effModel->passVeto( kinematics ) );
1554  }
1555  }
1556  }
1557  }
1558  } else {
1559  if (embeddedData) {
1560  embeddedData->getEmbeddedEvent(0);
1561  genSCF = this->storeSignalMCMatch( embeddedData );
1562  } else if ( useSCF_ ) {
1563  Double_t frac = scfFrac_.genValue();
1564  if ( frac < LauRandom::randomFun()->Rndm() ) {
1565  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1566  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1567  genSCF = kFALSE;
1568  } else {
1569  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1570  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1571  genSCF = kTRUE;
1572  }
1573  }
1574  }
1575  if (genOK) {
1576  if ( useSCF_ ) {
1577  if ( genSCF ) {
1578  this->generateExtraPdfValues(scfPdfs, embeddedData);
1579  } else {
1580  this->generateExtraPdfValues(sigPdfs, embeddedData);
1581  }
1582  } else {
1583  this->generateExtraPdfValues(sigPdfs, embeddedData);
1584  }
1585  }
1586  // Check for problems with the embedding
1587  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1588  cerr << "WARNING in LauCPFitModel::generateSignalEvent : Source of embedded signal events used up, clearing the list of used events." << endl;
1589  embeddedData->clearUsedList();
1590  }
1591 
1592  return genOK;
1593 }
1594 
1596 {
1597  // Generate Bkgnd event
1598  Bool_t genOK(kTRUE);
1599 
1600  LauAbsBkgndDPModel* model(0);
1601  LauEmbeddedData* embeddedData(0);
1602  LauPdfList* extraPdfs(0);
1603  LauKinematics* kinematics(0);
1604 
1605  if (curEvtCharge_<0) {
1606  model = negBkgndDPModels_[bkgndID];
1607  if (this->enableEmbedding()) {
1608  embeddedData = negBkgndTree_[bkgndID];
1609  }
1610  extraPdfs = &negBkgndPdfs_[bkgndID];
1611  kinematics = negKinematics_;
1612  } else {
1613  model = posBkgndDPModels_[bkgndID];
1614  if (this->enableEmbedding()) {
1615  embeddedData = posBkgndTree_[bkgndID];
1616  }
1617  if ( tagged_ ) {
1618  extraPdfs = &posBkgndPdfs_[bkgndID];
1619  } else {
1620  extraPdfs = &negBkgndPdfs_[bkgndID];
1621  }
1622  kinematics = posKinematics_;
1623  }
1624 
1625  if (this->useDP()) {
1626  if (embeddedData) {
1627  embeddedData->getEmbeddedEvent(kinematics);
1628  } else {
1629  if (model == 0) {
1630  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1631  cerr << "ERROR in LauCPFitModel::generateBkgndEvent : Can't find the DP model for background class \"" << bkgndClass << "\"." << endl;
1632  gSystem->Exit(EXIT_FAILURE);
1633  }
1634  genOK = model->generate();
1635  }
1636  } else {
1637  if (embeddedData) {
1638  embeddedData->getEmbeddedEvent(0);
1639  }
1640  }
1641  if (genOK) {
1642  this->generateExtraPdfValues(extraPdfs, embeddedData);
1643  }
1644 
1645  // Check for problems with the embedding
1646  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1647  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1648  cerr << "WARNING in LauCPFitModel::generateBkgndEvent : Source of embedded " << bkgndClass << " events used up, clearing the list of used events." << endl;
1649  embeddedData->clearUsedList();
1650  }
1651 
1652  return genOK;
1653 }
1654 
1656 {
1657  // Setup the required ntuple branches
1658  this->addGenNtupleDoubleBranch("evtWeight");
1659  this->addGenNtupleIntegerBranch("genSig");
1660  if ( useSCF_ || ( this->enableEmbedding() &&
1661  negSignalTree_ != 0 && negSignalTree_->haveBranch("mcMatch") &&
1662  posSignalTree_ != 0 && posSignalTree_->haveBranch("mcMatch") ) ) {
1663  this->addGenNtupleIntegerBranch("genTMSig");
1664  this->addGenNtupleIntegerBranch("genSCFSig");
1665  }
1666  const UInt_t nBkgnds = this->nBkgndClasses();
1667  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1668  TString name( this->bkgndClassName(iBkgnd) );
1669  name.Prepend("gen");
1670  this->addGenNtupleIntegerBranch(name);
1671  }
1672  this->addGenNtupleIntegerBranch("charge");
1673  if (this->useDP() == kTRUE) {
1674  this->addGenNtupleDoubleBranch("m12");
1675  this->addGenNtupleDoubleBranch("m23");
1676  this->addGenNtupleDoubleBranch("m13");
1677  this->addGenNtupleDoubleBranch("m12Sq");
1678  this->addGenNtupleDoubleBranch("m23Sq");
1679  this->addGenNtupleDoubleBranch("m13Sq");
1680  this->addGenNtupleDoubleBranch("cosHel12");
1681  this->addGenNtupleDoubleBranch("cosHel23");
1682  this->addGenNtupleDoubleBranch("cosHel13");
1684  this->addGenNtupleDoubleBranch("mPrime");
1685  this->addGenNtupleDoubleBranch("thPrime");
1686  }
1687  }
1688  for (LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter) {
1689  std::vector<TString> varNames = (*pdf_iter)->varNames();
1690  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1691  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1692  this->addGenNtupleDoubleBranch( (*var_iter) );
1693  }
1694  }
1695  }
1696 }
1697 
1699 {
1700  LauKinematics* kinematics(0);
1701  if (curEvtCharge_<0) {
1702  kinematics = negKinematics_;
1703  } else {
1704  kinematics = posKinematics_;
1705  }
1706 
1707  // Store all the DP information
1708  this->setGenNtupleDoubleBranchValue("m12", kinematics->getm12());
1709  this->setGenNtupleDoubleBranchValue("m23", kinematics->getm23());
1710  this->setGenNtupleDoubleBranchValue("m13", kinematics->getm13());
1711  this->setGenNtupleDoubleBranchValue("m12Sq", kinematics->getm12Sq());
1712  this->setGenNtupleDoubleBranchValue("m23Sq", kinematics->getm23Sq());
1713  this->setGenNtupleDoubleBranchValue("m13Sq", kinematics->getm13Sq());
1714  this->setGenNtupleDoubleBranchValue("cosHel12", kinematics->getc12());
1715  this->setGenNtupleDoubleBranchValue("cosHel23", kinematics->getc23());
1716  this->setGenNtupleDoubleBranchValue("cosHel13", kinematics->getc13());
1717  if (kinematics->squareDP()) {
1718  this->setGenNtupleDoubleBranchValue("mPrime", kinematics->getmPrime());
1719  this->setGenNtupleDoubleBranchValue("thPrime", kinematics->getThetaPrime());
1720  }
1721 }
1722 
1724 {
1725  LauKinematics* kinematics(0);
1726  if (curEvtCharge_<0) {
1727  kinematics = negKinematics_;
1728  } else {
1729  kinematics = posKinematics_;
1730  }
1731 
1732  if (!extraPdfs) {
1733  std::cerr << "ERROR in LauCPFitModel::generateExtraPdfValues : Null pointer to PDF list." << std::endl;
1734  gSystem->Exit(EXIT_FAILURE);
1735  }
1736 
1737  if (extraPdfs->empty()) {
1738  //std::cerr << "WARNING in LauCPFitModel::generateExtraPdfValues : PDF list is empty." << std::endl;
1739  return;
1740  }
1741 
1742  // Generate from the extra PDFs
1743  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1744  LauFitData genValues;
1745  if (embeddedData) {
1746  genValues = embeddedData->getValues( (*pdf_iter)->varNames() );
1747  } else {
1748  genValues = (*pdf_iter)->generate(kinematics);
1749  }
1750  for ( LauFitData::const_iterator var_iter = genValues.begin(); var_iter != genValues.end(); ++var_iter ) {
1751  TString varName = var_iter->first;
1752  if ( varName != "m13Sq" && varName != "m23Sq" ) {
1753  Double_t value = var_iter->second;
1754  this->setGenNtupleDoubleBranchValue(varName,value);
1755  }
1756  }
1757  }
1758 }
1759 
1761 {
1762  // Default to TM
1763  Bool_t genSCF(kFALSE);
1764  Int_t match(1);
1765 
1766  // Check that we have a valid pointer and that embedded data has
1767  // the mcMatch branch. If so then get the match value.
1768  if ( embeddedData && embeddedData->haveBranch("mcMatch") ) {
1769  match = TMath::Nint( embeddedData->getValue("mcMatch") );
1770  }
1771 
1772  // Set the variables accordingly.
1773  if (match) {
1774  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1775  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1776  genSCF = kFALSE;
1777  } else {
1778  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1779  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1780  genSCF = kTRUE;
1781  }
1782 
1783  return genSCF;
1784 }
1785 
1787 {
1788  // Update the signal parameters and then the total normalisation for the signal likelihood
1789  if (this->useDP() == kTRUE) {
1790  this->updateCoeffs();
1793  }
1794 
1795  // Update the signal fraction from the background fractions if not doing an extended fit
1796  if ( !this->doEMLFit() && !signalEvents_->fixed() ) {
1797  this->updateSigEvents();
1798  }
1799 }
1800 
1802 {
1803  // The background parameters will have been set from Minuit.
1804  // We need to update the signal events using these.
1805  Double_t nTotEvts = this->eventsPerExpt();
1806 
1807  signalEvents_->range(-2.0*nTotEvts,2.0*nTotEvts);
1808  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1809  (*iter)->range(-2.0*nTotEvts,2.0*nTotEvts);
1810  }
1811 
1812  if (signalEvents_->fixed()) {
1813  return;
1814  }
1815 
1816  // Subtract background events (if any) from signal.
1817  Double_t signalEvents = nTotEvts;
1818  if (usingBkgnd_ == kTRUE) {
1819  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1820  signalEvents -= (*iter)->value();
1821  }
1822  }
1823  signalEvents_->value(signalEvents);
1824 }
1825 
1827 {
1828  // Fill the internal data trees of the signal and background models.
1829  // Note that we store the events of both charges in both the
1830  // negative and the positive models. It's only later, at the stage
1831  // when the likelihood is being calculated, that we separate them.
1832 
1833  LauFitDataTree* inputFitData = this->fitData();
1834 
1835  // First the Dalitz plot variables (m_ij^2)
1836  if (this->useDP() == kTRUE) {
1837 
1838  // need to append SCF smearing bins before caching DP amplitudes
1839  if ( scfMap_ != 0 ) {
1840  this->appendBinCentres( inputFitData );
1841  }
1842 
1843  negSigModel_->fillDataTree(*inputFitData);
1844  posSigModel_->fillDataTree(*inputFitData);
1845 
1846  if (usingBkgnd_ == kTRUE) {
1847  for (LauBkgndDPModelList::iterator iter = negBkgndDPModels_.begin(); iter != negBkgndDPModels_.end(); ++iter) {
1848  (*iter)->fillDataTree(*inputFitData);
1849  }
1850  for (LauBkgndDPModelList::iterator iter = posBkgndDPModels_.begin(); iter != posBkgndDPModels_.end(); ++iter) {
1851  (*iter)->fillDataTree(*inputFitData);
1852  }
1853  }
1854  }
1855 
1856  // ...and then the extra PDFs
1857  this->cacheInfo(negSignalPdfs_, *inputFitData);
1858  this->cacheInfo(negScfPdfs_, *inputFitData);
1859  for (LauBkgndPdfsList::iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
1860  this->cacheInfo((*iter), *inputFitData);
1861  }
1862  if ( tagged_ ) {
1863  this->cacheInfo(posSignalPdfs_, *inputFitData);
1864  this->cacheInfo(posScfPdfs_, *inputFitData);
1865  for (LauBkgndPdfsList::iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
1866  this->cacheInfo((*iter), *inputFitData);
1867  }
1868  }
1869 
1870  // the SCF fractions and jacobians
1871  if ( useSCF_ && useSCFHist_ ) {
1872  if ( !inputFitData->haveBranch( "m13Sq" ) || !inputFitData->haveBranch( "m23Sq" ) ) {
1873  cerr << "ERROR in LauCPFitModel::cacheInputFitVars : Input data does not contain DP branches and so can't cache the SCF fraction." << endl;
1874  gSystem->Exit(EXIT_FAILURE);
1875  }
1876 
1877  UInt_t nEvents = inputFitData->nEvents();
1878  recoSCFFracs_.clear();
1879  recoSCFFracs_.reserve( nEvents );
1880  if ( negKinematics_->squareDP() ) {
1881  recoJacobians_.clear();
1882  recoJacobians_.reserve( nEvents );
1883  }
1884  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1885  const LauFitData& dataValues = inputFitData->getData(iEvt);
1886  LauFitData::const_iterator m13_iter = dataValues.find("m13Sq");
1887  LauFitData::const_iterator m23_iter = dataValues.find("m23Sq");
1888  negKinematics_->updateKinematics( m13_iter->second, m23_iter->second );
1889  Double_t scfFrac = scfFracHist_->calcEfficiency( negKinematics_ );
1890  recoSCFFracs_.push_back( scfFrac );
1891  if ( negKinematics_->squareDP() ) {
1893  }
1894  }
1895  }
1896 
1897  // finally cache the event charge
1898  evtCharges_.clear();
1899  if ( tagged_ ) {
1900  if ( !inputFitData->haveBranch( tagVarName_ ) ) {
1901  cerr << "ERROR in LauCPFitModel::cacheInputFitVars : Input data does not contain branch \"" << tagVarName_ << "\"." << endl;
1902  gSystem->Exit(EXIT_FAILURE);
1903  }
1904  UInt_t nEvents = inputFitData->nEvents();
1905  evtCharges_.reserve( nEvents );
1906  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1907  const LauFitData& dataValues = inputFitData->getData(iEvt);
1908  LauFitData::const_iterator iter = dataValues.find( tagVarName_ );
1909  curEvtCharge_ = static_cast<Int_t>( iter->second );
1910  evtCharges_.push_back( curEvtCharge_ );
1911  }
1912  }
1913 }
1914 
1916 {
1917  // We'll be caching the DP amplitudes and efficiencies of the centres of the true bins.
1918  // To do so, we attach some fake points at the end of inputData, the number of the entry
1919  // minus the total number of events corresponding to the number of the histogram for that
1920  // given true bin in the LauScfMap object. (What this means is that when Laura is provided with
1921  // the LauScfMap object by the user, it's the latter who has to make sure that it contains the
1922  // right number of histograms and in exactly the right order!)
1923 
1924  // Get the x and y co-ordinates of the bin centres
1925  std::vector<Double_t> binCentresXCoords;
1926  std::vector<Double_t> binCentresYCoords;
1927  scfMap_->listBinCentres(binCentresXCoords, binCentresYCoords);
1928 
1929  // The SCF histograms could be in square Dalitz plot histograms.
1930  // The dynamics takes normal Dalitz plot coords, so we might have to convert them back.
1931  Bool_t sqDP = negKinematics_->squareDP();
1932  UInt_t nBins = binCentresXCoords.size();
1933  fakeSCFFracs_.clear();
1934  fakeSCFFracs_.reserve( nBins );
1935  if ( sqDP ) {
1936  fakeJacobians_.clear();
1937  fakeJacobians_.reserve( nBins );
1938  }
1939 
1940  for (UInt_t iBin = 0; iBin < nBins; ++iBin) {
1941 
1942  if ( sqDP ) {
1943  negKinematics_->updateSqDPKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1944  binCentresXCoords[iBin] = negKinematics_->getm13Sq();
1945  binCentresYCoords[iBin] = negKinematics_->getm23Sq();
1947  } else {
1948  negKinematics_->updateKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1949  }
1950 
1952  }
1953 
1954  // Set up inputFitVars_ object to hold the fake events
1955  inputData->appendFakePoints(binCentresXCoords,binCentresYCoords);
1956 }
1957 
1959 {
1960  // Find out whether we have B- or B+
1961  if ( tagged_ ) {
1962  curEvtCharge_ = evtCharges_[iEvt];
1963 
1964  // check that the charge is either +1 or -1
1965  if (TMath::Abs(curEvtCharge_)!=1) {
1966  cerr << "ERROR in LauCPFitModel::getTotEvtLikelihood : Charge/tag not accepted value: " << curEvtCharge_ << endl;
1967  if (curEvtCharge_>0) {
1968  curEvtCharge_ = +1;
1969  } else {
1970  curEvtCharge_ = -1;
1971  }
1972  cerr << " : Making it: " << curEvtCharge_ << "." << endl;
1973  }
1974  }
1975 
1976  // Get the DP likelihood for signal and backgrounds
1977  this->getEvtDPLikelihood(iEvt);
1978 
1979  // Get the combined extra PDFs likelihood for signal and backgrounds
1980  this->getEvtExtraLikelihoods(iEvt);
1981 
1982  // If appropriate, combine the TM and SCF likelihoods
1983  Double_t sigLike = sigDPLike_ * sigExtraLike_;
1984  if ( useSCF_ ) {
1985  Double_t scfFrac(0.0);
1986  if (useSCFHist_) {
1987  scfFrac = recoSCFFracs_[iEvt];
1988  } else {
1989  scfFrac = scfFrac_.value();
1990  }
1991  sigLike *= (1.0 - scfFrac);
1992  if ( (scfMap_ != 0) && (this->useDP() == kTRUE) ) {
1993  // if we're smearing the SCF DP PDF then the SCF frac
1994  // is already included in the SCF DP likelihood
1995  sigLike += (scfDPLike_ * scfExtraLike_);
1996  } else {
1997  sigLike += (scfFrac * scfDPLike_ * scfExtraLike_);
1998  }
1999  }
2000 
2001  // Get the correct event fractions depending on the charge
2002  // Signal asymmetry is built into the DP model... but when the DP
2003  // isn't in the fit we need an explicit parameter
2004  Double_t signalEvents = signalEvents_->value() * 0.5;
2005  if (this->useDP() == kFALSE) {
2006  signalEvents *= (1.0 - curEvtCharge_ * signalAsym_->value());
2007  }
2008 
2009  // Construct the total event likelihood
2010  Double_t likelihood(0.0);
2011  if (usingBkgnd_) {
2012  likelihood = sigLike*signalEvents;
2013  const UInt_t nBkgnds = this->nBkgndClasses();
2014  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2015  Double_t bkgndEvents = bkgndEvents_[bkgndID]->value() * 0.5 * (1.0 - curEvtCharge_ * bkgndAsym_[bkgndID]->value());
2016  likelihood += bkgndEvents*bkgndDPLike_[bkgndID]*bkgndExtraLike_[bkgndID];
2017  }
2018  } else {
2019  likelihood = sigLike*0.5;
2020  }
2021  return likelihood;
2022 }
2023 
2025 {
2026  Double_t eventSum(0.0);
2027  eventSum += signalEvents_->value();
2028  if (usingBkgnd_) {
2029  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
2030  eventSum += (*iter)->value();
2031  }
2032  }
2033  return eventSum;
2034 }
2035 
2037 {
2038  // Function to return the signal and background likelihoods for the
2039  // Dalitz plot for the given event evtNo.
2040 
2041  if ( ! this->useDP() ) {
2042  // There's always going to be a term in the likelihood for the
2043  // signal, so we'd better not zero it.
2044  sigDPLike_ = 1.0;
2045  scfDPLike_ = 1.0;
2046 
2047  const UInt_t nBkgnds = this->nBkgndClasses();
2048  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2049  if (usingBkgnd_ == kTRUE) {
2050  bkgndDPLike_[bkgndID] = 1.0;
2051  } else {
2052  bkgndDPLike_[bkgndID] = 0.0;
2053  }
2054  }
2055 
2056  return;
2057  }
2058 
2059  const UInt_t nBkgnds = this->nBkgndClasses();
2060  if ( tagged_ ) {
2061  if (curEvtCharge_==+1) {
2065 
2066  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2067  if (usingBkgnd_ == kTRUE) {
2068  bkgndDPLike_[bkgndID] = posBkgndDPModels_[bkgndID]->getLikelihood(iEvt);
2069  } else {
2070  bkgndDPLike_[bkgndID] = 0.0;
2071  }
2072  }
2073  } else {
2077 
2078  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2079  if (usingBkgnd_ == kTRUE) {
2080  bkgndDPLike_[bkgndID] = negBkgndDPModels_[bkgndID]->getLikelihood(iEvt);
2081  } else {
2082  bkgndDPLike_[bkgndID] = 0.0;
2083  }
2084  }
2085  }
2086  } else {
2089 
2092 
2093  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2094  if (usingBkgnd_ == kTRUE) {
2095  bkgndDPLike_[bkgndID] = 0.5 * ( posBkgndDPModels_[bkgndID]->getLikelihood(iEvt) +
2096  negBkgndDPModels_[bkgndID]->getLikelihood(iEvt) );
2097  } else {
2098  bkgndDPLike_[bkgndID] = 0.0;
2099  }
2100  }
2101  }
2102 
2103  if ( useSCF_ == kTRUE ) {
2104  if ( scfMap_ == 0 ) {
2105  // we're not smearing the SCF DP position
2106  // so the likelihood is the same as the TM
2108  } else {
2109  // calculate the smeared SCF DP likelihood
2110  scfDPLike_ = this->getEvtSCFDPLikelihood(iEvt);
2111  }
2112  }
2113 
2114  // Calculate the signal normalisation
2115  // NB the 2.0 is there so that the 0.5 factor is applied to
2116  // signal and background in the same place otherwise you get
2117  // normalisation problems when you switch off the DP in the fit
2118  Double_t norm = negSigModel_->getDPNorm() + posSigModel_->getDPNorm();
2119  sigDPLike_ *= 2.0/norm;
2120  scfDPLike_ *= 2.0/norm;
2121 }
2122 
2124 {
2125  Double_t scfDPLike(0.0);
2126 
2127  Double_t recoJacobian(1.0);
2128  Double_t xCoord(0.0);
2129  Double_t yCoord(0.0);
2130 
2131  Bool_t squareDP = negKinematics_->squareDP();
2132  if ( squareDP ) {
2133  xCoord = negSigModel_->getEvtmPrime();
2134  yCoord = negSigModel_->getEvtthPrime();
2135  recoJacobian = recoJacobians_[iEvt];
2136  } else {
2137  xCoord = negSigModel_->getEvtm13Sq();
2138  yCoord = negSigModel_->getEvtm23Sq();
2139  }
2140 
2141  // Find the bin that our reco event falls in
2142  Int_t recoBin = scfMap_->binNumber( xCoord, yCoord );
2143 
2144  // Find out which true Bins contribute to the given reco bin
2145  const std::vector<Int_t>* trueBins = scfMap_->trueBins(recoBin);
2146 
2147  const Int_t nDataEvents = this->eventsPerExpt();
2148 
2149  // Loop over the true bins
2150  for (std::vector<Int_t>::const_iterator iter = trueBins->begin(); iter != trueBins->end(); ++iter)
2151  {
2152  Int_t trueBin = (*iter);
2153 
2154  // prob of a true event in the given true bin migrating to the reco bin
2155  Double_t pRecoGivenTrue = scfMap_->prob( recoBin, trueBin );
2156  Double_t pTrue(0.0);
2157 
2158  // We've cached the DP amplitudes and the efficiency for the
2159  // true bin centres, just after the data points
2160  if ( tagged_ ) {
2161  LauAbsDPDynamics* sigModel(0);
2162  if (curEvtCharge_<0) {
2163  sigModel = negSigModel_;
2164  } else {
2165  sigModel = posSigModel_;
2166  }
2167 
2168  sigModel->calcLikelihoodInfo( nDataEvents + trueBin );
2169 
2170  pTrue = sigModel->getEvtDPAmp().abs2() * sigModel->getEvtEff();
2171  } else {
2172  posSigModel_->calcLikelihoodInfo( nDataEvents + trueBin );
2173  negSigModel_->calcLikelihoodInfo( nDataEvents + trueBin );
2174 
2175  pTrue = 0.5 * ( posSigModel_->getEvtDPAmp().abs2() * posSigModel_->getEvtEff() +
2177  }
2178 
2179  // Get the cached SCF fraction (and jacobian if we're using the square DP)
2180  Double_t scfFraction = fakeSCFFracs_[ trueBin ];
2181  Double_t jacobian(1.0);
2182  if ( squareDP ) {
2183  jacobian = fakeJacobians_[ trueBin ];
2184  }
2185 
2186  scfDPLike += pTrue * jacobian * scfFraction * pRecoGivenTrue;
2187  }
2188 
2189  // Divide by the reco jacobian
2190  scfDPLike /= recoJacobian;
2191 
2192  return scfDPLike;
2193 }
2194 
2196 {
2197  // Function to return the signal and background likelihoods for the
2198  // extra variables for the given event evtNo.
2199 
2200  sigExtraLike_ = 1.0;
2201 
2202  const UInt_t nBkgnds = this->nBkgndClasses();
2203 
2204  if ( ! tagged_ || curEvtCharge_ < 0 ) {
2205  sigExtraLike_ = this->prodPdfValue( negSignalPdfs_, iEvt );
2206 
2207  if (useSCF_) {
2208  scfExtraLike_ = this->prodPdfValue( negScfPdfs_, iEvt );
2209  }
2210 
2211  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2212  if (usingBkgnd_) {
2213  bkgndExtraLike_[bkgndID] = this->prodPdfValue( negBkgndPdfs_[bkgndID], iEvt );
2214  } else {
2215  bkgndExtraLike_[bkgndID] = 0.0;
2216  }
2217  }
2218  } else {
2219  sigExtraLike_ = this->prodPdfValue( posSignalPdfs_, iEvt );
2220 
2221  if (useSCF_) {
2222  scfExtraLike_ = this->prodPdfValue( posScfPdfs_, iEvt );
2223  }
2224 
2225  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2226  if (usingBkgnd_) {
2227  bkgndExtraLike_[bkgndID] = this->prodPdfValue( posBkgndPdfs_[bkgndID], iEvt );
2228  } else {
2229  bkgndExtraLike_[bkgndID] = 0.0;
2230  }
2231  }
2232  }
2233 }
2234 
2236 {
2237  negCoeffs_.clear(); posCoeffs_.clear();
2238  negCoeffs_.reserve(nSigComp_); posCoeffs_.reserve(nSigComp_);
2239  for (UInt_t i = 0; i < nSigComp_; i++) {
2240  negCoeffs_.push_back(coeffPars_[i]->antiparticleCoeff());
2241  posCoeffs_.push_back(coeffPars_[i]->particleCoeff());
2242  }
2243 }
2244 
2246 {
2247  // add branches for storing the experiment number and the number of
2248  // the event within the current experiment
2249  this->addSPlotNtupleIntegerBranch("iExpt");
2250  this->addSPlotNtupleIntegerBranch("iEvtWithinExpt");
2251 
2252  // Store the efficiency of the event (for inclusive BF calculations).
2253  if (this->storeDPEff()) {
2254  this->addSPlotNtupleDoubleBranch("efficiency");
2256  this->addSPlotNtupleDoubleBranch("scffraction");
2257  }
2258  }
2259 
2260  // Store the total event likelihood for each species.
2261  if (useSCF_) {
2262  this->addSPlotNtupleDoubleBranch("sigTMTotalLike");
2263  this->addSPlotNtupleDoubleBranch("sigSCFTotalLike");
2264  this->addSPlotNtupleDoubleBranch("sigSCFFrac");
2265  } else {
2266  this->addSPlotNtupleDoubleBranch("sigTotalLike");
2267  }
2268  if (usingBkgnd_) {
2269  const UInt_t nBkgnds = this->nBkgndClasses();
2270  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2271  TString name( this->bkgndClassName(iBkgnd) );
2272  name += "TotalLike";
2273  this->addSPlotNtupleDoubleBranch(name);
2274  }
2275  }
2276 
2277  // Store the DP likelihoods
2278  if (this->useDP()) {
2279  if (useSCF_) {
2280  this->addSPlotNtupleDoubleBranch("sigTMDPLike");
2281  this->addSPlotNtupleDoubleBranch("sigSCFDPLike");
2282  } else {
2283  this->addSPlotNtupleDoubleBranch("sigDPLike");
2284  }
2285  if (usingBkgnd_) {
2286  const UInt_t nBkgnds = this->nBkgndClasses();
2287  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2288  TString name( this->bkgndClassName(iBkgnd) );
2289  name += "DPLike";
2290  this->addSPlotNtupleDoubleBranch(name);
2291  }
2292  }
2293  }
2294 
2295  // Store the likelihoods for each extra PDF
2296  if (useSCF_) {
2297  this->addSPlotNtupleBranches(&negSignalPdfs_, "sigTM");
2298  this->addSPlotNtupleBranches(&negScfPdfs_, "sigSCF");
2299  } else {
2300  this->addSPlotNtupleBranches(&negSignalPdfs_, "sig");
2301  }
2302  if (usingBkgnd_) {
2303  const UInt_t nBkgnds = this->nBkgndClasses();
2304  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2305  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2306  const LauPdfList* pdfList = &(negBkgndPdfs_[iBkgnd]);
2307  this->addSPlotNtupleBranches(pdfList, bkgndClass);
2308  }
2309  }
2310 }
2311 
2312 void LauCPFitModel::addSPlotNtupleBranches(const LauPdfList* extraPdfs, const TString& prefix)
2313 {
2314  if (extraPdfs) {
2315  // Loop through each of the PDFs
2316  for (LauPdfList::const_iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
2317 
2318  // Count the number of input variables that are not
2319  // DP variables (used in the case where there is DP
2320  // dependence for e.g. MVA)
2321  UInt_t nVars(0);
2322  std::vector<TString> varNames = (*pdf_iter)->varNames();
2323  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2324  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2325  ++nVars;
2326  }
2327  }
2328 
2329  if ( nVars == 1 ) {
2330  // If the PDF only has one variable then
2331  // simply add one branch for that variable
2332  TString varName = (*pdf_iter)->varName();
2333  TString name(prefix);
2334  name += varName;
2335  name += "Like";
2336  this->addSPlotNtupleDoubleBranch(name);
2337  } else if ( nVars == 2 ) {
2338  // If the PDF has two variables then we
2339  // need a branch for them both together and
2340  // branches for each
2341  TString allVars("");
2342  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2343  allVars += (*var_iter);
2344  TString name(prefix);
2345  name += (*var_iter);
2346  name += "Like";
2347  this->addSPlotNtupleDoubleBranch(name);
2348  }
2349  TString name(prefix);
2350  name += allVars;
2351  name += "Like";
2352  this->addSPlotNtupleDoubleBranch(name);
2353  } else {
2354  cerr << "WARNING in LauCPFitModel::addSPlotNtupleBranches : Can't yet deal with 3D PDFs." << endl;
2355  }
2356  }
2357  }
2358 }
2359 
2360 Double_t LauCPFitModel::setSPlotNtupleBranchValues(LauPdfList* extraPdfs, const TString& prefix, UInt_t iEvt)
2361 {
2362  // Store the PDF value for each variable in the list
2363  Double_t totalLike(1.0);
2364  Double_t extraLike(0.0);
2365  if (extraPdfs) {
2366  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
2367 
2368  // calculate the likelihood for this event
2369  (*pdf_iter)->calcLikelihoodInfo(iEvt);
2370  extraLike = (*pdf_iter)->getLikelihood();
2371  totalLike *= extraLike;
2372 
2373  // Count the number of input variables that are not
2374  // DP variables (used in the case where there is DP
2375  // dependence for e.g. MVA)
2376  UInt_t nVars(0);
2377  std::vector<TString> varNames = (*pdf_iter)->varNames();
2378  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2379  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2380  ++nVars;
2381  }
2382  }
2383 
2384  if ( nVars == 1 ) {
2385  // If the PDF only has one variable then
2386  // simply store the value for that variable
2387  TString varName = (*pdf_iter)->varName();
2388  TString name(prefix);
2389  name += varName;
2390  name += "Like";
2391  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
2392  } else if ( nVars == 2 ) {
2393  // If the PDF has two variables then we
2394  // store the value for them both together
2395  // and for each on their own
2396  TString allVars("");
2397  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2398  allVars += (*var_iter);
2399  TString name(prefix);
2400  name += (*var_iter);
2401  name += "Like";
2402  Double_t indivLike = (*pdf_iter)->getLikelihood( (*var_iter) );
2403  this->setSPlotNtupleDoubleBranchValue(name, indivLike);
2404  }
2405  TString name(prefix);
2406  name += allVars;
2407  name += "Like";
2408  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
2409  } else {
2410  cerr << "WARNING in LauCPFitModel::setSPlotNtupleBranchValues : Can't yet deal with 3D PDFs." << endl;
2411  }
2412  }
2413  }
2414  return totalLike;
2415 }
2416 
2418 {
2419  LauSPlot::NameSet nameSet;
2420  if (this->useDP()) {
2421  nameSet.insert("DP");
2422  }
2423  // Loop through all the signal PDFs
2424  for (LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter) {
2425  // Loop over the variables involved in each PDF
2426  std::vector<TString> varNames = (*pdf_iter)->varNames();
2427  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2428  // If they are not DP coordinates then add them
2429  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2430  nameSet.insert( (*var_iter) );
2431  }
2432  }
2433  }
2434  return nameSet;
2435 }
2436 
2438 {
2439  LauSPlot::NumbMap numbMap;
2440  if (!signalEvents_->fixed() && this->doEMLFit()) {
2441  numbMap["sig"] = signalEvents_->genValue();
2442  }
2443  if ( usingBkgnd_ ) {
2444  const UInt_t nBkgnds = this->nBkgndClasses();
2445  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2446  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2447  const LauParameter* par = bkgndEvents_[iBkgnd];
2448  if (!par->fixed()) {
2449  numbMap[bkgndClass] = par->genValue();
2450  }
2451  }
2452  }
2453  return numbMap;
2454 }
2455 
2457 {
2458  LauSPlot::NumbMap numbMap;
2459  if (signalEvents_->fixed() && this->doEMLFit()) {
2460  numbMap["sig"] = signalEvents_->genValue();
2461  }
2462  if ( usingBkgnd_ ) {
2463  const UInt_t nBkgnds = this->nBkgndClasses();
2464  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2465  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2466  const LauParameter* par = bkgndEvents_[iBkgnd];
2467  if (par->fixed()) {
2468  numbMap[bkgndClass] = par->genValue();
2469  }
2470  }
2471  }
2472  return numbMap;
2473 }
2474 
2476 {
2477  // This makes the assumption that the form of the positive and
2478  // negative PDFs are the same, which seems reasonable to me
2479 
2480  LauSPlot::TwoDMap twodimMap;
2481 
2482  for (LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter) {
2483  // Count the number of input variables that are not DP variables
2484  UInt_t nVars(0);
2485  std::vector<TString> varNames = (*pdf_iter)->varNames();
2486  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2487  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2488  ++nVars;
2489  }
2490  }
2491  if ( nVars == 2 ) {
2492  if (useSCF_) {
2493  twodimMap.insert( std::make_pair( "sigTM", std::make_pair( varNames[0], varNames[1] ) ) );
2494  } else {
2495  twodimMap.insert( std::make_pair( "sig", std::make_pair( varNames[0], varNames[1] ) ) );
2496  }
2497  }
2498  }
2499 
2500  if ( useSCF_ ) {
2501  for (LauPdfList::const_iterator pdf_iter = negScfPdfs_.begin(); pdf_iter != negScfPdfs_.end(); ++pdf_iter) {
2502  // Count the number of input variables that are not DP variables
2503  UInt_t nVars(0);
2504  std::vector<TString> varNames = (*pdf_iter)->varNames();
2505  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2506  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2507  ++nVars;
2508  }
2509  }
2510  if ( nVars == 2 ) {
2511  twodimMap.insert( std::make_pair( "sigSCF", std::make_pair( varNames[0], varNames[1] ) ) );
2512  }
2513  }
2514  }
2515 
2516  if (usingBkgnd_) {
2517  const UInt_t nBkgnds = this->nBkgndClasses();
2518  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2519  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2520  const LauPdfList& pdfList = negBkgndPdfs_[iBkgnd];
2521  for (LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
2522  // Count the number of input variables that are not DP variables
2523  UInt_t nVars(0);
2524  std::vector<TString> varNames = (*pdf_iter)->varNames();
2525  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2526  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2527  ++nVars;
2528  }
2529  }
2530  if ( nVars == 2 ) {
2531  twodimMap.insert( std::make_pair( bkgndClass, std::make_pair( varNames[0], varNames[1] ) ) );
2532  }
2533  }
2534  }
2535  }
2536 
2537  return twodimMap;
2538 }
2539 
2541 {
2542  cout << "Storing per-event likelihood values..." << endl;
2543 
2544  // if we've not been using the DP model then we need to cache all
2545  // the info here so that we can get the efficiency from it
2546 
2547  LauFitDataTree* inputFitData = this->fitData();
2548 
2549  if (!this->useDP() && this->storeDPEff()) {
2552  negSigModel_->fillDataTree(*inputFitData);
2553  posSigModel_->fillDataTree(*inputFitData);
2554  }
2555 
2556  UInt_t evtsPerExpt(this->eventsPerExpt());
2557 
2558  LauAbsDPDynamics* sigModel(0);
2559  LauPdfList* sigPdfs(0);
2560  LauPdfList* scfPdfs(0);
2561  LauBkgndPdfsList* bkgndPdfs(0);
2562 
2563  for (UInt_t iEvt = 0; iEvt < evtsPerExpt; ++iEvt) {
2564 
2565  this->setSPlotNtupleIntegerBranchValue("iExpt",this->iExpt());
2566  this->setSPlotNtupleIntegerBranchValue("iEvtWithinExpt",iEvt);
2567 
2568  // Find out whether we have B- or B+
2569  if ( tagged_ ) {
2570  const LauFitData& dataValues = inputFitData->getData(iEvt);
2571  LauFitData::const_iterator iter = dataValues.find("charge");
2572  curEvtCharge_ = static_cast<Int_t>(iter->second);
2573 
2574  if (curEvtCharge_==+1) {
2575  sigModel = posSigModel_;
2576  sigPdfs = &posSignalPdfs_;
2577  scfPdfs = &posScfPdfs_;
2578  bkgndPdfs = &posBkgndPdfs_;
2579  } else {
2580  sigModel = negSigModel_;
2581  sigPdfs = &negSignalPdfs_;
2582  scfPdfs = &negScfPdfs_;
2583  bkgndPdfs = &negBkgndPdfs_;
2584  }
2585  } else {
2586  sigPdfs = &negSignalPdfs_;
2587  scfPdfs = &negScfPdfs_;
2588  bkgndPdfs = &negBkgndPdfs_;
2589  }
2590 
2591  // the DP information
2592  this->getEvtDPLikelihood(iEvt);
2593  if (this->storeDPEff()) {
2594  if (!this->useDP()) {
2597  }
2598  if ( tagged_ ) {
2599  this->setSPlotNtupleDoubleBranchValue("efficiency",sigModel->getEvtEff());
2601  this->setSPlotNtupleDoubleBranchValue("scffraction",sigModel->getEvtScfFraction());
2602  }
2603  } else {
2607  }
2608  }
2609  }
2610  if (this->useDP()) {
2612  if (useSCF_) {
2614  this->setSPlotNtupleDoubleBranchValue("sigTMDPLike",sigDPLike_);
2615  this->setSPlotNtupleDoubleBranchValue("sigSCFDPLike",scfDPLike_);
2616  } else {
2617  this->setSPlotNtupleDoubleBranchValue("sigDPLike",sigDPLike_);
2618  }
2619  if (usingBkgnd_) {
2620  const UInt_t nBkgnds = this->nBkgndClasses();
2621  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2622  TString name = this->bkgndClassName(iBkgnd);
2623  name += "DPLike";
2624  this->setSPlotNtupleDoubleBranchValue(name,bkgndDPLike_[iBkgnd]);
2625  }
2626  }
2627  } else {
2628  sigTotalLike_ = 1.0;
2629  if (useSCF_) {
2630  scfTotalLike_ = 1.0;
2631  }
2632  if (usingBkgnd_) {
2633  const UInt_t nBkgnds = this->nBkgndClasses();
2634  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2635  bkgndTotalLike_[iBkgnd] = 1.0;
2636  }
2637  }
2638  }
2639 
2640  // the signal PDF values
2641  if ( useSCF_ ) {
2642  sigTotalLike_ *= this->setSPlotNtupleBranchValues(sigPdfs, "sigTM", iEvt);
2643  scfTotalLike_ *= this->setSPlotNtupleBranchValues(scfPdfs, "sigSCF", iEvt);
2644  } else {
2645  sigTotalLike_ *= this->setSPlotNtupleBranchValues(sigPdfs, "sig", iEvt);
2646  }
2647 
2648  // the background PDF values
2649  if (usingBkgnd_) {
2650  const UInt_t nBkgnds = this->nBkgndClasses();
2651  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2652  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2653  LauPdfList& pdfs = (*bkgndPdfs)[iBkgnd];
2654  bkgndTotalLike_[iBkgnd] *= this->setSPlotNtupleBranchValues(&(pdfs), bkgndClass, iEvt);
2655  }
2656  }
2657 
2658  // the total likelihoods
2659  if (useSCF_) {
2660  Double_t scfFrac(0.0);
2661  if ( useSCFHist_ ) {
2662  scfFrac = recoSCFFracs_[iEvt];
2663  } else {
2664  scfFrac = scfFrac_.value();
2665  }
2666  this->setSPlotNtupleDoubleBranchValue("sigSCFFrac",scfFrac);
2667  sigTotalLike_ *= ( 1.0 - scfFrac );
2668  if ( scfMap_ == 0 ) {
2669  scfTotalLike_ *= scfFrac;
2670  }
2671  this->setSPlotNtupleDoubleBranchValue("sigTMTotalLike",sigTotalLike_);
2672  this->setSPlotNtupleDoubleBranchValue("sigSCFTotalLike",scfTotalLike_);
2673  } else {
2674  this->setSPlotNtupleDoubleBranchValue("sigTotalLike",sigTotalLike_);
2675  }
2676  if (usingBkgnd_) {
2677  const UInt_t nBkgnds = this->nBkgndClasses();
2678  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2679  TString name = this->bkgndClassName(iBkgnd);
2680  name += "TotalLike";
2682  }
2683  }
2684  // fill the tree
2685  this->fillSPlotNtupleBranches();
2686  }
2687  cout << "Finished storing per-event likelihood values." << endl;
2688 }
2689 
2690 void LauCPFitModel::embedNegSignal(const TString& fileName, const TString& treeName,
2691  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
2692  Bool_t useReweighting)
2693 {
2694  if (negSignalTree_) {
2695  cerr << "ERROR in LauCPFitModel::embedNegSignal : Already embedding signal from a file." << endl;
2696  return;
2697  }
2698 
2699  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2700  cerr << "WARNING in LauCPFitModel::embedNegSignal : Conflicting options provided, will not reuse events at all." << endl;
2701  reuseEventsWithinExperiment = kFALSE;
2702  }
2703 
2704  negSignalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2705  Bool_t dataOK = negSignalTree_->findBranches();
2706  if (!dataOK) {
2707  delete negSignalTree_; negSignalTree_ = 0;
2708  cerr << "ERROR in LauCPFitModel::embedNegSignal : Problem creating data tree for embedding." << endl;
2709  return;
2710  }
2711  reuseSignal_ = reuseEventsWithinEnsemble;
2712  useNegReweighting_ = useReweighting;
2713  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2714 }
2715 
2716 void LauCPFitModel::embedNegBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
2717  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
2718 {
2719  if ( ! this->validBkgndClass( bkgndClass ) ) {
2720  cerr << "ERROR in LauCPFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
2721  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
2722  return;
2723  }
2724 
2725  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
2726 
2727  if (negBkgndTree_[bkgndID]) {
2728  cerr << "ERROR in LauCPFitModel::embedNegBkgnd : Already embedding background from a file." << endl;
2729  return;
2730  }
2731 
2732  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2733  cerr << "WARNING in LauCPFitModel::embedNegBkgnd : Conflicting options provided, will not reuse events at all." << endl;
2734  reuseEventsWithinExperiment = kFALSE;
2735  }
2736 
2737  negBkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2738  Bool_t dataOK = negBkgndTree_[bkgndID]->findBranches();
2739  if (!dataOK) {
2740  delete negBkgndTree_[bkgndID]; negBkgndTree_[bkgndID] = 0;
2741  cerr << "ERROR in LauCPFitModel::embedNegBkgnd : Problem creating data tree for embedding." << endl;
2742  return;
2743  }
2744  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
2745  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2746 }
2747 
2748 void LauCPFitModel::embedPosSignal(const TString& fileName, const TString& treeName,
2749  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
2750  Bool_t useReweighting)
2751 {
2752  if (posSignalTree_) {
2753  cerr << "ERROR in LauCPFitModel::embedPosSignal : Already embedding signal from a file." << endl;
2754  return;
2755  }
2756 
2757  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2758  cerr << "WARNING in LauCPFitModel::embedPosSignal : Conflicting options provided, will not reuse events at all." << endl;
2759  reuseEventsWithinExperiment = kFALSE;
2760  }
2761 
2762  posSignalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2763  Bool_t dataOK = posSignalTree_->findBranches();
2764  if (!dataOK) {
2765  delete posSignalTree_; posSignalTree_ = 0;
2766  cerr << "ERROR in LauCPFitModel::embedPosSignal : Problem creating data tree for embedding." << endl;
2767  return;
2768  }
2769  reuseSignal_ = reuseEventsWithinEnsemble;
2770  usePosReweighting_ = useReweighting;
2771  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2772 }
2773 
2774 void LauCPFitModel::embedPosBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
2775  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
2776 {
2777  if ( ! this->validBkgndClass( bkgndClass ) ) {
2778  cerr << "ERROR in LauCPFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
2779  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
2780  return;
2781  }
2782 
2783  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
2784 
2785  if (posBkgndTree_[bkgndID]) {
2786  cerr << "ERROR in LauCPFitModel::embedPosBkgnd : Already embedding background from a file." << endl;
2787  return;
2788  }
2789 
2790  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2791  cerr << "WARNING in LauCPFitModel::embedPosBkgnd : Conflicting options provided, will not reuse events at all." << endl;
2792  reuseEventsWithinExperiment = kFALSE;
2793  }
2794 
2795  posBkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2796  Bool_t dataOK = posBkgndTree_[bkgndID]->findBranches();
2797  if (!dataOK) {
2798  delete posBkgndTree_[bkgndID]; posBkgndTree_[bkgndID] = 0;
2799  cerr << "ERROR in LauCPFitModel::embedPosBkgnd : Problem creating data tree for embedding." << endl;
2800  return;
2801  }
2802  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
2803  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2804 }
2805 
2806 void LauCPFitModel::weightEvents( const TString& /*dataFileName*/, const TString& /*dataTreeName*/ )
2807 {
2808  cerr << "ERROR in LauCPFitModel::weightEvents : Method not available for this fit model." << endl;
2809  return;
2810 }
2811 
Double_t range() const
The range allowed for the parameter.
void generateExtraPdfValues(LauPdfList *extraPdfs, LauEmbeddedData *embeddedData)
Generate from the extra PDFs.
void splitSignalComponent(const TH2 *dpHisto, Bool_t upperHalf=kFALSE, LauScfMap *scfMap=0)
Split the signal component into well-reconstructed and mis-reconstructed parts.
Bool_t generateSignalEvent()
Generate signal event.
virtual void addSPlotNtupleIntegerBranch(const TString &name)
Add a branch to the sPlot tree for storing an integer.
Class for defining the abstract interface for signal Dalitz plot dynamics.
The abstract interface for a background Dalitz plot model.
void storeCorrMatrix(UInt_t iExpt, Double_t NLL, Int_t fitStatus)
Store the correlation matrix and other fit information.
Definition: LauFitNtuple.cc:61
Bool_t squareDP() const
Are the square Dalitz plot co-ordinates being calculated?
std::vector< LauParameter > acp_
A_CP parameter.
virtual LauSPlot::NameSet variableNames() const
Returns the names of all variables in the fit.
Double_t getc23() const
Get the cosine of the helicity angle theta23.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
std::vector< LauParameter > LauParameterList
List of parameters.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
LauKinematics * negKinematics_
The B- Dalitz plot kinematics object.
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.
virtual void printAsymmetries(std::ostream &output)
Print the asymmetries.
Bool_t writeLatexTable() const
Determine whether writing out of the latex table is enabled.
LauPdfList negSignalPdfs_
The B- signal PDFs.
void addSPlotNtupleBranches(const LauPdfList *extraPdfs, const TString &prefix)
Add sPlot branches for the extra PDFs.
Bool_t forceAsym_
Option to force an asymmetry.
Double_t scfTotalLike_
Total SCF likelihood.
void setExtraPdfParameters()
Set the fit parameters for the extra PDFs.
void updateSqDPKinematics(Double_t mPrime, Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m&#39; and theta&#39;.
LauBkgndDPModelList posBkgndDPModels_
The B+ background Dalitz plot models.
void setEffHisto(const TH2 *effHisto, Bool_t useInterpolation=kTRUE, Bool_t fluctuateBins=kFALSE, Double_t avEff=-1.0, Double_t absError=-1.0, Bool_t useUpperHalfOnly=kFALSE, Bool_t squareDP=kFALSE)
Set the efficiency variation across the phase space using a predetermined 2D histogram.
Definition: LauEffModel.cc:63
void cacheInfo(LauPdfList &pdfList, const LauFitDataTree &theData)
Have all PDFs in the list cache the data.
virtual void setGenNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the gen tree.
File containing declaration of LauAsymmCalc class.
UInt_t nNormPar_
Number of parameters.
LauEmbeddedData * posSignalTree_
The B+ signal event tree.
virtual TString baseName() const
Retrieve the base name of the coefficient set.
LauParameter()
Default constructor.
Definition: LauParameter.cc:30
Double_t getm23() const
Get the m23 invariant mass.
Int_t curEvtCharge_
Current event charge.
void calcExtraFractions(Bool_t initValues=kFALSE)
Calculate the CP-conserving and CP-violating fit fractions.
const TString & name() const
The parameter name.
virtual LauParameter acp()=0
Calculate the CP asymmetry.
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.
Class for calculating the asymmetry between two variables.
Definition: LauAsymmCalc.hh:25
LauEmbeddedData * negSignalTree_
The B- signal event tree.
Bool_t useDP() const
Is the Dalitz plot term in the likelihood.
LauDaughters * getDaughters()
Retrieve the daughters.
Double_t getc13() const
Get the cosine of the helicity angle theta13.
File containing declaration of LauAbsCoeffSet class.
void clearUsedList()
Clear the list of used events.
Bool_t useNegReweighting_
Boolean to use reweighting for B-.
TH2 * trueHist(Int_t trueBin)
Retrieve the migration histogram for trueBin.
Definition: LauScfMap.cc:176
std::vector< LauComplex > negCoeffs_
The complex coefficients for B-.
LauBkgndPdfsList posBkgndPdfs_
The B+ background PDFs.
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
LauParArray negFitFrac_
The B- fit fractions.
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
virtual void getEvtExtraLikelihoods(UInt_t iEvt)
Determine the signal and background likelihood for the extra variables for a given event...
Double_t sigDPLike_
Signal DP likelihood value.
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
LauPdfList negScfPdfs_
The B- SCF PDFs.
Bool_t useSCF_
Is the signal split into TM and SCF.
virtual Bool_t genExpt()
Toy MC generation and fitting overloaded functions.
virtual void setGenNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the gen tree.
TString negParent_
Name of the parent particle.
TString posParent_
Name of the parent particle.
void embedNegSignal(const TString &fileName, const TString &treeName, Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment=kFALSE, Bool_t useReweighting=kFALSE)
Embed full simulation events for the B- signal, rather than generating toy from the PDFs...
File containing declaration of LauScfMap class.
std::vector< Double_t > fakeSCFFracs_
The cached values of the SCF fraction for each bin centre.
std::vector< std::vector< LauParameter > > LauParArray
Type to define an array of parameters.
LauParameter getDPRate() const
Retrieve the overall Dalitz plot rate.
Int_t fitStatus() const
Access the fit status information.
Bool_t storeSignalMCMatch(LauEmbeddedData *embeddedData)
Store the MC truth info on the TM/SCF nature of the embedded signal event.
Bool_t usingBkgnd_
Background boolean.
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...
LauGenInfo eventsToGenerate()
Determine the number of events to generate for each hypothesis.
Double_t prodPdfValue(LauPdfList &pdfList, UInt_t iEvt)
Calculate the product of the per-event likelihoods of the PDFs in the list.
std::vector< LauAbsPdf * > LauPdfList
List of Pdfs.
std::map< TString, Double_t > LauFitData
Type for holding event data.
LauBkgndYieldList bkgndAsym_
Background asymmetries(s)
Double_t getm23Sq() const
Get the m23 invariant mass square.
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.
std::vector< LauParameter > getExtraParameters()
Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions)
LauParameter negDPRate_
The average DP rate for B-.
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< Double_t > bkgndTotalLike_
Total background likelihood(s)
virtual void fillGenNtupleBranches()
Fill the gen tuple branches.
std::vector< Double_t > recoSCFFracs_
The cached values of the SCF fraction for each event.
LauParameter posDPRate_
The average DP rate for B+.
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.
std::vector< Double_t > bkgndExtraLike_
Background likelihood value(s) from extra PDFs.
void setupGenNtupleBranches()
Setup the required ntuple branches.
void printFitParameters(const LauPdfList &pdfList, std::ostream &fout) const
Print the fit parameters for all PDFs in the list.
LauAbsDPDynamics * posSigModel_
The B+ signal Dalitz plot model.
Double_t scfExtraLike_
SCF likelihood from extra PDFs.
Bool_t haveBranch(const TString &name) const
Boolean to determine whether branch exists.
File containing declaration of LauCPFitModel class.
Bool_t reuseSignal_
Boolean to reuse signal events.
LauParArray posFitFrac_
The B+ fit fractions.
virtual void cacheInputFitVars()
Read in the input fit data variables, e.g. m13Sq and m23Sq.
void setSignalDPParameters()
Set the fit parameters for the DP model.
void randomiseInitFitPars()
Randomise the initial fit parameters.
void appendBinCentres(LauFitDataTree *inputData)
virtual UInt_t index() const
Retrieve the index number of the coefficient set.
void updateKinematics(Double_t m13Sq, Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
virtual void checkInitFitParams()
Check the initial fit parameters.
const LauParameterList & extraPars() const
Access the extra variables.
Double_t abs2() const
Obtain the square of the absolute value of the complex number.
Definition: LauComplex.hh:229
std::vector< LauParameter * > LauParameterPList
List of parameter pointers.
virtual void fillDataTree(const LauFitDataTree &inputFitTree)=0
Obtain data from a fit tree.
const Bool_t tagged_
IS the analysis tagged?
UInt_t nExtraPdfPar_
Number of extra PDF parameters.
virtual void finaliseFitResults(const TString &tablePrefixName)
Get the fit results and store them.
Bool_t passVeto(const LauKinematics *kinematics) const
Determine whether the given DP position is outside the vetoes.
Definition: LauEffModel.cc:183
virtual void printFitFractions(std::ostream &output)
Print the fit fractions, total DP rate and mean efficiency.
Double_t calcEfficiency(const LauKinematics *kinematics) const
Determine the efficiency for a given point in the Dalitz plot.
Definition: LauEffModel.cc:130
virtual void calcExtraInfo(Bool_t init=kFALSE)=0
Calculate the fit fractions, mean efficiency and total DP rate.
Class for defining a CP fit model.
Abstract interface to the fitting and toy MC model.
UInt_t iExpt() const
Obtain the number of the current experiment.
const TString tagVarName_
Event charge.
LauParameter * signalEvents_
Signal yield.
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.
std::vector< Int_t > evtCharges_
Vector to store event charges.
const LauFitNtuple * fitNtuple() const
Access the fit ntuple.
Int_t binNumber(Double_t xCoord, Double_t yCoord) const
Find the global bin number for the given co-ordinates.
Definition: LauScfMap.cc:144
Bool_t doPoissonSmearing() const
Determine whether Poisson smearing is enabled for the toy MC generation.
Double_t getc12() const
Get the cosine of the helicity angle theta12.
Double_t getm13() const
Get the m13 invariant mass.
std::vector< LauPdfList > LauBkgndPdfsList
Typedef for a vector of background PDFs.
void calcAsymmetries(Bool_t initValues=kFALSE)
Calculate the CP asymmetries.
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.
LauScfMap * scfMap_
The smearing matrix for the SCF DP PDF.
File containing declaration of LauComplex class.
Double_t sigExtraLike_
Signal likelihood from extra PDFs.
virtual TString getConjResName(const TString &resName) const
Retrieve the name of the charge conjugate of a named resonance.
LauFitData getValues(const std::vector< TString > &names) const
Get values of specified branches.
virtual void setNBkgndEvents(LauParameter *nBkgndEvents)
Set the background event yield(s)
virtual void initialise()
Initialise the fit.
virtual LauSPlot::NumbMap fixdSpeciesNames() const
Returns the names and yields of species that are fixed in the fit.
Class for defining the abstract interface for complex coefficient classes.
Double_t nll() const
Access the current NLL value.
Double_t sigTotalLike_
Total signal likelihood.
void clearFitParVectors()
Clear the vectors containing fit parameters.
void setAmpCoeffSet(LauAbsCoeffSet *coeffSet)
Set the DP amplitude coefficients.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:31
Class to store the data for embedding in toy experiments.
virtual Double_t getEvtmPrime() const =0
Retrieve the square Dalitz plot coordinate, m&#39;, for the current event.
virtual void addSPlotNtupleDoubleBranch(const TString &name)
Add a branch to the sPlot tree for storing a double.
UInt_t nSigComp_
Number of signal components.
virtual void setupSPlotNtupleBranches()
Add branches to store experiment number and the event number within the experiment.
void clearExtraVarVectors()
Clear the vectors containing extra ntuple variables.
LauParameter * signalAsym_
Signal asymmetry.
void getEmbeddedEvent(LauKinematics *kinematics)
Retrieve an event from the data sample.
std::vector< Double_t > bkgndDPLike_
Background DP likelihood value(s)
void embedNegBkgnd(const TString &bgClass, 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...
LauBkgndEmbDataList negBkgndTree_
The B- background event tree.
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.
LauBkgndReuseEventsList reuseBkgnd_
Vector of booleans to reuse background events.
std::vector< LauAbsCoeffSet * > coeffPars_
Magnitudes and Phases.
LauParameter posMeanEff_
The mean efficiency for B+ model.
LauAbsDPDynamics * negSigModel_
The B- signal Dalitz plot model.
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.
LauParArray CPCFitFrac_
The CP conserving fit fraction.
virtual void setNSigEvents(LauParameter *nSigEvents)
Set the signal event yield.
virtual void weightEvents(const TString &dataFileName, const TString &dataTreeName)
Weight events based on the DP model.
Bool_t validBkgndClass(const TString &className) const
Check if the given background class is in the list.
Bool_t usePosReweighting_
Boolean to use reweighting for B+.
Bool_t generateBkgndEvent(UInt_t bgID)
Generate background event.
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.
Double_t getDPNorm() const
Retrieve the normalisation factor for the log-likelihood function.
Class that implements the efficiency description across the signal Dalitz plot.
Definition: LauEffModel.hh:37
TString getNameParent() const
Get name of the parent particle.
UInt_t nUsedEvents() const
Retrieve the number of events that have already been sampled.
virtual void updateCoeffs(const std::vector< LauComplex > &coeffs)
Update the complex coefficients for the resonances.
Double_t getm13Sq() const
Get the m13 invariant mass square.
void storeParsAndErrors(const std::vector< LauParameter * > &fitVars, const std::vector< LauParameter > &extraVars)
Store parameters and their errors.
virtual void updateCoeffs()
Update the coefficients.
virtual void setSPlotNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the sPlot tree.
Double_t setSPlotNtupleBranchValues(LauPdfList *extraPdfs, const TString &prefix, UInt_t iEvt)
Set the branches for the sPlot ntuple with extra PDFs.
Double_t initValue() const
The initial value of the parameter.
virtual ~LauCPFitModel()
Destructor.
std::vector< LauComplex > posCoeffs_
The complex coefficients for B+.
virtual void initialiseDPModels()
Initialise the signal DP models.
void updateSigEvents()
Update the signal events after Minuit sets background parameters.
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.
virtual void fillSPlotNtupleBranches()
Fill the sPlot tuple.
virtual void setupBkgndVectors()
Define the length of the background vectors.
File containing LauConstants namespace.
void setSignalPdfs(LauAbsPdf *negPdf, LauAbsPdf *posPdf)
Set the signal PDFs.
void setFitNEvents()
Set the initial yields.
virtual LauSPlot::TwoDMap twodimPDFs() const
Returns the species and variables for all 2D PDFs in the fit.
virtual LauSPlot::NumbMap freeSpeciesNames() const
Returns the names and yields of species that are free in the fit.
LauEffModel * scfFracHist_
The histogram giving the DP-dependence of the SCF fraction.
void printFormat(std::ostream &stream, Double_t value) const
Method to choose the printing format to a specified level of precision.
Definition: LauPrint.cc:32
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:41
UInt_t nSigDPPar_
Number of signal DP parameters.
Bool_t haveBranch(const TString &name) const
Check if the named branch is stored.
virtual void propagateParUpdates()
Class for representing the 4D smearing matrix for mis-reconstructed signal (self cross feed) ...
Definition: LauScfMap.hh:29
virtual void storePerEvtLlhds()
Store the per event likelihood values.
Double_t scfDPLike_
SCF DP likelihood value.
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.
std::vector< Double_t > recoJacobians_
The cached values of the sqDP jacobians for each event.
Bool_t enableEmbedding() const
Determine whether embedding of events is enabled in the generation.
virtual Double_t getEventSum() const
Get the total number of events.
virtual TString name() const
Retrieve the name of the coefficient set.
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:40
const LauParameterPList & fitPars() const
Access the fit variables.
LauPdfList posScfPdfs_
The B+ SCF PDFs.
virtual const LauComplex & getEvtDPAmp() const =0
Retrieve the total amplitude of all amplitude components at the current point in the Dalitz plot...
LauBkgndYieldList bkgndEvents_
Background yield(s)
Class for calculating 3-body kinematic quantities.
Double_t value() const
The value of the parameter.
LauBkgndDPModelList negBkgndDPModels_
The B- background Dalitz plot models.
const LauParArray & getFitFractions() const
Retrieve the fit fractions for the amplitude components.
UInt_t nEvents() const
Retrieve the number of events.
LauParArray CPVFitFrac_
The CP violating fit fraction.
virtual void addGenNtupleDoubleBranch(const TString &name)
Add a branch to the gen tree for storing a double.
LauPdfList posSignalPdfs_
The B+ signal PDFs.
Double_t calcSqDPJacobian()
Calculate the Jacobian for the transformation m23^2, m13^2 -&gt; m&#39;, theta&#39; (square DP) ...
void updateFitParameters(LauPdfList &pdfList)
Update the fit parameters for all PDFs in the list.
void embedPosBkgnd(const TString &bgClass, 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...
std::vector< Double_t > fakeJacobians_
The cached values of the sqDP jacobians for each true bin.
virtual Double_t getEvtEff() const =0
Retrieve the efficiency for the current event.
Double_t getAsymmetry()
Retrieve the asymmetry.
Definition: LauAsymmCalc.hh:42
virtual Double_t getEvtSCFDPLikelihood(UInt_t iEvt)
Calculate the SCF likelihood for the DP for a given event.
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 setSCFPdfs(LauAbsPdf *negPdf, LauAbsPdf *posPdf)
Set the SCF PDF for a given variable.
std::map< std::pair< TString, Int_t >, std::pair< Int_t, Double_t > > LauGenInfo
Define a map to be used to store a category name and numbers.
Bool_t useSCFHist_
Is the SCF fraction DP-dependent.
LauParameter negMeanEff_
The mean efficiency for B- model.
void setBkgndPdfs(const TString &bkgndClass, LauAbsPdf *negPdf, LauAbsPdf *posPdf)
Set the background PDFs.
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
virtual Double_t getTotEvtLikelihood(UInt_t iEvt)
Get the total likelihood for each event.
virtual void getEvtDPLikelihood(UInt_t iEvt)
Calculate the signal and background likelihoods for the DP for a given event.
LauEffModel * getEffModel()
Retrieve the model for the efficiency across the Dalitz plot.
Double_t genValue() const
The value generated for the parameter.
LauBkgndEmbDataList posBkgndTree_
The B+ background event tree.
void setBkgndDPModels(const TString &bkgndClass, LauAbsBkgndDPModel *negModel, LauAbsBkgndDPModel *posModel)
Set the background DP models.
LauBkgndPdfsList negBkgndPdfs_
The B- background PDFs.
virtual void initialise(const std::vector< LauComplex > &coeffs)=0
Initialise the Dalitz plot dynamics.
Class to store the input fit variables.
LauParameter scfFrac_
The (global) SCF fraction parameter.
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.
void setDPBranchValues()
Store all of the DP information.
LauKinematics * posKinematics_
The B+ Dalitz plot kinematics object.
void embedPosSignal(const TString &fileName, const TString &treeName, Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment=kFALSE, Bool_t useReweighting=kFALSE)
Embed full simulation events for the B+ signal, rather than generating toy from the PDFs...
std::vector< LauParameter > fitFracAsymm_
The fit fraction asymmetries.
virtual Bool_t generate()=0
Generate a toy MC event from the model.
void setExtraNtupleVars()
Set-up other parameters that are derived from the fit results, e.g. fit fractions.
virtual Bool_t generate()=0
Generate a toy MC signal event.