laura is hosted by Hepforge, IPPP Durham
Laura++  v3r2
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauFitNtuple.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 
17 #include "TFile.h"
18 #include "TMath.h"
19 #include "TMatrixD.h"
20 #include "TSystem.h"
21 #include "TTree.h"
22 
23 #include "LauAbsFitter.hh"
24 #include "LauFitNtuple.hh"
25 #include "LauFitter.hh"
26 #include "LauParameter.hh"
27 #include "LauParamFixed.hh"
28 
30 
31 
32 LauFitNtuple::LauFitNtuple(const TString& fileName, Bool_t storeAsymErrors) :
33  rootFileName_(fileName),
34  rootFile_(0),
35  fitResults_(0),
36  definedFitTree_(kFALSE),
37  storeAsymErrors_(storeAsymErrors),
38  fitStatus_({-1,0.0,0.0}),
39  nFitPars_(0),
40  nFreePars_(0),
41  nExtraPars_(0),
42  iExpt_(0)
43 {
44  rootFile_ = TFile::Open(rootFileName_, "recreate");
45  rootFile_->cd();
46  fitResults_ = new TTree("fitResults", "fitResults");
47  fitResults_->SetDirectory(rootFile_);
48 
49  fitVars_.clear(); extraVars_.clear();
50 }
51 
53 {
54  if (rootFile_ && rootFile_->IsOpen()) {
55  delete fitResults_; fitResults_ = 0;
56  }
57  delete rootFile_; rootFile_ = 0;
58 }
59 
60 void LauFitNtuple::storeCorrMatrix(const UInt_t iExpt, const LauAbsFitter::FitStatus& fitStatus, const TMatrixD& covMatrix)
61 {
62  // store the minimised NLL value, correlation matrix status and experiment number
63  iExpt_ = iExpt;
64  fitStatus_ = fitStatus;
65 
66  // make the correlation matrix the correct dimensions
67  if (definedFitTree_ == kFALSE) {
68  corrMatrix_.clear();
69  corrMatrix_.resize(nFitPars_);
70  for (UInt_t i = 0; i < nFitPars_; ++i) {corrMatrix_[i].resize(nFitPars_);}
71  }
72 
73  // under certain circumstances, e.g. if the fit has failed in the first
74  // stage of a two-stage fit, the covariance matrix might not have the
75  // expected dimensions, or it might even be empty
76  Bool_t needsPadding = kFALSE;
77  const UInt_t nElements = covMatrix.GetNoElements();
78  if ( nElements == 0 ) {
79 
80  // if it's empty we can just make a diagonal matrix
81  std::cerr << "WARNING in LauFitNtuple::storeCorrMatrix : received empty covariance matrix - will store diagonal correlation matrix" << std::endl;
82  for (UInt_t i = 0; i < nFitPars_; ++i) {
83  for (UInt_t j = 0; j < nFitPars_; ++j) {
84  if (i == j) {
85  corrMatrix_[i][j] = 1.0;
86  } else {
87  corrMatrix_[i][j] = 0.0;
88  }
89  }
90  }
91  return;
92 
93  } else if ( nElements != nFreePars_*nFreePars_ ) {
94 
95  UInt_t dimension = covMatrix.GetNrows();
96  UInt_t nSecondStage = 0;
97  for (UInt_t i = 0; i < nFitPars_; ++i) {
98  if ( fitVars_[i]->secondStage() ) {
99  ++nSecondStage;
100  }
101  }
102 
103  if ( (dimension + nSecondStage) == nFreePars_ ) {
104  needsPadding = kTRUE;
105  std::cerr << "WARNING in LauFitNtuple::storeCorrMatrix : received smaller covariance matrix than expected, likely due to failure of first stage fit - will pad the correlation matrix" << std::endl;
106  std::cerr << " : nFitPars_ = " << nFitPars_ << std::endl;
107  std::cerr << " : nFreePars_ = " << nFreePars_ << std::endl;
108  std::cerr << " : nSecondStage = " << nSecondStage << std::endl;
109  std::cerr << " : covMatrix size = " << dimension << std::endl;
110  } else {
111  std::cerr << "WARNING in LauFitNtuple::storeCorrMatrix : received smaller covariance matrix than expected, for unknown reasons - will store diagonal correlation matrix" << std::endl;
112  std::cerr << " : nFitPars_ = " << nFitPars_ << std::endl;
113  std::cerr << " : nFreePars_ = " << nFreePars_ << std::endl;
114  std::cerr << " : nSecondStage = " << nSecondStage << std::endl;
115  std::cerr << " : covMatrix size = " << dimension << std::endl;
116  for (UInt_t i = 0; i < nFitPars_; ++i) {
117  for (UInt_t j = 0; j < nFitPars_; ++j) {
118  if (i == j) {
119  corrMatrix_[i][j] = 1.0;
120  } else {
121  corrMatrix_[i][j] = 0.0;
122  }
123  }
124  }
125  return;
126  }
127 
128  }
129 
130  // calculate the correlation matrix information from the fitter covariance matrix
131  Bool_t iFixed(kFALSE);
132  Bool_t jFixed(kFALSE);
133  Bool_t iSecondStage(kFALSE);
134  Bool_t jSecondStage(kFALSE);
135  UInt_t iFree(0);
136  UInt_t jFree(0);
137 
138  for (UInt_t i = 0; i < nFitPars_; ++i) {
139 
140  iFixed = fitVars_[i]->fixed();
141  iSecondStage = fitVars_[i]->secondStage();
142 
143  // reset the "j" free parameter counter
144  jFree = 0;
145 
146  // NB the supplied covariance matrix is of dimension nFreePars_ x nFreePars_
147  for (UInt_t j = 0; j < nFitPars_; ++j) {
148 
149  jFixed = fitVars_[j]->fixed();
150  jSecondStage = fitVars_[j]->secondStage();
151 
152  if (i == j) {
153  corrMatrix_[i][j] = 1.0;
154  } else if (iFixed == kTRUE || jFixed == kTRUE) {
155  corrMatrix_[i][j] = 0.0;
156  } else if ( needsPadding && ( iSecondStage || jSecondStage ) ) {
157  corrMatrix_[i][j] = 0.0;
158  } else {
159  Double_t r_ij = covMatrix(iFree,jFree);
160  Double_t r_ii = covMatrix(iFree,iFree);
161  Double_t r_jj = covMatrix(jFree,jFree);
162  Double_t denom = r_ii * r_jj;
163  if (denom < 0.0) {
164  r_ij = 0.0;
165  denom = 1.0;
166  }
167  denom = TMath::Sqrt(denom);
168  if (denom < 1e-30) {
169  r_ij = 0.0;
170  denom = 1.0;
171  }
172  corrMatrix_[i][j] = r_ij / denom;
173  }
174 
175  if ( !jFixed && !(needsPadding && jSecondStage) ) {
176  ++jFree;
177  }
178  }
179 
180  if ( !iFixed && !(needsPadding && iSecondStage) ) {
181  ++iFree;
182  }
183  }
184 }
185 
186 void LauFitNtuple::storeParsAndErrors(const std::vector<LauParameter*>& fitVars, const std::vector<LauParameter>& extraVars)
187 {
188  fitVars_ = fitVars;
189  UInt_t nFitPars = fitVars_.size();
190 
191  // the number of parameters being given to us should be the same as the number from the last fit
192  // OR it's the first time so the "last" number is zero
193  if (nFitPars_ != 0 && nFitPars_ != nFitPars) {
194  std::cerr << "ERROR in LauFitNtuple::storeParsAndErrors : expected total number of parameters (" << nFitPars_
195  << ") not the same as the number provided (" << nFitPars << ")." << std::endl;
196  gSystem->Exit(EXIT_FAILURE);
197  }
198 
199  LauParamFixed pred;
200  UInt_t nFreePars = nFitPars - std::count_if(fitVars_.begin(),fitVars_.end(),pred);
201 
202  // the number of free parameters being given to us should be the same as the number from the last fit
203  // OR it's the first time so the "last" number is zero
204  // (NB we check whether nFitPars_ is zero for this since it is possible to have zero free parameters, albeit rather daft)
205  if (nFitPars_ != 0 && nFreePars_ != nFreePars) {
206  std::cerr << "ERROR in LauFitNtuple::storeParsAndErrors : expected number of free parameters (" << nFreePars_
207  << ") not the same as the number provided (" << nFreePars << ")." << std::endl;
208  gSystem->Exit(EXIT_FAILURE);
209  }
210 
211  nFitPars_ = nFitPars;
212  nFreePars_ = nFreePars;
213 
214  extraVars_ = extraVars;
215  nExtraPars_ = extraVars_.size();
216 }
217 
219 {
220  // Now create and fill the stored fit results into an ntuple (TTree)
221  if (definedFitTree_ == kFALSE) {
222 
223  std::cout << "INFO in LauFitNtuple::updateFitNtuple : totNoPars = " << nFitPars_ << std::endl;
224 
225  // Add experiment number as a branch
226  fitResults_->Branch("iExpt", &iExpt_, "iExpt/I");
227  fitResults_->Branch("fitStatus", &fitStatus_.status, "fitStatus/I");
228 
229  // Add NLL (negative log-likelihood) and EDM values from fit
230  fitResults_->Branch("NLL", &fitStatus_.NLL, "NLL/D");
231  fitResults_->Branch("EDM", &fitStatus_.EDM, "EDM/D");
232 
233  for (UInt_t i = 0; i < nFitPars_; i++) {
234 
235  TString parName = fitVars_[i]->name();
236  TString parNameD(parName); parNameD += "/D";
237  fitResults_->Branch(parName.Data(), &fitVars_[i]->value_, parNameD.Data());
238 
239  TString parInitName(parName); parInitName += "_True";
240  TString parInitNameD(parInitName); parInitNameD += "/D";
241  fitResults_->Branch(parInitName.Data(), &fitVars_[i]->genValue_, parInitNameD.Data());
242 
243  if (!fitVars_[i]->fixed()) {
244  TString parErrName(parName); parErrName += "_Error";
245  TString parErrNameD(parErrName); parErrNameD += "/D";
246  fitResults_->Branch(parErrName.Data(), &fitVars_[i]->error_, parErrNameD.Data());
247 
248  if ( storeAsymErrors_ ) {
249  TString parNegErrName(parName); parNegErrName += "_NegError";
250  TString parNegErrNameD(parNegErrName); parNegErrNameD += "/D";
251  fitResults_->Branch(parNegErrName.Data(), &fitVars_[i]->negError_, parNegErrNameD.Data());
252 
253  TString parPosErrName(parName); parPosErrName += "_PosError";
254  TString parPosErrNameD(parPosErrName); parPosErrNameD += "/D";
255  fitResults_->Branch(parPosErrName.Data(), &fitVars_[i]->posError_, parPosErrNameD.Data());
256  }
257 
258  TString parPullName(parName); parPullName += "_Pull";
259  TString parPullNameD(parPullName); parPullNameD += "/D";
260  fitResults_->Branch(parPullName.Data(), &fitVars_[i]->pull_, parPullNameD.Data());
261  }
262 
263  // Now add in the correlation matrix values (only for floating parameters)
264  if (!fitVars_[i]->fixed()) {
265  // First the global correlation coeffs
266  TString parGCCName(parName); parGCCName += "_GCC";
267  TString parGCCNameD(parGCCName); parGCCNameD += "/D";
268  fitResults_->Branch(parGCCName.Data(), &fitVars_[i]->gcc_, parGCCNameD.Data());
269 
270  if ( ! corrMatrix_.empty() ) {
271  // Then the rest
272  for (UInt_t j = 0; j < nFitPars_; j++) {
273  if (!fitVars_[j]->fixed() && i!=j) {
274  TString parName2 = fitVars_[j]->name();
275  TString corrName("corr__");
276  corrName += parName; corrName += "__"; corrName += parName2;
277 
278  TString corrNameD(corrName); corrNameD += "/D";
279  fitResults_->Branch(corrName.Data(), &corrMatrix_[i][j], corrNameD.Data());
280  }
281  }
282  }
283  }
284  }
285 
286  // Update extra parameter values...
287  for (UInt_t i = 0; i < nExtraPars_; i++) {
288 
289  TString parName = extraVars_[i].name();
290  TString parNameD(parName); parNameD += "/D";
291  fitResults_->Branch(parName.Data(), &extraVars_[i].value_, parNameD.Data());
292 
293  TString parInitName(parName); parInitName += "_True";
294  TString parInitNameD(parInitName); parInitNameD += "/D";
295  fitResults_->Branch(parInitName.Data(), &extraVars_[i].genValue_, parInitNameD.Data());
296 
297  //TString parErrName(parName); parErrName += "_Error";
298  //TString parErrNameD(parErrName); parErrNameD += "/D";
299  //fitResults_->Branch(parErrName.Data(), &extraVars_[i].error_, parErrNameD.Data());
300 
301  // Also find the fit fraction pull and store it
302  //TString pullName(parName); pullName += "_Pull";
303  //TString pullNameD(pullName); pullNameD += "/D";
304  //fitResults_->Branch(pullName.Data(), &extraVars_[i].pull_, pullNameD.Data());
305 
306  }
307 
308  definedFitTree_ = kTRUE;
309 
310  }
311 
312  std::cout << "INFO in LauFitNtuple::updateFitNtuple : Stored fit values for experiment " << iExpt_ << "\n";
313  std::cout << " : fitStatus = " << fitStatus_.status << std::endl;
314  std::cout << " : NLL = " << fitStatus_.NLL << std::endl;
315  std::cout << " : EDM = " << fitStatus_.EDM << std::endl;
316 
317  fitResults_->Fill();
318 }
319 
321 {
322  // Write out the fit ntuple to the appropriate root file
323  rootFile_->cd();
324  fitResults_->Write("",TObject::kOverwrite);
325  rootFile_->Close();
326  delete rootFile_; rootFile_ = 0;
327 }
328 
Bool_t definedFitTree_
Flags whether the fit tree has been defined.
std::vector< LauParameter * > fitVars_
Fit variables.
Definition: LauFitNtuple.hh:93
Double_t EDM
The estimated distance to the minimum.
Definition: LauAbsFitter.hh:47
Bool_t fixed() const
Check whether the parameter is fixed or floated.
File containing declaration of LauParamFixed class.
nExtraPars_(0)
Int_t iExpt_
Experiment number.
ClassImp(LauAbsCoeffSet)
UInt_t nExtraPars_
Number of extra parameters.
File containing declaration of LauAbsFitter class.
LauAbsFitter::FitStatus fitStatus_
Status of fit.
nFitPars_(0)
std::vector< std::vector< Double_t > > corrMatrix_
Correlation matrix.
Definition: LauFitNtuple.hh:98
std::vector< LauParameter > extraVars_
Extra variables.
Definition: LauFitNtuple.hh:95
void updateFitNtuple()
Update the fit ntuple.
TFile * rootFile_
Root file.
Definition: LauFitNtuple.hh:88
File containing declaration of LauFitter class.
UInt_t nFitPars_
Number of fit parameters.
void storeCorrMatrix(const UInt_t iExpt, const LauAbsFitter::FitStatus &fitStatus, const TMatrixD &covMatrix)
Store the correlation matrix and other fit information.
Definition: LauFitNtuple.cc:60
void writeOutFitResults()
Write out fit results.
nFreePars_(0)
Double_t NLL
The negative log-likelihood.
Definition: LauAbsFitter.hh:45
Int_t status
The status code of the fit.
Definition: LauAbsFitter.hh:43
File containing declaration of LauParameter class.
Struct to store fit status information.
Definition: LauAbsFitter.hh:41
Predicate to allow counting of the number of fixed parameters.
Bool_t secondStage() const
Check whether the parameter should be floated only in the second stage of a two stage fit...
virtual ~LauFitNtuple()
Destructor.
Definition: LauFitNtuple.cc:52
void storeParsAndErrors(const std::vector< LauParameter * > &fitVars, const std::vector< LauParameter > &extraVars)
Store parameters and their errors.
TTree * fitResults_
Fit results.
Definition: LauFitNtuple.hh:90
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:44
iExpt_(0)
Definition: LauFitNtuple.cc:42
UInt_t nFreePars_
Number of free parameters.
File containing declaration of LauFitNtuple class.
Bool_t storeAsymErrors_
Flags whether or not to store the asymmetric error information.