laura is hosted by Hepforge, IPPP Durham
Laura++  v3r2
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauFitObject.hh
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2013 - 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 
19 #ifndef LAU_FIT_OBJECT
20 #define LAU_FIT_OBJECT
21 
22 #include <vector>
23 
24 #include "TMatrixD.h"
25 #include "TObject.h"
26 #include "TString.h"
27 
28 #include "LauAbsFitter.hh"
29 
30 class LauFitObject : public TObject {
31 
32  public:
34  virtual ~LauFitObject() {}
35 
37 
40  void useAsymmFitErrors(Bool_t useAsymmErrors) {useAsymmFitErrors_ = useAsymmErrors;}
41 
43  Bool_t useAsymmFitErrors() const {return useAsymmFitErrors_;}
44 
46 
56  void twoStageFit(Bool_t doTwoStageFit) {twoStageFit_ = doTwoStageFit;}
57 
59  Bool_t twoStageFit() const {return twoStageFit_;}
60 
62 
68  virtual void withinAsymErrorCalc(const Bool_t inAsymErrCalc) {withinAsymErrorCalc_ = inAsymErrCalc;}
69 
71 
74  Bool_t withinAsymErrorCalc() const {return withinAsymErrorCalc_;}
75 
77 
81  void setNExpts(UInt_t nExperiments, UInt_t firstExperiment = 0) {
82  nExpt_ = nExperiments;
83  firstExpt_ = firstExperiment;
84  }
85 
87  UInt_t eventsPerExpt() const {return evtsPerExpt_;}
88 
90  UInt_t nExpt() const {return nExpt_;}
91 
93  UInt_t firstExpt() const {return firstExpt_;}
94 
96  UInt_t iExpt() const {return iExpt_;}
97 
99 
106  virtual void setParsFromMinuit(Double_t* par, Int_t npar) = 0;
107 
109 
113  virtual Double_t getTotNegLogLikelihood() = 0;
114 
116 
122  virtual void addConstraint(const TString& formula, const std::vector<TString>& pars, const Double_t mean, const Double_t width);
123 
124  protected:
126  LauFitObject();
127 
128  // Setup a struct to store information on constrained fit parameters
135  TString formula_;
137  std::vector<TString> conPars_;
139  Double_t mean_;
141  Double_t width_;
142  };
143 
145  const std::vector<StoreConstraints>& constraintsStore() const {return storeCon_;}
146 
148  std::vector<StoreConstraints>& constraintsStore() {return storeCon_;}
149 
151  void resetFitCounters();
152 
154 
157  void setCurrentExperiment( const UInt_t curExpt ) { iExpt_ = curExpt; }
158 
160 
164  void startNewFit( const UInt_t nPars, const UInt_t nFreePars );
165 
167  void eventsPerExpt(UInt_t nEvents) {evtsPerExpt_ = nEvents;}
168 
170  Double_t worstLogLike() const {return worstLogLike_;}
171 
173 
176  void worstLogLike( const Double_t newWorstLogLike ) {worstLogLike_ = newWorstLogLike;}
177 
179 
183  void storeFitStatus( const LauAbsFitter::FitStatus& status, const TMatrixD& covMatrix );
184 
186  UInt_t nTotParams() const {return nParams_;}
187 
189  UInt_t nFreeParams() const {return nFreeParams_;}
190 
193 
195  Double_t nll() const {return fitStatus_.NLL;}
196 
198  Double_t edm() const {return fitStatus_.EDM;}
199 
201  Int_t statusCode() const {return fitStatus_.status;}
202 
204  const TMatrixD& covarianceMatrix() const {return covMatrix_;}
205 
207  UInt_t numberOKFits() const {return numberOKFits_;}
208 
210  UInt_t numberBadFits() const {return numberBadFits_;}
211 
212  private:
214  LauFitObject(const LauFitObject& rhs);
215 
217  LauFitObject& operator=(const LauFitObject& rhs);
218 
220  std::vector<StoreConstraints> storeCon_;
221 
223  Bool_t twoStageFit_;
224 
227 
229  UInt_t nParams_;
230 
232  UInt_t nFreeParams_;
233 
236 
238  UInt_t firstExpt_;
239 
241  UInt_t nExpt_;
242 
244  UInt_t iExpt_;
245 
247  UInt_t evtsPerExpt_;
248 
251 
253  Double_t worstLogLike_;
254 
256  TMatrixD covMatrix_;
257 
260 
263 
264  ClassDef(LauFitObject,0)
265 };
266 
267 #endif
268 
UInt_t iExpt_
The current experiment number.
UInt_t nTotParams() const
Access the total number of fit parameters.
Int_t statusCode() const
Access the fit status code.
Double_t width_
The width of the Gaussian constraint to be applied.
Double_t EDM
The estimated distance to the minimum.
Definition: LauAbsFitter.hh:47
Double_t mean_
The mean value of the Gaussian constraint to be applied.
virtual Double_t getTotNegLogLikelihood()=0
Calculate the new value of the negative log likelihood.
void eventsPerExpt(UInt_t nEvents)
Set the number of events in the current experiment.
Bool_t withinAsymErrorCalc() const
Query whether the fit is calculating the asymmetric errors.
Definition: LauFitObject.hh:74
void setNExpts(UInt_t nExperiments, UInt_t firstExperiment=0)
Set the number of experiments and the first experiment.
Definition: LauFitObject.hh:81
Double_t worstLogLike() const
Access the worst log likelihood found so far.
LauFitObject()
Constructor.
Definition: LauFitObject.cc:20
UInt_t evtsPerExpt_
The number of events in the current experiment.
std::vector< StoreConstraints > storeCon_
Store the constraints for fit parameters until initialisation is complete.
virtual void setParsFromMinuit(Double_t *par, Int_t npar)=0
This function sets the parameter values from Minuit.
File containing declaration of LauAbsFitter class.
virtual ~LauFitObject()
Destructor.
Definition: LauFitObject.hh:34
void storeFitStatus(const LauAbsFitter::FitStatus &status, const TMatrixD &covMatrix)
Store fit status information.
Definition: LauFitObject.cc:56
Bool_t twoStageFit() const
Report whether the two-stage fit is enabled.
Definition: LauFitObject.hh:59
TString formula_
The formula to be used in the LauFormulaPar.
Bool_t twoStageFit_
Option to perform a two stage fit.
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.
Definition: LauFitObject.cc:77
void useAsymmFitErrors(Bool_t useAsymmErrors)
Turn on or off the computation of asymmetric errors (e.g. MINOS routine in Minuit) ...
Definition: LauFitObject.hh:40
TMatrixD covMatrix_
The fit covariance matrix.
LauAbsFitter::FitStatus fitStatus_
The status of the current fit.
UInt_t nFreeParams() const
Access the total number of fit parameters.
void resetFitCounters()
Reset the good/bad fit counters.
Definition: LauFitObject.cc:39
UInt_t firstExpt_
The number of the first experiment to consider.
Double_t edm() const
Access the current EDM value.
UInt_t eventsPerExpt() const
Obtain the total number of events in the current experiment.
Definition: LauFitObject.hh:87
void setCurrentExperiment(const UInt_t curExpt)
Set the ID of the current experiment.
UInt_t numberOKFits_
The number of successful fits.
const std::vector< StoreConstraints > & constraintsStore() const
Const access to the constraints store.
Double_t NLL
The negative log-likelihood.
Definition: LauAbsFitter.hh:45
Int_t status
The status code of the fit.
Definition: LauAbsFitter.hh:43
Struct to store fit status information.
Definition: LauAbsFitter.hh:41
UInt_t numberBadFits_
The number of fit failures.
UInt_t nExpt() const
Obtain the number of experiments.
Definition: LauFitObject.hh:90
const TMatrixD & covarianceMatrix() const
Access the fit covariance matrix.
virtual void withinAsymErrorCalc(const Bool_t inAsymErrCalc)
Mark that the fit is calculating asymmetric errors.
Definition: LauFitObject.hh:68
UInt_t nParams_
The number of fit parameters.
Bool_t useAsymmFitErrors_
Option to use asymmetric errors.
Bool_t useAsymmFitErrors() const
Report whether or not calculation of asymmetric errors is enabled.
Definition: LauFitObject.hh:43
void twoStageFit(Bool_t doTwoStageFit)
Turn on or off the two stage fit.
Definition: LauFitObject.hh:56
UInt_t numberOKFits() const
Access the number of successful fits.
UInt_t numberBadFits() const
Access the number of failed fits.
Bool_t withinAsymErrorCalc_
Flag to indicate if the asymmetric error calculation (e.g. MINOS) is currently running.
Double_t nll() const
Access the current NLL value.
std::vector< StoreConstraints > & constraintsStore()
Access to the constraints store.
Struct to store constraint information until the fit is run.
void startNewFit(const UInt_t nPars, const UInt_t nFreePars)
Indicate the start of a new fit.
Definition: LauFitObject.cc:46
const LauAbsFitter::FitStatus & fitStatus() const
Access the fit status information.
void worstLogLike(const Double_t newWorstLogLike)
Set a new value for the worst log likelihood.
UInt_t nExpt_
The number of experiments to consider.
UInt_t iExpt() const
Obtain the number of the current experiment.
Definition: LauFitObject.hh:96
Double_t worstLogLike_
The worst log likelihood value found so far.
LauFitObject & operator=(const LauFitObject &rhs)
Copy assignment operator (not implemented)
UInt_t firstExpt() const
Obtain the number of the first experiment.
Definition: LauFitObject.hh:93
std::vector< TString > conPars_
The list of LauParameter names to be used in the LauFormulaPar.
The abstract interface for the objects that control the calculation of the likelihood.
Definition: LauFitObject.hh:30
UInt_t nFreeParams_
The number of free fit parameters.