laura is hosted by Hepforge, IPPP Durham
Laura++  v3r1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauAbsFitModel.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2014.
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 <limits>
17 #include <vector>
18 
19 #include "TMessage.h"
20 #include "TMonitor.h"
21 #include "TServerSocket.h"
22 #include "TSocket.h"
23 #include "TSystem.h"
24 #include "TVectorD.h"
25 #include "TVirtualFitter.h"
26 
27 #include "LauAbsFitModel.hh"
28 #include "LauAbsFitter.hh"
29 #include "LauAbsPdf.hh"
30 #include "LauComplex.hh"
31 #include "LauFitter.hh"
32 #include "LauFitDataTree.hh"
33 #include "LauFitNtuple.hh"
34 #include "LauGenNtuple.hh"
35 #include "LauParameter.hh"
36 #include "LauParamFixed.hh"
37 #include "LauPrint.hh"
38 #include "LauSPlot.hh"
39 
41 
42 
44  storeCon_(0),
45  twoStageFit_(kFALSE),
46  useAsymmFitErrors_(kFALSE),
47  compareFitData_(kFALSE),
48  writeLatexTable_(kFALSE),
49  writeSPlotData_(kFALSE),
50  storeDPEff_(kFALSE),
51  randomFit_(kFALSE),
52  emlFit_(kFALSE),
53  poissonSmear_(kFALSE),
54  enableEmbedding_(kFALSE),
55  usingDP_(kTRUE),
56  pdfsDependOnDP_(kFALSE),
57  firstExpt_(0),
58  nExpt_(0),
59  evtsPerExpt_(0),
60  iExpt_(0),
61  inputFitData_(0),
62  fitNtuple_(0),
63  genNtuple_(0),
64  sPlotNtuple_(0),
65  fitStatus_(0),
66  NLL_(0),
67  numberOKFits_(0),
68  numberBadFits_(0),
69  nParams_(0),
70  nFreeParams_(0),
71  worstLogLike_(std::numeric_limits<Double_t>::max()),
72  withinAsymErrorCalc_(kFALSE),
73  nullString_(""),
74  doSFit_(kFALSE),
75  sWeightBranchName_(""),
76  sWeightScaleFactor_(1.0),
77  outputTableName_(""),
78  fitToyMCFileName_("fitToyMC.root"),
79  fitToyMCTableName_("fitToyMCTable"),
80  fitToyMCScale_(10),
81  fitToyMCPoissonSmear_(kFALSE),
82  sPlotFileName_(""),
83  sPlotTreeName_(""),
84  sPlotVerbosity_("")
85 {
86 }
87 
88 
90 {
91  delete inputFitData_; inputFitData_ = 0;
92  delete fitNtuple_; fitNtuple_ = 0;
93  delete genNtuple_; genNtuple_ = 0;
94  delete sPlotNtuple_; sPlotNtuple_ = 0;
95 
96  // Remove the components created to apply constraints to fit parameters
97  for (std::vector<LauAbsRValue*>::iterator iter = conVars_.begin(); iter != conVars_.end(); ++iter){
98  if ( !(*iter)->isLValue() ){
99  delete (*iter);
100  (*iter) = 0;
101  }
102  }
103 }
104 
105 void LauAbsFitModel::run(const TString& applicationCode, const TString& dataFileName, const TString& dataTreeName,
106  const TString& histFileName, const TString& tableFileName)
107 {
108  // Chose whether you want to generate or fit events in the Dalitz plot.
109  // To generate events choose applicationCode = "gen", to fit events choose
110  // applicationCode = "fit".
111 
112  TString runCode(applicationCode);
113  runCode.ToLower();
114 
115  TString histFileNameCopy(histFileName);
116  TString tableFileNameCopy(tableFileName);
117  TString dataFileNameCopy(dataFileName);
118  TString dataTreeNameCopy(dataTreeName);
119 
120  // Initialise the fit par vectors. Each class that inherits from this one
121  // must implement this sensibly for all vectors specified in clearFitParVectors,
122  // i.e. specify parameter names, initial, min, max and fixed values
123  this->initialise();
124 
125  // Add variables to Gaussian constrain to a list
126  this->addConParameters();
127 
128  if (dataFileNameCopy == "") {dataFileNameCopy = "data.root";}
129  if (dataTreeNameCopy == "") {dataTreeNameCopy = "genResults";}
130 
131  if (runCode.Contains("gen")) {
132 
133  if (histFileNameCopy == "") {histFileNameCopy = "parInfo.root";}
134  if (tableFileNameCopy == "") {tableFileNameCopy = "genResults";}
135 
136  this->setGenValues();
137  this->generate(dataFileNameCopy, dataTreeNameCopy, histFileNameCopy, tableFileNameCopy);
138 
139  } else if (runCode.Contains("fit")) {
140 
141  if (histFileNameCopy == "") {histFileNameCopy = "parInfo.root";}
142  if (tableFileNameCopy == "") {tableFileNameCopy = "fitResults";}
143 
144  this->fit(dataFileNameCopy, dataTreeNameCopy, histFileNameCopy, tableFileNameCopy);
145 
146  } else if (runCode.Contains("plot")) {
147 
148  this->savePDFPlots("plot");
149 
150  } else if (runCode.Contains("reweight")) {
151 
152  this->weightEvents(dataFileNameCopy, dataTreeNameCopy);
153  }
154 }
155 
156 void LauAbsFitModel::runSlave(const TString& dataFileName, const TString& dataTreeName,
157  const TString& histFileName, const TString& tableFileName,
158  const TString& addressMaster, const UInt_t portMaster)
159 {
160  // Establish the connection to the master process
161  this->connectToMaster( addressMaster, portMaster );
162 
163  // Initialise the fit par vectors. Each class that inherits from this one
164  // must implement this sensibly for all vectors specified in clearFitParVectors,
165  // i.e. specify parameter names, initial, min, max and fixed values
166  this->initialise();
167 
168  // NB call to addConParameters() is intentionally not included here cf.
169  // run() since this has to be dealt with by the master to avoid
170  // multiple inclusions of each penalty term
171 
172  TString dataFileNameCopy(dataFileName);
173  TString dataTreeNameCopy(dataTreeName);
174  TString histFileNameCopy(histFileName);
175  TString tableFileNameCopy(tableFileName);
176 
177  if (dataFileNameCopy == "") {dataFileNameCopy = "data.root";}
178  if (dataTreeNameCopy == "") {dataTreeNameCopy = "genResults";}
179  if (histFileNameCopy == "") {histFileNameCopy = "parInfo.root";}
180  if (tableFileNameCopy == "") {tableFileNameCopy = "fitResults";}
181 
182  this->fitSlave(dataFileNameCopy, dataTreeNameCopy, histFileNameCopy, tableFileNameCopy);
183 
184  std::cout << "INFO in LauAbsFitModel::runSlave : Fit slave " << this->slaveId() << " has finished successfully" << std::endl;
185 }
186 
187 void LauAbsFitModel::doSFit( const TString& sWeightBranchName, Double_t scaleFactor )
188 {
189  if ( sWeightBranchName == "" ) {
190  std::cerr << "WARNING in LauAbsFitModel::doSFit : sWeight branch name is empty string, not setting-up sFit." << std::endl;
191  return;
192  }
193 
194  doSFit_ = kTRUE;
195  sWeightBranchName_ = sWeightBranchName;
196  sWeightScaleFactor_ = scaleFactor;
197 }
198 
199 void LauAbsFitModel::setBkgndClassNames( const std::vector<TString>& names )
200 {
201  if ( !bkgndClassNames_.empty() ) {
202  std::cerr << "WARNING in LauAbsFitModel::setBkgndClassNames : Names already stored, not changing them." << std::endl;
203  return;
204  }
205 
206  UInt_t nBkgnds = names.size();
207  for ( UInt_t i(0); i < nBkgnds; ++i ) {
208  bkgndClassNames_.insert( std::make_pair( i, names[i] ) );
209  }
210 
211  this->setupBkgndVectors();
212 }
213 
214 Bool_t LauAbsFitModel::validBkgndClass( const TString& className ) const
215 {
216  if ( bkgndClassNames_.empty() ) {
217  return kFALSE;
218  }
219 
220  Bool_t found(kFALSE);
221  for ( LauBkgndClassMap::const_iterator iter = bkgndClassNames_.begin(); iter != bkgndClassNames_.end(); ++iter ) {
222  if ( iter->second == className ) {
223  found = kTRUE;
224  break;
225  }
226  }
227 
228  return found;
229 }
230 
231 UInt_t LauAbsFitModel::bkgndClassID( const TString& className ) const
232 {
233  if ( ! this->validBkgndClass( className ) ) {
234  std::cerr << "ERROR in LauAbsFitModel::bkgndClassID : Request for ID for invalid background class \"" << className << "\"." << std::endl;
235  return (bkgndClassNames_.size() + 1);
236  }
237 
238  UInt_t bgID(0);
239  for ( LauBkgndClassMap::const_iterator iter = bkgndClassNames_.begin(); iter != bkgndClassNames_.end(); ++iter ) {
240  if ( iter->second == className ) {
241  bgID = iter->first;
242  break;
243  }
244  }
245 
246  return bgID;
247 }
248 
249 const TString& LauAbsFitModel::bkgndClassName( UInt_t classID ) const
250 {
251  LauBkgndClassMap::const_iterator iter = bkgndClassNames_.find( classID );
252 
253  if ( iter == bkgndClassNames_.end() ) {
254  std::cerr << "ERROR in LauAbsFitModel::bkgndClassName : Request for name of invalid background class ID " << classID << "." << std::endl;
255  return nullString_;
256  }
257 
258  return iter->second;
259 }
260 
262 {
263  std::cout << "INFO in LauAbsFitModel::clearFitParVectors : Clearing fit variable vectors" << std::endl;
264 
265  // Remove the components created to apply constraints to fit parameters
266  for (std::vector<LauAbsRValue*>::iterator iter = conVars_.begin(); iter != conVars_.end(); ++iter){
267  if ( !(*iter)->isLValue() ){
268  delete (*iter);
269  (*iter) = 0;
270  }
271  }
272  conVars_.clear();
273  fitVars_.clear();
274 }
275 
277 {
278  std::cout << "INFO in LauAbsFitModel::clearExtraVarVectors : Clearing extra ntuple variable vectors" << std::endl;
279  extraVars_.clear();
280 }
281 
283 {
284  // makes sure each parameter holds its genValue as its current value
285  for (LauParameterPList::iterator iter = fitVars_.begin(); iter != fitVars_.end(); ++iter) {
286  (*iter)->value((*iter)->genValue());
287  }
288  this->propagateParUpdates();
289 }
290 
291 void LauAbsFitModel::writeSPlotData(const TString& fileName, const TString& treeName, Bool_t storeDPEfficiency, const TString& verbosity)
292 {
293  if (this->writeSPlotData()) {
294  std::cerr << "ERROR in LauAbsFitModel::writeSPlotData : Already have an sPlot ntuple setup, not doing it again." << std::endl;
295  return;
296  }
297  writeSPlotData_ = kTRUE;
298  sPlotFileName_ = fileName;
299  sPlotTreeName_ = treeName;
300  sPlotVerbosity_ = verbosity;
301  storeDPEff_ = storeDPEfficiency;
302 }
303 
304 // TODO : histFileName isn't used here at the moment but could be used for
305 // storing the values of the parameters used in the generation.
306 // These could then be read and used for setting the "true" values
307 // in a subsequent fit.
308 void LauAbsFitModel::generate(const TString& dataFileName, const TString& dataTreeName, const TString& /*histFileName*/, const TString& tableFileNameBase)
309 {
310  // Create the ntuple for storing the results
311  std::cout << "INFO in LauAbsFitModel::generate : Creating generation ntuple." << std::endl;
312  if (genNtuple_ != 0) {delete genNtuple_; genNtuple_ = 0;}
313  genNtuple_ = new LauGenNtuple(dataFileName,dataTreeName);
314 
315  // add branches for storing the experiment number and the number of
316  // the event within the current experiment
317  this->addGenNtupleIntegerBranch("iExpt");
318  this->addGenNtupleIntegerBranch("iEvtWithinExpt");
319  this->setupGenNtupleBranches();
320 
321  // Start the cumulative timer
322  cumulTimer_.Start();
323 
324  Bool_t genOK(kTRUE);
325  do {
326  // Loop over the number of experiments
327  for (iExpt_ = firstExpt_; iExpt_ < (firstExpt_+nExpt_); ++iExpt_) {
328 
329  // Start the timer to see how long each experiment takes to generate
330  timer_.Start();
331 
332  // Store the experiment number in the ntuple
333  this->setGenNtupleIntegerBranchValue("iExpt",iExpt_);
334 
335  // Do the generation for this experiment
336  std::cout << "INFO in LauAbsFitModel::generate : Generating experiment number " << iExpt_ << std::endl;
337  genOK = this->genExpt();
338 
339  // Stop the timer and see how long the program took so far
340  timer_.Stop();
341  timer_.Print();
342 
343  if (!genOK) {
344  // delete and recreate an empty tree
346 
347  // then break out of the experiment loop
348  std::cerr << "WARNING in LauAbsFitModel::generate : Problem in toy MC generation. Starting again with updated parameters..." << std::endl;
349  break;
350  }
351 
352  if (this->writeLatexTable()) {
353  TString tableFileName(tableFileNameBase);
354  tableFileName += "_";
355  tableFileName += iExpt_;
356  tableFileName += ".tex";
357  this->writeOutTable(tableFileName);
358  }
359 
360  } // Loop over number of experiments
361  } while (!genOK);
362 
363  // Print out total timing info.
364  cumulTimer_.Stop();
365  std::cout << "INFO in LauAbsFitModel::generate : Finished generating all experiments." << std::endl;
366  std::cout << "INFO in LauAbsFitModel::generate : Cumulative timing:" << std::endl;
367  cumulTimer_.Print();
368 
369  // Build the event index
370  std::cout << "INFO in LauAbsFitModel::generate : Building experiment:event index." << std::endl;
371  // TODO - can test this return value?
372  //Int_t nIndexEntries =
373  genNtuple_->buildIndex("iExpt","iEvtWithinExpt");
374 
375  // Write out toy MC ntuple
376  std::cout << "INFO in LauAbsFitModel::generate : Writing data to file " << dataFileName << "." << std::endl;
378 }
379 
381 {
383 }
384 
386 {
388 }
389 
391 {
392  genNtuple_->setIntegerBranchValue(name,value);
393 }
394 
396 {
397  genNtuple_->setDoubleBranchValue(name,value);
398 }
399 
401 {
402  return genNtuple_->getIntegerBranchValue(name);
403 }
404 
406 {
407  return genNtuple_->getDoubleBranchValue(name);
408 }
409 
411 {
413 }
414 
416 {
418 }
419 
421 {
423 }
424 
426 {
427  sPlotNtuple_->setIntegerBranchValue(name,value);
428 }
429 
431 {
432  sPlotNtuple_->setDoubleBranchValue(name,value);
433 }
434 
436 {
438 }
439 
440 void LauAbsFitModel::fit(const TString& dataFileName, const TString& dataTreeName, const TString& histFileName, const TString& tableFileNameBase)
441 {
442  // Routine to perform the total fit.
443 
444  std::cout << "INFO in LauAbsFitModel::fit : First experiment = " << firstExpt_ << std::endl;
445  std::cout << "INFO in LauAbsFitModel::fit : Number of experiments = " << nExpt_ << std::endl;
446 
447  // Start the cumulative timer
448  cumulTimer_.Start();
449 
450  numberOKFits_ = 0, numberBadFits_ = 0;
451  fitStatus_ = -1;
452 
453  // Create and setup the fit results ntuple
454  std::cout << "INFO in LauAbsFitModel::fit : Creating fit ntuple." << std::endl;
455  if (fitNtuple_ != 0) {delete fitNtuple_; fitNtuple_ = 0;}
456  fitNtuple_ = new LauFitNtuple(histFileName, this->useAsymmFitErrors());
457 
458  // Create and setup the sPlot ntuple
459  if (this->writeSPlotData()) {
460  std::cout << "INFO in LauAbsFitModel::fit : Creating sPlot ntuple." << std::endl;
461  if (sPlotNtuple_ != 0) {delete sPlotNtuple_; sPlotNtuple_ = 0;}
463  this->setupSPlotNtupleBranches();
464  }
465 
466  // This reads in the given dataFile and creates an input
467  // fit data tree that stores them for all events and experiments.
468  Bool_t dataOK = this->cacheFitData(dataFileName,dataTreeName);
469  if (!dataOK) {
470  std::cerr << "ERROR in LauAbsFitModel::fit : Problem caching the fit data." << std::endl;
471  gSystem->Exit(EXIT_FAILURE);
472  }
473 
474  // Loop over the number of experiments
475  for (iExpt_ = firstExpt_; iExpt_ < (firstExpt_+nExpt_); ++iExpt_) {
476 
477  // Start the timer to see how long each fit takes
478  timer_.Start();
479 
480  UInt_t nEvents = this->readExperimentData( iExpt_ );
481  if (nEvents < 1) {
482  std::cerr << "WARNING in LauAbsFitModel::fit : Zero events in experiment " << iExpt_ << ", skipping..." << std::endl;
483  timer_.Stop();
484  continue;
485  }
486 
487  // Now the sub-classes must implement whatever they need to do
488  // to cache any more input fit data they need in order to
489  // calculate the likelihoods during the fit.
490  // They need to use the inputFitData_ tree as input. For example,
491  // inputFitData_ contains m13Sq and m23Sq. The appropriate fit model
492  // then caches the resonance dynamics for the signal model, as
493  // well as the background likelihood values in the Dalitz plot
494  this->cacheInputFitVars();
495  if ( this->doSFit() ) {
496  this->cacheInputSWeights();
497  }
498 
499  // Do the fit for this experiment
500  this->fitExpt();
501 
502  // Write the results into the ntuple
503  this->finaliseFitResults(tableFileNameBase);
504 
505  // Stop the timer and see how long the program took so far
506  timer_.Stop();
507  timer_.Print();
508 
509  // Store the per-event likelihood values
510  if ( this->writeSPlotData() ) {
511  this->storePerEvtLlhds();
512  }
513 
514  // Create a toy MC sample using the fitted parameters so that
515  // the user can compare the fit to the data.
516  if (compareFitData_ == kTRUE && fitStatus_ == 3) {
518  }
519 
520  // Keep track of how many fits worked or failed
521  // NB values of fitStatus_ now indicate the status of the error matrix:
522  // 0= not calculated at all
523  // 1= approximation only, not accurate
524  // 2= full matrix, but forced positive-definite
525  // 3= full accurate covariance matrix
526  if (fitStatus_ == 3) {
527  numberOKFits_++;
528  } else {
529  numberBadFits_++;
530  }
531 
532  } // Loop over number of experiments
533 
534  // Print out total timing info.
535  cumulTimer_.Stop();
536  std::cout << "INFO in LauAbsFitModel::fit : Cumulative timing:" << std::endl;
537  cumulTimer_.Print();
538 
539  // Print out stats on OK fits.
540  std::cout << "INFO in LauAbsFitModel::fit : Number of OK Fits = " << numberOKFits_ << std::endl;
541  std::cout << "INFO in LauAbsFitModel::fit : Number of Failed Fits = " << numberBadFits_ << std::endl;
542  Double_t fitEff(0.0);
543  if (nExpt_ != 0) {fitEff = numberOKFits_/(1.0*nExpt_);}
544  std::cout << "INFO in LauAbsFitModel::fit : Fit efficiency = " << fitEff*100.0 << "%." << std::endl;
545 
546  // Write out any fit results (ntuples etc...).
547  this->writeOutAllFitResults();
548  if ( this->writeSPlotData() ) {
549  this->calculateSPlotData();
550  }
551 }
552 
553 void LauAbsFitModel::fitSlave(const TString& dataFileName, const TString& dataTreeName, const TString& histFileName, const TString& tableFileNameBase)
554 {
555  outputTableName_ = tableFileNameBase;
556 
557  // Create and setup the fit results ntuple
558  std::cout << "INFO in LauAbsFitModel::fitSlave : Creating fit ntuple." << std::endl;
559  if (fitNtuple_ != 0) {delete fitNtuple_; fitNtuple_ = 0;}
560  fitNtuple_ = new LauFitNtuple(histFileName, this->useAsymmFitErrors());
561 
562  // This reads in the given dataFile and creates an input
563  // fit data tree that stores them for all events and experiments.
564  Bool_t dataOK = this->cacheFitData(dataFileName,dataTreeName);
565  if (!dataOK) {
566  std::cerr << "ERROR in LauAbsFitModel::fitSlave : Problem caching the fit data." << std::endl;
567  gSystem->Exit(EXIT_FAILURE);
568  }
569 
570  // Now process the various requests from the master
571  this->processMasterRequests();
572 }
573 
574 Bool_t LauAbsFitModel::cacheFitData(const TString& dataFileName, const TString& dataTreeName)
575 {
576  // From the input data stream, store the variables into the
577  // internal tree inputFitData_ that can be used by the sub-classes
578  // in calculating their likelihood functions for the fit
579  delete inputFitData_;
580  inputFitData_ = new LauFitDataTree(dataFileName,dataTreeName);
581  Bool_t dataOK = inputFitData_->findBranches();
582 
583  if (!dataOK) {
584  delete inputFitData_; inputFitData_ = 0;
585  }
586 
587  return dataOK;
588 }
589 
591 {
592  Bool_t hasBranch = inputFitData_->haveBranch( sWeightBranchName_ );
593  if ( ! hasBranch ) {
594  std::cerr << "ERROR in LauAbsFitModel::cacheInputSWeights : Input data does not contain variable \"" << sWeightBranchName_ << "\".\n";
595  std::cerr << " : Turning off sFit!" << std::endl;
596  doSFit_ = kFALSE;
597  return;
598  }
599 
600  UInt_t nEvents = this->eventsPerExpt();
601  sWeights_.clear();
602  sWeights_.reserve( nEvents );
603 
604  for (UInt_t iEvt = 0; iEvt < nEvents; ++iEvt) {
605 
606  const LauFitData& dataValues = inputFitData_->getData(iEvt);
607 
608  LauFitData::const_iterator iter = dataValues.find( sWeightBranchName_ );
609 
610  sWeights_.push_back( iter->second * sWeightScaleFactor_ );
611  }
612 }
613 
615 {
616  // Routine to perform the actual fit for the given experiment
617 
618  // Reset the worst likelihood found to its catch-all value
619  worstLogLike_ = std::numeric_limits<Double_t>::max();
620 
621  // Update initial fit parameters if required (e.g. if using random numbers).
622  this->checkInitFitParams();
623 
624  // Initialise the fitter
628 
631 
632  // Now ready for minimisation step
633  std::cout << "\nINFO in LauAbsFitModel::fitExpt : Start minimisation...\n";
634  std::pair<Int_t,Double_t> fitResult = LauFitter::fitter()->minimise();
635 
636  fitStatus_ = fitResult.first;
637  NLL_ = fitResult.second;
638 
639  // If we're doing a two stage fit we can now release (i.e. float)
640  // the 2nd stage parameters and re-fit
641  if (this->twoStageFit()) {
642 
643  if ( fitStatus_ != 3 ) {
644  std::cerr << "WARNING in LauAbsFitModel:fitExpt : Not running second stage fit since first stage failed." << std::endl;
646  } else {
650  fitResult = LauFitter::fitter()->minimise();
651  }
652  }
653 
654  fitStatus_ = fitResult.first;
655  NLL_ = fitResult.second;
656  const TMatrixD& covMat = LauFitter::fitter()->covarianceMatrix();
657  covMatrix_.Clear();
658  covMatrix_.ResizeTo( covMat.GetNrows(), covMat.GetNcols() );
659  covMatrix_.SetMatrixArray( covMat.GetMatrixArray() );
660 
661  // Store the final fit results and errors into protected internal vectors that
662  // all sub-classes can use within their own finalFitResults implementation
663  // used below (e.g. putting them into an ntuple in a root file)
665 }
666 
668 {
669  // Write out histograms at end
670  if (fitNtuple_ != 0) {
672  }
673 }
674 
676 {
677  if (sPlotNtuple_ != 0) {
680  LauSPlot splot(sPlotNtuple_->fileName(), sPlotNtuple_->treeName(), this->firstExpt(), this->nExpt(),
681  this->variableNames(), this->freeSpeciesNames(), this->fixdSpeciesNames(), this->twodimPDFs(),
682  this->splitSignal(), this->scfDPSmear());
684  splot.writeOutResults();
685  }
686 }
687 
688 void LauAbsFitModel::compareFitData(UInt_t toyMCScale, const TString& mcFileName, const TString& tableFileName, Bool_t poissonSmearing)
689 {
690  compareFitData_ = kTRUE;
691  fitToyMCScale_ = toyMCScale;
692  fitToyMCFileName_ = mcFileName;
693  fitToyMCTableName_ = tableFileName;
694  fitToyMCPoissonSmear_ = poissonSmearing;
695 }
696 
697 void LauAbsFitModel::createFitToyMC(const TString& mcFileName, const TString& tableFileName)
698 {
699  // Create a toy MC sample so that the user can compare the fitted
700  // result with the data.
701  // Generate more toy MC to reduce statistical fluctuations:
702  // - use the rescaling value fitToyMCScale_
703 
704  // Store the info on the number of experiments, first expt and current expt
705  const UInt_t oldNExpt(this->nExpt());
706  const UInt_t oldFirstExpt(this->firstExpt());
707  const UInt_t oldIExpt(iExpt_);
708 
709  // Turn off Poisson smearing if required
710  const Bool_t poissonSmearing(this->doPoissonSmearing());
712 
713  // Turn off embedding, since we need toy MC, not reco'ed events
714  const Bool_t enableEmbeddingOrig(this->enableEmbedding());
715  this->enableEmbedding(kFALSE);
716 
717  // Need to make sure that the generation of the DP co-ordinates is
718  // switched on if any of our PDFs depend on it
719  const Bool_t origUseDP = this->useDP();
720  if ( !origUseDP && this->pdfsDependOnDP() ) {
721  this->useDP( kTRUE );
722  this->initialiseDPModels();
723  }
724 
725  // Construct a unique filename for this experiment
726  TString exptString("_expt");
727  exptString += oldIExpt;
728  TString fileName( mcFileName );
729  fileName.Insert( fileName.Last('.'), exptString );
730 
731  // Generate the toy MC
732  std::cout << "INFO in LauAbsFitModel::createFitToyMC : Generating toy MC in " << fileName << " to compare fit with data..." << std::endl;
733  std::cout << " : Number of experiments to generate = " << fitToyMCScale_ << "." << std::endl;
734  std::cout << " : This is to allow the toy MC to be made with reduced statistical fluctuations." << std::endl;
735 
736  // Set the genValue of each parameter to its current (fitted) value
737  // but first store the original genValues for restoring later
738  std::vector<Double_t> origGenValues; origGenValues.reserve(nParams_);
739  Bool_t blind(kFALSE);
740  for (LauParameterPList::iterator iter = fitVars_.begin(); iter != fitVars_.end(); ++iter) {
741  origGenValues.push_back((*iter)->genValue());
742  (*iter)->genValue((*iter)->unblindValue());
743  if ( (*iter)->blind() ) {
744  blind = kTRUE;
745  }
746  }
747  if ( blind ) {
748  std::cerr << "WARNING in LauAbsFitModel::createFitToyMC : One or more parameters are blind but the toy will be created using the unblind values - use with caution!!" << std::endl;
749  }
750 
751  // If we're asked to generate more than 100 experiments then split it
752  // up into multiple files since otherwise can run into memory issues
753  // when building the index
754 
755  UInt_t totalExpts = fitToyMCScale_;
756  if ( totalExpts > 100 ) {
757  UInt_t nFiles = totalExpts/100;
758  if ( totalExpts%100 ) {
759  nFiles += 1;
760  }
761  for ( UInt_t iFile(0); iFile < nFiles; ++iFile ) {
762 
763  UInt_t firstExp( iFile*100 );
764 
765  // Set number of experiments and first experiment to generate
766  UInt_t nExp = ((firstExp + 100)>totalExpts) ? totalExpts-firstExp : 100;
767  this->setNExpts(nExp, firstExp);
768 
769  // Create a unique filename and generate the events
770  TString extraname = "_file";
771  extraname += iFile;
772  fileName.Insert( fileName.Last('.'), extraname );
773  this->generate(fileName, "genResults", "dummy.root", tableFileName);
774  }
775  } else {
776  // Set number of experiments to new value
777  this->setNExpts(fitToyMCScale_, 0);
778  // Generate the toy
779  this->generate(fileName, "genResults", "dummy.root", tableFileName);
780  }
781 
782  // Reset number of experiments to original value
783  iExpt_ = oldIExpt;
784  this->setNExpts(oldNExpt, oldFirstExpt);
785 
786  // Restore the Poisson smearing to its former value
787  this->doPoissonSmearing(poissonSmearing);
788 
789  // Restore the embedding status
790  this->enableEmbedding(enableEmbeddingOrig);
791 
792  // Restore "useDP" to its former status
793  this->useDP( origUseDP );
794 
795  // Restore the original genValue to each parameter
796  for (UInt_t i(0); i<nParams_; ++i) {
797  fitVars_[i]->genValue(origGenValues[i]);
798  }
799 
800  std::cout << "INFO in LauAbsFitModel::createFitToyMC : Finished in createFitToyMC." << std::endl;
801 }
802 
804 {
805  // Calculate the total negative log-likelihood over all events.
806  // This function assumes that the fit parameters and data tree have
807  // already been set-up correctly.
808 
809  // Loop over the data points to calculate the log likelihood
810  Double_t logLike = this->getLogLikelihood( 0, this->eventsPerExpt() );
811 
812  // Include the Poisson term in the extended likelihood if required
813  if (this->doEMLFit()) {
814  logLike -= this->getEventSum();
815  }
816 
817  // Calculate any penalty terms from Gaussian constrained variables
818  if ( ! conVars_.empty() ){
819  logLike -= this->getLogLikelihoodPenalty();
820  }
821 
822  Double_t totNegLogLike = -logLike;
823  return totNegLogLike;
824 }
825 
827 {
828  Double_t penalty(0.0);
829 
830  for ( LauAbsRValuePList::const_iterator iter = conVars_.begin(); iter != conVars_.end(); ++iter ) {
831  Double_t val = (*iter)->unblindValue();
832  Double_t mean = (*iter)->constraintMean();
833  Double_t width = (*iter)->constraintWidth();
834 
835  Double_t term = ( val - mean )*( val - mean );
836  penalty += term/( 2*width*width );
837  }
838 
839  return penalty;
840 }
841 
842 Double_t LauAbsFitModel::getLogLikelihood( UInt_t iStart, UInt_t iEnd )
843 {
844  // Calculate the total negative log-likelihood over all events.
845  // This function assumes that the fit parameters and data tree have
846  // already been set-up correctly.
847 
848  // Loop over the data points to calculate the log likelihood
849  Double_t logLike(0.0);
850 
851  // Loop over the number of events in this experiment
852  Bool_t ok(kTRUE);
853  for (UInt_t iEvt = iStart; iEvt < iEnd; ++iEvt) {
854 
855  Double_t likelihood = this->getTotEvtLikelihood(iEvt);
856 
857  if (likelihood > DBL_MIN) { // Is the likelihood zero?
858  Double_t evtLogLike = TMath::Log(likelihood);
859  if ( doSFit_ ) {
860  evtLogLike *= sWeights_[iEvt];
861  }
862  logLike += evtLogLike;
863  } else {
864  ok = kFALSE;
865  std::cerr << "WARNING in LauAbsFitModel::getLogLikelihood : Strange likelihood value for event " << iEvt << ": " << likelihood << "\n";
866  this->printEventInfo(iEvt);
867  this->printVarsInfo(); //Write the values of the floated variables for which the likelihood is zero
868  break;
869  }
870  }
871 
872  if (!ok) {
873  std::cerr << " : Returning worst NLL found so far to force MINUIT out of this region." << std::endl;
874  logLike = worstLogLike_;
875  } else if (logLike < worstLogLike_) {
876  worstLogLike_ = logLike;
877  }
878 
879  return logLike;
880 }
881 
882 void LauAbsFitModel::setParsFromMinuit(Double_t* par, Int_t npar)
883 {
884  // This function sets the internal parameters based on the values
885  // that Minuit is using when trying to minimise the total likelihood function.
886 
887  // MINOS reports different numbers of free parameters depending on the
888  // situation, so disable this check
889  if ( ! withinAsymErrorCalc_ ) {
890  if (static_cast<UInt_t>(npar) != nFreeParams_) {
891  std::cerr << "ERROR in LauAbsFitModel::setParsFromMinuit : Unexpected number of free parameters: " << npar << ".\n";
892  std::cerr << " Expected: " << nFreeParams_ << ".\n" << std::endl;
893  gSystem->Exit(EXIT_FAILURE);
894  }
895  }
896 
897  // Despite npar being the number of free parameters
898  // the par array actually contains all the parameters,
899  // free and floating...
900 
901  // Update all the floating ones with their new values
902  // Also check if we have any parameters on which the DP integrals depend
903  // and whether they have changed since the last iteration
904  Bool_t recalcNorm(kFALSE);
905  const LauParameterPSet::const_iterator resVarsEnd = resVars_.end();
906  for (UInt_t i(0); i<nParams_; ++i) {
907  if (!fitVars_[i]->fixed()) {
908  if ( resVars_.find( fitVars_[i] ) != resVarsEnd ) {
909  if ( fitVars_[i]->value() != par[i] ) {
910  recalcNorm = kTRUE;
911  }
912  }
913  fitVars_[i]->value(par[i]);
914  }
915  }
916 
917  // If so, then recalculate the normalisation
918  if (recalcNorm) {
919  this->recalculateNormalisation();
920  }
921 
922  this->propagateParUpdates();
923 }
924 
926 {
927  UInt_t nParsAdded(0);
928  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
929  LauAbsPdf* pdf = (*pdf_iter);
930  if ( pdf->isDPDependent() ) {
931  this->pdfsDependOnDP( kTRUE );
932  }
933  LauAbsRValuePList& pars = pdf->getParameters();
934  for (LauAbsRValuePList::iterator pars_iter = pars.begin(); pars_iter != pars.end(); ++pars_iter) {
935  LauParameterPList params = (*pars_iter)->getPars();
936  for (LauParameterPList::iterator params_iter = params.begin(); params_iter != params.end(); ++params_iter) {
937  if ( !(*params_iter)->clone() && ( !(*params_iter)->fixed() || ( this->twoStageFit() && (*params_iter)->secondStage() ) ) ) {
938  fitVars_.push_back(*params_iter);
939  ++nParsAdded;
940  }
941  }
942  }
943  }
944  return nParsAdded;
945 }
946 
947 void LauAbsFitModel::addConstraint(const TString& formula, const std::vector<TString>& pars, const Double_t mean, const Double_t width)
948 {
949  StoreConstraints newCon;
950  newCon.formula_ = formula;
951  newCon.conPars_ = pars;
952  newCon.mean_ = mean;
953  newCon.width_ = width;
954  storeCon_.push_back(newCon);
955 }
956 
958 {
959  for ( LauParameterPList::const_iterator iter = fitVars_.begin(); iter != fitVars_.end(); ++iter ) {
960  if ( (*iter)->gaussConstraint() ) {
961  conVars_.push_back( *iter );
962  std::cout << "INFO in LauAbsFitModel::addConParameters : Added Gaussian constraint to parameter "<< (*iter)->name() << std::endl;
963  }
964  }
965 
966  // Add penalties from the constraints to fit parameters
967  for ( std::vector<StoreConstraints>::iterator iter = storeCon_.begin(); iter != storeCon_.end(); ++iter ) {
968  std::vector<TString> names = (*iter).conPars_;
969  std::vector<LauParameter*> params;
970  for ( std::vector<TString>::iterator iternames = names.begin(); iternames != names.end(); ++iternames ) {
971  for ( LauParameterPList::const_iterator iterfit = fitVars_.begin(); iterfit != fitVars_.end(); ++iterfit ) {
972  if ( (*iternames) == (*iterfit)->name() ){
973  params.push_back(*iterfit);
974  }
975  }
976  }
977 
978  // If the parameters are not found, skip it
979  if ( params.size() != (*iter).conPars_.size() ) {
980  std::cerr << "WARNING in LauAbsFitModel::addConParameters: Could not find parameters to constrain in the formula... skipping" << std::endl;
981  continue;
982  }
983 
984  LauFormulaPar* formPar = new LauFormulaPar( (*iter).formula_, (*iter).formula_, params );
985  formPar->addGaussianConstraint( (*iter).mean_, (*iter).width_ );
986  conVars_.push_back(formPar);
987 
988  std::cout << "INFO in LauAbsFitModel::addConParameters : Added Gaussian constraint to formula\n";
989  std::cout << " : Formula: " << (*iter).formula_ << std::endl;
990  for ( std::vector<LauParameter*>::iterator iterparam = params.begin(); iterparam != params.end(); ++iterparam ) {
991  std::cout << " : Parameter: " << (*iterparam)->name() << std::endl;
992  }
993  }
994 
995 }
996 
998 {
999  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1000  (*pdf_iter)->updatePulls();
1001  }
1002 }
1003 
1004 void LauAbsFitModel::printFitParameters(const LauPdfList& pdfList, std::ostream& fout) const
1005 {
1006  LauPrint print;
1007  for (LauPdfList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1008  const LauAbsRValuePList& pars = (*pdf_iter)->getParameters();
1009  for (LauAbsRValuePList::const_iterator pars_iter = pars.begin(); pars_iter != pars.end(); ++pars_iter) {
1010  LauParameterPList params = (*pars_iter)->getPars();
1011  for (LauParameterPList::iterator params_iter = params.begin(); params_iter != params.end(); ++params_iter) {
1012  if (!(*params_iter)->clone()) {
1013  fout << (*params_iter)->name() << " & $";
1014  print.printFormat(fout, (*params_iter)->value());
1015  if ((*params_iter)->fixed() == kTRUE) {
1016  fout << "$ (fixed) \\\\";
1017  } else {
1018  fout << " \\pm ";
1019  print.printFormat(fout, (*params_iter)->error());
1020  fout << "$ \\\\" << std::endl;
1021  }
1022  }
1023  }
1024  }
1025  }
1026 }
1027 
1029 {
1030  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1031  (*pdf_iter)->cacheInfo(theData);
1032  }
1033 }
1034 
1035 Double_t LauAbsFitModel::prodPdfValue(LauPdfList& pdfList, UInt_t iEvt)
1036 {
1037  Double_t pdfVal = 1.0;
1038  for (LauPdfList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end(); ++pdf_iter) {
1039  (*pdf_iter)->calcLikelihoodInfo(iEvt);
1040  pdfVal *= (*pdf_iter)->getLikelihood();
1041  }
1042  return pdfVal;
1043 }
1044 
1045 void LauAbsFitModel::printEventInfo(UInt_t iEvt) const
1046 {
1047  const LauFitData& data = inputFitData_->getData(iEvt);
1048  std::cerr << " : Input data values for this event:" << std::endl;
1049  for (LauFitData::const_iterator iter = data.begin(); iter != data.end(); ++iter) {
1050  std::cerr << " " << iter->first << " = " << iter->second << std::endl;
1051  }
1052 }
1053 
1055 {
1056  std::cerr << " : Current values of fit parameters:" << std::endl;
1057  for (UInt_t i(0); i<nParams_; ++i) {
1058  std::cerr << " " << (fitVars_[i]->name()).Data() << " = " << fitVars_[i]->value() << std::endl;
1059  }
1060 }
1061 
1063 {
1064  // Update initial fit parameters if required (e.g. if using random numbers).
1065  this->checkInitFitParams();
1066 
1067  // Store the total number of parameters and the number of free parameters
1068  nParams_ = fitVars_.size();
1069  nFreeParams_ = 0;
1070 
1071  // Send the fit parameters
1072  for ( LauParameterPList::iterator iter = fitVars_.begin(); iter != fitVars_.end(); ++iter ) {
1073  if ( ! (*iter)->fixed() ) {
1074  ++nFreeParams_;
1075  }
1076  array.Add( *iter );
1077  }
1078 }
1079 
1080 void LauAbsFitModel::finaliseResults( const Int_t fitStat, const Double_t NLL, const TObjArray* parsFromMaster, const TMatrixD* covMat, TObjArray& parsToMaster )
1081 {
1082  // Copy the fit status information
1083  fitStatus_ = fitStat;
1084  NLL_ = NLL;
1085 
1086  // Copy the contents of the covariance matrix
1087  covMatrix_.Clear();
1088  covMatrix_.ResizeTo( covMat->GetNrows(), covMat->GetNcols() );
1089  covMatrix_.SetMatrixArray( covMat->GetMatrixArray() );
1090 
1091  // Now process the parameters
1092  UInt_t nPars = parsFromMaster->GetEntries();
1093  if ( nPars != nParams_ ) {
1094  std::cerr << "ERROR in LauAbsFitModel::finaliseResults : Unexpected number of parameters received from master" << std::endl;
1095  std::cerr << " : Received " << nPars << " when expecting " << nParams_ << std::endl;
1096  gSystem->Exit( EXIT_FAILURE );
1097  }
1098 
1099  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
1100  LauParameter* parameter = dynamic_cast<LauParameter*>( (*parsFromMaster)[iPar] );
1101  if ( ! parameter ) {
1102  std::cerr << "ERROR in LauAbsFitModel::finaliseResults : Error reading parameter from master" << std::endl;
1103  gSystem->Exit( EXIT_FAILURE );
1104  }
1105 
1106  if ( parameter->name() != fitVars_[iPar]->name() ) {
1107  std::cerr << "ERROR in LauAbsFitModel::finaliseResults : Error reading parameter from master" << std::endl;
1108  gSystem->Exit( EXIT_FAILURE );
1109  }
1110 
1111  *(fitVars_[iPar]) = *parameter;
1112  }
1113 
1114  // Write the results into the ntuple
1116 
1117  // Store the per-event likelihood values
1118  if ( this->writeSPlotData() ) {
1119  this->storePerEvtLlhds();
1120  }
1121 
1122  // Create a toy MC sample using the fitted parameters so that
1123  // the user can compare the fit to the data.
1124  if (compareFitData_ == kTRUE && fitStatus_ == 3) {
1126  }
1127 
1128  // Send the finalised fit parameters
1129  for ( LauParameterPList::iterator iter = fitVars_.begin(); iter != fitVars_.end(); ++iter ) {
1130  parsToMaster.Add( *iter );
1131  }
1132 }
1133 
1134 UInt_t LauAbsFitModel::readExperimentData( const UInt_t exptIndex )
1135 {
1136  // set our record of which experiment we're examining
1137  iExpt_ = exptIndex;
1138 
1139  // retrieve the data and find out how many events have been read
1141  UInt_t nEvent = inputFitData_->nEvents();
1142  this->eventsPerExpt( nEvent );
1143 
1144  return nEvent;
1145 }
1146 
Class for defining the SPlot technique.
Definition: LauSPlot.hh:55
TString fitToyMCTableName_
The output table name for Toy MC.
virtual void addSPlotNtupleIntegerBranch(const TString &name)
Add a branch to the sPlot tree for storing an integer.
Bool_t findBranches()
Find all of the branches in the tree.
virtual Bool_t twoStageFit() const
Determine whether the two-stage fit is enabled.
virtual Bool_t genExpt()=0
The method that actually generates the toy MC events for the given experiment.
void createFitToyMC(const TString &mcFileName, const TString &tableFileName)
Create a toy MC sample from the fitted parameters.
UInt_t firstExpt() const
Obtain the number of the first experiment.
Double_t width_
The width of the Gaussian constraint to be applied.
Definition: LauFitObject.hh:84
File containing declaration of LauFitDataTree class.
virtual Double_t getEventSum() const =0
Returns the sum of the expected events over all hypotheses; used in the EML fit scenario.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
virtual void updateParameters()=0
Update the values and errors of the parameters based on the fit minimum.
UInt_t eventsPerExpt() const
Obtain the total number of events in the current experiment.
void setGenValues()
Make sure all parameters hold their genValue as the current value.
TString sPlotFileName_
The name of the sPlot file.
Double_t mean_
The mean value of the Gaussian constraint to be applied.
Definition: LauFitObject.hh:82
virtual void writeOutTable(const TString &outputFile)=0
Write the latex table.
Bool_t writeLatexTable() const
Determine whether writing out of the latex table is enabled.
Bool_t doSFit_
Option to perfom the sFit.
Bool_t writeSPlotData() const
Determine whether the sPlot data is to be written out.
TStopwatch timer_
The fit timer.
UInt_t nParams_
The number of fit parameters.
File containing declaration of LauParamFixed class.
virtual Double_t getTotNegLogLikelihood()
Calculates the total negative log-likelihood.
LauAbsRValuePList conVars_
Internal vectors of Gaussian parameters.
Bool_t storeDPEff_
Option to store DP efficiencies in the sPlot ntuple.
void cacheInfo(LauPdfList &pdfList, const LauFitDataTree &theData)
Have all PDFs in the list cache the data.
void runCalculations(const TString &option="q")
Method to calculate the sWeights and cN coeffs.
Definition: LauSPlot.cc:626
virtual void setGenNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the gen tree.
ClassImp(LauAbsCoeffSet)
virtual void setupBkgndVectors()=0
Method to set up the storage for background-related quantities called by setBkgndClassNames.
void processMasterRequests()
Listen for requests from the master and act accordingly.
LauGenNtuple * genNtuple_
The generated ntuple.
File containing declaration of LauAbsFitModel class.
Bool_t withinAsymErrorCalc_
Flag to indicate if the asymmetric error calculation (e.g. MINOS) is currently running.
Bool_t useAsymmFitErrors() const
Determine whether calculation of asymmetric errors is enabled.
const TString & name() const
The parameter name.
Class for defining combinations of fit parameter objects.
virtual void setupGenNtupleBranches()=0
Setup the generation ntuple branches.
UInt_t addFitParameters(LauPdfList &pdfList)
Add parameters of the PDFs in the list to the list of all fit parameters.
virtual void generate(const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileNameBase)
Generate toy MC.
virtual void setBkgndClassNames(const std::vector< TString > &names)
Setup the background class names.
Bool_t useDP() const
Is the Dalitz plot term in the likelihood.
Int_t buildIndex(const TString &majorName, const TString &minorName="0")
Create an index table using leaves of the tree.
void compareFitData(UInt_t toyMCScale=10, const TString &mcFileName="fitToyMC.root", const TString &tableFileName="fitToyMCTable.tex", Bool_t poissonSmearing=kTRUE)
Specify that a toy MC sample should be created for a successful fit to an experiment.
const TString & treeName() const
Ntuple tree name.
Definition: LauGenNtuple.hh:56
Int_t getIntegerBranchValue(const TString &name) const
Get value of an integer branch.
virtual void cacheInputSWeights()
Cache the value of the sWeights to be used in the sFit.
virtual Bool_t twoStageFit() const =0
Determine whether the two-stage fit is enabled.
virtual void setGenNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the gen tree.
void deleteAndRecreateTree()
Delete and recreate tree.
File containing declaration of LauAbsFitter class.
void addFriendTree(const TString &rootFileName, const TString &rootTreeName)
Add a friend tree.
virtual void storePerEvtLlhds()=0
Store the per-event likelihood values.
virtual void printEventInfo(UInt_t iEvt) const
Prints the values of all the fit variables for the specified event - useful for diagnostics.
Bool_t pdfsDependOnDP() const
Do any of the PDFs have a dependence on the DP?
virtual ~LauAbsFitModel()
Destructor.
Bool_t compareFitData_
Option to make toy from 1st successful experiment.
File containing declaration of LauPrint class.
Double_t prodPdfValue(LauPdfList &pdfList, UInt_t iEvt)
Calculate the product of the per-event likelihoods of the PDFs in the list.
virtual void addConstraint(const TString &formula, const std::vector< TString > &pars, const Double_t mean, const Double_t width)
Store constraint information for fit parameters.
Double_t getDoubleBranchValue(const TString &name) const
Get value of a double branch.
std::vector< LauAbsPdf * > LauPdfList
List of Pdfs.
TString formula_
The formula to be used in the LauFormulaPar.
Definition: LauFitObject.hh:78
std::map< TString, Double_t > LauFitData
Type for holding event data.
void fillBranches()
Fill branches in the ntuple.
virtual Int_t getGenNtupleIntegerBranchValue(const TString &name) const
Get the value of an integer branch in the gen tree.
virtual LauSPlot::NumbMap fixdSpeciesNames() const =0
Returns the names and yields of species that are fixed in the fit.
LauBkgndClassMap bkgndClassNames_
The background class names.
Class to define various output print commands.
Definition: LauPrint.hh:29
virtual void setupSPlotNtupleBranches()=0
Setup the branches of the sPlot tuple.
virtual void fillGenNtupleBranches()
Fill the gen tuple branches.
void fitSlave(const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileNameBase)
Slaves required when performing a simultaneous fit.
virtual void setSPlotNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the sPlot tree.
void printFitParameters(const LauPdfList &pdfList, std::ostream &fout) const
Print the fit parameters for all PDFs in the list.
virtual UInt_t readExperimentData(const UInt_t exptIndex)
Read in the data for the specified experiment.
void addGaussianConstraint(Double_t newGaussMean, Double_t newGaussWidth)
Add a Gaussian constraint (or modify an existing one)
void connectToMaster(const TString &addressMaster, const UInt_t portMaster)
Establish the connection to the master process.
virtual void initialise(LauFitObject *fitObj, const std::vector< LauParameter * > &parameters)=0
Initialise the fitter, setting the information on the parameters.
void readExperimentData(UInt_t iExpt)
Read events only for the given experiment.
virtual std::pair< Int_t, Double_t > minimise()=0
Perform the minimisation of the fit function.
std::vector< LauAbsRValue * > LauAbsRValuePList
List of parameter pointers.
virtual Double_t getTotEvtLikelihood(UInt_t iEvt)=0
Calculates the likelihood for a given event.
TStopwatch cumulTimer_
The total fit timer.
File containing declaration of LauFitter class.
virtual void propagateParUpdates()=0
This function (specific to each model) calculates anything that depends on the fit parameter values...
void addConParameters()
Add parameters to the list of Gaussian constrained parameters.
virtual void recalculateNormalisation()=0
Recalculate normalisation the signal DP model(s)
LauFitDataTree * inputFitData_
The input data.
Bool_t doSFit() const
Return the flag to store the status of using an sFit or not.
void run(const TString &applicationCode, const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileName="")
Start the toy generation / fitting.
std::vector< LauParameter * > LauParameterPList
List of parameter pointers.
const TString & treeName() const
Retrieve the tree name.
LauParameterList extraVars_
Extra variables that aren&#39;t in the fit but are stored in the ntuple.
Bool_t blind() const
The blinding state.
std::vector< StoreConstraints > storeCon_
Store the constraints for fit parameters until initialisation is complete.
TMatrixD covMatrix_
The fit covariance matrix.
void writeOutFitResults()
Write out fit results.
void addIntegerBranch(const TString &name)
Add integer branch to tree.
Definition: LauGenNtuple.cc:73
const TString & fileName() const
Ntuple file name.
Definition: LauGenNtuple.hh:50
Abstract interface to the fitting and toy MC model.
UInt_t numberBadFits_
The number of bad fits.
virtual void addGenNtupleIntegerBranch(const TString &name)
Add a branch to the gen tree for storing an integer.
virtual void prepareInitialParArray(TObjArray &array)
Package the initial fit parameters for transmission to the master.
UInt_t nExpt() const
Obtain the number of experiments.
File containing declaration of LauParameter class.
Bool_t doPoissonSmearing() const
Determine whether Poisson smearing is enabled for the toy MC generation.
virtual Double_t getGenNtupleDoubleBranchValue(const TString &name) const
Get the value of a double branch in the gen tree.
void setNExpts(UInt_t nExperiments, UInt_t firstExperiment=0)
Set the number of experiments and the first experiment.
virtual void weightEvents(const TString &dataFileName, const TString &dataTreeName)=0
Reweighting - allows e.g. MC events to be weighted by the DP model.
Bool_t doEMLFit() const
Determine whether an extended maximum likelihood fit it being performed.
File containing declaration of LauComplex class.
UInt_t numberOKFits_
The number of good fits.
const TString & fileName() const
Retrieve the file name.
virtual Bool_t isDPDependent() const
Specifies whether or not the PDF is DP dependent.
Definition: LauAbsPdf.hh:111
virtual UInt_t nParameters() const =0
Get the total number of fit parameters.
virtual LauSPlot::NumbMap freeSpeciesNames() const =0
Returns the names and yields of species that are free in the fit.
Class to store the results from the toy MC generation into an ntuple.
Definition: LauGenNtuple.hh:32
File containing declaration of LauSPlot class.
virtual Bool_t splitSignal() const =0
Check if the signal is split into well-reconstructed and mis-reconstructed types. ...
Double_t worstLogLike_
The worst LL value found so far.
virtual void releaseSecondStageParameters()=0
Release parameters marked as &quot;second stage&quot;.
TString sPlotTreeName_
The name of the sPlot tree.
virtual void initialiseDPModels()=0
Initialise the DP models.
LauGenNtuple * sPlotNtuple_
The sPlot ntuple.
TString sPlotVerbosity_
Control the verbosity of the sFit.
const TString nullString_
An empty string.
void clearFitParVectors()
Clear the vectors containing fit parameters.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:35
virtual void addSPlotNtupleDoubleBranch(const TString &name)
Add a branch to the sPlot tree for storing a double.
void clearExtraVarVectors()
Clear the vectors containing extra ntuple variables.
virtual void cacheInputFitVars()=0
Cache the input data values to calculate the likelihood during the fit.
void fitExpt()
Routine to perform the actual fit for a given experiment.
virtual const std::vector< LauAbsRValue * > & getParameters() const
Retrieve the parameters of the PDF, e.g. so that they can be loaded into a fit.
Definition: LauAbsPdf.hh:239
Double_t sWeightScaleFactor_
The sWeight scaling factor.
File containing declaration of LauAbsPdf class.
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.
virtual void writeOutAllFitResults()
Write out any fit results.
UInt_t slaveId() const
Obtain the ID number of this slave.
virtual void finaliseFitResults(const TString &tableFileName)=0
Write the results of the fit into the ntuple.
UInt_t firstExpt_
The number of the first experiment to consider.
static LauAbsFitter * fitter()
Method that provides access to the singleton fitter.
Definition: LauFitter.cc:34
void setIntegerBranchValue(const TString &name, Int_t value)
Set value of an integer branch.
Definition: LauGenNtuple.cc:91
virtual void checkInitFitParams()=0
Update initial fit parameters if required.
virtual void printVarsInfo() const
Same as printEventInfo, but printing out the values of the variables in the fit.
Bool_t writeSPlotData_
Option to write sPlot data.
Bool_t validBkgndClass(const TString &className) const
Check if the given background class is in the list.
Bool_t fitToyMCPoissonSmear_
Option to perform Poisson smearing.
TString sWeightBranchName_
The name of the sWeight branch.
virtual UInt_t nFreeParameters() const =0
Get the number of floating fit parameters.
UInt_t iExpt_
The number of the current experiment.
LauParameterPList fitVars_
Internal vector of fit parameters.
virtual Bool_t scfDPSmear() const =0
Check if the mis-reconstructed signal is to be smeared in the DP.
Bool_t cacheFitData(const TString &dataFileName, const TString &dataTreeName)
Store variables from the input file into the internal data storage.
virtual void finaliseResults(const Int_t fitStat, const Double_t NLL, const TObjArray *parsFromMaster, const TMatrixD *covMat, TObjArray &parsToMaster)
Perform all finalisation actions.
void fit(const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileNameBase)
Perform the total fit.
virtual void setSPlotNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the sPlot tree.
UInt_t fitToyMCScale_
The scaling factor (toy vs data statistics)
virtual void fillSPlotNtupleBranches()
Fill the sPlot tuple.
TString fitToyMCFileName_
The output file name for Toy MC.
void writeOutGenResults()
Write out the results from the generation.
void printFormat(std::ostream &stream, Double_t value) const
Method to choose the printing format to a specified level of precision.
Definition: LauPrint.cc:32
void runSlave(const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileName="", const TString &addressMaster="localhost", const UInt_t portMaster=9090)
Start the slave process for simultaneous fitting.
virtual LauSPlot::TwoDMap twodimPDFs() const =0
Returns the species and variables for all 2D PDFs in the fit.
Double_t getLogLikelihoodPenalty()
Calculate the penalty terms to the log likelihood from Gaussian constraints.
Bool_t haveBranch(const TString &name) const
Check if the named branch is stored.
virtual void savePDFPlots(const TString &label)=0
Save the pdf Plots for all the resonances of experiment number fitExp.
const TString & bkgndClassName(UInt_t classID) const
Get the name of a background class from the number.
Struct to store constraint information until the fit is run.
Definition: LauFitObject.hh:76
LauParameterPSet resVars_
Internal set of fit parameters upon which the DP normalisation depends.
virtual Bool_t useAsymmFitErrors() const =0
Determine whether calculation of asymmetric errors is enabled.
Bool_t enableEmbedding() const
Determine whether embedding of events is enabled in the generation.
Double_t NLL_
The negative log-likelihood.
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:41
Double_t getLogLikelihood(UInt_t iStart, UInt_t iEnd)
Calculate the sum of the log-likelihood over the specified events.
UInt_t nExpt_
The number of experiments to consider.
LauFitNtuple * fitNtuple_
The fit ntuple.
virtual LauSPlot::NameSet variableNames() const =0
Returns the names of all variables in the fit.
Double_t value() const
The value of the parameter.
UInt_t nEvents() const
Retrieve the number of events.
virtual void addGenNtupleDoubleBranch(const TString &name)
Add a branch to the gen tree for storing a double.
friend class LauFitNtuple
LauFitNtuple is a friend class.
void updateFitParameters(LauPdfList &pdfList)
Update the fit parameters for all PDFs in the list.
void addDoubleBranch(const TString &name)
Add double branch to tree.
Definition: LauGenNtuple.cc:82
TString outputTableName_
The output table name.
virtual void setParsFromMinuit(Double_t *par, Int_t npar)
This function sets the parameter values from Minuit.
Int_t fitStatus_
The status of the fit.
virtual void initialise()=0
Initialise the fit par vectors.
std::vector< TString > conPars_
The list of LauParameter names to be used in the LauFormulaPar.
Definition: LauFitObject.hh:80
Class to store the input fit variables.
File containing declaration of LauFitNtuple class.
virtual const TMatrixD & covarianceMatrix() const =0
Retrieve the fit covariance matrix.
virtual void calculateSPlotData()
Calculate the sPlot data.
File containing declaration of LauGenNtuple class.
void setDoubleBranchValue(const TString &name, Double_t value)
Set value of a double branch.
Definition: LauGenNtuple.cc:96
std::vector< Double_t > sWeights_
The vector of sWeights.
UInt_t nFreeParams_
The number of free fit parameters.