laura is hosted by Hepforge, IPPP Durham
Laura++  v1r1p1
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 != nSigComp_ || nPosAmp != nSigComp_) {
569  cerr << "ERROR in LauCPFitModel::setSignalDPParameters : Number of signal DP components with magnitude and phase set not right." << endl;
570  gSystem->Exit(EXIT_FAILURE);
571  }
572 
573 
574  // Place signal model parameters in vector of fit variables
575  LauParameterPList& fitVars = this->fitPars();
576  for (UInt_t i = 0; i < nSigComp_; i++) {
577  LauParameterPList pars = coeffPars_[i]->getParameters();
578  for (LauParameterPList::iterator iter = pars.begin(); iter != pars.end(); ++iter) {
579  if ( !(*iter)->clone() ) {
580  fitVars.push_back(*iter);
581  ++nSigDPPar_;
582  }
583  }
584  }
585 }
586 
588 {
589  // Include all the parameters of the PDF in the fit
590  // NB all of them are passed to the fit, even though some have been fixed through parameter.fixed(kTRUE)
591  // With the new "cloned parameter" scheme only "original" parameters are passed to the fit.
592  // Their clones are updated automatically when the originals are updated.
593 
594  nExtraPdfPar_ = 0;
595 
597  if ( tagged_ ) {
599  }
600 
601  if (useSCF_ == kTRUE) {
603  if ( tagged_ ) {
605  }
606  }
607 
608  if (usingBkgnd_ == kTRUE) {
609  for (LauBkgndPdfsList::iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
610  nExtraPdfPar_ += this->addFitParameters(*iter);
611  }
612  if ( tagged_ ) {
613  for (LauBkgndPdfsList::iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
614  nExtraPdfPar_ += this->addFitParameters(*iter);
615  }
616  }
617  }
618 }
619 
621 {
622  if ( signalEvents_ == 0 ) {
623  cerr << "ERROR in LauCPFitModel::setFitNEvents : Signal yield not defined." << endl;
624  return;
625  }
626  nNormPar_ = 0;
627 
628  // initialise the total number of events to be the sum of all the hypotheses
629  Double_t nTotEvts = signalEvents_->value();
630  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
631  nTotEvts += (*iter)->value();
632  if ( (*iter) == 0 ) {
633  cerr << "ERROR in LauCPFitModel::setFitNEvents : Background yield not defined." << endl;
634  return;
635  }
636  }
637  this->eventsPerExpt(TMath::FloorNint(nTotEvts));
638 
639  LauParameterPList& fitVars = this->fitPars();
640 
641  // if doing an extended ML fit add the number of signal events into the fit parameters
642  if (this->doEMLFit()) {
643  cout << "INFO in LauCPFitModel::setFitNEvents : Initialising number of events for signal and background components..." << endl;
644  // add the signal fraction to the list of fit parameters
645  fitVars.push_back(signalEvents_);
646  ++nNormPar_;
647  } else {
648  cout << "INFO in LauCPFitModel::setFitNEvents : Initialising number of events for background components (and hence signal)..." << endl;
649  }
650 
651  // if not using the DP in the model we need an explicit signal asymmetry parameter
652  if (this->useDP() == kFALSE) {
653  fitVars.push_back(signalAsym_);
654  ++nNormPar_;
655  }
656 
657  if (useSCF_ && !useSCFHist_) {
658  fitVars.push_back(&scfFrac_);
659  ++nNormPar_;
660  }
661 
662  if (usingBkgnd_ == kTRUE) {
663  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
664  LauParameter* parameter = (*iter);
665  fitVars.push_back(parameter);
666  ++nNormPar_;
667  }
668  for (LauBkgndYieldList::iterator iter = bkgndAsym_.begin(); iter != bkgndAsym_.end(); ++iter) {
669  LauParameter* parameter = (*iter);
670  fitVars.push_back(parameter);
671  ++nNormPar_;
672  }
673  }
674 }
675 
677 {
678  // Set-up other parameters derived from the fit results, e.g. fit fractions.
679 
680  if (this->useDP() != kTRUE) {
681  return;
682  }
683 
684  // First clear the vectors so we start from scratch
685  this->clearExtraVarVectors();
686 
687  // Add the positive and negative fit fractions for each signal component
689  if (negFitFrac_.size() != nSigComp_) {
690  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << negFitFrac_.size() << endl;
691  gSystem->Exit(EXIT_FAILURE);
692  }
693  for (UInt_t i(0); i<nSigComp_; ++i) {
694  if (negFitFrac_[i].size() != nSigComp_) {
695  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << negFitFrac_[i].size() << endl;
696  gSystem->Exit(EXIT_FAILURE);
697  }
698  }
699 
700  LauParameterList& extraVars = this->extraPars();
701 
702  for (UInt_t i(0); i<nSigComp_; ++i) {
703  for (UInt_t j(i); j<nSigComp_; ++j) {
704  TString name = negFitFrac_[i][j].name();
705  name.Insert( name.Index("FitFrac"), "Neg" );
706  negFitFrac_[i][j].name(name);
707  extraVars.push_back(negFitFrac_[i][j]);
708  }
709  }
710 
712  if (posFitFrac_.size() != nSigComp_) {
713  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << posFitFrac_.size() << endl;
714  gSystem->Exit(EXIT_FAILURE);
715  }
716  for (UInt_t i(0); i<nSigComp_; ++i) {
717  if (posFitFrac_[i].size() != nSigComp_) {
718  cerr << "ERROR in LauCPFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: " << posFitFrac_[i].size() << endl;
719  gSystem->Exit(EXIT_FAILURE);
720  }
721  }
722  for (UInt_t i(0); i<nSigComp_; ++i) {
723  for (UInt_t j(i); j<nSigComp_; ++j) {
724  TString name = posFitFrac_[i][j].name();
725  name.Insert( name.Index("FitFrac"), "Pos" );
726  posFitFrac_[i][j].name(name);
727  extraVars.push_back(posFitFrac_[i][j]);
728  }
729  }
730 
731  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
732  std::vector<LauParameter> negExtraPars = negSigModel_->getExtraParameters();
733  std::vector<LauParameter>::iterator negExtraIter;
734  for (negExtraIter = negExtraPars.begin(); negExtraIter != negExtraPars.end(); ++negExtraIter) {
735  LauParameter negExtraParameter = (*negExtraIter);
736  extraVars.push_back(negExtraParameter);
737  }
738 
739  std::vector<LauParameter> posExtraPars = posSigModel_->getExtraParameters();
740  std::vector<LauParameter>::iterator posExtraIter;
741  for (posExtraIter = posExtraPars.begin(); posExtraIter != posExtraPars.end(); ++posExtraIter) {
742  LauParameter posExtraParameter = (*posExtraIter);
743  extraVars.push_back(posExtraParameter);
744  }
745 
746  // Now add in the DP efficiency value
747  Double_t initMeanEff = negSigModel_->getMeanEff().initValue();
748  negMeanEff_.value(initMeanEff);
749  negMeanEff_.genValue(initMeanEff);
750  negMeanEff_.initValue(initMeanEff);
751  extraVars.push_back(negMeanEff_);
752 
753  initMeanEff = posSigModel_->getMeanEff().initValue();
754  posMeanEff_.value(initMeanEff);
755  posMeanEff_.genValue(initMeanEff);
756  posMeanEff_.initValue(initMeanEff);
757  extraVars.push_back(posMeanEff_);
758 
759  // Also add in the DP rates
760  Double_t initDPRate = negSigModel_->getDPRate().initValue();
761  negDPRate_.value(initDPRate);
762  negDPRate_.genValue(initDPRate);
763  negDPRate_.initValue(initDPRate);
764  extraVars.push_back(negDPRate_);
765 
766  initDPRate = posSigModel_->getDPRate().initValue();
767  posDPRate_.value(initDPRate);
768  posDPRate_.genValue(initDPRate);
769  posDPRate_.initValue(initDPRate);
770  extraVars.push_back(posDPRate_);
771 
772  // Calculate the CPC and CPV Fit Fractions, ACPs and FitFrac asymmetries
773  this->calcExtraFractions(kTRUE);
774  this->calcAsymmetries(kTRUE);
775 
776  // Add the CP violating and CP conserving fit fractions for each signal component
777  for (UInt_t i = 0; i < nSigComp_; i++) {
778  for (UInt_t j = i; j < nSigComp_; j++) {
779  extraVars.push_back(CPVFitFrac_[i][j]);
780  }
781  }
782  for (UInt_t i = 0; i < nSigComp_; i++) {
783  for (UInt_t j = i; j < nSigComp_; j++) {
784  extraVars.push_back(CPCFitFrac_[i][j]);
785  }
786  }
787 
788  // Add the Fit Fraction asymmetry for each signal component
789  for (UInt_t i = 0; i < nSigComp_; i++) {
790  extraVars.push_back(fitFracAsymm_[i]);
791  }
792 
793  // Add the calculated CP asymmetry for each signal component
794  for (UInt_t i = 0; i < nSigComp_; i++) {
795  extraVars.push_back(acp_[i]);
796  }
797 }
798 
799 void LauCPFitModel::calcExtraFractions(Bool_t initValues)
800 {
801  // Calculate the CP-conserving and CP-violating fit fractions
802 
803  if (initValues) {
804  // create the structure
805  CPCFitFrac_.clear();
806  CPVFitFrac_.clear();
807  CPCFitFrac_.resize(nSigComp_);
808  CPVFitFrac_.resize(nSigComp_);
809  for (UInt_t i(0); i<nSigComp_; ++i) {
810  CPCFitFrac_[i].resize(nSigComp_);
811  CPVFitFrac_[i].resize(nSigComp_);
812 
813  for (UInt_t j(i); j<nSigComp_; ++j) {
814  TString name = negFitFrac_[i][j].name();
815  name.Replace( name.Index("Neg"), 3, "CPC" );
816  CPCFitFrac_[i][j].name( name );
817  CPCFitFrac_[i][j].valueAndRange( 0.0, -100.0, 100.0 );
818 
819  name = negFitFrac_[i][j].name();
820  name.Replace( name.Index("Neg"), 3, "CPV" );
821  CPVFitFrac_[i][j].name( name );
822  CPVFitFrac_[i][j].valueAndRange( 0.0, -100.0, 100.0 );
823  }
824  }
825  }
826 
827  Double_t denom = negDPRate_ + posDPRate_;
828 
829  for (UInt_t i(0); i<nSigComp_; ++i) {
830  for (UInt_t j(i); j<nSigComp_; ++j) {
831 
832  Double_t negTerm = negFitFrac_[i][j]*negDPRate_;
833  Double_t posTerm = posFitFrac_[i][j]*posDPRate_;
834 
835  Double_t cpcFitFrac = (negTerm + posTerm)/denom;
836  Double_t cpvFitFrac = (negTerm - posTerm)/denom;
837 
838  CPCFitFrac_[i][j] = cpcFitFrac;
839  CPVFitFrac_[i][j] = cpvFitFrac;
840 
841  if (initValues) {
842  CPCFitFrac_[i][j].genValue(cpcFitFrac);
843  CPCFitFrac_[i][j].initValue(cpcFitFrac);
844 
845  CPVFitFrac_[i][j].genValue(cpvFitFrac);
846  CPVFitFrac_[i][j].initValue(cpvFitFrac);
847  }
848  }
849  }
850 }
851 
852 void LauCPFitModel::calcAsymmetries(Bool_t initValues)
853 {
854  // Calculate the CP asymmetries
855  // Also calculate the fit fraction asymmetries
856 
857  for (UInt_t i = 0; i < nSigComp_; i++) {
858 
859  acp_[i] = coeffPars_[i]->acp();
860 
861  LauAsymmCalc asymmCalc(negFitFrac_[i][i].value(), posFitFrac_[i][i].value());
862  Double_t asym = asymmCalc.getAsymmetry();
863  fitFracAsymm_[i] = asym;
864  if (initValues) {
865  fitFracAsymm_[i].genValue(asym);
866  fitFracAsymm_[i].initValue(asym);
867  }
868  }
869 }
870 
871 void LauCPFitModel::finaliseFitResults(const TString& tablePrefixName)
872 {
873  // Retrieve parameters from the fit results for calculations and toy generation
874  // and eventually store these in output root ntuples/text files
875 
876  // Now take the fit parameters and update them as necessary
877  // i.e. to make mag > 0.0, phase in the right range.
878  // This function will also calculate any other values, such as the
879  // fit fractions, using any errors provided by fitParErrors as appropriate.
880  // Also obtain the pull values: (measured - generated)/(average error)
881 
882  if (this->useDP() == kTRUE) {
883  for (UInt_t i = 0; i < nSigComp_; ++i) {
884  // Check whether we have "a/b > 0.0", and phases in the right range
885  coeffPars_[i]->finaliseValues();
886  }
887  }
888 
889  // update the pulls on the event fractions and asymmetries
890  if (this->doEMLFit()) {
892  }
893  if (this->useDP() == kFALSE) {
895  }
896  if (useSCF_ && !useSCFHist_) {
898  }
899  if (usingBkgnd_ == kTRUE) {
900  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
901  (*iter)->updatePull();
902  }
903  for (LauBkgndYieldList::iterator iter = bkgndAsym_.begin(); iter != bkgndAsym_.end(); ++iter) {
904  (*iter)->updatePull();
905  }
906  }
907 
908  // Update the pulls on all the extra PDFs' parameters
911  if (useSCF_ == kTRUE) {
914  }
915  if (usingBkgnd_ == kTRUE) {
916  for (LauBkgndPdfsList::iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
917  this->updateFitParameters(*iter);
918  }
919  for (LauBkgndPdfsList::iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
920  this->updateFitParameters(*iter);
921  }
922  }
923 
924  // Fill the fit results to the ntuple
925 
926  // update the coefficients and then calculate the fit fractions and ACP's
927  if (this->useDP() == kTRUE) {
928  this->updateCoeffs();
931 
932  LauParArray negFitFrac = negSigModel_->getFitFractions();
933  if (negFitFrac.size() != nSigComp_) {
934  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << negFitFrac.size() << endl;
935  gSystem->Exit(EXIT_FAILURE);
936  }
937  for (UInt_t i(0); i<nSigComp_; ++i) {
938  if (negFitFrac[i].size() != nSigComp_) {
939  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << negFitFrac[i].size() << endl;
940  gSystem->Exit(EXIT_FAILURE);
941  }
942  }
943 
944  LauParArray posFitFrac = posSigModel_->getFitFractions();
945  if (posFitFrac.size() != nSigComp_) {
946  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << posFitFrac.size() << endl;
947  gSystem->Exit(EXIT_FAILURE);
948  }
949  for (UInt_t i(0); i<nSigComp_; ++i) {
950  if (posFitFrac[i].size() != nSigComp_) {
951  cerr << "ERROR in LauCPFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: " << posFitFrac[i].size() << endl;
952  gSystem->Exit(EXIT_FAILURE);
953  }
954  }
955 
956  for (UInt_t i(0); i<nSigComp_; ++i) {
957  for (UInt_t j(i); j<nSigComp_; ++j) {
958  negFitFrac_[i][j].value(negFitFrac[i][j].value());
959  posFitFrac_[i][j].value(posFitFrac[i][j].value());
960  }
961  }
962 
967 
968  this->calcExtraFractions();
969  this->calcAsymmetries();
970 
971  // Then store the final fit parameters, and any extra parameters for
972  // the signal model (e.g. fit fractions, FF asymmetries, ACPs, mean efficiency and DP rate)
973 
974  this->clearExtraVarVectors();
975  LauParameterList& extraVars = this->extraPars();
976 
977  // Add the positive and negative fit fractions for each signal component
978  for (UInt_t i(0); i<nSigComp_; ++i) {
979  for (UInt_t j(i); j<nSigComp_; ++j) {
980  extraVars.push_back(negFitFrac_[i][j]);
981  }
982  }
983  for (UInt_t i(0); i<nSigComp_; ++i) {
984  for (UInt_t j(i); j<nSigComp_; ++j) {
985  extraVars.push_back(posFitFrac_[i][j]);
986  }
987  }
988 
989  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
990  std::vector<LauParameter> negExtraPars = negSigModel_->getExtraParameters();
991  std::vector<LauParameter>::iterator negExtraIter;
992  for (negExtraIter = negExtraPars.begin(); negExtraIter != negExtraPars.end(); ++negExtraIter) {
993  LauParameter negExtraParameter = (*negExtraIter);
994  extraVars.push_back(negExtraParameter);
995  }
996 
997  std::vector<LauParameter> posExtraPars = posSigModel_->getExtraParameters();
998  std::vector<LauParameter>::iterator posExtraIter;
999  for (posExtraIter = posExtraPars.begin(); posExtraIter != posExtraPars.end(); ++posExtraIter) {
1000  LauParameter posExtraParameter = (*posExtraIter);
1001  extraVars.push_back(posExtraParameter);
1002  }
1003 
1004  extraVars.push_back(negMeanEff_);
1005  extraVars.push_back(posMeanEff_);
1006  extraVars.push_back(negDPRate_);
1007  extraVars.push_back(posDPRate_);
1008 
1009  for (UInt_t i = 0; i < nSigComp_; i++) {
1010  for (UInt_t j(i); j<nSigComp_; ++j) {
1011  extraVars.push_back(CPVFitFrac_[i][j]);
1012  }
1013  }
1014  for (UInt_t i = 0; i < nSigComp_; i++) {
1015  for (UInt_t j(i); j<nSigComp_; ++j) {
1016  extraVars.push_back(CPCFitFrac_[i][j]);
1017  }
1018  }
1019  for (UInt_t i(0); i<nSigComp_; ++i) {
1020  extraVars.push_back(fitFracAsymm_[i]);
1021  }
1022  for (UInt_t i(0); i<nSigComp_; ++i) {
1023  extraVars.push_back(acp_[i]);
1024  }
1025 
1026  this->printFitFractions(cout);
1027  this->printAsymmetries(cout);
1028  }
1029 
1030  const LauParameterPList& fitVars = this->fitPars();
1031  const LauParameterList& extraVars = this->extraPars();
1032  LauFitNtuple* ntuple = this->fitNtuple();
1033  ntuple->storeParsAndErrors(fitVars, extraVars);
1034 
1035  // find out the correlation matrix for the parameters
1036  ntuple->storeCorrMatrix(this->iExpt(), this->nll(), this->fitStatus());
1037 
1038  // Fill the data into ntuple
1039  ntuple->updateFitNtuple();
1040 
1041  // Print out the partial fit fractions, phases and the
1042  // averaged efficiency, reweighted by the dynamics (and anything else)
1043  if (this->writeLatexTable()) {
1044  TString sigOutFileName(tablePrefixName);
1045  sigOutFileName += "_"; sigOutFileName += this->iExpt(); sigOutFileName += "Expt.tex";
1046  this->writeOutTable(sigOutFileName);
1047  }
1048 }
1049 
1050 void LauCPFitModel::printFitFractions(std::ostream& output)
1051 {
1052  // Print out Fit Fractions, total DP rate and mean efficiency
1053  // First for the B- events
1054  for (UInt_t i = 0; i < nSigComp_; i++) {
1055  const TString compName(coeffPars_[i]->name());
1056  output << negParent_ << " FitFraction for component " << i << " (" << compName << ") = " << negFitFrac_[i][i] << endl;
1057  }
1058  output << negParent_ << " overall DP rate (integral of matrix element squared) = " << negDPRate_ << endl;
1059  output << negParent_ << " average efficiency weighted by whole DP dynamics = " << negMeanEff_ << endl;
1060 
1061  // Then for the positive sample
1062  for (UInt_t i = 0; i < nSigComp_; i++) {
1063  const TString compName(coeffPars_[i]->name());
1064  const TString conjName(negSigModel_->getConjResName(compName));
1065  output << posParent_ << " FitFraction for component " << i << " (" << conjName << ") = " << posFitFrac_[i][i] << endl;
1066  }
1067  output << posParent_ << " overall DP rate (integral of matrix element squared) = " << posDPRate_ << endl;
1068  output << posParent_ << " average efficiency weighted by whole DP dynamics = " << posMeanEff_ << endl;
1069 }
1070 
1071 void LauCPFitModel::printAsymmetries(std::ostream& output)
1072 {
1073  for (UInt_t i = 0; i < nSigComp_; i++) {
1074  const TString compName(coeffPars_[i]->name());
1075  output << "Fit Fraction asymmetry for component " << i << " (" << compName << ") = " << fitFracAsymm_[i] << endl;
1076  }
1077  for (UInt_t i = 0; i < nSigComp_; i++) {
1078  const TString compName(coeffPars_[i]->name());
1079  output << "ACP for component " << i << " (" << compName << ") = " << acp_[i] << endl;
1080  }
1081 }
1082 
1083 void LauCPFitModel::writeOutTable(const TString& outputFile)
1084 {
1085  // Write out the results of the fit to a tex-readable table
1086  // TODO - need to include the yields in this table
1087  std::ofstream fout(outputFile);
1088  LauPrint print;
1089 
1090  cout << "Writing out results of the fit to the tex file " << outputFile << endl;
1091 
1092  if (this->useDP() == kTRUE) {
1093  // print the fit coefficients in one table
1094  coeffPars_.front()->printTableHeading(fout);
1095  for (UInt_t i = 0; i < nSigComp_; i++) {
1096  coeffPars_[i]->printTableRow(fout);
1097  }
1098  fout << "\\hline" << endl;
1099  fout << "\\end{tabular}" << endl << endl;
1100 
1101  // print the fit fractions and asymmetries in another
1102  fout << "\\begin{tabular}{|l|c|c|c|c|}" << endl;
1103  fout << "\\hline" << endl;
1104  fout << "Component & " << negParent_ << " Fit Fraction & " << posParent_ << " Fit Fraction & Fit Fraction Asymmetry & ACP \\\\" << endl;
1105  fout << "\\hline" << endl;
1106 
1107  Double_t negFitFracSum(0.0);
1108  Double_t posFitFracSum(0.0);
1109 
1110  for (UInt_t i = 0; i < nSigComp_; i++) {
1111  TString resName = coeffPars_[i]->name();
1112  resName = resName.ReplaceAll("_", "\\_");
1113 
1114  Double_t negFitFrac = negFitFrac_[i][i].value();
1115  Double_t posFitFrac = posFitFrac_[i][i].value();
1116  negFitFracSum += negFitFrac;
1117  posFitFracSum += posFitFrac;
1118 
1119  Double_t fitFracAsymm = fitFracAsymm_[i].value();
1120 
1121  Double_t acp = acp_[i].value();
1122  Double_t acpErr = acp_[i].error();
1123 
1124  fout << resName << " & $";
1125 
1126  print.printFormat(fout, negFitFrac);
1127  fout << "$ & $";
1128  print.printFormat(fout, posFitFrac);
1129  fout << "$ & $";
1130 
1131  print.printFormat(fout, fitFracAsymm);
1132  fout << "$ & $";
1133 
1134  print.printFormat(fout, acp);
1135  fout << " \\pm ";
1136  print.printFormat(fout, acpErr);
1137  fout << "$ \\\\" << endl;
1138  }
1139  fout << "\\hline" << endl;
1140 
1141  // Also print out sum of fit fractions
1142  fout << "Fit Fraction Sum & $";
1143  print.printFormat(fout, negFitFracSum);
1144  fout << "$ & $";
1145  print.printFormat(fout, posFitFracSum);
1146  fout << "$ & & \\\\" << endl;
1147  fout << "\\hline" << endl;
1148 
1149  fout << "DP rate & $";
1150  print.printFormat(fout, negDPRate_.value());
1151  fout << "$ & $";
1152  print.printFormat(fout, posDPRate_.value());
1153  fout << "$ & & \\\\" << endl;
1154 
1155  fout << "$< \\varepsilon > $ & $";
1156  print.printFormat(fout, negMeanEff_.value());
1157  fout << "$ & $";
1158  print.printFormat(fout, posMeanEff_.value());
1159  fout << "$ & & \\\\" << endl;
1160  fout << "\\hline" << endl;
1161  fout << "\\end{tabular}" << endl << endl;
1162  }
1163 
1164  if (!negSignalPdfs_.empty()) {
1165  fout << "\\begin{tabular}{|l|c|}" << endl;
1166  fout << "\\hline" << endl;
1167  if (useSCF_ == kTRUE) {
1168  fout << "\\Extra TM Signal PDFs' Parameters: & \\\\" << endl;
1169  } else {
1170  fout << "\\Extra Signal PDFs' Parameters: & \\\\" << endl;
1171  }
1172  this->printFitParameters(negSignalPdfs_, fout);
1173  if ( tagged_ ) {
1174  this->printFitParameters(posSignalPdfs_, fout);
1175  }
1176  if (useSCF_ == kTRUE && !negScfPdfs_.empty()) {
1177  fout << "\\hline" << endl;
1178  fout << "\\Extra SCF Signal PDFs' Parameters: & \\\\" << endl;
1179  this->printFitParameters(negScfPdfs_, fout);
1180  if ( tagged_ ) {
1181  this->printFitParameters(posScfPdfs_, fout);
1182  }
1183  }
1184  if (usingBkgnd_ == kTRUE && !negBkgndPdfs_.empty()) {
1185  fout << "\\hline" << endl;
1186  fout << "\\Extra Background PDFs' Parameters: & \\\\" << endl;
1187  for (LauBkgndPdfsList::const_iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
1188  this->printFitParameters(*iter, fout);
1189  }
1190  if ( tagged_ ) {
1191  for (LauBkgndPdfsList::const_iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
1192  this->printFitParameters(*iter, fout);
1193  }
1194  }
1195  }
1196  fout << "\\hline \n\\end{tabular}" << endl << endl;
1197  }
1198 }
1199 
1201 {
1202  // Update the number of signal events to be total-sum(background events)
1203  this->updateSigEvents();
1204 
1205  // Check whether we want to have randomised initial fit parameters for the signal model
1206  if (this->useRandomInitFitPars() == kTRUE) {
1207  cout << "Setting random parameters for the signal model" << endl;
1208  this->randomiseInitFitPars();
1209  }
1210 }
1211 
1213 {
1214  // Only randomise those parameters that are not fixed!
1215  cout << "Randomising the initial fit magnitudes and phases of the components..." << endl;
1216 
1217  for (UInt_t i = 0; i < nSigComp_; i++) {
1218  coeffPars_[i]->randomiseInitValues();
1219  }
1220 }
1221 
1223 {
1224  // Determine the number of events to generate for each hypothesis
1225  // If we're smearing then smear each one individually
1226 
1227  LauGenInfo nEvtsGen;
1228 
1229  // Signal
1230  Double_t evtWeight(1.0);
1231  Double_t nEvts = signalEvents_->genValue();
1232  if ( nEvts < 0.0 ) {
1233  evtWeight = -1.0;
1234  nEvts = TMath::Abs( nEvts );
1235  }
1236  Double_t asym(0.0);
1237  Double_t sigAsym(0.0);
1238  // need to include this as an alternative in case the DP isn't in the model
1239  if ( !this->useDP() || forceAsym_ ) {
1240  sigAsym = signalAsym_->genValue();
1241  } else {
1242  Double_t negRate = negSigModel_->getDPNorm();
1243  Double_t posRate = posSigModel_->getDPNorm();
1244  if (negRate+posRate>1e-30) {
1245  sigAsym = (negRate-posRate)/(negRate+posRate);
1246  }
1247  }
1248  asym = sigAsym;
1249  Int_t nPosEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 - asym)) + 0.5);
1250  Int_t nNegEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 + asym)) + 0.5);
1251  if (this->doPoissonSmearing()) {
1252  nNegEvts = LauRandom::randomFun()->Poisson(nNegEvts);
1253  nPosEvts = LauRandom::randomFun()->Poisson(nPosEvts);
1254  }
1255  nEvtsGen[std::make_pair("signal",-1)] = std::make_pair(nNegEvts,evtWeight);
1256  nEvtsGen[std::make_pair("signal",+1)] = std::make_pair(nPosEvts,evtWeight);
1257 
1258  // backgrounds
1259  const UInt_t nBkgnds = this->nBkgndClasses();
1260  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1261  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1262  const LauParameter* evtsPar = bkgndEvents_[bkgndID];
1263  const LauParameter* asymPar = bkgndAsym_[bkgndID];
1264  evtWeight = 1.0;
1265  nEvts = TMath::FloorNint( evtsPar->genValue() );
1266  if ( nEvts < 0 ) {
1267  evtWeight = -1.0;
1268  nEvts = TMath::Abs( nEvts );
1269  }
1270  asym = asymPar->genValue();
1271  nPosEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 - asym)) + 0.5);
1272  nNegEvts = static_cast<Int_t>((nEvts/2.0 * (1.0 + asym)) + 0.5);
1273  if (this->doPoissonSmearing()) {
1274  nNegEvts = LauRandom::randomFun()->Poisson(nNegEvts);
1275  nPosEvts = LauRandom::randomFun()->Poisson(nPosEvts);
1276  }
1277  nEvtsGen[std::make_pair(bkgndClass,-1)] = std::make_pair(nNegEvts,evtWeight);
1278  nEvtsGen[std::make_pair(bkgndClass,+1)] = std::make_pair(nPosEvts,evtWeight);
1279  }
1280 
1281  cout << "Generating toy MC with:" << endl;
1282  cout << "Signal asymmetry = " << sigAsym << " and number of signal events = " << signalEvents_->genValue() << endl;
1283  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1284  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1285  const LauParameter* evtsPar = bkgndEvents_[bkgndID];
1286  const LauParameter* asymPar = bkgndAsym_[bkgndID];
1287  cout << bkgndClass << " asymmetry = " << asymPar->genValue() << " and number of " << bkgndClass << " events = " << evtsPar->genValue() << endl;
1288  }
1289 
1290  return nEvtsGen;
1291 }
1292 
1294 {
1295  // Routine to generate toy Monte Carlo events according to the various models we have defined.
1296 
1297  // Determine the number of events to generate for each hypothesis
1298  LauGenInfo nEvts = this->eventsToGenerate();
1299 
1300  Bool_t genOK(kTRUE);
1301  Int_t evtNum(0);
1302 
1303  const UInt_t nBkgnds = this->nBkgndClasses();
1304  std::vector<TString> bkgndClassNames(nBkgnds);
1305  std::vector<TString> bkgndClassNamesGen(nBkgnds);
1306  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1307  TString name( this->bkgndClassName(iBkgnd) );
1308  bkgndClassNames[iBkgnd] = name;
1309  bkgndClassNamesGen[iBkgnd] = "gen"+name;
1310  }
1311 
1312  const Bool_t storeSCFTruthInfo = ( useSCF_ || ( this->enableEmbedding() &&
1313  negSignalTree_ != 0 && negSignalTree_->haveBranch("mcMatch") &&
1314  posSignalTree_ != 0 && posSignalTree_->haveBranch("mcMatch") ) );
1315 
1316  // Loop over the hypotheses and generate the requested number of events for each one
1317  for (LauGenInfo::const_iterator iter = nEvts.begin(); iter != nEvts.end(); ++iter) {
1318 
1319  const TString& type(iter->first.first);
1320  curEvtCharge_ = iter->first.second;
1321  Double_t evtWeight( iter->second.second );
1322  Int_t nEvtsGen( iter->second.first );
1323 
1324  for (Int_t iEvt(0); iEvt<nEvtsGen; ++iEvt) {
1325 
1326  this->setGenNtupleDoubleBranchValue( "evtWeight", evtWeight );
1327 
1328  if (type == "signal") {
1329  this->setGenNtupleIntegerBranchValue("genSig",1);
1330  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1331  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], 0 );
1332  }
1333  genOK = this->generateSignalEvent();
1334  } else {
1335  this->setGenNtupleIntegerBranchValue("genSig",0);
1336  if ( storeSCFTruthInfo ) {
1337  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1338  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1339  }
1340  UInt_t bkgndID(0);
1341  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1342  Int_t gen(0);
1343  if ( bkgndClassNames[iBkgnd] == type ) {
1344  gen = 1;
1345  bkgndID = iBkgnd;
1346  }
1347  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], gen );
1348  }
1349  genOK = this->generateBkgndEvent(bkgndID);
1350  }
1351  if (!genOK) {
1352  // If there was a problem with the generation then break out and return.
1353  // The problem model will have adjusted itself so that all should be OK next time.
1354  break;
1355  }
1356  if (this->useDP() == kTRUE) {
1357  this->setDPBranchValues();
1358  }
1359 
1360  // Store the event charge
1362 
1363  // Store the event number (within this experiment)
1364  // and then increment it
1365  this->setGenNtupleIntegerBranchValue("iEvtWithinExpt",evtNum);
1366  ++evtNum;
1367 
1368  this->fillGenNtupleBranches();
1369  if (iEvt%500 == 0) {cout << "Generated event number " << iEvt << " out of " << nEvtsGen << " " << type << " events." << endl;}
1370  }
1371 
1372  if (!genOK) {
1373  break;
1374  }
1375  }
1376 
1377  if (this->useDP() && genOK) {
1378 
1379  negSigModel_->checkToyMC(kTRUE,kTRUE);
1380  posSigModel_->checkToyMC(kTRUE,kTRUE);
1381 
1382  // Get the fit fractions if they're to be written into the latex table
1383  if (this->writeLatexTable()) {
1384  LauParArray negFitFrac = negSigModel_->getFitFractions();
1385  if (negFitFrac.size() != nSigComp_) {
1386  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << negFitFrac.size() << endl;
1387  gSystem->Exit(EXIT_FAILURE);
1388  }
1389  for (UInt_t i(0); i<nSigComp_; ++i) {
1390  if (negFitFrac[i].size() != nSigComp_) {
1391  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << negFitFrac[i].size() << endl;
1392  gSystem->Exit(EXIT_FAILURE);
1393  }
1394  }
1395  LauParArray posFitFrac = posSigModel_->getFitFractions();
1396  if (posFitFrac.size() != nSigComp_) {
1397  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << posFitFrac.size() << endl;
1398  gSystem->Exit(EXIT_FAILURE);
1399  }
1400  for (UInt_t i(0); i<nSigComp_; ++i) {
1401  if (posFitFrac[i].size() != nSigComp_) {
1402  cerr << "ERROR in LauCPFitModel::genExpt : Fit Fraction array of unexpected dimension: " << posFitFrac[i].size() << endl;
1403  gSystem->Exit(EXIT_FAILURE);
1404  }
1405  }
1406 
1407  for (UInt_t i(0); i<nSigComp_; ++i) {
1408  for (UInt_t j(i); j<nSigComp_; ++j) {
1409  negFitFrac_[i][j].value(negFitFrac[i][j].value());
1410  posFitFrac_[i][j].value(posFitFrac[i][j].value());
1411  }
1412  }
1417  }
1418  }
1419 
1420  // If we're reusing embedded events or if the generation is being
1421  // reset then clear the lists of used events
1422  if (reuseSignal_ || !genOK) {
1423  if (negSignalTree_) {
1425  }
1426  if (posSignalTree_) {
1428  }
1429  }
1430  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1431  LauEmbeddedData* data = negBkgndTree_[bkgndID];
1432  if (reuseBkgnd_[bkgndID] || !genOK) {
1433  if (data) {
1434  data->clearUsedList();
1435  }
1436  }
1437  }
1438  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
1439  LauEmbeddedData* data = posBkgndTree_[bkgndID];
1440  if (reuseBkgnd_[bkgndID] || !genOK) {
1441  if (data) {
1442  data->clearUsedList();
1443  }
1444  }
1445  }
1446 
1447  return genOK;
1448 }
1449 
1451 {
1452  // Generate signal event
1453  Bool_t genOK(kTRUE);
1454  Bool_t genSCF(kFALSE);
1455 
1456  LauAbsDPDynamics* model(0);
1457  LauKinematics* kinematics(0);
1458  LauEmbeddedData* embeddedData(0);
1459  LauPdfList* sigPdfs(0);
1460  LauPdfList* scfPdfs(0);
1461 
1462  Bool_t doReweighting(kFALSE);
1463 
1464  if (curEvtCharge_<0) {
1465  model = negSigModel_;
1466  kinematics = negKinematics_;
1467  sigPdfs = &negSignalPdfs_;
1468  scfPdfs = &negScfPdfs_;
1469  if (this->enableEmbedding()) {
1470  embeddedData = negSignalTree_;
1471  doReweighting = useNegReweighting_;
1472  }
1473  } else {
1474  model = posSigModel_;
1475  kinematics = posKinematics_;
1476  if ( tagged_ ) {
1477  sigPdfs = &posSignalPdfs_;
1478  scfPdfs = &posScfPdfs_;
1479  } else {
1480  sigPdfs = &negSignalPdfs_;
1481  scfPdfs = &negScfPdfs_;
1482  }
1483  if (this->enableEmbedding()) {
1484  embeddedData = posSignalTree_;
1485  doReweighting = usePosReweighting_;
1486  }
1487  }
1488 
1489  if (this->useDP()) {
1490  if (embeddedData) {
1491  if (doReweighting) {
1492  // Select a (random) event from the generated data. Then store the
1493  // reconstructed DP co-ords, together with other pdf information,
1494  // as the embedded data.
1495  genOK = embeddedData->getReweightedEvent(model);
1496  } else {
1497  // Just get the information of a (randomly) selected event in the
1498  // embedded data
1499  embeddedData->getEmbeddedEvent(kinematics);
1500  }
1501  genSCF = this->storeSignalMCMatch( embeddedData );
1502  } else {
1503  genOK = model->generate();
1504  if ( genOK && useSCF_ ) {
1505  Double_t frac(0.0);
1506  if ( useSCFHist_ ) {
1507  frac = scfFracHist_->calcEfficiency( kinematics );
1508  } else {
1509  frac = scfFrac_.genValue();
1510  }
1511  if ( frac < LauRandom::randomFun()->Rndm() ) {
1512  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1513  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1514  genSCF = kFALSE;
1515  } else {
1516  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1517  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1518  genSCF = kTRUE;
1519 
1520  // Optionally smear the DP position
1521  // of the SCF event
1522  if ( scfMap_ != 0 ) {
1523  Double_t xCoord(0.0), yCoord(0.0);
1524  if ( kinematics->squareDP() ) {
1525  xCoord = kinematics->getmPrime();
1526  yCoord = kinematics->getThetaPrime();
1527  } else {
1528  xCoord = kinematics->getm13Sq();
1529  yCoord = kinematics->getm23Sq();
1530  }
1531 
1532  // Find the bin number where this event is generated
1533  Int_t binNo = scfMap_->binNumber( xCoord, yCoord );
1534 
1535  // Retrieve the migration histogram
1536  TH2* histo = scfMap_->trueHist( binNo );
1537 
1538  LauEffModel * effModel = model->getEffModel();
1539  do {
1540  // Get a random point from the histogram
1541  histo->GetRandom2( xCoord, yCoord );
1542 
1543  // Update the kinematics
1544  if ( kinematics->squareDP() ) {
1545  kinematics->updateSqDPKinematics( xCoord, yCoord );
1546  } else {
1547  kinematics->updateKinematics( xCoord, yCoord );
1548  }
1549  } while ( ! effModel->passVeto( kinematics ) );
1550  }
1551  }
1552  }
1553  }
1554  } else {
1555  if (embeddedData) {
1556  embeddedData->getEmbeddedEvent(0);
1557  genSCF = this->storeSignalMCMatch( embeddedData );
1558  } else if ( useSCF_ ) {
1559  Double_t frac = scfFrac_.genValue();
1560  if ( frac < LauRandom::randomFun()->Rndm() ) {
1561  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1562  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1563  genSCF = kFALSE;
1564  } else {
1565  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1566  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1567  genSCF = kTRUE;
1568  }
1569  }
1570  }
1571  if (genOK) {
1572  if ( useSCF_ ) {
1573  if ( genSCF ) {
1574  this->generateExtraPdfValues(scfPdfs, embeddedData);
1575  } else {
1576  this->generateExtraPdfValues(sigPdfs, embeddedData);
1577  }
1578  } else {
1579  this->generateExtraPdfValues(sigPdfs, embeddedData);
1580  }
1581  }
1582  // Check for problems with the embedding
1583  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1584  cerr << "WARNING in LauCPFitModel::generateSignalEvent : Source of embedded signal events used up, clearing the list of used events." << endl;
1585  embeddedData->clearUsedList();
1586  }
1587 
1588  return genOK;
1589 }
1590 
1592 {
1593  // Generate Bkgnd event
1594  Bool_t genOK(kTRUE);
1595 
1596  LauAbsBkgndDPModel* model(0);
1597  LauEmbeddedData* embeddedData(0);
1598  LauPdfList* extraPdfs(0);
1599  LauKinematics* kinematics(0);
1600 
1601  if (curEvtCharge_<0) {
1602  model = negBkgndDPModels_[bkgndID];
1603  if (this->enableEmbedding()) {
1604  embeddedData = negBkgndTree_[bkgndID];
1605  }
1606  extraPdfs = &negBkgndPdfs_[bkgndID];
1607  kinematics = negKinematics_;
1608  } else {
1609  model = posBkgndDPModels_[bkgndID];
1610  if (this->enableEmbedding()) {
1611  embeddedData = posBkgndTree_[bkgndID];
1612  }
1613  if ( tagged_ ) {
1614  extraPdfs = &posBkgndPdfs_[bkgndID];
1615  } else {
1616  extraPdfs = &negBkgndPdfs_[bkgndID];
1617  }
1618  kinematics = posKinematics_;
1619  }
1620 
1621  if (this->useDP()) {
1622  if (embeddedData) {
1623  embeddedData->getEmbeddedEvent(kinematics);
1624  } else {
1625  if (model == 0) {
1626  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1627  cerr << "ERROR in LauCPFitModel::generateBkgndEvent : Can't find the DP model for background class \"" << bkgndClass << "\"." << endl;
1628  gSystem->Exit(EXIT_FAILURE);
1629  }
1630  genOK = model->generate();
1631  }
1632  } else {
1633  if (embeddedData) {
1634  embeddedData->getEmbeddedEvent(0);
1635  }
1636  }
1637  if (genOK) {
1638  this->generateExtraPdfValues(extraPdfs, embeddedData);
1639  }
1640 
1641  // Check for problems with the embedding
1642  if (embeddedData && (embeddedData->nEvents() == embeddedData->nUsedEvents())) {
1643  const TString& bkgndClass = this->bkgndClassName(bkgndID);
1644  cerr << "WARNING in LauCPFitModel::generateBkgndEvent : Source of embedded " << bkgndClass << " events used up, clearing the list of used events." << endl;
1645  embeddedData->clearUsedList();
1646  }
1647 
1648  return genOK;
1649 }
1650 
1652 {
1653  // Setup the required ntuple branches
1654  this->addGenNtupleDoubleBranch("evtWeight");
1655  this->addGenNtupleIntegerBranch("genSig");
1656  if ( useSCF_ || ( this->enableEmbedding() &&
1657  negSignalTree_ != 0 && negSignalTree_->haveBranch("mcMatch") &&
1658  posSignalTree_ != 0 && posSignalTree_->haveBranch("mcMatch") ) ) {
1659  this->addGenNtupleIntegerBranch("genTMSig");
1660  this->addGenNtupleIntegerBranch("genSCFSig");
1661  }
1662  const UInt_t nBkgnds = this->nBkgndClasses();
1663  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
1664  TString name( this->bkgndClassName(iBkgnd) );
1665  name.Prepend("gen");
1666  this->addGenNtupleIntegerBranch(name);
1667  }
1668  this->addGenNtupleIntegerBranch("charge");
1669  if (this->useDP() == kTRUE) {
1670  this->addGenNtupleDoubleBranch("m12");
1671  this->addGenNtupleDoubleBranch("m23");
1672  this->addGenNtupleDoubleBranch("m13");
1673  this->addGenNtupleDoubleBranch("m12Sq");
1674  this->addGenNtupleDoubleBranch("m23Sq");
1675  this->addGenNtupleDoubleBranch("m13Sq");
1676  this->addGenNtupleDoubleBranch("cosHel12");
1677  this->addGenNtupleDoubleBranch("cosHel23");
1678  this->addGenNtupleDoubleBranch("cosHel13");
1680  this->addGenNtupleDoubleBranch("mPrime");
1681  this->addGenNtupleDoubleBranch("thPrime");
1682  }
1683  }
1684  for (LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter) {
1685  std::vector<TString> varNames = (*pdf_iter)->varNames();
1686  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
1687  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
1688  this->addGenNtupleDoubleBranch( (*var_iter) );
1689  }
1690  }
1691  }
1692 }
1693 
1695 {
1696  LauKinematics* kinematics(0);
1697  if (curEvtCharge_<0) {
1698  kinematics = negKinematics_;
1699  } else {
1700  kinematics = posKinematics_;
1701  }
1702 
1703  // Store all the DP information
1704  this->setGenNtupleDoubleBranchValue("m12", kinematics->getm12());
1705  this->setGenNtupleDoubleBranchValue("m23", kinematics->getm23());
1706  this->setGenNtupleDoubleBranchValue("m13", kinematics->getm13());
1707  this->setGenNtupleDoubleBranchValue("m12Sq", kinematics->getm12Sq());
1708  this->setGenNtupleDoubleBranchValue("m23Sq", kinematics->getm23Sq());
1709  this->setGenNtupleDoubleBranchValue("m13Sq", kinematics->getm13Sq());
1710  this->setGenNtupleDoubleBranchValue("cosHel12", kinematics->getc12());
1711  this->setGenNtupleDoubleBranchValue("cosHel23", kinematics->getc23());
1712  this->setGenNtupleDoubleBranchValue("cosHel13", kinematics->getc13());
1713  if (kinematics->squareDP()) {
1714  this->setGenNtupleDoubleBranchValue("mPrime", kinematics->getmPrime());
1715  this->setGenNtupleDoubleBranchValue("thPrime", kinematics->getThetaPrime());
1716  }
1717 }
1718 
1720 {
1721  LauKinematics* kinematics(0);
1722  if (curEvtCharge_<0) {
1723  kinematics = negKinematics_;
1724  } else {
1725  kinematics = posKinematics_;
1726  }
1727 
1728  if (!extraPdfs) {
1729  std::cerr << "ERROR in LauCPFitModel::generateExtraPdfValues : Null pointer to PDF list." << std::endl;
1730  gSystem->Exit(EXIT_FAILURE);
1731  }
1732 
1733  if (extraPdfs->empty()) {
1734  //std::cerr << "WARNING in LauCPFitModel::generateExtraPdfValues : PDF list is empty." << std::endl;
1735  return;
1736  }
1737 
1738  // Generate from the extra PDFs
1739  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
1740  LauFitData genValues;
1741  if (embeddedData) {
1742  genValues = embeddedData->getValues( (*pdf_iter)->varNames() );
1743  } else {
1744  genValues = (*pdf_iter)->generate(kinematics);
1745  }
1746  for ( LauFitData::const_iterator var_iter = genValues.begin(); var_iter != genValues.end(); ++var_iter ) {
1747  TString varName = var_iter->first;
1748  if ( varName != "m13Sq" && varName != "m23Sq" ) {
1749  Double_t value = var_iter->second;
1750  this->setGenNtupleDoubleBranchValue(varName,value);
1751  }
1752  }
1753  }
1754 }
1755 
1757 {
1758  // Default to TM
1759  Bool_t genSCF(kFALSE);
1760  Int_t match(1);
1761 
1762  // Check that we have a valid pointer and that embedded data has
1763  // the mcMatch branch. If so then get the match value.
1764  if ( embeddedData && embeddedData->haveBranch("mcMatch") ) {
1765  match = TMath::Nint( embeddedData->getValue("mcMatch") );
1766  }
1767 
1768  // Set the variables accordingly.
1769  if (match) {
1770  this->setGenNtupleIntegerBranchValue("genTMSig",1);
1771  this->setGenNtupleIntegerBranchValue("genSCFSig",0);
1772  genSCF = kFALSE;
1773  } else {
1774  this->setGenNtupleIntegerBranchValue("genTMSig",0);
1775  this->setGenNtupleIntegerBranchValue("genSCFSig",1);
1776  genSCF = kTRUE;
1777  }
1778 
1779  return genSCF;
1780 }
1781 
1783 {
1784  // Update the signal parameters and then the total normalisation for the signal likelihood
1785  if (this->useDP() == kTRUE) {
1786  this->updateCoeffs();
1789  }
1790 
1791  // Update the signal fraction from the background fractions if not doing an extended fit
1792  if ( !this->doEMLFit() && !signalEvents_->fixed() ) {
1793  this->updateSigEvents();
1794  }
1795 }
1796 
1798 {
1799  // The background parameters will have been set from Minuit.
1800  // We need to update the signal events using these.
1801  Double_t nTotEvts = this->eventsPerExpt();
1802 
1803  signalEvents_->range(-2.0*nTotEvts,2.0*nTotEvts);
1804  for (LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1805  (*iter)->range(-2.0*nTotEvts,2.0*nTotEvts);
1806  }
1807 
1808  if (signalEvents_->fixed()) {
1809  return;
1810  }
1811 
1812  // Subtract background events (if any) from signal.
1813  Double_t signalEvents = nTotEvts;
1814  if (usingBkgnd_ == kTRUE) {
1815  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
1816  signalEvents -= (*iter)->value();
1817  }
1818  }
1819  signalEvents_->value(signalEvents);
1820 }
1821 
1823 {
1824  // Fill the internal data trees of the signal and background models.
1825  // Note that we store the events of both charges in both the
1826  // negative and the positive models. It's only later, at the stage
1827  // when the likelihood is being calculated, that we separate them.
1828 
1829  LauFitDataTree* inputFitData = this->fitData();
1830 
1831  // First the Dalitz plot variables (m_ij^2)
1832  if (this->useDP() == kTRUE) {
1833 
1834  // need to append SCF smearing bins before caching DP amplitudes
1835  if ( scfMap_ != 0 ) {
1836  this->appendBinCentres( inputFitData );
1837  }
1838 
1839  negSigModel_->fillDataTree(*inputFitData);
1840  posSigModel_->fillDataTree(*inputFitData);
1841 
1842  if (usingBkgnd_ == kTRUE) {
1843  for (LauBkgndDPModelList::iterator iter = negBkgndDPModels_.begin(); iter != negBkgndDPModels_.end(); ++iter) {
1844  (*iter)->fillDataTree(*inputFitData);
1845  }
1846  for (LauBkgndDPModelList::iterator iter = posBkgndDPModels_.begin(); iter != posBkgndDPModels_.end(); ++iter) {
1847  (*iter)->fillDataTree(*inputFitData);
1848  }
1849  }
1850  }
1851 
1852  // ...and then the extra PDFs
1853  this->cacheInfo(negSignalPdfs_, *inputFitData);
1854  this->cacheInfo(negScfPdfs_, *inputFitData);
1855  for (LauBkgndPdfsList::iterator iter = negBkgndPdfs_.begin(); iter != negBkgndPdfs_.end(); ++iter) {
1856  this->cacheInfo((*iter), *inputFitData);
1857  }
1858  if ( tagged_ ) {
1859  this->cacheInfo(posSignalPdfs_, *inputFitData);
1860  this->cacheInfo(posScfPdfs_, *inputFitData);
1861  for (LauBkgndPdfsList::iterator iter = posBkgndPdfs_.begin(); iter != posBkgndPdfs_.end(); ++iter) {
1862  this->cacheInfo((*iter), *inputFitData);
1863  }
1864  }
1865 
1866  // the SCF fractions and jacobians
1867  if ( useSCF_ && useSCFHist_ ) {
1868  if ( !inputFitData->haveBranch( "m13Sq" ) || !inputFitData->haveBranch( "m23Sq" ) ) {
1869  cerr << "ERROR in LauCPFitModel::cacheInputFitVars : Input data does not contain DP branches and so can't cache the SCF fraction." << endl;
1870  gSystem->Exit(EXIT_FAILURE);
1871  }
1872 
1873  UInt_t nEvents = inputFitData->nEvents();
1874  recoSCFFracs_.clear();
1875  recoSCFFracs_.reserve( nEvents );
1876  if ( negKinematics_->squareDP() ) {
1877  recoJacobians_.clear();
1878  recoJacobians_.reserve( nEvents );
1879  }
1880  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1881  const LauFitData& dataValues = inputFitData->getData(iEvt);
1882  LauFitData::const_iterator m13_iter = dataValues.find("m13Sq");
1883  LauFitData::const_iterator m23_iter = dataValues.find("m23Sq");
1884  negKinematics_->updateKinematics( m13_iter->second, m23_iter->second );
1885  Double_t scfFrac = scfFracHist_->calcEfficiency( negKinematics_ );
1886  recoSCFFracs_.push_back( scfFrac );
1887  if ( negKinematics_->squareDP() ) {
1889  }
1890  }
1891  }
1892 
1893  // finally cache the event charge
1894  evtCharges_.clear();
1895  if ( tagged_ ) {
1896  if ( !inputFitData->haveBranch( tagVarName_ ) ) {
1897  cerr << "ERROR in LauCPFitModel::cacheInputFitVars : Input data does not contain branch \"" << tagVarName_ << "\"." << endl;
1898  gSystem->Exit(EXIT_FAILURE);
1899  }
1900  UInt_t nEvents = inputFitData->nEvents();
1901  evtCharges_.reserve( nEvents );
1902  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
1903  const LauFitData& dataValues = inputFitData->getData(iEvt);
1904  LauFitData::const_iterator iter = dataValues.find( tagVarName_ );
1905  curEvtCharge_ = static_cast<Int_t>( iter->second );
1906  evtCharges_.push_back( curEvtCharge_ );
1907  }
1908  }
1909 }
1910 
1912 {
1913  // We'll be caching the DP amplitudes and efficiencies of the centres of the true bins.
1914  // To do so, we attach some fake points at the end of inputData, the number of the entry
1915  // minus the total number of events corresponding to the number of the histogram for that
1916  // given true bin in the LauScfMap object. (What this means is that when Laura is provided with
1917  // the LauScfMap object by the user, it's the latter who has to make sure that it contains the
1918  // right number of histograms and in exactly the right order!)
1919 
1920  // Get the x and y co-ordinates of the bin centres
1921  std::vector<Double_t> binCentresXCoords;
1922  std::vector<Double_t> binCentresYCoords;
1923  scfMap_->listBinCentres(binCentresXCoords, binCentresYCoords);
1924 
1925  // The SCF histograms could be in square Dalitz plot histograms.
1926  // The dynamics takes normal Dalitz plot coords, so we might have to convert them back.
1927  Bool_t sqDP = negKinematics_->squareDP();
1928  UInt_t nBins = binCentresXCoords.size();
1929  fakeSCFFracs_.clear();
1930  fakeSCFFracs_.reserve( nBins );
1931  if ( sqDP ) {
1932  fakeJacobians_.clear();
1933  fakeJacobians_.reserve( nBins );
1934  }
1935 
1936  for (UInt_t iBin = 0; iBin < nBins; ++iBin) {
1937 
1938  if ( sqDP ) {
1939  negKinematics_->updateSqDPKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1940  binCentresXCoords[iBin] = negKinematics_->getm13Sq();
1941  binCentresYCoords[iBin] = negKinematics_->getm23Sq();
1943  } else {
1944  negKinematics_->updateKinematics(binCentresXCoords[iBin],binCentresYCoords[iBin]);
1945  }
1946 
1948  }
1949 
1950  // Set up inputFitVars_ object to hold the fake events
1951  inputData->appendFakePoints(binCentresXCoords,binCentresYCoords);
1952 }
1953 
1955 {
1956  // Find out whether we have B- or B+
1957  if ( tagged_ ) {
1958  curEvtCharge_ = evtCharges_[iEvt];
1959 
1960  // check that the charge is either +1 or -1
1961  if (TMath::Abs(curEvtCharge_)!=1) {
1962  cerr << "ERROR in LauCPFitModel::getTotEvtLikelihood : Charge/tag not accepted value: " << curEvtCharge_ << endl;
1963  if (curEvtCharge_>0) {
1964  curEvtCharge_ = +1;
1965  } else {
1966  curEvtCharge_ = -1;
1967  }
1968  cerr << " : Making it: " << curEvtCharge_ << "." << endl;
1969  }
1970  }
1971 
1972  // Get the DP likelihood for signal and backgrounds
1973  this->getEvtDPLikelihood(iEvt);
1974 
1975  // Get the combined extra PDFs likelihood for signal and backgrounds
1976  this->getEvtExtraLikelihoods(iEvt);
1977 
1978  // If appropriate, combine the TM and SCF likelihoods
1979  Double_t sigLike = sigDPLike_ * sigExtraLike_;
1980  if ( useSCF_ ) {
1981  Double_t scfFrac(0.0);
1982  if (useSCFHist_) {
1983  scfFrac = recoSCFFracs_[iEvt];
1984  } else {
1985  scfFrac = scfFrac_.value();
1986  }
1987  sigLike *= (1.0 - scfFrac);
1988  if ( (scfMap_ != 0) && (this->useDP() == kTRUE) ) {
1989  // if we're smearing the SCF DP PDF then the SCF frac
1990  // is already included in the SCF DP likelihood
1991  sigLike += (scfDPLike_ * scfExtraLike_);
1992  } else {
1993  sigLike += (scfFrac * scfDPLike_ * scfExtraLike_);
1994  }
1995  }
1996 
1997  // Get the correct event fractions depending on the charge
1998  // Signal asymmetry is built into the DP model... but when the DP
1999  // isn't in the fit we need an explicit parameter
2000  Double_t signalEvents = signalEvents_->value() * 0.5;
2001  if (this->useDP() == kFALSE) {
2002  signalEvents *= (1.0 - curEvtCharge_ * signalAsym_->value());
2003  }
2004 
2005  // Construct the total event likelihood
2006  Double_t likelihood(0.0);
2007  if (usingBkgnd_) {
2008  likelihood = sigLike*signalEvents;
2009  const UInt_t nBkgnds = this->nBkgndClasses();
2010  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2011  Double_t bkgndEvents = bkgndEvents_[bkgndID]->value() * 0.5 * (1.0 - curEvtCharge_ * bkgndAsym_[bkgndID]->value());
2012  likelihood += bkgndEvents*bkgndDPLike_[bkgndID]*bkgndExtraLike_[bkgndID];
2013  }
2014  } else {
2015  likelihood = sigLike*0.5;
2016  }
2017  return likelihood;
2018 }
2019 
2021 {
2022  Double_t eventSum(0.0);
2023  eventSum += signalEvents_->value();
2024  if (usingBkgnd_) {
2025  for (LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end(); ++iter) {
2026  eventSum += (*iter)->value();
2027  }
2028  }
2029  return eventSum;
2030 }
2031 
2033 {
2034  // Function to return the signal and background likelihoods for the
2035  // Dalitz plot for the given event evtNo.
2036 
2037  if ( ! this->useDP() ) {
2038  // There's always going to be a term in the likelihood for the
2039  // signal, so we'd better not zero it.
2040  sigDPLike_ = 1.0;
2041  scfDPLike_ = 1.0;
2042 
2043  const UInt_t nBkgnds = this->nBkgndClasses();
2044  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2045  if (usingBkgnd_ == kTRUE) {
2046  bkgndDPLike_[bkgndID] = 1.0;
2047  } else {
2048  bkgndDPLike_[bkgndID] = 0.0;
2049  }
2050  }
2051 
2052  return;
2053  }
2054 
2055  const UInt_t nBkgnds = this->nBkgndClasses();
2056  if ( tagged_ ) {
2057  if (curEvtCharge_==+1) {
2061 
2062  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2063  if (usingBkgnd_ == kTRUE) {
2064  bkgndDPLike_[bkgndID] = posBkgndDPModels_[bkgndID]->getLikelihood(iEvt);
2065  } else {
2066  bkgndDPLike_[bkgndID] = 0.0;
2067  }
2068  }
2069  } else {
2073 
2074  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2075  if (usingBkgnd_ == kTRUE) {
2076  bkgndDPLike_[bkgndID] = negBkgndDPModels_[bkgndID]->getLikelihood(iEvt);
2077  } else {
2078  bkgndDPLike_[bkgndID] = 0.0;
2079  }
2080  }
2081  }
2082  } else {
2085 
2088 
2089  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2090  if (usingBkgnd_ == kTRUE) {
2091  bkgndDPLike_[bkgndID] = 0.5 * ( posBkgndDPModels_[bkgndID]->getLikelihood(iEvt) +
2092  negBkgndDPModels_[bkgndID]->getLikelihood(iEvt) );
2093  } else {
2094  bkgndDPLike_[bkgndID] = 0.0;
2095  }
2096  }
2097  }
2098 
2099  if ( useSCF_ == kTRUE ) {
2100  if ( scfMap_ == 0 ) {
2101  // we're not smearing the SCF DP position
2102  // so the likelihood is the same as the TM
2104  } else {
2105  // calculate the smeared SCF DP likelihood
2106  scfDPLike_ = this->getEvtSCFDPLikelihood(iEvt);
2107  }
2108  }
2109 
2110  // Calculate the signal normalisation
2111  // NB the 2.0 is there so that the 0.5 factor is applied to
2112  // signal and background in the same place otherwise you get
2113  // normalisation problems when you switch off the DP in the fit
2114  Double_t norm = negSigModel_->getDPNorm() + posSigModel_->getDPNorm();
2115  sigDPLike_ *= 2.0/norm;
2116  scfDPLike_ *= 2.0/norm;
2117 }
2118 
2120 {
2121  Double_t scfDPLike(0.0);
2122 
2123  Double_t recoJacobian(1.0);
2124  Double_t xCoord(0.0);
2125  Double_t yCoord(0.0);
2126 
2127  Bool_t squareDP = negKinematics_->squareDP();
2128  if ( squareDP ) {
2129  xCoord = negSigModel_->getEvtmPrime();
2130  yCoord = negSigModel_->getEvtthPrime();
2131  recoJacobian = recoJacobians_[iEvt];
2132  } else {
2133  xCoord = negSigModel_->getEvtm13Sq();
2134  yCoord = negSigModel_->getEvtm23Sq();
2135  }
2136 
2137  // Find the bin that our reco event falls in
2138  Int_t recoBin = scfMap_->binNumber( xCoord, yCoord );
2139 
2140  // Find out which true Bins contribute to the given reco bin
2141  const std::vector<Int_t>* trueBins = scfMap_->trueBins(recoBin);
2142 
2143  const Int_t nDataEvents = this->eventsPerExpt();
2144 
2145  // Loop over the true bins
2146  for (std::vector<Int_t>::const_iterator iter = trueBins->begin(); iter != trueBins->end(); ++iter)
2147  {
2148  Int_t trueBin = (*iter);
2149 
2150  // prob of a true event in the given true bin migrating to the reco bin
2151  Double_t pRecoGivenTrue = scfMap_->prob( recoBin, trueBin );
2152  Double_t pTrue(0.0);
2153 
2154  // We've cached the DP amplitudes and the efficiency for the
2155  // true bin centres, just after the data points
2156  if ( tagged_ ) {
2157  LauAbsDPDynamics* sigModel(0);
2158  if (curEvtCharge_<0) {
2159  sigModel = negSigModel_;
2160  } else {
2161  sigModel = posSigModel_;
2162  }
2163 
2164  sigModel->calcLikelihoodInfo( nDataEvents + trueBin );
2165 
2166  pTrue = sigModel->getEvtDPAmp().abs2() * sigModel->getEvtEff();
2167  } else {
2168  posSigModel_->calcLikelihoodInfo( nDataEvents + trueBin );
2169  negSigModel_->calcLikelihoodInfo( nDataEvents + trueBin );
2170 
2171  pTrue = 0.5 * ( posSigModel_->getEvtDPAmp().abs2() * posSigModel_->getEvtEff() +
2173  }
2174 
2175  // Get the cached SCF fraction (and jacobian if we're using the square DP)
2176  Double_t scfFraction = fakeSCFFracs_[ trueBin ];
2177  Double_t jacobian(1.0);
2178  if ( squareDP ) {
2179  jacobian = fakeJacobians_[ trueBin ];
2180  }
2181 
2182  scfDPLike += pTrue * jacobian * scfFraction * pRecoGivenTrue;
2183  }
2184 
2185  // Divide by the reco jacobian
2186  scfDPLike /= recoJacobian;
2187 
2188  return scfDPLike;
2189 }
2190 
2192 {
2193  // Function to return the signal and background likelihoods for the
2194  // extra variables for the given event evtNo.
2195 
2196  sigExtraLike_ = 1.0;
2197 
2198  const UInt_t nBkgnds = this->nBkgndClasses();
2199 
2200  if ( ! tagged_ || curEvtCharge_ < 0 ) {
2201  sigExtraLike_ = this->prodPdfValue( negSignalPdfs_, iEvt );
2202 
2203  if (useSCF_) {
2204  scfExtraLike_ = this->prodPdfValue( negScfPdfs_, iEvt );
2205  }
2206 
2207  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2208  if (usingBkgnd_) {
2209  bkgndExtraLike_[bkgndID] = this->prodPdfValue( negBkgndPdfs_[bkgndID], iEvt );
2210  } else {
2211  bkgndExtraLike_[bkgndID] = 0.0;
2212  }
2213  }
2214  } else {
2215  sigExtraLike_ = this->prodPdfValue( posSignalPdfs_, iEvt );
2216 
2217  if (useSCF_) {
2218  scfExtraLike_ = this->prodPdfValue( posScfPdfs_, iEvt );
2219  }
2220 
2221  for ( UInt_t bkgndID(0); bkgndID < nBkgnds; ++bkgndID ) {
2222  if (usingBkgnd_) {
2223  bkgndExtraLike_[bkgndID] = this->prodPdfValue( posBkgndPdfs_[bkgndID], iEvt );
2224  } else {
2225  bkgndExtraLike_[bkgndID] = 0.0;
2226  }
2227  }
2228  }
2229 }
2230 
2232 {
2233  negCoeffs_.clear(); posCoeffs_.clear();
2234  negCoeffs_.reserve(nSigComp_); posCoeffs_.reserve(nSigComp_);
2235  for (UInt_t i = 0; i < nSigComp_; i++) {
2236  negCoeffs_.push_back(coeffPars_[i]->antiparticleCoeff());
2237  posCoeffs_.push_back(coeffPars_[i]->particleCoeff());
2238  }
2239 }
2240 
2242 {
2243  // add branches for storing the experiment number and the number of
2244  // the event within the current experiment
2245  this->addSPlotNtupleIntegerBranch("iExpt");
2246  this->addSPlotNtupleIntegerBranch("iEvtWithinExpt");
2247 
2248  // Store the efficiency of the event (for inclusive BF calculations).
2249  if (this->storeDPEff()) {
2250  this->addSPlotNtupleDoubleBranch("efficiency");
2252  this->addSPlotNtupleDoubleBranch("scffraction");
2253  }
2254  }
2255 
2256  // Store the total event likelihood for each species.
2257  if (useSCF_) {
2258  this->addSPlotNtupleDoubleBranch("sigTMTotalLike");
2259  this->addSPlotNtupleDoubleBranch("sigSCFTotalLike");
2260  this->addSPlotNtupleDoubleBranch("sigSCFFrac");
2261  } else {
2262  this->addSPlotNtupleDoubleBranch("sigTotalLike");
2263  }
2264  if (usingBkgnd_) {
2265  const UInt_t nBkgnds = this->nBkgndClasses();
2266  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2267  TString name( this->bkgndClassName(iBkgnd) );
2268  name += "TotalLike";
2269  this->addSPlotNtupleDoubleBranch(name);
2270  }
2271  }
2272 
2273  // Store the DP likelihoods
2274  if (this->useDP()) {
2275  if (useSCF_) {
2276  this->addSPlotNtupleDoubleBranch("sigTMDPLike");
2277  this->addSPlotNtupleDoubleBranch("sigSCFDPLike");
2278  } else {
2279  this->addSPlotNtupleDoubleBranch("sigDPLike");
2280  }
2281  if (usingBkgnd_) {
2282  const UInt_t nBkgnds = this->nBkgndClasses();
2283  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2284  TString name( this->bkgndClassName(iBkgnd) );
2285  name += "DPLike";
2286  this->addSPlotNtupleDoubleBranch(name);
2287  }
2288  }
2289  }
2290 
2291  // Store the likelihoods for each extra PDF
2292  if (useSCF_) {
2293  this->addSPlotNtupleBranches(&negSignalPdfs_, "sigTM");
2294  this->addSPlotNtupleBranches(&negScfPdfs_, "sigSCF");
2295  } else {
2296  this->addSPlotNtupleBranches(&negSignalPdfs_, "sig");
2297  }
2298  if (usingBkgnd_) {
2299  const UInt_t nBkgnds = this->nBkgndClasses();
2300  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2301  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2302  const LauPdfList* pdfList = &(negBkgndPdfs_[iBkgnd]);
2303  this->addSPlotNtupleBranches(pdfList, bkgndClass);
2304  }
2305  }
2306 }
2307 
2308 void LauCPFitModel::addSPlotNtupleBranches(const LauPdfList* extraPdfs, const TString& prefix)
2309 {
2310  if (extraPdfs) {
2311  // Loop through each of the PDFs
2312  for (LauPdfList::const_iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
2313 
2314  // Count the number of input variables that are not
2315  // DP variables (used in the case where there is DP
2316  // dependence for e.g. MVA)
2317  UInt_t nVars(0);
2318  std::vector<TString> varNames = (*pdf_iter)->varNames();
2319  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2320  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2321  ++nVars;
2322  }
2323  }
2324 
2325  if ( nVars == 1 ) {
2326  // If the PDF only has one variable then
2327  // simply add one branch for that variable
2328  TString varName = (*pdf_iter)->varName();
2329  TString name(prefix);
2330  name += varName;
2331  name += "Like";
2332  this->addSPlotNtupleDoubleBranch(name);
2333  } else if ( nVars == 2 ) {
2334  // If the PDF has two variables then we
2335  // need a branch for them both together and
2336  // branches for each
2337  TString allVars("");
2338  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2339  allVars += (*var_iter);
2340  TString name(prefix);
2341  name += (*var_iter);
2342  name += "Like";
2343  this->addSPlotNtupleDoubleBranch(name);
2344  }
2345  TString name(prefix);
2346  name += allVars;
2347  name += "Like";
2348  this->addSPlotNtupleDoubleBranch(name);
2349  } else {
2350  cerr << "WARNING in LauCPFitModel::addSPlotNtupleBranches : Can't yet deal with 3D PDFs." << endl;
2351  }
2352  }
2353  }
2354 }
2355 
2356 Double_t LauCPFitModel::setSPlotNtupleBranchValues(LauPdfList* extraPdfs, const TString& prefix, UInt_t iEvt)
2357 {
2358  // Store the PDF value for each variable in the list
2359  Double_t totalLike(1.0);
2360  Double_t extraLike(0.0);
2361  if (extraPdfs) {
2362  for (LauPdfList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end(); ++pdf_iter) {
2363 
2364  // calculate the likelihood for this event
2365  (*pdf_iter)->calcLikelihoodInfo(iEvt);
2366  extraLike = (*pdf_iter)->getLikelihood();
2367  totalLike *= extraLike;
2368 
2369  // Count the number of input variables that are not
2370  // DP variables (used in the case where there is DP
2371  // dependence for e.g. MVA)
2372  UInt_t nVars(0);
2373  std::vector<TString> varNames = (*pdf_iter)->varNames();
2374  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2375  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2376  ++nVars;
2377  }
2378  }
2379 
2380  if ( nVars == 1 ) {
2381  // If the PDF only has one variable then
2382  // simply store the value for that variable
2383  TString varName = (*pdf_iter)->varName();
2384  TString name(prefix);
2385  name += varName;
2386  name += "Like";
2387  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
2388  } else if ( nVars == 2 ) {
2389  // If the PDF has two variables then we
2390  // store the value for them both together
2391  // and for each on their own
2392  TString allVars("");
2393  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2394  allVars += (*var_iter);
2395  TString name(prefix);
2396  name += (*var_iter);
2397  name += "Like";
2398  Double_t indivLike = (*pdf_iter)->getLikelihood( (*var_iter) );
2399  this->setSPlotNtupleDoubleBranchValue(name, indivLike);
2400  }
2401  TString name(prefix);
2402  name += allVars;
2403  name += "Like";
2404  this->setSPlotNtupleDoubleBranchValue(name, extraLike);
2405  } else {
2406  cerr << "WARNING in LauCPFitModel::setSPlotNtupleBranchValues : Can't yet deal with 3D PDFs." << endl;
2407  }
2408  }
2409  }
2410  return totalLike;
2411 }
2412 
2414 {
2415  LauSPlot::NameSet nameSet;
2416  if (this->useDP()) {
2417  nameSet.insert("DP");
2418  }
2419  // Loop through all the signal PDFs
2420  for (LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter) {
2421  // Loop over the variables involved in each PDF
2422  std::vector<TString> varNames = (*pdf_iter)->varNames();
2423  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2424  // If they are not DP coordinates then add them
2425  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2426  nameSet.insert( (*var_iter) );
2427  }
2428  }
2429  }
2430  return nameSet;
2431 }
2432 
2434 {
2435  LauSPlot::NumbMap numbMap;
2436  if (!signalEvents_->fixed() && this->doEMLFit()) {
2437  numbMap["sig"] = signalEvents_->genValue();
2438  }
2439  if ( usingBkgnd_ ) {
2440  const UInt_t nBkgnds = this->nBkgndClasses();
2441  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2442  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2443  const LauParameter* par = bkgndEvents_[iBkgnd];
2444  if (!par->fixed()) {
2445  numbMap[bkgndClass] = par->genValue();
2446  }
2447  }
2448  }
2449  return numbMap;
2450 }
2451 
2453 {
2454  LauSPlot::NumbMap numbMap;
2455  if (signalEvents_->fixed() && this->doEMLFit()) {
2456  numbMap["sig"] = signalEvents_->genValue();
2457  }
2458  if ( usingBkgnd_ ) {
2459  const UInt_t nBkgnds = this->nBkgndClasses();
2460  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2461  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2462  const LauParameter* par = bkgndEvents_[iBkgnd];
2463  if (par->fixed()) {
2464  numbMap[bkgndClass] = par->genValue();
2465  }
2466  }
2467  }
2468  return numbMap;
2469 }
2470 
2472 {
2473  // This makes the assumption that the form of the positive and
2474  // negative PDFs are the same, which seems reasonable to me
2475 
2476  LauSPlot::TwoDMap twodimMap;
2477 
2478  for (LauPdfList::const_iterator pdf_iter = negSignalPdfs_.begin(); pdf_iter != negSignalPdfs_.end(); ++pdf_iter) {
2479  // Count the number of input variables that are not DP variables
2480  UInt_t nVars(0);
2481  std::vector<TString> varNames = (*pdf_iter)->varNames();
2482  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2483  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2484  ++nVars;
2485  }
2486  }
2487  if ( nVars == 2 ) {
2488  if (useSCF_) {
2489  twodimMap.insert( std::make_pair( "sigTM", std::make_pair( varNames[0], varNames[1] ) ) );
2490  } else {
2491  twodimMap.insert( std::make_pair( "sig", std::make_pair( varNames[0], varNames[1] ) ) );
2492  }
2493  }
2494  }
2495 
2496  if ( useSCF_ ) {
2497  for (LauPdfList::const_iterator pdf_iter = negScfPdfs_.begin(); pdf_iter != negScfPdfs_.end(); ++pdf_iter) {
2498  // Count the number of input variables that are not DP variables
2499  UInt_t nVars(0);
2500  std::vector<TString> varNames = (*pdf_iter)->varNames();
2501  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2502  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2503  ++nVars;
2504  }
2505  }
2506  if ( nVars == 2 ) {
2507  twodimMap.insert( std::make_pair( "sigSCF", std::make_pair( varNames[0], varNames[1] ) ) );
2508  }
2509  }
2510  }
2511 
2512  if (usingBkgnd_) {
2513  const UInt_t nBkgnds = this->nBkgndClasses();
2514  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2515  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2516  const LauPdfList& pdfList = negBkgndPdfs_[iBkgnd];
2517  for (LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
2518  // Count the number of input variables that are not DP variables
2519  UInt_t nVars(0);
2520  std::vector<TString> varNames = (*pdf_iter)->varNames();
2521  for ( std::vector<TString>::const_iterator var_iter = varNames.begin(); var_iter != varNames.end(); ++var_iter ) {
2522  if ( (*var_iter) != "m13Sq" && (*var_iter) != "m23Sq" ) {
2523  ++nVars;
2524  }
2525  }
2526  if ( nVars == 2 ) {
2527  twodimMap.insert( std::make_pair( bkgndClass, std::make_pair( varNames[0], varNames[1] ) ) );
2528  }
2529  }
2530  }
2531  }
2532 
2533  return twodimMap;
2534 }
2535 
2537 {
2538  cout << "Storing per-event likelihood values..." << endl;
2539 
2540  // if we've not been using the DP model then we need to cache all
2541  // the info here so that we can get the efficiency from it
2542 
2543  LauFitDataTree* inputFitData = this->fitData();
2544 
2545  if (!this->useDP() && this->storeDPEff()) {
2548  negSigModel_->fillDataTree(*inputFitData);
2549  posSigModel_->fillDataTree(*inputFitData);
2550  }
2551 
2552  UInt_t evtsPerExpt(this->eventsPerExpt());
2553 
2554  LauAbsDPDynamics* sigModel(0);
2555  LauPdfList* sigPdfs(0);
2556  LauPdfList* scfPdfs(0);
2557  LauBkgndPdfsList* bkgndPdfs(0);
2558 
2559  for (UInt_t iEvt = 0; iEvt < evtsPerExpt; ++iEvt) {
2560 
2561  this->setSPlotNtupleIntegerBranchValue("iExpt",this->iExpt());
2562  this->setSPlotNtupleIntegerBranchValue("iEvtWithinExpt",iEvt);
2563 
2564  // Find out whether we have B- or B+
2565  if ( tagged_ ) {
2566  const LauFitData& dataValues = inputFitData->getData(iEvt);
2567  LauFitData::const_iterator iter = dataValues.find("charge");
2568  curEvtCharge_ = static_cast<Int_t>(iter->second);
2569 
2570  if (curEvtCharge_==+1) {
2571  sigModel = posSigModel_;
2572  sigPdfs = &posSignalPdfs_;
2573  scfPdfs = &posScfPdfs_;
2574  bkgndPdfs = &posBkgndPdfs_;
2575  } else {
2576  sigModel = negSigModel_;
2577  sigPdfs = &negSignalPdfs_;
2578  scfPdfs = &negScfPdfs_;
2579  bkgndPdfs = &negBkgndPdfs_;
2580  }
2581  } else {
2582  sigPdfs = &negSignalPdfs_;
2583  scfPdfs = &negScfPdfs_;
2584  bkgndPdfs = &negBkgndPdfs_;
2585  }
2586 
2587  // the DP information
2588  this->getEvtDPLikelihood(iEvt);
2589  if (this->storeDPEff()) {
2590  if (!this->useDP()) {
2593  }
2594  if ( tagged_ ) {
2595  this->setSPlotNtupleDoubleBranchValue("efficiency",sigModel->getEvtEff());
2597  this->setSPlotNtupleDoubleBranchValue("scffraction",sigModel->getEvtScfFraction());
2598  }
2599  } else {
2603  }
2604  }
2605  }
2606  if (this->useDP()) {
2608  if (useSCF_) {
2610  this->setSPlotNtupleDoubleBranchValue("sigTMDPLike",sigDPLike_);
2611  this->setSPlotNtupleDoubleBranchValue("sigSCFDPLike",scfDPLike_);
2612  } else {
2613  this->setSPlotNtupleDoubleBranchValue("sigDPLike",sigDPLike_);
2614  }
2615  if (usingBkgnd_) {
2616  const UInt_t nBkgnds = this->nBkgndClasses();
2617  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2618  TString name = this->bkgndClassName(iBkgnd);
2619  name += "DPLike";
2620  this->setSPlotNtupleDoubleBranchValue(name,bkgndDPLike_[iBkgnd]);
2621  }
2622  }
2623  } else {
2624  sigTotalLike_ = 1.0;
2625  if (useSCF_) {
2626  scfTotalLike_ = 1.0;
2627  }
2628  if (usingBkgnd_) {
2629  const UInt_t nBkgnds = this->nBkgndClasses();
2630  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2631  bkgndTotalLike_[iBkgnd] = 1.0;
2632  }
2633  }
2634  }
2635 
2636  // the signal PDF values
2637  if ( useSCF_ ) {
2638  sigTotalLike_ *= this->setSPlotNtupleBranchValues(sigPdfs, "sigTM", iEvt);
2639  scfTotalLike_ *= this->setSPlotNtupleBranchValues(scfPdfs, "sigSCF", iEvt);
2640  } else {
2641  sigTotalLike_ *= this->setSPlotNtupleBranchValues(sigPdfs, "sig", iEvt);
2642  }
2643 
2644  // the background PDF values
2645  if (usingBkgnd_) {
2646  const UInt_t nBkgnds = this->nBkgndClasses();
2647  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2648  const TString& bkgndClass = this->bkgndClassName(iBkgnd);
2649  LauPdfList& pdfs = (*bkgndPdfs)[iBkgnd];
2650  bkgndTotalLike_[iBkgnd] *= this->setSPlotNtupleBranchValues(&(pdfs), bkgndClass, iEvt);
2651  }
2652  }
2653 
2654  // the total likelihoods
2655  if (useSCF_) {
2656  Double_t scfFrac(0.0);
2657  if ( useSCFHist_ ) {
2658  scfFrac = recoSCFFracs_[iEvt];
2659  } else {
2660  scfFrac = scfFrac_.value();
2661  }
2662  this->setSPlotNtupleDoubleBranchValue("sigSCFFrac",scfFrac);
2663  sigTotalLike_ *= ( 1.0 - scfFrac );
2664  if ( scfMap_ == 0 ) {
2665  scfTotalLike_ *= scfFrac;
2666  }
2667  this->setSPlotNtupleDoubleBranchValue("sigTMTotalLike",sigTotalLike_);
2668  this->setSPlotNtupleDoubleBranchValue("sigSCFTotalLike",scfTotalLike_);
2669  } else {
2670  this->setSPlotNtupleDoubleBranchValue("sigTotalLike",sigTotalLike_);
2671  }
2672  if (usingBkgnd_) {
2673  const UInt_t nBkgnds = this->nBkgndClasses();
2674  for ( UInt_t iBkgnd(0); iBkgnd < nBkgnds; ++iBkgnd ) {
2675  TString name = this->bkgndClassName(iBkgnd);
2676  name += "TotalLike";
2678  }
2679  }
2680  // fill the tree
2681  this->fillSPlotNtupleBranches();
2682  }
2683  cout << "Finished storing per-event likelihood values." << endl;
2684 }
2685 
2686 void LauCPFitModel::embedNegSignal(const TString& fileName, const TString& treeName,
2687  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
2688  Bool_t useReweighting)
2689 {
2690  if (negSignalTree_) {
2691  cerr << "ERROR in LauCPFitModel::embedNegSignal : Already embedding signal from a file." << endl;
2692  return;
2693  }
2694 
2695  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2696  cerr << "WARNING in LauCPFitModel::embedNegSignal : Conflicting options provided, will not reuse events at all." << endl;
2697  reuseEventsWithinExperiment = kFALSE;
2698  }
2699 
2700  negSignalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2701  Bool_t dataOK = negSignalTree_->findBranches();
2702  if (!dataOK) {
2703  delete negSignalTree_; negSignalTree_ = 0;
2704  cerr << "ERROR in LauCPFitModel::embedNegSignal : Problem creating data tree for embedding." << endl;
2705  return;
2706  }
2707  reuseSignal_ = reuseEventsWithinEnsemble;
2708  useNegReweighting_ = useReweighting;
2709  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2710 }
2711 
2712 void LauCPFitModel::embedNegBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
2713  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
2714 {
2715  if ( ! this->validBkgndClass( bkgndClass ) ) {
2716  cerr << "ERROR in LauCPFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
2717  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
2718  return;
2719  }
2720 
2721  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
2722 
2723  if (negBkgndTree_[bkgndID]) {
2724  cerr << "ERROR in LauCPFitModel::embedNegBkgnd : Already embedding background from a file." << endl;
2725  return;
2726  }
2727 
2728  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2729  cerr << "WARNING in LauCPFitModel::embedNegBkgnd : Conflicting options provided, will not reuse events at all." << endl;
2730  reuseEventsWithinExperiment = kFALSE;
2731  }
2732 
2733  negBkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2734  Bool_t dataOK = negBkgndTree_[bkgndID]->findBranches();
2735  if (!dataOK) {
2736  delete negBkgndTree_[bkgndID]; negBkgndTree_[bkgndID] = 0;
2737  cerr << "ERROR in LauCPFitModel::embedNegBkgnd : Problem creating data tree for embedding." << endl;
2738  return;
2739  }
2740  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
2741  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2742 }
2743 
2744 void LauCPFitModel::embedPosSignal(const TString& fileName, const TString& treeName,
2745  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment,
2746  Bool_t useReweighting)
2747 {
2748  if (posSignalTree_) {
2749  cerr << "ERROR in LauCPFitModel::embedPosSignal : Already embedding signal from a file." << endl;
2750  return;
2751  }
2752 
2753  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2754  cerr << "WARNING in LauCPFitModel::embedPosSignal : Conflicting options provided, will not reuse events at all." << endl;
2755  reuseEventsWithinExperiment = kFALSE;
2756  }
2757 
2758  posSignalTree_ = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2759  Bool_t dataOK = posSignalTree_->findBranches();
2760  if (!dataOK) {
2761  delete posSignalTree_; posSignalTree_ = 0;
2762  cerr << "ERROR in LauCPFitModel::embedPosSignal : Problem creating data tree for embedding." << endl;
2763  return;
2764  }
2765  reuseSignal_ = reuseEventsWithinEnsemble;
2766  usePosReweighting_ = useReweighting;
2767  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2768 }
2769 
2770 void LauCPFitModel::embedPosBkgnd(const TString& bkgndClass, const TString& fileName, const TString& treeName,
2771  Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment)
2772 {
2773  if ( ! this->validBkgndClass( bkgndClass ) ) {
2774  cerr << "ERROR in LauCPFitModel::embedBkgnd : Invalid background class \"" << bkgndClass << "\"." << std::endl;
2775  cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed." << endl;
2776  return;
2777  }
2778 
2779  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
2780 
2781  if (posBkgndTree_[bkgndID]) {
2782  cerr << "ERROR in LauCPFitModel::embedPosBkgnd : Already embedding background from a file." << endl;
2783  return;
2784  }
2785 
2786  if (!reuseEventsWithinEnsemble && reuseEventsWithinExperiment) {
2787  cerr << "WARNING in LauCPFitModel::embedPosBkgnd : Conflicting options provided, will not reuse events at all." << endl;
2788  reuseEventsWithinExperiment = kFALSE;
2789  }
2790 
2791  posBkgndTree_[bkgndID] = new LauEmbeddedData(fileName,treeName,reuseEventsWithinExperiment);
2792  Bool_t dataOK = posBkgndTree_[bkgndID]->findBranches();
2793  if (!dataOK) {
2794  delete posBkgndTree_[bkgndID]; posBkgndTree_[bkgndID] = 0;
2795  cerr << "ERROR in LauCPFitModel::embedPosBkgnd : Problem creating data tree for embedding." << endl;
2796  return;
2797  }
2798  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
2799  if (this->enableEmbedding() == kFALSE) {this->enableEmbedding(kTRUE);}
2800 }
2801 
2802 void LauCPFitModel::weightEvents( const TString& /*dataFileName*/, const TString& /*dataTreeName*/ )
2803 {
2804  cerr << "ERROR in LauCPFitModel::weightEvents : Method not available for this fit model." << endl;
2805  return;
2806 }
2807 
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:60
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:134
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:95
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.