laura is hosted by Hepforge, IPPP Durham
Laura++  v1r0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauBkgndDPModel.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <cstdlib>
16 #include <iostream>
17 using std::cout;
18 using std::cerr;
19 using std::endl;
20 
21 #include "TRandom.h"
22 #include "TSystem.h"
23 
24 #include "Lau2DHistDPPdf.hh"
25 #include "LauBkgndDPModel.hh"
26 #include "LauDaughters.hh"
27 #include "LauFitDataTree.hh"
28 #include "LauKinematics.hh"
29 #include "LauRandom.hh"
30 #include "LauVetoes.hh"
31 
32 ClassImp(LauBkgndDPModel)
33 
34 
36  LauAbsBkgndDPModel(daughters, vetoes),
37  symmetricalDP_(kFALSE),
38  squareDP_(kFALSE),
39  bgHistDPPdf_(0),
40  curEvtHistVal_(0.0),
41  maxPdfHeight_(1.0),
42  pdfNorm_(1.0),
43  doneGenWarning_(kFALSE)
44 {
45  if (daughters != 0) {
46  symmetricalDP_ = daughters->gotSymmetricalDP();
47  }
48 }
49 
51 {
52  if (bgHistDPPdf_ != 0) {
53  delete bgHistDPPdf_; bgHistDPPdf_ = 0;
54  }
55 }
56 
57 void LauBkgndDPModel::setBkgndHisto(const TH2* histo, Bool_t useInterpolation, Bool_t fluctuateBins, Bool_t useUpperHalfOnly, Bool_t squareDP)
58 {
59  Bool_t upperHalf = kFALSE;
60  if (symmetricalDP_ == kTRUE && useUpperHalfOnly == kTRUE) {upperHalf = kTRUE;}
61  cout<<"Bg histogram has upperHalf = "<<static_cast<Int_t>(upperHalf)<<endl;
62 
63  squareDP_ = squareDP;
64 
65  LauKinematics* kinematics = this->getKinematics();
66  const LauVetoes* vetoes = this->getVetoes();
67  bgHistDPPdf_ = new Lau2DHistDPPdf(histo, kinematics, vetoes, useInterpolation, fluctuateBins, upperHalf, squareDP);
68 
71 }
72 
73 Double_t LauBkgndDPModel::calcHistValue(Double_t xVal, Double_t yVal)
74 {
75  // Get the likelihood value of the background in the Dalitz plot.
76 
77  // Check that we have a valid histogram PDF
78  if (bgHistDPPdf_ == 0) {
79  cerr << "WARNING in LauBkgndDPModel::calcHistValue : We don't have a histogram so assuming the likelihood is flat in the Dalitz plot." << endl;
80  this->setBkgndHisto( 0, kFALSE, kFALSE, kFALSE, kFALSE );
81  }
82 
83  // Find out the un-normalised PDF value
84  Double_t value = bgHistDPPdf_->interpolateXY(xVal, yVal);
85 
86  LauKinematics* kinematics = this->getKinematics();
87 
88  // For square DP co-ordinates, need to divide by Jacobian
89  if (squareDP_ == kTRUE) {
90  // Make sure that square DP kinematics are correct, then get Jacobian
91  kinematics->updateSqDPKinematics(xVal, yVal);
92  Double_t jacobian = kinematics->calcSqDPJacobian();
93  value /= jacobian;
94  }
95 
96  return value;
97 }
98 
100 {
101 }
102 
104 {
105  // Routine to generate the background, using data provided by an
106  // already defined histogram.
107 
108  LauKinematics* kinematics = this->getKinematics();
109 
110  Bool_t gotBG(kFALSE);
111 
112  while (gotBG == kFALSE) {
113 
114  if (squareDP_ == kTRUE) {
115  // Generate a point in m', theta' space. By construction, this point
116  // is already within the DP region.
117  Double_t mPrime(0.0), thetaPrime(0.0);
118  kinematics->genFlatSqDP(mPrime, thetaPrime);
119 
120  // If we're in a symmetrical DP then we should only generate events in one half
121  if ( symmetricalDP_ && thetaPrime > 0.5 ) {
122  thetaPrime = 1.0 - thetaPrime;
123  }
124 
125  // Calculate histogram height for DP point and
126  // compare with the maximum height
127  if ( bgHistDPPdf_ != 0 ) {
128  Double_t bgContDP = bgHistDPPdf_->interpolateXY(mPrime, thetaPrime)/maxPdfHeight_;
129  if (LauRandom::randomFun()->Rndm() < bgContDP) {
130  kinematics->updateSqDPKinematics(mPrime, thetaPrime);
131  gotBG = kTRUE;
132  }
133  } else {
134  if ( !doneGenWarning_ ) {
135  cerr << "WARNING in LauBkgndDPModel::generate : We don't have a histogram so generating flat in the square DP, which won't be flat in the conventional DP!" << endl;
136  cerr << "WARNING in LauBkgndDPModel::generate : This should never happen!! What have you done?!" << endl;
137  doneGenWarning_ = kTRUE;
138  }
139  kinematics->updateSqDPKinematics(mPrime, thetaPrime);
140  gotBG = kTRUE;
141  }
142  } else {
143  // Generate a point in the Dalitz plot (phase-space).
144  Double_t m13Sq(0.0), m23Sq(0.0);
145  kinematics->genFlatPhaseSpace(m13Sq, m23Sq);
146 
147  // If we're in a symmetrical DP then we should only generate events in one half
148  if ( symmetricalDP_ && m13Sq > m23Sq ) {
149  Double_t tmpSq = m13Sq;
150  m13Sq = m23Sq;
151  m23Sq = tmpSq;
152  }
153 
154  // Calculate histogram height for DP point and
155  // compare with the maximum height
156  if ( bgHistDPPdf_ != 0 ) {
157  Double_t bgContDP = bgHistDPPdf_->interpolateXY(m13Sq, m23Sq)/maxPdfHeight_;
158  if (LauRandom::randomFun()->Rndm() < bgContDP) {
159  kinematics->updateKinematics(m13Sq, m23Sq);
160  gotBG = kTRUE;
161  }
162  } else {
163  if ( !doneGenWarning_ ) {
164  cerr << "WARNING in LauBkgndDPModel::generate : We don't have a histogram so generating flat in the DP." << endl;
165  doneGenWarning_ = kTRUE;
166  }
167  kinematics->updateKinematics(m13Sq, m23Sq);
168  gotBG = kTRUE;
169  }
170  }
171  }
172 
173  // Implement veto
174  Bool_t vetoOK(kTRUE);
175  const LauVetoes* vetoes = this->getVetoes();
176  if (vetoes) {
177  vetoOK = vetoes->passVeto(kinematics);
178  }
179  // Call this function recusively until we pass the veto.
180  if (vetoOK == kFALSE) {
181  this->generate();
182  }
183  return kTRUE;
184 }
185 
187 {
188  // In LauFitDataTree, the first two variables should always be
189  // m13^2 and m23^2. Other variables follow thus: charge.
190 
191  Int_t nBranches = inputFitTree.nBranches();
192  if (nBranches < 2) {
193  cout<<"Error in LauBkgndDPModel::fillDataTree. Expecting at least 2 variables "
194  <<"in input data tree, but there are "<<nBranches<<"! Make sure you have "
195  <<"the right number of variables in your input data file!"<<endl;
196  return;
197  }
198 
199  Double_t m13Sq(0.0), m23Sq(0.0), mPrime(0.0), thetaPrime(0.0);
200 
201  UInt_t nEvents = inputFitTree.nEvents();
202  LauKinematics* kinematics = this->getKinematics();
203 
204  // clear and resize the data vector
205  bgData_.clear();
206  bgData_.resize(nEvents);
207 
208  for (UInt_t iEvt = 0; iEvt < nEvents; iEvt++) {
209 
210  const LauFitData& dataValues = inputFitTree.getData(iEvt);
211  LauFitData::const_iterator iter = dataValues.find("m13Sq");
212  m13Sq = iter->second;
213  iter = dataValues.find("m23Sq");
214  m23Sq = iter->second;
215 
216  // Update the kinematics. This will also update m' and theta' if squareDP = true
217  kinematics->updateKinematics(m13Sq, m23Sq);
218 
219  if (squareDP_ == kTRUE) {
220  mPrime = kinematics->getmPrime();
221  thetaPrime = kinematics->getThetaPrime();
222  curEvtHistVal_ = this->calcHistValue(mPrime, thetaPrime);
223  } else {
224  curEvtHistVal_ = this->calcHistValue(m13Sq, m23Sq);
225  }
226 
227  bgData_[iEvt] = curEvtHistVal_;
228  }
229 }
230 
231 Double_t LauBkgndDPModel::getUnNormValue(UInt_t iEvt)
232 {
233  // Retrieve the likelihood for the given event
234  this->setDataEventNo(iEvt);
235  return curEvtHistVal_;
236 }
237 
238 Double_t LauBkgndDPModel::getLikelihood(UInt_t iEvt)
239 {
240  // Retrieve the likelihood for the given event
241  this->setDataEventNo(iEvt);
242  Double_t llhd = curEvtHistVal_ / this->getPdfNorm();
243  return llhd;
244 }
245 
247 {
248  // Retrieve the data for event iEvt
249  if (bgData_.size() > iEvt) {
250  curEvtHistVal_ = bgData_[iEvt];
251  } else {
252  cerr<<"ERROR in LauBkgndDPModel::setDataEventNo : Event index too large: "<<iEvt<<" >= "<<bgData_.size()<<"."<<endl;
253  }
254 }
255 
The abstract interface for a background Dalitz plot model.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
File containing declaration of LauFitDataTree class.
virtual void initialise()
Initialise the model.
Double_t maxPdfHeight_
Maximum height of PDF.
void updateSqDPKinematics(Double_t mPrime, Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m&#39; and theta&#39;.
File containing declaration of LauBkgndDPModel class.
std::vector< Double_t > bgData_
Cached histogram values for each event.
Double_t getmPrime() const
Get m&#39; value.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:33
File containing declaration of LauDaughters class.
virtual Bool_t generate()
Generate a toy MC event from the model.
Class for defining a 2D DP histogram PDF.
const LauKinematics * getKinematics() const
Get the Dalitz plot kinematics.
Double_t getHistNorm() const
Retrieve PDF normalisation.
std::map< TString, Double_t > LauFitData
Type for holding event data.
UInt_t nBranches() const
Obtain the number of branches in the tree.
Bool_t symmetricalDP_
Flags whether the DP is symmetrical or not.
File containing declaration of Lau2DHistDPPdf class.
File containing declaration of LauKinematics class.
void updateKinematics(Double_t m13Sq, Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
Bool_t doneGenWarning_
Boolean to indicate if the warning that there is no histogram has already been issued.
void setBkgndHisto(const TH2 *histo, Bool_t useInterpolation, Bool_t fluctuateBins, Bool_t useUpperHalfOnly, Bool_t squareDP=kFALSE)
Set background histogram.
virtual Double_t getLikelihood(UInt_t iEvt)
Get likelihood for a given event.
void genFlatSqDP(Double_t &mPrime, Double_t &thetaPrime) const
Routine to generate events flat in the square Dalitz plot.
const LauVetoes * getVetoes() const
Get vetoes in the Dalitz plot.
Double_t calcHistValue(Double_t xVal, Double_t yVal)
Calulate histogram value at a given point.
Double_t curEvtHistVal_
Histogram value for the current event.
Bool_t squareDP_
Flags whether or not to work in square DP coordinates.
Class for defining a histogram-based background Dalitz plot model.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
File containing LauRandom namespace.
virtual Double_t getUnNormValue(UInt_t iEvt)
Get unnormalised likelihood for a given event.
virtual ~LauBkgndDPModel()
Destructor.
Bool_t passVeto(Double_t &m12Sq, Double_t &m23Sq, Double_t &m13Sq) const
Check whether the specified Dalitz plot point passes the vetoes.
Definition: LauVetoes.cc:109
Double_t getMaxHeight() const
Retrieve maximum height.
Double_t getThetaPrime() const
Get theta&#39; value.
virtual void fillDataTree(const LauFitDataTree &fitDataTree)
Cache the input data and (if appropriate) the per-event likelihood values.
Lau2DHistDPPdf * bgHistDPPdf_
PDF of Dalitz plot background, from a 2D histogram.
virtual void setDataEventNo(UInt_t iEvt)
Set data event number.
void genFlatPhaseSpace(Double_t &m13Sq, Double_t &m23Sq) const
Routine to generate events flat in phase-space.
Class for calculating 3-body kinematic quantities.
Double_t value() const
The value of the parameter.
Double_t pdfNorm_
Normalisation of PDF.
UInt_t nEvents() const
Retrieve the number of events.
File containing declaration of LauVetoes class.
Double_t calcSqDPJacobian()
Calculate the Jacobian for the transformation m23^2, m13^2 -&gt; m&#39;, theta&#39; (square DP) ...
Double_t interpolateXY(Double_t x, Double_t y) const
Perform the interpolation (unnormalised)
Class for defining vetoes within the Dalitz plot.
Definition: LauVetoes.hh:33
Class to store the input fit variables.
virtual Double_t getPdfNorm() const
Get PDF normalisation constant.