laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauBkgndDPModel.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 "LauBkgndDPModel.hh"
30 
31 #include "Lau2DHistDPPdf.hh"
32 #include "Lau2DSplineDPPdf.hh"
33 #include "LauDaughters.hh"
34 #include "LauFitDataTree.hh"
35 #include "LauKinematics.hh"
36 #include "LauRandom.hh"
37 #include "LauVetoes.hh"
38 
39 #include "TRandom.h"
40 #include "TSystem.h"
41 
42 #include <cstdlib>
43 #include <iostream>
44 
46  LauAbsBkgndDPModel( daughters, vetoes ),
47  symmetricalDP_( kFALSE ),
48  squareDP_( kFALSE ),
49  bgHistDPPdf_( 0 ),
50  curEvtHistVal_( 0.0 ),
51  maxPdfHeight_( 1.0 ),
52  pdfNorm_( 1.0 ),
53  doneGenWarning_( kFALSE ),
54  lowBinWarningIssued_( kFALSE )
55 {
56  if ( daughters != 0 ) {
57  symmetricalDP_ = daughters->gotSymmetricalDP();
58  }
59 }
60 
62 {
63  if ( bgHistDPPdf_ != 0 ) {
64  delete bgHistDPPdf_;
65  bgHistDPPdf_ = 0;
66  }
67 }
68 
69 void LauBkgndDPModel::setBkgndHisto( const TH2* histo,
70  Bool_t useInterpolation,
71  Bool_t fluctuateBins,
72  Bool_t useUpperHalfOnly,
73  Bool_t squareDP )
74 {
75  if ( ! histo ) {
76  std::cerr << "WARNING in LauBkgndDPModel::setBkgndHisto : Supplied background histogram pointer is null, likelihood for this component will be flat in the Dalitz plot"
77  << std::endl;
78  }
79 
80  Bool_t upperHalf = kFALSE;
81  if ( symmetricalDP_ == kTRUE && useUpperHalfOnly == kTRUE ) {
82  upperHalf = kTRUE;
83  }
84  std::cout << "INFO in LauBkgndDPModel::setBkgndHisto : Background histogram has upperHalf = "
85  << static_cast<Int_t>( upperHalf ) << std::endl;
86 
87  squareDP_ = squareDP;
88 
89  LauKinematics* kinematics = this->getKinematics();
90  const LauVetoes* vetoes = this->getVetoes();
91  bgHistDPPdf_ = new Lau2DHistDPPdf( histo,
92  kinematics,
93  vetoes,
94  useInterpolation,
95  fluctuateBins,
96  upperHalf,
97  squareDP );
98 
101 }
102 
103 void LauBkgndDPModel::setBkgndSpline( const TH2* histo,
104  Bool_t fluctuateBins,
105  Bool_t useUpperHalfOnly,
106  Bool_t squareDP )
107 {
108  if ( ! histo ) {
109  std::cerr << "WARNING in LauBkgndDPModel::setBkgndSpline : Supplied background histogram pointer is null, construction of spline will fail"
110  << std::endl;
111  }
112 
113  Bool_t upperHalf = kFALSE;
114  if ( symmetricalDP_ == kTRUE && useUpperHalfOnly == kTRUE ) {
115  upperHalf = kTRUE;
116  }
117  std::cout << "INFO in LauBkgndDPModel::setBkgndSpline : Background histogram has upperHalf = "
118  << static_cast<Int_t>( upperHalf ) << std::endl;
119 
120  squareDP_ = squareDP;
121 
122  LauKinematics* kinematics = this->getKinematics();
123  const LauVetoes* vetoes = this->getVetoes();
124  bgHistDPPdf_ =
125  new Lau2DSplineDPPdf( histo, kinematics, vetoes, fluctuateBins, upperHalf, squareDP );
126 
129 }
130 
131 Double_t LauBkgndDPModel::calcHistValue( Double_t xVal, Double_t yVal )
132 {
133  // Get the likelihood value of the background in the Dalitz plot.
134 
135  // Check that we have a valid histogram PDF
136  if ( bgHistDPPdf_ == 0 ) {
137  std::cerr << "WARNING in LauBkgndDPModel::calcHistValue : We don't have a histogram so assuming the likelihood is flat in the Dalitz plot."
138  << std::endl;
139  this->setBkgndHisto( 0, kFALSE, kFALSE, kFALSE, kFALSE );
140  }
141 
142  // Find out the un-normalised PDF value
143  Double_t value = bgHistDPPdf_->interpolateXY( xVal, yVal );
144 
145  // Check that the value is greater than zero
146  // If we're using a spline then negative values can be caused by adjacent bins that all contain a value of zero.
147  // The spline requires the value, its first derivatives and the mixed second derivative to be continuous and to match the input histogram
148  // at the bin centres. Derivatives are calculated using a finite difference approximation taking the difference between the neighbouring bins.
149  // If two bins are zero but the third is not then the second bin will have a positive first derivative causing the spline to dip below zero
150  // between the two zero bins to remain smooth. Such dips are unavoidable but are correctly removed here.
151  if ( value < 0.0 ) {
152  if ( ! lowBinWarningIssued_ ) {
153  std::cerr << "WARNING in LauBkgndDPModel::calcHistValue : Value " << value
154  << " is less than 0 - setting to 0. You may want to check your histogram!"
155  << std::endl
156  << " : If you are using a spline then this could be caused by adjacent empty bins. Further warnings will be suppressed."
157  << std::endl;
158  lowBinWarningIssued_ = kTRUE;
159  }
160  return 0.0;
161  }
162 
163  LauKinematics* kinematics = this->getKinematics();
164 
165  // For square DP co-ordinates, need to divide by Jacobian
166  if ( squareDP_ == kTRUE ) {
167  // Make sure that square DP kinematics are correct, then get Jacobian
168  kinematics->updateSqDPKinematics( xVal, yVal );
169  Double_t jacobian = kinematics->calcSqDPJacobian();
170  value /= jacobian;
171  }
172 
173  return value;
174 }
175 
177 {
178 }
179 
181 {
182  // Routine to generate the background, using data provided by an
183  // already defined histogram.
184 
185  LauKinematics* kinematics = this->getKinematics();
186 
187  Bool_t gotBG( kFALSE );
188 
189  while ( gotBG == kFALSE ) {
190 
191  if ( squareDP_ == kTRUE ) {
192  // Generate a point in m', theta' space. By construction, this point
193  // is already within the DP region.
194  Double_t mPrime( 0.0 ), thetaPrime( 0.0 );
195  kinematics->genFlatSqDP( mPrime, thetaPrime );
196 
197  // If we're in a symmetrical DP then we should only generate events in one half
198  if ( symmetricalDP_ && thetaPrime > 0.5 ) {
199  thetaPrime = 1.0 - thetaPrime;
200  }
201 
202  // Calculate histogram height for DP point and
203  // compare with the maximum height
204  if ( bgHistDPPdf_ != 0 ) {
205  Double_t bgContDP = bgHistDPPdf_->interpolateXY( mPrime, thetaPrime ) / maxPdfHeight_;
206  if ( LauRandom::randomFun()->Rndm() < bgContDP ) {
207  kinematics->updateSqDPKinematics( mPrime, thetaPrime );
208  gotBG = kTRUE;
209  }
210  } else {
211  if ( ! doneGenWarning_ ) {
212  std::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!"
213  << std::endl;
214  std::cerr << "WARNING in LauBkgndDPModel::generate : This should never happen!! What have you done?!"
215  << std::endl;
216  doneGenWarning_ = kTRUE;
217  }
218  kinematics->updateSqDPKinematics( mPrime, thetaPrime );
219  gotBG = kTRUE;
220  }
221  } else {
222  // Generate a point in the Dalitz plot (phase-space).
223  Double_t m13Sq( 0.0 ), m23Sq( 0.0 );
224  kinematics->genFlatPhaseSpace( m13Sq, m23Sq );
225 
226  // If we're in a symmetrical DP then we should only generate events in one half
227  if ( symmetricalDP_ && m13Sq > m23Sq ) {
228  Double_t tmpSq = m13Sq;
229  m13Sq = m23Sq;
230  m23Sq = tmpSq;
231  }
232 
233  // Calculate histogram height for DP point and
234  // compare with the maximum height
235  if ( bgHistDPPdf_ != 0 ) {
236  Double_t bgContDP = bgHistDPPdf_->interpolateXY( m13Sq, m23Sq ) / maxPdfHeight_;
237  if ( LauRandom::randomFun()->Rndm() < bgContDP ) {
238  kinematics->updateKinematics( m13Sq, m23Sq );
239  gotBG = kTRUE;
240  }
241  } else {
242  if ( ! doneGenWarning_ ) {
243  std::cerr << "WARNING in LauBkgndDPModel::generate : We don't have a histogram so generating flat in the DP."
244  << std::endl;
245  doneGenWarning_ = kTRUE;
246  }
247  kinematics->updateKinematics( m13Sq, m23Sq );
248  gotBG = kTRUE;
249  }
250  }
251  }
252 
253  // Implement veto
254  Bool_t vetoOK( kTRUE );
255  const LauVetoes* vetoes = this->getVetoes();
256  if ( vetoes ) {
257  vetoOK = vetoes->passVeto( kinematics );
258  }
259  // Call this function recusively until we pass the veto.
260  if ( vetoOK == kFALSE ) {
261  this->generate();
262  }
263  return kTRUE;
264 }
265 
267 {
268  // In LauFitDataTree, the first two variables should always be
269  // m13^2 and m23^2. Other variables follow thus: charge.
270 
271  Int_t nBranches = inputFitTree.nBranches();
272  if ( nBranches < 2 ) {
273  std::cerr << "ERROR in LauBkgndDPModel::fillDataTree : Expecting at least 2 variables "
274  << "in input data tree, but there are " << nBranches << "!\n";
275  std::cerr << " : Make sure you have the right number of variables in your input data file!"
276  << std::endl;
277  return;
278  }
279 
280  const UInt_t nEvents { inputFitTree.nEvents() };
281 
282  // clear and resize the data vector
283  bgData_.clear();
284  bgData_.resize( nEvents );
285 
286  Double_t m13Sq { 0.0 }, m23Sq { 0.0 };
287 
288  for ( UInt_t iEvt = 0; iEvt < nEvents; iEvt++ ) {
289 
290  const LauFitData& dataValues = inputFitTree.getData( iEvt );
291  m13Sq = dataValues.at( "m13Sq" );
292  m23Sq = dataValues.at( "m23Sq" );
293 
294  bgData_[iEvt] = this->getUnNormValue( m13Sq, m23Sq );
295  }
296 }
297 
298 Double_t LauBkgndDPModel::getUnNormValue( const Double_t m13Sq, const Double_t m23Sq )
299 {
300  // Update the kinematics. This will also update m' and theta' if squareDP = true
301  LauKinematics* kinematics = this->getKinematics();
302  kinematics->updateKinematics( m13Sq, m23Sq );
303 
304  if ( squareDP_ ) {
305  const Double_t mPrime { kinematics->getmPrime() };
306  const Double_t thetaPrime { kinematics->getThetaPrime() };
307  curEvtHistVal_ = this->calcHistValue( mPrime, thetaPrime );
308  } else {
309  curEvtHistVal_ = this->calcHistValue( m13Sq, m23Sq );
310  }
311 
312  return curEvtHistVal_;
313 }
314 
315 Double_t LauBkgndDPModel::getLikelihood( const Double_t m13Sq, const Double_t m23Sq )
316 {
317  this->getUnNormValue( m13Sq, m23Sq );
318  return curEvtHistVal_ / this->getPdfNorm();
319 }
320 
321 Double_t LauBkgndDPModel::getUnNormValue( UInt_t iEvt )
322 {
323  // Retrieve the likelihood for the given event
324  this->setDataEventNo( iEvt );
325  return curEvtHistVal_;
326 }
327 
328 Double_t LauBkgndDPModel::getLikelihood( UInt_t iEvt )
329 {
330  // Retrieve the likelihood for the given event
331  this->setDataEventNo( iEvt );
332  Double_t llhd = curEvtHistVal_ / this->getPdfNorm();
333  return llhd;
334 }
335 
337 {
338  // Retrieve the data for event iEvt
339  if ( bgData_.size() > iEvt ) {
340  curEvtHistVal_ = bgData_[iEvt];
341  } else {
342  std::cerr << "ERROR in LauBkgndDPModel::setDataEventNo : Event index too large: " << iEvt
343  << " >= " << bgData_.size() << "." << std::endl;
344  }
345 }
Double_t pdfNorm_
Normalisation of PDF.
File containing LauRandom namespace.
Bool_t squareDP_
Flags whether or not to work in square DP coordinates.
virtual Double_t getUnNormValue(const Double_t m13Sq, const Double_t m23Sq)
Get unnormalised likelihood for a given DP position.
Double_t value() const
The value of the parameter.
UInt_t nBranches() const
Obtain the number of branches in the tree.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
Double_t getThetaPrime() const
Get theta' value.
virtual ~LauBkgndDPModel()
Destructor.
The abstract interface for a background Dalitz plot model.
Double_t getmPrime() const
Get m' value.
void setBkgndSpline(const TH2 *histo, Bool_t fluctuateBins, Bool_t useUpperHalfOnly, Bool_t squareDP)
Set the background histogram and generate a spline.
Double_t calcHistValue(Double_t xVal, Double_t yVal)
Calulate histogram value at a given point.
Bool_t passVeto(const LauKinematics *kinematics) const
Check whether the specified Dalitz plot point passes the vetoes.
Definition: LauVetoes.cc:124
void genFlatSqDP(Double_t &mPrime, Double_t &thetaPrime) const
Routine to generate events flat in the square Dalitz plot.
LauBkgndDPModel(LauDaughters *daughters, LauVetoes *vetoes)
Constructor.
void updateKinematics(const Double_t m13Sq, const Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
File containing declaration of LauFitDataTree class.
void setBkgndHisto(const TH2 *histo, Bool_t useInterpolation, Bool_t fluctuateBins, Bool_t useUpperHalfOnly, Bool_t squareDP=kFALSE)
Set background histogram.
Class to store the input fit variables.
Class for defining a 2D DP spline PDF.
Bool_t lowBinWarningIssued_
Flag to track whether a warning has been issued for bin values less than zero.
Lau2DAbsDPPdf * bgHistDPPdf_
PDF of Dalitz plot background, from a 2D histogram.
virtual Double_t getHistNorm() const =0
Retrieve PDF normalisation.
Bool_t doneGenWarning_
Boolean to indicate if the warning that there is no histogram has already been issued.
File containing declaration of Lau2DHistDPPdf class.
virtual Double_t getMaxHeight() const =0
Retrieve maximum height.
virtual void fillDataTree(const LauFitDataTree &fitDataTree)
Cache the input data and (if appropriate) the per-event likelihood values.
File containing declaration of Lau2DSplineDPPdf class.
File containing declaration of LauDaughters class.
const LauVetoes * getVetoes() const
Get vetoes in the Dalitz plot.
virtual Double_t getLikelihood(const Double_t m13Sq, const Double_t m23Sq)
Get likelihood for a given DP position.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:33
const LauKinematics * getKinematics() const
Get the Dalitz plot kinematics.
virtual Double_t interpolateXY(Double_t x, Double_t y) const =0
Perform the interpolation (unnormalised)
void genFlatPhaseSpace(Double_t &m13Sq, Double_t &m23Sq) const
Routine to generate events flat in phase-space.
virtual Double_t getPdfNorm() const
Get PDF normalisation constant.
std::map< TString, Double_t > LauFitData
Type for holding event data.
void updateSqDPKinematics(const Double_t mPrime, const Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m' and theta'.
virtual void setDataEventNo(UInt_t iEvt)
Set data event number.
virtual Bool_t generate()
Generate a toy MC event from the model.
Class for defining vetoes within the Dalitz plot.
Definition: LauVetoes.hh:49
Double_t maxPdfHeight_
Maximum height of PDF.
Class for calculating 3-body kinematic quantities.
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...
Class for defining a 2D DP histogram PDF.
std::vector< Double_t > bgData_
Cached histogram values for each event.
File containing declaration of LauBkgndDPModel class.
Double_t curEvtHistVal_
Histogram value for the current event.
File containing declaration of LauVetoes class.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:84
virtual void initialise()
Initialise the model.
UInt_t nEvents() const
Retrieve the number of events.
Bool_t symmetricalDP_
Flags whether the DP is symmetrical or not.
File containing declaration of LauKinematics class.