laura is hosted by Hepforge, IPPP Durham
Laura++  v2r0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauAbsFitModel.hh
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 
47 #ifndef LAU_ABS_FIT_MODEL
48 #define LAU_ABS_FIT_MODEL
49 
50 #include "TMatrixD.h"
51 #include "TString.h"
52 #include "TStopwatch.h"
53 #include "TVectorDfwd.h"
54 
55 #include <vector>
56 #include <iosfwd>
57 
58 #include "LauFitObject.hh"
59 // LauSPlot included to get LauSPlot::NameSet typedef
60 #include "LauSPlot.hh"
61 
62 class TMessage;
63 class TMonitor;
64 class TSocket;
65 class TTree;
66 class LauAbsCoeffSet;
67 class LauAbsPdf;
68 class LauComplex;
69 class LauFitDataTree;
70 class LauFitNtuple;
71 class LauGenNtuple;
72 class LauParameter;
73 
74 class LauAbsFitModel : public LauFitObject {
75 
76  public:
79 
81  virtual ~LauAbsFitModel();
82 
84  Bool_t useDP() const { return usingDP_; }
85 
87 
90  void useDP(Bool_t usingDP) { usingDP_ = usingDP; }
91 
93  Bool_t doSFit() const { return doSFit_; }
94 
96 
100  void doSFit( const TString& sWeightBranchName, Double_t scaleFactor = 1.0 );
101 
103  Bool_t doEMLFit() const {return emlFit_;}
104 
106 
109  void doEMLFit(Bool_t emlFit) {emlFit_ = emlFit;}
110 
112  virtual Bool_t twoStageFit() const {return twoStageFit_;}
113 
115 
125  virtual void twoStageFit(Bool_t doTwoStageFit) {twoStageFit_ = doTwoStageFit;}
126 
128  Bool_t useAsymmFitErrors() const {return useAsymmFitErrors_;}
129 
131 
134  void useAsymmFitErrors(Bool_t useAsymmErrors) {useAsymmFitErrors_ = useAsymmErrors;}
135 
137  Bool_t doPoissonSmearing() const {return poissonSmear_;}
138 
140 
143  void doPoissonSmearing(Bool_t poissonSmear) {poissonSmear_ = poissonSmear;}
144 
146  Bool_t enableEmbedding() const {return enableEmbedding_;}
147 
149 
152  void enableEmbedding(Bool_t enable) {enableEmbedding_ = enable;}
153 
155 
161  virtual void withinAsymErrorCalc(Bool_t inAsymErrCalc) {withinAsymErrorCalc_ = inAsymErrCalc;}
162 
164  Bool_t writeLatexTable() const {return writeLatexTable_;}
165 
167 
170  void writeLatexTable(Bool_t writeTable) {writeLatexTable_ = writeTable;}
171 
173 
179  void writeSPlotData(const TString& fileName, const TString& treeName, Bool_t storeDPEfficiency, const TString& verbosity = "q");
180 
182  Bool_t writeSPlotData() const {return writeSPlotData_;}
183 
185  Bool_t storeDPEff() const {return storeDPEff_;}
186 
188  Bool_t useRandomInitFitPars() const {return randomFit_;}
189 
191  void useRandomInitFitPars(Bool_t boolean) {randomFit_ = boolean;}
192 
194 
198  void setNExpts(UInt_t nExperiments, UInt_t firstExperiment = 0) {
199  nExpt_ = nExperiments;
200  firstExpt_ = firstExperiment;
201  }
202 
204  UInt_t eventsPerExpt() const {return evtsPerExpt_;}
205 
207  UInt_t nExpt() const {return nExpt_;}
208 
210  UInt_t firstExpt() const {return firstExpt_;}
211 
213  UInt_t iExpt() const {return iExpt_;}
214 
216 
219  virtual void setBkgndClassNames( const std::vector<TString>& names );
220 
222  inline UInt_t nBkgndClasses() const {return bkgndClassNames_.size();}
223 
225 
228  virtual void setNSigEvents(LauParameter* nSigEvents) = 0;
229 
231 
237  virtual void setNBkgndEvents(LauParameter* nBkgndEvents) = 0;
238 
240 
243  virtual void setAmpCoeffSet(LauAbsCoeffSet* coeffSet) = 0;
244 
246 
254  void compareFitData(UInt_t toyMCScale = 10, const TString& mcFileName = "fitToyMC.root",
255  const TString& tableFileName = "fitToyMCTable.tex", Bool_t poissonSmearing = kTRUE);
256 
258 
262  virtual void weightEvents( const TString& dataFileName, const TString& dataTreeName ) = 0;
263 
265 
272  void run(const TString& applicationCode, const TString& dataFileName, const TString& dataTreeName,
273  const TString& histFileName, const TString& tableFileName = "");
274 
276 
284  void runSlave(const TString& dataFileName, const TString& dataTreeName,
285  const TString& histFileName, const TString& tableFileName = "",
286  const TString& addressMaster = "localhost", const UInt_t portMaster = 9090);
287 
289 
296  virtual void setParsFromMinuit(Double_t* par, Int_t npar);
297 
299 
303  Double_t getTotNegLogLikelihood();
304 
305  protected:
306 
307  // Some typedefs
308 
310  typedef std::vector<LauAbsPdf*> LauPdfList;
312  typedef std::vector<LauParameter*> LauParameterPList;
314  typedef std::vector<LauParameter> LauParameterList;
316  typedef std::map<UInt_t,TString> LauBkgndClassMap;
317 
319  void initSockets();
320 
322  void clearFitParVectors();
323 
325  void clearExtraVarVectors();
326 
328 
334  virtual void generate(const TString& dataFileName, const TString& dataTreeName, const TString& histFileName, const TString& tableFileNameBase);
335 
337 
340  virtual Bool_t genExpt() = 0;
341 
343 
349  void fit(const TString& dataFileName, const TString& dataTreeName, const TString& histFileName, const TString& tableFileNameBase);
350 
352 
358  void fitSlave(const TString& dataFileName, const TString& dataTreeName, const TString& histFileName, const TString& tableFileNameBase);
359 
361  void fitExpt();
362 
364 
367  Bool_t runMinimisation();
368 
370 
374  void createFitToyMC(const TString& mcFileName, const TString& tableFileName);
375 
377 
381  Bool_t cacheFitData(const TString& dataFileName, const TString& dataTreeName);
382 
384  virtual void cacheInputFitVars() = 0;
385 
387  virtual void cacheInputSWeights();
388 
390 
396  virtual void initialise() = 0;
397 
399  virtual void initialiseDPModels() = 0;
400 
408  virtual void updateCoeffs() = 0;
409 
411  virtual void propagateParUpdates() = 0;
412 
414 
418  Double_t getLogLikelihood( UInt_t iStart, UInt_t iEnd );
419 
421  Double_t getLogLikelihoodPenalty();
422 
424 
427  virtual Double_t getTotEvtLikelihood(UInt_t iEvt) = 0;
428 
430  virtual Double_t getEventSum() const = 0;
431 
433 
436  virtual void printEventInfo(UInt_t iEvt) const;
437 
439  virtual void printVarsInfo() const;
440 
442  virtual void checkInitFitParams() = 0;
443 
445 
448  virtual void finaliseFitResults(const TString& tableFileName) = 0;
449 
451 
454  virtual void writeOutTable(const TString& outputFile) = 0;
455 
457  virtual void storePerEvtLlhds() = 0;
458 
460  virtual void writeOutAllFitResults();
461 
463  virtual void calculateSPlotData();
464 
466  void setGenValues();
467 
469  virtual void setupBkgndVectors() = 0;
470 
472 
476  Bool_t validBkgndClass( const TString& className ) const;
477 
479 
483  UInt_t bkgndClassID( const TString& className ) const;
484 
486 
490  const TString& bkgndClassName( UInt_t classID ) const;
491 
493  void eventsPerExpt(UInt_t nEvents) {evtsPerExpt_ = nEvents;}
494 
496  virtual void setupGenNtupleBranches() = 0;
497 
499 
502  virtual void addGenNtupleIntegerBranch(const TString& name);
503 
505 
508  virtual void addGenNtupleDoubleBranch(const TString& name);
509 
511 
515  virtual void setGenNtupleIntegerBranchValue(const TString& name, Int_t value);
516 
518 
522  virtual void setGenNtupleDoubleBranchValue(const TString& name, Double_t value);
523 
525 
529  virtual Int_t getGenNtupleIntegerBranchValue(const TString& name) const;
530 
532 
536  virtual Double_t getGenNtupleDoubleBranchValue(const TString& name) const;
537 
539  virtual void fillGenNtupleBranches();
540 
542  virtual void setupSPlotNtupleBranches() = 0;
543 
545 
548  virtual void addSPlotNtupleIntegerBranch(const TString& name);
549 
551 
554  virtual void addSPlotNtupleDoubleBranch(const TString& name);
555 
557 
561  virtual void setSPlotNtupleIntegerBranchValue(const TString& name, Int_t value);
562 
564 
568  virtual void setSPlotNtupleDoubleBranchValue(const TString& name, Double_t value);
569 
571  virtual void fillSPlotNtupleBranches();
572 
574  virtual LauSPlot::NameSet variableNames() const = 0;
575 
577  virtual LauSPlot::NumbMap freeSpeciesNames() const = 0;
578 
580  virtual LauSPlot::NumbMap fixdSpeciesNames() const = 0;
581 
583  virtual LauSPlot::TwoDMap twodimPDFs() const = 0;
584 
586  virtual Bool_t splitSignal() const = 0;
587 
589  virtual Bool_t scfDPSmear() const = 0;
590 
592 
596  UInt_t addFitParameters(LauPdfList& pdfList);
597 
599  void addConParameters();
600 
602 
606  void printFitParameters(const LauPdfList& pdfList, std::ostream& fout) const;
607 
609 
612  void updateFitParameters(LauPdfList& pdfList);
613 
615 
619  void cacheInfo(LauPdfList& pdfList, const LauFitDataTree& theData);
620 
622 
626  Double_t prodPdfValue(LauPdfList& pdfList, UInt_t iEvt);
627 
629 
632  Bool_t pdfsDependOnDP() const {return pdfsDependOnDP_;}
633 
635 
638  void pdfsDependOnDP(Bool_t dependOnDP) { pdfsDependOnDP_ = dependOnDP; }
639 
641  const LauParameterPList& fitPars() const {return fitVars_;}
643 
645  const LauParameterList& extraPars() const {return extraVars_;}
647 
649  const LauParameterPList& conPars() const {return conVars_;}
651 
653  const LauFitNtuple* fitNtuple() const {return fitNtuple_;}
655 
657  const LauGenNtuple* genNtuple() const {return genNtuple_;}
659 
661  const LauGenNtuple* sPlotNtuple() const {return sPlotNtuple_;}
663 
665  const LauFitDataTree* fitData() const {return inputFitData_;}
667 
669  Double_t nll() const {return NLL_;}
670 
672  Int_t fitStatus() const {return fitStatus_;}
673 
675  const TMatrixD& covarianceMatrix() const {return covMatrix_;}
676 
677  private:
678  // Various control booleans
679 
681  Bool_t twoStageFit_;
691  Bool_t storeDPEff_;
693  Bool_t randomFit_;
695  Bool_t emlFit_;
697  Bool_t poissonSmear_;
701  Bool_t usingDP_;
704 
705  // Info on number of experiments and number of events
706 
708  UInt_t firstExpt_;
710  UInt_t nExpt_;
712  UInt_t evtsPerExpt_;
714  UInt_t iExpt_;
715 
718 
721 
724 
725  // Input data and output ntuple
726 
735 
736  // Fit bookeeping variables
737 
739  Int_t fitStatus_;
741  Double_t NLL_;
743  TMatrixD covMatrix_;
745  UInt_t numberOKFits_;
747  UInt_t numberBadFits_;
749  UInt_t nParams_;
751  UInt_t nFreeParams_;
753  Double_t worstLogLike_;
756 
757  // Background class names
761  const TString nullString_;
762 
763  // sFit related variables
764 
766  Bool_t doSFit_;
770  std::vector<Double_t> sWeights_;
773 
774  // Fit timers
775 
777  TStopwatch timer_;
779  TStopwatch cumulTimer_;
780 
781  // Comparison toy MC related variables
782 
791 
792  // sPlot related variables
793 
795  TString sPlotFileName_;
797  TString sPlotTreeName_;
800 
801  // Parallel fitting variables
802 
804  TSocket* sMaster_;
808  UInt_t slaveId_;
810  UInt_t nSlaves_;
812  Double_t* parValues_;
813 
814  ClassDef(LauAbsFitModel,0) // Abstract interface to fit/toyMC model
815 };
816 
817 #endif
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.
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.
std::vector< LauParameter > LauParameterList
List of parameters.
virtual Double_t getEventSum() const =0
Returns the sum of the expected events over all hypotheses; used in the EML fit scenario.
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 * parValues_
Parameter values array (for reading from the master)
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.
Double_t getTotNegLogLikelihood()
Calculates the total negative log-likelihood.
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.
virtual void setGenNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the gen tree.
void writeLatexTable(Bool_t writeTable)
Turn on or off the writing out of the latex table.
virtual void setupBkgndVectors()=0
Method to set up the storage for background-related quantities called by setBkgndClassNames.
LauGenNtuple * genNtuple_
The generated ntuple.
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.
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.
void doPoissonSmearing(Bool_t poissonSmear)
Turn Poisson smearing (for the toy MC generation) on or off.
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
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.
Bool_t randomFit_
Option to randomise the initial values of the fit parameters.
virtual void cacheInputSWeights()
Cache the value of the sWeights to be used in the sFit.
virtual void setGenNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the gen tree.
virtual void setNSigEvents(LauParameter *nSigEvents)=0
Set the number of signal events.
void eventsPerExpt(UInt_t nEvents)
Set the number of events in the current experiment.
TSocket * sMaster_
A socket to enable parallel setup.
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.
Int_t fitStatus() const
Access the fit status information.
Bool_t compareFitData_
Option to make toy from 1st successful experiment.
Bool_t writeLatexTable_
Option to output a Latex format table.
Double_t prodPdfValue(LauPdfList &pdfList, UInt_t iEvt)
Calculate the product of the per-event likelihoods of the PDFs in the list.
Bool_t emlFit_
Option to perform an extended ML fit.
LauFitDataTree * fitData()
std::vector< LauAbsPdf * > LauPdfList
List of Pdfs.
virtual Int_t getGenNtupleIntegerBranchValue(const TString &name) const
Get the value of an integer branch in the gen tree.
UInt_t nSlaves_
The total number of slaves.
Bool_t poissonSmear_
Option to perform Poisson smearing.
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.
TMessage * messageFromMaster_
Message from master to the slaves.
Bool_t twoStageFit_
Option to perform a two stage fit.
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.
Bool_t useAsymmFitErrors_
Option to use asymmetric errors.
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.
Bool_t enableEmbedding_
Option to enable embedding.
UInt_t slaveId_
Slave id number.
LauParameterList & extraPars()
LauFitNtuple * fitNtuple()
virtual Double_t getTotEvtLikelihood(UInt_t iEvt)=0
Calculates the likelihood for a given event.
TStopwatch cumulTimer_
The total fit timer.
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.
void initSockets()
Initialise socket connections for the slaves.
LauFitDataTree * inputFitData_
The input data.
const LauParameterList & extraPars() const
Access the extra variables.
void doEMLFit(Bool_t emlFit)
Choice to perform an extended maximum likelihood fit.
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.
LauParameterList extraVars_
Extra variables that aren&#39;t in the fit but are stored in the ntuple.
TMatrixD covMatrix_
The fit covariance matrix.
Abstract interface to the fitting and toy MC model.
UInt_t iExpt() const
Obtain the number of the current experiment.
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.
const LauFitNtuple * fitNtuple() const
Access the fit ntuple.
UInt_t nExpt() const
Obtain the number of experiments.
virtual void withinAsymErrorCalc(Bool_t inAsymErrCalc)
Mark that the fit is calculating asymmetric errors.
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 setAmpCoeffSet(LauAbsCoeffSet *coeffSet)=0
Set the DP amplitude coefficients.
virtual void weightEvents(const TString &dataFileName, const TString &dataTreeName)=0
Reweighting - allows e.g. MC events to be weighted by the DP model.
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.
std::map< UInt_t, TString > LauBkgndClassMap
A type to store background classes.
virtual void twoStageFit(Bool_t doTwoStageFit)
Turn on or off the two stage fit.
UInt_t numberOKFits_
The number of good fits.
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. ...
LauParameterPList & fitPars()
Double_t worstLogLike_
The worst LL value found so far.
TString sPlotTreeName_
The name of the sPlot tree.
virtual void initialiseDPModels()=0
Initialise the DP models.
Class for defining the abstract interface for complex coefficient classes.
Double_t nll() const
Access the current NLL value.
LauAbsFitModel()
Constructor.
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:32
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 setNBkgndEvents(LauParameter *nBkgndEvents)=0
Set the number of background events.
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.
const LauGenNtuple * genNtuple() const
Access the gen ntuple.
LauGenNtuple * sPlotNtuple()
LauGenNtuple * genNtuple()
Double_t sWeightScaleFactor_
The sWeight scaling factor.
Bool_t storeDPEff() const
Determine whether the efficiency information should be stored in the sPlot ntuple.
const LauGenNtuple * sPlotNtuple() const
Access the sPlot ntuple.
UInt_t bkgndClassID(const TString &className) const
The number assigned to a background class.
virtual void writeOutAllFitResults()
Write out any fit results.
Bool_t useRandomInitFitPars() const
Determine whether the initial values of the fit parameters, in particular the isobar coefficient para...
virtual void finaliseFitResults(const TString &tableFileName)=0
Write the results of the fit into the ntuple.
const LauFitDataTree * fitData() const
Access the data store.
UInt_t firstExpt_
The number of the first experiment to consider.
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.
UInt_t iExpt_
The number of the current experiment.
LauParameterPList fitVars_
Internal vectors 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.
void fit(const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileNameBase)
Perform the total fit.
void enableEmbedding(Bool_t enable)
Turn on or off embedding of events in the generation.
void pdfsDependOnDP(Bool_t dependOnDP)
Do any of the PDFs have a dependence on the DP?
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)
Bool_t runMinimisation()
Routine to perform the minimisation.
File containing declaration of LauFitObject class.
virtual void fillSPlotNtupleBranches()
Fill the sPlot tuple.
virtual void updateCoeffs()=0
TString fitToyMCFileName_
The output file name for Toy MC.
Bool_t usingDP_
Option to include the DP as part of the fit.
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:42
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.
LauParameterPList & conPars()
Double_t getLogLikelihoodPenalty()
Calculate the penalty terms to the log likelihood from Gaussian constraints.
Class for defining a complex number.
Definition: LauComplex.hh:47
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.
void useRandomInitFitPars(Bool_t boolean)
Randomise the initial values of the fit parameters, in particular the isobar coefficient parameters...
void useDP(Bool_t usingDP)
Switch on/off the Dalitz plot term in the Likelihood (allows fits to other quantities, e.g. B mass)
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:40
Double_t getLogLikelihood(UInt_t iStart, UInt_t iEnd)
Calculate the sum of the log-likelihood over the specified events.
const TMatrixD & covarianceMatrix() const
Access the fit covariance matrix.
const LauParameterPList & fitPars() const
Access the fit variables.
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.
UInt_t evtsPerExpt_
The number of events per experiment.
virtual void addGenNtupleDoubleBranch(const TString &name)
Add a branch to the gen tree for storing a double.
void useAsymmFitErrors(Bool_t useAsymmErrors)
Turn on or off the computation of asymmetric errors (e.g. MINOS routine in Minuit) ...
void updateFitParameters(LauPdfList &pdfList)
Update the fit parameters for all PDFs in the list.
LauParameterPList conVars_
Internal vectors of Gaussian parameters.
const LauParameterPList & conPars() const
Access the Gaussian constrained variables.
Bool_t pdfsDependOnDP_
Option to state if pdfs depend on DP position.
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.
Class to store the input fit variables.
std::set< TString > NameSet
Type to store names, e.g. of the discriminating/control variables.
Definition: LauSPlot.hh:59
The abstract interface for the objects that control the calculation of the likelihood.
Definition: LauFitObject.hh:26
virtual void calculateSPlotData()
Calculate the sPlot data.
std::vector< Double_t > sWeights_
The vector of sWeights.
UInt_t nFreeParams_
The number of free fit parameters.