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