ROOT logo
/*****************************************************************************
 * Project: RooFit                                                           *
 * Package: RooFitModels                                                     *
 * @(#)root/roofit:$Id$
 * Authors:                                                                  *
 *   AB, Adrian Bevan, Liverpool University, bevan@slac.stanford.edu         *
 *                                                                           *
 * Copyright (c) 2000-2005, Regents of the University of California,         *
 *                          Liverpool University,                            *
 *                          and Stanford University. All rights reserved.    *
 *                                                                           *
 * Redistribution and use in source and binary forms,                        *
 * with or without modification, are permitted according to the terms        *
 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
 *****************************************************************************/

//////////////////////////////////////////////////////////////////////////////
//
// BEGIN_HTML
// Two-dimensional kernel estimation p.d.f. 
// 
// <b>This function has been superceded by the more general RooNDKeysPdf</b>
// 
// END_HTML
//

#include "RooFit.h"

#include "Roo2DKeysPdf.h"
#include "Roo2DKeysPdf.h"
#include "RooRealVar.h"
#include "TTree.h"
#include "TH2.h"
#include "TFile.h"
#include "TBranch.h"
#include "TMath.h"

//#include <math.h>

using namespace std;

ClassImp(Roo2DKeysPdf)


//_____________________________________________________________________________
Roo2DKeysPdf::Roo2DKeysPdf(const char *name, const char *title,
                       RooAbsReal& xx, RooAbsReal & yy, RooDataSet& data,  TString options, Double_t widthScaleFactor):
  RooAbsPdf(name,title),
  x("x", "x dimension",this, xx),
  y("y", "y dimension",this, yy)
{
  setWidthScaleFactor(widthScaleFactor);
  loadDataSet(data, options);
}


//_____________________________________________________________________________
Roo2DKeysPdf::Roo2DKeysPdf(const Roo2DKeysPdf & other, const char* name) :
  RooAbsPdf(other,name),
  x("x", this, other.x),
  y("y", this, other.y)
{
  if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2DKeysPdf copy ctor" << endl; }

  _xMean   = other._xMean;
  _xSigma  = other._xSigma;
  _yMean   = other._yMean;
  _ySigma  = other._ySigma;
  _n       = other._n;

  _BandWidthType    = other._BandWidthType;
  _MirrorAtBoundary = other._MirrorAtBoundary;
  _widthScaleFactor = other._widthScaleFactor;

  _2pi     = other._2pi;
  _sqrt2pi = other._sqrt2pi;
  _nEvents = other._nEvents;
  _n16     = other._n16;
  _debug   = other._debug;
  _verbosedebug   = other._verbosedebug;
  _vverbosedebug  = other._vverbosedebug;

  _lox       = other._lox;
  _hix       = other._hix;
  _loy       = other._loy;
  _hiy       = other._hiy;
  _xoffset   = other._xoffset;
  _yoffset   = other._yoffset;

  _x  = new Double_t[_nEvents];
  _y  = new Double_t[_nEvents];
  _hx = new Double_t[_nEvents];
  _hy = new Double_t[_nEvents];

  //copy the data and bandwidths
  for(Int_t iEvt = 0; iEvt< _nEvents; iEvt++)
  {
    _x[iEvt]  = other._x[iEvt];
    _y[iEvt]  = other._y[iEvt];
    _hx[iEvt] = other._hx[iEvt];
    _hy[iEvt] = other._hy[iEvt];
  }
}


//_____________________________________________________________________________
Roo2DKeysPdf::~Roo2DKeysPdf() {
  if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2KeysPdf dtor" << endl; }
    delete[] _x;
    delete[] _hx;
    delete[] _y;
    delete[] _hy;
}

/////////////////////////////////////////////////////////////////////////////
// Load a new data set into the class instance.  If the calculation fails, //
//    return 1                                                             //
//    return 0 indicates a success                                         //
/////////////////////////////////////////////////////////////////////////////

//_____________________________________________________________________________
Int_t Roo2DKeysPdf::loadDataSet(RooDataSet& data, TString options)
{
  if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet" << endl; }

  setOptions(options);

  if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet(RooDataSet& data, TString options)" << endl; }

  _2pi       = 2.0*TMath::Pi() ;   //use pi from math.h
  _sqrt2pi   = sqrt(_2pi);
  _nEvents   = (Int_t)data.numEntries();
  if(_nEvents == 0) 
  {
    cout << "ERROR:  Roo2DKeysPdf::loadDataSet The input data set is empty.  Unable to begin generating the PDF" << endl;
    return 1;
  }
  _n16       =  TMath::Power(_nEvents, -0.166666666); // = (4/[n(dim(R) + 2)])^1/(dim(R)+4); dim(R) = 2

  _lox       = x.min();
  _hix       = x.max();
  _loy       = y.min();
  _hiy       = y.max();

  _x  = new Double_t[_nEvents];
  _y  = new Double_t[_nEvents];
  _hx = new Double_t[_nEvents];
  _hy = new Double_t[_nEvents];

  Double_t x0 = 0.0;
  Double_t x1 = 0.0;
  Double_t x_2 = 0.0;
  Double_t y0 = 0.0;
  Double_t y1 = 0.0;
  Double_t y_2 = 0.0;

  //check that the data contain the variable we are interested in  
  Int_t bad = 0;
  const RooAbsReal & xx = x.arg();
  const RooAbsReal & yy = y.arg();
  if(! (RooRealVar*)( (RooArgSet*)data.get(0) )->find( xx.GetName() ) )
  {
    cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<xx.GetName()<<" not in the data set" <<endl;
    bad = 1;
  }
  if(! (RooRealVar*)( (RooArgSet*)data.get(0) )->find( yy.GetName() ) )
  {
    cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<yy.GetName()<<" not in the data set" << endl;
    bad = 1;
  }
  if(bad)
  {
    cout << "Roo2DKeysPdf::Roo2DKeysPdf Unable to initilize object; incompatible RooDataSet doesn't contain"<<endl;
    cout << "                           all of the RooAbsReal arguments"<<endl;
    return 1;
  }

  //copy the data into local arrays
  const RooArgSet * values = data.get();
  const RooRealVar* X = ((RooRealVar*)(values->find(xx.GetName())) ) ;
  const RooRealVar* Y = ((RooRealVar*)(values->find(yy.GetName())) ) ;

  for (Int_t j=0;j<_nEvents;++j) 
  {
    data.get(j) ;

    _x[j] = X->getVal() ;
    _y[j] = Y->getVal() ;

    x0+=1; x1+=_x[j]; x_2+=_x[j]*_x[j];
    y0+=1; y1+=_y[j]; y_2+=_y[j]*_y[j];
  }

  //==========================================//
  //calculate the mean and sigma for the data //
  //==========================================//
  if(_nEvents == 0) 
  {
    cout << "Roo2DKeysPdf::Roo2DKeysPdf Empty data set was used; can't generate a PDF"<<endl;
  }

  _xMean  = x1/x0;
  _xSigma = sqrt(x_2/_nEvents-_xMean*_xMean);
  
  _yMean  = y1/y0;
  _ySigma = sqrt(y_2/_nEvents-_yMean*_yMean);

  _n = Double_t(1)/(_2pi*_nEvents*_xSigma*_ySigma);

  //calculate the PDF
  return calculateBandWidth(_BandWidthType);
}


//_____________________________________________________________________________
void Roo2DKeysPdf::setOptions(TString options)
{
  if(_verbosedebug) { cout << "Roo2DKeysPdf::setOptions" << endl; }

  options.ToLower(); 
  if( options.Contains("a") )   _BandWidthType    = 0;
  else                          _BandWidthType    = 1;
  if( options.Contains("n") )   _BandWidthType    = 1;
  else                          _BandWidthType    = 0;
  if( options.Contains("m") )   _MirrorAtBoundary = 1;
  else                          _MirrorAtBoundary = 0;
  if( options.Contains("d") )   _debug            = 1;
  else                          _debug            = 0;
  if( options.Contains("v") )   { _debug         = 1; _verbosedebug = 1; }
  else                            _verbosedebug  = 0;
  if( options.Contains("vv") )  { _vverbosedebug = 1; }
  else                            _vverbosedebug = 0;

  if( _debug )
  {
    cout << "Roo2DKeysPdf::setOptions(TString options)    options = "<< options << endl;
    cout << "\t_BandWidthType    = " << _BandWidthType    << endl;
    cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
    cout << "\t_debug            = " << _debug            << endl;
    cout << "\t_verbosedebug     = " << _verbosedebug     << endl;
    cout << "\t_vverbosedebug    = " << _vverbosedebug    << endl;
  }
}


//_____________________________________________________________________________
void Roo2DKeysPdf::getOptions(void) const
{
  cout << "Roo2DKeysPdf::getOptions(void)" << endl;
  cout << "\t_BandWidthType                           = " << _BandWidthType    << endl;
  cout << "\t_MirrorAtBoundary                        = " << _MirrorAtBoundary << endl;
  cout << "\t_debug                                   = " << _debug            << endl;
  cout << "\t_verbosedebug                            = " << _verbosedebug     << endl;
  cout << "\t_vverbosedebug                           = " << _vverbosedebug    << endl;
}

//=====================================================//
// calculate the kernal bandwith for x & y             //
// & Calculate the probability look up table _p[i][j]  //
//=====================================================//

//_____________________________________________________________________________
Int_t Roo2DKeysPdf::calculateBandWidth(Int_t kernel)
{
  if(_verbosedebug) { cout << "Roo2DKeysPdf::calculateBandWidth(Int_t kernel)" << endl; }
  if(kernel != -999)
  {
    _BandWidthType = kernel;
  }

  Double_t h = 0.0;

  Double_t sigSum       = _xSigma*_xSigma + _ySigma*_ySigma;
  Double_t sqrtSum      = sqrt( sigSum );
  Double_t sigProd      = _ySigma*_xSigma;
  if(sigProd != 0.0)  h = _n16*sqrt( sigSum/sigProd );
  if(sqrtSum == 0)
  { 
    cout << "Roo2DKeysPdf::calculateBandWidth The sqr(variance sum) == 0.0. " << " Your dataset represents a delta function."<<endl;
    return 1;
  }

  Double_t hXSigma = h * _xSigma; 
  Double_t hYSigma = h * _ySigma; 
  Double_t xhmin   = hXSigma * sqrt(2.)/10;  //smallest anticipated bandwidth
  Double_t yhmin   = hYSigma * sqrt(2.)/10;

  //////////////////////////////////////
  //calculate bandwidths from the data//
  //////////////////////////////////////
  if(_BandWidthType == 1)  //calculate a trivial bandwith
  {
    cout << "Roo2DKeysPdf::calculateBandWidth Using a normal bandwith (same for a given dimension) based on"<<endl;
    cout << "                                 h_j = n^{-1/6}*sigma_j for the j^th dimension and n events * "<<_widthScaleFactor<<endl;
    Double_t hxGaussian = _n16 * _xSigma * _widthScaleFactor;
    Double_t hyGaussian = _n16 * _ySigma * _widthScaleFactor;
    for(Int_t j=0;j<_nEvents;++j) 
    {
      _hx[j] = hxGaussian;
      _hy[j] = hyGaussian;
      if(_hx[j]<xhmin) _hx[j] = xhmin;
      if(_hy[j]<yhmin) _hy[j] = yhmin;
     }
  }
  else //use an adaptive bandwith to reduce the dependance on global data distribution
  {
    cout << "Roo2DKeysPdf::calculateBandWidth Using an adaptive bandwith (in general different for all events) [default]"<<endl;
    cout << "                                 scaled by a factor of "<<_widthScaleFactor<<endl;
    Double_t xnorm   = h * TMath::Power(_xSigma/sqrtSum, 1.5) * _widthScaleFactor;
    Double_t ynorm   = h * TMath::Power(_ySigma/sqrtSum, 1.5) * _widthScaleFactor;
    for(Int_t j=0;j<_nEvents;++j) 
    {
      Double_t f_ti =  TMath::Power( g(_x[j], _x, hXSigma, _y[j], _y, hYSigma), -0.25 ) ;
      _hx[j] = xnorm * f_ti;
      _hy[j] = ynorm * f_ti;
      if(_hx[j]<xhmin) _hx[j] = xhmin;
      if(_hy[j]<yhmin) _hy[j] = yhmin;
    }
  }

  return 0;
}

//=======================================================================================//
// evaluate the kernal estimation for x,y, interpolating between the points if necessary //
//=======================================================================================//

//_____________________________________________________________________________
Double_t Roo2DKeysPdf::evaluate() const
{
  // use the cacheing intrinsic in RFC to bypass the grid and remove
  // the grid and extrapolation approximation in the kernel estimation method 
  //implementation - cheers Wouter :)
  if(_vverbosedebug) { cout << "Roo2DKeysPdf::evaluate()" << endl; }
  return evaluateFull(x,y);
}

/////////////////////////////////////////////////////////
// Evaluate the sum of the product of the 2D kernels   //
// for use in calculating the fixed kernel estimate, f //
// given the bandwiths _hx[j] and _hy[j]               //
/////////////////////////////////////////////////////////
// _n is calculated once in the constructor

//_____________________________________________________________________________
Double_t Roo2DKeysPdf::evaluateFull(Double_t thisX, Double_t thisY) const
{
  if( _vverbosedebug ) { cout << "Roo2DKeysPdf::evaluateFull()" << endl; }

  Double_t f=0.0;

  Double_t rx2, ry2, zx, zy;
  if( _MirrorAtBoundary )
  {
    for (Int_t j = 0; j < _nEvents; ++j) 
    {
      rx2 = 0.0; ry2 = 0.0; zx = 0.0; zy = 0.0;
      if(_hx[j] != 0.0) rx2 = (thisX - _x[j])/_hx[j];
      if(_hy[j] != 0.0) ry2 = (thisY - _y[j])/_hy[j];

      if(_hx[j] != 0.0) zx = exp(-0.5*rx2*rx2)/_hx[j];
      if(_hy[j] != 0.0) zy = exp(-0.5*ry2*ry2)/_hy[j];

      zx += highBoundaryCorrection(thisX, _hx[j], x.max(), _x[j])
 	 +   lowBoundaryCorrection(thisX, _hx[j], x.min(), _x[j]);
      zy += highBoundaryCorrection(thisY, _hy[j], y.max(), _y[j])
 	 +   lowBoundaryCorrection(thisY, _hy[j], y.min(), _y[j]);
      f += zy * zx;
      //      f += _n * zy * zx; // ooops this is a normalisation factor :(
    }
  }
  else
  {
    for (Int_t j = 0; j < _nEvents; ++j) 
    {
      rx2 = 0.0; ry2 = 0.0; zx = 0.0; zy = 0.0;
      if(_hx[j] != 0.0) rx2 = (thisX - _x[j])/_hx[j];
      if(_hy[j] != 0.0) ry2 = (thisY - _y[j])/_hy[j];

      if(_hx[j] != 0.0) zx = exp(-0.5*rx2*rx2)/_hx[j];
      if(_hy[j] != 0.0) zy = exp(-0.5*ry2*ry2)/_hy[j];
      f += zy * zx;
      //      f += _n * zy * zx; // ooops this is a normalisation factor  :(
    }
  }
  return f;
}

// Apply the mirror at boundary correction to a dimension given the space position to evaluate 
// at (thisVar), the bandwidth at this position (thisH), the boundary (high/low) and the
// value of the data kernal that this correction is being applied to  tVar (i.e. the _x[ix] etc.)

//_____________________________________________________________________________
Double_t Roo2DKeysPdf::highBoundaryCorrection(Double_t thisVar, Double_t thisH, Double_t high, Double_t tVar) const
{
  if(_vverbosedebug) { cout << "Roo2DKeysPdf::highBoundaryCorrection" << endl; }

  if(thisH == 0.0) return 0.0;
  Double_t correction = (thisVar + tVar - 2.0* high )/thisH;
  return exp(-0.5*correction*correction)/thisH;
}


//_____________________________________________________________________________
Double_t Roo2DKeysPdf::lowBoundaryCorrection(Double_t thisVar, Double_t thisH, Double_t low, Double_t tVar) const
{
  if(_vverbosedebug) { cout << "Roo2DKeysPdf::lowBoundaryCorrection" << endl; }

  if(thisH == 0.0) return 0.0;
  Double_t correction = (thisVar + tVar - 2.0* low )/thisH;
  return exp(-0.5*correction*correction)/thisH;
}

//==========================================================================================//
// calculate f(t_i) for the bandwidths                                                      //
//                                                                                          //
// g = 1/(Nevt * sigma_j * sqrt2pi)*sum_{all evts}{prod d K[ exp{-(xd - ti)/sigma_jd^2} ]}  //
//                                                                                          //
//==========================================================================================//

//_____________________________________________________________________________
Double_t Roo2DKeysPdf::g(Double_t varMean1, Double_t * _var1, Double_t sigma1, Double_t varMean2, Double_t * _var2, Double_t sigma2) const
{
  if((_nEvents == 0.0) || (sigma1 == 0.0) || (sigma2 == 0)) return 0.0;

  Double_t c1 = -1.0/(2.0*sigma1*sigma1);
  Double_t c2 = -1.0/(2.0*sigma2*sigma2);
  Double_t d  = 4.0*c1*c2  /(_sqrt2pi*_nEvents);
  Double_t z  = 0.0;

  for (Int_t i = 0; i < _nEvents; ++i) 
  {
    Double_t r1 =  _var1[i] - varMean1; 
    Double_t r2 =  _var2[i] - varMean2; 
    z          += exp( c1 * r1*r1 ) * exp( c2 * r2*r2 );
  }
  z = z*d;
  return z;
}


//_____________________________________________________________________________
Int_t Roo2DKeysPdf::getBandWidthType() const
{
  if(_BandWidthType == 1)  cout << "The Bandwidth Type selected is Trivial" << endl;
  else                     cout << "The Bandwidth Type selected is Adaptive" << endl;

  return _BandWidthType;
}


//_____________________________________________________________________________
Double_t Roo2DKeysPdf::getMean(const char * axis) const
{
  if(!strcmp(axis,x.GetName()) || !strcmp(axis,"x") || !strcmp(axis,"X"))      return _xMean;
  else if(!strcmp(axis,y.GetName()) || !strcmp(axis,"y") || !strcmp(axis,"Y")) return _yMean;
  else 
  {
    cout << "Roo2DKeysPdf::getMean unknown axis "<<axis<<endl;
  }
  return 0.0;
}


//_____________________________________________________________________________
Double_t Roo2DKeysPdf::getSigma(const char * axis) const
{
  if(!strcmp(axis,x.GetName()) || !strcmp(axis,"x") || !strcmp(axis,"X"))      return _xSigma;
  else if(!strcmp(axis,y.GetName()) || !strcmp(axis,"y") || !strcmp(axis,"Y")) return _ySigma;
  else 
  {
    cout << "Roo2DKeysPdf::getSigma unknown axis "<<axis<<endl;
  }
  return 0.0;
}



//_____________________________________________________________________________
void Roo2DKeysPdf::writeToFile(char * outputFile, const char * name) const
{
  TString histName = name;
  histName        += "_hist";
  TString nName    = name;
  nName           += "_Ntuple";
  writeHistToFile( outputFile,    histName);
  writeNTupleToFile( outputFile,  nName);
}

// plot the PDf as a histogram and save to file
// so that it can be loaded in as a Roo2DHist Pdf in the future to 
// save on calculation time

//_____________________________________________________________________________
void Roo2DKeysPdf::writeHistToFile(char * outputFile, const char * histName) const
{
  TFile * file = 0;
  cout << "Roo2DKeysPdf::writeHistToFile This member function is temporarily disabled" <<endl;
  //make sure that any existing file is not over written
  file = new TFile(outputFile, "UPDATE"); 
  if (!file)
  {
    cout << "Roo2DKeysPdf::writeHistToFile unable to open file "<< outputFile <<endl;
    return;
  }


  const RooAbsReal & xx = x.arg();
  const RooAbsReal & yy = y.arg();
  RooArgSet values( RooArgList( xx, yy ));
  RooRealVar * xArg = ((RooRealVar*)(values.find(xx.GetName())) ) ;
  RooRealVar * yArg = ((RooRealVar*)(values.find(yy.GetName())) ) ;

  TH2F * hist = (TH2F*)xArg->createHistogram("hist", *yArg);
  hist = (TH2F*)this->fillHistogram(hist, RooArgList(*xArg, *yArg) ); 
  hist->SetName(histName);

  file->Write();
  file->Close();
}

// save the data and calculated bandwidths to file
// as a record of what produced the PDF and to give a reduced
// data set in order to facilitate re-calculation in the future

//_____________________________________________________________________________
void Roo2DKeysPdf::writeNTupleToFile(char * outputFile, const char * name) const
{
  TFile * file = 0;

  //make sure that any existing file is not over written
  file = new TFile(outputFile, "UPDATE"); 
  if (!file)
  {
    cout << "Roo2DKeysPdf::writeNTupleToFile unable to open file "<< outputFile <<endl;
    return;
  }
  RooAbsReal & xArg = (RooAbsReal&)x.arg();
  RooAbsReal & yArg = (RooAbsReal&)y.arg();

  Double_t theX, theY, hx/*, hy*/;
  TString label = name;
  label += " the source data for 2D Keys PDF";
  TTree * _theTree =  new TTree(name, label);
  if(!_theTree) { cout << "Unable to get a TTree for output" << endl; return; }
  _theTree->SetAutoSave(1000000000);  // autosave when 1 Gbyte written

  //name the TBranches the same as the RooAbsReal's
  const char * xname = xArg.GetName();
  const char * yname = yArg.GetName();
  if (!strcmp(xname,"")) xname = "x";
  if (!strcmp(yname,"")) yname = "y";

  _theTree->Branch(xname, &theX, " x/D");
  _theTree->Branch(yname, &theY, " y/D");
  _theTree->Branch("hx",  &hx,   " hx/D");
  _theTree->Branch("hy",  &hx,  " hy/D");

  for(Int_t iEvt = 0; iEvt < _nEvents; iEvt++)
  {
    theX = _x[iEvt];
    theY = _y[iEvt];
    hx   = _hx[iEvt];
    hx   = _hy[iEvt];
    _theTree->Fill();
  }
  file->Write();
  file->Close();
}

/////////////////////////////////////////////////////
// print out _p[_nPoints][_nPoints] indicating the //
// domain limits                                   //
/////////////////////////////////////////////////////

//_____________________________________________________________________________
void Roo2DKeysPdf::PrintInfo(ostream & out) const
{
  out << "Roo2DKeysPDF instance domain information:"<<endl;
  out << "\tX_min          = " << _lox <<endl;
  out << "\tX_max          = " << _hix <<endl;
  out << "\tY_min          = " << _loy <<endl;
  out << "\tY_max          = " << _hiy <<endl;

  out << "Data information:" << endl;
  out << "\t<x>             = " << _xMean <<endl;
  out << "\tsigma(x)       = " << _xSigma <<endl;
  out << "\t<y>             = " << _yMean <<endl;
  out << "\tsigma(y)       = " << _ySigma <<endl;

  out << "END of info for Roo2DKeys pdf instance"<< endl;
}


 Roo2DKeysPdf.cxx:1
 Roo2DKeysPdf.cxx:2
 Roo2DKeysPdf.cxx:3
 Roo2DKeysPdf.cxx:4
 Roo2DKeysPdf.cxx:5
 Roo2DKeysPdf.cxx:6
 Roo2DKeysPdf.cxx:7
 Roo2DKeysPdf.cxx:8
 Roo2DKeysPdf.cxx:9
 Roo2DKeysPdf.cxx:10
 Roo2DKeysPdf.cxx:11
 Roo2DKeysPdf.cxx:12
 Roo2DKeysPdf.cxx:13
 Roo2DKeysPdf.cxx:14
 Roo2DKeysPdf.cxx:15
 Roo2DKeysPdf.cxx:16
 Roo2DKeysPdf.cxx:17
 Roo2DKeysPdf.cxx:18
 Roo2DKeysPdf.cxx:19
 Roo2DKeysPdf.cxx:20
 Roo2DKeysPdf.cxx:21
 Roo2DKeysPdf.cxx:22
 Roo2DKeysPdf.cxx:23
 Roo2DKeysPdf.cxx:24
 Roo2DKeysPdf.cxx:25
 Roo2DKeysPdf.cxx:26
 Roo2DKeysPdf.cxx:27
 Roo2DKeysPdf.cxx:28
 Roo2DKeysPdf.cxx:29
 Roo2DKeysPdf.cxx:30
 Roo2DKeysPdf.cxx:31
 Roo2DKeysPdf.cxx:32
 Roo2DKeysPdf.cxx:33
 Roo2DKeysPdf.cxx:34
 Roo2DKeysPdf.cxx:35
 Roo2DKeysPdf.cxx:36
 Roo2DKeysPdf.cxx:37
 Roo2DKeysPdf.cxx:38
 Roo2DKeysPdf.cxx:39
 Roo2DKeysPdf.cxx:40
 Roo2DKeysPdf.cxx:41
 Roo2DKeysPdf.cxx:42
 Roo2DKeysPdf.cxx:43
 Roo2DKeysPdf.cxx:44
 Roo2DKeysPdf.cxx:45
 Roo2DKeysPdf.cxx:46
 Roo2DKeysPdf.cxx:47
 Roo2DKeysPdf.cxx:48
 Roo2DKeysPdf.cxx:49
 Roo2DKeysPdf.cxx:50
 Roo2DKeysPdf.cxx:51
 Roo2DKeysPdf.cxx:52
 Roo2DKeysPdf.cxx:53
 Roo2DKeysPdf.cxx:54
 Roo2DKeysPdf.cxx:55
 Roo2DKeysPdf.cxx:56
 Roo2DKeysPdf.cxx:57
 Roo2DKeysPdf.cxx:58
 Roo2DKeysPdf.cxx:59
 Roo2DKeysPdf.cxx:60
 Roo2DKeysPdf.cxx:61
 Roo2DKeysPdf.cxx:62
 Roo2DKeysPdf.cxx:63
 Roo2DKeysPdf.cxx:64
 Roo2DKeysPdf.cxx:65
 Roo2DKeysPdf.cxx:66
 Roo2DKeysPdf.cxx:67
 Roo2DKeysPdf.cxx:68
 Roo2DKeysPdf.cxx:69
 Roo2DKeysPdf.cxx:70
 Roo2DKeysPdf.cxx:71
 Roo2DKeysPdf.cxx:72
 Roo2DKeysPdf.cxx:73
 Roo2DKeysPdf.cxx:74
 Roo2DKeysPdf.cxx:75
 Roo2DKeysPdf.cxx:76
 Roo2DKeysPdf.cxx:77
 Roo2DKeysPdf.cxx:78
 Roo2DKeysPdf.cxx:79
 Roo2DKeysPdf.cxx:80
 Roo2DKeysPdf.cxx:81
 Roo2DKeysPdf.cxx:82
 Roo2DKeysPdf.cxx:83
 Roo2DKeysPdf.cxx:84
 Roo2DKeysPdf.cxx:85
 Roo2DKeysPdf.cxx:86
 Roo2DKeysPdf.cxx:87
 Roo2DKeysPdf.cxx:88
 Roo2DKeysPdf.cxx:89
 Roo2DKeysPdf.cxx:90
 Roo2DKeysPdf.cxx:91
 Roo2DKeysPdf.cxx:92
 Roo2DKeysPdf.cxx:93
 Roo2DKeysPdf.cxx:94
 Roo2DKeysPdf.cxx:95
 Roo2DKeysPdf.cxx:96
 Roo2DKeysPdf.cxx:97
 Roo2DKeysPdf.cxx:98
 Roo2DKeysPdf.cxx:99
 Roo2DKeysPdf.cxx:100
 Roo2DKeysPdf.cxx:101
 Roo2DKeysPdf.cxx:102
 Roo2DKeysPdf.cxx:103
 Roo2DKeysPdf.cxx:104
 Roo2DKeysPdf.cxx:105
 Roo2DKeysPdf.cxx:106
 Roo2DKeysPdf.cxx:107
 Roo2DKeysPdf.cxx:108
 Roo2DKeysPdf.cxx:109
 Roo2DKeysPdf.cxx:110
 Roo2DKeysPdf.cxx:111
 Roo2DKeysPdf.cxx:112
 Roo2DKeysPdf.cxx:113
 Roo2DKeysPdf.cxx:114
 Roo2DKeysPdf.cxx:115
 Roo2DKeysPdf.cxx:116
 Roo2DKeysPdf.cxx:117
 Roo2DKeysPdf.cxx:118
 Roo2DKeysPdf.cxx:119
 Roo2DKeysPdf.cxx:120
 Roo2DKeysPdf.cxx:121
 Roo2DKeysPdf.cxx:122
 Roo2DKeysPdf.cxx:123
 Roo2DKeysPdf.cxx:124
 Roo2DKeysPdf.cxx:125
 Roo2DKeysPdf.cxx:126
 Roo2DKeysPdf.cxx:127
 Roo2DKeysPdf.cxx:128
 Roo2DKeysPdf.cxx:129
 Roo2DKeysPdf.cxx:130
 Roo2DKeysPdf.cxx:131
 Roo2DKeysPdf.cxx:132
 Roo2DKeysPdf.cxx:133
 Roo2DKeysPdf.cxx:134
 Roo2DKeysPdf.cxx:135
 Roo2DKeysPdf.cxx:136
 Roo2DKeysPdf.cxx:137
 Roo2DKeysPdf.cxx:138
 Roo2DKeysPdf.cxx:139
 Roo2DKeysPdf.cxx:140
 Roo2DKeysPdf.cxx:141
 Roo2DKeysPdf.cxx:142
 Roo2DKeysPdf.cxx:143
 Roo2DKeysPdf.cxx:144
 Roo2DKeysPdf.cxx:145
 Roo2DKeysPdf.cxx:146
 Roo2DKeysPdf.cxx:147
 Roo2DKeysPdf.cxx:148
 Roo2DKeysPdf.cxx:149
 Roo2DKeysPdf.cxx:150
 Roo2DKeysPdf.cxx:151
 Roo2DKeysPdf.cxx:152
 Roo2DKeysPdf.cxx:153
 Roo2DKeysPdf.cxx:154
 Roo2DKeysPdf.cxx:155
 Roo2DKeysPdf.cxx:156
 Roo2DKeysPdf.cxx:157
 Roo2DKeysPdf.cxx:158
 Roo2DKeysPdf.cxx:159
 Roo2DKeysPdf.cxx:160
 Roo2DKeysPdf.cxx:161
 Roo2DKeysPdf.cxx:162
 Roo2DKeysPdf.cxx:163
 Roo2DKeysPdf.cxx:164
 Roo2DKeysPdf.cxx:165
 Roo2DKeysPdf.cxx:166
 Roo2DKeysPdf.cxx:167
 Roo2DKeysPdf.cxx:168
 Roo2DKeysPdf.cxx:169
 Roo2DKeysPdf.cxx:170
 Roo2DKeysPdf.cxx:171
 Roo2DKeysPdf.cxx:172
 Roo2DKeysPdf.cxx:173
 Roo2DKeysPdf.cxx:174
 Roo2DKeysPdf.cxx:175
 Roo2DKeysPdf.cxx:176
 Roo2DKeysPdf.cxx:177
 Roo2DKeysPdf.cxx:178
 Roo2DKeysPdf.cxx:179
 Roo2DKeysPdf.cxx:180
 Roo2DKeysPdf.cxx:181
 Roo2DKeysPdf.cxx:182
 Roo2DKeysPdf.cxx:183
 Roo2DKeysPdf.cxx:184
 Roo2DKeysPdf.cxx:185
 Roo2DKeysPdf.cxx:186
 Roo2DKeysPdf.cxx:187
 Roo2DKeysPdf.cxx:188
 Roo2DKeysPdf.cxx:189
 Roo2DKeysPdf.cxx:190
 Roo2DKeysPdf.cxx:191
 Roo2DKeysPdf.cxx:192
 Roo2DKeysPdf.cxx:193
 Roo2DKeysPdf.cxx:194
 Roo2DKeysPdf.cxx:195
 Roo2DKeysPdf.cxx:196
 Roo2DKeysPdf.cxx:197
 Roo2DKeysPdf.cxx:198
 Roo2DKeysPdf.cxx:199
 Roo2DKeysPdf.cxx:200
 Roo2DKeysPdf.cxx:201
 Roo2DKeysPdf.cxx:202
 Roo2DKeysPdf.cxx:203
 Roo2DKeysPdf.cxx:204
 Roo2DKeysPdf.cxx:205
 Roo2DKeysPdf.cxx:206
 Roo2DKeysPdf.cxx:207
 Roo2DKeysPdf.cxx:208
 Roo2DKeysPdf.cxx:209
 Roo2DKeysPdf.cxx:210
 Roo2DKeysPdf.cxx:211
 Roo2DKeysPdf.cxx:212
 Roo2DKeysPdf.cxx:213
 Roo2DKeysPdf.cxx:214
 Roo2DKeysPdf.cxx:215
 Roo2DKeysPdf.cxx:216
 Roo2DKeysPdf.cxx:217
 Roo2DKeysPdf.cxx:218
 Roo2DKeysPdf.cxx:219
 Roo2DKeysPdf.cxx:220
 Roo2DKeysPdf.cxx:221
 Roo2DKeysPdf.cxx:222
 Roo2DKeysPdf.cxx:223
 Roo2DKeysPdf.cxx:224
 Roo2DKeysPdf.cxx:225
 Roo2DKeysPdf.cxx:226
 Roo2DKeysPdf.cxx:227
 Roo2DKeysPdf.cxx:228
 Roo2DKeysPdf.cxx:229
 Roo2DKeysPdf.cxx:230
 Roo2DKeysPdf.cxx:231
 Roo2DKeysPdf.cxx:232
 Roo2DKeysPdf.cxx:233
 Roo2DKeysPdf.cxx:234
 Roo2DKeysPdf.cxx:235
 Roo2DKeysPdf.cxx:236
 Roo2DKeysPdf.cxx:237
 Roo2DKeysPdf.cxx:238
 Roo2DKeysPdf.cxx:239
 Roo2DKeysPdf.cxx:240
 Roo2DKeysPdf.cxx:241
 Roo2DKeysPdf.cxx:242
 Roo2DKeysPdf.cxx:243
 Roo2DKeysPdf.cxx:244
 Roo2DKeysPdf.cxx:245
 Roo2DKeysPdf.cxx:246
 Roo2DKeysPdf.cxx:247
 Roo2DKeysPdf.cxx:248
 Roo2DKeysPdf.cxx:249
 Roo2DKeysPdf.cxx:250
 Roo2DKeysPdf.cxx:251
 Roo2DKeysPdf.cxx:252
 Roo2DKeysPdf.cxx:253
 Roo2DKeysPdf.cxx:254
 Roo2DKeysPdf.cxx:255
 Roo2DKeysPdf.cxx:256
 Roo2DKeysPdf.cxx:257
 Roo2DKeysPdf.cxx:258
 Roo2DKeysPdf.cxx:259
 Roo2DKeysPdf.cxx:260
 Roo2DKeysPdf.cxx:261
 Roo2DKeysPdf.cxx:262
 Roo2DKeysPdf.cxx:263
 Roo2DKeysPdf.cxx:264
 Roo2DKeysPdf.cxx:265
 Roo2DKeysPdf.cxx:266
 Roo2DKeysPdf.cxx:267
 Roo2DKeysPdf.cxx:268
 Roo2DKeysPdf.cxx:269
 Roo2DKeysPdf.cxx:270
 Roo2DKeysPdf.cxx:271
 Roo2DKeysPdf.cxx:272
 Roo2DKeysPdf.cxx:273
 Roo2DKeysPdf.cxx:274
 Roo2DKeysPdf.cxx:275
 Roo2DKeysPdf.cxx:276
 Roo2DKeysPdf.cxx:277
 Roo2DKeysPdf.cxx:278
 Roo2DKeysPdf.cxx:279
 Roo2DKeysPdf.cxx:280
 Roo2DKeysPdf.cxx:281
 Roo2DKeysPdf.cxx:282
 Roo2DKeysPdf.cxx:283
 Roo2DKeysPdf.cxx:284
 Roo2DKeysPdf.cxx:285
 Roo2DKeysPdf.cxx:286
 Roo2DKeysPdf.cxx:287
 Roo2DKeysPdf.cxx:288
 Roo2DKeysPdf.cxx:289
 Roo2DKeysPdf.cxx:290
 Roo2DKeysPdf.cxx:291
 Roo2DKeysPdf.cxx:292
 Roo2DKeysPdf.cxx:293
 Roo2DKeysPdf.cxx:294
 Roo2DKeysPdf.cxx:295
 Roo2DKeysPdf.cxx:296
 Roo2DKeysPdf.cxx:297
 Roo2DKeysPdf.cxx:298
 Roo2DKeysPdf.cxx:299
 Roo2DKeysPdf.cxx:300
 Roo2DKeysPdf.cxx:301
 Roo2DKeysPdf.cxx:302
 Roo2DKeysPdf.cxx:303
 Roo2DKeysPdf.cxx:304
 Roo2DKeysPdf.cxx:305
 Roo2DKeysPdf.cxx:306
 Roo2DKeysPdf.cxx:307
 Roo2DKeysPdf.cxx:308
 Roo2DKeysPdf.cxx:309
 Roo2DKeysPdf.cxx:310
 Roo2DKeysPdf.cxx:311
 Roo2DKeysPdf.cxx:312
 Roo2DKeysPdf.cxx:313
 Roo2DKeysPdf.cxx:314
 Roo2DKeysPdf.cxx:315
 Roo2DKeysPdf.cxx:316
 Roo2DKeysPdf.cxx:317
 Roo2DKeysPdf.cxx:318
 Roo2DKeysPdf.cxx:319
 Roo2DKeysPdf.cxx:320
 Roo2DKeysPdf.cxx:321
 Roo2DKeysPdf.cxx:322
 Roo2DKeysPdf.cxx:323
 Roo2DKeysPdf.cxx:324
 Roo2DKeysPdf.cxx:325
 Roo2DKeysPdf.cxx:326
 Roo2DKeysPdf.cxx:327
 Roo2DKeysPdf.cxx:328
 Roo2DKeysPdf.cxx:329
 Roo2DKeysPdf.cxx:330
 Roo2DKeysPdf.cxx:331
 Roo2DKeysPdf.cxx:332
 Roo2DKeysPdf.cxx:333
 Roo2DKeysPdf.cxx:334
 Roo2DKeysPdf.cxx:335
 Roo2DKeysPdf.cxx:336
 Roo2DKeysPdf.cxx:337
 Roo2DKeysPdf.cxx:338
 Roo2DKeysPdf.cxx:339
 Roo2DKeysPdf.cxx:340
 Roo2DKeysPdf.cxx:341
 Roo2DKeysPdf.cxx:342
 Roo2DKeysPdf.cxx:343
 Roo2DKeysPdf.cxx:344
 Roo2DKeysPdf.cxx:345
 Roo2DKeysPdf.cxx:346
 Roo2DKeysPdf.cxx:347
 Roo2DKeysPdf.cxx:348
 Roo2DKeysPdf.cxx:349
 Roo2DKeysPdf.cxx:350
 Roo2DKeysPdf.cxx:351
 Roo2DKeysPdf.cxx:352
 Roo2DKeysPdf.cxx:353
 Roo2DKeysPdf.cxx:354
 Roo2DKeysPdf.cxx:355
 Roo2DKeysPdf.cxx:356
 Roo2DKeysPdf.cxx:357
 Roo2DKeysPdf.cxx:358
 Roo2DKeysPdf.cxx:359
 Roo2DKeysPdf.cxx:360
 Roo2DKeysPdf.cxx:361
 Roo2DKeysPdf.cxx:362
 Roo2DKeysPdf.cxx:363
 Roo2DKeysPdf.cxx:364
 Roo2DKeysPdf.cxx:365
 Roo2DKeysPdf.cxx:366
 Roo2DKeysPdf.cxx:367
 Roo2DKeysPdf.cxx:368
 Roo2DKeysPdf.cxx:369
 Roo2DKeysPdf.cxx:370
 Roo2DKeysPdf.cxx:371
 Roo2DKeysPdf.cxx:372
 Roo2DKeysPdf.cxx:373
 Roo2DKeysPdf.cxx:374
 Roo2DKeysPdf.cxx:375
 Roo2DKeysPdf.cxx:376
 Roo2DKeysPdf.cxx:377
 Roo2DKeysPdf.cxx:378
 Roo2DKeysPdf.cxx:379
 Roo2DKeysPdf.cxx:380
 Roo2DKeysPdf.cxx:381
 Roo2DKeysPdf.cxx:382
 Roo2DKeysPdf.cxx:383
 Roo2DKeysPdf.cxx:384
 Roo2DKeysPdf.cxx:385
 Roo2DKeysPdf.cxx:386
 Roo2DKeysPdf.cxx:387
 Roo2DKeysPdf.cxx:388
 Roo2DKeysPdf.cxx:389
 Roo2DKeysPdf.cxx:390
 Roo2DKeysPdf.cxx:391
 Roo2DKeysPdf.cxx:392
 Roo2DKeysPdf.cxx:393
 Roo2DKeysPdf.cxx:394
 Roo2DKeysPdf.cxx:395
 Roo2DKeysPdf.cxx:396
 Roo2DKeysPdf.cxx:397
 Roo2DKeysPdf.cxx:398
 Roo2DKeysPdf.cxx:399
 Roo2DKeysPdf.cxx:400
 Roo2DKeysPdf.cxx:401
 Roo2DKeysPdf.cxx:402
 Roo2DKeysPdf.cxx:403
 Roo2DKeysPdf.cxx:404
 Roo2DKeysPdf.cxx:405
 Roo2DKeysPdf.cxx:406
 Roo2DKeysPdf.cxx:407
 Roo2DKeysPdf.cxx:408
 Roo2DKeysPdf.cxx:409
 Roo2DKeysPdf.cxx:410
 Roo2DKeysPdf.cxx:411
 Roo2DKeysPdf.cxx:412
 Roo2DKeysPdf.cxx:413
 Roo2DKeysPdf.cxx:414
 Roo2DKeysPdf.cxx:415
 Roo2DKeysPdf.cxx:416
 Roo2DKeysPdf.cxx:417
 Roo2DKeysPdf.cxx:418
 Roo2DKeysPdf.cxx:419
 Roo2DKeysPdf.cxx:420
 Roo2DKeysPdf.cxx:421
 Roo2DKeysPdf.cxx:422
 Roo2DKeysPdf.cxx:423
 Roo2DKeysPdf.cxx:424
 Roo2DKeysPdf.cxx:425
 Roo2DKeysPdf.cxx:426
 Roo2DKeysPdf.cxx:427
 Roo2DKeysPdf.cxx:428
 Roo2DKeysPdf.cxx:429
 Roo2DKeysPdf.cxx:430
 Roo2DKeysPdf.cxx:431
 Roo2DKeysPdf.cxx:432
 Roo2DKeysPdf.cxx:433
 Roo2DKeysPdf.cxx:434
 Roo2DKeysPdf.cxx:435
 Roo2DKeysPdf.cxx:436
 Roo2DKeysPdf.cxx:437
 Roo2DKeysPdf.cxx:438
 Roo2DKeysPdf.cxx:439
 Roo2DKeysPdf.cxx:440
 Roo2DKeysPdf.cxx:441
 Roo2DKeysPdf.cxx:442
 Roo2DKeysPdf.cxx:443
 Roo2DKeysPdf.cxx:444
 Roo2DKeysPdf.cxx:445
 Roo2DKeysPdf.cxx:446
 Roo2DKeysPdf.cxx:447
 Roo2DKeysPdf.cxx:448
 Roo2DKeysPdf.cxx:449
 Roo2DKeysPdf.cxx:450
 Roo2DKeysPdf.cxx:451
 Roo2DKeysPdf.cxx:452
 Roo2DKeysPdf.cxx:453
 Roo2DKeysPdf.cxx:454
 Roo2DKeysPdf.cxx:455
 Roo2DKeysPdf.cxx:456
 Roo2DKeysPdf.cxx:457
 Roo2DKeysPdf.cxx:458
 Roo2DKeysPdf.cxx:459
 Roo2DKeysPdf.cxx:460
 Roo2DKeysPdf.cxx:461
 Roo2DKeysPdf.cxx:462
 Roo2DKeysPdf.cxx:463
 Roo2DKeysPdf.cxx:464
 Roo2DKeysPdf.cxx:465
 Roo2DKeysPdf.cxx:466
 Roo2DKeysPdf.cxx:467
 Roo2DKeysPdf.cxx:468
 Roo2DKeysPdf.cxx:469
 Roo2DKeysPdf.cxx:470
 Roo2DKeysPdf.cxx:471
 Roo2DKeysPdf.cxx:472
 Roo2DKeysPdf.cxx:473
 Roo2DKeysPdf.cxx:474
 Roo2DKeysPdf.cxx:475
 Roo2DKeysPdf.cxx:476
 Roo2DKeysPdf.cxx:477
 Roo2DKeysPdf.cxx:478
 Roo2DKeysPdf.cxx:479
 Roo2DKeysPdf.cxx:480
 Roo2DKeysPdf.cxx:481
 Roo2DKeysPdf.cxx:482
 Roo2DKeysPdf.cxx:483
 Roo2DKeysPdf.cxx:484
 Roo2DKeysPdf.cxx:485
 Roo2DKeysPdf.cxx:486
 Roo2DKeysPdf.cxx:487
 Roo2DKeysPdf.cxx:488
 Roo2DKeysPdf.cxx:489
 Roo2DKeysPdf.cxx:490
 Roo2DKeysPdf.cxx:491
 Roo2DKeysPdf.cxx:492
 Roo2DKeysPdf.cxx:493
 Roo2DKeysPdf.cxx:494
 Roo2DKeysPdf.cxx:495
 Roo2DKeysPdf.cxx:496
 Roo2DKeysPdf.cxx:497
 Roo2DKeysPdf.cxx:498
 Roo2DKeysPdf.cxx:499
 Roo2DKeysPdf.cxx:500
 Roo2DKeysPdf.cxx:501
 Roo2DKeysPdf.cxx:502
 Roo2DKeysPdf.cxx:503
 Roo2DKeysPdf.cxx:504
 Roo2DKeysPdf.cxx:505
 Roo2DKeysPdf.cxx:506
 Roo2DKeysPdf.cxx:507
 Roo2DKeysPdf.cxx:508
 Roo2DKeysPdf.cxx:509
 Roo2DKeysPdf.cxx:510
 Roo2DKeysPdf.cxx:511
 Roo2DKeysPdf.cxx:512
 Roo2DKeysPdf.cxx:513
 Roo2DKeysPdf.cxx:514
 Roo2DKeysPdf.cxx:515
 Roo2DKeysPdf.cxx:516
 Roo2DKeysPdf.cxx:517
 Roo2DKeysPdf.cxx:518
 Roo2DKeysPdf.cxx:519
 Roo2DKeysPdf.cxx:520
 Roo2DKeysPdf.cxx:521
 Roo2DKeysPdf.cxx:522
 Roo2DKeysPdf.cxx:523
 Roo2DKeysPdf.cxx:524
 Roo2DKeysPdf.cxx:525
 Roo2DKeysPdf.cxx:526
 Roo2DKeysPdf.cxx:527
 Roo2DKeysPdf.cxx:528
 Roo2DKeysPdf.cxx:529
 Roo2DKeysPdf.cxx:530
 Roo2DKeysPdf.cxx:531
 Roo2DKeysPdf.cxx:532
 Roo2DKeysPdf.cxx:533
 Roo2DKeysPdf.cxx:534
 Roo2DKeysPdf.cxx:535
 Roo2DKeysPdf.cxx:536
 Roo2DKeysPdf.cxx:537
 Roo2DKeysPdf.cxx:538
 Roo2DKeysPdf.cxx:539
 Roo2DKeysPdf.cxx:540
 Roo2DKeysPdf.cxx:541
 Roo2DKeysPdf.cxx:542
 Roo2DKeysPdf.cxx:543
 Roo2DKeysPdf.cxx:544
 Roo2DKeysPdf.cxx:545
 Roo2DKeysPdf.cxx:546
 Roo2DKeysPdf.cxx:547
 Roo2DKeysPdf.cxx:548
 Roo2DKeysPdf.cxx:549
 Roo2DKeysPdf.cxx:550
 Roo2DKeysPdf.cxx:551
 Roo2DKeysPdf.cxx:552
 Roo2DKeysPdf.cxx:553
 Roo2DKeysPdf.cxx:554
 Roo2DKeysPdf.cxx:555
 Roo2DKeysPdf.cxx:556
 Roo2DKeysPdf.cxx:557
 Roo2DKeysPdf.cxx:558
 Roo2DKeysPdf.cxx:559
 Roo2DKeysPdf.cxx:560
 Roo2DKeysPdf.cxx:561
 Roo2DKeysPdf.cxx:562
 Roo2DKeysPdf.cxx:563
 Roo2DKeysPdf.cxx:564
 Roo2DKeysPdf.cxx:565
 Roo2DKeysPdf.cxx:566
 Roo2DKeysPdf.cxx:567
 Roo2DKeysPdf.cxx:568
 Roo2DKeysPdf.cxx:569
 Roo2DKeysPdf.cxx:570
 Roo2DKeysPdf.cxx:571
 Roo2DKeysPdf.cxx:572
 Roo2DKeysPdf.cxx:573
 Roo2DKeysPdf.cxx:574
 Roo2DKeysPdf.cxx:575
 Roo2DKeysPdf.cxx:576
 Roo2DKeysPdf.cxx:577
 Roo2DKeysPdf.cxx:578
 Roo2DKeysPdf.cxx:579
 Roo2DKeysPdf.cxx:580
 Roo2DKeysPdf.cxx:581
 Roo2DKeysPdf.cxx:582
 Roo2DKeysPdf.cxx:583
 Roo2DKeysPdf.cxx:584
 Roo2DKeysPdf.cxx:585
 Roo2DKeysPdf.cxx:586
 Roo2DKeysPdf.cxx:587
 Roo2DKeysPdf.cxx:588
 Roo2DKeysPdf.cxx:589
 Roo2DKeysPdf.cxx:590
 Roo2DKeysPdf.cxx:591
 Roo2DKeysPdf.cxx:592
 Roo2DKeysPdf.cxx:593
 Roo2DKeysPdf.cxx:594