laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauDPDepSumPdf.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2009 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 "LauConstants.hh"
36 #include "LauDPDepSumPdf.hh"
37 #include "LauDaughters.hh"
38 #include "LauEffModel.hh"
39 #include "LauParameter.hh"
40 
41 #include "TMath.h"
42 #include "TSystem.h"
43 
45  LauAbsPdf* pdf2,
46  const LauDaughters* daughters,
47  const TH2* dpHisto,
48  Bool_t upperHalf,
49  Bool_t useSpline ) :
50  LauAbsPdf( pdf1 ? pdf1->varNames() : std::vector<TString>(),
51  std::vector<LauAbsRValue*>(),
52  pdf1 ? pdf1->getMinAbscissas() : LauFitData(),
53  pdf1 ? pdf1->getMaxAbscissas() : LauFitData() ),
54  daughters_( new LauDaughters( *daughters ) ),
55  pdf1_( pdf1 ),
56  pdf2_( pdf2 ),
57  frac_( 0 ),
58  fracVal_( 0.5 ),
59  dpDependence_( new LauEffModel( daughters, 0 ) ),
60  dpAxis_( CentreDist )
61 {
62  // Constructor for the sum PDF.
63  // We are defining the sum as:
64  // f x (PDF1/S(PDF1)) + (1-f) x (PDF2/S(PDF2))
65  // where f is the fraction, x is multiplication, PDFi is the i'th PDF,
66  // and S(PDFi) is the integral of the i'th PDF.
67  // The value of the fraction is read from the DP histogram.
68 
69  if ( useSpline ) {
70  dpDependence_->setEffSpline( dpHisto, kFALSE, 0.0, 0.0, upperHalf, daughters->squareDP() );
71  } else {
72  dpDependence_->setEffHisto( dpHisto, kTRUE, kFALSE, 0.0, 0.0, upperHalf, daughters->squareDP() );
73  }
74 
75  // So the first thing we have to do is check the pointers are all valid.
76  if ( ! pdf1 || ! pdf2 ) {
77  cerr << "ERROR in LauDPDepSumPdf constructor: one of the 2 PDF pointers is null." << endl;
78  gSystem->Exit( EXIT_FAILURE );
79  }
80 
81  // Next check that the abscissa ranges are the same for each PDF
82  if ( pdf1->getMinAbscissa() != pdf2->getMinAbscissa() ) {
83  cerr << "ERROR in LauDPDepSumPdf constructor: minimum abscissa values not the same for the two PDFs."
84  << endl;
85  gSystem->Exit( EXIT_FAILURE );
86  }
87  if ( pdf1->getMaxAbscissa() != pdf2->getMaxAbscissa() ) {
88  cerr << "ERROR in LauDPDepSumPdf constructor: maximum abscissa values not the same for the two PDFs."
89  << endl;
90  gSystem->Exit( EXIT_FAILURE );
91  }
92 
93  // Also check that both PDFs are expecting the same number of input variables
94  if ( pdf1->nInputVars() != pdf2->nInputVars() ) {
95  cerr << "ERROR in LauDPDepSumPdf constructor: number of input variables not the same for the two PDFs."
96  << endl;
97  gSystem->Exit( EXIT_FAILURE );
98  }
99 
100  // Also check that both PDFs are expecting the same variable name(s)
101  if ( pdf1->varNames() != pdf2->varNames() ) {
102  cerr << "ERROR in LauDPDepSumPdf constructor: variable name(s) not the same for the two PDFs."
103  << endl;
104  gSystem->Exit( EXIT_FAILURE );
105  }
106 
107  // Then we need to grab all the parameters and pass them to the base class.
108  // This is so that when we are asked for them they can be put into the fit.
109  // The number of parameters is the number in PDF1 + the number in PDF2.
110  UInt_t nPar( pdf1->nParameters() + pdf2->nParameters() );
111  std::vector<LauAbsRValue*> params;
112  params.reserve( nPar );
113  std::vector<LauAbsRValue*>& pdf1pars = pdf1->getParameters();
114  std::vector<LauAbsRValue*>& pdf2pars = pdf2->getParameters();
115  for ( std::vector<LauAbsRValue*>::iterator iter = pdf1pars.begin(); iter != pdf1pars.end();
116  ++iter ) {
117  params.push_back( *iter );
118  }
119  for ( std::vector<LauAbsRValue*>::iterator iter = pdf2pars.begin(); iter != pdf2pars.end();
120  ++iter ) {
121  params.push_back( *iter );
122  }
123  this->addParameters( params );
124 
125  // Cache the normalisation factor
126  this->calcNorm();
127 }
128 
130  LauAbsPdf* pdf2,
131  LauParameter* frac,
132  const LauDaughters* daughters,
133  const std::vector<Double_t>& fracCoeffs,
134  DPAxis dpAxis ) :
135  LauAbsPdf( pdf1 ? pdf1->varNames() : std::vector<TString>(),
136  std::vector<LauAbsRValue*>(),
137  pdf1 ? pdf1->getMinAbscissas() : LauFitData(),
138  pdf1 ? pdf1->getMaxAbscissas() : LauFitData() ),
139  daughters_( new LauDaughters( *daughters ) ),
140  pdf1_( pdf1 ),
141  pdf2_( pdf2 ),
142  frac_( frac ),
143  fracVal_( frac ? frac->unblindValue() : 0.0 ),
144  dpDependence_( 0 ),
145  fracCoeffs_( fracCoeffs ),
146  dpAxis_( dpAxis )
147 {
148  // Constructor for the sum PDF.
149  // We are defining the sum as:
150  // f x (PDF1/S(PDF1)) + (1-f) x (PDF2/S(PDF2))
151  // where f is the fraction, x is multiplication, PDFi is the i'th PDF,
152  // and S(PDFi) is the integral of the i'th PDF.
153  // The value of the fraction has a polynomial dependence on one of
154  // the DP axes.
155 
156  // So the first thing we have to do is check the pointers are all valid.
157  if ( ! pdf1 || ! pdf2 ) {
158  cerr << "ERROR in LauDPDepSumPdf constructor: one of the 2 PDF pointers is null." << endl;
159  gSystem->Exit( EXIT_FAILURE );
160  }
161  if ( ! frac ) {
162  cerr << "ERROR in LauDPDepSumPdf constructor: the fraction parameter pointer is null."
163  << endl;
164  gSystem->Exit( EXIT_FAILURE );
165  }
166 
167  // Next check that the abscissa ranges are the same for each PDF
168  if ( pdf1->getMinAbscissa() != pdf2->getMinAbscissa() ) {
169  cerr << "ERROR in LauDPDepSumPdf constructor: minimum abscissa values not the same for the two PDFs."
170  << endl;
171  gSystem->Exit( EXIT_FAILURE );
172  }
173  if ( pdf1->getMaxAbscissa() != pdf2->getMaxAbscissa() ) {
174  cerr << "ERROR in LauDPDepSumPdf constructor: maximum abscissa values not the same for the two PDFs."
175  << endl;
176  gSystem->Exit( EXIT_FAILURE );
177  }
178 
179  // Also check that both PDFs are expecting the same number of input variables
180  if ( pdf1->nInputVars() != pdf2->nInputVars() ) {
181  cerr << "ERROR in LauDPDepSumPdf constructor: number of input variables not the same for the two PDFs."
182  << endl;
183  gSystem->Exit( EXIT_FAILURE );
184  }
185 
186  // Also check that both PDFs are expecting the same variable name(s)
187  if ( pdf1->varNames() != pdf2->varNames() ) {
188  cerr << "ERROR in LauDPDepSumPdf constructor: variable name(s) not the same for the two PDFs."
189  << endl;
190  gSystem->Exit( EXIT_FAILURE );
191  }
192 
193  // Then we need to grab all the parameters and pass them to the base class.
194  // This is so that when we are asked for them they can be put into the fit.
195  // The number of parameters is the number in PDF1 + the number in PDF2.
196  UInt_t nPar( pdf1->nParameters() + pdf2->nParameters() + 1 );
197  std::vector<LauAbsRValue*> params;
198  params.reserve( nPar );
199  params.push_back( frac );
200  std::vector<LauAbsRValue*>& pdf1pars = pdf1->getParameters();
201  std::vector<LauAbsRValue*>& pdf2pars = pdf2->getParameters();
202  for ( std::vector<LauAbsRValue*>::iterator iter = pdf1pars.begin(); iter != pdf1pars.end();
203  ++iter ) {
204  params.push_back( *iter );
205  }
206  for ( std::vector<LauAbsRValue*>::iterator iter = pdf2pars.begin(); iter != pdf2pars.end();
207  ++iter ) {
208  params.push_back( *iter );
209  }
210  this->addParameters( params );
211 
212  // Now check that we can find the fraction parameter ok
213  frac_ = this->findParameter( "frac" );
214  if ( frac_ == 0 ) {
215  cerr << "ERROR in LauDPDepSumPdf constructor: parameter \"frac\" not found." << endl;
216  gSystem->Exit( EXIT_FAILURE );
217  }
218 
219  // Cache the normalisation factor
220  this->calcNorm();
221 }
222 
224 {
225  // Destructor
226  delete daughters_;
227  daughters_ = 0;
228  delete dpDependence_;
229  dpDependence_ = 0;
230 }
231 
233 {
234  // Check that the given abscissa is within the allowed range
235  if ( ! this->checkRange( abscissas ) ) {
236  gSystem->Exit( EXIT_FAILURE );
237  }
238 
239  LauAbscissas noDPVars( 1 );
240  noDPVars[0] = abscissas[0];
241 
242  // Evaluate the normalised PDF values
243  if ( pdf1_->isDPDependent() ) {
244  pdf1_->calcLikelihoodInfo( abscissas );
245  } else {
246  pdf1_->calcLikelihoodInfo( noDPVars );
247  }
248  if ( pdf2_->isDPDependent() ) {
249  pdf2_->calcLikelihoodInfo( abscissas );
250  } else {
251  pdf2_->calcLikelihoodInfo( noDPVars );
252  }
253  Double_t result1 = pdf1_->getLikelihood();
254  Double_t result2 = pdf2_->getLikelihood();
255 
256  // Determine the fraction
257  // The DP variables will be abscissas[nInputVars] and
258  // abscissas[nInputVars+1] (if present).
259  UInt_t nVars = this->nInputVars();
260  if ( abscissas.size() == nVars + 2 ) {
261  Double_t m13Sq = abscissas[nVars];
262  Double_t m23Sq = abscissas[nVars + 1];
263  LauKinematics* kinematics = daughters_->getKinematics();
264  if ( dpDependence_ ) {
265  kinematics->updateKinematics( m13Sq, m23Sq );
266  fracVal_ = dpDependence_->calcEfficiency( kinematics );
267  } else {
269  Double_t dpPos( 0.0 );
270  if ( dpAxis_ == M12 ) {
271  dpPos = kinematics->calcThirdMassSq( m13Sq, m23Sq );
272  } else if ( dpAxis_ == M13 ) {
273  dpPos = m13Sq;
274  } else if ( dpAxis_ == M23 ) {
275  dpPos = m23Sq;
276  } else if ( dpAxis_ == MMIN ) {
277  dpPos = TMath::Min( m13Sq, m23Sq );
278  } else if ( dpAxis_ == MMAX ) {
279  dpPos = TMath::Max( m13Sq, m23Sq );
280  } else {
281  dpPos = kinematics->distanceFromDPCentre( m13Sq, m23Sq );
282  }
283  this->scaleFrac( dpPos );
284  }
285  }
286 
287  // Add them together
288  Double_t result = fracVal_ * result1 + ( 1.0 - fracVal_ ) * result2;
289  this->setUnNormPDFVal( result );
290 }
291 
292 void LauDPDepSumPdf::scaleFrac( Double_t dpPos )
293 {
294  Int_t power = 1;
295  for ( std::vector<Double_t>::const_iterator iter = fracCoeffs_.begin(); iter != fracCoeffs_.end();
296  ++iter ) {
297  Double_t coeff = ( *iter );
298  fracVal_ += coeff * TMath::Power( dpPos, power );
299  ++power;
300  }
301 }
302 
304 {
305  // Nothing to do here, since it is already normalized
306  this->setNorm( 1.0 );
307 }
308 
310 {
311  // This method gives you the maximum possible height of the PDF.
312  // It combines the maximum heights of the two individual PDFs.
313  // So it would give the true maximum if the two individual maxima coincided.
314  // It is guaranteed to always be >= the true maximum.
315 
316  // Update the heights of the individual PDFs
317  pdf1_->calcPDFHeight( kinematics );
318  pdf2_->calcPDFHeight( kinematics );
319 
320  // Get the up to date parameter values
321  if ( dpDependence_ ) {
322  fracVal_ = dpDependence_->calcEfficiency( kinematics );
323  } else {
325  Double_t dpPos( 0.0 );
326  if ( dpAxis_ == M12 ) {
327  dpPos = kinematics->getm12Sq();
328  } else if ( dpAxis_ == M13 ) {
329  dpPos = kinematics->getm13Sq();
330  } else if ( dpAxis_ == M23 ) {
331  dpPos = kinematics->getm23Sq();
332  } else if ( dpAxis_ == MMIN ) {
333  Double_t m13Sq = kinematics->getm13Sq();
334  Double_t m23Sq = kinematics->getm23Sq();
335  dpPos = TMath::Min( m13Sq, m23Sq );
336  } else if ( dpAxis_ == MMAX ) {
337  Double_t m13Sq = kinematics->getm13Sq();
338  Double_t m23Sq = kinematics->getm23Sq();
339  dpPos = TMath::Max( m13Sq, m23Sq );
340  } else {
341  dpPos = kinematics->distanceFromDPCentre();
342  }
343  this->scaleFrac( dpPos );
344  }
345 
346  // Find the (un-normalised) individual PDF maxima
347  Double_t height1 = pdf1_->getMaxHeight();
348  Double_t height2 = pdf2_->getMaxHeight();
349 
350  // Get the individual PDF normalisation factors
351  Double_t norm1 = pdf1_->getNorm();
352  Double_t norm2 = pdf2_->getNorm();
353 
354  // Calculate the normalised individual PDF maxima
355  height1 /= norm1;
356  height2 /= norm2;
357 
358  // Combine these heights together
359  Double_t height = fracVal_ * height1 + ( 1 - fracVal_ ) * height2;
360  this->setMaxHeight( height );
361 }
362 
363 // Override the base class methods for cacheInfo and calcLikelihoodInfo(UInt_t iEvt).
364 // For both of these we delegate to the two constituent PDFs.
365 
367 {
368  // delegate to the two sub-PDFs to cache their information
369  pdf1_->cacheInfo( inputData );
370  pdf2_->cacheInfo( inputData );
371 
372  // now check that the DP variables are included in the data
373  Bool_t hasBranch = inputData.haveBranch( "m13Sq" );
374  hasBranch &= inputData.haveBranch( "m23Sq" );
375  if ( ! hasBranch ) {
376  cerr << "ERROR in LauDPDepSumPdf::cacheInfo : Input data does not contain Dalitz plot variables m13Sq and/or m23Sq."
377  << endl;
378  return;
379  }
380 
381  // determine whether we are caching our PDF value
382  Bool_t doCaching( this->nFixedParameters() == this->nParameters() );
383  this->cachePDF( doCaching );
384 
385  if ( ! doCaching ) {
386  // in this case we seem to be doing a fit where the parameters are floating
387  // so need to mark that the PDF height is no longer up to date
388  this->heightUpToDate( kFALSE );
389  }
390 
391  // clear the vector and reserve enough space
392  UInt_t nEvents = inputData.nEvents();
393  fractions_.clear();
394  fractions_.reserve( nEvents );
395 
396  // loop through the events, determine the fraction and store
397  for ( UInt_t iEvt = 0; iEvt < nEvents; iEvt++ ) {
398 
399  const LauFitData& dataValues = inputData.getData( iEvt );
400 
401  Double_t m13Sq = dataValues.find( "m13Sq" )->second;
402  Double_t m23Sq = dataValues.find( "m23Sq" )->second;
403 
404  LauKinematics* kinematics = daughters_->getKinematics();
405  if ( dpDependence_ ) {
406  // if we're using the histogram then just
407  // determine the fraction and store
408  kinematics->updateKinematics( m13Sq, m23Sq );
409  fracVal_ = dpDependence_->calcEfficiency( kinematics );
410  } else {
411  // if we're scaling the fraction parameter then we
412  // just store the scaling info since the parameter
413  // might be floating
415  Double_t dpPos( 0.0 );
416  if ( dpAxis_ == M12 ) {
417  dpPos = kinematics->calcThirdMassSq( m13Sq, m23Sq );
418  } else if ( dpAxis_ == M13 ) {
419  dpPos = m13Sq;
420  } else if ( dpAxis_ == M23 ) {
421  dpPos = m23Sq;
422  } else if ( dpAxis_ == MMIN ) {
423  dpPos = TMath::Min( m13Sq, m23Sq );
424  } else if ( dpAxis_ == MMAX ) {
425  dpPos = TMath::Max( m13Sq, m23Sq );
426  } else {
427  dpPos = kinematics->distanceFromDPCentre( m13Sq, m23Sq );
428  }
429  this->scaleFrac( dpPos );
431  }
432 
433  fractions_.push_back( fracVal_ );
434  }
435 }
436 
438 {
439  // Get the fraction value for this event
440  fracVal_ = fractions_[iEvt];
441  if ( frac_ ) {
442  // if we're scaling the parameter then need to add the
443  // current value of the parameter
445  }
446 
447  // Evaluate the normalised PDF values
448  pdf1_->calcLikelihoodInfo( iEvt );
449  pdf2_->calcLikelihoodInfo( iEvt );
450  Double_t result1 = pdf1_->getLikelihood();
451  Double_t result2 = pdf2_->getLikelihood();
452 
453  // Add them together
454  Double_t result = fracVal_ * result1 + ( 1 - fracVal_ ) * result2;
455  this->setUnNormPDFVal( result );
456 }
virtual ~LauDPDepSumPdf()
Destructor.
virtual void calcPDFHeight(const LauKinematics *kinematics)
Calculate the PDF height.
LauDPDepSumPdf(LauAbsPdf *pdf1, LauAbsPdf *pdf2, const LauDaughters *daughters, const TH2 *dpHisto, Bool_t upperHalf=kFALSE, Bool_t useSpline=kFALSE)
Constructor - fraction determined by 2D histogram.
File containing declaration of LauEffModel class.
virtual UInt_t nFixedParameters() const
Retrieve the number of fixed PDF parameters.
Definition: LauAbsPdf.cc:133
LauAbsPdf * pdf1_
First PDF.
virtual void calcNorm()
Calculate the normalisation.
DPAxis dpAxis_
The DP axis we depend on.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
virtual Double_t getMaxAbscissa() const
Retrieve the maximum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:141
Double_t unblindValue() const
The unblinded value of the parameter.
virtual void cacheInfo(const LauFitDataTree &inputData)
Cache information from data.
File containing declaration of LauDPDepSumPdf class.
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)
Calculate the likelihood (and intermediate info) for a given abscissa.
virtual Double_t getMaxHeight() const
Retrieve the maximum height.
Definition: LauAbsPdf.hh:244
virtual Double_t getMinAbscissa() const
Retrieve the minimum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:135
Double_t getm12Sq() const
Get the m12 invariant mass square.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
Double_t fracVal_
Fractional contribution of first PDF to the final PDF.
LauAbsRValue * frac_
Fractional contribution of first PDF to the final PDF.
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
virtual Bool_t cachePDF() const
Check if the PDF is to be cached.
Definition: LauAbsPdf.hh:299
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.
Double_t calcEfficiency(const LauKinematics *kinematics) const
Determine the efficiency for a given point in the Dalitz plot.
Definition: LauEffModel.cc:367
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.
void scaleFrac(Double_t dpPos)
Scale fraction according to DP position.
LauKinematics * getKinematics()
Retrieve the Dalitz plot kinematics.
virtual void calcPDFHeight(const LauKinematics *kinematics)=0
Calculate the maximum height of the PDF.
virtual LauAbsRValue * findParameter(const TString &parName)
Retrieve the specified parameter.
Definition: LauAbsPdf.cc:431
virtual void cacheInfo(const LauFitDataTree &inputData)
Cache information from data.
Definition: LauAbsPdf.cc:275
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)
Class to store the input fit variables.
virtual UInt_t nParameters() const
Retrieve the number of PDF parameters.
Definition: LauAbsPdf.hh:109
Class that implements the efficiency description across the signal Dalitz plot.
Definition: LauEffModel.hh:50
virtual void setNorm(Double_t norm)
Set the normalisation factor.
Definition: LauAbsPdf.hh:348
virtual Double_t unblindValue() const =0
The unblinded value of the parameter.
virtual std::vector< TString > varNames() const
Retrieve the names of the abscissas.
Definition: LauAbsPdf.cc:122
LauAbsPdf * pdf2_
Second PDF.
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
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
virtual Double_t getLikelihood() const
Retrieve the normalised likelihood value.
Definition: LauAbsPdf.cc:403
Pure abstract base class for defining a parameter containing an R value.
Definition: LauAbsRValue.hh:45
void setEffSpline(const TH2 *effHisto, 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 spline based on a predetermined 2D histog...
Definition: LauEffModel.cc:154
File containing LauConstants namespace.
const std::vector< Double_t > fracCoeffs_
Polynomial used to scale fraction.
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
LauEffModel * dpDependence_
DP dependence.
Class for calculating 3-body kinematic quantities.
LauDaughters * daughters_
Daughter particles.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getm23Sq() const
Get the m23 invariant mass square.
DPAxis
Define possibilties for the DP axes.
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
std::vector< Double_t > fractions_
Cached values of fractions.
UInt_t nEvents() const
Retrieve the number of events.
virtual Bool_t checkRange(const LauAbscissas &abscissas) const
Check that all abscissas are within their allowed ranges.
Definition: LauAbsPdf.cc:243