laura is hosted by Hepforge, IPPP Durham
Laura++  v3r3
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauDPDepMapPdf.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2009 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 #include <vector>
17 using std::cout;
18 using std::cerr;
19 using std::endl;
20 
21 #include "TH1.h"
22 #include "TMath.h"
23 #include "TSystem.h"
24 
25 #include "LauConstants.hh"
26 #include "LauDaughters.hh"
27 #include "Lau2DHistDP.hh"
28 #include "LauParameter.hh"
29 #include "LauDPDepMapPdf.hh"
30 
32 
33 
34 LauDPDepMapPdf::LauDPDepMapPdf(const std::vector<LauAbsPdf*>& pdfs,
35  const LauDaughters* daughters,
36  const TH2* dpHisto, Bool_t upperHalf) :
37  LauAbsPdf((!pdfs.empty() && pdfs[0]) ? pdfs[0]->varNames() : std::vector<TString>(), std::vector<LauAbsRValue*>(), (!pdfs.empty() && pdfs[0]) ? pdfs[0]->getMinAbscissas() : LauFitData(), (!pdfs.empty() && pdfs[0]) ? pdfs[0]->getMaxAbscissas() : LauFitData()),
38  daughters_( new LauDaughters(*daughters) ),
39  pdfs_(pdfs),
40  dpDependence_( new Lau2DHistDP(dpHisto, daughters, kFALSE, kFALSE, 0, 0, upperHalf, daughters->squareDP()) ),
41  dpAxisDependence_(0),
42  dpAxis_( CentreDist ),
43  indices_(0)
44 {
45  // Constructor for the PDF map.
46  // The index into the PDF collection is read from the DP histogram.
47 
48 
49  // So the first thing we have to do is check the pointers are all valid.
50  for ( std::vector<LauAbsPdf*>::const_iterator iter = pdfs_.begin(); iter != pdfs_.end(); ++iter){
51 
52  if (!(*iter)) {
53  cerr<<"ERROR in LauDPDepMapPdf constructor: one of the PDF pointers is null."<<endl;
54  gSystem->Exit(EXIT_FAILURE);
55  }
56 
57  // Next check that the abscissa ranges are the same for each PDF
58  if (pdfs_[0]->getMinAbscissa() != (*iter)->getMinAbscissa()) {
59  cerr<<"ERROR in LauDPDepMapPdf constructor: minimum abscissa values not the same for two PDFs."<<endl;
60  gSystem->Exit(EXIT_FAILURE);
61  }
62  if (pdfs_[0]->getMaxAbscissa() != (*iter)->getMaxAbscissa()) {
63  cerr<<"ERROR in LauDPDepMapPdf constructor: maximum abscissa values not the same for two PDFs."<<endl;
64  gSystem->Exit(EXIT_FAILURE);
65  }
66 
67  // Also check that both PDFs are expecting the same number of input variables
68  if (pdfs_[0]->nInputVars() != (*iter)->nInputVars()) {
69  cerr<<"ERROR in LauDPDepMapPdf constructor: number of input variables not the same for two PDFs."<<endl;
70  gSystem->Exit(EXIT_FAILURE);
71  }
72 
73  // Also check that both PDFs are expecting the same variable name(s)
74  if (pdfs_[0]->varNames() != (*iter)->varNames()) {
75  cerr<<"ERROR in LauDPDepMapPdf constructor: variable name(s) not the same for two PDFs."<<endl;
76  gSystem->Exit(EXIT_FAILURE);
77  }
78  }
79 
80  // Then we need to grab all the parameters and pass them to the base class.
81  // This is so that when we are asked for them they can be put into the fit.
82 
83 
84  for ( std::vector<LauAbsPdf*>::const_iterator iter = pdfs_.begin(); iter != pdfs_.end(); ++iter){
85  std::vector<LauAbsRValue*>& pdfpars = (*iter)->getParameters();
86  this->addParameters(pdfpars);
87 
88  }
89 
90  // Cache the normalisation factor
91  this->calcNorm();
92 }
93 
94 LauDPDepMapPdf::LauDPDepMapPdf(const std::vector<LauAbsPdf*>& pdfs,
95  const LauDaughters* daughters,
96  const TH1* dpAxisHisto,
97  DPAxis dpAxis) :
98  LauAbsPdf((!pdfs.empty() && pdfs[0]) ? pdfs[0]->varNames() : std::vector<TString>(), std::vector<LauAbsRValue*>(), (!pdfs.empty() && pdfs[0]) ? pdfs[0]->getMinAbscissas() : LauFitData(), (!pdfs.empty() && pdfs[0]) ? pdfs[0]->getMaxAbscissas() : LauFitData()),
99  daughters_( new LauDaughters(*daughters) ),
100  pdfs_(pdfs),
101  dpDependence_( 0 ),
102  dpAxisDependence_(dpAxisHisto ? dynamic_cast<TH1*>(dpAxisHisto->Clone()) : 0),
103  dpAxis_( dpAxis ),
104  indices_( 0 )
105 {
106  // Constructor for the PDFs map.
107  // The index into the PDF collection is read from one of
108  // the DP axes.
109 
110  if ( dpAxisDependence_ ) {
111  dpAxisDependence_->SetDirectory(0);
112  }
113 
114  // So the first thing we have to do is check the pointers are all valid.
115  for ( std::vector<LauAbsPdf*>::const_iterator iter = pdfs_.begin(); iter != pdfs_.end(); ++iter){
116 
117  if (!(*iter)) {
118  cerr<<"ERROR in LauDPDepMapPdf constructor: one of the PDF pointers is null."<<endl;
119  gSystem->Exit(EXIT_FAILURE);
120  }
121 
122  // Next check that the abscissa ranges are the same for each PDF
123  if (pdfs_[0]->getMinAbscissa() != (*iter)->getMinAbscissa()) {
124  cerr<<"ERROR in LauDPDepMapPdf constructor: minimum abscissa values not the same for two PDFs."<<endl;
125  gSystem->Exit(EXIT_FAILURE);
126  }
127  if (pdfs_[0]->getMaxAbscissa() != (*iter)->getMaxAbscissa()) {
128  cerr<<"ERROR in LauDPDepMapPdf constructor: maximum abscissa values not the same for two PDFs."<<endl;
129  gSystem->Exit(EXIT_FAILURE);
130  }
131 
132  // Also check that both PDFs are expecting the same number of input variables
133  if (pdfs_[0]->nInputVars() != (*iter)->nInputVars()) {
134  cerr<<"ERROR in LauDPDepMapPdf constructor: number of input variables not the same for two PDFs."<<endl;
135  gSystem->Exit(EXIT_FAILURE);
136  }
137 
138  // Also check that both PDFs are expecting the same variable name(s)
139  if (pdfs_[0]->varNames() != (*iter)->varNames()) {
140  cerr<<"ERROR in LauDPDepMapPdf constructor: variable name(s) not the same for two PDFs."<<endl;
141  gSystem->Exit(EXIT_FAILURE);
142  }
143  }
144 
145 
146  // Then we need to grab all the parameters and pass them to the base class.
147  // This is so that when we are asked for them they can be put into the fit.
148 
149  for ( std::vector<LauAbsPdf*>::const_iterator iter = pdfs_.begin(); iter != pdfs_.end(); ++iter){
150  std::vector<LauAbsRValue*>& pdfpars = (*iter)->getParameters();
151  this->addParameters(pdfpars);
152 
153  }
154 
155  // Cache the normalisation factor
156  this->calcNorm();
157 }
158 
160 {
161  // Destructor
162  delete daughters_; daughters_ = 0;
163  delete dpDependence_; dpDependence_ = 0;
165 }
166 
167 UInt_t LauDPDepMapPdf::determineDPRegion( Double_t m13Sq, Double_t m23Sq ) const
168 {
169  UInt_t regionIndex(0);
170 
171  LauKinematics* kinematics = daughters_->getKinematics();
172  kinematics->updateKinematics( m13Sq, m23Sq );
173 
174  if ( dpDependence_ ) {
175  if (daughters_->squareDP()){
176  Double_t mprime = kinematics->getmPrime();
177  Double_t thprime = kinematics->getThetaPrime();
178  regionIndex = static_cast<UInt_t>( dpDependence_->interpolateXY( mprime, thprime ) );
179  } else {
180  regionIndex = static_cast<UInt_t>( dpDependence_->interpolateXY( m13Sq, m23Sq ) );
181  }
182  } else if ( dpAxisDependence_ ) {
183  Double_t dpPos(0.0);
184  if ( dpAxis_ == M12 ) {
185  dpPos = kinematics->calcThirdMassSq(m13Sq,m23Sq);
186  } else if ( dpAxis_ == M13 ) {
187  dpPos = m13Sq;
188  } else if ( dpAxis_ == M23 ) {
189  dpPos = m23Sq;
190  } else if ( dpAxis_ == MMIN ) {
191  dpPos = TMath::Min( m13Sq, m23Sq );
192  } else if ( dpAxis_ == MMAX ) {
193  dpPos = TMath::Max( m13Sq, m23Sq );
194  } else {
195  dpPos = kinematics->distanceFromDPCentre(m13Sq,m23Sq);
196  }
197  Int_t bin = dpAxisDependence_->FindFixBin( dpPos );
198  regionIndex = static_cast<UInt_t>( dpAxisDependence_->GetBinContent( bin ) );
199  } else {
200  // This should never happen
201  cerr << "ERROR in LauDPDepMapPdf::determineDPRegion : No means of determining the region!" << endl;
202  gSystem->Exit(EXIT_FAILURE);
203  }
204 
205  return regionIndex;
206 }
207 
209 {
210  // Check that the given abscissa is within the allowed range
211  if (!this->checkRange(abscissas)) {
212  gSystem->Exit(EXIT_FAILURE);
213  }
214 
215  // Create a set of absissas that doesn't contain the DP variables
216  UInt_t nVars = this->nInputVars();
217  LauAbscissas noDPVars( nVars );
218  for ( UInt_t i(0); i < nVars; ++i ) {
219  noDPVars[i] = abscissas[i];
220  }
221 
222  // Find which region of the DP we're in
223  // The DP variables will be abscissas[nInputVars] and
224  // abscissas[nInputVars+1] (if present).
225  UInt_t regionIndex(0);
226  if ( abscissas.size() == nVars+2 ) {
227  Double_t m13Sq = abscissas[nVars];
228  Double_t m23Sq = abscissas[nVars+1];
229  regionIndex = this->determineDPRegion( m13Sq, m23Sq );
230  } else {
231  // This should never happen
232  cerr << "ERROR in LauDPDepMapPdf::calcLikelihoodInfo : DP vars not supplied, no means of determining the region!" << endl;
233  gSystem->Exit(EXIT_FAILURE);
234  }
235 
236  // Check that the region index is valid
237  if ( regionIndex >= pdfs_.size() ) {
238  cerr << "ERROR in LauDPDepMapPdf::calcLikelihoodInfo : No PDF supplied for region " << regionIndex << endl;
239  gSystem->Exit(EXIT_FAILURE);
240  }
241 
242  // Evaluate the normalised PDF values
243  LauAbsPdf* pdf = pdfs_[regionIndex];
244 
245  if ( pdf->isDPDependent() ) {
246  pdf->calcLikelihoodInfo(abscissas);
247  } else {
248  pdf->calcLikelihoodInfo(noDPVars);
249  }
250 
251  Double_t result = pdf->getUnNormLikelihood();
252  this->setUnNormPDFVal(result);
253 
254  Double_t norm = pdf->getNorm();
255  this->setNorm(norm);
256 }
257 
259 {
260  // Nothing to do here, since it is already normalized
261  this->setNorm(1.0);
262 }
263 
265 {
266  // This method gives you the maximum possible height of the PDF.
267  // It combines the maximum heights of the two individual PDFs.
268  // So it would give the true maximum if the two individual maxima coincided.
269  // It is guaranteed to always be >= the true maximum.
270 
271  Double_t m13Sq = kinematics->getm13Sq();
272  Double_t m23Sq = kinematics->getm23Sq();
273 
274  // Find which region of the DP we're in
275  UInt_t regionIndex = this->determineDPRegion( m13Sq, m23Sq );
276 
277  // Check that the region index is valid
278  if ( regionIndex >= pdfs_.size() ) {
279  cerr << "ERROR in LauDPDepMapPdf::calcPDFHeight : No PDF supplied for region " << regionIndex << endl;
280  gSystem->Exit(EXIT_FAILURE);
281  }
282 
283  // Evaluate the normalised PDF values
284  LauAbsPdf* pdf = pdfs_[regionIndex];
285 
286  // Update the heights of the individual PDFs
287  pdf->calcPDFHeight( kinematics );
288 
289  // Find the PDF maximum
290  Double_t height = pdf->getMaxHeight();
291  this->setMaxHeight(height);
292 }
293 
294 // Override the base class methods for cacheInfo and calcLikelihoodInfo(UInt_t iEvt).
295 // For both of these we delegate some functionality to the constituent PDFs.
296 
298 {
299  // delegate to the sub-PDFs to cache their information
300  for ( std::vector<LauAbsPdf*>::const_iterator iter = pdfs_.begin(); iter != pdfs_.end(); ++iter){
301  (*iter)->cacheInfo( inputData );
302  }
303 
304  // now check that the DP variables are included in the data
305  Bool_t hasBranch = inputData.haveBranch( "m13Sq" );
306  hasBranch &= inputData.haveBranch( "m23Sq" );
307  if (!hasBranch) {
308  cerr<<"ERROR in LauDPDepMapPdf::cacheInfo : Input data does not contain Dalitz plot variables m13Sq and/or m23Sq."<<endl;
309  return;
310  }
311 
312  // determine whether we are caching our PDF value
313  Bool_t doCaching( this->nFixedParameters() == this->nParameters() );
314  this->cachePDF( doCaching );
315 
316  if ( !doCaching ) {
317  // in this case we seem to be doing a fit where the parameters are floating
318  // so need to mark that the PDF height is no longer up to date
319  this->heightUpToDate(kFALSE);
320  }
321 
322  // clear the vector and reserve enough space
323  UInt_t nEvents = inputData.nEvents();
324  indices_.clear();
325  indices_.reserve(nEvents);
326 
327  // loop through the events, determine the fraction and store
328  for (UInt_t iEvt = 0; iEvt < nEvents; ++iEvt) {
329 
330  const LauFitData& dataValues = inputData.getData(iEvt);
331 
332  Double_t m13Sq = dataValues.find("m13Sq")->second;
333  Double_t m23Sq = dataValues.find("m23Sq")->second;
334 
335  UInt_t regionIndex = this->determineDPRegion( m13Sq, m23Sq );
336  indices_.push_back( regionIndex );
337  }
338 }
339 
341 {
342  // Get the fraction value for this event
343  UInt_t regionIndex = indices_[iEvt];
344 
345  // Evaluate the normalised PDF value
346  LauAbsPdf* pdf = pdfs_[regionIndex];
347  pdf->calcLikelihoodInfo(iEvt);
348 
349  Double_t result = pdf->getUnNormLikelihood();
350  this->setUnNormPDFVal(result);
351 
352  Double_t norm = pdf->getNorm();
353  this->setNorm(norm);
354 }
355 
virtual void setUnNormPDFVal(Double_t unNormPDFVal)
Set the unnormalised likelihood.
Definition: LauAbsPdf.hh:369
File containing declaration of LauDPDepMapPdf class.
virtual Double_t getMinAbscissa() const
Retrieve the minimum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:117
LauDaughters * daughters_
Daughter particles.
virtual Bool_t heightUpToDate() const
Check if the maximum height of the PDF is up to date.
Definition: LauAbsPdf.hh:264
ClassImp(LauAbsCoeffSet)
Lau2DAbsDP * dpDependence_
2D histogram - DP
Double_t getmPrime() const
Get m&#39; value.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:33
virtual UInt_t nParameters() const
Retrieve the number of PDF parameters.
Definition: LauAbsPdf.hh:92
virtual Double_t getUnNormLikelihood() const
Retrieve the unnormalised likelihood value.
Definition: LauAbsPdf.hh:196
File containing declaration of LauDaughters class.
virtual void calcPDFHeight(const LauKinematics *kinematics)=0
Calculate the maximum height of the PDF.
UInt_t determineDPRegion(Double_t m13Sq, Double_t m23Sq) const
Determine the DP region.
virtual void calcPDFHeight(const LauKinematics *kinematics)
Calculate the PDF height.
virtual void setNorm(Double_t norm)
Set the normalisation factor.
Definition: LauAbsPdf.hh:325
virtual Bool_t checkRange(const LauAbscissas &abscissas) const
Check that all abscissas are within their allowed ranges.
Definition: LauAbsPdf.cc:213
TH1 * dpAxisDependence_
1D histogram - DP axis
virtual ~LauDPDepMapPdf()
Destructor.
virtual Double_t interpolateXY(Double_t x, Double_t y) const =0
Perform the interpolation.
std::map< TString, Double_t > LauFitData
Type for holding event data.
Double_t getm23Sq() const
Get the m23 invariant mass square.
Bool_t squareDP() const
Determine to use or not the square Dalitz plot.
Definition: LauDaughters.hh:84
std::vector< LauAbsPdf * > pdfs_
The PDFs.
Double_t distanceFromDPCentre() const
Calculate the distance from the currently set (m13Sq, m23Sq) point to the centre of the Dalitz plot (...
LauDPDepMapPdf(const std::vector< LauAbsPdf * > &pdfs, const LauDaughters *daughters, const TH2 *dpHisto, Bool_t upperHalf=kFALSE)
Constructor - map described by 2D histogram.
virtual Double_t getMaxHeight() const
Retrieve the maximum height.
Definition: LauAbsPdf.hh:221
std::vector< UInt_t > indices_
Cached indices values.
virtual void cacheInfo(const LauFitDataTree &inputData)
Cache information from data.
Double_t calcThirdMassSq(const Double_t firstMassSq, const Double_t secondMassSq) const
Calculate the third invariant mass square from the two provided (e.g. mjkSq from mijSq and mikSq) ...
File containing declaration of LauParameter class.
LauKinematics * getKinematics()
Retrieve the Dalitz plot kinematics.
virtual std::vector< TString > varNames() const
Retrieve the names of the abscissas.
Definition: LauAbsPdf.cc:103
virtual Double_t getMaxAbscissa() const
Retrieve the maximum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:123
virtual Bool_t isDPDependent() const
Specifies whether or not the PDF is DP dependent.
Definition: LauAbsPdf.hh:111
virtual void setMaxHeight(Double_t maxHeight)
Set the maximum height.
Definition: LauAbsPdf.hh:331
File containing declaration of Lau2DHistDP class.
virtual void addParameters(std::vector< LauAbsRValue * > &params)
Add parameters to the PDF.
Definition: LauAbsPdf.cc:423
virtual const std::vector< LauAbsRValue * > & getParameters() const
Retrieve the parameters of the PDF, e.g. so that they can be loaded into a fit.
Definition: LauAbsPdf.hh:239
void updateKinematics(const Double_t m13Sq, const Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
virtual Bool_t cachePDF() const
Check if the PDF is to be cached.
Definition: LauAbsPdf.hh:276
Double_t getm13Sq() const
Get the m13 invariant mass square.
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)=0
Calculate the likelihood (and all associated information) given value(s) of the abscissa(s) ...
Class for defining a 2D DP histogram.
Definition: Lau2DHistDP.hh:34
File containing LauConstants namespace.
Bool_t haveBranch(const TString &name) const
Check if the named branch is stored.
Double_t getThetaPrime() const
Get theta&#39; value.
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:41
Class for calculating 3-body kinematic quantities.
virtual UInt_t nFixedParameters() const
Retrieve the number of fixed PDF parameters.
Definition: LauAbsPdf.cc:113
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)
Calculate the likelihood (and intermediate info) for a given abscissa.
UInt_t nEvents() const
Retrieve the number of events.
virtual void calcNorm()
Calculate the normalisation.
virtual UInt_t nInputVars() const
Retrieve the number of abscissas.
Definition: LauAbsPdf.hh:104
Pure abstract base class for defining a parameter containing an R value.
Definition: LauAbsRValue.hh:29
virtual Double_t getNorm() const
Retrieve the normalisation factor.
Definition: LauAbsPdf.hh:202
Class to allow having different PDFs in different regions of the DP.
Class to store the input fit variables.
DPAxis dpAxis_
The DP axis we depend on.
std::vector< Double_t > LauAbscissas
The type used for containing multiple abscissa values.
Definition: LauAbsPdf.hh:45
DPAxis
Define possibilties for the DP axes.