laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauSimpleFitModel.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2004 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 "LauSimpleFitModel.hh"
30 
31 #include "LauAbsBkgndDPModel.hh"
32 #include "LauAbsCoeffSet.hh"
33 #include "LauAbsPdf.hh"
34 #include "LauComplex.hh"
35 #include "LauConstants.hh"
36 #include "LauDaughters.hh"
37 #include "LauEffModel.hh"
38 #include "LauEmbeddedData.hh"
39 #include "LauFitNtuple.hh"
40 #include "LauGenNtuple.hh"
41 #include "LauIsobarDynamics.hh"
42 #include "LauKinematics.hh"
43 #include "LauPrint.hh"
44 #include "LauRandom.hh"
45 #include "LauScfMap.hh"
46 
47 #include "TCanvas.h"
48 #include "TFile.h"
49 #include "TGraph.h"
50 #include "TGraph2D.h"
51 #include "TH2.h"
52 #include "TMinuit.h"
53 #include "TRandom.h"
54 #include "TStyle.h"
55 #include "TSystem.h"
56 #include "TVirtualFitter.h"
57 
58 #include <fstream>
59 #include <iomanip>
60 #include <iostream>
61 #include <typeinfo>
62 
65  sigDPModel_( sigModel ),
66  kinematics_( sigModel ? sigModel->getKinematics() : 0 ),
67  usingBkgnd_( kFALSE ),
68  nSigComp_( 0 ),
69  nSigDPPar_( 0 ),
70  nExtraPdfPar_( 0 ),
71  nNormPar_( 0 ),
72  meanEff_( "meanEff", 0.0, 0.0, 1.0 ),
73  dpRate_( "dpRate", 0.0, 0.0, 100.0 ),
74  //signalEvents_("signalEvents",1.0,0.0,1.0),
75  signalEvents_( 0 ),
76  useSCF_( kFALSE ),
77  useSCFHist_( kFALSE ),
78  scfFrac_( "scfFrac", 0.0, 0.0, 1.0 ),
79  scfFracHist_( 0 ),
80  scfMap_( 0 ),
81  compareFitData_( kFALSE ),
82  signalTree_( 0 ),
83  reuseSignal_( kFALSE ),
84  useReweighting_( kFALSE ),
85  sigDPLike_( 0.0 ),
86  scfDPLike_( 0.0 ),
87  sigExtraLike_( 0.0 ),
88  scfExtraLike_( 0.0 ),
89  sigTotalLike_( 0.0 ),
90  scfTotalLike_( 0.0 )
91 {
92 }
93 
95 {
96  delete signalTree_;
97  for ( LauBkgndEmbDataList::iterator iter = bkgndTree_.begin(); iter != bkgndTree_.end(); ++iter ) {
98  delete ( *iter );
99  }
100  delete scfFracHist_;
101 }
102 
104 {
105  UInt_t nBkgnds = this->nBkgndClasses();
106  bkgndDPModels_.resize( nBkgnds );
107  bkgndPdfs_.resize( nBkgnds );
108  bkgndEvents_.resize( nBkgnds );
109  bkgndTree_.resize( nBkgnds );
110  reuseBkgnd_.resize( nBkgnds );
111  bkgndDPLike_.resize( nBkgnds );
112  bkgndExtraLike_.resize( nBkgnds );
113  bkgndTotalLike_.resize( nBkgnds );
114 }
115 
117 {
118  if ( nSigEvents == 0 ) {
119  std::cerr << "ERROR in LauSimpleFitModel::setNSigEvents : The signal yield LauParameter pointer is null."
120  << std::endl;
121  gSystem->Exit( EXIT_FAILURE );
122  }
123  if ( signalEvents_ != 0 ) {
124  std::cerr << "ERROR in LauSimpleFitModel::setNSigEvents : You are trying to overwrite the signal yield."
125  << std::endl;
126  return;
127  }
128 
129  signalEvents_ = nSigEvents;
130  TString name = signalEvents_->name();
131  if ( ! name.Contains( "signalEvents" ) &&
132  ! ( name.BeginsWith( "signal" ) && name.EndsWith( "Events" ) ) ) {
133  signalEvents_->name( "signalEvents" );
134  }
135  Double_t value = nSigEvents->value();
136  signalEvents_->range( -2.0 * ( TMath::Abs( value ) + 1.0 ), 2.0 * ( TMath::Abs( value ) + 1.0 ) );
137 }
138 
140 {
141  if ( nBkgndEvents == 0 ) {
142  std::cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : The background yield LauParameter pointer is null."
143  << std::endl;
144  gSystem->Exit( EXIT_FAILURE );
145  }
146 
147  if ( ! this->validBkgndClass( nBkgndEvents->name() ) ) {
148  std::cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : Invalid background class \""
149  << nBkgndEvents->name() << "\"." << std::endl;
150  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed."
151  << std::endl;
152  gSystem->Exit( EXIT_FAILURE );
153  }
154 
155  UInt_t bkgndID = this->bkgndClassID( nBkgndEvents->name() );
156 
157  if ( bkgndEvents_[bkgndID] != 0 ) {
158  std::cerr << "ERROR in LauSimpleFitModel::setNBkgndEvents : You are trying to overwrite the background yield."
159  << std::endl;
160  return;
161  }
162 
163  nBkgndEvents->name( nBkgndEvents->name() + "Events" );
164  if ( nBkgndEvents->isLValue() ) {
165  Double_t value = nBkgndEvents->value();
166  LauParameter* yield = dynamic_cast<LauParameter*>( nBkgndEvents );
167  yield->range( -2.0 * ( TMath::Abs( value ) + 1.0 ), 2.0 * ( TMath::Abs( value ) + 1.0 ) );
168  }
169  bkgndEvents_[bkgndID] = nBkgndEvents;
170 }
171 
173  const Bool_t upperHalf,
174  const Bool_t fluctuateBins,
175  LauScfMap* scfMap )
176 {
177  if ( useSCF_ == kTRUE ) {
178  std::cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : Have already setup SCF."
179  << std::endl;
180  return;
181  }
182 
183  if ( dpHisto == 0 ) {
184  std::cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : The histogram pointer is null."
185  << std::endl;
186  return;
187  }
188 
189  const LauDaughters* daughters = sigDPModel_->getDaughters();
190  scfFracHist_ = new LauEffModel( daughters, 0 );
192  ->setEffHisto( dpHisto, kTRUE, fluctuateBins, 0.0, 0.0, upperHalf, daughters->squareDP() );
193 
194  scfMap_ = scfMap;
195 
196  useSCF_ = kTRUE;
197  useSCFHist_ = kTRUE;
198 }
199 
200 void LauSimpleFitModel::splitSignalComponent( const Double_t scfFrac, const Bool_t fixed )
201 {
202  if ( useSCF_ == kTRUE ) {
203  std::cerr << "ERROR in LauSimpleFitModel::splitSignalComponent : Have already setup SCF."
204  << std::endl;
205  return;
206  }
207 
208  scfFrac_.range( 0.0, 1.0 );
209  scfFrac_.value( scfFrac );
210  scfFrac_.initValue( scfFrac );
211  scfFrac_.genValue( scfFrac );
212  scfFrac_.fixed( fixed );
213 
214  useSCF_ = kTRUE;
215  useSCFHist_ = kFALSE;
216 }
217 
218 void LauSimpleFitModel::setBkgndDPModel( const TString& bkgndClass, LauAbsBkgndDPModel* bkgndDPModel )
219 {
220  if ( bkgndDPModel == 0 ) {
221  std::cerr << "ERROR in LauSimpleFitModel::setBkgndDPModel : The model pointer is null."
222  << std::endl;
223  return;
224  }
225 
226  // check that this background name is valid
227  if ( ! this->validBkgndClass( bkgndClass ) ) {
228  std::cerr << "ERROR in LauSimpleFitModel::setBkgndDPModel : Invalid background class \""
229  << bkgndClass << "\"." << std::endl;
230  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed."
231  << std::endl;
232  return;
233  }
234 
235  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
236  bkgndDPModels_[bkgndID] = bkgndDPModel;
237 
238  usingBkgnd_ = kTRUE;
239 }
240 
242 {
243  if ( pdf == 0 ) {
244  std::cerr << "ERROR in LauSimpleFitModel::setSignalPdf : The PDF pointer is null."
245  << std::endl;
246  return;
247  }
248  signalPdfs_.push_back( pdf );
249 }
250 
252 {
253  if ( pdf == 0 ) {
254  std::cerr << "ERROR in LauSimpleFitModel::setSCFPdf : The PDF pointer is null." << std::endl;
255  return;
256  }
257  scfPdfs_.push_back( pdf );
258 }
259 
260 void LauSimpleFitModel::setBkgndPdf( const TString& bkgndClass, LauAbsPdf* pdf )
261 {
262  if ( pdf == 0 ) {
263  std::cerr << "ERROR in LauSimpleFitModel::setBkgndPdf : The PDF pointer is null."
264  << std::endl;
265  return;
266  }
267 
268  // check that this background name is valid
269  if ( ! this->validBkgndClass( bkgndClass ) ) {
270  std::cerr << "ERROR in LauSimpleFitModel::setBkgndPdf : Invalid background class \""
271  << bkgndClass << "\"." << std::endl;
272  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed."
273  << std::endl;
274  return;
275  }
276 
277  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
278  bkgndPdfs_[bkgndID].push_back( pdf );
279 
280  usingBkgnd_ = kTRUE;
281 }
282 
284 {
285  // Resize the coeffPars vector if not already done
286  if ( coeffPars_.empty() ) {
287  coeffPars_.resize( sigDPModel_->getnTotAmp() );
288  for ( std::vector<LauAbsCoeffSet*>::iterator iter = coeffPars_.begin();
289  iter != coeffPars_.end();
290  ++iter ) {
291  ( *iter ) = 0;
292  }
293  }
294 
295  // Is there a component called compName in the signal model?
296  TString compName( coeffSet->name() );
297  const Int_t index = sigDPModel_->resonanceIndex( compName );
298  if ( index < 0 ) {
299  std::cerr << "ERROR in LauSimpleFitModel::setAmpCoeffSet : Signal DP model doesn't contain component \""
300  << compName << "\"." << std::endl;
301  return;
302  }
303 
304  // Do we already have it in our list of names?
305  if ( coeffPars_[index] != 0 && coeffPars_[index]->name() == compName ) {
306  std::cerr << "ERROR in LauSimpleFitModel::setAmpCoeffSet : Have already set coefficients for \""
307  << compName << "\"." << std::endl;
308  return;
309  }
310 
311  coeffSet->index( index );
312  coeffPars_[index] = coeffSet;
313 
314  ++nSigComp_;
315 
316  std::cout << "INFO in LauSimpleFitModel::setAmpCoeffSet : Added coefficients for component A"
317  << index << ": \"" << compName << "\" to the fit model." << std::endl;
318  coeffSet->printParValues();
319 }
320 
322 {
323  // Initialisation
324  if ( ! this->useDP() && signalPdfs_.empty() ) {
325  std::cerr << "ERROR in LauSimpleFitModel::initialise : Signal model doesn't exist for any variable."
326  << std::endl;
327  gSystem->Exit( EXIT_FAILURE );
328  }
329 
330  if ( this->useDP() ) {
331  // Check that we have all the Dalitz-plot models
332  if ( sigDPModel_ == 0 ) {
333  std::cerr << "ERROR in LauSimpleFitModel::initialise : The pointer to the signal DP model is null.\n";
334  std::cerr << " : Removing the Dalitz Plot from the model."
335  << std::endl;
336  this->useDP( kFALSE );
337  }
338  if ( usingBkgnd_ ) {
339  for ( LauBkgndDPModelList::const_iterator dpmodel_iter = bkgndDPModels_.begin();
340  dpmodel_iter != bkgndDPModels_.end();
341  ++dpmodel_iter ) {
342  if ( ( *dpmodel_iter ) == 0 ) {
343  std::cerr << "ERROR in LauSimpleFitModel::initialise : The pointer to one of the background DP models is null.\n";
344  std::cerr << " : Removing the Dalitz Plot from the model."
345  << std::endl;
346  this->useDP( kFALSE );
347  break;
348  }
349  }
350  }
351 
352  // Need to check that the number of components we have and that the dynamics has matches up
353  UInt_t nAmp = sigDPModel_->getnTotAmp();
354  if ( nAmp != nSigComp_ ) {
355  std::cerr << "ERROR in LauSimpleFitModel::initialise : Number of signal DP components with magnitude and phase set not right."
356  << std::endl;
357  gSystem->Exit( EXIT_FAILURE );
358  }
359 
360  // From the initial parameter values calculate the coefficients
361  // so they can be passed to the signal model
362  this->updateCoeffs();
363 
364  // If all is well, go ahead and initialise them
365  this->initialiseDPModels();
366  }
367 
368  // Next check that, if a given component is being used we've got the
369  // right number of PDFs for all the variables involved
370  // TODO - should probably check variable names and so on as well
371 
372  UInt_t nsigpdfvars( 0 );
373  for ( LauPdfPList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end();
374  ++pdf_iter ) {
375  std::vector<TString> varNames = ( *pdf_iter )->varNames();
376  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
377  var_iter != varNames.end();
378  ++var_iter ) {
379  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
380  ++nsigpdfvars;
381  }
382  }
383  }
384  if ( useSCF_ ) {
385  UInt_t nscfpdfvars( 0 );
386  for ( LauPdfPList::const_iterator pdf_iter = scfPdfs_.begin(); pdf_iter != scfPdfs_.end();
387  ++pdf_iter ) {
388  std::vector<TString> varNames = ( *pdf_iter )->varNames();
389  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
390  var_iter != varNames.end();
391  ++var_iter ) {
392  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
393  ++nscfpdfvars;
394  }
395  }
396  }
397  if ( nscfpdfvars != nsigpdfvars ) {
398  std::cerr << "ERROR in LauSimpleFitModel::initialise : There are " << nsigpdfvars
399  << " TM signal PDF variables but " << nscfpdfvars
400  << " SCF signal PDF variables." << std::endl;
401  gSystem->Exit( EXIT_FAILURE );
402  }
403  }
404  if ( usingBkgnd_ ) {
405  for ( LauBkgndPdfsList::const_iterator bgclass_iter = bkgndPdfs_.begin();
406  bgclass_iter != bkgndPdfs_.end();
407  ++bgclass_iter ) {
408  UInt_t nbkgndpdfvars( 0 );
409  const LauPdfPList& pdfList = ( *bgclass_iter );
410  for ( LauPdfPList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end();
411  ++pdf_iter ) {
412  std::vector<TString> varNames = ( *pdf_iter )->varNames();
413  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
414  var_iter != varNames.end();
415  ++var_iter ) {
416  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
417  ++nbkgndpdfvars;
418  }
419  }
420  }
421  if ( nbkgndpdfvars != nsigpdfvars ) {
422  std::cerr << "ERROR in LauSimpleFitModel::initialise : There are " << nsigpdfvars
423  << " signal PDF variables but " << nbkgndpdfvars
424  << " bkgnd PDF variables." << std::endl;
425  gSystem->Exit( EXIT_FAILURE );
426  }
427  }
428  }
429 
430  // Clear the vectors of parameter information so we can start from scratch
431  this->clearFitParVectors();
432 
433  // Set the fit parameters for signal and background models
434  this->setSignalDPParameters();
435 
436  // Set the fit parameters for the various extra PDFs
437  this->setExtraPdfParameters();
438 
439  // Set the initial bg and signal events
440  this->setFitNEvents();
441 
442  // Check that we have the expected number of fit variables
443  const LauParameterPList& fitVars = this->fitPars();
444  if ( fitVars.size() != ( nSigDPPar_ + nExtraPdfPar_ + nNormPar_ ) ) {
445  std::cerr << "ERROR in LauSimpleFitModel::initialise : Number of fit parameters not of expected size."
446  << std::endl;
447  gSystem->Exit( EXIT_FAILURE );
448  }
449 
450  this->setExtraNtupleVars();
451 }
452 
454 {
455  std::cout << "INFO in LauSimpleFitModel::initialiseDPModels : Initialising DP models"
456  << std::endl;
457 
459 
460  if ( usingBkgnd_ == kTRUE ) {
461  for ( LauBkgndDPModelList::iterator iter = bkgndDPModels_.begin();
462  iter != bkgndDPModels_.end();
463  ++iter ) {
464  ( *iter )->initialise();
465  }
466  }
467 }
468 
470 {
471  //std::cout << "INFO in LauSimpleFitModel::recalculateNormalizationInDPModels : Recalc Norm in DP model" << std::endl;
474 }
475 
477 {
478  // Set the fit parameters for the signal model.
479 
480  nSigDPPar_ = 0;
481  if ( ! this->useDP() ) {
482  return;
483  }
484 
485  std::cout << "INFO in LauSimpleFitModel::setSignalDPParameters : Setting the initial fit parameters for the signal DP model."
486  << std::endl;
487 
488  // Place isobar coefficient parameters in vector of fit variables
489  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
490  LauParameterPList pars = coeffPars_[i]->getParameters();
491  nSigDPPar_ += this->addFitParameters( pars, kTRUE );
492  }
493 
494  // Obtain the resonance parameters and place them in the vector of fit variables and in a separate vector
496  nSigDPPar_ += this->addResonanceParameters( sigDPPars );
497 }
498 
500 {
501  // Include all the parameters of the various PDFs in the fit
502  // NB all of them are passed to the fit, even though some have been fixed through parameter.fixed(kTRUE)
503  // With the new "cloned parameter" scheme only "original" parameters are passed to the fit.
504  // Their clones are updated automatically when the originals are updated.
505 
506  nExtraPdfPar_ = 0;
507 
509 
510  if ( useSCF_ == kTRUE ) {
512  }
513 
514  if ( usingBkgnd_ == kTRUE ) {
515  for ( LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end();
516  ++iter ) {
517  nExtraPdfPar_ += this->addFitParameters( *iter );
518  }
519  }
520 }
521 
523 {
524  if ( signalEvents_ == 0 ) {
525  std::cerr << "ERROR in LauSimpleFitModel::setFitNEvents : Signal yield not defined."
526  << std::endl;
527  return;
528  }
529  nNormPar_ = 0;
530 
531  // initialise the total number of events to be the sum of all the hypotheses
532  Double_t nTotEvts = signalEvents_->value();
533  for ( LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end();
534  ++iter ) {
535  nTotEvts += ( *iter )->value();
536  if ( ( *iter ) == 0 ) {
537  std::cerr << "ERROR in LauSimpleFitModel::setFitNEvents : Background yield not defined."
538  << std::endl;
539  return;
540  }
541  }
542  this->eventsPerExpt( TMath::FloorNint( nTotEvts ) );
543 
544  // if doing an extended ML fit add the number of signal events into the fit parameters
545  if ( this->doEMLFit() ) {
546  std::cout << "INFO in LauSimpleFitModel::setFitNEvents : Initialising number of events for signal and background components..."
547  << std::endl;
548  // add the signal events to the list of fit parameters
550  } else {
551  std::cout << "INFO in LauSimpleFitModel::setFitNEvents : Initialising number of events for background components (and hence signal)..."
552  << std::endl;
553  }
554 
555  if ( useSCF_ && ! useSCFHist_ ) {
556  nNormPar_ += this->addFitParameters( &scfFrac_ );
557  }
558 
559  if ( usingBkgnd_ == kTRUE ) {
561  }
562 }
563 
565 {
566  // Set-up other parameters derived from the fit results, e.g. fit fractions.
567  if ( this->useDP() != kTRUE ) {
568  return;
569  }
570 
571  // First clear the vectors so we start from scratch
572  this->clearExtraVarVectors();
573 
574  LauParameterList& extraVars = this->extraPars();
575 
576  // Add a fit fraction for each signal component
578  if ( fitFrac_.size() != nSigComp_ ) {
579  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: "
580  << fitFrac_.size() << std::endl;
581  gSystem->Exit( EXIT_FAILURE );
582  }
583  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
584  if ( fitFrac_[i].size() != nSigComp_ ) {
585  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: "
586  << fitFrac_[i].size() << std::endl;
587  gSystem->Exit( EXIT_FAILURE );
588  }
589  }
590 
591  // Add the fit fraction that has not been corrected for the efficiency for each signal component
593  if ( fitFracEffUnCorr_.size() != nSigComp_ ) {
594  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: "
595  << fitFracEffUnCorr_.size() << std::endl;
596  gSystem->Exit( EXIT_FAILURE );
597  }
598  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
599  if ( fitFracEffUnCorr_[i].size() != nSigComp_ ) {
600  std::cerr << "ERROR in LauSimpleFitModel::setExtraNtupleVars : Initial Fit Fraction array of unexpected dimension: "
601  << fitFracEffUnCorr_[i].size() << std::endl;
602  gSystem->Exit( EXIT_FAILURE );
603  }
604  }
605 
606  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
607  for ( UInt_t j = i; j < nSigComp_; j++ ) {
608  extraVars.push_back( fitFrac_[i][j] );
609  extraVars.push_back( fitFracEffUnCorr_[i][j] );
610  }
611  }
612 
613  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
614  std::vector<LauParameter> extraParams = sigDPModel_->getExtraParameters();
615  std::vector<LauParameter>::iterator extraIter;
616  for ( extraIter = extraParams.begin(); extraIter != extraParams.end(); ++extraIter ) {
617  LauParameter extraParameter = ( *extraIter );
618  extraVars.push_back( extraParameter );
619  }
620 
621  // Now add in the DP efficiency value
622  Double_t initMeanEff = sigDPModel_->getMeanEff().initValue();
623  meanEff_.value( initMeanEff );
624  meanEff_.initValue( initMeanEff );
625  meanEff_.genValue( initMeanEff );
626  extraVars.push_back( meanEff_ );
627 
628  // Also add in the DP rate
629  Double_t initDPRate = sigDPModel_->getDPRate().initValue();
630  dpRate_.value( initDPRate );
631  dpRate_.initValue( initDPRate );
632  dpRate_.genValue( initDPRate );
633  extraVars.push_back( dpRate_ );
634 }
635 
636 void LauSimpleFitModel::finaliseFitResults( const TString& tablePrefixName )
637 {
638  // Retrieve parameters from the fit results for calculations and toy generation
639  // and eventually store these in output root ntuples/text files
640 
641  // Now take the fit parameters and update them as necessary
642  // e.g. to make mag > 0.0 and phase in the right range.
643  // This function will also calculate any other values, such as the
644  // fit fractions, using any errors provided by fitParErrors as appropriate.
645  // Also obtain the pull values: (measured - generated)/(average error)
646 
647  if ( this->useDP() == kTRUE ) {
648  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
649  // Check whether we have mag > 0.0, and phase in the right range
650  coeffPars_[i]->finaliseValues();
651  }
652  }
653 
654  // update the pulls on the events
655  if ( this->doEMLFit() ) {
657  }
658  if ( useSCF_ && ! useSCFHist_ ) {
660  }
661  if ( usingBkgnd_ == kTRUE ) {
662  for ( LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end();
663  ++iter ) {
664  std::vector<LauParameter*> parameters = ( *iter )->getPars();
665  for ( LauParameter* parameter : parameters ) {
666  parameter->updatePull();
667  }
668  }
669  }
670 
671  // Update the pulls on all the extra PDFs' parameters
673  if ( useSCF_ == kTRUE ) {
674  this->updateFitParameters( scfPdfs_ );
675  }
676  if ( usingBkgnd_ == kTRUE ) {
677  for ( LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end();
678  ++iter ) {
679  this->updateFitParameters( *iter );
680  }
681  }
682 
683  // Fill the fit results to the ntuple for current experiment
684 
685  // update the coefficients and then calculate the fit fractions
686  if ( this->useDP() == kTRUE ) {
687  this->updateCoeffs();
690 
692  if ( fitFrac.size() != nSigComp_ ) {
693  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: "
694  << fitFrac.size() << std::endl;
695  gSystem->Exit( EXIT_FAILURE );
696  }
697  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
698  if ( fitFrac[i].size() != nSigComp_ ) {
699  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: "
700  << fitFrac[i].size() << std::endl;
701  gSystem->Exit( EXIT_FAILURE );
702  }
703  }
704 
706  if ( fitFracEffUnCorr.size() != nSigComp_ ) {
707  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: "
708  << fitFracEffUnCorr.size() << std::endl;
709  gSystem->Exit( EXIT_FAILURE );
710  }
711  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
712  if ( fitFracEffUnCorr[i].size() != nSigComp_ ) {
713  std::cerr << "ERROR in LauSimpleFitModel::finaliseFitResults : Fit Fraction array of unexpected dimension: "
714  << fitFracEffUnCorr[i].size() << std::endl;
715  gSystem->Exit( EXIT_FAILURE );
716  }
717  }
718 
719  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
720  for ( UInt_t j( i ); j < nSigComp_; ++j ) {
721  fitFrac_[i][j].value( fitFrac[i][j].value() );
722  fitFracEffUnCorr_[i][j].value( fitFracEffUnCorr[i][j].value() );
723  }
724  }
725 
728 
729  this->clearExtraVarVectors();
730  LauParameterList& extraVars = this->extraPars();
731 
732  // Then store the final fit parameters, and any extra parameters for
733  // the signal model (e.g. fit fractions)
734  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
735  for ( UInt_t j( i ); j < nSigComp_; ++j ) {
736  extraVars.push_back( fitFrac_[i][j] );
737  extraVars.push_back( fitFracEffUnCorr_[i][j] );
738  }
739  }
740 
741  // Store any extra parameters/quantities from the DP model (e.g. K-matrix total fit fractions)
742  std::vector<LauParameter> extraParams = sigDPModel_->getExtraParameters();
743  std::vector<LauParameter>::iterator extraIter;
744  for ( extraIter = extraParams.begin(); extraIter != extraParams.end(); ++extraIter ) {
745  LauParameter extraParameter = ( *extraIter );
746  extraVars.push_back( extraParameter );
747  }
748 
749  // Now add in the DP efficiency value
750  extraVars.push_back( meanEff_ );
751 
752  // Also add in the DP rate
753  extraVars.push_back( dpRate_ );
754 
755  this->printFitFractions( std::cout );
756  }
757 
758  LauFitNtuple* ntuple = this->fitNtuple();
759  ntuple->storeParsAndErrors( this->fitPars(), this->multiDimConstrainedPars(), this->extraPars() );
760 
761  // find out the correlation matrix for the parameters
762  ntuple->storeCorrMatrix( this->iExpt(), this->fitStatus(), this->covarianceMatrix() );
763 
764  // Fill the data into ntuple
765  ntuple->updateFitNtuple();
766 
767  // Print out the partial fit fractions, phases and the
768  // averaged efficiency, reweighted by the dynamics (and anything else)
769  if ( this->writeLatexTable() ) {
770  TString sigOutFileName( tablePrefixName );
771  sigOutFileName += "_";
772  sigOutFileName += this->iExpt();
773  sigOutFileName += "Expt.tex";
774  this->writeOutTable( sigOutFileName );
775  }
776 }
777 
778 void LauSimpleFitModel::printFitFractions( std::ostream& output )
779 {
780  // Print out Fit Fractions, total DP rate and mean efficiency
781  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
782  output << "FitFraction for component " << i << " (" << coeffPars_[i]->name()
783  << ") = " << fitFrac_[i][i] << std::endl;
784  }
785  output << "Overall DP rate (integral of matrix element squared) = " << dpRate_ << std::endl;
786  output << "Average efficiency weighted by whole DP dynamics = " << meanEff_ << std::endl;
787 }
788 
789 void LauSimpleFitModel::writeOutTable( const TString& outputFile )
790 {
791  // Write out the results of the fit to a tex-readable table
792  // TODO - need to include the yields in this table
793  std::ofstream fout( outputFile );
794  LauPrint print;
795 
796  std::cout << "INFO in LauSimpleFitModel::writeOutTable : Writing out results of the fit to the tex file "
797  << outputFile << std::endl;
798 
799  if ( this->useDP() == kTRUE ) {
800  // print the fit coefficients in one table
801  coeffPars_.front()->printTableHeading( fout );
802  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
803  coeffPars_[i]->printTableRow( fout );
804  }
805  fout << "\\hline" << std::endl;
806  fout << "\\end{tabular}" << std::endl << std::endl;
807 
808  // print the fit fractions in another
809  fout << "\\begin{tabular}{|l|c|}" << std::endl;
810  fout << "\\hline" << std::endl;
811  fout << "Component & FitFraction \\\\" << std::endl;
812  fout << "\\hline" << std::endl;
813  Double_t fitFracSum( 0.0 );
814  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
815  TString resName = coeffPars_[i]->name();
816  resName = resName.ReplaceAll( "_", "\\_" );
817  fout << resName << " & $";
818  Double_t fitFrac = fitFrac_[i][i].value();
819  fitFracSum += fitFrac;
820  print.printFormat( fout, fitFrac );
821  fout << "$ \\\\" << std::endl;
822  }
823  fout << "\\hline" << std::endl;
824 
825  // Also print out sum of fit fractions
826  fout << "Fit Fraction Sum & $";
827  print.printFormat( fout, fitFracSum );
828  fout << "$ \\\\" << std::endl;
829  fout << "\\hline" << std::endl;
830 
831  fout << "DP rate & $";
832  print.printFormat( fout, dpRate_.value() );
833  fout << "$ \\\\" << std::endl;
834  fout << "\\hline" << std::endl;
835 
836  fout << "$< \\varepsilon >$ & $";
837  print.printFormat( fout, meanEff_.value() );
838  fout << "$ \\\\" << std::endl;
839  fout << "\\hline" << std::endl;
840  fout << "\\end{tabular}" << std::endl << std::endl;
841  }
842 
843  if ( ! signalPdfs_.empty() ) {
844  fout << "\\begin{tabular}{|l|c|}" << std::endl;
845  fout << "\\hline" << std::endl;
846  if ( useSCF_ == kTRUE ) {
847  fout << "\\Extra TM Signal PDFs' Parameters: & \\\\" << std::endl;
848  } else {
849  fout << "\\Extra Signal PDFs' Parameters: & \\\\" << std::endl;
850  }
851  this->printFitParameters( signalPdfs_, fout );
852  if ( useSCF_ == kTRUE && ! scfPdfs_.empty() ) {
853  fout << "\\hline" << std::endl;
854  fout << "\\Extra SCF Signal PDFs' Parameters: & \\\\" << std::endl;
855  this->printFitParameters( scfPdfs_, fout );
856  }
857  if ( usingBkgnd_ == kTRUE && ! bkgndPdfs_.empty() ) {
858  fout << "\\hline" << std::endl;
859  fout << "\\Extra Background PDFs' Parameters: & \\\\" << std::endl;
860  for ( LauBkgndPdfsList::const_iterator iter = bkgndPdfs_.begin();
861  iter != bkgndPdfs_.end();
862  ++iter ) {
863  this->printFitParameters( *iter, fout );
864  }
865  }
866  fout << "\\hline \n\\end{tabular}" << std::endl << std::endl;
867  }
868 }
869 
871 {
872  // Update the number of signal events to be total-sum(background events)
873  this->updateSigEvents();
874 
875  // Check whether we want to have randomised initial fit parameters for the signal model
876  if ( this->useRandomInitFitPars() == kTRUE ) {
877  std::cout << "INFO in LauSimpleFitModel::checkInitFitParams : Setting random parameters for the signal DP model"
878  << std::endl;
879  this->randomiseInitFitPars();
880  }
881 }
882 
884 {
885  // Only randomise those parameters that are not fixed!
886  std::cout << "INFO in LauSimpleFitModel::randomiseInitFitPars : Randomising the initial fit magnitudes and phases of the resonances..."
887  << std::endl;
888 
889  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
890  coeffPars_[i]->randomiseInitValues();
891  }
892 }
893 
894 std::pair<LauSimpleFitModel::LauGenInfo, Bool_t> LauSimpleFitModel::eventsToGenerate()
895 {
896  // Determine the number of events to generate for each hypothesis
897  // If we're smearing then smear each one individually
898 
899  LauGenInfo nEvtsGen;
900 
901  // Keep track of whether any yield or asymmetry parameters are blinded
902  Bool_t blind = kFALSE;
903 
904  // Signal
905  if ( signalEvents_->blind() ) {
906  blind = kTRUE;
907  }
908 
909  Double_t evtWeight( 1.0 );
910  Double_t nEvts = signalEvents_->genValue();
911  if ( nEvts < 0 ) {
912  evtWeight = -1.0;
913  nEvts = TMath::Abs( nEvts );
914  }
915 
916  Int_t nEvtsToGen { static_cast<Int_t>( nEvts ) };
917  if ( this->doPoissonSmearing() ) {
918  nEvtsToGen = LauRandom::randomFun()->Poisson( nEvts );
919  }
920 
921  nEvtsGen["signal"] = std::make_pair( nEvtsToGen, evtWeight );
922 
923  // Backgrounds
924  const UInt_t nBkgnds = this->nBkgndClasses();
925  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
926  const LauAbsRValue* evtsPar = bkgndEvents_[bkgndID];
927  if ( evtsPar->blind() ) {
928  blind = kTRUE;
929  }
930 
931  evtWeight = 1.0;
932  nEvts = evtsPar->genValue();
933  if ( nEvts < 0 ) {
934  evtWeight = -1.0;
935  nEvts = TMath::Abs( nEvts );
936  }
937 
938  nEvtsToGen = static_cast<Int_t>( nEvts );
939  if ( this->doPoissonSmearing() ) {
940  nEvtsToGen = LauRandom::randomFun()->Poisson( nEvts );
941  }
942 
943  const TString& bkgndClass = this->bkgndClassName( bkgndID );
944  nEvtsGen[bkgndClass] = std::make_pair( nEvtsToGen, evtWeight );
945  }
946 
947  return std::make_pair( nEvtsGen, blind );
948 }
949 
951 {
952  // Routine to generate toy Monte Carlo events according to the various models we have defined.
953 
954  // Determine the number of events to generate for each hypothesis
955  std::pair<LauGenInfo, Bool_t> info = this->eventsToGenerate();
956  LauGenInfo nEvts = info.first;
957  const Bool_t blind = info.second;
958 
959  Bool_t genOK( kTRUE );
960  Int_t evtNum( 0 );
961 
962  const UInt_t nBkgnds = this->nBkgndClasses();
963  std::vector<TString> bkgndClassNames( nBkgnds );
964  std::vector<TString> bkgndClassNamesGen( nBkgnds );
965  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
966  TString name( this->bkgndClassName( iBkgnd ) );
967  bkgndClassNames[iBkgnd] = name;
968  bkgndClassNamesGen[iBkgnd] = "gen" + name;
969  }
970 
971  const Bool_t storeSCFTruthInfo = ( useSCF_ || ( this->enableEmbedding() && signalTree_ != 0 &&
972  signalTree_->haveBranch( "mcMatch" ) ) );
973 
974  // Loop over the hypotheses and generate the requested number of events for each one
975  for ( LauGenInfo::const_iterator iter = nEvts.begin(); iter != nEvts.end(); ++iter ) {
976 
977  const TString& type( iter->first );
978  Int_t nEvtsGen( iter->second.first );
979  Double_t evtWeight( iter->second.second );
980 
981  for ( Int_t iEvt( 0 ); iEvt < nEvtsGen; ++iEvt ) {
982 
983  this->setGenNtupleDoubleBranchValue( "evtWeight", evtWeight );
984  // Add efficiency information
985  this->setGenNtupleDoubleBranchValue( "efficiency", 1 );
986 
987  if ( type == "signal" ) {
988  this->setGenNtupleIntegerBranchValue( "genSig", 1 );
989  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
990  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], 0 );
991  }
992  genOK = this->generateSignalEvent();
993  this->setGenNtupleDoubleBranchValue( "efficiency", sigDPModel_->getEvtEff() );
994  } else {
995  this->setGenNtupleIntegerBranchValue( "genSig", 0 );
996  if ( storeSCFTruthInfo ) {
997  this->setGenNtupleIntegerBranchValue( "genTMSig", 0 );
998  this->setGenNtupleIntegerBranchValue( "genSCFSig", 0 );
999  }
1000  UInt_t bkgndID( 0 );
1001  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1002  Int_t gen( 0 );
1003  if ( bkgndClassNames[iBkgnd] == type ) {
1004  gen = 1;
1005  bkgndID = iBkgnd;
1006  }
1007  this->setGenNtupleIntegerBranchValue( bkgndClassNamesGen[iBkgnd], gen );
1008  }
1009  genOK = this->generateBkgndEvent( bkgndID );
1010  }
1011 
1012  if ( ! genOK ) {
1013  // If there was a problem with the generation then break out and return.
1014  // The problem model will have adjusted itself so that all should be OK next time.
1015  break;
1016  }
1017 
1018  if ( this->useDP() == kTRUE ) {
1019  this->setDPBranchValues();
1020  }
1021 
1022  // Store the event number (within this experiment)
1023  this->setGenNtupleIntegerBranchValue( "iEvtWithinExpt", evtNum );
1024  // and then increment it
1025  ++evtNum;
1026 
1027  this->fillGenNtupleBranches();
1028  if ( ! blind && ( iEvt % 500 == 0 ) ) {
1029  std::cout << "INFO in LauSimpleFitModel::genExpt : Generated event number " << iEvt
1030  << " out of " << nEvtsGen << " " << type << " events." << std::endl;
1031  }
1032  }
1033 
1034  if ( ! genOK ) {
1035  break;
1036  }
1037  }
1038 
1039  if ( this->useDP() && genOK ) {
1040 
1041  sigDPModel_->checkToyMC( kTRUE, kTRUE );
1042 
1043  // Get the fit fractions if they're to be written into the latex table
1044  if ( this->writeLatexTable() ) {
1046  if ( fitFrac.size() != nSigComp_ ) {
1047  std::cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: "
1048  << fitFrac.size() << std::endl;
1049  gSystem->Exit( EXIT_FAILURE );
1050  }
1051  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
1052  if ( fitFrac[i].size() != nSigComp_ ) {
1053  std::cerr << "ERROR in LauSimpleFitModel::genExpt : Fit Fraction array of unexpected dimension: "
1054  << fitFrac[i].size() << std::endl;
1055  gSystem->Exit( EXIT_FAILURE );
1056  }
1057  }
1058 
1059  for ( UInt_t i( 0 ); i < nSigComp_; ++i ) {
1060  for ( UInt_t j = i; j < nSigComp_; j++ ) {
1061  fitFrac_[i][j].value( fitFrac[i][j].value() );
1062  }
1063  }
1066  }
1067  }
1068 
1069  // If we're reusing embedded events or if the generation is being
1070  // reset then clear the lists of used events
1071  if ( reuseSignal_ || ! genOK ) {
1072  if ( signalTree_ ) {
1074  }
1075  }
1076  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
1077  LauEmbeddedData* data = bkgndTree_[bkgndID];
1078  if ( reuseBkgnd_[bkgndID] || ! genOK ) {
1079  if ( data ) {
1080  data->clearUsedList();
1081  }
1082  }
1083  }
1084 
1085  return genOK;
1086 }
1087 
1089 {
1090  // Generate signal event
1091  Bool_t genOK( kTRUE );
1092  Bool_t genSCF( kFALSE );
1093 
1094  if ( this->useDP() ) {
1095  if ( this->enableEmbedding() && signalTree_ ) {
1096  if ( useReweighting_ ) {
1097  // Select a (random) event from the generated data. Then store the
1098  // reconstructed DP co-ords, together with other pdf information,
1099  // as the embedded data.
1101  } else {
1102 
1104  }
1105  if ( signalTree_->haveBranch( "mcMatch" ) ) {
1106  Int_t match = static_cast<Int_t>( signalTree_->getValue( "mcMatch" ) );
1107  if ( match ) {
1108  this->setGenNtupleIntegerBranchValue( "genTMSig", 1 );
1109  this->setGenNtupleIntegerBranchValue( "genSCFSig", 0 );
1110  genSCF = kFALSE;
1111  } else {
1112  this->setGenNtupleIntegerBranchValue( "genTMSig", 0 );
1113  this->setGenNtupleIntegerBranchValue( "genSCFSig", 1 );
1114  genSCF = kTRUE;
1115  }
1116  }
1117  } else {
1118  genOK = sigDPModel_->generate();
1119  if ( genOK && useSCF_ ) {
1120  Double_t frac( 0.0 );
1121  if ( useSCFHist_ ) {
1123  } else {
1124  frac = scfFrac_.genValue();
1125  }
1126  if ( frac < LauRandom::randomFun()->Rndm() ) {
1127  this->setGenNtupleIntegerBranchValue( "genTMSig", 1 );
1128  this->setGenNtupleIntegerBranchValue( "genSCFSig", 0 );
1129  genSCF = kFALSE;
1130  } else {
1131  this->setGenNtupleIntegerBranchValue( "genTMSig", 0 );
1132  this->setGenNtupleIntegerBranchValue( "genSCFSig", 1 );
1133  genSCF = kTRUE;
1134 
1135  // Optionally smear the DP position
1136  // of the SCF event
1137  if ( scfMap_ != 0 ) {
1138  Double_t xCoord( 0.0 ), yCoord( 0.0 );
1139  if ( kinematics_->squareDP() ) {
1140  xCoord = kinematics_->getmPrime();
1141  yCoord = kinematics_->getThetaPrime();
1142  } else {
1143  xCoord = kinematics_->getm13Sq();
1144  yCoord = kinematics_->getm23Sq();
1145  }
1146 
1147  // Find the bin number where this event is generated
1148  Int_t binNo = scfMap_->binNumber( xCoord, yCoord );
1149 
1150  // Retrieve the migration histogram
1151  TH2* histo = scfMap_->trueHist( binNo );
1152 
1153  const LauAbsEffModel* effModel = sigDPModel_->getEffModel();
1154  do {
1155  // Get a random point from the histogram
1156  histo->GetRandom2( xCoord, yCoord );
1157 
1158  // Update the kinematics
1159  if ( kinematics_->squareDP() ) {
1160  kinematics_->updateSqDPKinematics( xCoord, yCoord );
1161  } else {
1162  kinematics_->updateKinematics( xCoord, yCoord );
1163  }
1164  } while ( ! effModel->passVeto( kinematics_ ) );
1165  }
1166  }
1167  }
1168  }
1169  } else {
1170  if ( this->enableEmbedding() && signalTree_ ) {
1172  if ( signalTree_->haveBranch( "mcMatch" ) ) {
1173  Int_t match = static_cast<Int_t>( signalTree_->getValue( "mcMatch" ) );
1174  if ( match ) {
1175  this->setGenNtupleIntegerBranchValue( "genTMSig", 1 );
1176  this->setGenNtupleIntegerBranchValue( "genSCFSig", 0 );
1177  genSCF = kFALSE;
1178  } else {
1179  this->setGenNtupleIntegerBranchValue( "genTMSig", 0 );
1180  this->setGenNtupleIntegerBranchValue( "genSCFSig", 1 );
1181  genSCF = kTRUE;
1182  }
1183  }
1184  } else if ( useSCF_ ) {
1185  Double_t frac = scfFrac_.genValue();
1186  if ( frac < LauRandom::randomFun()->Rndm() ) {
1187  this->setGenNtupleIntegerBranchValue( "genTMSig", 1 );
1188  this->setGenNtupleIntegerBranchValue( "genSCFSig", 0 );
1189  genSCF = kFALSE;
1190  } else {
1191  this->setGenNtupleIntegerBranchValue( "genTMSig", 0 );
1192  this->setGenNtupleIntegerBranchValue( "genSCFSig", 1 );
1193  genSCF = kTRUE;
1194  }
1195  }
1196  }
1197  if ( genOK ) {
1198  if ( useSCF_ ) {
1199  if ( genSCF ) {
1201  } else {
1203  }
1204  } else {
1206  }
1207  }
1208  // Check for problems with the embedding
1209  if ( this->enableEmbedding() && signalTree_ &&
1210  ( signalTree_->nEvents() == signalTree_->nUsedEvents() ) ) {
1211  std::cerr << "WARNING in LauSimpleFitModel::generateSignalEvent : Source of embedded signal events used up, clearing the list of used events."
1212  << std::endl;
1214  }
1215  return genOK;
1216 }
1217 
1219 {
1220  // Generate background event
1221  Bool_t genOK( kTRUE );
1222 
1223  LauAbsBkgndDPModel* model = bkgndDPModels_[bkgndID];
1224  LauEmbeddedData* embeddedData( 0 );
1225  if ( this->enableEmbedding() ) {
1226  embeddedData = bkgndTree_[bkgndID];
1227  }
1228  LauPdfPList* extraPdfs = &bkgndPdfs_[bkgndID];
1229  if ( this->useDP() ) {
1230  if ( embeddedData ) {
1231  embeddedData->getEmbeddedEvent( kinematics_ );
1232  } else {
1233  if ( model == 0 ) {
1234  std::cerr << "ERROR in LauSimpleFitModel::generateBkgndEvent : Can't find the DP model for background class \""
1235  << this->bkgndClassName( bkgndID ) << "\"." << std::endl;
1236  gSystem->Exit( EXIT_FAILURE );
1237  }
1238  genOK = model->generate();
1239  }
1240  } else {
1241  if ( embeddedData ) {
1242  embeddedData->getEmbeddedEvent( 0 );
1243  }
1244  }
1245  if ( genOK ) {
1246  this->generateExtraPdfValues( extraPdfs, embeddedData );
1247  }
1248  // Check for problems with the embedding
1249  if ( embeddedData && ( embeddedData->nEvents() == embeddedData->nUsedEvents() ) ) {
1250  std::cerr << "WARNING in LauSimpleFitModel::generateBkgndEvent : Source of embedded "
1251  << this->bkgndClassName( bkgndID )
1252  << " events used up, clearing the list of used events." << std::endl;
1253  embeddedData->clearUsedList();
1254  }
1255  return genOK;
1256 }
1257 
1259 {
1260  // Setup the required ntuple branches
1261  this->addGenNtupleDoubleBranch( "evtWeight" );
1262  this->addGenNtupleIntegerBranch( "genSig" );
1263  this->addGenNtupleDoubleBranch( "efficiency" );
1264  if ( useSCF_ ||
1265  ( this->enableEmbedding() && signalTree_ != 0 && signalTree_->haveBranch( "mcMatch" ) ) ) {
1266  this->addGenNtupleIntegerBranch( "genTMSig" );
1267  this->addGenNtupleIntegerBranch( "genSCFSig" );
1268  }
1269  const UInt_t nBkgnds = this->nBkgndClasses();
1270  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1271  TString name( this->bkgndClassName( iBkgnd ) );
1272  name.Prepend( "gen" );
1273  this->addGenNtupleIntegerBranch( name );
1274  }
1275  if ( this->useDP() == kTRUE ) {
1276  this->addGenNtupleDoubleBranch( "m12" );
1277  this->addGenNtupleDoubleBranch( "m23" );
1278  this->addGenNtupleDoubleBranch( "m13" );
1279  this->addGenNtupleDoubleBranch( "m12Sq" );
1280  this->addGenNtupleDoubleBranch( "m23Sq" );
1281  this->addGenNtupleDoubleBranch( "m13Sq" );
1282  this->addGenNtupleDoubleBranch( "cosHel12" );
1283  this->addGenNtupleDoubleBranch( "cosHel23" );
1284  this->addGenNtupleDoubleBranch( "cosHel13" );
1285  if ( kinematics_->squareDP() ) {
1286  this->addGenNtupleDoubleBranch( "mPrime" );
1287  this->addGenNtupleDoubleBranch( "thPrime" );
1288  }
1289  }
1290  for ( LauPdfPList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end();
1291  ++pdf_iter ) {
1292  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1293  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1294  var_iter != varNames.end();
1295  ++var_iter ) {
1296  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
1297  this->addGenNtupleDoubleBranch( ( *var_iter ) );
1298  }
1299  }
1300  }
1301 }
1302 
1304 {
1305  // Store all the DP information
1309  this->setGenNtupleDoubleBranchValue( "m12Sq", kinematics_->getm12Sq() );
1310  this->setGenNtupleDoubleBranchValue( "m23Sq", kinematics_->getm23Sq() );
1311  this->setGenNtupleDoubleBranchValue( "m13Sq", kinematics_->getm13Sq() );
1312  this->setGenNtupleDoubleBranchValue( "cosHel12", kinematics_->getc12() );
1313  this->setGenNtupleDoubleBranchValue( "cosHel23", kinematics_->getc23() );
1314  this->setGenNtupleDoubleBranchValue( "cosHel13", kinematics_->getc13() );
1315  if ( kinematics_->squareDP() ) {
1316  this->setGenNtupleDoubleBranchValue( "mPrime", kinematics_->getmPrime() );
1318  }
1319 }
1320 
1322 {
1323  // Generate from the extra PDFs
1324  if ( extraPdfs ) {
1325  for ( LauPdfPList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end();
1326  ++pdf_iter ) {
1327  LauFitData genValues;
1328  if ( embeddedData ) {
1329  genValues = embeddedData->getValues( ( *pdf_iter )->varNames() );
1330  } else {
1331  genValues = ( *pdf_iter )->generate( kinematics_ );
1332  }
1333  for ( LauFitData::const_iterator var_iter = genValues.begin();
1334  var_iter != genValues.end();
1335  ++var_iter ) {
1336  TString varName = var_iter->first;
1337  if ( varName != "m13Sq" && varName != "m23Sq" ) {
1338  Double_t value = var_iter->second;
1339  this->setGenNtupleDoubleBranchValue( varName, value );
1340  }
1341  }
1342  }
1343  }
1344 }
1345 
1347 {
1348  // Update the total normalisation for the signal likelihood
1349  if ( this->useDP() == kTRUE ) {
1350  this->updateCoeffs();
1352  }
1353 
1354  // Update the signal events from the background events if not doing an extended fit
1355  if ( ! this->doEMLFit() && ! signalEvents_->fixed() ) {
1356  this->updateSigEvents();
1357  }
1358 }
1359 
1361 {
1362  // The background parameters will have been set from Minuit.
1363  // We need to update the signal events using these.
1364  Double_t nTotEvts = this->eventsPerExpt();
1365 
1366  signalEvents_->range( -2.0 * nTotEvts, 2.0 * nTotEvts );
1367  for ( LauBkgndYieldList::iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end();
1368  ++iter ) {
1369  LauAbsRValue* nBkgndEvents = ( *iter );
1370  if ( nBkgndEvents->isLValue() ) {
1371  LauParameter* yield = dynamic_cast<LauParameter*>( nBkgndEvents );
1372  yield->range( -2.0 * nTotEvts, 2.0 * nTotEvts );
1373  }
1374  }
1375 
1376  if ( signalEvents_->fixed() ) {
1377  return;
1378  }
1379 
1380  // Subtract background events (if any) from signal.
1381  Double_t signalEvents = nTotEvts;
1382  if ( usingBkgnd_ == kTRUE ) {
1383  for ( LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin();
1384  iter != bkgndEvents_.end();
1385  ++iter ) {
1386  signalEvents -= ( *iter )->value();
1387  }
1388  }
1389  signalEvents_->value( signalEvents );
1390 }
1391 
1393 {
1394  // Fill the internal data trees of the signal and backgrond models.
1395 
1396  LauFitDataTree* inputFitData = this->fitData();
1397 
1398  // First the Dalitz plot variables (m_ij^2)
1399  if ( this->useDP() == kTRUE ) {
1400 
1401  // need to append SCF smearing bins before caching DP amplitudes
1402  if ( scfMap_ != 0 ) {
1403  this->appendBinCentres( inputFitData );
1404  }
1405  sigDPModel_->fillDataTree( *inputFitData );
1406 
1407  if ( usingBkgnd_ == kTRUE ) {
1408  for ( LauBkgndDPModelList::iterator iter = bkgndDPModels_.begin();
1409  iter != bkgndDPModels_.end();
1410  ++iter ) {
1411  ( *iter )->fillDataTree( *inputFitData );
1412  }
1413  }
1414  }
1415 
1416  // ...and then the extra PDFs
1417  this->cacheInfo( signalPdfs_, *inputFitData );
1418  this->cacheInfo( scfPdfs_, *inputFitData );
1419  for ( LauBkgndPdfsList::iterator iter = bkgndPdfs_.begin(); iter != bkgndPdfs_.end(); ++iter ) {
1420  this->cacheInfo( ( *iter ), *inputFitData );
1421  }
1422 
1423  // and finally the SCF fractions and jacobians
1424  if ( useSCF_ && useSCFHist_ ) {
1425  if ( ! inputFitData->haveBranch( "m13Sq" ) || ! inputFitData->haveBranch( "m23Sq" ) ) {
1426  std::cerr << "ERROR in LauSimpleFitModel::cacheInputFitVars : Input data does not contain DP branches and so can't cache the SCF fraction."
1427  << std::endl;
1428  gSystem->Exit( EXIT_FAILURE );
1429  }
1430 
1431  UInt_t nEvents = inputFitData->nEvents();
1432  recoSCFFracs_.clear();
1433  recoSCFFracs_.reserve( nEvents );
1434  if ( kinematics_->squareDP() ) {
1435  recoJacobians_.clear();
1436  recoJacobians_.reserve( nEvents );
1437  }
1438  for ( UInt_t iEvt = 0; iEvt < nEvents; iEvt++ ) {
1439  const LauFitData& dataValues = inputFitData->getData( iEvt );
1440  LauFitData::const_iterator m13_iter = dataValues.find( "m13Sq" );
1441  LauFitData::const_iterator m23_iter = dataValues.find( "m23Sq" );
1442  kinematics_->updateKinematics( m13_iter->second, m23_iter->second );
1443  Double_t scfFrac = scfFracHist_->calcEfficiency( kinematics_ );
1444  recoSCFFracs_.push_back( scfFrac );
1445  if ( kinematics_->squareDP() ) {
1447  }
1448  }
1449  }
1450 }
1451 
1453 {
1454  // We'll be caching the DP amplitudes and efficiencies of the centres of the true bins.
1455  // To do so, we attach some fake points at the end of inputData, the number of the entry
1456  // minus the total number of events corresponding to the number of the histogram for that
1457  // given true bin in the LauScfMap object. (What this means is that when Laura is provided with
1458  // the LauScfMap object by the user, it's the latter who has to make sure that it contains the
1459  // right number of histograms and in exactly the right order!)
1460 
1461  // Get the x and y co-ordinates of the bin centres
1462  std::vector<Double_t> binCentresXCoords;
1463  std::vector<Double_t> binCentresYCoords;
1464  scfMap_->listBinCentres( binCentresXCoords, binCentresYCoords );
1465 
1466  // The SCF histograms could be in square Dalitz plot histograms.
1467  // The dynamics takes normal Dalitz plot coords, so we might have to convert them back.
1468  Bool_t sqDP = kinematics_->squareDP();
1469  UInt_t nBins = binCentresXCoords.size();
1470  fakeSCFFracs_.clear();
1471  fakeSCFFracs_.reserve( nBins );
1472  if ( sqDP ) {
1473  fakeJacobians_.clear();
1474  fakeJacobians_.reserve( nBins );
1475  }
1476 
1477  for ( UInt_t iBin = 0; iBin < nBins; ++iBin ) {
1478 
1479  if ( sqDP ) {
1480  kinematics_->updateSqDPKinematics( binCentresXCoords[iBin], binCentresYCoords[iBin] );
1481  binCentresXCoords[iBin] = kinematics_->getm13Sq();
1482  binCentresYCoords[iBin] = kinematics_->getm23Sq();
1484  } else {
1485  kinematics_->updateKinematics( binCentresXCoords[iBin], binCentresYCoords[iBin] );
1486  }
1487 
1489  }
1490 
1491  // Set up inputFitVars_ object to hold the fake events
1492  inputData->appendFakePoints( binCentresXCoords, binCentresYCoords );
1493 }
1494 
1496 {
1497  // Get the DP likelihood for signal and backgrounds
1498  this->getEvtDPLikelihood( iEvt );
1499 
1500  // Get the combined extra PDFs likelihood for signal and backgrounds
1501  this->getEvtExtraLikelihoods( iEvt );
1502 
1503  // If appropriate, combine the TM and SCF likelihoods
1504  Double_t sigLike = sigDPLike_ * sigExtraLike_;
1505  if ( useSCF_ ) {
1506  Double_t scfFrac( 0.0 );
1507  if ( useSCFHist_ ) {
1508  scfFrac = recoSCFFracs_[iEvt];
1509  } else {
1510  scfFrac = scfFrac_.unblindValue();
1511  }
1512  sigLike *= ( 1.0 - scfFrac );
1513  if ( ( scfMap_ != 0 ) && ( this->useDP() == kTRUE ) ) {
1514  // if we're smearing the SCF DP PDF then the SCF frac
1515  // is already included in the SCF DP likelihood
1516  sigLike += ( scfDPLike_ * scfExtraLike_ );
1517  } else {
1518  sigLike += ( scfFrac * scfDPLike_ * scfExtraLike_ );
1519  }
1520  }
1521 
1522  // Construct the total event likelihood
1523  Double_t likelihood = signalEvents_->unblindValue() * sigLike;
1524  const UInt_t nBkgnds = this->nBkgndClasses();
1525  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
1526  likelihood += ( bkgndEvents_[bkgndID]->unblindValue() * bkgndDPLike_[bkgndID] *
1527  bkgndExtraLike_[bkgndID] );
1528  }
1529 
1530  return likelihood;
1531 }
1532 
1534 {
1535  Double_t eventSum( 0.0 );
1536  eventSum += signalEvents_->unblindValue();
1537  for ( LauBkgndYieldList::const_iterator iter = bkgndEvents_.begin(); iter != bkgndEvents_.end();
1538  ++iter ) {
1539  eventSum += ( *iter )->unblindValue();
1540  }
1541  return eventSum;
1542 }
1543 
1545 {
1546  // Function to return the signal and background likelihoods for the
1547  // Dalitz plot for the given event evtNo.
1548 
1549  if ( this->useDP() == kTRUE ) {
1552 
1553  if ( useSCF_ == kTRUE ) {
1554  if ( scfMap_ == 0 ) {
1555  // we're not smearing the SCF DP position
1556  // so the likelihood is the same as the TM
1558  } else {
1559  // calculate the smeared SCF DP likelihood
1560  scfDPLike_ = this->getEvtSCFDPLikelihood( iEvt );
1561  }
1562  }
1563 
1564  const UInt_t nBkgnds = this->nBkgndClasses();
1565  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
1566  if ( usingBkgnd_ == kTRUE ) {
1567  bkgndDPLike_[bkgndID] = bkgndDPModels_[bkgndID]->getLikelihood( iEvt );
1568  } else {
1569  bkgndDPLike_[bkgndID] = 0.0;
1570  }
1571  }
1572  } else {
1573  // There's always going to be a term in the likelihood for the
1574  // signal, so we'd better not zero it.
1575  sigDPLike_ = 1.0;
1576  scfDPLike_ = 1.0;
1577 
1578  const UInt_t nBkgnds = this->nBkgndClasses();
1579  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
1580  if ( usingBkgnd_ == kTRUE ) {
1581  bkgndDPLike_[bkgndID] = 1.0;
1582  } else {
1583  bkgndDPLike_[bkgndID] = 0.0;
1584  }
1585  }
1586  }
1587 }
1588 
1590 {
1591  Double_t scfDPLike( 0.0 );
1592 
1593  Double_t recoJacobian( 1.0 );
1594  Double_t xCoord( 0.0 );
1595  Double_t yCoord( 0.0 );
1596  Bool_t squareDP = kinematics_->squareDP();
1597  if ( squareDP ) {
1598  xCoord = sigDPModel_->getEvtmPrime();
1599  yCoord = sigDPModel_->getEvtthPrime();
1600  recoJacobian = recoJacobians_[iEvt];
1601  } else {
1602  xCoord = sigDPModel_->getEvtm13Sq();
1603  yCoord = sigDPModel_->getEvtm23Sq();
1604  }
1605 
1606  // Find the bin that our reco event falls in
1607  Int_t recoBin = scfMap_->binNumber( xCoord, yCoord );
1608 
1609  // Find out which true Bins contribute to the given reco bin
1610  const std::vector<Int_t>* trueBins = scfMap_->trueBins( recoBin );
1611 
1612  Int_t nDataEvents = this->eventsPerExpt();
1613 
1614  // Loop over the true bins
1615  for ( std::vector<Int_t>::const_iterator iter = trueBins->begin(); iter != trueBins->end();
1616  ++iter ) {
1617  Int_t trueBin = ( *iter );
1618 
1619  // prob of a true event in the given true bin migrating to the reco bin
1620  Double_t pRecoGivenTrue = scfMap_->prob( recoBin, trueBin );
1621 
1622  // We've cached the DP amplitudes and the efficiency for the
1623  // true bin centres, just after the data points
1624  sigDPModel_->calcLikelihoodInfo( nDataEvents + trueBin );
1625  Double_t pTrue = sigDPModel_->getEvtLikelihood();
1626 
1627  // Get the cached SCF fraction (and jacobian if we're using the square DP)
1628  Double_t scfFraction = fakeSCFFracs_[trueBin];
1629  Double_t jacobian( 1.0 );
1630  if ( squareDP ) {
1631  jacobian = fakeJacobians_[trueBin];
1632  }
1633 
1634  scfDPLike += pTrue * jacobian * scfFraction * pRecoGivenTrue;
1635  }
1636 
1637  // Divide by the reco jacobian
1638  scfDPLike /= recoJacobian;
1639 
1640  return scfDPLike;
1641 }
1642 
1644 {
1645  // Function to return the signal and background likelihoods for the
1646  // extra variables for the given event evtNo.
1647 
1648  sigExtraLike_ = 1.0;
1649 
1650  for ( LauPdfPList::iterator iter = signalPdfs_.begin(); iter != signalPdfs_.end(); ++iter ) {
1651  ( *iter )->calcLikelihoodInfo( iEvt );
1652  sigExtraLike_ *= ( *iter )->getLikelihood();
1653  }
1654  if ( useSCF_ ) {
1655  scfExtraLike_ = 1.0;
1656  for ( LauPdfPList::iterator iter = scfPdfs_.begin(); iter != scfPdfs_.end(); ++iter ) {
1657  ( *iter )->calcLikelihoodInfo( iEvt );
1658  scfExtraLike_ *= ( *iter )->getLikelihood();
1659  }
1660  }
1661 
1662  const UInt_t nBkgnds = this->nBkgndClasses();
1663  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
1664  if ( usingBkgnd_ ) {
1665  bkgndExtraLike_[bkgndID] = 1.0;
1666  LauPdfPList& pdfList = bkgndPdfs_[bkgndID];
1667  for ( LauPdfPList::iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end();
1668  ++pdf_iter ) {
1669  ( *pdf_iter )->calcLikelihoodInfo( iEvt );
1670  bkgndExtraLike_[bkgndID] *= ( *pdf_iter )->getLikelihood();
1671  }
1672  } else {
1673  bkgndExtraLike_[bkgndID] = 0.0;
1674  }
1675  }
1676 }
1677 
1679 {
1680  coeffs_.clear();
1681  coeffs_.reserve( nSigComp_ );
1682  for ( UInt_t i = 0; i < nSigComp_; i++ ) {
1683  coeffs_.push_back( coeffPars_[i]->particleCoeff() );
1684  }
1685 }
1686 
1688 {
1689  // add branches for storing the experiment number and the number of
1690  // the event within the current experiment
1691  this->addSPlotNtupleIntegerBranch( "iExpt" );
1692  this->addSPlotNtupleIntegerBranch( "iEvtWithinExpt" );
1693 
1694  // Store the efficiency of the event (for inclusive BF calculations).
1695  if ( this->storeDPEff() ) {
1696  this->addSPlotNtupleDoubleBranch( "efficiency" );
1697  if ( sigDPModel_->usingScfModel() ) {
1698  this->addSPlotNtupleDoubleBranch( "scffraction" );
1699  }
1700  }
1701 
1702  // Store the total event likelihood for each species.
1703  if ( useSCF_ ) {
1704  this->addSPlotNtupleDoubleBranch( "sigTMTotalLike" );
1705  this->addSPlotNtupleDoubleBranch( "sigSCFTotalLike" );
1706  this->addSPlotNtupleDoubleBranch( "sigSCFFrac" );
1707  } else {
1708  this->addSPlotNtupleDoubleBranch( "sigTotalLike" );
1709  }
1710  if ( usingBkgnd_ ) {
1711  const UInt_t nBkgnds = this->nBkgndClasses();
1712  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1713  TString name( this->bkgndClassName( iBkgnd ) );
1714  name += "TotalLike";
1715  this->addSPlotNtupleDoubleBranch( name );
1716  }
1717  }
1718 
1719  // Store the DP likelihoods
1720  if ( this->useDP() ) {
1721  if ( useSCF_ ) {
1722  this->addSPlotNtupleDoubleBranch( "sigTMDPLike" );
1723  this->addSPlotNtupleDoubleBranch( "sigSCFDPLike" );
1724  } else {
1725  this->addSPlotNtupleDoubleBranch( "sigDPLike" );
1726  }
1727  if ( usingBkgnd_ ) {
1728  const UInt_t nBkgnds = this->nBkgndClasses();
1729  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1730  TString name( this->bkgndClassName( iBkgnd ) );
1731  name += "DPLike";
1732  this->addSPlotNtupleDoubleBranch( name );
1733  }
1734  }
1735  }
1736 
1737  // Store the likelihoods for each extra PDF
1738  if ( useSCF_ ) {
1739  this->addSPlotNtupleBranches( &signalPdfs_, "sigTM" );
1740  this->addSPlotNtupleBranches( &scfPdfs_, "sigSCF" );
1741  } else {
1742  this->addSPlotNtupleBranches( &signalPdfs_, "sig" );
1743  }
1744  if ( usingBkgnd_ ) {
1745  const UInt_t nBkgnds = this->nBkgndClasses();
1746  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1747  const TString& bkgndClass = this->bkgndClassName( iBkgnd );
1748  const LauPdfPList* pdfList = &( bkgndPdfs_[iBkgnd] );
1749  this->addSPlotNtupleBranches( pdfList, bkgndClass );
1750  }
1751  }
1752 }
1753 
1754 void LauSimpleFitModel::addSPlotNtupleBranches( const LauPdfPList* extraPdfs, const TString& prefix )
1755 {
1756  if ( extraPdfs ) {
1757  // Loop through each of the PDFs
1758  for ( LauPdfPList::const_iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end();
1759  ++pdf_iter ) {
1760 
1761  // Count the number of input variables that are not
1762  // DP variables (used in the case where there is DP
1763  // dependence for e.g. MVA)
1764  UInt_t nVars( 0 );
1765  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1766  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1767  var_iter != varNames.end();
1768  ++var_iter ) {
1769  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
1770  ++nVars;
1771  }
1772  }
1773 
1774  if ( nVars == 1 ) {
1775  // If the PDF only has one variable then
1776  // simply add one branch for that variable
1777  TString varName = ( *pdf_iter )->varName();
1778  TString name( prefix );
1779  name += varName;
1780  name += "Like";
1781  this->addSPlotNtupleDoubleBranch( name );
1782  } else if ( nVars == 2 ) {
1783  // If the PDF has two variables then we
1784  // need a branch for them both together and
1785  // branches for each
1786  TString allVars( "" );
1787  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1788  var_iter != varNames.end();
1789  ++var_iter ) {
1790  allVars += ( *var_iter );
1791  TString name( prefix );
1792  name += ( *var_iter );
1793  name += "Like";
1794  this->addSPlotNtupleDoubleBranch( name );
1795  }
1796  TString name( prefix );
1797  name += allVars;
1798  name += "Like";
1799  this->addSPlotNtupleDoubleBranch( name );
1800  } else {
1801  std::cerr << "WARNING in LauSimpleFitModel::addSPlotNtupleBranches : Can't yet deal with 3D PDFs."
1802  << std::endl;
1803  }
1804  }
1805  }
1806 }
1807 
1809  const TString& prefix,
1810  UInt_t iEvt )
1811 {
1812  // Store the PDF value for each variable in the list
1813  Double_t totalLike( 1.0 );
1814  Double_t extraLike( 0.0 );
1815  if ( extraPdfs ) {
1816  for ( LauPdfPList::iterator pdf_iter = extraPdfs->begin(); pdf_iter != extraPdfs->end();
1817  ++pdf_iter ) {
1818 
1819  // calculate the likelihood for this event
1820  ( *pdf_iter )->calcLikelihoodInfo( iEvt );
1821  extraLike = ( *pdf_iter )->getLikelihood();
1822  totalLike *= extraLike;
1823 
1824  // Count the number of input variables that are not
1825  // DP variables (used in the case where there is DP
1826  // dependence for e.g. MVA)
1827  UInt_t nVars( 0 );
1828  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1829  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1830  var_iter != varNames.end();
1831  ++var_iter ) {
1832  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
1833  ++nVars;
1834  }
1835  }
1836 
1837  if ( nVars == 1 ) {
1838  // If the PDF only has one variable then
1839  // simply store the value for that variable
1840  TString varName = ( *pdf_iter )->varName();
1841  TString name( prefix );
1842  name += varName;
1843  name += "Like";
1844  this->setSPlotNtupleDoubleBranchValue( name, extraLike );
1845  } else if ( nVars == 2 ) {
1846  // If the PDF has two variables then we
1847  // store the value for them both together
1848  // and for each on their own
1849  TString allVars( "" );
1850  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1851  var_iter != varNames.end();
1852  ++var_iter ) {
1853  allVars += ( *var_iter );
1854  TString name( prefix );
1855  name += ( *var_iter );
1856  name += "Like";
1857  Double_t indivLike = ( *pdf_iter )->getLikelihood( ( *var_iter ) );
1858  this->setSPlotNtupleDoubleBranchValue( name, indivLike );
1859  }
1860  TString name( prefix );
1861  name += allVars;
1862  name += "Like";
1863  this->setSPlotNtupleDoubleBranchValue( name, extraLike );
1864  } else {
1865  std::cerr << "WARNING in LauSimpleFitModel::setSPlotNtupleBranchValues : Can't yet deal with 3D PDFs."
1866  << std::endl;
1867  }
1868  }
1869  }
1870  return totalLike;
1871 }
1872 
1874 {
1875  LauSPlot::NameSet nameSet;
1876  if ( this->useDP() ) {
1877  nameSet.insert( "DP" );
1878  }
1879  // Loop through all the signal PDFs
1880  for ( LauPdfPList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end();
1881  ++pdf_iter ) {
1882  // Loop over the variables involved in each PDF
1883  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1884  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1885  var_iter != varNames.end();
1886  ++var_iter ) {
1887  // If they are not DP coordinates then add them
1888  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
1889  nameSet.insert( ( *var_iter ) );
1890  }
1891  }
1892  }
1893  return nameSet;
1894 }
1895 
1897 {
1898  LauSPlot::NumbMap numbMap;
1899  if ( ! signalEvents_->fixed() && this->doEMLFit() ) {
1900  numbMap["sig"] = signalEvents_->genValue();
1901  }
1902  if ( usingBkgnd_ ) {
1903  const UInt_t nBkgnds = this->nBkgndClasses();
1904  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1905  const TString& bkgndClass = this->bkgndClassName( iBkgnd );
1906  const LauAbsRValue* par = bkgndEvents_[iBkgnd];
1907  if ( ! par->fixed() ) {
1908  numbMap[bkgndClass] = par->genValue();
1909  if ( ! par->isLValue() ) {
1910  std::cerr << "WARNING in LauSimpleFitModel::freeSpeciesNames : \"" << par->name()
1911  << "\" is a LauFormulaPar, which implies it is perhaps not entirely free to float in the fit, so the sWeight calculation may not be reliable"
1912  << std::endl;
1913  }
1914  }
1915  }
1916  }
1917  return numbMap;
1918 }
1919 
1921 {
1922  LauSPlot::NumbMap numbMap;
1923  if ( signalEvents_->fixed() && this->doEMLFit() ) {
1924  numbMap["sig"] = signalEvents_->genValue();
1925  }
1926  if ( usingBkgnd_ ) {
1927  const UInt_t nBkgnds = this->nBkgndClasses();
1928  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1929  const TString& bkgndClass = this->bkgndClassName( iBkgnd );
1930  const LauAbsRValue* par = bkgndEvents_[iBkgnd];
1931  if ( par->fixed() ) {
1932  numbMap[bkgndClass] = par->genValue();
1933  }
1934  }
1935  }
1936  return numbMap;
1937 }
1938 
1940 {
1941  LauSPlot::TwoDMap twodimMap;
1942 
1943  for ( LauPdfPList::const_iterator pdf_iter = signalPdfs_.begin(); pdf_iter != signalPdfs_.end();
1944  ++pdf_iter ) {
1945  // Count the number of input variables that are not DP variables
1946  UInt_t nVars( 0 );
1947  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1948  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1949  var_iter != varNames.end();
1950  ++var_iter ) {
1951  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
1952  ++nVars;
1953  }
1954  }
1955  if ( nVars == 2 ) {
1956  if ( useSCF_ ) {
1957  twodimMap.insert(
1958  std::make_pair( "sigTM", std::make_pair( varNames[0], varNames[1] ) ) );
1959  } else {
1960  twodimMap.insert(
1961  std::make_pair( "sig", std::make_pair( varNames[0], varNames[1] ) ) );
1962  }
1963  }
1964  }
1965 
1966  if ( useSCF_ ) {
1967  for ( LauPdfPList::const_iterator pdf_iter = scfPdfs_.begin(); pdf_iter != scfPdfs_.end();
1968  ++pdf_iter ) {
1969  // Count the number of input variables that are not DP variables
1970  UInt_t nVars( 0 );
1971  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1972  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1973  var_iter != varNames.end();
1974  ++var_iter ) {
1975  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
1976  ++nVars;
1977  }
1978  }
1979  if ( nVars == 2 ) {
1980  twodimMap.insert(
1981  std::make_pair( "sigSCF", std::make_pair( varNames[0], varNames[1] ) ) );
1982  }
1983  }
1984  }
1985 
1986  if ( usingBkgnd_ ) {
1987  const UInt_t nBkgnds = this->nBkgndClasses();
1988  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
1989  const TString& bkgndClass = this->bkgndClassName( iBkgnd );
1990  const LauPdfPList& pdfList = bkgndPdfs_[iBkgnd];
1991  for ( LauPdfPList::const_iterator pdf_iter = pdfList.begin(); pdf_iter != pdfList.end();
1992  ++pdf_iter ) {
1993  // Count the number of input variables that are not DP variables
1994  UInt_t nVars( 0 );
1995  std::vector<TString> varNames = ( *pdf_iter )->varNames();
1996  for ( std::vector<TString>::const_iterator var_iter = varNames.begin();
1997  var_iter != varNames.end();
1998  ++var_iter ) {
1999  if ( ( *var_iter ) != "m13Sq" && ( *var_iter ) != "m23Sq" ) {
2000  ++nVars;
2001  }
2002  }
2003  if ( nVars == 2 ) {
2004  twodimMap.insert(
2005  std::make_pair( bkgndClass, std::make_pair( varNames[0], varNames[1] ) ) );
2006  }
2007  }
2008  }
2009  }
2010 
2011  return twodimMap;
2012 }
2013 
2015 {
2016  std::cout << "INFO in LauSimpleFitModel::storePerEvtLlhds : Storing per-event likelihood values..."
2017  << std::endl;
2018 
2019  // if we've not been using the DP model then we need to cache all
2020  // the info here so that we can get the efficiency from it
2021 
2022  LauFitDataTree* inputFitData = this->fitData();
2023 
2024  if ( ! this->useDP() && this->storeDPEff() ) {
2026  sigDPModel_->fillDataTree( *inputFitData );
2027  }
2028 
2029  UInt_t evtsPerExpt( this->eventsPerExpt() );
2030 
2031  for ( UInt_t iEvt = 0; iEvt < evtsPerExpt; ++iEvt ) {
2032 
2033  this->setSPlotNtupleIntegerBranchValue( "iExpt", this->iExpt() );
2034  this->setSPlotNtupleIntegerBranchValue( "iEvtWithinExpt", iEvt );
2035 
2036  // the DP information
2037  this->getEvtDPLikelihood( iEvt );
2038  if ( this->storeDPEff() ) {
2039  if ( ! this->useDP() ) {
2041  }
2042  this->setSPlotNtupleDoubleBranchValue( "efficiency", sigDPModel_->getEvtEff() );
2043  if ( sigDPModel_->usingScfModel() ) {
2044  this->setSPlotNtupleDoubleBranchValue( "scffraction",
2046  }
2047  }
2048  if ( this->useDP() ) {
2050  if ( useSCF_ ) {
2052  this->setSPlotNtupleDoubleBranchValue( "sigTMDPLike", sigDPLike_ );
2053  this->setSPlotNtupleDoubleBranchValue( "sigSCFDPLike", scfDPLike_ );
2054  } else {
2055  this->setSPlotNtupleDoubleBranchValue( "sigDPLike", sigDPLike_ );
2056  }
2057  if ( usingBkgnd_ ) {
2058  const UInt_t nBkgnds = this->nBkgndClasses();
2059  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
2060  TString name = this->bkgndClassName( iBkgnd );
2061  name += "DPLike";
2062  this->setSPlotNtupleDoubleBranchValue( name, bkgndDPLike_[iBkgnd] );
2063  }
2064  }
2065  } else {
2066  sigTotalLike_ = 1.0;
2067  if ( useSCF_ ) {
2068  scfTotalLike_ = 1.0;
2069  }
2070  if ( usingBkgnd_ ) {
2071  const UInt_t nBkgnds = this->nBkgndClasses();
2072  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
2073  bkgndTotalLike_[iBkgnd] = 1.0;
2074  }
2075  }
2076  }
2077 
2078  // the signal PDF values
2079  if ( useSCF_ ) {
2080  sigTotalLike_ *= this->setSPlotNtupleBranchValues( &signalPdfs_, "sigTM", iEvt );
2081  scfTotalLike_ *= this->setSPlotNtupleBranchValues( &scfPdfs_, "sigSCF", iEvt );
2082  } else {
2083  sigTotalLike_ *= this->setSPlotNtupleBranchValues( &signalPdfs_, "sig", iEvt );
2084  }
2085 
2086  // the background PDF values
2087  if ( usingBkgnd_ ) {
2088  const UInt_t nBkgnds = this->nBkgndClasses();
2089  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
2090  const TString& bkgndClass = this->bkgndClassName( iBkgnd );
2091  LauPdfPList& pdfs = bkgndPdfs_[iBkgnd];
2092  bkgndTotalLike_[iBkgnd] *=
2093  this->setSPlotNtupleBranchValues( &( pdfs ), bkgndClass, iEvt );
2094  }
2095  }
2096 
2097  // the total likelihoods
2098  if ( useSCF_ ) {
2099  Double_t scfFrac( 0.0 );
2100  if ( useSCFHist_ ) {
2101  scfFrac = recoSCFFracs_[iEvt];
2102  } else {
2103  scfFrac = scfFrac_.unblindValue();
2104  }
2105  this->setSPlotNtupleDoubleBranchValue( "sigSCFFrac", scfFrac );
2106  sigTotalLike_ *= ( 1.0 - scfFrac );
2107  if ( scfMap_ == 0 ) {
2108  scfTotalLike_ *= scfFrac;
2109  }
2110  this->setSPlotNtupleDoubleBranchValue( "sigTMTotalLike", sigTotalLike_ );
2111  this->setSPlotNtupleDoubleBranchValue( "sigSCFTotalLike", scfTotalLike_ );
2112  } else {
2113  this->setSPlotNtupleDoubleBranchValue( "sigTotalLike", sigTotalLike_ );
2114  }
2115  if ( usingBkgnd_ ) {
2116  const UInt_t nBkgnds = this->nBkgndClasses();
2117  for ( UInt_t iBkgnd( 0 ); iBkgnd < nBkgnds; ++iBkgnd ) {
2118  TString name = this->bkgndClassName( iBkgnd );
2119  name += "TotalLike";
2120  this->setSPlotNtupleDoubleBranchValue( name, bkgndTotalLike_[iBkgnd] );
2121  }
2122  }
2123  // fill the tree
2124  this->fillSPlotNtupleBranches();
2125  }
2126  std::cout << "INFO in LauSimpleFitModel::storePerEvtLlhds : Finished storing per-event likelihood values."
2127  << std::endl;
2128 }
2129 
2130 std::map<TString, LauComplex> LauSimpleFitModel::getDPAmps( const Double_t m13Sq,
2131  const Double_t m23Sq )
2132 {
2133  // Initialise the DP model, if not already done
2134  if ( coeffs_.empty() ) {
2135  this->updateCoeffs();
2136  this->initialiseDPModels();
2137  }
2138 
2139  if ( ! kinematics_->withinDPLimits( m13Sq, m23Sq ) ) {
2140  return {
2141  {"signal", {}}
2142  };
2143  }
2144 
2145  sigDPModel_->calcLikelihoodInfo( m13Sq, m23Sq );
2146 
2147  return {
2148  {"signal", sigDPModel_->getEvtDPAmp()}
2149  };
2150 }
2151 
2152 std::map<TString, Double_t> LauSimpleFitModel::getDPLikelihoods( const Double_t m13Sq,
2153  const Double_t m23Sq )
2154 {
2155  // Initialise the DP model, if not already done
2156  if ( coeffs_.empty() ) {
2157  this->updateCoeffs();
2158  this->initialiseDPModels();
2159  }
2160 
2161  std::map<TString, Double_t> likelihoods;
2162 
2163  if ( ! kinematics_->withinDPLimits( m13Sq, m23Sq ) ) {
2164  likelihoods.emplace( "signal", 0.0 );
2165  if ( usingBkgnd_ ) {
2166  const UInt_t nBkgnds { this->nBkgndClasses() };
2167  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
2168  likelihoods.emplace( this->bkgndClassName( bkgndID ), 0.0 );
2169  }
2170  }
2171  return likelihoods;
2172  }
2173 
2174  sigDPModel_->calcLikelihoodInfo( m13Sq, m23Sq );
2175  likelihoods.emplace( "signal", sigDPModel_->getEvtLikelihood() );
2176 
2177  // TODO - SCF signal
2178  static bool warningIssued { false };
2179  if ( useSCF_ && ! warningIssued ) {
2180  warningIssued = true;
2181  std::cerr << "WARNING in LauSimpleFitModel::getDPLikelihoods : calculation of SCF likelihoods not currently implemented in this function\n";
2182  std::cerr << " : signal likelihood will just be the truth-matched value";
2183  std::cerr << std::endl;
2184  }
2185 
2186  if ( usingBkgnd_ ) {
2187  const UInt_t nBkgnds { this->nBkgndClasses() };
2188  for ( UInt_t bkgndID( 0 ); bkgndID < nBkgnds; ++bkgndID ) {
2189  likelihoods.emplace( this->bkgndClassName( bkgndID ),
2190  bkgndDPModels_[bkgndID]->getLikelihood( m13Sq, m23Sq ) );
2191  }
2192  }
2193 
2194  return likelihoods;
2195 }
2196 
2197 void LauSimpleFitModel::embedSignal( const TString& fileName,
2198  const TString& treeName,
2199  Bool_t reuseEventsWithinEnsemble,
2200  Bool_t reuseEventsWithinExperiment,
2201  Bool_t useReweighting )
2202 {
2203  if ( signalTree_ ) {
2204  std::cerr << "ERROR in LauSimpleFitModel::embedSignal : Already embedding signal from a file."
2205  << std::endl;
2206  return;
2207  }
2208 
2209  signalTree_ = new LauEmbeddedData( fileName, treeName, reuseEventsWithinExperiment );
2210 
2211  Bool_t dataOK = signalTree_->findBranches();
2212  if ( ! dataOK ) {
2213  delete signalTree_;
2214  signalTree_ = 0;
2215  std::cerr << "ERROR in LauSimpleFitModel::embedSignal : Problem creating data tree for embedding."
2216  << std::endl;
2217  return;
2218  }
2219 
2220  reuseSignal_ = reuseEventsWithinEnsemble;
2221  useReweighting_ = useReweighting;
2222 
2223  if ( this->enableEmbedding() == kFALSE ) {
2224  this->enableEmbedding( kTRUE );
2225  }
2226 }
2227 
2228 void LauSimpleFitModel::embedBkgnd( const TString& bkgndClass,
2229  const TString& fileName,
2230  const TString& treeName,
2231  Bool_t reuseEventsWithinEnsemble,
2232  Bool_t reuseEventsWithinExperiment )
2233 {
2234  if ( ! this->validBkgndClass( bkgndClass ) ) {
2235  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Invalid background class \""
2236  << bkgndClass << "\"." << std::endl;
2237  std::cerr << " : Background class names must be provided in \"setBkgndClassNames\" before any other background-related actions can be performed."
2238  << std::endl;
2239  return;
2240  }
2241 
2242  UInt_t bkgndID = this->bkgndClassID( bkgndClass );
2243 
2244  if ( bkgndTree_[bkgndID] ) {
2245  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Already embedding background from a file."
2246  << std::endl;
2247  return;
2248  }
2249 
2250  bkgndTree_[bkgndID] = new LauEmbeddedData( fileName, treeName, reuseEventsWithinExperiment );
2251 
2252  Bool_t dataOK = bkgndTree_[bkgndID]->findBranches();
2253  if ( ! dataOK ) {
2254  delete bkgndTree_[bkgndID];
2255  bkgndTree_[bkgndID] = 0;
2256  std::cerr << "ERROR in LauSimpleFitModel::embedBkgnd : Problem creating data tree for embedding."
2257  << std::endl;
2258  return;
2259  }
2260 
2261  reuseBkgnd_[bkgndID] = reuseEventsWithinEnsemble;
2262 
2263  if ( this->enableEmbedding() == kFALSE ) {
2264  this->enableEmbedding( kTRUE );
2265  }
2266 }
2267 
2268 void LauSimpleFitModel::weightEvents( const TString& dataFileName, const TString& dataTreeName )
2269 {
2270  // Routine to provide weights for events that are uniformly distributed
2271  // in the DP (or square DP) so as to reproduce the given DP model
2272 
2273  const Bool_t squareDP { kinematics_->squareDP() };
2274  if ( squareDP ) {
2275  std::cout << "INFO in LauSimpleFitModel::weightEvents : will store DP model weights and the square DP jacobian\n";
2276  std::cout << " : the DP model weights can be used on their own to weight events that were generated flat in phase space\n";
2277  std::cout << " : or they can be multiplied by the jacobian to weight events that were generated flat in the square DP\n";
2278  std::cout << " : or they can be multiplied by max(1.0, jacobian) to weight events that were generated quasi-flat in the square DP"
2279  << std::endl;
2280  } else {
2281  std::cout << "INFO in LauSimpleFitModel::weightEvents : will store DP model weights suitable for weighting events that were generated flat in phase space"
2282  << std::endl;
2283  }
2284 
2285  // This reads in the given dataFile and creates an input
2286  // fit data tree that stores them for all events and experiments.
2287  const Bool_t dataOK { this->verifyFitData( dataFileName, dataTreeName ) };
2288  if ( ! dataOK ) {
2289  std::cerr << "ERROR in LauSimpleFitModel::weightEvents : Problem caching the data."
2290  << std::endl;
2291  return;
2292  }
2293 
2294  LauFitDataTree* inputFitData { this->fitData() };
2295 
2296  if ( ! inputFitData->haveBranch( "m13Sq_MC" ) || ! inputFitData->haveBranch( "m23Sq_MC" ) ) {
2297  std::cerr << "ERROR in LauSimpleFitModel::weightEvents : Cannot find MC truth DP coordinate branches in supplied data, aborting."
2298  << std::endl;
2299  return;
2300  }
2301 
2302  // Create the ntuple to hold the DP weights
2303  TString weightsFileName { dataFileName };
2304  const Ssiz_t index { weightsFileName.Last( '.' ) };
2305  weightsFileName.Insert( index, "_DPweights" );
2306 
2307  LauGenNtuple weightsTuple { weightsFileName, dataTreeName };
2308  weightsTuple.addIntegerBranch( "iExpt" );
2309  weightsTuple.addIntegerBranch( "iEvtWithinExpt" );
2310  weightsTuple.addDoubleBranch( "dpModelWeight" );
2311  if ( squareDP ) {
2312  weightsTuple.addDoubleBranch( "sqDPJacobian" );
2313  }
2314 
2315  const UInt_t nExpmt { this->nExpt() };
2316  const UInt_t firstExpmt { this->firstExpt() };
2317  for ( UInt_t iExpmt { firstExpmt }; iExpmt < ( firstExpmt + nExpmt ); ++iExpmt ) {
2318 
2319  inputFitData->readExperimentData( iExpmt );
2320 
2321  const UInt_t nEvents { inputFitData->nEvents() };
2322  if ( nEvents < 1 ) {
2323  std::cerr << "WARNING in LauSimpleFitModel::weightEvents : Zero events in experiment "
2324  << iExpmt << ", skipping..." << std::endl;
2325  continue;
2326  }
2327 
2328  weightsTuple.setIntegerBranchValue( "iExpt", iExpmt );
2329 
2330  // Calculate and store the weights for the events in this experiment
2331  for ( UInt_t iEvent { 0 }; iEvent < nEvents; ++iEvent ) {
2332 
2333  weightsTuple.setIntegerBranchValue( "iEvtWithinExpt", iEvent );
2334 
2335  const LauFitData& evtData { inputFitData->getData( iEvent ) };
2336 
2337  const Double_t m13Sq_MC { evtData.at( "m13Sq_MC" ) };
2338  const Double_t m23Sq_MC { evtData.at( "m23Sq_MC" ) };
2339 
2340  Double_t dpModelWeight { 0.0 };
2341  Double_t jacobian { 1.0 };
2342 
2343  if ( kinematics_->withinDPLimits( m13Sq_MC, m23Sq_MC ) ) {
2344 
2345  kinematics_->updateKinematics( m13Sq_MC, m23Sq_MC );
2346  dpModelWeight = sigDPModel_->getEventWeight();
2347 
2348  if ( squareDP ) {
2349  jacobian = kinematics_->calcSqDPJacobian();
2350  }
2351 
2352  dpModelWeight /= sigDPModel_->getDPNorm();
2353  }
2354 
2355  weightsTuple.setDoubleBranchValue( "dpModelWeight", dpModelWeight );
2356  if ( squareDP ) {
2357  weightsTuple.setDoubleBranchValue( "sqDPJacobian", jacobian );
2358  }
2359  weightsTuple.fillBranches();
2360  }
2361  }
2362 
2363  weightsTuple.buildIndex( "iExpt", "iEvtWithinExpt" );
2364  weightsTuple.addFriendTree( dataFileName, dataTreeName );
2365  weightsTuple.writeOutGenResults();
2366 }
2367 
2368 void LauSimpleFitModel::savePDFPlots( const TString& label )
2369 {
2370  savePDFPlotsWave( label, 0 );
2371  savePDFPlotsWave( label, 1 );
2372  savePDFPlotsWave( label, 2 );
2373 
2374  std::cout << "LauCPFitModel::plot" << std::endl;
2375  // ((LauIsobarDynamics*)sigDPModel_)->plot();
2376 
2377  //Double_t minm13 = sigDPModel_->getKinematics()->getm13Min();
2378  Double_t minm13 = 0.0;
2379  Double_t maxm13 = sigDPModel_->getKinematics()->getm13Max();
2380  //Double_t minm23 = sigDPModel_->getKinematics()->getm23Min();
2381  Double_t minm23 = 0.0;
2382  Double_t maxm23 = sigDPModel_->getKinematics()->getm23Max();
2383 
2384  Double_t mins13 = minm13 * minm13;
2385  Double_t maxs13 = maxm13 * maxm13;
2386  Double_t mins23 = minm23 * minm23;
2387  Double_t maxs23 = maxm23 * maxm23;
2388 
2389  Double_t s13, s23, chPdf;
2390 
2391  Int_t n13 = 200.00, n23 = 200.00;
2392  Double_t delta13, delta23;
2393  delta13 = ( maxs13 - mins13 ) / n13;
2394  delta23 = ( maxs23 - mins23 ) / n23;
2395  UInt_t nAmp = sigDPModel_->getnCohAmp();
2396  for ( UInt_t resID = 0; resID <= nAmp; ++resID ) {
2397  TGraph2D* dt = new TGraph2D();
2398  TString resName = "TotalAmp";
2399  if ( resID != nAmp ) {
2400  TString tStrResID = Form( "%d", resID );
2401  const LauIsobarDynamics* model = sigDPModel_;
2402  const LauAbsResonance* resonance = model->getResonance( resID );
2403  resName = resonance->getResonanceName();
2404  std::cout << "resName = " << resName << std::endl;
2405  }
2406 
2407  resName.ReplaceAll( "(", "" );
2408  resName.ReplaceAll( ")", "" );
2409  resName.ReplaceAll( "*", "Star" );
2410 
2411  TCanvas* c = new TCanvas( "c" + resName + label, resName + " (" + label + ")", 0, 0, 600, 400 );
2412  dt->SetName( resName + label );
2413  dt->SetTitle( resName + " (" + label + ")" );
2414  Int_t count = 0;
2415  for ( Int_t i = 0; i < n13; i++ ) {
2416  s13 = mins13 + i * delta13;
2417  for ( Int_t j = 0; j < n23; j++ ) {
2418  s23 = mins23 + j * delta23;
2419  if ( sigDPModel_->getKinematics()->withinDPLimits2( s23, s13 ) ) {
2420  //if (s13 > 4) continue;
2421  sigDPModel_->calcLikelihoodInfo( s13, s23 );
2422  LauComplex chAmp = sigDPModel_->getEvtDPAmp();
2423  if ( resID != nAmp ) {
2424  chAmp = sigDPModel_->getFullAmplitude( resID );
2425  }
2426  chPdf = chAmp.abs2();
2427  //if ((z > 0.04)||(z < -0.03)) continue;
2428  //z = TMath::Log(z);
2430  Double_t sLow = s13;
2431  Double_t sHigh = s23;
2432  if ( sLow > sHigh ) {
2433  continue;
2434  //sLow = s23;
2435  //sHigh = s13;
2436  }
2437 
2438  //if (sLow > 3.5) continue;
2439  //if (i<10) std::cout << "SymmetricalDP" << std::endl;
2440  //if (z>0.02) z = 0.02;
2441  //if (z<-0.02) z = -0.02;
2442 
2443  dt->SetPoint( count, sHigh, sLow, chPdf );
2444  count++;
2445  } else {
2446  dt->SetPoint( count, s13, s23, chPdf );
2447  count++;
2448  }
2449  }
2450  }
2451  }
2452  gStyle->SetPalette( 1 );
2453  dt->GetXaxis()->SetTitle( "m_{K#pi}(low)" );
2454  dt->GetYaxis()->SetTitle( "m_{K#pi}(high)" );
2455  dt->Draw( "SURF1" );
2456  dt->GetXaxis()->SetTitle( "m_{K#pi}(low)" );
2457  dt->GetYaxis()->SetTitle( "m_{K#pi}(high)" );
2458  c->SaveAs( "plot_2D_" + resName + "_" + label + ".C" );
2459  }
2460 }
2461 
2462 void LauSimpleFitModel::savePDFPlotsWave( const TString& label, const Int_t& spin )
2463 {
2464 
2465  std::cout << "label = " << label << ", spin = " << spin << std::endl;
2466 
2467  TString tStrResID = "S_Wave";
2468  if ( spin == 1 )
2469  tStrResID = "P_Wave";
2470  if ( spin == 2 )
2471  tStrResID = "D_Wave";
2472 
2473  std::cout << "LauSimpleFitModel::savePDFPlotsWave: " << tStrResID << std::endl;
2474 
2475  TCanvas* c = new TCanvas( "c" + tStrResID + label, tStrResID + " (" + label + ")", 0, 0, 600, 400 );
2476 
2477  Double_t minm13 = 0.0;
2478  Double_t maxm13 = sigDPModel_->getKinematics()->getm13Max();
2479  Double_t minm23 = 0.0;
2480  Double_t maxm23 = sigDPModel_->getKinematics()->getm23Max();
2481 
2482  Double_t mins13 = minm13 * minm13;
2483  Double_t maxs13 = maxm13 * maxm13;
2484  Double_t mins23 = minm23 * minm23;
2485  Double_t maxs23 = maxm23 * maxm23;
2486 
2487  Double_t s13, s23, chPdf;
2488  TGraph2D* dt = new TGraph2D();
2489  dt->SetName( tStrResID + label );
2490  dt->SetTitle( tStrResID + " (" + label + ")" );
2491 
2492  Int_t n13 = 200.00, n23 = 200.00;
2493  Double_t delta13, delta23;
2494  delta13 = ( maxs13 - mins13 ) / n13;
2495  delta23 = ( maxs23 - mins23 ) / n23;
2496  UInt_t nAmp = sigDPModel_->getnCohAmp();
2497 
2498  Int_t count = 0;
2499  for ( Int_t i = 0; i < n13; i++ ) {
2500  s13 = mins13 + i * delta13;
2501  for ( Int_t j = 0; j < n23; j++ ) {
2502  s23 = mins23 + j * delta23;
2503  if ( sigDPModel_->getKinematics()->withinDPLimits2( s23, s13 ) ) {
2504  //if (s13 > 4) continue;
2505  LauComplex chAmp( 0, 0 );
2506  Bool_t noWaveRes = kTRUE;
2507  sigDPModel_->calcLikelihoodInfo( s13, s23 );
2508  for ( UInt_t resID = 0; resID < nAmp; ++resID ) {
2509  const LauIsobarDynamics* model = sigDPModel_;
2510  const LauAbsResonance* resonance = model->getResonance( resID );
2511  Int_t spin_res = resonance->getSpin();
2512  if ( spin != spin_res )
2513  continue;
2514  noWaveRes = kFALSE;
2515  chAmp += sigDPModel_->getFullAmplitude( resID );
2516  }
2517 
2518  if ( noWaveRes )
2519  return;
2520 
2521  chPdf = chAmp.abs2();
2523  Double_t sLow = s13;
2524  Double_t sHigh = s23;
2525  if ( sLow > sHigh ) {
2526  continue;
2527  //sLow = s23;
2528  //sHigh = s13;
2529  }
2530  dt->SetPoint( count, sHigh, sLow, chPdf );
2531  count++;
2532  } else {
2533  dt->SetPoint( count, s13, s23, chPdf );
2534  count++;
2535  }
2536  }
2537  }
2538  }
2539  gStyle->SetPalette( 1 );
2540  dt->GetXaxis()->SetTitle( "pipi" );
2541  dt->GetYaxis()->SetTitle( "pipi" );
2542  dt->Draw( "SURF1" );
2543  dt->GetXaxis()->SetTitle( "pipi" );
2544  dt->GetYaxis()->SetTitle( "pipi" );
2545  c->SaveAs( "plot_2D_" + tStrResID + "_" + label + ".C" );
2546 }
std::map< TString, std::pair< Int_t, Double_t > > LauGenInfo
Define a map to be used to store a category name and numbers.
virtual void printParValues() const =0
Print the current values of the parameters.
File containing declaration of LauEffModel class.
File containing declaration of LauAbsPdf class.
void splitSignalComponent(const TH2 *dpHisto, const Bool_t upperHalf=kFALSE, const Bool_t fluctuateBins=kFALSE, LauScfMap *scfMap=0)
Split the signal component into well-reconstructed and mis-reconstructed parts.
Double_t getEventWeight()
Calculate the acceptance rate, for events with the current kinematics, when generating events accordi...
virtual void weightEvents(const TString &dataFileName, const TString &dataTreeName)
Weight events based on the DP model.
Double_t scfExtraLike_
SCF likelihood from extra PDFs.
File containing LauRandom namespace.
Bool_t doEMLFit() const
Determine whether an extended maximum likelihood fit it being performed.
virtual void updateCoeffs()
Update the coefficients.
virtual void getEvtDPLikelihood(UInt_t iEvt)
Calculate the signal and background likelihoods for the DP for a given event.
Int_t resonanceIndex(const TString &resName) const
Retrieve the index for the given resonance.
const LauAbsResonance * getResonance(const UInt_t resIndex) const
Retrieve a resonance by its index.
std::vector< Double_t > fakeSCFFracs_
The cached values of the SCF fraction for each bin centre.
Double_t scfDPLike_
SCF likelihood value.
UInt_t addResonanceParameters(LauParameter *param)
Add the given parameter to the list of resonance parameters and the list of all fit parameters.
void setSCFPdf(LauAbsPdf *pdf)
Set the SCF PDF for a given variable.
Class for defining the abstract interface for complex coefficient classes.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
const TString & getResonanceName() const
Get the name of the resonance.
TH2 * trueHist(Int_t trueBin)
Retrieve the migration histogram for trueBin.
Definition: LauScfMap.cc:190
virtual void setNSigEvents(LauParameter *nSigEvents)
Set the signal event yield.
Double_t unblindValue() const
The unblinded value of the parameter.
Double_t getc13() const
Get the cosine of the helicity angle theta13.
virtual Double_t getEventSum() const
Get the total number of events.
virtual void finaliseFitResults(const TString &tablePrefixName)
Get the fit results and store them.
const LauParArray & getFitFractionsEfficiencyUncorrected() const
Retrieve the fit fractions for the amplitude components.
Bool_t doPoissonSmearing() const
Determine whether Poisson smearing is enabled for the toy MC generation.
Double_t getc23() const
Get the cosine of the helicity angle theta23.
void listBinCentres(std::vector< Double_t > &xCoords, std::vector< Double_t > &yCoords) const
Create lists of the co-ordinates of the centres of true bins.
Definition: LauScfMap.cc:128
void calcLikelihoodInfo(const UInt_t iEvt)
Calculate the likelihood (and all associated information) for the given event number.
Double_t value() const
The value of the parameter.
Double_t getDPNorm() const
Retrieve the normalisation factor for the log-likelihood function.
virtual void setGenNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the gen tree.
virtual void initialiseDPModels()
Initialise the signal DP model.
Double_t scfTotalLike_
Total SCF likelihood.
void setBkgndPdf(const TString &bkgndClass, LauAbsPdf *pdf)
Set the background PDF.
Double_t getValue(const TString &name) const
Get the value of a specified branch.
Double_t getm12Sq() const
Get the m12 invariant mass square.
Double_t getEvtEff() const
Retrieve the efficiency for the current event.
virtual Bool_t verifyFitData(const TString &dataFileName, const TString &dataTreeName)
Open the input file and verify that all required variables are present.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
virtual void storePerEvtLlhds()
Store the per event likelihood values.
Double_t getThetaPrime() const
Get theta' value.
virtual void printFitFractions(std::ostream &output)
Print the fit fractions, total DP rate and mean efficiency.
UInt_t firstExpt() const
Obtain the number of the first experiment.
LauEmbeddedData * signalTree_
The signal event tree.
void recalculateNormalisation()
recalculate Normalization
The abstract interface for a background Dalitz plot model.
Double_t getmPrime() const
Get m' value.
LauParameter dpRate_
The average DP rate.
Double_t sigExtraLike_
Signal likelihood from extra PDFs.
LauIsobarDynamics * sigDPModel_
The signal Dalitz plot model.
virtual TString name() const
Retrieve the name of the coefficient set.
UInt_t nNormPar_
Number of normalisation parameters (i.e. yields)
Bool_t withinDPLimits2(const Double_t m13Sq, const Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot (alternat...
std::vector< Double_t > bkgndExtraLike_
Background likelihood value(s) from extra PDFs.
void appendFakePoints(const std::vector< Double_t > &xCoords, const std::vector< Double_t > &yCoords)
Add fake events to the data.
Bool_t storeDPEff() const
Determine whether the efficiency information should be stored in the sPlot ntuple.
File containing declaration of LauScfMap class.
void setExtraPdfParameters()
Set the fit parameters for the extra PDFs.
Bool_t squareDP() const
Determine to use or not the square Dalitz plot.
virtual void cacheInputFitVars()
Read in the input fit data variables, e.g. m13Sq and m23Sq.
Double_t calcEfficiency(const LauKinematics *kinematics) const
Determine the efficiency for a given point in the Dalitz plot.
Definition: LauEffModel.cc:367
Double_t getm13Max() const
Get the m13 maximum defined as (mParent - m2)
LauEffModel * scfFracHist_
The histogram giving the DP-dependence of the SCF fraction.
void cacheInfo(LauPdfPList &pdfList, const LauFitDataTree &theData)
Have all PDFs in the list cache the data.
void randomiseInitFitPars()
Randomise the initial fit parameters.
void setSignalPdf(LauAbsPdf *pdf)
Set the signal PDF for a given variable.
File containing declaration of LauFitNtuple class.
void clearUsedList()
Clear the list of used events.
Bool_t enableEmbedding() const
Determine whether embedding of events is enabled in the generation.
const LauFitDataTree * fitData() const
Const access the data store.
virtual Double_t getTotEvtLikelihood(UInt_t iEvt)
Get the total likelihood for each event.
void fillDataTree(const LauFitDataTree &fitDataTree)
Fill the internal data structure that caches the resonance dynamics.
virtual LauSPlot::NameSet variableNames() const
Returns the names of all variables in the fit.
void setDPBranchValues()
Store all of the DP information.
void setEffHisto(const TH2 *effHisto, Bool_t useInterpolation=kTRUE, Bool_t fluctuateBins=kFALSE, Double_t avEff=-1.0, Double_t absError=-1.0, Bool_t useUpperHalfOnly=kFALSE, Bool_t squareDP=kFALSE)
Set the efficiency variation across the phase space using a predetermined 2D histogram.
Definition: LauEffModel.cc:68
void updateKinematics(const Double_t m13Sq, const Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
virtual void addGenNtupleDoubleBranch(const TString &name)
Add a branch to the gen tree for storing a double.
virtual Bool_t blind() const =0
The blinding state.
File containing declaration of LauAbsCoeffSet class.
const TString & bkgndClassName(UInt_t classID) const
Get the name of a background class from the number.
virtual std::map< TString, Double_t > getDPLikelihoods(const Double_t m13Sq, const Double_t m23Sq)
Calculate the DP likelihood(s) for a given DP position.
virtual Double_t value() const =0
Return the value of the parameter.
Bool_t haveBranch(const TString &name) const
Boolean to determine whether branch exists.
Double_t getEvtmPrime() const
Retrieve the square Dalitz plot coordinate, m', for the current event.
Bool_t withinDPLimits(const Double_t m13Sq, const Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot.
const std::set< TString > & multiDimConstrainedPars() const
Const access to the parameter names used in ND constraints.
std::vector< LauParameter > LauParameterList
List of parameters.
virtual void savePDFPlots(const TString &label)
Save the pdf Plots for all the resonances of experiment number fitExp.
UInt_t bkgndClassID(const TString &className) const
The number assigned to a background class.
std::vector< LauParameter * > & getFloatingParameters()
Retrieve the floating parameters of the resonance models.
std::vector< LauParameter * > LauParameterPList
List of parameter pointers.
Double_t sigDPLike_
Signal likelihood value.
std::vector< LauAbsCoeffSet * > coeffPars_
Magnitudes and Phases.
LauParArray fitFracEffUnCorr_
Fit fractions (uncorrected for the efficiency)
Bool_t haveBranch(const TString &name) const
Check if the named branch is stored.
std::set< TString > NameSet
Type to store names, e.g. of the discriminating/control variables.
Definition: LauSPlot.hh:72
Int_t getSpin() const
Get the spin of the resonance.
LauBkgndPdfsList bkgndPdfs_
The background PDFs.
virtual void recalculateNormalisation()
Recalculate Normalization the signal DP models.
UInt_t nSigComp_
Number of signal components.
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
Class to store the input fit variables.
const LauAbsFitter::FitStatus & fitStatus() const
Access the fit status information.
Double_t getm23Max() const
Get the m23 maximum defined as (mParent - m1)
Bool_t generate()
Generate a toy MC signal event.
LauBkgndReuseEventsList reuseBkgnd_
Vector of booleans to reuse background events.
Bool_t useSCF_
Is the signal split into TM and SCF.
Class for defining a complex number.
Definition: LauComplex.hh:61
virtual void initialise()
Initialise the fit.
virtual void propagateParUpdates()
Calculate things that depend on the fit parameters after they have been updated by Minuit.
Bool_t useSCFHist_
Is the SCF fraction DP-dependent.
Double_t getc12() const
Get the cosine of the helicity angle theta12.
UInt_t getnTotAmp() const
Retrieve the total number of amplitude components.
virtual void savePDFPlotsWave(const TString &label, const Int_t &spin)
Save the pdf Plots for the sum of ressonances correspondint to "sin" of experiment number fitExp.
Class that implements the efficiency description across the signal Dalitz plot.
Definition: LauEffModel.hh:50
virtual Bool_t passVeto(const LauKinematics *kinematics) const =0
Determine whether the given DP position is outside the vetoes.
std::vector< Double_t > fakeJacobians_
The cached values of the sqDP jacobians for each true bin.
std::vector< Double_t > bkgndTotalLike_
Total background likelihood(s)
Bool_t blind() const
The blinding state.
Bool_t useRandomInitFitPars() const
Determine whether the initial values of the fit parameters, in particular the isobar coefficient para...
virtual LauSPlot::TwoDMap twodimPDFs() const
Returns the species and variables for all 2D PDFs in the fit.
Bool_t useReweighting_
Boolean to use reweighting.
virtual Bool_t generate()=0
Generate a toy MC event from the model.
UInt_t addFitParameters(LauParameter *param, const Bool_t addFixed=kFALSE)
Add the given parameter to the list of all fit parameters.
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:75
void squareDP(const Bool_t calcSquareDPCoords)
Enable/disable the calculation of square Dalitz plot co-ordinates.
void generateExtraPdfValues(LauPdfPList *extraPdfs, LauEmbeddedData *embeddedData)
Generate from the extra PDFs.
void embedBkgnd(const TString &bkgndClass, const TString &fileName, const TString &treeName, Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment=kFALSE)
Embed full simulation events for the given background class, rather than generating toy from the PDFs...
std::vector< Double_t > recoSCFFracs_
The cached values of the SCF fraction for each event.
std::vector< LauComplex > coeffs_
The complex coefficients.
virtual Double_t getEvtSCFDPLikelihood(UInt_t iEvt)
Calculate the SCF likelihood for the DP for a given event.
Double_t getm23() const
Get the m23 invariant mass.
Double_t setSPlotNtupleBranchValues(LauPdfPList *extraPdfs, const TString &prefix, UInt_t iEvt)
Set the branches for the sPlot ntuple with extra PDFs.
virtual void setSPlotNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the sPlot tree.
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:58
File containing declaration of LauComplex class.
UInt_t iExpt() const
Obtain the number of the current experiment.
virtual void fillSPlotNtupleBranches()
Fill the sPlot tuple.
virtual LauSPlot::NumbMap freeSpeciesNames() const
Returns the names and yields of species that are free in the fit.
File containing declaration of LauDaughters class.
Bool_t usingScfModel() const
Check whether a self cross feed fraction model is being used.
Bool_t generateBkgndEvent(UInt_t bkgndID)
Generate background event.
ToyMCStatus checkToyMC(Bool_t printErrorMessages=kTRUE, Bool_t printInfoMessages=kFALSE)
Check the status of the toy MC generation.
File containing declaration of LauGenNtuple class.
UInt_t nUsedEvents() const
Retrieve the number of events that have already been sampled.
LauParameter scfFrac_
The (global) SCF fraction parameter.
virtual void addSPlotNtupleIntegerBranch(const TString &name)
Add a branch to the sPlot tree for storing an integer.
UInt_t getnCohAmp() const
Retrieve the number of coherent amplitude components.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:33
const std::vector< LauParameter > & getExtraParameters() const
Retrieve any extra parameters/quantities (e.g. K-matrix total fit fractions)
std::vector< std::vector< LauParameter > > LauParArray
Type to define an array of parameters.
Double_t getEvtLikelihood() const
Retrieve the likelihood for the current event.
std::vector< Double_t > bkgndDPLike_
Background likelihood value(s)
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:54
virtual void setGenNtupleDoubleBranchValue(const TString &name, Double_t value)
Set the value of a double branch in the gen tree.
Double_t abs2() const
Obtain the square of the absolute value of the complex number.
Definition: LauComplex.hh:260
const LauFitNtuple * fitNtuple() const
Const access to the fit ntuple.
Double_t getm12() const
Get the m12 invariant mass.
LauSimpleFitModel(LauIsobarDynamics *sigDPModel)
Constructor.
Int_t binNumber(Double_t xCoord, Double_t yCoord) const
Find the global bin number for the given co-ordinates.
Definition: LauScfMap.cc:158
LauBkgndYieldList bkgndEvents_
Background yield(s)
Bool_t fixed() const
Check whether the parameter is fixed or floated.
virtual LauSPlot::NumbMap fixdSpeciesNames() const
Returns the names and yields of species that are fixed in the fit.
const LauDaughters * getDaughters() const
Retrieve the daughters.
void setupGenNtupleBranches()
Setup the required ntuple branches.
Bool_t writeLatexTable() const
Determine whether writing out of the latex table is enabled.
Double_t getm13() const
Get the m13 invariant mass.
Double_t getEvtm23Sq() const
Retrieve the invariant mass squared of the second and third daughters in the current event.
void setExtraNtupleVars()
Set-up other parameters that are derived from the fit results, e.g. fit fractions.
Pure abstract base class for defining a parameter containing an R value.
Definition: LauAbsRValue.hh:45
Bool_t getReweightedEvent(LauIsobarDynamics *dynamics)
Retrieve an event from the data sample, applying an accept/reject based on the given DP model.
const LauParameterList & extraPars() const
Const access the extra variables.
Double_t sigTotalLike_
Total signal likelihood.
virtual ~LauSimpleFitModel()
Destructor.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
Class for representing the 4D smearing matrix for mis-reconstructed signal (self cross feed)
Definition: LauScfMap.hh:45
virtual Bool_t genExpt()
Toy MC generation and fitting overloaded functions.
File containing declaration of LauAbsBkgndDPModel class.
virtual void addSPlotNtupleDoubleBranch(const TString &name)
Add a branch to the sPlot tree for storing a double.
void setFitNEvents()
Set the initial yields.
void updateSigEvents()
Update the signal events after Minuit sets background parameters.
void initialise(const std::vector< LauComplex > &coeffs)
Initialise the Dalitz plot dynamics.
virtual void getEvtExtraLikelihoods(UInt_t iEvt)
Determine the signal and background likelihood for the extra variables for a given event.
void addIntegerBranch(const TString &name)
Add integer branch to tree.
Definition: LauGenNtuple.cc:88
UInt_t nBkgndClasses() const
Returns the number of background classes.
const std::vector< Int_t > * trueBins(Int_t recoBin) const
Find which true bins contribute to the given reco bin.
Definition: LauScfMap.cc:169
UInt_t nEvents() const
Retrieve the number of events.
LauKinematics * kinematics_
The Dalitz plot kinematics object.
const LauParameterPList & fitPars() const
Const access to the fit variables.
void clearFitParVectors()
Clear the vectors containing fit parameters.
std::pair< LauGenInfo, Bool_t > eventsToGenerate()
Determine the number of events to generate for each hypothesis.
LauPdfPList scfPdfs_
The SCF PDFs.
LauBkgndEmbDataList bkgndTree_
The background event tree.
const LauComplex & getEvtDPAmp() const
Retrieve the total amplitude for the current event.
Double_t getEvtScfFraction() const
Retrieve the fraction of events that are poorly reconstructed (the self cross feed fraction) for the ...
void appendBinCentres(LauFitDataTree *inputData)
Append fake data points to the inputData for each bin in the SCF smearing matrix.
virtual void setupBkgndVectors()
Define the length of the background vectors.
void storeParsAndErrors(const std::vector< LauParameter * > &fitVars, const std::set< TString > &constrainedVars, const std::vector< LauParameter > &extraVars)
Store parameters and their errors.
virtual void checkInitFitParams()
Check the initial fit parameters.
virtual void setupSPlotNtupleBranches()
Add branches to store experiment number and the event number within the experiment.
void setSignalDPParameters()
Set the fit parameters for the DP model.
virtual void setNBkgndEvents(LauAbsRValue *nBkgndEvents)
Set the background event yield(s)
std::map< TString, Double_t > LauFitData
Type for holding event data.
virtual UInt_t index() const
Retrieve the index number of the coefficient set.
void modifyDataTree()
Recache the amplitude values for those that have changed.
void embedSignal(const TString &fileName, const TString &treeName, Bool_t reuseEventsWithinEnsemble, Bool_t reuseEventsWithinExperiment=kFALSE, Bool_t useReweighting=kFALSE)
Embed full simulation events for the signal, rather than generating toy from the PDFs.
void setBkgndDPModel(const TString &bkgndClass, LauAbsBkgndDPModel *bkgndModel)
Set the background DP models.
void updateSqDPKinematics(const Double_t mPrime, const Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m' and theta'.
Bool_t reuseSignal_
Boolean to reuse signal events.
void addSPlotNtupleBranches(const LauPdfPList *extraPdfs, const TString &prefix)
Add sPlot branches for the extra PDFs.
virtual void addGenNtupleIntegerBranch(const TString &name)
Add a branch to the gen tree for storing an integer.
virtual Bool_t fixed() const =0
Check is the parameter is fixed or floated.
LauBkgndDPModelList bkgndDPModels_
The background Dalitz Plot model.
Double_t range() const
The range allowed for the parameter.
void updatePull()
Call to update the bias and pull values.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
File containing declaration of LauPrint class and LauOutputLevel enum.
Class to define various output print commands.
Definition: LauPrint.hh:54
void getEmbeddedEvent(LauKinematics *kinematics)
Retrieve an event from the data sample.
LauScfMap * scfMap_
The smearing matrix for the SCF DP PDF.
Bool_t useDP() const
Is the Dalitz plot term in the likelihood.
const TMatrixD & covarianceMatrix() const
Access the fit covariance matrix.
LauFitData getValues(const std::vector< TString > &names) const
Get values of specified branches.
const LauParameter & getDPRate() const
Retrieve the overall Dalitz plot rate.
virtual void setSPlotNtupleIntegerBranchValue(const TString &name, Int_t value)
Set the value of an integer branch in the sPlot tree.
void clearExtraVarVectors()
Clear the vectors containing extra ntuple variables.
const LauParameter & getMeanEff() const
Retrieve the mean efficiency across the Dalitz plot.
Double_t calcSqDPJacobian(const Double_t mPrime, const Double_t thPrime) const
Calculate the Jacobian for the transformation m23^2, m13^2 -> m', theta' (square DP) at the given poi...
Double_t genValue() const
The value generated for the parameter.
UInt_t nExpt() const
Obtain the number of experiments.
void printFitParameters(const LauPdfPList &pdfList, std::ostream &fout) const
Print the fit parameters for all PDFs in the list.
Pure abstract base class for defining the efficiency description across the signal Dalitz plot.
Bool_t generateSignalEvent()
Generate signal event.
File containing declaration of LauSimpleFitModel class.
UInt_t eventsPerExpt() const
Obtain the total number of events in the current experiment.
void printFormat(std::ostream &stream, Double_t value) const
Method to choose the printing format to a specified level of precision.
Definition: LauPrint.cc:43
const LauAbsEffModel * getEffModel() const
Retrieve the model for the efficiency across the Dalitz plot.
Bool_t findBranches()
Find and read the branches in data tree.
void updateCoeffs(const std::vector< LauComplex > &coeffs)
Update the complex coefficients for the resonances.
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:81
LauParameter meanEff_
The mean efficiency.
virtual Double_t genValue() const =0
The value generated for the parameter.
Class to store the data for embedding in toy experiments.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
void updateFitNtuple()
Update the fit ntuple.
Double_t getm23Sq() const
Get the m23 invariant mass square.
Double_t getEvtthPrime() const
Retrieve the square Dalitz plot coordinate, theta', for the current event.
Bool_t usingBkgnd_
Background boolean.
LauComplex getFullAmplitude(const Int_t resID) const
Retrieve the Amplitude of resonance resID.
Class for defining signal dynamics using the isobar model.
virtual void fillGenNtupleBranches()
Fill the gen tuple branches.
File containing declaration of LauEmbeddedData class.
Double_t getm13Sq() const
Get the m13 invariant mass square.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:84
void updateFitParameters(LauPdfPList &pdfList)
Update the fit parameters for all PDFs in the list.
virtual Bool_t isLValue() const =0
Is the parameter also an L value or not.
UInt_t nSigDPPar_
Number of signal DP parameters.
void calcExtraInfo(const Bool_t init=kFALSE)
Calculate the fit fractions, mean efficiency and total DP rate.
const LauKinematics * getKinematics() const
Retrieve the Dalitz plot kinematics.
virtual const TString & name() const =0
Return the name of the parameter.
Double_t initValue() const
The initial value of the parameter.
UInt_t nEvents() const
Retrieve the number of events.
Abstract interface to the fitting and toy MC model.
File containing declaration of LauIsobarDynamics class.
Bool_t validBkgndClass(const TString &className) const
Check if the given background class is in the list.
LauParArray fitFrac_
Fit fractions.
LauPdfPList signalPdfs_
The signal PDFs.
Double_t getEvtm13Sq() const
Retrieve the invariant mass squared of the first and third daughters in the current event.
UInt_t nExtraPdfPar_
Number of extra PDF parameters.
Double_t prob(Int_t recoBin, Int_t trueBin) const
Probability of a true event in the given true bin migrating to the reco bin.
Definition: LauScfMap.cc:179
virtual void setAmpCoeffSet(LauAbsCoeffSet *coeffSet)
Set the DP amplitude coefficients.
virtual void writeOutTable(const TString &outputFile)
Write the fit results in latex table format.
File containing declaration of LauKinematics class.
Class to store the results from the toy MC generation into an ntuple.
Definition: LauGenNtuple.hh:45
std::vector< LauAbsPdf * > LauPdfPList
List of Pdfs.
virtual std::map< TString, LauComplex > getDPAmps(const Double_t m13Sq, const Double_t m23Sq)
Calculate the DP amplitude(s) for a given DP position.
const LauParArray & getFitFractions() const
Retrieve the fit fractions for the amplitude components.
std::vector< Double_t > recoJacobians_
The cached values of the sqDP jacobians for each event.
LauParameter * signalEvents_
Signal yield.