laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauRooFitTask.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2017 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
29 #include "LauRooFitTask.hh"
30 
31 #include "LauFitNtuple.hh"
32 #include "LauParameter.hh"
33 #include "LauSimFitTask.hh"
34 
35 #include "RooAbsPdf.h"
36 #include "RooDataSet.h"
37 #include "RooFormulaVar.h"
38 #include "RooNLLVar.h"
39 #include "RooRealVar.h"
40 #include "TFile.h"
41 #include "TString.h"
42 #include "TSystem.h"
43 #include "TTree.h"
44 
45 #include <iostream>
46 #include <vector>
47 
48 LauRooFitTask::LauRooFitTask( RooAbsPdf& model,
49  const Bool_t extended,
50  const RooArgSet& vars,
51  const TString& weightVarName ) :
52  LauSimFitTask(),
53  model_( model ),
54  dataVars_( vars ),
55  weightVarName_( weightVarName ),
56  dataFile_( 0 ),
57  dataTree_( 0 ),
58  exptData_( 0 ),
59  extended_( extended ),
60  iExptSet_(),
61  nllVar_( 0 )
62 {
63 }
64 
66 {
67  delete nllVar_;
68  nllVar_ = 0;
69  this->cleanData();
70 }
71 
73 {
74  if ( dataFile_ != 0 ) {
75  dataFile_->Close();
76  delete dataFile_;
77  dataTree_ = 0;
78  dataFile_ = 0;
79  }
80  delete exptData_;
81  exptData_ = 0;
82 }
83 
85 {
86  if ( weightVarName_ != "" ) {
87  Bool_t weightVarFound = kFALSE;
88  RooFIter argset_iter = dataVars_.fwdIterator();
89  RooAbsArg* param( 0 );
90  while ( ( param = argset_iter.next() ) ) {
91  TString name = param->GetName();
92  if ( name == weightVarName_ ) {
93  weightVarFound = kTRUE;
94  break;
95  }
96  }
97  if ( ! weightVarFound ) {
98  std::cerr << "ERROR in LauRooFitTask::initialise : The set of data variables does not contain the weighting variable \""
99  << weightVarName_ << std::endl;
100  std::cerr << " : Weighting will be disabled."
101  << std::endl;
102  weightVarName_ = "";
103  }
104  }
105 }
106 
107 Bool_t LauRooFitTask::verifyFitData( const TString& dataFileName, const TString& dataTreeName )
108 {
109  // Clean-up from any previous runs
110  if ( dataFile_ != 0 ) {
111  this->cleanData();
112  }
113 
114  // Open the data file
115  dataFile_ = TFile::Open( dataFileName );
116  if ( ! dataFile_ ) {
117  std::cerr << "ERROR in LauRooFitTask::verifyFitData : Problem opening data file \""
118  << dataFileName << "\"" << std::endl;
119  return kFALSE;
120  }
121 
122  // Retrieve the tree
123  dataTree_ = dynamic_cast<TTree*>( dataFile_->Get( dataTreeName ) );
124  if ( ! dataTree_ ) {
125  std::cerr << "ERROR in LauRooFitTask::verifyFitData : Problem retrieving tree \""
126  << dataTreeName << "\" from data file \"" << dataFileName << "\"" << std::endl;
127  dataFile_->Close();
128  delete dataFile_;
129  dataFile_ = 0;
130  return kFALSE;
131  }
132 
133  // Check that the tree contains branches for all the fit variables
134  RooFIter argset_iter = dataVars_.fwdIterator();
135  RooAbsArg* param( 0 );
136  Bool_t allOK( kTRUE );
137  while ( ( param = argset_iter.next() ) ) {
138  TString name = param->GetName();
139  TBranch* branch = dataTree_->GetBranch( name );
140  if ( branch == 0 ) {
141  std::cerr << "ERROR in LauRooFitTask::verifyFitData : The data tree does not contain a branch for fit variable \""
142  << name << std::endl;
143  allOK = kFALSE;
144  }
145  }
146  if ( ! allOK ) {
147  return kFALSE;
148  }
149 
150  // Check whether the tree has the branch iExpt
151  TBranch* branch = dataTree_->GetBranch( "iExpt" );
152  if ( branch == 0 ) {
153  std::cout << "WARNING in LauRooFitTask::verifyFitData : Cannot find branch \"iExpt\" in the tree, will treat all data as being from a single experiment"
154  << std::endl;
155  } else {
156  // Define the valid values for iExpt
157  iExptSet_.clear();
158  const UInt_t firstExp = dataTree_->GetMinimum( "iExpt" );
159  const UInt_t lastExp = dataTree_->GetMaximum( "iExpt" );
160  for ( UInt_t iExp = firstExp; iExp <= lastExp; ++iExp ) {
161  iExptSet_.insert( iExp );
162  }
163  }
164 
165  return kTRUE;
166 }
167 
169 {
170  // Check that the NLL variable has been initialised
171  if ( ! nllVar_ ) {
172  std::cerr << "ERROR in LauRooFitTask::prepareInitialParArray : NLL var not initialised"
173  << std::endl;
174  return;
175  }
176 
177  // If we already prepared the entries in the fitPars_ vector then we only need to add the contents to the array
178  if ( ! fitPars_.empty() ) {
179  for ( std::vector<LauParameter*>::iterator iter = fitPars_.begin(); iter != fitPars_.end();
180  ++iter ) {
181  array.Add( *iter );
182  }
183  return;
184  }
185 
186  // Store the set of parameters and the total number of parameters
187  RooArgSet* varSet = nllVar_->getParameters( exptData_ );
188  UInt_t nFreePars( 0 );
189 
190  // Loop through the fit parameters
191  RooFIter argset_iter = varSet->fwdIterator();
192  RooAbsArg* param( 0 );
193  while ( ( param = argset_iter.next() ) ) {
194  // Only consider the free parameters
195  if ( ! param->isConstant() ) {
196  // Add the parameter
197  RooRealVar* rrvar = dynamic_cast<RooRealVar*>( param );
198  if ( rrvar != 0 ) {
199  // Count the number of free parameters
200  ++nFreePars;
201  // Do the conversion and add it to the array
202  LauParameter* lpar = this->convertToLauParameter( rrvar );
203  fitVars_.push_back( rrvar );
204  fitPars_.push_back( lpar );
205  array.Add( lpar );
206  } else {
207  RooFormulaVar* rfvar = dynamic_cast<RooFormulaVar*>( param );
208  if ( rfvar == 0 ) {
209  std::cerr << "ERROR in LauRooFitTask::prepareInitialParArray : The parameter is neither a RooRealVar nor a RooFormulaVar, don't know what to do"
210  << std::endl;
211  continue;
212  }
213  std::vector<std::pair<RooRealVar*, LauParameter*>> lpars =
214  this->convertToLauParameters( rfvar );
215  for ( std::vector<std::pair<RooRealVar*, LauParameter*>>::iterator iter =
216  lpars.begin();
217  iter != lpars.end();
218  ++iter ) {
219  RooRealVar* rrv = iter->first;
220  LauParameter* lpar = iter->second;
221  if ( ! rrv->isConstant() ) {
222  continue;
223  }
224 
225  // Count the number of free parameters
226  ++nFreePars;
227  // Add the parameter to the array
228  fitVars_.push_back( rrvar );
229  fitPars_.push_back( lpar );
230  array.Add( lpar );
231  }
232  }
233  }
234  }
235  delete varSet;
236 
237  this->startNewFit( nFreePars, nFreePars );
238 }
239 
240 LauParameter* LauRooFitTask::convertToLauParameter( const RooRealVar* rooParameter ) const
241 {
242  return new LauParameter( rooParameter->GetName(),
243  rooParameter->getVal(),
244  rooParameter->getMin(),
245  rooParameter->getMax(),
246  rooParameter->isConstant() );
247 }
248 
249 std::vector<std::pair<RooRealVar*, LauParameter*>> LauRooFitTask::convertToLauParameters(
250  const RooFormulaVar* rooFormula ) const
251 {
252  // Create the empty vector
253  std::vector<std::pair<RooRealVar*, LauParameter*>> lauParameters;
254 
255  Int_t parIndex( 0 );
256  RooAbsArg* rabsarg( 0 );
257  RooRealVar* rrvar( 0 );
258  RooFormulaVar* rfvar( 0 );
259  // Loop through all the parameters of the formula
260  while ( ( rabsarg = rooFormula->getParameter( parIndex ) ) ) {
261  // First try converting to a RooRealVar
262  rrvar = dynamic_cast<RooRealVar*>( rabsarg );
263  if ( rrvar ) {
264  // Do the conversion and add it to the array
265  LauParameter* lpar = this->convertToLauParameter( rrvar );
266  lauParameters.push_back( std::make_pair( rrvar, lpar ) );
267  continue;
268  }
269 
270  // If that didn't work, try converting to a RooFormulaVar
271  rfvar = dynamic_cast<RooFormulaVar*>( rabsarg );
272  if ( rfvar ) {
273  // Do the conversion and add these to the array
274  std::vector<std::pair<RooRealVar*, LauParameter*>> lpars = this->convertToLauParameters(
275  rfvar );
276  for ( std::vector<std::pair<RooRealVar*, LauParameter*>>::iterator iter = lpars.begin();
277  iter != lpars.end();
278  ++iter ) {
279  lauParameters.push_back( *iter );
280  }
281  continue;
282  }
283 
284  // If neither of those worked we don't know what to do, so print an error message and continue
285  std::cerr << "ERROR in LauRooFitTask::convertToLauParameters : One of the parameters is not a RooRealVar nor a RooFormulaVar, it is a: "
286  << rabsarg->ClassName() << std::endl;
287  std::cerr << " : Do not know how to process that - it will be skipped."
288  << std::endl;
289  }
290 
291  return lauParameters;
292 }
293 
295 {
296  Double_t nLL = ( nllVar_ != 0 ) ? nllVar_->getVal() : 0.0;
297  return nLL;
298 }
299 
300 void LauRooFitTask::setParsFromMinuit( Double_t* par, Int_t npar )
301 {
302  // This function sets the internal parameters based on the values
303  // that Minuit is using when trying to minimise the total likelihood function.
304 
305  // MINOS reports different numbers of free parameters depending on the
306  // situation, so disable this check
307  const UInt_t nFreePars = this->nFreeParams();
308  if ( ! this->withinAsymErrorCalc() ) {
309  if ( static_cast<UInt_t>( npar ) != nFreePars ) {
310  std::cerr << "ERROR in LauRooFitTask::setParsFromMinuit : Unexpected number of free parameters: "
311  << npar << ".\n";
312  std::cerr << " Expected: " << nFreePars
313  << ".\n"
314  << std::endl;
315  gSystem->Exit( EXIT_FAILURE );
316  }
317  }
318 
319  // Despite npar being the number of free parameters
320  // the par array actually contains all the parameters,
321  // free and floating...
322 
323  // Update all the floating ones with their new values
324  for ( UInt_t i( 0 ); i < nFreePars; ++i ) {
325  if ( ! fitPars_[i]->fixed() ) {
326  // Set both the RooRealVars and the LauParameters
327  fitPars_[i]->value( par[i] );
328  fitVars_[i]->setVal( par[i] );
329  }
330  }
331 }
332 
334 {
335  // check that we're being asked to read a valid index
336  const UInt_t exptIndex = this->iExpt();
337  if ( iExptSet_.empty() && exptIndex != 0 ) {
338  std::cerr << "ERROR in LauRooFitTask::readExperimentData : Invalid experiment number "
339  << exptIndex << ", data contains only one experiment" << std::endl;
340  return 0;
341  } else if ( iExptSet_.find( exptIndex ) == iExptSet_.end() ) {
342  std::cerr << "ERROR in LauRooFitTask::readExperimentData : Invalid experiment number "
343  << exptIndex << std::endl;
344  return 0;
345  }
346 
347  // cleanup the data from any previous experiment
348  delete exptData_;
349 
350  // retrieve the data and find out how many events have been read
351  if ( iExptSet_.empty() ) {
352  exptData_ = new RooDataSet( TString::Format( "expt%dData", exptIndex ),
353  "",
354  dataTree_,
355  dataVars_,
356  "",
357  ( weightVarName_ != "" ) ? weightVarName_.Data() : 0 );
358  } else {
359  const TString selectionString { TString::Format( "iExpt==%d", exptIndex ) };
360  TTree* exptTree { dataTree_->CopyTree( selectionString ) };
361  exptData_ = new RooDataSet( TString::Format( "expt%dData", exptIndex ),
362  "",
363  exptTree,
364  dataVars_,
365  "",
366  ( weightVarName_ != "" ) ? weightVarName_.Data() : 0 );
367  delete exptTree;
368  }
369 
370  const UInt_t nEvent = exptData_->numEntries();
371  this->eventsPerExpt( nEvent );
372  return nEvent;
373 }
374 
376 {
377  // cleanup the old NLL info
378  delete nllVar_;
379 
380  // construct the new NLL variable for this dataset
381  nllVar_ = new RooNLLVar( "nllVar", "", model_, *exptData_, extended_ );
382 }
383 
385  const TObjArray* parsFromCoordinator,
386  const TMatrixD* covMat,
387  TObjArray& parsToCoordinator )
388 {
389  // Copy the fit status information
390  this->storeFitStatus( fitStat, *covMat );
391 
392  // Now process the parameters
393  const UInt_t nFreePars = this->nFreeParams();
394  UInt_t nPars = parsFromCoordinator->GetEntries();
395  if ( nPars != nFreePars ) {
396  std::cerr << "ERROR in LauRooFitTask::finaliseExperiment : Unexpected number of parameters received from coordinator"
397  << std::endl;
398  std::cerr << " : Received " << nPars
399  << " when expecting " << nFreePars << std::endl;
400  gSystem->Exit( EXIT_FAILURE );
401  }
402 
403  for ( UInt_t iPar( 0 ); iPar < nPars; ++iPar ) {
404  LauParameter* parameter = dynamic_cast<LauParameter*>( ( *parsFromCoordinator )[iPar] );
405  if ( ! parameter ) {
406  std::cerr << "ERROR in LauRooFitTask::finaliseExperiment : Error reading parameter from coordinator"
407  << std::endl;
408  gSystem->Exit( EXIT_FAILURE );
409  }
410 
411  if ( parameter->name() != fitPars_[iPar]->name() ) {
412  std::cerr << "ERROR in LauRooFitTask::finaliseExperiment : Error reading parameter from coordinator"
413  << std::endl;
414  gSystem->Exit( EXIT_FAILURE );
415  }
416 
417  *( fitPars_[iPar] ) = *parameter;
418 
419  RooRealVar* rrv = fitVars_[iPar];
420  rrv->setVal( parameter->value() );
421  rrv->setError( parameter->error() );
422  rrv->setAsymError( parameter->negError(), parameter->posError() );
423  }
424 
425  // Update the pulls and add each finalised fit parameter to the list to
426  // send back to the coordinator
427  for ( std::vector<LauParameter*>::iterator iter = fitPars_.begin(); iter != fitPars_.end();
428  ++iter ) {
429  ( *iter )->updatePull();
430  parsToCoordinator.Add( *iter );
431  }
432 
433  // Write the results into the ntuple
434  LauFitNtuple* ntuple = this->fitNtuple();
435  ntuple->storeParsAndErrors( fitPars_, this->multiDimConstrainedPars(), {} );
436 
437  // find out the correlation matrix for the parameters
438  ntuple->storeCorrMatrix( this->iExpt(), this->fitStatus(), this->covarianceMatrix() );
439 
440  // Fill the data into ntuple
441  ntuple->updateFitNtuple();
442 }
LauParameter * convertToLauParameter(const RooRealVar *rooParameter) const
Convert a RooRealVar into a LauParameter.
const Bool_t extended_
Is the PDF extended?
virtual Bool_t verifyFitData(const TString &dataFileName, const TString &dataTreeName)
Open the input file and verify that all required variables are present.
RooNLLVar * nllVar_
The NLL variable.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
Struct to store fit status information.
Definition: LauAbsFitter.hh:54
Double_t value() const
The value of the parameter.
virtual ~LauRooFitTask()
Destructor.
File containing declaration of LauParameter class.
File containing declaration of LauFitNtuple class.
RooAbsPdf & model_
The model.
The base class for any task process for simultaneous/combined fits.
void storeFitStatus(const LauAbsFitter::FitStatus &status, const TMatrixD &covMatrix)
Store fit status information.
Definition: LauFitObject.cc:93
const std::set< TString > & multiDimConstrainedPars() const
Const access to the parameter names used in ND constraints.
virtual UInt_t readExperimentData()
Read in the data for the current experiment.
void storeCorrMatrix(const UInt_t iExpt, const LauAbsFitter::FitStatus &fitStatus, const TMatrixD &covMatrix)
Store the correlation matrix and other fit information.
Definition: LauFitNtuple.cc:76
virtual void initialise()
Initialise the fit model.
const LauAbsFitter::FitStatus & fitStatus() const
Access the fit status information.
virtual void finaliseExperiment(const LauAbsFitter::FitStatus &fitStat, const TObjArray *parsFromCoordinator, const TMatrixD *covMat, TObjArray &parsToCoordinator)
Perform all finalisation actions.
std::set< UInt_t > iExptSet_
The experiment category variable.
LauRooFitTask(RooAbsPdf &model, const Bool_t extended, const RooArgSet &vars, const TString &weightVarName="")
Constructor.
virtual void setParsFromMinuit(Double_t *par, Int_t npar)
This function sets the parameter values from Minuit.
Bool_t withinAsymErrorCalc() const
Query whether the fit is calculating the asymmetric errors.
Definition: LauFitObject.hh:94
std::vector< std::pair< RooRealVar *, LauParameter * > > convertToLauParameters(const RooFormulaVar *rooFormula) const
Convert a RooFormulaVar into LauParameters.
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:58
UInt_t iExpt() const
Obtain the number of the current experiment.
std::vector< LauParameter * > fitPars_
The fit parameters (as LauParameter's)
Double_t posError() const
The upper error on the parameter.
const LauFitNtuple * fitNtuple() const
Const access to the fit ntuple.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
TString weightVarName_
The name of the (optional) weight variable in the dataset.
virtual void prepareInitialParArray(TObjArray &array)
Package the initial fit parameters for transmission to the coordinator.
const TString & name() const
The parameter name.
void startNewFit(const UInt_t nPars, const UInt_t nFreePars)
Indicate the start of a new fit.
Definition: LauFitObject.cc:83
Double_t error() const
The error on the parameter.
void storeParsAndErrors(const std::vector< LauParameter * > &fitVars, const std::set< TString > &constrainedVars, const std::vector< LauParameter > &extraVars)
Store parameters and their errors.
std::vector< RooRealVar * > fitVars_
The fit parameters (as RooRealVar's)
TTree * dataTree_
The data tree.
void cleanData()
Cleanup the data.
const TMatrixD & covarianceMatrix() const
Access the fit covariance matrix.
RooAbsData * exptData_
The data for the current experiment.
File containing declaration of LauSimFitTask class.
Double_t negError() const
The lower error on the parameter.
TFile * dataFile_
The data file.
UInt_t eventsPerExpt() const
Obtain the total number of events in the current experiment.
void updateFitNtuple()
Update the fit ntuple.
File containing declaration of LauRooFitTask class.
RooArgSet dataVars_
The dataset variables.
virtual void cacheInputFitVars()
Cache the input data values to calculate the likelihood during the fit.
UInt_t nFreeParams() const
Access the total number of fit parameters.
virtual Double_t getTotNegLogLikelihood()
Calculates the total negative log-likelihood.