laura is hosted by Hepforge, IPPP Durham
Laura++  v2r2
A maximum likelihood fitting package for performing Dalitz-plot analysis.
Lau2DHistPdf.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2008 - 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 #include <vector>
18 
19 #include "TAxis.h"
20 #include "TH2.h"
21 #include "TRandom.h"
22 #include "TSystem.h"
23 
24 #include "Lau1DHistPdf.hh"
25 #include "Lau2DHistPdf.hh"
26 #include "LauRandom.hh"
27 
28 class LauParameter;
29 
31 
32 
33 Lau2DHistPdf::Lau2DHistPdf(const std::vector<TString>& theVarNames, const TH2* hist,
34  const LauFitData& minVals, const LauFitData& maxVals,
35  Bool_t useInterpolation, Bool_t fluctuateBins) :
36  LauAbsPdf(theVarNames, std::vector<LauAbsRValue*>(), minVals, maxVals),
37  hist_(hist ? dynamic_cast<TH2*>(hist->Clone()) : 0),
38  xProj_(0),
39  yProj_(0),
40  xVarPdf_(0),
41  yVarPdf_(0),
42  xName_(""),
43  yName_(""),
44  nBinsX_(0),
45  nBinsY_(0),
46  minX_(0.0),
47  maxX_(0.0),
48  minY_(0.0),
49  maxY_(0.0),
50  rangeX_(0.0),
51  rangeY_(0.0),
52  binXWidth_(0.0),
53  binYWidth_(0.0),
54  invBinXWidth_(0.0),
55  invBinYWidth_(0.0),
56  useInterpolation_(useInterpolation),
57  fluctuateBins_(fluctuateBins)
58 {
59  // Constructor
60 
61  // Check that we have two variables
62  if ( this->nInputVars() != 2 ) {
63  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Have not been provided with exactly two variables." << std::endl;
64  gSystem->Exit(EXIT_FAILURE);
65  }
66 
67  // Set the variable names from the abstract class
68  xName_ = this->varNames()[0];
69  yName_ = this->varNames()[1];
70 
71  // Set the variable limits from the abstract class
72  minX_ = this->getMinAbscissa(xName_);
73  maxX_ = this->getMaxAbscissa(xName_);
74  minY_ = this->getMinAbscissa(yName_);
75  maxY_ = this->getMaxAbscissa(yName_);
76  rangeX_ = maxX_ - minX_;
77  rangeY_ = maxY_ - minY_;
78 
79  // Have we got a valid histogram
80  if ( hist_ == 0 ) {
81  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Histogram pointer is null." << std::endl;
82  gSystem->Exit(EXIT_FAILURE);
83  }
84 
85  // Set the directory for the histogram
86  hist_->SetDirectory(0);
87 
88  // Save various attributes of the histogram
89  nBinsX_ = hist_->GetNbinsX();
90  nBinsY_ = hist_->GetNbinsY();
91 
92  TAxis* xAxis = hist_->GetXaxis();
93  Double_t xAxisMin = xAxis->GetXmin();
94  Double_t xAxisMax = xAxis->GetXmax();
95 
96  TAxis* yAxis = hist_->GetYaxis();
97  Double_t yAxisMin = yAxis->GetXmin();
98  Double_t yAxisMax = yAxis->GetXmax();
99 
100  // Check that axis ranges corresponds to ranges of abscissas
101  if (TMath::Abs(minX_ - xAxisMin)>1e-6) {
102  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Histogram x-axis minimum: " << xAxisMin <<
103  " does not correspond to " << xName_ << " minimum: " << minX_ << "." << std::endl;
104  gSystem->Exit(EXIT_FAILURE);
105  }
106  if (TMath::Abs(maxX_ - xAxisMax)>1e-6) {
107  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Histogram x-axis maximum: " << xAxisMax <<
108  " does not correspond to " << xName_ << " maximum: " << maxX_ << "." << std::endl;
109  gSystem->Exit(EXIT_FAILURE);
110  }
111  if (TMath::Abs(minY_ - yAxisMin)>1e-6) {
112  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Histogram y-axis minimum: " << yAxisMin <<
113  " does not correspond to " << yName_ << " minimum: " << minY_ << "." << std::endl;
114  gSystem->Exit(EXIT_FAILURE);
115  }
116  if (TMath::Abs(maxY_ - yAxisMax)>1e-6) {
117  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Histogram y-axis maximum: " << yAxisMax <<
118  " does not correspond to " << yName_ << " maximum: " << maxY_ << "." << std::endl;
119  gSystem->Exit(EXIT_FAILURE);
120  }
121 
122  // Calculate the bin widths and inverse bin widths
123  binXWidth_ = static_cast<Double_t>(TMath::Abs(rangeX_)/(nBinsX_*1.0));
124  binYWidth_ = static_cast<Double_t>(TMath::Abs(rangeY_)/(nBinsY_*1.0));
125  if (binXWidth_ > 1e-10) {invBinXWidth_ = 1.0/binXWidth_;}
126  if (binYWidth_ > 1e-10) {invBinYWidth_ = 1.0/binYWidth_;}
127 
128  // If the bins are to be fluctuated then do so now before
129  // calculating anything that depends on the bin content.
130  if (fluctuateBins) {
131  this->doBinFluctuation();
132  }
133 
134  // Create the projections and 1D PDFs from those
135  xProj_ = hist_->ProjectionX();
136  yProj_ = hist_->ProjectionY();
137  if ( !xProj_ || !yProj_ ) {
138  std::cerr << "ERROR in Lau2DHistPdf::Lau2DHistPdf : Can't get X or Y projection from 2D histogram." << std::endl;
139  gSystem->Exit(EXIT_FAILURE);
140  }
141  xVarPdf_ = new Lau1DHistPdf(xName_, xProj_, minX_, maxX_, useInterpolation_, fluctuateBins_);
142  yVarPdf_ = new Lau1DHistPdf(yName_, yProj_, minY_, maxY_, useInterpolation_, fluctuateBins_);
143 
144  // Calculate the PDF normalisation.
145  this->calcNorm();
146 
147  // And check it is OK.
148  this->checkNormalisation();
149 }
150 
152 {
153  // Destructor
154  delete hist_; hist_ = 0;
155  delete xProj_; xProj_ = 0;
156  delete yProj_; yProj_ = 0;
157  delete xVarPdf_; xVarPdf_ = 0;
158  delete yVarPdf_; yVarPdf_ = 0;
159 }
160 
161 Lau2DHistPdf::Lau2DHistPdf(const Lau2DHistPdf& other) : LauAbsPdf(other.varName(), other.getParameters(), other.getMinAbscissa(), other.getMaxAbscissa())
162 {
163  hist_ = dynamic_cast<TH2*>(other.hist_->Clone());
164  xName_ = other.xName_;
165  yName_ = other.yName_;
166  nBinsX_ = other.nBinsX_;
167  nBinsY_ = other.nBinsY_;
168  minX_ = other.minX_;
169  maxX_ = other.maxX_;
170  minY_ = other.minY_;
171  maxY_ = other.maxY_;
172  rangeX_ = other.rangeX_;
173  rangeY_ = other.rangeY_;
174  binXWidth_ = other.binXWidth_;
175  binYWidth_ = other.binYWidth_;
180 
181  xProj_ = hist_->ProjectionX();
182  yProj_ = hist_->ProjectionY();
185 
186  this->setRandomFun(other.getRandomFun());
187  this->calcNorm();
188 }
189 
190 void Lau2DHistPdf::calcPDFHeight( const LauKinematics* /*kinematics*/ )
191 {
192  if (this->heightUpToDate()) {
193  return;
194  }
195 
196  // Get the maximum height of the histogram
197  Int_t maxBin = hist_->GetMaximumBin();
198  Double_t height = hist_->GetBinContent(maxBin);
199  this->setMaxHeight(height);
200 }
201 
203 {
204  // Calculate the histogram normalisation.
205 
206  // Loop over the total x and y range to get the total area under x
207  // and y. Just sum the contributions up using 1e-3 increments of
208  // the range in x and y. Multiply the end result by dx and dy.
209 
210  Double_t dx(1e-3*rangeX_), dy(1e-3*rangeY_);
211  Double_t area(0.0);
212 
213  Double_t x(minX_ + dx/2.0);
214  while (x > minX_ && x < maxX_) {
215  Double_t y(minY_ + dy/2.0);
216  while (y > minY_ && y < maxY_) {
217  area += this->interpolateXY(x,y);
218  y += dy;
219  } // y while loop
220  x += dx;
221  } // x while loop
222 
223  Double_t norm = area*dx*dy;
224  this->setNorm(norm);
225 }
226 
228 {
229  Double_t dx(1e-3*rangeX_), dy(1e-3*rangeY_);
230  Double_t area(0.0);
231  Double_t areaNoNorm(0.0);
232 
233  // Preserve the value of a variable we change temporarily
234  Bool_t interpolate = useInterpolation_;
235 
236  Double_t x(minX_ + dx/2.0);
237  while (x > minX_ && x < maxX_) {
238  Double_t y(minY_ + dy/2.0);
239  while (y > minY_ && y < maxY_) {
240  area += this->interpolateXYNorm(x,y);
241  useInterpolation_ = kFALSE;
242  areaNoNorm += this->interpolateXY(x,y);
243  useInterpolation_ = kTRUE;
244  y += dy;
245  } // y while loop
246  x += dx;
247  } // x while loop
248  Double_t norm = area*dx*dy;
249 
250  useInterpolation_ = interpolate; //Restore old value
251 
252  std::cout << "INFO in Lau2DHistPdf::checkNormalisation : Area = " << area << ", dx = " << dx << ", dy = " << dy << ", dx*dy = " << dx*dy << std::endl;
253  std::cout << " : Area with no norm = " << areaNoNorm << "*dx*dy = " << areaNoNorm*dx*dy << std::endl;
254  std::cout << " : The total area of the normalised histogram PDF is " << norm << std::endl;
255 }
256 
257 Double_t Lau2DHistPdf::getBinHistValue(Int_t i, Int_t j) const
258 {
259  // Check that x co-ord is in range [1 , nBinsX_]
260  if ((i < 1) || (i > nBinsX_)) {
261  return 0.0;
262  }
263  // Check that y co-ord is in range [1 , nBinsY_]
264  if ((j < 1) || (j > nBinsY_)) {
265  return 0.0;
266  }
267  Double_t value = hist_->GetBinContent(i,j);
268  // protect against negative values
269  if ( value < 0.0 ) {
270  std::cerr << "WARNING in Lau2DHistPdf::getBinHistValue : Negative bin content set to zero!" << std::endl;
271  value = 0.0;
272  }
273  return value;
274 }
275 
276 Double_t Lau2DHistPdf::interpolateXYNorm(Double_t x, Double_t y) const
277 {
278  // Get the normalised interpolated value.
279  Double_t value = this->interpolateXY(x,y);
280  Double_t norm = this->getNorm();
281  return value/norm;
282 }
283 
284 Double_t Lau2DHistPdf::interpolateXY(Double_t x, Double_t y) const
285 {
286  // This function returns the interpolated value of the histogram
287  // function for the given values of x and y by finding the adjacent
288  // bins and extrapolating using weights based on the inverse
289  // distance of the point from the adajcent bin centres.
290 
291  // Find the 2D histogram bin for x and y
292  Int_t i = static_cast<Int_t>((x - minX_)*invBinXWidth_)+1;
293  Int_t j = static_cast<Int_t>((y - minY_)*invBinYWidth_)+1;
294 
295  if (i > nBinsX_) {i = nBinsX_;}
296  if (j > nBinsY_) {j = nBinsY_;}
297 
298  // Ask whether we want to do the interpolation, since this method is
299  // not reliable for low statistics histograms.
300  if (useInterpolation_ == kFALSE) {
301  return this->getBinHistValue(i,j);
302  }
303 
304  // Find the bin centres (actual co-ordinate positions, not histogram indices)
305  Double_t cbinx = static_cast<Double_t>(i-0.5)*binXWidth_ + minX_;
306  Double_t cbiny = static_cast<Double_t>(j-0.5)*binYWidth_ + minY_;
307 
308  // Find the adjacent bins
309  Double_t deltax = x - cbinx;
310  Double_t deltay = y - cbiny;
311 
312  Int_t i_adj(0), j_adj(0);
313  if (deltax > 0.0) {
314  i_adj = i + 1;
315  } else {
316  i_adj = i - 1;
317  }
318  if (deltay > 0.0) {
319  j_adj = j + 1;
320  } else {
321  j_adj = j - 1;
322  }
323 
324  Bool_t isXBoundary(kFALSE), isYBoundary(kFALSE);
325 
326  Double_t value(0.0);
327 
328  if (i_adj > nBinsX_ || i_adj < 1) {isYBoundary = kTRUE;}
329  if (j_adj > nBinsY_ || j_adj < 1) {isXBoundary = kTRUE;}
330 
331  // In the corners, do no interpolation. Use entry in bin (i,j)
332  if (isXBoundary == kTRUE && isYBoundary == kTRUE) {
333 
334  value = this->getBinHistValue(i,j);
335 
336  } else if (isXBoundary == kTRUE && isYBoundary == kFALSE) {
337 
338  // Find the adjacent x bin centre
339  Double_t cbinx_adj = static_cast<Double_t>(i_adj-0.5)*binXWidth_ + minX_;
340 
341  Double_t dx0 = TMath::Abs(x - cbinx);
342  Double_t dx1 = TMath::Abs(cbinx_adj - x);
343  Double_t inter_denom = dx0 + dx1;
344 
345  Double_t value1 = this->getBinHistValue(i,j);
346  Double_t value2 = this->getBinHistValue(i_adj,j);
347 
348  value = (value1*dx1 + value2*dx0)/inter_denom;
349 
350  } else if (isYBoundary == kTRUE && isXBoundary == kFALSE) {
351 
352  // Find the adjacent y bin centre
353  Double_t cbiny_adj = static_cast<Double_t>(j_adj-0.5)*binYWidth_ + minY_;
354 
355  Double_t dy0 = TMath::Abs(y - cbiny);
356  Double_t dy1 = TMath::Abs(cbiny_adj - y);
357  Double_t inter_denom = dy0 + dy1;
358 
359  Double_t value1 = this->getBinHistValue(i,j);
360  Double_t value2 = this->getBinHistValue(i,j_adj);
361 
362  value = (value1*dy1 + value2*dy0)/inter_denom;
363 
364  } else {
365 
366  // 2D linear interpolation using inverse distance as weights.
367  // Find the adjacent x and y bin centres
368  Double_t cbinx_adj = static_cast<Double_t>(i_adj-0.5)*binXWidth_ + minX_;
369  Double_t cbiny_adj = static_cast<Double_t>(j_adj-0.5)*binYWidth_ + minY_;
370 
371  Double_t dx0 = TMath::Abs(x - cbinx);
372  Double_t dx1 = TMath::Abs(cbinx_adj - x);
373  Double_t dy0 = TMath::Abs(y - cbiny);
374  Double_t dy1 = TMath::Abs(cbiny_adj - y);
375 
376  Double_t inter_denom = (dx0 + dx1)*(dy0 + dy1);
377 
378  Double_t value1 = this->getBinHistValue(i,j);
379  Double_t value2 = this->getBinHistValue(i_adj,j);
380  Double_t value3 = this->getBinHistValue(i,j_adj);
381  Double_t value4 = this->getBinHistValue(i_adj,j_adj);
382 
383  value = value1*dx1*dy1 + value2*dx0*dy1 + value3*dx1*dy0 + value4*dx0*dy0;
384  value /= inter_denom;
385  }
386 
387  return value;
388 }
389 
391 {
392  // Check that the given abscissa is within the allowed range
393  if (!this->checkRange(abscissas)) {
394  gSystem->Exit(EXIT_FAILURE);
395  }
396 
397  // Calculate the interpolated value
398  Double_t x = abscissas[0];
399  Double_t y = abscissas[1];
400  Double_t value = this->interpolateXY(x,y);
401  this->setUnNormPDFVal(value);
402 
403  // TODO - do we need this? I think probably not.
404  // Have the 1D PDFs calculate their likelihoods as well
405  //xVarPdf_->calcLikelihoodInfo(x);
406  //yVarPdf_->calcLikelihoodInfo(y);
407 }
408 
410 {
411  // Call the base class method
412  this->LauAbsPdf::calcLikelihoodInfo(iEvt);
413 
414  // Have the 1D PDFs retrieve their values as well
417 }
418 
419 Double_t Lau2DHistPdf::getLikelihood( const TString& theVarName ) const
420 {
421  if ( theVarName == xName_ ) {
422  return xVarPdf_->getLikelihood();
423  } else if ( theVarName == yName_ ) {
424  return yVarPdf_->getLikelihood();
425  } else {
426  std::cerr << "ERROR in Lau2DHistPdf::getLikelihood : Unrecognised variable name \"" << theVarName << "\", cannot determine likelihood." << std::endl;
427  return 0.0;
428  }
429 }
430 
432 {
433  // Call the base class method
434  this->LauAbsPdf::cacheInfo(inputData);
435 
436  // Have the 1D PDFs cache their info as well
437  xVarPdf_->cacheInfo(inputData);
438  yVarPdf_->cacheInfo(inputData);
439 }
440 
442 {
443  TRandom* random(LauRandom::randomFun());
444  for (Int_t i(0); i<nBinsX_; i++) {
445  for (Int_t j(0); j<nBinsY_; j++) {
446  Double_t currentContent = hist_->GetBinContent(i+1,j+1);
447  Double_t currentError = hist_->GetBinError(i+1,j+1);
448  Double_t newContent = random->Gaus(currentContent,currentError);
449  if (newContent<0.0) {
450  hist_->SetBinContent(i+1,j+1,0.0);
451  } else {
452  hist_->SetBinContent(i+1,j+1,newContent);
453  }
454  }
455  }
456 }
457 
459 {
460  this->withinGeneration(kTRUE);
461 
462  // Check that the PDF height is up to date
463  this->calcPDFHeight( kinematics );
464 
465  // Generate the value of the abscissa.
466  Bool_t gotAbscissa(kFALSE);
467  if (this->getRandomFun() == 0) {
468  std::cerr << "ERROR in Lau2DHistPdf::generate : Please set the random number generator for this PDF by using the setRandomFun(TRandom*) function." << std::endl;
469  this->withinGeneration(kFALSE);
470  gSystem->Exit(EXIT_FAILURE);
471  }
472 
473  Double_t genPDFVal(0.0);
474  LauAbscissas genAbscissas(2);
475 
476  Double_t PDFheight = this->getMaxHeight()*(1.0+1e-11);
477  while (!gotAbscissa) {
478 
479  genAbscissas[0] = this->getRandomFun()->Rndm()*this->getRange(xName_) + this->getMinAbscissa(xName_);
480  genAbscissas[1] = this->getRandomFun()->Rndm()*this->getRange(yName_) + this->getMinAbscissa(yName_);
481 
482  this->calcLikelihoodInfo(genAbscissas);
483  genPDFVal = this->getUnNormLikelihood();
484 
485  if (this->getRandomFun()->Rndm() <= genPDFVal/PDFheight) {gotAbscissa = kTRUE;}
486 
487  if (genPDFVal > PDFheight) {
488  std::cerr << "WARNING in LauAbsPdf::generate : genPDFVal = " << genPDFVal << " is larger than the specified PDF height " << this->getMaxHeight() << " for (x,y) = (" << genAbscissas[0] << "," << genAbscissas[1] << ")." << std::endl;
489  std::cerr << " : Need to reset height to be larger than " << genPDFVal << " by using the setMaxHeight(Double_t) function and re-run the Monte Carlo generation!" << std::endl;
490  }
491  }
492 
493  LauFitData genData;
494  genData[ xName_ ] = genAbscissas[0];
495  genData[ yName_ ] = genAbscissas[1];
496 
497  this->withinGeneration(kFALSE);
498 
499  return genData;
500 }
501 
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)
Calculate the likelihood (and intermediate info) for a given abscissa.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
virtual void setUnNormPDFVal(Double_t unNormPDFVal)
Set the unnormalised likelihood.
Definition: LauAbsPdf.hh:369
virtual Double_t getMinAbscissa() const
Retrieve the minimum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:117
virtual Bool_t heightUpToDate() const
Check if the maximum height of the PDF is up to date.
Definition: LauAbsPdf.hh:264
Lau1DHistPdf * xVarPdf_
1D PDF for x variable
Double_t interpolateXY(Double_t x, Double_t y) const
Perform the interpolation (unnormalised)
ClassImp(LauAbsCoeffSet)
Bool_t useInterpolation_
Control boolean for using the linear interpolation.
Double_t interpolateXYNorm(Double_t x, Double_t y) const
Perform the interpolation and divide by the normalisation.
void checkNormalisation()
Check the normalisation calculation.
virtual void cacheInfo(const LauFitDataTree &inputData)
Cache information from data.
virtual Double_t getUnNormLikelihood() const
Retrieve the unnormalised likelihood value.
Definition: LauAbsPdf.hh:196
File containing declaration of Lau2DHistPdf class.
TH2 * hist_
The underlying histogram.
Double_t invBinXWidth_
The histogram x-axis inverse bin width.
Class for defining a 1D histogram PDF.
Definition: Lau1DHistPdf.hh:31
virtual void setNorm(Double_t norm)
Set the normalisation factor.
Definition: LauAbsPdf.hh:325
Class for defining a 2D histogram PDF.
Definition: Lau2DHistPdf.hh:33
virtual Bool_t checkRange(const LauAbscissas &abscissas) const
Check that all abscissas are within their allowed ranges.
Definition: LauAbsPdf.cc:213
Double_t minY_
The histogram y-axis minimum.
std::map< TString, Double_t > LauFitData
Type for holding event data.
virtual LauFitData generate(const LauKinematics *kinematics)
Generate an event from the PDF.
virtual TRandom * getRandomFun() const
Retrieve the random function used for MC generation.
Definition: LauAbsPdf.hh:387
void doBinFluctuation()
Fluctuate the histogram bin contents in accordance with their errors.
TString xName_
x variable name
virtual Double_t getMaxHeight() const
Retrieve the maximum height.
Definition: LauAbsPdf.hh:221
virtual ~Lau2DHistPdf()
Destructor.
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)
Calculate the likelihood (and intermediate info) for a given value of the abscissas.
TH1 * yProj_
Projection of histogram y-axis.
Double_t binYWidth_
The histogram y-axis bin width.
Bool_t fluctuateBins_
Control boolean for performing the fluctuation of the histogram bin contents.
Double_t getBinHistValue(Int_t i, Int_t j) const
Get the bin content from the histogram.
TH1 * xProj_
Projection of histogram x-axis.
virtual void setMaxHeight(Double_t maxHeight)
Set the maximum height.
Definition: LauAbsPdf.hh:331
Lau1DHistPdf * yVarPdf_
1D PDF for y variable
Class for defining the fit parameter objects.
Definition: LauParameter.hh:33
virtual void calcPDFHeight(const LauKinematics *kinematics)
Calculate the PDF height.
Double_t binXWidth_
The histogram x-axis bin width.
virtual Bool_t withinGeneration() const
Check whether the generate method is running.
Definition: LauAbsPdf.hh:435
File containing LauRandom namespace.
virtual void calcNorm()
Calculate the normalisation.
Double_t maxX_
The histogram x-axis maximum.
virtual Double_t getLikelihood() const
Retrieve the normalised likelihood value.
Definition: LauAbsPdf.cc:354
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)=0
Calculate the likelihood (and all associated information) given value(s) of the abscissa(s) ...
Double_t invBinYWidth_
The histogram y-axis inverse bin width.
Double_t minX_
The histogram x-axis minimum.
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:41
Class for calculating 3-body kinematic quantities.
TString yName_
y variable name
Double_t maxY_
The histogram y-axis maximum.
virtual void cacheInfo(const LauFitDataTree &inputData)
Cache information from data.
Definition: LauAbsPdf.cc:241
Double_t value() const
The value of the parameter.
virtual Double_t getRange() const
Retrieve the range of the (primary) abscissa.
Definition: LauAbsPdf.hh:129
Lau2DHistPdf(const std::vector< TString > &theVarNames, const TH2 *hist, const LauFitData &minVals, const LauFitData &maxVals, Bool_t useInterpolation=kTRUE, Bool_t fluctuateBins=kFALSE)
Constructor.
Definition: Lau2DHistPdf.cc:33
virtual void setRandomFun(TRandom *randomFun)
Set the random function used for toy MC generation.
Definition: LauAbsPdf.hh:233
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
Double_t rangeX_
The histogram x-axis range.
Class to store the input fit variables.
std::vector< Double_t > LauAbscissas
The type used for containing multiple abscissa values.
Definition: LauAbsPdf.hh:45
Int_t nBinsX_
The number of bins on the x-axis of the histogram.
Int_t nBinsY_
The number of bins on the y-axis of the histogram.
Double_t rangeY_
The histogram y-axis range.
File containing declaration of Lau1DHistPdf class.