laura is hosted by Hepforge, IPPP Durham
Laura++  v3r0
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 void Lau2DHistPdf::calcPDFHeight( const LauKinematics* /*kinematics*/ )
162 {
163  if (this->heightUpToDate()) {
164  return;
165  }
166 
167  // Get the maximum height of the histogram
168  Int_t maxBin = hist_->GetMaximumBin();
169  Double_t height = hist_->GetBinContent(maxBin);
170  this->setMaxHeight(height);
171 }
172 
174 {
175  // Calculate the histogram normalisation.
176 
177  // Loop over the total x and y range to get the total area under x
178  // and y. Just sum the contributions up using 1e-3 increments of
179  // the range in x and y. Multiply the end result by dx and dy.
180 
181  Double_t dx(1e-3*rangeX_), dy(1e-3*rangeY_);
182  Double_t area(0.0);
183 
184  Double_t x(minX_ + dx/2.0);
185  while (x > minX_ && x < maxX_) {
186  Double_t y(minY_ + dy/2.0);
187  while (y > minY_ && y < maxY_) {
188  area += this->interpolateXY(x,y);
189  y += dy;
190  } // y while loop
191  x += dx;
192  } // x while loop
193 
194  Double_t norm = area*dx*dy;
195  this->setNorm(norm);
196 }
197 
199 {
200  Double_t dx(1e-3*rangeX_), dy(1e-3*rangeY_);
201  Double_t area(0.0);
202  Double_t areaNoNorm(0.0);
203 
204  // Preserve the value of a variable we change temporarily
205  Bool_t interpolate = useInterpolation_;
206 
207  Double_t x(minX_ + dx/2.0);
208  while (x > minX_ && x < maxX_) {
209  Double_t y(minY_ + dy/2.0);
210  while (y > minY_ && y < maxY_) {
211  area += this->interpolateXYNorm(x,y);
212  useInterpolation_ = kFALSE;
213  areaNoNorm += this->interpolateXY(x,y);
214  useInterpolation_ = kTRUE;
215  y += dy;
216  } // y while loop
217  x += dx;
218  } // x while loop
219  Double_t norm = area*dx*dy;
220 
221  useInterpolation_ = interpolate; //Restore old value
222 
223  std::cout << "INFO in Lau2DHistPdf::checkNormalisation : Area = " << area << ", dx = " << dx << ", dy = " << dy << ", dx*dy = " << dx*dy << std::endl;
224  std::cout << " : Area with no norm = " << areaNoNorm << "*dx*dy = " << areaNoNorm*dx*dy << std::endl;
225  std::cout << " : The total area of the normalised histogram PDF is " << norm << std::endl;
226 }
227 
228 Double_t Lau2DHistPdf::getBinHistValue(Int_t i, Int_t j) const
229 {
230  // Check that x co-ord is in range [1 , nBinsX_]
231  if ((i < 1) || (i > nBinsX_)) {
232  return 0.0;
233  }
234  // Check that y co-ord is in range [1 , nBinsY_]
235  if ((j < 1) || (j > nBinsY_)) {
236  return 0.0;
237  }
238  Double_t value = hist_->GetBinContent(i,j);
239  // protect against negative values
240  if ( value < 0.0 ) {
241  std::cerr << "WARNING in Lau2DHistPdf::getBinHistValue : Negative bin content set to zero!" << std::endl;
242  value = 0.0;
243  }
244  return value;
245 }
246 
247 Double_t Lau2DHistPdf::interpolateXYNorm(Double_t x, Double_t y) const
248 {
249  // Get the normalised interpolated value.
250  Double_t value = this->interpolateXY(x,y);
251  Double_t norm = this->getNorm();
252  return value/norm;
253 }
254 
255 Double_t Lau2DHistPdf::interpolateXY(Double_t x, Double_t y) const
256 {
257  // This function returns the interpolated value of the histogram
258  // function for the given values of x and y by finding the adjacent
259  // bins and extrapolating using weights based on the inverse
260  // distance of the point from the adajcent bin centres.
261 
262  // Find the 2D histogram bin for x and y
263  Int_t i = static_cast<Int_t>((x - minX_)*invBinXWidth_)+1;
264  Int_t j = static_cast<Int_t>((y - minY_)*invBinYWidth_)+1;
265 
266  if (i > nBinsX_) {i = nBinsX_;}
267  if (j > nBinsY_) {j = nBinsY_;}
268 
269  // Ask whether we want to do the interpolation, since this method is
270  // not reliable for low statistics histograms.
271  if (useInterpolation_ == kFALSE) {
272  return this->getBinHistValue(i,j);
273  }
274 
275  // Find the bin centres (actual co-ordinate positions, not histogram indices)
276  Double_t cbinx = static_cast<Double_t>(i-0.5)*binXWidth_ + minX_;
277  Double_t cbiny = static_cast<Double_t>(j-0.5)*binYWidth_ + minY_;
278 
279  // Find the adjacent bins
280  Double_t deltax = x - cbinx;
281  Double_t deltay = y - cbiny;
282 
283  Int_t i_adj(0), j_adj(0);
284  if (deltax > 0.0) {
285  i_adj = i + 1;
286  } else {
287  i_adj = i - 1;
288  }
289  if (deltay > 0.0) {
290  j_adj = j + 1;
291  } else {
292  j_adj = j - 1;
293  }
294 
295  Bool_t isXBoundary(kFALSE), isYBoundary(kFALSE);
296 
297  Double_t value(0.0);
298 
299  if (i_adj > nBinsX_ || i_adj < 1) {isYBoundary = kTRUE;}
300  if (j_adj > nBinsY_ || j_adj < 1) {isXBoundary = kTRUE;}
301 
302  // In the corners, do no interpolation. Use entry in bin (i,j)
303  if (isXBoundary == kTRUE && isYBoundary == kTRUE) {
304 
305  value = this->getBinHistValue(i,j);
306 
307  } else if (isXBoundary == kTRUE && isYBoundary == kFALSE) {
308 
309  // Find the adjacent x bin centre
310  Double_t cbinx_adj = static_cast<Double_t>(i_adj-0.5)*binXWidth_ + minX_;
311 
312  Double_t dx0 = TMath::Abs(x - cbinx);
313  Double_t dx1 = TMath::Abs(cbinx_adj - x);
314  Double_t inter_denom = dx0 + dx1;
315 
316  Double_t value1 = this->getBinHistValue(i,j);
317  Double_t value2 = this->getBinHistValue(i_adj,j);
318 
319  value = (value1*dx1 + value2*dx0)/inter_denom;
320 
321  } else if (isYBoundary == kTRUE && isXBoundary == kFALSE) {
322 
323  // Find the adjacent y bin centre
324  Double_t cbiny_adj = static_cast<Double_t>(j_adj-0.5)*binYWidth_ + minY_;
325 
326  Double_t dy0 = TMath::Abs(y - cbiny);
327  Double_t dy1 = TMath::Abs(cbiny_adj - y);
328  Double_t inter_denom = dy0 + dy1;
329 
330  Double_t value1 = this->getBinHistValue(i,j);
331  Double_t value2 = this->getBinHistValue(i,j_adj);
332 
333  value = (value1*dy1 + value2*dy0)/inter_denom;
334 
335  } else {
336 
337  // 2D linear interpolation using inverse distance as weights.
338  // Find the adjacent x and y bin centres
339  Double_t cbinx_adj = static_cast<Double_t>(i_adj-0.5)*binXWidth_ + minX_;
340  Double_t cbiny_adj = static_cast<Double_t>(j_adj-0.5)*binYWidth_ + minY_;
341 
342  Double_t dx0 = TMath::Abs(x - cbinx);
343  Double_t dx1 = TMath::Abs(cbinx_adj - x);
344  Double_t dy0 = TMath::Abs(y - cbiny);
345  Double_t dy1 = TMath::Abs(cbiny_adj - y);
346 
347  Double_t inter_denom = (dx0 + dx1)*(dy0 + dy1);
348 
349  Double_t value1 = this->getBinHistValue(i,j);
350  Double_t value2 = this->getBinHistValue(i_adj,j);
351  Double_t value3 = this->getBinHistValue(i,j_adj);
352  Double_t value4 = this->getBinHistValue(i_adj,j_adj);
353 
354  value = value1*dx1*dy1 + value2*dx0*dy1 + value3*dx1*dy0 + value4*dx0*dy0;
355  value /= inter_denom;
356  }
357 
358  return value;
359 }
360 
362 {
363  // Check that the given abscissa is within the allowed range
364  if (!this->checkRange(abscissas)) {
365  gSystem->Exit(EXIT_FAILURE);
366  }
367 
368  // Calculate the interpolated value
369  Double_t x = abscissas[0];
370  Double_t y = abscissas[1];
371  Double_t value = this->interpolateXY(x,y);
372  this->setUnNormPDFVal(value);
373 
374  // TODO - do we need this? I think probably not.
375  // Have the 1D PDFs calculate their likelihoods as well
376  //xVarPdf_->calcLikelihoodInfo(x);
377  //yVarPdf_->calcLikelihoodInfo(y);
378 }
379 
381 {
382  // Call the base class method
383  this->LauAbsPdf::calcLikelihoodInfo(iEvt);
384 
385  // Have the 1D PDFs retrieve their values as well
388 }
389 
390 Double_t Lau2DHistPdf::getLikelihood( const TString& theVarName ) const
391 {
392  if ( theVarName == xName_ ) {
393  return xVarPdf_->getLikelihood();
394  } else if ( theVarName == yName_ ) {
395  return yVarPdf_->getLikelihood();
396  } else {
397  std::cerr << "ERROR in Lau2DHistPdf::getLikelihood : Unrecognised variable name \"" << theVarName << "\", cannot determine likelihood." << std::endl;
398  return 0.0;
399  }
400 }
401 
403 {
404  // Call the base class method
405  this->LauAbsPdf::cacheInfo(inputData);
406 
407  // Have the 1D PDFs cache their info as well
408  xVarPdf_->cacheInfo(inputData);
409  yVarPdf_->cacheInfo(inputData);
410 }
411 
413 {
414  TRandom* random(LauRandom::randomFun());
415  for (Int_t i(0); i<nBinsX_; i++) {
416  for (Int_t j(0); j<nBinsY_; j++) {
417  Double_t currentContent = hist_->GetBinContent(i+1,j+1);
418  Double_t currentError = hist_->GetBinError(i+1,j+1);
419  Double_t newContent = random->Gaus(currentContent,currentError);
420  if (newContent<0.0) {
421  hist_->SetBinContent(i+1,j+1,0.0);
422  } else {
423  hist_->SetBinContent(i+1,j+1,newContent);
424  }
425  }
426  }
427 }
428 
430 {
431  this->withinGeneration(kTRUE);
432 
433  // Check that the PDF height is up to date
434  this->calcPDFHeight( kinematics );
435 
436  // Generate the value of the abscissa.
437  Bool_t gotAbscissa(kFALSE);
438  if (this->getRandomFun() == 0) {
439  std::cerr << "ERROR in Lau2DHistPdf::generate : Please set the random number generator for this PDF by using the setRandomFun(TRandom*) function." << std::endl;
440  this->withinGeneration(kFALSE);
441  gSystem->Exit(EXIT_FAILURE);
442  }
443 
444  Double_t genPDFVal(0.0);
445  LauAbscissas genAbscissas(2);
446 
447  Double_t PDFheight = this->getMaxHeight()*(1.0+1e-11);
448  while (!gotAbscissa) {
449 
450  genAbscissas[0] = this->getRandomFun()->Rndm()*this->getRange(xName_) + this->getMinAbscissa(xName_);
451  genAbscissas[1] = this->getRandomFun()->Rndm()*this->getRange(yName_) + this->getMinAbscissa(yName_);
452 
453  this->calcLikelihoodInfo(genAbscissas);
454  genPDFVal = this->getUnNormLikelihood();
455 
456  if (this->getRandomFun()->Rndm() <= genPDFVal/PDFheight) {gotAbscissa = kTRUE;}
457 
458  if (genPDFVal > PDFheight) {
459  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;
460  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;
461  }
462  }
463 
464  LauFitData genData;
465  genData[ xName_ ] = genAbscissas[0];
466  genData[ yName_ ] = genAbscissas[1];
467 
468  this->withinGeneration(kFALSE);
469 
470  return genData;
471 }
472 
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.
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:34
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
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.