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