/*****************************************************************************
 * Project: RooFit                                                           *
 * Package: RooFitCore                                                       *
 * @(#)root/roofitcore:$Id$
 * Authors:                                                                  *
 *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
 *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
 *                                                                           *
 * Copyright (c) 2000-2005, Regents of the University of California          *
 *                          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
// RooDataSet is a container class to hold N-dimensional binned data. Each bins central 
// coordinates in N-dimensional space are represented by a RooArgSet of RooRealVar, RooCategory 
// or RooStringVar objects, thus data can be binned in real and/or discrete dimensions
// END_HTML
//
//

#include "RooFit.h"
#include "Riostream.h"

#include "TH1.h"
#include "TH1.h"
#include "TDirectory.h"
#include "TMath.h"
#include "RooMsgService.h"
#include "RooDataHist.h"
#include "RooDataHistSliceIter.h"
#include "RooAbsLValue.h"
#include "RooArgList.h"
#include "RooRealVar.h"
#include "RooMath.h"
#include "RooBinning.h"
#include "RooPlot.h"
#include "RooHistError.h"
#include "RooCategory.h"
#include "RooCmdConfig.h"
#include "RooLinkedListIter.h"
#include "RooTreeDataStore.h"
#include "RooVectorDataStore.h"
#include "TTree.h"
#include "RooTrace.h"
#include "RooTreeData.h"

using namespace std ;

ClassImp(RooDataHist) 
;



//_____________________________________________________________________________
RooDataHist::RooDataHist() : _pbinvCacheMgr(0,10)
{
  // Default constructor
  _arrSize = 0 ;
  _wgt = 0 ;
  _errLo = 0 ;
  _errHi = 0 ;
  _sumw2 = 0 ;
  _binv = 0 ;
  _pbinv = 0 ;
  _curWeight = 0 ;
  _curIndex = -1 ;
  _realIter = _realVars.createIterator() ;
  _binValid = 0 ;
  _curSumW2 = 0 ;
  _curVolume = 1 ;
  _curWgtErrHi = 0 ;
  _curWgtErrLo = 0 ;
  _cache_sum_valid = 0 ;
  TRACE_CREATE
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char *name, const char *title, const RooArgSet& vars, const char* binningName) : 
  RooAbsData(name,title,vars), _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of an empty data hist from a RooArgSet defining the dimensions
  // of the data space. The range and number of bins in each dimensions are taken
  // from getMin()getMax(),getBins() of each RooAbsArg representing that
  // dimension.
  //
  // For real dimensions, the fit range and number of bins can be set independently
  // of the plot range and number of bins, but it is advisable to keep the
  // ratio of the plot bin width and the fit bin width an integer value.
  // For category dimensions, the fit ranges always comprises all defined states
  // and each state is always has its individual bin
  //
  // To effective achive binning of real dimensions with variable bins sizes,
  // construct a RooThresholdCategory of the real dimension to be binned variably.
  // Set the thresholds at the desired bin boundaries, and construct the
  // data hist as function of the threshold category instead of the real variable.

  // Initialize datastore
  _dstore = (defaultStorageType==Tree) ? ((RooAbsDataStore*) new RooTreeDataStore(name,title,_vars)) : 
                                         ((RooAbsDataStore*) new RooVectorDataStore(name,title,_vars)) ;
  
  initialize(binningName) ;

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;

  appendToDir(this,kTRUE) ;
  TRACE_CREATE
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char *name, const char *title, const RooArgSet& vars, const RooAbsData& data, Double_t wgt) :
  RooAbsData(name,title,vars), _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of a data hist from an existing data collection (binned or unbinned)
  // The RooArgSet 'vars' defines the dimensions of the histogram. 
  // The range and number of bins in each dimensions are taken
  // from getMin()getMax(),getBins() of each RooAbsArg representing that
  // dimension.
  //
  // For real dimensions, the fit range and number of bins can be set independently
  // of the plot range and number of bins, but it is advisable to keep the
  // ratio of the plot bin width and the fit bin width an integer value.
  // For category dimensions, the fit ranges always comprises all defined states
  // and each state is always has its individual bin
  //
  // To effective achive binning of real dimensions with variable bins sizes,
  // construct a RooThresholdCategory of the real dimension to be binned variably.
  // Set the thresholds at the desired bin boundaries, and construct the
  // data hist as function of the threshold category instead of the real variable.
  //
  // If the constructed data hist has less dimensions that in source data collection,
  // all missing dimensions will be projected.

  // Initialize datastore
  _dstore = (defaultStorageType==Tree) ? ((RooAbsDataStore*) new RooTreeDataStore(name,title,_vars)) : 
                                         ((RooAbsDataStore*) new RooVectorDataStore(name,title,_vars)) ;

  initialize() ;
  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;

  add(data,(const RooFormulaVar*)0,wgt) ;
  appendToDir(this,kTRUE) ;
  TRACE_CREATE
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char *name, const char *title, const RooArgList& vars, RooCategory& indexCat, 
			 map<string,TH1*> histMap, Double_t wgt) :
  RooAbsData(name,title,RooArgSet(vars,&indexCat)), 
  _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of a data hist from a map of TH1,TH2 or TH3 that are collated into a x+1 dimensional
  // RooDataHist where the added dimension is a category that labels the input source as defined
  // in the histMap argument. The state names used in histMap must correspond to predefined states
  // 'indexCat'
  //
  // The RooArgList 'vars' defines the dimensions of the histogram. 
  // The ranges and number of bins are taken from the input histogram and must be the same in all histograms

  // Initialize datastore
  _dstore = (defaultStorageType==Tree) ? ((RooAbsDataStore*) new RooTreeDataStore(name,title,_vars)) : 
                                         ((RooAbsDataStore*) new RooVectorDataStore(name,title,_vars)) ;
  
  importTH1Set(vars, indexCat, histMap, wgt, kFALSE) ;

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;
  TRACE_CREATE
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char *name, const char *title, const RooArgList& vars, RooCategory& indexCat, 
			 map<string,RooDataHist*> dhistMap, Double_t wgt) :
  RooAbsData(name,title,RooArgSet(vars,&indexCat)), 
  _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of a data hist from a map of RooDataHists that are collated into a x+1 dimensional
  // RooDataHist where the added dimension is a category that labels the input source as defined
  // in the histMap argument. The state names used in histMap must correspond to predefined states
  // 'indexCat'
  //
  // The RooArgList 'vars' defines the dimensions of the histogram. 
  // The ranges and number of bins are taken from the input histogram and must be the same in all histograms

  // Initialize datastore
  _dstore = (defaultStorageType==Tree) ? ((RooAbsDataStore*) new RooTreeDataStore(name,title,_vars)) : 
                                         ((RooAbsDataStore*) new RooVectorDataStore(name,title,_vars)) ;
  
  importDHistSet(vars, indexCat, dhistMap, wgt) ;

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;
  TRACE_CREATE
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char *name, const char *title, const RooArgList& vars, const TH1* hist, Double_t wgt) :
  RooAbsData(name,title,vars), _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of a data hist from an TH1,TH2 or TH3
  // The RooArgSet 'vars' defines the dimensions of the histogram. The ranges
  // and number of bins are taken from the input histogram, and the corresponding
  // values are set accordingly on the arguments in 'vars'

  // Initialize datastore
  _dstore = (defaultStorageType==Tree) ? ((RooAbsDataStore*) new RooTreeDataStore(name,title,_vars)) : 
                                         ((RooAbsDataStore*) new RooVectorDataStore(name,title,_vars)) ;

  // Check consistency in number of dimensions
  if (vars.getSize() != hist->GetDimension()) {
    coutE(InputArguments) << "RooDataHist::ctor(" << GetName() << ") ERROR: dimension of input histogram must match "
			  << "number of dimension variables" << endl ;
    assert(0) ; 
  }

  importTH1(vars,*const_cast<TH1*>(hist),wgt, kFALSE) ;

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;
  TRACE_CREATE
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char *name, const char *title, const RooArgList& vars, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
			 const RooCmdArg& arg4,const RooCmdArg& arg5,const RooCmdArg& arg6,const RooCmdArg& arg7,const RooCmdArg& arg8) :
  RooAbsData(name,title,RooArgSet(vars,(RooAbsArg*)RooCmdConfig::decodeObjOnTheFly("RooDataHist::RooDataHist", "IndexCat",0,0,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8))), 
  _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of a binned dataset from a RooArgSet defining the dimensions
  // of the data space. The range and number of bins in each dimensions are taken
  // from getMin() getMax(),getBins() of each RooAbsArg representing that
  // dimension.
  //
  // This constructor takes the following optional arguments
  //
  // Import(TH1&, Bool_t impDens) -- Import contents of the given TH1/2/3 into this binned dataset. The 
  //                                 ranges and binning of the binned dataset are automatically adjusted to
  //                                 match those of the imported histogram. 
  //
  //                                 Please note: for TH1& with unequal binning _only_,
  //                                 you should decide if you want to import the absolute bin content,
  //                                 or the bin content expressed as density. The latter is default and will
  //                                 result in the same histogram as the original TH1. For certain type of
  //                                 bin contents (containing efficiencies, asymmetries, or ratio is general)
  //                                 you should import the absolute value and set impDens to kFALSE
  //                                 
  //
  // Weight(Double_t)          -- Apply given weight factor when importing histograms
  //
  // Index(RooCategory&)       -- Prepare import of multiple TH1/1/2/3 into a N+1 dimensional RooDataHist
  //                              where the extra discrete dimension labels the source of the imported histogram
  //                              If the index category defines states for which no histogram is be imported
  //                              the corresponding bins will be left empty.
  //                              
  // Import(const char*, TH1&) -- Import a THx to be associated with the given state name of the index category
  //                              specified in Index(). If the given state name is not yet defined in the index
  //                              category it will be added on the fly. The import command can be specified
  //                              multiple times. 
  // Import(map<string,TH1*>&) -- As above, but allows specification of many imports in a single operation
  //                              

  // Initialize datastore
  _dstore = (defaultStorageType==Tree) ? ((RooAbsDataStore*) new RooTreeDataStore(name,title,_vars)) : 
                                         ((RooAbsDataStore*) new RooVectorDataStore(name,title,_vars)) ;

  // Define configuration for this method
  RooCmdConfig pc(Form("RooDataHist::ctor(%s)",GetName())) ;
  pc.defineObject("impHist","ImportHisto",0) ;
  pc.defineInt("impDens","ImportHisto",0) ;
  pc.defineObject("indexCat","IndexCat",0) ;
  pc.defineObject("impSliceHist","ImportHistoSlice",0,0,kTRUE) ; // array
  pc.defineString("impSliceState","ImportHistoSlice",0,"",kTRUE) ; // array
  pc.defineObject("impSliceDHist","ImportDataHistSlice",0,0,kTRUE) ; // array
  pc.defineString("impSliceDState","ImportDataHistSlice",0,"",kTRUE) ; // array
  pc.defineDouble("weight","Weight",0,1) ; 
  pc.defineObject("dummy1","ImportDataHistSliceMany",0) ;
  pc.defineObject("dummy2","ImportHistoSliceMany",0) ;
  pc.defineMutex("ImportHisto","ImportHistoSlice","ImportDataHistSlice") ;
  pc.defineDependency("ImportHistoSlice","IndexCat") ;
  pc.defineDependency("ImportDataHistSlice","IndexCat") ;

  RooLinkedList l ;
  l.Add((TObject*)&arg1) ;  l.Add((TObject*)&arg2) ;  
  l.Add((TObject*)&arg3) ;  l.Add((TObject*)&arg4) ;
  l.Add((TObject*)&arg5) ;  l.Add((TObject*)&arg6) ;  
  l.Add((TObject*)&arg7) ;  l.Add((TObject*)&arg8) ;

  // Process & check varargs 
  pc.process(l) ;
  if (!pc.ok(kTRUE)) {
    assert(0) ;
    return ;
  }

  TH1* impHist = static_cast<TH1*>(pc.getObject("impHist")) ;
  Bool_t impDens = pc.getInt("impDens") ;
  Double_t initWgt = pc.getDouble("weight") ;
  const char* impSliceNames = pc.getString("impSliceState","",kTRUE) ;
  const RooLinkedList& impSliceHistos = pc.getObjectList("impSliceHist") ;
  RooCategory* indexCat = static_cast<RooCategory*>(pc.getObject("indexCat")) ;
  const char* impSliceDNames = pc.getString("impSliceDState","",kTRUE) ;
  const RooLinkedList& impSliceDHistos = pc.getObjectList("impSliceDHist") ;


  if (impHist) {
    
    // Initialize importing contents from TH1
    importTH1(vars,*impHist,initWgt, impDens) ;

  } else if (indexCat) {


    if (impSliceHistos.GetSize()>0) {

      // Initialize importing mapped set of TH1s
      map<string,TH1*> hmap ;
      char tmp[1024] ;
      strlcpy(tmp,impSliceNames,1024) ;
      char* token = strtok(tmp,",") ;
      TIterator* hiter = impSliceHistos.MakeIterator() ;
      while(token) {
	hmap[token] = (TH1*) hiter->Next() ;
	token = strtok(0,",") ;
      }
      importTH1Set(vars,*indexCat,hmap,initWgt,kFALSE) ;
    } else {

      // Initialize importing mapped set of RooDataHists
      map<string,RooDataHist*> dmap ;
      char tmp[1024] ;
      strlcpy(tmp,impSliceDNames,1024) ;
      char* token = strtok(tmp,",") ;
      TIterator* hiter = impSliceDHistos.MakeIterator() ;
      while(token) {
	dmap[token] = (RooDataHist*) hiter->Next() ;
	token = strtok(0,",") ;
      }
      importDHistSet(vars,*indexCat,dmap,initWgt) ;

    }


  } else {

    // Initialize empty
    initialize() ;
    appendToDir(this,kTRUE) ;    

  }

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;
  TRACE_CREATE

}




//_____________________________________________________________________________
void RooDataHist::importTH1(const RooArgList& vars, TH1& histo, Double_t wgt, Bool_t doDensityCorrection) 
{
  // Import data from given TH1/2/3 into this RooDataHist

  // Adjust binning of internal observables to match that of input THx
  Int_t offset[3] ;
  adjustBinning(vars,histo,offset) ;
  
  // Initialize internal data structure
  initialize() ;
  appendToDir(this,kTRUE) ;

  // Define x,y,z as 1st, 2nd and 3rd observable
  RooRealVar* xvar = (RooRealVar*) _vars.find(vars.at(0)->GetName()) ;
  RooRealVar* yvar = (RooRealVar*) (vars.at(1) ? _vars.find(vars.at(1)->GetName()) : 0 ) ;
  RooRealVar* zvar = (RooRealVar*) (vars.at(2) ? _vars.find(vars.at(2)->GetName()) : 0 ) ;

  // Transfer contents
  Int_t xmin(0),ymin(0),zmin(0) ;
  RooArgSet vset(*xvar) ;
  Double_t volume = xvar->getMax()-xvar->getMin() ;
  xmin = offset[0] ;
  if (yvar) {
    vset.add(*yvar) ;
    ymin = offset[1] ;
    volume *= (yvar->getMax()-yvar->getMin()) ;
  }
  if (zvar) {
    vset.add(*zvar) ;
    zmin = offset[2] ;
    volume *= (zvar->getMax()-zvar->getMin()) ;
  }
  //Double_t avgBV = volume / numEntries() ;
//   cout << "average bin volume = " << avgBV << endl ;

  Int_t ix(0),iy(0),iz(0) ;
  for (ix=0 ; ix < xvar->getBins() ; ix++) {
    xvar->setBin(ix) ;
    if (yvar) {
      for (iy=0 ; iy < yvar->getBins() ; iy++) {
	yvar->setBin(iy) ;
	if (zvar) {
	  for (iz=0 ; iz < zvar->getBins() ; iz++) {
	    zvar->setBin(iz) ;
	    Double_t bv = doDensityCorrection ? binVolume(vset) : 1;
	    add(vset,bv*histo.GetBinContent(ix+1+xmin,iy+1+ymin,iz+1+zmin)*wgt,bv*TMath::Power(histo.GetBinError(ix+1+xmin,iy+1+ymin,iz+1+zmin)*wgt,2)) ;
	  }
	} else {
	  Double_t bv = doDensityCorrection ? binVolume(vset) : 1;
	  add(vset,bv*histo.GetBinContent(ix+1+xmin,iy+1+ymin)*wgt,bv*TMath::Power(histo.GetBinError(ix+1+xmin,iy+1+ymin)*wgt,2)) ;
	}
      }
    } else {
      Double_t bv = doDensityCorrection ? binVolume(vset) : 1 ;
      add(vset,bv*histo.GetBinContent(ix+1+xmin)*wgt,bv*TMath::Power(histo.GetBinError(ix+1+xmin)*wgt,2)) ;	    
    }
  }  
  
}





//_____________________________________________________________________________
void RooDataHist::importTH1Set(const RooArgList& vars, RooCategory& indexCat, map<string,TH1*> hmap, Double_t wgt, Bool_t doDensityCorrection) 
{
  // Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
  // in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
  // and the import source

  RooCategory* icat = (RooCategory*) _vars.find(indexCat.GetName()) ;

  TH1* histo(0) ;  
  Bool_t init(kFALSE) ;
  for (map<string,TH1*>::iterator hiter = hmap.begin() ; hiter!=hmap.end() ; ++hiter) {
    // Store pointer to first histogram from which binning specification will be taken
    if (!histo) {
      histo = hiter->second ;
    }
    // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
    if (!indexCat.lookupType(hiter->first.c_str())) {
      indexCat.defineType(hiter->first.c_str()) ;
      coutI(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << ") defining state \"" << hiter->first << "\" in index category " << indexCat.GetName() << endl ;
    }
    if (!icat->lookupType(hiter->first.c_str())) {	
      icat->defineType(hiter->first.c_str()) ;
    }
  }

  // Check consistency in number of dimensions
  if (histo && (vars.getSize() != histo->GetDimension())) {
    coutE(InputArguments) << "RooDataHist::ctor(" << GetName() << ") ERROR: dimension of input histogram must match "
			  << "number of continuous variables" << endl ;
    assert(0) ; 
  }
  
  // Copy bins and ranges from THx to dimension observables
  Int_t offset[3] ;
  adjustBinning(vars,*histo,offset) ;
  
  // Initialize internal data structure
  if (!init) {
    initialize() ;
    appendToDir(this,kTRUE) ;
    init = kTRUE ;
  }

  // Define x,y,z as 1st, 2nd and 3rd observable
  RooRealVar* xvar = (RooRealVar*) _vars.find(vars.at(0)->GetName()) ;
  RooRealVar* yvar = (RooRealVar*) (vars.at(1) ? _vars.find(vars.at(1)->GetName()) : 0 ) ;
  RooRealVar* zvar = (RooRealVar*) (vars.at(2) ? _vars.find(vars.at(2)->GetName()) : 0 ) ;

  // Transfer contents
  Int_t xmin(0),ymin(0),zmin(0) ;
  RooArgSet vset(*xvar) ;
  Double_t volume = xvar->getMax()-xvar->getMin() ;
  xmin = offset[0] ;
  if (yvar) {
    vset.add(*yvar) ;
    ymin = offset[1] ;
    volume *= (yvar->getMax()-yvar->getMin()) ;
  }
  if (zvar) {
    vset.add(*zvar) ;
    zmin = offset[2] ;
    volume *= (zvar->getMax()-zvar->getMin()) ;
  }
  Double_t avgBV = volume / numEntries() ;
  
  Int_t ic(0),ix(0),iy(0),iz(0) ;
  for (ic=0 ; ic < icat->numBins(0) ; ic++) {
    icat->setBin(ic) ;
    histo = hmap[icat->getLabel()] ;
    for (ix=0 ; ix < xvar->getBins() ; ix++) {
      xvar->setBin(ix) ;
      if (yvar) {
	for (iy=0 ; iy < yvar->getBins() ; iy++) {
	  yvar->setBin(iy) ;
	  if (zvar) {
	    for (iz=0 ; iz < zvar->getBins() ; iz++) {
	      zvar->setBin(iz) ;
	      Double_t bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
	      add(vset,bv*histo->GetBinContent(ix+1+xmin,iy+1+ymin,iz+1+zmin)*wgt,bv*TMath::Power(histo->GetBinError(ix+1+xmin,iy+1+ymin,iz+1+zmin)*wgt,2)) ;
	    }
	  } else {
	    Double_t bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
	    add(vset,bv*histo->GetBinContent(ix+1+xmin,iy+1+ymin)*wgt,bv*TMath::Power(histo->GetBinError(ix+1+xmin,iy+1+ymin)*wgt,2)) ;
	  }
	}
      } else {
	Double_t bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
	add(vset,bv*histo->GetBinContent(ix+1+xmin)*wgt,bv*TMath::Power(histo->GetBinError(ix+1+xmin)*wgt,2)) ;	    
      }
    }  
  }

}



//_____________________________________________________________________________
void RooDataHist::importDHistSet(const RooArgList& /*vars*/, RooCategory& indexCat, std::map<std::string,RooDataHist*> dmap, Double_t initWgt) 
{
  // Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
  // in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
  // and the import source

  RooCategory* icat = (RooCategory*) _vars.find(indexCat.GetName()) ;

  for (map<string,RooDataHist*>::iterator diter = dmap.begin() ; diter!=dmap.end() ; ++diter) {

    // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
    if (!indexCat.lookupType(diter->first.c_str())) {
      indexCat.defineType(diter->first.c_str()) ;
      coutI(InputArguments) << "RooDataHist::importDHistSet(" << GetName() << ") defining state \"" << diter->first << "\" in index category " << indexCat.GetName() << endl ;
    }
    if (!icat->lookupType(diter->first.c_str())) {	
      icat->defineType(diter->first.c_str()) ;
    }
  }

  initialize() ;
  appendToDir(this,kTRUE) ;  


  for (map<string,RooDataHist*>::iterator diter = dmap.begin() ; diter!=dmap.end() ; ++diter) {

    RooDataHist* dhist = diter->second ;

    icat->setLabel(diter->first.c_str()) ;

    // Transfer contents
    for (Int_t i=0 ; i<dhist->numEntries() ; i++) {
      _vars = *dhist->get(i) ;
      add(_vars,dhist->weight()*initWgt, pow(dhist->weightError(SumW2),2) ) ;
    }

  }
}

//_____________________________________________________________________________
void RooDataHist::adjustBinning(const RooArgList& vars, TH1& href, Int_t* offset) 
{
  // Adjust binning specification on first and optionally second and third
  // observable to binning in given reference TH1. Used by constructors
  // that import data from an external TH1

  // X
  RooRealVar* xvar = (RooRealVar*) _vars.find(*vars.at(0)) ;
  if (!dynamic_cast<RooRealVar*>(xvar)) {
    coutE(InputArguments) << "RooDataHist::adjustBinning(" << GetName() << ") ERROR: dimension " << xvar->GetName() << " must be real" << endl ;
    assert(0) ;
  }

  Double_t xlo = ((RooRealVar*)vars.at(0))->getMin() ;
  Double_t xhi = ((RooRealVar*)vars.at(0))->getMax() ;
  Int_t xmin(0) ;
  if (href.GetXaxis()->GetXbins()->GetArray()) {

    RooBinning xbins(href.GetNbinsX(),href.GetXaxis()->GetXbins()->GetArray()) ;

    Double_t tolerance = 1e-6*xbins.averageBinWidth() ;
    
    // Adjust xlo/xhi to nearest boundary
    Double_t xloAdj = xbins.binLow(xbins.binNumber(xlo+tolerance)) ;
    Double_t xhiAdj = xbins.binHigh(xbins.binNumber(xhi-tolerance)) ;
    xbins.setRange(xloAdj,xhiAdj) ;

    ((RooRealVar*)vars.at(0))->setBinning(xbins) ;
    if (fabs(xloAdj-xlo)>tolerance||fabs(xhiAdj-xhi)<tolerance) {
      coutI(DataHandling) << "RooDataHist::adjustBinning(" << GetName() << "): fit range of variable " << xvar->GetName() << " expanded to nearest bin boundaries: [" 
			  << xlo << "," << xhi << "] --> [" << xloAdj << "," << xhiAdj << "]" << endl ;
    }

    xvar->setBinning(xbins) ;
    xmin = xbins.rawBinNumber(xloAdj+tolerance) ;
    if (offset) {
      offset[0] = xmin ;
    }

  } else {

    RooBinning xbins(href.GetXaxis()->GetXmin(),href.GetXaxis()->GetXmax()) ;
    xbins.addUniform(href.GetNbinsX(),href.GetXaxis()->GetXmin(),href.GetXaxis()->GetXmax()) ;

    Double_t tolerance = 1e-6*xbins.averageBinWidth() ;

    // Adjust xlo/xhi to nearest boundary
    Double_t xloAdj = xbins.binLow(xbins.binNumber(xlo+tolerance)) ;
    Double_t xhiAdj = xbins.binHigh(xbins.binNumber(xhi-tolerance)) ;
    xbins.setRange(xloAdj,xhiAdj) ;
    ((RooRealVar*)vars.at(0))->setRange(xloAdj,xhiAdj) ;
    if (fabs(xloAdj-xlo)>tolerance||fabs(xhiAdj-xhi)<tolerance) {
      coutI(DataHandling) << "RooDataHist::adjustBinning(" << GetName() << "): fit range of variable " << xvar->GetName() << " expanded to nearest bin boundaries: [" 
			  << xlo << "," << xhi << "] --> [" << xloAdj << "," << xhiAdj << "]" << endl ;
    }


    RooUniformBinning xbins2(xloAdj,xhiAdj,xbins.numBins()) ;
    xvar->setBinning(xbins2) ;
    xmin = xbins.rawBinNumber(xloAdj+tolerance) ;
    if (offset) {
      offset[0] = xmin ;
    }
  }



  // Y
  RooRealVar* yvar = (RooRealVar*) (vars.at(1) ? _vars.find(*vars.at(1)) : 0 ) ;
  Int_t ymin(0) ;
  if (yvar) {
    Double_t ylo = ((RooRealVar*)vars.at(1))->getMin() ;
    Double_t yhi = ((RooRealVar*)vars.at(1))->getMax() ;

    if (!dynamic_cast<RooRealVar*>(yvar)) {
      coutE(InputArguments) << "RooDataHist::adjustBinning(" << GetName() << ") ERROR: dimension " << yvar->GetName() << " must be real" << endl ;
      assert(0) ;
    }

    if (href.GetYaxis()->GetXbins()->GetArray()) {

      RooBinning ybins(href.GetNbinsY(),href.GetYaxis()->GetXbins()->GetArray()) ;

      Double_t tolerance = 1e-6*ybins.averageBinWidth() ;
      
      // Adjust ylo/yhi to nearest boundary
      Double_t yloAdj = ybins.binLow(ybins.binNumber(ylo+tolerance)) ;
      Double_t yhiAdj = ybins.binHigh(ybins.binNumber(yhi-tolerance)) ;
      ybins.setRange(yloAdj,yhiAdj) ;
      ((RooRealVar*)vars.at(1))->setBinning(ybins) ;
      if (fabs(yloAdj-ylo)>tolerance||fabs(yhiAdj-yhi)<tolerance) {
	coutI(DataHandling) << "RooDataHist::adjustBinning(" << GetName() << "): fit range of variable " << yvar->GetName() << " expanded to nearest bin boundaries: [" 
			    << ylo << "," << yhi << "] --> [" << yloAdj << "," << yhiAdj << "]" << endl ;
      }

      yvar->setBinning(ybins) ;
      ymin = ybins.rawBinNumber(yloAdj+tolerance) ;
      if (offset) {
	offset[1] = ymin ;
      }

    } else {

      RooBinning ybins(href.GetYaxis()->GetXmin(),href.GetYaxis()->GetXmax()) ;
      ybins.addUniform(href.GetNbinsY(),href.GetYaxis()->GetXmin(),href.GetYaxis()->GetXmax()) ;

      Double_t tolerance = 1e-6*ybins.averageBinWidth() ;
      
      // Adjust ylo/yhi to nearest boundary
      Double_t yloAdj = ybins.binLow(ybins.binNumber(ylo+tolerance)) ;
      Double_t yhiAdj = ybins.binHigh(ybins.binNumber(yhi-tolerance)) ;
      ybins.setRange(yloAdj,yhiAdj) ;
      ((RooRealVar*)vars.at(1))->setRange(yloAdj,yhiAdj) ;
      if (fabs(yloAdj-ylo)>tolerance||fabs(yhiAdj-yhi)<tolerance) {
	coutI(DataHandling) << "RooDataHist::adjustBinning(" << GetName() << "): fit range of variable " << yvar->GetName() << " expanded to nearest bin boundaries: [" 
			    << ylo << "," << yhi << "] --> [" << yloAdj << "," << yhiAdj << "]" << endl ;
      }
      
      RooUniformBinning ybins2(yloAdj,yhiAdj,ybins.numBins()) ;
      yvar->setBinning(ybins2) ;
      ymin = ybins.rawBinNumber(yloAdj+tolerance) ;
      if (offset) {
	offset[1] = ymin ;
      }

    }    
  }
  
  // Z
  RooRealVar* zvar = (RooRealVar*) (vars.at(2) ? _vars.find(*vars.at(2)) : 0 ) ;
  Int_t zmin(0) ;
  if (zvar) {
    Double_t zlo = ((RooRealVar*)vars.at(2))->getMin() ;
    Double_t zhi = ((RooRealVar*)vars.at(2))->getMax() ;

    if (!dynamic_cast<RooRealVar*>(zvar)) {
      coutE(InputArguments) << "RooDataHist::adjustBinning(" << GetName() << ") ERROR: dimension " << zvar->GetName() << " must be real" << endl ;
      assert(0) ;
    }

    if (href.GetZaxis()->GetXbins()->GetArray()) {

      RooBinning zbins(href.GetNbinsZ(),href.GetZaxis()->GetXbins()->GetArray()) ;

      Double_t tolerance = 1e-6*zbins.averageBinWidth() ;
      
      // Adjust zlo/zhi to nearest boundary
      Double_t zloAdj = zbins.binLow(zbins.binNumber(zlo+tolerance)) ;
      Double_t zhiAdj = zbins.binHigh(zbins.binNumber(zhi-tolerance)) ;
      zbins.setRange(zloAdj,zhiAdj) ;
      ((RooRealVar*)vars.at(2))->setBinning(zbins) ;
      if (fabs(zloAdj-zlo)>tolerance||fabs(zhiAdj-zhi)<tolerance) {
	coutI(DataHandling) << "RooDataHist::adjustBinning(" << GetName() << "): fit range of variable " << zvar->GetName() << " expanded to nearest bin boundaries: [" 
			    << zlo << "," << zhi << "] --> [" << zloAdj << "," << zhiAdj << "]" << endl ;
      }
      
      zvar->setBinning(zbins) ;
      zmin = zbins.rawBinNumber(zloAdj+tolerance) ;
      if (offset) {
	offset[2] = zmin ;
      }
      
    } else {

      RooBinning zbins(href.GetZaxis()->GetXmin(),href.GetZaxis()->GetXmax()) ;
      zbins.addUniform(href.GetNbinsZ(),href.GetZaxis()->GetXmin(),href.GetZaxis()->GetXmax()) ;

      Double_t tolerance = 1e-6*zbins.averageBinWidth() ;
      
      // Adjust zlo/zhi to nearest boundary
      Double_t zloAdj = zbins.binLow(zbins.binNumber(zlo+tolerance)) ;
      Double_t zhiAdj = zbins.binHigh(zbins.binNumber(zhi-tolerance)) ;
      zbins.setRange(zloAdj,zhiAdj) ;
      ((RooRealVar*)vars.at(2))->setRange(zloAdj,zhiAdj) ;
      if (fabs(zloAdj-zlo)>tolerance||fabs(zhiAdj-zhi)<tolerance) {
	coutI(DataHandling) << "RooDataHist::adjustBinning(" << GetName() << "): fit range of variable " << zvar->GetName() << " expanded to nearest bin boundaries: [" 
			    << zlo << "," << zhi << "] --> [" << zloAdj << "," << zhiAdj << "]" << endl ;
      }
      
      RooUniformBinning zbins2(zloAdj,zhiAdj,zbins.numBins()) ;
      zvar->setBinning(zbins2) ;
      zmin = zbins.rawBinNumber(zloAdj+tolerance) ;
      if (offset) {
	offset[2] = zmin ;
      }
    }
  }

}





//_____________________________________________________________________________
void RooDataHist::initialize(const char* binningName, Bool_t fillTree)
{
  // Initialization procedure: allocate weights array, calculate
  // multipliers needed for N-space to 1-dim array jump table,
  // and fill the internal tree with all bin center coordinates


  // Save real dimensions of dataset separately
  RooAbsArg* real ;
  _iterator->Reset() ;
  while((real=(RooAbsArg*)_iterator->Next())) {
    if (dynamic_cast<RooAbsReal*>(real)) _realVars.add(*real);
  }
  _realIter = _realVars.createIterator() ;

  // Fill array of LValue pointers to variables
  _iterator->Reset();
  RooAbsArg* rvarg;
  while((rvarg=(RooAbsArg*)_iterator->Next())) {
    if (binningName) {
      RooRealVar* rrv = dynamic_cast<RooRealVar*>(rvarg); 
      if (rrv) {
	rrv->setBinning(rrv->getBinning(binningName));
      }
    }
    // coverity[FORWARD_NULL]
    _lvvars.push_back(dynamic_cast<RooAbsLValue*>(rvarg));    
    // coverity[FORWARD_NULL]
    const RooAbsBinning* binning = dynamic_cast<RooAbsLValue*>(rvarg)->getBinningPtr(0);
    _lvbins.push_back(binning ? binning->clone() : 0);
  }

  
  // Allocate coefficients array
  _idxMult.resize(_vars.getSize()) ;

  _arrSize = 1 ;
  _iterator->Reset() ;
  RooAbsLValue* arg ;
  Int_t n(0), i ;
  while((arg=dynamic_cast<RooAbsLValue*>(_iterator->Next()))) {
    
    // Calculate sub-index multipliers for master index
    for (i=0 ; i<n ; i++) {
      _idxMult[i] *= arg->numBins() ;
    }
    _idxMult[n++] = 1 ;

    // Calculate dimension of weight array
    _arrSize *= arg->numBins() ;
  }  

  // Allocate and initialize weight array if necessary
  if (!_wgt) {
    _wgt = new Double_t[_arrSize] ;
    _errLo = new Double_t[_arrSize] ;
    _errHi = new Double_t[_arrSize] ;
    _sumw2 = new Double_t[_arrSize] ;
    _binv = new Double_t[_arrSize] ;

    // Refill array pointers in data store when reading
    // from Streamer
    if (fillTree==kFALSE) {
      _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;      
    }
    
    for (i=0 ; i<_arrSize ; i++) {
      _wgt[i] = 0 ;
      _errLo[i] = -1 ;
      _errHi[i] = -1 ;
      _sumw2[i] = 0 ;
    }
  }

  if (!fillTree) return ;

  // Fill TTree with bin center coordinates
  // Calculate plot bins of components from master index

  Int_t ibin ;
  for (ibin=0 ; ibin<_arrSize ; ibin++) {
    _iterator->Reset() ;
    RooAbsLValue* arg2 ;
    Int_t j(0), idx(0), tmp(ibin) ;
    Double_t theBinVolume(1) ;
    while((arg2=dynamic_cast<RooAbsLValue*>(_iterator->Next()))) {
      idx  = tmp / _idxMult[j] ;
      tmp -= idx*_idxMult[j++] ;
      RooAbsLValue* arglv = dynamic_cast<RooAbsLValue*>(arg2) ;
      arglv->setBin(idx) ;
      theBinVolume *= arglv->getBinWidth(idx) ;
//       cout << "init: bin width at idx=" << idx << " = " << arglv->getBinWidth(idx) << " binv[" << idx << "] = " << theBinVolume << endl ;
    }
    _binv[ibin] = theBinVolume ;
//     cout << "binv[" << ibin << "] = " << theBinVolume << endl ;
    fill() ;
  }


}


//_____________________________________________________________________________
void RooDataHist::checkBinBounds() const
{
  if (!_binbounds.empty()) return;
  for (std::vector<const RooAbsBinning*>::const_iterator it = _lvbins.begin();
      _lvbins.end() != it; ++it) {
    _binbounds.push_back(std::vector<Double_t>());
    if (*it) {
      std::vector<Double_t>& bounds = _binbounds.back();
      bounds.reserve(2 * (*it)->numBins());
      for (Int_t i = 0; i < (*it)->numBins(); ++i) {
	bounds.push_back((*it)->binLow(i));
	bounds.push_back((*it)->binHigh(i));
      }
    }
  }
}

//_____________________________________________________________________________
RooDataHist::RooDataHist(const RooDataHist& other, const char* newname) :
  RooAbsData(other,newname), RooDirItem(), _idxMult(other._idxMult), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(other._pbinvCacheMgr,0), _cache_sum_valid(0)
{
  // Copy constructor

  Int_t i ;

  // Allocate and initialize weight array 
  _arrSize = other._arrSize ;
  _wgt = new Double_t[_arrSize] ;
  _errLo = new Double_t[_arrSize] ;
  _errHi = new Double_t[_arrSize] ;
  _binv = new Double_t[_arrSize] ;
  _sumw2 = new Double_t[_arrSize] ;
  for (i=0 ; i<_arrSize ; i++) {
    _wgt[i] = other._wgt[i] ;
    _errLo[i] = other._errLo[i] ;
    _errHi[i] = other._errHi[i] ;
    _sumw2[i] = other._sumw2[i] ;
    _binv[i] = other._binv[i] ;
  }  

  // Save real dimensions of dataset separately
  RooAbsArg* arg ;
  _iterator->Reset() ;
  while((arg=(RooAbsArg*)_iterator->Next())) {
    if (dynamic_cast<RooAbsReal*>(arg)) _realVars.add(*arg) ;
  }
  _realIter = _realVars.createIterator() ;

  // Fill array of LValue pointers to variables
  _iterator->Reset() ;
  RooAbsArg* rvarg ;
  while((rvarg=(RooAbsArg*)_iterator->Next())) {
    // coverity[FORWARD_NULL]
    _lvvars.push_back(dynamic_cast<RooAbsLValue*>(rvarg)) ;
    // coverity[FORWARD_NULL]
    const RooAbsBinning* binning = dynamic_cast<RooAbsLValue*>(rvarg)->getBinningPtr(0) ;
    _lvbins.push_back(binning ? binning->clone() : 0) ;    
  }

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;

 appendToDir(this,kTRUE) ;
}



//_____________________________________________________________________________
RooDataHist::RooDataHist(const char* name, const char* title, RooDataHist* h, const RooArgSet& varSubset, 
			 const RooFormulaVar* cutVar, const char* cutRange, Int_t nStart, Int_t nStop, Bool_t copyCache) :
  RooAbsData(name,title,varSubset),
  _wgt(0), _binValid(0), _curWeight(0), _curVolume(1), _pbinv(0), _pbinvCacheMgr(0,10), _cache_sum_valid(0)
{
  // Constructor of a data hist from (part of) an existing data hist. The dimensions
  // of the data set are defined by the 'vars' RooArgSet, which can be identical
  // to 'dset' dimensions, or a subset thereof. Reduced dimensions will be projected
  // in the output data hist. The optional 'cutVar' formula variable can used to 
  // select the subset of bins to be copied.
  //
  // For most uses the RooAbsData::reduce() wrapper function, which uses this constructor, 
  // is the most convenient way to create a subset of an existing data  

  // Initialize datastore
  _dstore = new RooTreeDataStore(name,title,*h->_dstore,_vars,cutVar,cutRange,nStart,nStop,copyCache) ;
  
  initialize(0,kFALSE) ;

  _dstore->setExternalWeightArray(_wgt,_errLo,_errHi,_sumw2) ;

  // Copy weight array etc
  Int_t i ;
  for (i=0 ; i<_arrSize ; i++) {
    _wgt[i] = h->_wgt[i] ;
    _errLo[i] = h->_errLo[i] ;
    _errHi[i] = h->_errHi[i] ;
    _sumw2[i] = h->_sumw2[i] ;
    _binv[i] = h->_binv[i] ;
  }  


  appendToDir(this,kTRUE) ;
  TRACE_CREATE
}


//_____________________________________________________________________________
RooAbsData* RooDataHist::cacheClone(const RooAbsArg* newCacheOwner, const RooArgSet* newCacheVars, const char* newName) 
{
  // Construct a clone of this dataset that contains only the cached variables
  checkInit() ;

  RooDataHist* dhist = new RooDataHist(newName?newName:GetName(),GetTitle(),this,*get(),0,0,0,2000000000,kTRUE) ; 

  RooArgSet* selCacheVars = (RooArgSet*) newCacheVars->selectCommon(dhist->_cachedVars) ;
  dhist->attachCache(newCacheOwner, *selCacheVars) ;
  delete selCacheVars ;

  return dhist ;
}



//_____________________________________________________________________________
RooAbsData* RooDataHist::reduceEng(const RooArgSet& varSubset, const RooFormulaVar* cutVar, const char* cutRange, 
				   Int_t nStart, Int_t nStop, Bool_t /*copyCache*/)
{
  // Implementation of RooAbsData virtual method that drives the RooAbsData::reduce() methods

  checkInit() ;
  RooArgSet* myVarSubset = (RooArgSet*) _vars.selectCommon(varSubset) ;
  RooDataHist *rdh = new RooDataHist(GetName(), GetTitle(), *myVarSubset) ;
  delete myVarSubset ;

  RooFormulaVar* cloneVar = 0;
  RooArgSet* tmp(0) ;
  if (cutVar) {
    // Deep clone cutVar and attach clone to this dataset
    tmp = (RooArgSet*) RooArgSet(*cutVar).snapshot() ;
    if (!tmp) {
      coutE(DataHandling) << "RooDataHist::reduceEng(" << GetName() << ") Couldn't deep-clone cut variable, abort," << endl ;
      return 0 ;
    }
    cloneVar = (RooFormulaVar*) tmp->find(*cutVar) ;
    cloneVar->attachDataSet(*this) ;
  }

  Int_t i ;
  Double_t lo,hi ;
  Int_t nevt = nStop < numEntries() ? nStop : numEntries() ;
  TIterator* vIter = get()->createIterator() ;
  for (i=nStart ; i<nevt ; i++) {
    const RooArgSet* row = get(i) ;

    Bool_t doSelect(kTRUE) ;
    if (cutRange) {
      RooAbsArg* arg ;
      vIter->Reset() ;
      while((arg=(RooAbsArg*)vIter->Next())) {	
	if (!arg->inRange(cutRange)) {
	  doSelect = kFALSE ;
	  break ;
	}
      }
    }
    if (!doSelect) continue ;

    if (!cloneVar || cloneVar->getVal()) {
      weightError(lo,hi,SumW2) ;
      rdh->add(*row,weight(),lo*lo) ;
    }
  }
  delete vIter ;

  if (cloneVar) {
    delete tmp ;
  } 
  
    return rdh ;
  }



//_____________________________________________________________________________
RooDataHist::~RooDataHist() 
{
  // Destructor

  if (_wgt) delete[] _wgt ;
  if (_errLo) delete[] _errLo ;
  if (_errHi) delete[] _errHi ;
  if (_sumw2) delete[] _sumw2 ;
  if (_binv) delete[] _binv ;
  if (_realIter) delete _realIter ;
  if (_binValid) delete[] _binValid ;
  vector<const RooAbsBinning*>::iterator iter = _lvbins.begin() ;
  while(iter!=_lvbins.end()) {
    delete *iter ;
    iter++ ;
  }

   removeFromDir(this) ;
  TRACE_DESTROY
}




//_____________________________________________________________________________
Int_t RooDataHist::getIndex(const RooArgSet& coord, Bool_t fast)
{
  checkInit() ;
  if (fast) {
    _vars.assignFast(coord,kFALSE) ;
  } else {
    _vars.assignValueOnly(coord) ;
  }
  return calcTreeIndex() ;
}




//_____________________________________________________________________________
Int_t RooDataHist::calcTreeIndex() const 
{
  // Calculate the index for the weights array corresponding to 
  // to the bin enclosing the current coordinates of the internal argset

  Int_t masterIdx(0), i(0) ;
  vector<RooAbsLValue*>::const_iterator iter = _lvvars.begin() ;
  vector<const RooAbsBinning*>::const_iterator biter = _lvbins.begin() ;
  for (;iter!=_lvvars.end() ; ++iter) {
    const RooAbsBinning* binning = (*biter) ;
    masterIdx += _idxMult[i++]*(*iter)->getBin(binning) ;
    biter++ ;
  }
  return masterIdx ;
}



//_____________________________________________________________________________
void RooDataHist::dump2() 
{  
  // Debug stuff, should go...
  cout << "_arrSize = " << _arrSize << endl ;
  for (Int_t i=0 ; i<_arrSize ; i++) {
    cout << "wgt[" << i << "] = " << _wgt[i] << "sumw2[" << i << "] = " << _sumw2[i] << " vol[" << i << "] = " << _binv[i] << endl ;
  }
}



//_____________________________________________________________________________
RooPlot *RooDataHist::plotOn(RooPlot *frame, PlotOpt o) const 
{
  // Back end function to plotting functionality. Plot RooDataHist on given
  // frame in mode specified by plot options 'o'. The main purpose of
  // this function is to match the specified binning on 'o' to the
  // internal binning of the plot observable in this RooDataHist.
  checkInit() ;
  if (o.bins) return RooAbsData::plotOn(frame,o) ;

  if(0 == frame) {
    coutE(InputArguments) << ClassName() << "::" << GetName() << ":plotOn: frame is null" << endl;
    return 0;
  }
  RooAbsRealLValue *var= (RooAbsRealLValue*) frame->getPlotVar();
  if(0 == var) {
    coutE(InputArguments) << ClassName() << "::" << GetName()
	 << ":plotOn: frame does not specify a plot variable" << endl;
    return 0;
  }

  RooRealVar* dataVar = (RooRealVar*) _vars.find(*var) ;
  if (!dataVar) {
    coutE(InputArguments) << ClassName() << "::" << GetName()
	 << ":plotOn: dataset doesn't contain plot frame variable" << endl;
    return 0;
  }

  o.bins = &dataVar->getBinning() ;
  o.correctForBinWidth = kFALSE ;
  return RooAbsData::plotOn(frame,o) ;
}




//_____________________________________________________________________________
Double_t RooDataHist::weightSquared() const {
  return _curSumW2 ;
}



//_____________________________________________________________________________
Double_t RooDataHist::weight(const RooArgSet& bin, Int_t intOrder, Bool_t correctForBinSize, Bool_t cdfBoundaries, Bool_t oneSafe) 
{
  // Return the weight at given coordinates with optional
  // interpolation. If intOrder is zero, the weight
  // for the bin enclosing the coordinates
  // contained in 'bin' is returned. For higher values,
  // the result is interpolated in the real dimensions 
  // of the dataset
  
  //cout << "RooDataHist::weight(" << bin << "," << intOrder << "," << correctForBinSize << "," << cdfBoundaries << "," << oneSafe << ")" << endl ;

  checkInit() ;

  // Handle illegal intOrder values
  if (intOrder<0) {
    coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << endl ;
    return 0 ;
  }

  // Handle no-interpolation case
  if (intOrder==0) {    
    _vars.assignValueOnly(bin,oneSafe) ;
    Int_t idx = calcTreeIndex() ;
    //cout << "intOrder 0, idx = " << idx << endl ;
    if (correctForBinSize) {
      //calculatePartialBinVolume(*get()) ;
      //cout << "binw[" << idx << "] = " << _wgt[idx] <<  " / " << _binv[idx] << endl ;
      return _wgt[idx] / _binv[idx] ;
    } else {
      //cout << "binw[" << idx << "] = " << _wgt[idx] << endl ;
      return _wgt[idx] ;
    }
  }

  // Handle all interpolation cases
  _vars.assignValueOnly(bin) ;

  Double_t wInt(0) ;
  if (_realVars.getSize()==1) {

    // 1-dimensional interpolation
    RooFIter realIter = _realVars.fwdIterator() ;
    RooRealVar* real=(RooRealVar*)realIter.next() ;
    const RooAbsBinning* binning = real->getBinningPtr(0) ;
    wInt = interpolateDim(*real,binning,((RooAbsReal*)bin.find(*real))->getVal(), intOrder, correctForBinSize, cdfBoundaries) ;
    
  } else if (_realVars.getSize()==2) {

    // 2-dimensional interpolation
    RooFIter realIter = _realVars.fwdIterator() ;
    RooRealVar* realX=(RooRealVar*)realIter.next() ;
    RooRealVar* realY=(RooRealVar*)realIter.next() ;
    Double_t xval = ((RooAbsReal*)bin.find(*realX))->getVal() ;
    Double_t yval = ((RooAbsReal*)bin.find(*realY))->getVal() ;
    
    Int_t ybinC = realY->getBin() ;
    Int_t ybinLo = ybinC-intOrder/2 - ((yval<realY->getBinning().binCenter(ybinC))?1:0) ;
    Int_t ybinM = realY->numBins() ;
    
    Int_t i ;
    Double_t yarr[10] ;
    Double_t xarr[10] ;
    const RooAbsBinning* binning = realX->getBinningPtr(0) ;
    for (i=ybinLo ; i<=intOrder+ybinLo ; i++) {
      Int_t ibin ;
      if (i>=0 && i<ybinM) {
	// In range
	ibin = i ;
	realY->setBin(ibin) ;
	xarr[i-ybinLo] = realY->getVal() ;
      } else if (i>=ybinM) {
	// Overflow: mirror
	ibin = 2*ybinM-i-1 ;
	realY->setBin(ibin) ;
	xarr[i-ybinLo] = 2*realY->getMax()-realY->getVal() ;
      } else {
	// Underflow: mirror
	ibin = -i -1;
	realY->setBin(ibin) ;
	xarr[i-ybinLo] = 2*realY->getMin()-realY->getVal() ;
      }
      yarr[i-ybinLo] = interpolateDim(*realX,binning,xval,intOrder,correctForBinSize,kFALSE) ;	
    }

    if (gDebug>7) {
      cout << "RooDataHist interpolating data is" << endl ;
      cout << "xarr = " ;
      for (int q=0; q<=intOrder ; q++) cout << xarr[q] << " " ;
      cout << " yarr = " ;
      for (int q=0; q<=intOrder ; q++) cout << yarr[q] << " " ;
      cout << endl ;
    }
    wInt = RooMath::interpolate(xarr,yarr,intOrder+1,yval) ;
    
  } else {

    // Higher dimensional scenarios not yet implemented
    coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") interpolation in " 
	 << _realVars.getSize() << " dimensions not yet implemented" << endl ;
    return weight(bin,0) ;

  }

  // Cut off negative values
//   if (wInt<=0) {
//     wInt=0 ; 
//   }

  //cout << "RooDataHist wInt = " << wInt << endl ;
  return wInt ;
}




//_____________________________________________________________________________
void RooDataHist::weightError(Double_t& lo, Double_t& hi, ErrorType etype) const 
{ 
  // Return the error on current weight
  checkInit() ;

  switch (etype) {

  case Auto:
    throw string(Form("RooDataHist::weightError(%s) error type Auto not allowed here",GetName())) ;
    break ;

  case Expected:
    throw string(Form("RooDataHist::weightError(%s) error type Expected not allowed here",GetName())) ;
    break ;

  case Poisson:
    if (_curWgtErrLo>=0) {
      // Weight is preset or precalculated    
      lo = _curWgtErrLo ;
      hi = _curWgtErrHi ;
      return ;
    }
    
    // Calculate poisson errors
    Double_t ym,yp ;  
    RooHistError::instance().getPoissonInterval(Int_t(weight()+0.5),ym,yp,1) ;
    _curWgtErrLo = weight()-ym ;
    _curWgtErrHi = yp-weight() ;
    _errLo[_curIndex] = _curWgtErrLo ;
    _errHi[_curIndex] = _curWgtErrHi ;
    lo = _curWgtErrLo ;
    hi = _curWgtErrHi ;
    return ;

  case SumW2:
    lo = sqrt(_curSumW2) ;
    hi = sqrt(_curSumW2) ;
    return ;

  case None:
    lo = 0 ;
    hi = 0 ;
    return ;
  }
}


// wve adjust for variable bin sizes

//_____________________________________________________________________________
Double_t RooDataHist::interpolateDim(RooRealVar& dim, const RooAbsBinning* binning, Double_t xval, Int_t intOrder, Bool_t correctForBinSize, Bool_t cdfBoundaries) 
{
  // Perform boundary safe 'intOrder'-th interpolation of weights in dimension 'dim'
  // at current value 'xval'
 
  // Fill workspace arrays spanning interpolation area
  Int_t fbinC = dim.getBin(*binning) ;
  Int_t fbinLo = fbinC-intOrder/2 - ((xval<binning->binCenter(fbinC))?1:0) ;
  Int_t fbinM = dim.numBins(*binning) ;


  Int_t i ;
  Double_t yarr[10] ;
  Double_t xarr[10] ;
  for (i=fbinLo ; i<=intOrder+fbinLo ; i++) {
    Int_t ibin ;
    if (i>=0 && i<fbinM) {
      // In range
      ibin = i ;
      dim.setBinFast(ibin,*binning) ;
      //cout << "INRANGE: dim.getVal(ibin=" << ibin << ") = " << dim.getVal() << endl ;
      xarr[i-fbinLo] = dim.getVal() ;
      Int_t idx = calcTreeIndex() ;      
      yarr[i-fbinLo] = _wgt[idx] ; 
      if (correctForBinSize) yarr[i-fbinLo] /=  _binv[idx] ;
    } else if (i>=fbinM) {
      // Overflow: mirror
      ibin = 2*fbinM-i-1 ;
      dim.setBinFast(ibin,*binning) ;
      //cout << "OVERFLOW: dim.getVal(ibin=" << ibin << ") = " << dim.getVal() << endl ;
      if (cdfBoundaries) {	
	xarr[i-fbinLo] = dim.getMax()+1e-10*(i-fbinM+1) ;
	yarr[i-fbinLo] = 1.0 ;
      } else {
	Int_t idx = calcTreeIndex() ;      
	xarr[i-fbinLo] = 2*dim.getMax()-dim.getVal() ;
	yarr[i-fbinLo] = _wgt[idx] ; 
	if (correctForBinSize) yarr[i-fbinLo] /=  _binv[idx] ;
      }
    } else {
      // Underflow: mirror
      ibin = -i - 1 ;
      dim.setBinFast(ibin,*binning) ;
      //cout << "UNDERFLOW: dim.getVal(ibin=" << ibin << ") = " << dim.getVal() << endl ;
      if (cdfBoundaries) {
	xarr[i-fbinLo] = dim.getMin()-ibin*(1e-10) ; ;
	yarr[i-fbinLo] = 0.0 ;
      } else {
	Int_t idx = calcTreeIndex() ;      
	xarr[i-fbinLo] = 2*dim.getMin()-dim.getVal() ;
	yarr[i-fbinLo] = _wgt[idx] ; 
	if (correctForBinSize) yarr[i-fbinLo] /=  _binv[idx] ;
      }
    }
    //cout << "ibin = " << ibin << endl ;
  }
//   for (int k=0 ; k<=intOrder ; k++) {
//     cout << "k=" << k << " x = " << xarr[k] << " y = " << yarr[k] << endl ;
//   }
  dim.setBinFast(fbinC,*binning) ;
  Double_t ret = RooMath::interpolate(xarr,yarr,intOrder+1,xval) ;
  return ret ;
}




//_____________________________________________________________________________
void RooDataHist::add(const RooArgSet& row, Double_t wgt, Double_t sumw2) 
{
  // Increment the weight of the bin enclosing the coordinates given
  // by 'row' by the specified amount. Add the sum of weights squared
  // for the bin by 'sumw2' rather than wgt^2

  checkInit() ;

//   cout << "RooDataHist::add() accepted coordinate is " << endl ;
//   _vars.Print("v") ;

  _vars = row ;
  Int_t idx = calcTreeIndex() ;
  _wgt[idx] += wgt ;  
  _sumw2[idx] += (sumw2>0?sumw2:wgt*wgt) ;
  _errLo[idx] = -1 ;
  _errHi[idx] = -1 ;

  _cache_sum_valid = kFALSE ;
}



//_____________________________________________________________________________
void RooDataHist::set(const RooArgSet& row, Double_t wgt, Double_t wgtErrLo, Double_t wgtErrHi) 
{
  // Increment the weight of the bin enclosing the coordinates
  // given by 'row' by the specified amount. Associate errors
  // [wgtErrLo,wgtErrHi] with the event weight on this bin.

  checkInit() ;

  _vars = row ;
  Int_t idx = calcTreeIndex() ;
  _wgt[idx] = wgt ;  
  _errLo[idx] = wgtErrLo ;  
  _errHi[idx] = wgtErrHi ;  

  _cache_sum_valid = kFALSE ;
}



//_____________________________________________________________________________
void RooDataHist::set(Double_t wgt, Double_t wgtErr) 
{
  // Increment the weight of the bin enclosing the coordinates
  // given by 'row' by the specified amount. Associate errors
  // [wgtErrLo,wgtErrHi] with the event weight on this bin.

  checkInit() ;
  
  if (_curIndex<0) {
    _curIndex = calcTreeIndex() ;
  }

  _wgt[_curIndex] = wgt ;  
  _errLo[_curIndex] = wgtErr ;  
  _errHi[_curIndex] = wgtErr ;  
  _sumw2[_curIndex] = wgtErr*wgtErr ;

  _cache_sum_valid = kFALSE ;
}



//_____________________________________________________________________________
void RooDataHist::set(const RooArgSet& row, Double_t wgt, Double_t wgtErr) 
{
  // Increment the weight of the bin enclosing the coordinates
  // given by 'row' by the specified amount. Associate errors
  // [wgtErrLo,wgtErrHi] with the event weight on this bin.

  checkInit() ;

  _vars = row ;
  Int_t idx = calcTreeIndex() ;
  _wgt[idx] = wgt ;  
  _errLo[idx] = wgtErr ;  
  _errHi[idx] = wgtErr ;  
  _sumw2[idx] = wgtErr*wgtErr ;

  _cache_sum_valid = kFALSE ;
}



//_____________________________________________________________________________
void RooDataHist::add(const RooAbsData& dset, const char* cut, Double_t wgt) 
{  
  // Add all data points contained in 'dset' to this data set with given weight.
  // Optional cut string expression selects the data points to be added and can
  // reference any variable contained in this data set

  RooFormulaVar cutVar("select",cut,*dset.get()) ;
  add(dset,&cutVar,wgt) ;
}



//_____________________________________________________________________________
void RooDataHist::add(const RooAbsData& dset, const RooFormulaVar* cutVar, Double_t wgt) 
{
  // Add all data points contained in 'dset' to this data set with given weight.
  // Optional RooFormulaVar pointer selects the data points to be added.

  checkInit() ;

  RooFormulaVar* cloneVar = 0;
  RooArgSet* tmp(0) ;
  if (cutVar) {
    // Deep clone cutVar and attach clone to this dataset
    tmp = (RooArgSet*) RooArgSet(*cutVar).snapshot() ;
    if (!tmp) {
      coutE(DataHandling) << "RooDataHist::add(" << GetName() << ") Couldn't deep-clone cut variable, abort," << endl ;
      return ;
    }

    cloneVar = (RooFormulaVar*) tmp->find(*cutVar) ;
    cloneVar->attachDataSet(dset) ;
  }


  Int_t i ;
  for (i=0 ; i<dset.numEntries() ; i++) {
    const RooArgSet* row = dset.get(i) ;
    if (!cloneVar || cloneVar->getVal()) {
       add(*row,wgt*dset.weight(), wgt*wgt*dset.weightSquared()) ;
    }
  }

  if (cloneVar) {
    delete tmp ;
  } 

  _cache_sum_valid = kFALSE ;
}



//_____________________________________________________________________________
Double_t RooDataHist::sum(Bool_t correctForBinSize, Bool_t inverseBinCor) const 
{
  // Return the sum of the weights of all hist bins.
  //
  // If correctForBinSize is specified, the sum of weights
  // is multiplied by the N-dimensional bin volume,
  // making the return value the integral over the function
  // represented by this histogram

  checkInit() ;

  // Check if result was cached
  Int_t cache_code = 1 + (correctForBinSize?1:0) + ((correctForBinSize&&inverseBinCor)?1:0) ;
  if (_cache_sum_valid==cache_code) {
    return _cache_sum ;
  }

  Int_t i ;
  Double_t total(0), carry(0);
  for (i=0 ; i<_arrSize ; i++) {
    
    Double_t theBinVolume = correctForBinSize ? (inverseBinCor ? 1/_binv[i] : _binv[i]) : 1.0 ;
    // cout << "total += " << _wgt[i] << "*" << theBinVolume << endl ;
    Double_t y = _wgt[i]*theBinVolume - carry;
    Double_t t = total + y;
    carry = (t - total) - y;
    total = t;
  }

  // Store result in cache
  _cache_sum_valid=cache_code ;
  _cache_sum = total ;

  return total ;
}



//_____________________________________________________________________________
Double_t RooDataHist::sum(const RooArgSet& sumSet, const RooArgSet& sliceSet, Bool_t correctForBinSize, Bool_t inverseBinCor)
{
  // Return the sum of the weights of a multi-dimensional slice of the histogram
  // by summing only over the dimensions specified in sumSet.
  //   
  // The coordinates of all other dimensions are fixed to those given in sliceSet
  //
  // If correctForBinSize is specified, the sum of weights
  // is multiplied by the M-dimensional bin volume, (M = N(sumSet)),
  // making the return value the integral over the function
  // represented by this histogram

  checkInit() ;

  RooArgSet varSave ;
  varSave.addClone(_vars) ;

  RooArgSet* sliceOnlySet = new RooArgSet(sliceSet) ;
  sliceOnlySet->remove(sumSet,kTRUE,kTRUE) ;

  _vars = *sliceOnlySet ;
  calculatePartialBinVolume(*sliceOnlySet) ;
  delete sliceOnlySet ;

  TIterator* ssIter = sumSet.createIterator() ;
  
  // Calculate mask and refence plot bins for non-iterating variables
  RooAbsArg* arg ;
  Bool_t* mask = new Bool_t[_vars.getSize()] ;
  Int_t*  refBin = new Int_t[_vars.getSize()] ;

  Int_t i(0) ;
  _iterator->Reset() ;
  while((arg=(RooAbsArg*)_iterator->Next())) {
    if (sumSet.find(*arg)) {
      mask[i] = kFALSE ;
    } else {
      mask[i] = kTRUE ;
      // coverity[FORWARD_NULL]
      refBin[i] = (dynamic_cast<RooAbsLValue*>(arg))->getBin() ;
    }
    i++ ;
  }
    
  // Loop over entire data set, skipping masked entries
  Double_t total(0), carry(0);
  Int_t ibin ;
  for (ibin=0 ; ibin<_arrSize ; ibin++) {

    Int_t idx(0), tmp(ibin), ivar(0) ;
    Bool_t skip(kFALSE) ;

    // Check if this bin belongs in selected slice
    _iterator->Reset() ;
    // coverity[UNUSED_VALUE]
    while((!skip && (arg=(RooAbsArg*)_iterator->Next()))) {
      idx  = tmp / _idxMult[ivar] ;
      tmp -= idx*_idxMult[ivar] ;
      if (mask[ivar] && idx!=refBin[ivar]) skip=kTRUE ;
      ivar++ ;
    }
    
    if (!skip) {
      Double_t theBinVolume = correctForBinSize ? (inverseBinCor ? 1/(*_pbinv)[i] : (*_pbinv)[i] ) : 1.0 ;
//       cout << "adding bin[" << ibin << "] to sum wgt = " << _wgt[ibin] << " binv = " << theBinVolume << endl ;
      Double_t y = _wgt[ibin]*theBinVolume - carry;
      Double_t t = total + y;
      carry = (t - total) - y;
      total = t;
    }
  }
  delete ssIter ;

  delete[] mask ;
  delete[] refBin ;

  _vars = varSave ;

  return total ;
}

//_____________________________________________________________________________
Double_t RooDataHist::sum(const RooArgSet& sumSet, const RooArgSet& sliceSet,
	Bool_t correctForBinSize, Bool_t inverseBinCor,
	const std::map<const RooAbsArg*, std::pair<Double_t, Double_t> >& ranges)
{
  // Return the sum of the weights of a multi-dimensional slice of the histogram
  // by summing only over the dimensions specified in sumSet.
  //   
  // The coordinates of all other dimensions are fixed to those given in sliceSet
  //
  // If correctForBinSize is specified, the sum of weights
  // is multiplied by the M-dimensional bin volume, (M = N(sumSet)),
  // or the fraction of it that falls inside the range rangeName,
  // making the return value the integral over the function
  // represented by this histogram
  //
  // If correctForBinSize is not specified, the weights are multiplied by the
  // fraction of the bin volume that falls inside the range, i.e. a factor or
  // binVolumeInRange/totalBinVolume.

  checkInit();
  checkBinBounds();
  RooArgSet varSave;
  varSave.addClone(_vars);
  {
    RooArgSet sliceOnlySet(sliceSet);
    sliceOnlySet.remove(sumSet,kTRUE,kTRUE);
    _vars = sliceOnlySet;
  }

  // Calculate mask and refence plot bins for non-iterating variables,
  // and get ranges for iterating variables
  std::vector<bool> mask(_vars.getSize());
  std::vector<Int_t> refBin(_vars.getSize());
  std::vector<Double_t> rangeLo(_vars.getSize(), -std::numeric_limits<Double_t>::infinity());
  std::vector<Double_t> rangeHi(_vars.getSize(), +std::numeric_limits<Double_t>::infinity());

  _iterator->Reset();
  RooAbsArg* arg;
  for (Int_t i = 0; (arg=(RooAbsArg*)_iterator->Next()); ++i) {
    RooAbsArg* sumsetv = sumSet.find(*arg);
    RooAbsArg* slicesetv = sliceSet.find(*arg);
    mask[i] = !sumsetv;
    if (mask[i]) {
      // coverity[FORWARD_NULL]
      refBin[i] = (dynamic_cast<RooAbsLValue*>(arg))->getBin();
    }
    std::map<const RooAbsArg*, std::pair<Double_t, Double_t> >::const_iterator
	it = ranges.find(sumsetv ? sumsetv : slicesetv);
    if (ranges.end() != it) {
      rangeLo[i] = it->second.first;
      rangeHi[i] = it->second.second;
    }
  }

  // Loop over entire data set, skipping masked entries
  Double_t total(0), carry(0);
  for (Int_t ibin = 0; ibin < _arrSize; ++ibin) {
    // Check if this bin belongs in selected slice
    _iterator->Reset();
    Bool_t skip(kFALSE);
    // coverity[UNUSED_VALUE]
    for (Int_t ivar = 0, tmp = ibin;
	(!skip && (arg=(RooAbsArg*)_iterator->Next())); ++ivar) {
      const Int_t idx = tmp / _idxMult[ivar];
      tmp -= idx*_idxMult[ivar];
      if (mask[ivar] && idx!=refBin[ivar]) skip=kTRUE;
    }
    if (skip) continue;
    _iterator->Reset();
    // work out bin volume
    Double_t theBinVolume = 1.;
    for (Int_t ivar = 0, tmp = ibin;
	(arg=(RooAbsArg*)_iterator->Next()); ++ivar) {
      const Int_t idx = tmp / _idxMult[ivar];
      tmp -= idx*_idxMult[ivar];
      if (_binbounds[ivar].empty()) continue;
      const Double_t binLo = _binbounds[ivar][2 * idx];
      const Double_t binHi = _binbounds[ivar][2 * idx + 1];
      if (binHi < rangeLo[ivar] || binLo > rangeHi[ivar]) {
	// bin is outside of allowed range - effective bin volume is zero
	theBinVolume = 0.;
	break;
      }
      theBinVolume *= 
	(std::min(rangeHi[ivar], binHi) - std::max(rangeLo[ivar], binLo));
    }
    const Double_t corrPartial = theBinVolume / _binv[ibin];
    if (0. == corrPartial) continue;
    const Double_t corr = correctForBinSize ? (inverseBinCor ? 1. / _binv[ibin] : _binv[ibin] ) : 1.0;
    //cout << "adding bin[" << ibin << "] to sum wgt = " << _wgt[ibin] << " binv = " << theBinVolume << " _binv[" << ibin << "] " << _binv[ibin] << endl;
    
    const Double_t y = _wgt[ibin] * corr * corrPartial - carry;
    const Double_t t = total + y;
    carry = (t - total) - y;
    total = t;
  }

  _vars = varSave;

  return total;
}



//_____________________________________________________________________________
void RooDataHist::calculatePartialBinVolume(const RooArgSet& dimSet) const 
{
  // Fill the transient cache with partial bin volumes with up-to-date
  // values for the partial volume specified by observables 'dimSet'

  // Allocate cache if not yet existing
  vector<Double_t> *pbinv = _pbinvCacheMgr.getObj(&dimSet) ;
  if (pbinv) {
    _pbinv = pbinv ;
    return ;
  }

  pbinv = new vector<Double_t>(_arrSize) ;

  // Calculate plot bins of components from master index
  Bool_t* selDim = new Bool_t[_vars.getSize()] ;
  _iterator->Reset() ;
  RooAbsArg* v ;
  Int_t i(0) ;
  while((v=(RooAbsArg*)_iterator->Next())) {
    selDim[i++] = dimSet.find(*v) ? kTRUE : kFALSE ;
  }

  // Recalculate partial bin volume cache
  Int_t ibin ;
  for (ibin=0 ; ibin<_arrSize ; ibin++) {
    _iterator->Reset() ;
    RooAbsLValue* arg ;
    Int_t j(0), idx(0), tmp(ibin) ;
    Double_t theBinVolume(1) ;
    while((arg=dynamic_cast<RooAbsLValue*>(_iterator->Next()))) {
      idx  = tmp / _idxMult[j] ;
      tmp -= idx*_idxMult[j++] ;
      if (selDim[j-1]) {
	RooAbsLValue* arglv = dynamic_cast<RooAbsLValue*>(arg) ;
	theBinVolume *= arglv->getBinWidth(idx) ;
      }
    }
    (*pbinv)[ibin] = theBinVolume ;
  }

  delete[] selDim ;

  // Put in cache (which takes ownership) 
  _pbinvCacheMgr.setObj(&dimSet,pbinv) ;

  // Publicize the array
  _pbinv = pbinv ;
}



//_____________________________________________________________________________
Int_t RooDataHist::numEntries() const 
{
  // Return the number of bins

  return RooAbsData::numEntries() ;
}



//_____________________________________________________________________________
Double_t RooDataHist::sumEntries() const 
{
  Int_t i ;
  Double_t n(0), carry(0);
  for (i=0 ; i<_arrSize ; i++) {
    if (!_binValid || _binValid[i]) {
      Double_t y = _wgt[i] - carry;
      Double_t t = n + y;
      carry = (t - n) - y;
      n = t;
    }
  }
  return n ;
}



//_____________________________________________________________________________
Double_t RooDataHist::sumEntries(const char* cutSpec, const char* cutRange) const
{
  // Return the sum of weights in all entries matching cutSpec (if specified)
  // and in named range cutRange (if specified)
  // Return the
  checkInit() ;

  if (cutSpec==0 && cutRange==0) {
    return sumEntries();
  } else {
    
    // Setup RooFormulaVar for cutSpec if it is present
    RooFormula* select = 0 ;
    if (cutSpec) {
      select = new RooFormula("select",cutSpec,*get()) ;
    }
    
    // Otherwise sum the weights in the event
    Double_t sumw(0), carry(0);
    Int_t i ;
    for (i=0 ; i<numEntries() ; i++) {
      get(i) ;
      if (select && select->eval()==0.) continue ;
      if (cutRange && !_vars.allInRange(cutRange)) continue ;

      if (!_binValid || _binValid[i]) {
	Double_t y = weight() - carry;
	Double_t t = sumw + y;
	carry = (t - sumw) - y;
	sumw = t;
      }
    }
    
    if (select) delete select ;
    
    return sumw ;
  }
}



//_____________________________________________________________________________
void RooDataHist::reset() 
{
  // Reset all bin weights to zero

  // WVE DO NOT CALL RooTreeData::reset() for binned
  // datasets as this will delete the bin definitions

  Int_t i ;
  for (i=0 ; i<_arrSize ; i++) {
    _wgt[i] = 0. ;
    _errLo[i] = -1 ;
    _errHi[i] = -1 ;
  }
  _curWeight = 0 ;
  _curWgtErrLo = -1 ;
  _curWgtErrHi = -1 ;
  _curVolume = 1 ;

  _cache_sum_valid = kFALSE ;

}



//_____________________________________________________________________________
const RooArgSet* RooDataHist::get(Int_t masterIdx) const  
{
  // Return an argset with the bin center coordinates for 
  // bin sequential number 'masterIdx'. For iterative use.
  checkInit() ;
  _curWeight = _wgt[masterIdx] ;
  _curWgtErrLo = _errLo[masterIdx] ;
  _curWgtErrHi = _errHi[masterIdx] ;
  _curSumW2 = _sumw2[masterIdx] ;
  _curVolume = _binv[masterIdx] ; 
  _curIndex  = masterIdx ;
  return RooAbsData::get(masterIdx) ;  
}



//_____________________________________________________________________________
const RooArgSet* RooDataHist::get(const RooArgSet& coord) const
{
  // Return a RooArgSet with center coordinates of the bin
  // enclosing the point 'coord'

  ((RooDataHist*)this)->_vars = coord ;
  return get(calcTreeIndex()) ;
}



//_____________________________________________________________________________
Double_t RooDataHist::binVolume(const RooArgSet& coord) 
{
  // Return the volume of the bin enclosing coordinates 'coord'
  checkInit() ;
  ((RooDataHist*)this)->_vars = coord ;
  return _binv[calcTreeIndex()] ;
}


//_____________________________________________________________________________
void RooDataHist::setAllWeights(Double_t value) 
{
  // Set all the event weight of all bins to the specified value
  for (Int_t i=0 ; i<_arrSize ; i++) {
    _wgt[i] = value ;
  }

  _cache_sum_valid = kFALSE ;
}



//_____________________________________________________________________________
TIterator* RooDataHist::sliceIterator(RooAbsArg& sliceArg, const RooArgSet& otherArgs) 
{
  // Create an iterator over all bins in a slice defined by the subset of observables
  // listed in sliceArg. The position of the slice is given by otherArgs

  // Update to current position
  _vars = otherArgs ;
  _curIndex = calcTreeIndex() ;
  
  RooAbsArg* intArg = _vars.find(sliceArg) ;
  if (!intArg) {
    coutE(InputArguments) << "RooDataHist::sliceIterator() variable " << sliceArg.GetName() << " is not part of this RooDataHist" << endl ;
    return 0 ;
  }
  return new RooDataHistSliceIter(*this,*intArg) ;
}


//_____________________________________________________________________________
void RooDataHist::SetName(const char *name) 
{
  // Change the name of the RooDataHist

  if (_dir) _dir->GetList()->Remove(this);
  TNamed::SetName(name) ;
  if (_dir) _dir->GetList()->Add(this);
}


//_____________________________________________________________________________
void RooDataHist::SetNameTitle(const char *name, const char* title) 
{
  // Change the title of this RooDataHist

  if (_dir) _dir->GetList()->Remove(this);
  TNamed::SetNameTitle(name,title) ;
  if (_dir) _dir->GetList()->Add(this);
}


//_____________________________________________________________________________
void RooDataHist::printValue(ostream& os) const 
{
  // Print value of the dataset, i.e. the sum of weights contained in the dataset
  os << numEntries() << " bins (" << sumEntries() << " weights)" ;
}




//_____________________________________________________________________________
void RooDataHist::printArgs(ostream& os) const 
{
  // Print argument of dataset, i.e. the observable names

  os << "[" ;    
  _iterator->Reset() ;
  RooAbsArg* arg ;
  Bool_t first(kTRUE) ;
  while((arg=(RooAbsArg*)_iterator->Next())) {
    if (first) {
      first=kFALSE ;
    } else {
      os << "," ;
    }
    os << arg->GetName() ;
  }
  os << "]" ;
}



//_____________________________________________________________________________
void RooDataHist::cacheValidEntries() 
{
  // Cache the datahist entries with bin centers that are inside/outside the
  // current observable definitio
  checkInit() ;

  if (!_binValid) {
    _binValid = new Bool_t[_arrSize] ;
  }
  TIterator* iter = _vars.createIterator() ;
  RooAbsArg* arg ;
  for (Int_t i=0 ; i<_arrSize ; i++) {
    get(i) ;
    _binValid[i] = kTRUE ;
    iter->Reset() ;
    while((arg=(RooAbsArg*)iter->Next())) {
      // coverity[CHECKED_RETURN]
      _binValid[i] &= arg->inRange(0) ;      
    }
  }
  delete iter ;

}


//_____________________________________________________________________________
Bool_t RooDataHist::valid() const 
{
  // Return true if currently loaded coordinate is considered valid within
  // the current range definitions of all observables

  // If caching is enabled, use the precached result
  if (_binValid) {
    return _binValid[_curIndex] ;
  }

  return kTRUE ;
}



//_____________________________________________________________________________
Bool_t RooDataHist::isNonPoissonWeighted() const
{
  // Returns true if datasets contains entries with a non-integer weight

  for (int i=0 ; i<numEntries() ; i++) {
    if (fabs(_wgt[i]-Int_t(_wgt[i]))>1e-10) return kTRUE ;
  }
  return kFALSE ;
}




//_____________________________________________________________________________
void RooDataHist::printMultiline(ostream& os, Int_t content, Bool_t verbose, TString indent) const 
{
  // Print the details on the dataset contents

  RooAbsData::printMultiline(os,content,verbose,indent) ;  

  os << indent << "Binned Dataset " << GetName() << " (" << GetTitle() << ")" << endl ;
  os << indent << "  Contains " << numEntries() << " bins with a total weight of " << sumEntries() << endl;
  
  if (!verbose) {
    os << indent << "  Observables " << _vars << endl ;
  } else {
    os << indent << "  Observables: " ;
    _vars.printStream(os,kName|kValue|kExtras|kTitle,kVerbose,indent+"  ") ;
  }

  if(verbose) {
    if (_cachedVars.getSize()>0) {
      os << indent << "  Caches " << _cachedVars << endl ;
    }
  }
}



//______________________________________________________________________________
void RooDataHist::Streamer(TBuffer &R__b)
{
   // Stream an object of class RooDataHist.

   if (R__b.IsReading()) {

     UInt_t R__s, R__c;
     Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
     
      if (R__v>2) {

	R__b.ReadClassBuffer(RooDataHist::Class(),this,R__v,R__s,R__c);
	initialize(0,kFALSE) ;

      } else {	

	// Legacy dataset conversion happens here. Legacy RooDataHist inherits from RooTreeData
	// which in turn inherits from RooAbsData. Manually stream RooTreeData contents on 
	// file here and convert it into a RooTreeDataStore which is installed in the 
	// new-style RooAbsData base class
	
	// --- This is the contents of the streamer code of RooTreeData version 2 ---
	UInt_t R__s1, R__c1;
	Version_t R__v1 = R__b.ReadVersion(&R__s1, &R__c1); if (R__v1) { }
	
	RooAbsData::Streamer(R__b);
	TTree* X_tree(0) ; R__b >> X_tree;
	RooArgSet X_truth ; X_truth.Streamer(R__b);
	TString X_blindString ; X_blindString.Streamer(R__b);
	R__b.CheckByteCount(R__s1, R__c1, RooTreeData::Class());
	// --- End of RooTreeData-v1 streamer
	
	// Construct RooTreeDataStore from X_tree and complete initialization of new-style RooAbsData
	_dstore = new RooTreeDataStore(X_tree,_vars) ;
	_dstore->SetName(GetName()) ;
	_dstore->SetTitle(GetTitle()) ;
	_dstore->checkInit() ;       
	
	RooDirItem::Streamer(R__b);
	R__b >> _arrSize;
	delete [] _wgt;
	_wgt = new Double_t[_arrSize];
	R__b.ReadFastArray(_wgt,_arrSize);
	delete [] _errLo;
	_errLo = new Double_t[_arrSize];
	R__b.ReadFastArray(_errLo,_arrSize);
	delete [] _errHi;
	_errHi = new Double_t[_arrSize];
	R__b.ReadFastArray(_errHi,_arrSize);
	delete [] _sumw2;
	_sumw2 = new Double_t[_arrSize];
	R__b.ReadFastArray(_sumw2,_arrSize);
	delete [] _binv;
	_binv = new Double_t[_arrSize];
	R__b.ReadFastArray(_binv,_arrSize);
	_realVars.Streamer(R__b);
	R__b >> _curWeight;
	R__b >> _curWgtErrLo;
	R__b >> _curWgtErrHi;
	R__b >> _curSumW2;
	R__b >> _curVolume;
	R__b >> _curIndex;
	R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());

      }
      
   } else {

      R__b.WriteClassBuffer(RooDataHist::Class(),this);
   }
}


 RooDataHist.cxx:1
 RooDataHist.cxx:2
 RooDataHist.cxx:3
 RooDataHist.cxx:4
 RooDataHist.cxx:5
 RooDataHist.cxx:6
 RooDataHist.cxx:7
 RooDataHist.cxx:8
 RooDataHist.cxx:9
 RooDataHist.cxx:10
 RooDataHist.cxx:11
 RooDataHist.cxx:12
 RooDataHist.cxx:13
 RooDataHist.cxx:14
 RooDataHist.cxx:15
 RooDataHist.cxx:16
 RooDataHist.cxx:17
 RooDataHist.cxx:18
 RooDataHist.cxx:19
 RooDataHist.cxx:20
 RooDataHist.cxx:21
 RooDataHist.cxx:22
 RooDataHist.cxx:23
 RooDataHist.cxx:24
 RooDataHist.cxx:25
 RooDataHist.cxx:26
 RooDataHist.cxx:27
 RooDataHist.cxx:28
 RooDataHist.cxx:29
 RooDataHist.cxx:30
 RooDataHist.cxx:31
 RooDataHist.cxx:32
 RooDataHist.cxx:33
 RooDataHist.cxx:34
 RooDataHist.cxx:35
 RooDataHist.cxx:36
 RooDataHist.cxx:37
 RooDataHist.cxx:38
 RooDataHist.cxx:39
 RooDataHist.cxx:40
 RooDataHist.cxx:41
 RooDataHist.cxx:42
 RooDataHist.cxx:43
 RooDataHist.cxx:44
 RooDataHist.cxx:45
 RooDataHist.cxx:46
 RooDataHist.cxx:47
 RooDataHist.cxx:48
 RooDataHist.cxx:49
 RooDataHist.cxx:50
 RooDataHist.cxx:51
 RooDataHist.cxx:52
 RooDataHist.cxx:53
 RooDataHist.cxx:54
 RooDataHist.cxx:55
 RooDataHist.cxx:56
 RooDataHist.cxx:57
 RooDataHist.cxx:58
 RooDataHist.cxx:59
 RooDataHist.cxx:60
 RooDataHist.cxx:61
 RooDataHist.cxx:62
 RooDataHist.cxx:63
 RooDataHist.cxx:64
 RooDataHist.cxx:65
 RooDataHist.cxx:66
 RooDataHist.cxx:67
 RooDataHist.cxx:68
 RooDataHist.cxx:69
 RooDataHist.cxx:70
 RooDataHist.cxx:71
 RooDataHist.cxx:72
 RooDataHist.cxx:73
 RooDataHist.cxx:74
 RooDataHist.cxx:75
 RooDataHist.cxx:76
 RooDataHist.cxx:77
 RooDataHist.cxx:78
 RooDataHist.cxx:79
 RooDataHist.cxx:80
 RooDataHist.cxx:81
 RooDataHist.cxx:82
 RooDataHist.cxx:83
 RooDataHist.cxx:84
 RooDataHist.cxx:85
 RooDataHist.cxx:86
 RooDataHist.cxx:87
 RooDataHist.cxx:88
 RooDataHist.cxx:89
 RooDataHist.cxx:90
 RooDataHist.cxx:91
 RooDataHist.cxx:92
 RooDataHist.cxx:93
 RooDataHist.cxx:94
 RooDataHist.cxx:95
 RooDataHist.cxx:96
 RooDataHist.cxx:97
 RooDataHist.cxx:98
 RooDataHist.cxx:99
 RooDataHist.cxx:100
 RooDataHist.cxx:101
 RooDataHist.cxx:102
 RooDataHist.cxx:103
 RooDataHist.cxx:104
 RooDataHist.cxx:105
 RooDataHist.cxx:106
 RooDataHist.cxx:107
 RooDataHist.cxx:108
 RooDataHist.cxx:109
 RooDataHist.cxx:110
 RooDataHist.cxx:111
 RooDataHist.cxx:112
 RooDataHist.cxx:113
 RooDataHist.cxx:114
 RooDataHist.cxx:115
 RooDataHist.cxx:116
 RooDataHist.cxx:117
 RooDataHist.cxx:118
 RooDataHist.cxx:119
 RooDataHist.cxx:120
 RooDataHist.cxx:121
 RooDataHist.cxx:122
 RooDataHist.cxx:123
 RooDataHist.cxx:124
 RooDataHist.cxx:125
 RooDataHist.cxx:126
 RooDataHist.cxx:127
 RooDataHist.cxx:128
 RooDataHist.cxx:129
 RooDataHist.cxx:130
 RooDataHist.cxx:131
 RooDataHist.cxx:132
 RooDataHist.cxx:133
 RooDataHist.cxx:134
 RooDataHist.cxx:135
 RooDataHist.cxx:136
 RooDataHist.cxx:137
 RooDataHist.cxx:138
 RooDataHist.cxx:139
 RooDataHist.cxx:140
 RooDataHist.cxx:141
 RooDataHist.cxx:142
 RooDataHist.cxx:143
 RooDataHist.cxx:144
 RooDataHist.cxx:145
 RooDataHist.cxx:146
 RooDataHist.cxx:147
 RooDataHist.cxx:148
 RooDataHist.cxx:149
 RooDataHist.cxx:150
 RooDataHist.cxx:151
 RooDataHist.cxx:152
 RooDataHist.cxx:153
 RooDataHist.cxx:154
 RooDataHist.cxx:155
 RooDataHist.cxx:156
 RooDataHist.cxx:157
 RooDataHist.cxx:158
 RooDataHist.cxx:159
 RooDataHist.cxx:160
 RooDataHist.cxx:161
 RooDataHist.cxx:162
 RooDataHist.cxx:163
 RooDataHist.cxx:164
 RooDataHist.cxx:165
 RooDataHist.cxx:166
 RooDataHist.cxx:167
 RooDataHist.cxx:168
 RooDataHist.cxx:169
 RooDataHist.cxx:170
 RooDataHist.cxx:171
 RooDataHist.cxx:172
 RooDataHist.cxx:173
 RooDataHist.cxx:174
 RooDataHist.cxx:175
 RooDataHist.cxx:176
 RooDataHist.cxx:177
 RooDataHist.cxx:178
 RooDataHist.cxx:179
 RooDataHist.cxx:180
 RooDataHist.cxx:181
 RooDataHist.cxx:182
 RooDataHist.cxx:183
 RooDataHist.cxx:184
 RooDataHist.cxx:185
 RooDataHist.cxx:186
 RooDataHist.cxx:187
 RooDataHist.cxx:188
 RooDataHist.cxx:189
 RooDataHist.cxx:190
 RooDataHist.cxx:191
 RooDataHist.cxx:192
 RooDataHist.cxx:193
 RooDataHist.cxx:194
 RooDataHist.cxx:195
 RooDataHist.cxx:196
 RooDataHist.cxx:197
 RooDataHist.cxx:198
 RooDataHist.cxx:199
 RooDataHist.cxx:200
 RooDataHist.cxx:201
 RooDataHist.cxx:202
 RooDataHist.cxx:203
 RooDataHist.cxx:204
 RooDataHist.cxx:205
 RooDataHist.cxx:206
 RooDataHist.cxx:207
 RooDataHist.cxx:208
 RooDataHist.cxx:209
 RooDataHist.cxx:210
 RooDataHist.cxx:211
 RooDataHist.cxx:212
 RooDataHist.cxx:213
 RooDataHist.cxx:214
 RooDataHist.cxx:215
 RooDataHist.cxx:216
 RooDataHist.cxx:217
 RooDataHist.cxx:218
 RooDataHist.cxx:219
 RooDataHist.cxx:220
 RooDataHist.cxx:221
 RooDataHist.cxx:222
 RooDataHist.cxx:223
 RooDataHist.cxx:224
 RooDataHist.cxx:225
 RooDataHist.cxx:226
 RooDataHist.cxx:227
 RooDataHist.cxx:228
 RooDataHist.cxx:229
 RooDataHist.cxx:230
 RooDataHist.cxx:231
 RooDataHist.cxx:232
 RooDataHist.cxx:233
 RooDataHist.cxx:234
 RooDataHist.cxx:235
 RooDataHist.cxx:236
 RooDataHist.cxx:237
 RooDataHist.cxx:238
 RooDataHist.cxx:239
 RooDataHist.cxx:240
 RooDataHist.cxx:241
 RooDataHist.cxx:242
 RooDataHist.cxx:243
 RooDataHist.cxx:244
 RooDataHist.cxx:245
 RooDataHist.cxx:246
 RooDataHist.cxx:247
 RooDataHist.cxx:248
 RooDataHist.cxx:249
 RooDataHist.cxx:250
 RooDataHist.cxx:251
 RooDataHist.cxx:252
 RooDataHist.cxx:253
 RooDataHist.cxx:254
 RooDataHist.cxx:255
 RooDataHist.cxx:256
 RooDataHist.cxx:257
 RooDataHist.cxx:258
 RooDataHist.cxx:259
 RooDataHist.cxx:260
 RooDataHist.cxx:261
 RooDataHist.cxx:262
 RooDataHist.cxx:263
 RooDataHist.cxx:264
 RooDataHist.cxx:265
 RooDataHist.cxx:266
 RooDataHist.cxx:267
 RooDataHist.cxx:268
 RooDataHist.cxx:269
 RooDataHist.cxx:270
 RooDataHist.cxx:271
 RooDataHist.cxx:272
 RooDataHist.cxx:273
 RooDataHist.cxx:274
 RooDataHist.cxx:275
 RooDataHist.cxx:276
 RooDataHist.cxx:277
 RooDataHist.cxx:278
 RooDataHist.cxx:279
 RooDataHist.cxx:280
 RooDataHist.cxx:281
 RooDataHist.cxx:282
 RooDataHist.cxx:283
 RooDataHist.cxx:284
 RooDataHist.cxx:285
 RooDataHist.cxx:286
 RooDataHist.cxx:287
 RooDataHist.cxx:288
 RooDataHist.cxx:289
 RooDataHist.cxx:290
 RooDataHist.cxx:291
 RooDataHist.cxx:292
 RooDataHist.cxx:293
 RooDataHist.cxx:294
 RooDataHist.cxx:295
 RooDataHist.cxx:296
 RooDataHist.cxx:297
 RooDataHist.cxx:298
 RooDataHist.cxx:299
 RooDataHist.cxx:300
 RooDataHist.cxx:301
 RooDataHist.cxx:302
 RooDataHist.cxx:303
 RooDataHist.cxx:304
 RooDataHist.cxx:305
 RooDataHist.cxx:306
 RooDataHist.cxx:307
 RooDataHist.cxx:308
 RooDataHist.cxx:309
 RooDataHist.cxx:310
 RooDataHist.cxx:311
 RooDataHist.cxx:312
 RooDataHist.cxx:313
 RooDataHist.cxx:314
 RooDataHist.cxx:315
 RooDataHist.cxx:316
 RooDataHist.cxx:317
 RooDataHist.cxx:318
 RooDataHist.cxx:319
 RooDataHist.cxx:320
 RooDataHist.cxx:321
 RooDataHist.cxx:322
 RooDataHist.cxx:323
 RooDataHist.cxx:324
 RooDataHist.cxx:325
 RooDataHist.cxx:326
 RooDataHist.cxx:327
 RooDataHist.cxx:328
 RooDataHist.cxx:329
 RooDataHist.cxx:330
 RooDataHist.cxx:331
 RooDataHist.cxx:332
 RooDataHist.cxx:333
 RooDataHist.cxx:334
 RooDataHist.cxx:335
 RooDataHist.cxx:336
 RooDataHist.cxx:337
 RooDataHist.cxx:338
 RooDataHist.cxx:339
 RooDataHist.cxx:340
 RooDataHist.cxx:341
 RooDataHist.cxx:342
 RooDataHist.cxx:343
 RooDataHist.cxx:344
 RooDataHist.cxx:345
 RooDataHist.cxx:346
 RooDataHist.cxx:347
 RooDataHist.cxx:348
 RooDataHist.cxx:349
 RooDataHist.cxx:350
 RooDataHist.cxx:351
 RooDataHist.cxx:352
 RooDataHist.cxx:353
 RooDataHist.cxx:354
 RooDataHist.cxx:355
 RooDataHist.cxx:356
 RooDataHist.cxx:357
 RooDataHist.cxx:358
 RooDataHist.cxx:359
 RooDataHist.cxx:360
 RooDataHist.cxx:361
 RooDataHist.cxx:362
 RooDataHist.cxx:363
 RooDataHist.cxx:364
 RooDataHist.cxx:365
 RooDataHist.cxx:366
 RooDataHist.cxx:367
 RooDataHist.cxx:368
 RooDataHist.cxx:369
 RooDataHist.cxx:370
 RooDataHist.cxx:371
 RooDataHist.cxx:372
 RooDataHist.cxx:373
 RooDataHist.cxx:374
 RooDataHist.cxx:375
 RooDataHist.cxx:376
 RooDataHist.cxx:377
 RooDataHist.cxx:378
 RooDataHist.cxx:379
 RooDataHist.cxx:380
 RooDataHist.cxx:381
 RooDataHist.cxx:382
 RooDataHist.cxx:383
 RooDataHist.cxx:384
 RooDataHist.cxx:385
 RooDataHist.cxx:386
 RooDataHist.cxx:387
 RooDataHist.cxx:388
 RooDataHist.cxx:389
 RooDataHist.cxx:390
 RooDataHist.cxx:391
 RooDataHist.cxx:392
 RooDataHist.cxx:393
 RooDataHist.cxx:394
 RooDataHist.cxx:395
 RooDataHist.cxx:396
 RooDataHist.cxx:397
 RooDataHist.cxx:398
 RooDataHist.cxx:399
 RooDataHist.cxx:400
 RooDataHist.cxx:401
 RooDataHist.cxx:402
 RooDataHist.cxx:403
 RooDataHist.cxx:404
 RooDataHist.cxx:405
 RooDataHist.cxx:406
 RooDataHist.cxx:407
 RooDataHist.cxx:408
 RooDataHist.cxx:409
 RooDataHist.cxx:410
 RooDataHist.cxx:411
 RooDataHist.cxx:412
 RooDataHist.cxx:413
 RooDataHist.cxx:414
 RooDataHist.cxx:415
 RooDataHist.cxx:416
 RooDataHist.cxx:417
 RooDataHist.cxx:418
 RooDataHist.cxx:419
 RooDataHist.cxx:420
 RooDataHist.cxx:421
 RooDataHist.cxx:422
 RooDataHist.cxx:423
 RooDataHist.cxx:424
 RooDataHist.cxx:425
 RooDataHist.cxx:426
 RooDataHist.cxx:427
 RooDataHist.cxx:428
 RooDataHist.cxx:429
 RooDataHist.cxx:430
 RooDataHist.cxx:431
 RooDataHist.cxx:432
 RooDataHist.cxx:433
 RooDataHist.cxx:434
 RooDataHist.cxx:435
 RooDataHist.cxx:436
 RooDataHist.cxx:437
 RooDataHist.cxx:438
 RooDataHist.cxx:439
 RooDataHist.cxx:440
 RooDataHist.cxx:441
 RooDataHist.cxx:442
 RooDataHist.cxx:443
 RooDataHist.cxx:444
 RooDataHist.cxx:445
 RooDataHist.cxx:446
 RooDataHist.cxx:447
 RooDataHist.cxx:448
 RooDataHist.cxx:449
 RooDataHist.cxx:450
 RooDataHist.cxx:451
 RooDataHist.cxx:452
 RooDataHist.cxx:453
 RooDataHist.cxx:454
 RooDataHist.cxx:455
 RooDataHist.cxx:456
 RooDataHist.cxx:457
 RooDataHist.cxx:458
 RooDataHist.cxx:459
 RooDataHist.cxx:460
 RooDataHist.cxx:461
 RooDataHist.cxx:462
 RooDataHist.cxx:463
 RooDataHist.cxx:464
 RooDataHist.cxx:465
 RooDataHist.cxx:466
 RooDataHist.cxx:467
 RooDataHist.cxx:468
 RooDataHist.cxx:469
 RooDataHist.cxx:470
 RooDataHist.cxx:471
 RooDataHist.cxx:472
 RooDataHist.cxx:473
 RooDataHist.cxx:474
 RooDataHist.cxx:475
 RooDataHist.cxx:476
 RooDataHist.cxx:477
 RooDataHist.cxx:478
 RooDataHist.cxx:479
 RooDataHist.cxx:480
 RooDataHist.cxx:481
 RooDataHist.cxx:482
 RooDataHist.cxx:483
 RooDataHist.cxx:484
 RooDataHist.cxx:485
 RooDataHist.cxx:486
 RooDataHist.cxx:487
 RooDataHist.cxx:488
 RooDataHist.cxx:489
 RooDataHist.cxx:490
 RooDataHist.cxx:491
 RooDataHist.cxx:492
 RooDataHist.cxx:493
 RooDataHist.cxx:494
 RooDataHist.cxx:495
 RooDataHist.cxx:496
 RooDataHist.cxx:497
 RooDataHist.cxx:498
 RooDataHist.cxx:499
 RooDataHist.cxx:500
 RooDataHist.cxx:501
 RooDataHist.cxx:502
 RooDataHist.cxx:503
 RooDataHist.cxx:504
 RooDataHist.cxx:505
 RooDataHist.cxx:506
 RooDataHist.cxx:507
 RooDataHist.cxx:508
 RooDataHist.cxx:509
 RooDataHist.cxx:510
 RooDataHist.cxx:511
 RooDataHist.cxx:512
 RooDataHist.cxx:513
 RooDataHist.cxx:514
 RooDataHist.cxx:515
 RooDataHist.cxx:516
 RooDataHist.cxx:517
 RooDataHist.cxx:518
 RooDataHist.cxx:519
 RooDataHist.cxx:520
 RooDataHist.cxx:521
 RooDataHist.cxx:522
 RooDataHist.cxx:523
 RooDataHist.cxx:524
 RooDataHist.cxx:525
 RooDataHist.cxx:526
 RooDataHist.cxx:527
 RooDataHist.cxx:528
 RooDataHist.cxx:529
 RooDataHist.cxx:530
 RooDataHist.cxx:531
 RooDataHist.cxx:532
 RooDataHist.cxx:533
 RooDataHist.cxx:534
 RooDataHist.cxx:535
 RooDataHist.cxx:536
 RooDataHist.cxx:537
 RooDataHist.cxx:538
 RooDataHist.cxx:539
 RooDataHist.cxx:540
 RooDataHist.cxx:541
 RooDataHist.cxx:542
 RooDataHist.cxx:543
 RooDataHist.cxx:544
 RooDataHist.cxx:545
 RooDataHist.cxx:546
 RooDataHist.cxx:547
 RooDataHist.cxx:548
 RooDataHist.cxx:549
 RooDataHist.cxx:550
 RooDataHist.cxx:551
 RooDataHist.cxx:552
 RooDataHist.cxx:553
 RooDataHist.cxx:554
 RooDataHist.cxx:555
 RooDataHist.cxx:556
 RooDataHist.cxx:557
 RooDataHist.cxx:558
 RooDataHist.cxx:559
 RooDataHist.cxx:560
 RooDataHist.cxx:561
 RooDataHist.cxx:562
 RooDataHist.cxx:563
 RooDataHist.cxx:564
 RooDataHist.cxx:565
 RooDataHist.cxx:566
 RooDataHist.cxx:567
 RooDataHist.cxx:568
 RooDataHist.cxx:569
 RooDataHist.cxx:570
 RooDataHist.cxx:571
 RooDataHist.cxx:572
 RooDataHist.cxx:573
 RooDataHist.cxx:574
 RooDataHist.cxx:575
 RooDataHist.cxx:576
 RooDataHist.cxx:577
 RooDataHist.cxx:578
 RooDataHist.cxx:579
 RooDataHist.cxx:580
 RooDataHist.cxx:581
 RooDataHist.cxx:582
 RooDataHist.cxx:583
 RooDataHist.cxx:584
 RooDataHist.cxx:585
 RooDataHist.cxx:586
 RooDataHist.cxx:587
 RooDataHist.cxx:588
 RooDataHist.cxx:589
 RooDataHist.cxx:590
 RooDataHist.cxx:591
 RooDataHist.cxx:592
 RooDataHist.cxx:593
 RooDataHist.cxx:594
 RooDataHist.cxx:595
 RooDataHist.cxx:596
 RooDataHist.cxx:597
 RooDataHist.cxx:598
 RooDataHist.cxx:599
 RooDataHist.cxx:600
 RooDataHist.cxx:601
 RooDataHist.cxx:602
 RooDataHist.cxx:603
 RooDataHist.cxx:604
 RooDataHist.cxx:605
 RooDataHist.cxx:606
 RooDataHist.cxx:607
 RooDataHist.cxx:608
 RooDataHist.cxx:609
 RooDataHist.cxx:610
 RooDataHist.cxx:611
 RooDataHist.cxx:612
 RooDataHist.cxx:613
 RooDataHist.cxx:614
 RooDataHist.cxx:615
 RooDataHist.cxx:616
 RooDataHist.cxx:617
 RooDataHist.cxx:618
 RooDataHist.cxx:619
 RooDataHist.cxx:620
 RooDataHist.cxx:621
 RooDataHist.cxx:622
 RooDataHist.cxx:623
 RooDataHist.cxx:624
 RooDataHist.cxx:625
 RooDataHist.cxx:626
 RooDataHist.cxx:627
 RooDataHist.cxx:628
 RooDataHist.cxx:629
 RooDataHist.cxx:630
 RooDataHist.cxx:631
 RooDataHist.cxx:632
 RooDataHist.cxx:633
 RooDataHist.cxx:634
 RooDataHist.cxx:635
 RooDataHist.cxx:636
 RooDataHist.cxx:637
 RooDataHist.cxx:638
 RooDataHist.cxx:639
 RooDataHist.cxx:640
 RooDataHist.cxx:641
 RooDataHist.cxx:642
 RooDataHist.cxx:643
 RooDataHist.cxx:644
 RooDataHist.cxx:645
 RooDataHist.cxx:646
 RooDataHist.cxx:647
 RooDataHist.cxx:648
 RooDataHist.cxx:649
 RooDataHist.cxx:650
 RooDataHist.cxx:651
 RooDataHist.cxx:652
 RooDataHist.cxx:653
 RooDataHist.cxx:654
 RooDataHist.cxx:655
 RooDataHist.cxx:656
 RooDataHist.cxx:657
 RooDataHist.cxx:658
 RooDataHist.cxx:659
 RooDataHist.cxx:660
 RooDataHist.cxx:661
 RooDataHist.cxx:662
 RooDataHist.cxx:663
 RooDataHist.cxx:664
 RooDataHist.cxx:665
 RooDataHist.cxx:666
 RooDataHist.cxx:667
 RooDataHist.cxx:668
 RooDataHist.cxx:669
 RooDataHist.cxx:670
 RooDataHist.cxx:671
 RooDataHist.cxx:672
 RooDataHist.cxx:673
 RooDataHist.cxx:674
 RooDataHist.cxx:675
 RooDataHist.cxx:676
 RooDataHist.cxx:677
 RooDataHist.cxx:678
 RooDataHist.cxx:679
 RooDataHist.cxx:680
 RooDataHist.cxx:681
 RooDataHist.cxx:682
 RooDataHist.cxx:683
 RooDataHist.cxx:684
 RooDataHist.cxx:685
 RooDataHist.cxx:686
 RooDataHist.cxx:687
 RooDataHist.cxx:688
 RooDataHist.cxx:689
 RooDataHist.cxx:690
 RooDataHist.cxx:691
 RooDataHist.cxx:692
 RooDataHist.cxx:693
 RooDataHist.cxx:694
 RooDataHist.cxx:695
 RooDataHist.cxx:696
 RooDataHist.cxx:697
 RooDataHist.cxx:698
 RooDataHist.cxx:699
 RooDataHist.cxx:700
 RooDataHist.cxx:701
 RooDataHist.cxx:702
 RooDataHist.cxx:703
 RooDataHist.cxx:704
 RooDataHist.cxx:705
 RooDataHist.cxx:706
 RooDataHist.cxx:707
 RooDataHist.cxx:708
 RooDataHist.cxx:709
 RooDataHist.cxx:710
 RooDataHist.cxx:711
 RooDataHist.cxx:712
 RooDataHist.cxx:713
 RooDataHist.cxx:714
 RooDataHist.cxx:715
 RooDataHist.cxx:716
 RooDataHist.cxx:717
 RooDataHist.cxx:718
 RooDataHist.cxx:719
 RooDataHist.cxx:720
 RooDataHist.cxx:721
 RooDataHist.cxx:722
 RooDataHist.cxx:723
 RooDataHist.cxx:724
 RooDataHist.cxx:725
 RooDataHist.cxx:726
 RooDataHist.cxx:727
 RooDataHist.cxx:728
 RooDataHist.cxx:729
 RooDataHist.cxx:730
 RooDataHist.cxx:731
 RooDataHist.cxx:732
 RooDataHist.cxx:733
 RooDataHist.cxx:734
 RooDataHist.cxx:735
 RooDataHist.cxx:736
 RooDataHist.cxx:737
 RooDataHist.cxx:738
 RooDataHist.cxx:739
 RooDataHist.cxx:740
 RooDataHist.cxx:741
 RooDataHist.cxx:742
 RooDataHist.cxx:743
 RooDataHist.cxx:744
 RooDataHist.cxx:745
 RooDataHist.cxx:746
 RooDataHist.cxx:747
 RooDataHist.cxx:748
 RooDataHist.cxx:749
 RooDataHist.cxx:750
 RooDataHist.cxx:751
 RooDataHist.cxx:752
 RooDataHist.cxx:753
 RooDataHist.cxx:754
 RooDataHist.cxx:755
 RooDataHist.cxx:756
 RooDataHist.cxx:757
 RooDataHist.cxx:758
 RooDataHist.cxx:759
 RooDataHist.cxx:760
 RooDataHist.cxx:761
 RooDataHist.cxx:762
 RooDataHist.cxx:763
 RooDataHist.cxx:764
 RooDataHist.cxx:765
 RooDataHist.cxx:766
 RooDataHist.cxx:767
 RooDataHist.cxx:768
 RooDataHist.cxx:769
 RooDataHist.cxx:770
 RooDataHist.cxx:771
 RooDataHist.cxx:772
 RooDataHist.cxx:773
 RooDataHist.cxx:774
 RooDataHist.cxx:775
 RooDataHist.cxx:776
 RooDataHist.cxx:777
 RooDataHist.cxx:778
 RooDataHist.cxx:779
 RooDataHist.cxx:780
 RooDataHist.cxx:781
 RooDataHist.cxx:782
 RooDataHist.cxx:783
 RooDataHist.cxx:784
 RooDataHist.cxx:785
 RooDataHist.cxx:786
 RooDataHist.cxx:787
 RooDataHist.cxx:788
 RooDataHist.cxx:789
 RooDataHist.cxx:790
 RooDataHist.cxx:791
 RooDataHist.cxx:792
 RooDataHist.cxx:793
 RooDataHist.cxx:794
 RooDataHist.cxx:795
 RooDataHist.cxx:796
 RooDataHist.cxx:797
 RooDataHist.cxx:798
 RooDataHist.cxx:799
 RooDataHist.cxx:800
 RooDataHist.cxx:801
 RooDataHist.cxx:802
 RooDataHist.cxx:803
 RooDataHist.cxx:804
 RooDataHist.cxx:805
 RooDataHist.cxx:806
 RooDataHist.cxx:807
 RooDataHist.cxx:808
 RooDataHist.cxx:809
 RooDataHist.cxx:810
 RooDataHist.cxx:811
 RooDataHist.cxx:812
 RooDataHist.cxx:813
 RooDataHist.cxx:814
 RooDataHist.cxx:815
 RooDataHist.cxx:816
 RooDataHist.cxx:817
 RooDataHist.cxx:818
 RooDataHist.cxx:819
 RooDataHist.cxx:820
 RooDataHist.cxx:821
 RooDataHist.cxx:822
 RooDataHist.cxx:823
 RooDataHist.cxx:824
 RooDataHist.cxx:825
 RooDataHist.cxx:826
 RooDataHist.cxx:827
 RooDataHist.cxx:828
 RooDataHist.cxx:829
 RooDataHist.cxx:830
 RooDataHist.cxx:831
 RooDataHist.cxx:832
 RooDataHist.cxx:833
 RooDataHist.cxx:834
 RooDataHist.cxx:835
 RooDataHist.cxx:836
 RooDataHist.cxx:837
 RooDataHist.cxx:838
 RooDataHist.cxx:839
 RooDataHist.cxx:840
 RooDataHist.cxx:841
 RooDataHist.cxx:842
 RooDataHist.cxx:843
 RooDataHist.cxx:844
 RooDataHist.cxx:845
 RooDataHist.cxx:846
 RooDataHist.cxx:847
 RooDataHist.cxx:848
 RooDataHist.cxx:849
 RooDataHist.cxx:850
 RooDataHist.cxx:851
 RooDataHist.cxx:852
 RooDataHist.cxx:853
 RooDataHist.cxx:854
 RooDataHist.cxx:855
 RooDataHist.cxx:856
 RooDataHist.cxx:857
 RooDataHist.cxx:858
 RooDataHist.cxx:859
 RooDataHist.cxx:860
 RooDataHist.cxx:861
 RooDataHist.cxx:862
 RooDataHist.cxx:863
 RooDataHist.cxx:864
 RooDataHist.cxx:865
 RooDataHist.cxx:866
 RooDataHist.cxx:867
 RooDataHist.cxx:868
 RooDataHist.cxx:869
 RooDataHist.cxx:870
 RooDataHist.cxx:871
 RooDataHist.cxx:872
 RooDataHist.cxx:873
 RooDataHist.cxx:874
 RooDataHist.cxx:875
 RooDataHist.cxx:876
 RooDataHist.cxx:877
 RooDataHist.cxx:878
 RooDataHist.cxx:879
 RooDataHist.cxx:880
 RooDataHist.cxx:881
 RooDataHist.cxx:882
 RooDataHist.cxx:883
 RooDataHist.cxx:884
 RooDataHist.cxx:885
 RooDataHist.cxx:886
 RooDataHist.cxx:887
 RooDataHist.cxx:888
 RooDataHist.cxx:889
 RooDataHist.cxx:890
 RooDataHist.cxx:891
 RooDataHist.cxx:892
 RooDataHist.cxx:893
 RooDataHist.cxx:894
 RooDataHist.cxx:895
 RooDataHist.cxx:896
 RooDataHist.cxx:897
 RooDataHist.cxx:898
 RooDataHist.cxx:899
 RooDataHist.cxx:900
 RooDataHist.cxx:901
 RooDataHist.cxx:902
 RooDataHist.cxx:903
 RooDataHist.cxx:904
 RooDataHist.cxx:905
 RooDataHist.cxx:906
 RooDataHist.cxx:907
 RooDataHist.cxx:908
 RooDataHist.cxx:909
 RooDataHist.cxx:910
 RooDataHist.cxx:911
 RooDataHist.cxx:912
 RooDataHist.cxx:913
 RooDataHist.cxx:914
 RooDataHist.cxx:915
 RooDataHist.cxx:916
 RooDataHist.cxx:917
 RooDataHist.cxx:918
 RooDataHist.cxx:919
 RooDataHist.cxx:920
 RooDataHist.cxx:921
 RooDataHist.cxx:922
 RooDataHist.cxx:923
 RooDataHist.cxx:924
 RooDataHist.cxx:925
 RooDataHist.cxx:926
 RooDataHist.cxx:927
 RooDataHist.cxx:928
 RooDataHist.cxx:929
 RooDataHist.cxx:930
 RooDataHist.cxx:931
 RooDataHist.cxx:932
 RooDataHist.cxx:933
 RooDataHist.cxx:934
 RooDataHist.cxx:935
 RooDataHist.cxx:936
 RooDataHist.cxx:937
 RooDataHist.cxx:938
 RooDataHist.cxx:939
 RooDataHist.cxx:940
 RooDataHist.cxx:941
 RooDataHist.cxx:942
 RooDataHist.cxx:943
 RooDataHist.cxx:944
 RooDataHist.cxx:945
 RooDataHist.cxx:946
 RooDataHist.cxx:947
 RooDataHist.cxx:948
 RooDataHist.cxx:949
 RooDataHist.cxx:950
 RooDataHist.cxx:951
 RooDataHist.cxx:952
 RooDataHist.cxx:953
 RooDataHist.cxx:954
 RooDataHist.cxx:955
 RooDataHist.cxx:956
 RooDataHist.cxx:957
 RooDataHist.cxx:958
 RooDataHist.cxx:959
 RooDataHist.cxx:960
 RooDataHist.cxx:961
 RooDataHist.cxx:962
 RooDataHist.cxx:963
 RooDataHist.cxx:964
 RooDataHist.cxx:965
 RooDataHist.cxx:966
 RooDataHist.cxx:967
 RooDataHist.cxx:968
 RooDataHist.cxx:969
 RooDataHist.cxx:970
 RooDataHist.cxx:971
 RooDataHist.cxx:972
 RooDataHist.cxx:973
 RooDataHist.cxx:974
 RooDataHist.cxx:975
 RooDataHist.cxx:976
 RooDataHist.cxx:977
 RooDataHist.cxx:978
 RooDataHist.cxx:979
 RooDataHist.cxx:980
 RooDataHist.cxx:981
 RooDataHist.cxx:982
 RooDataHist.cxx:983
 RooDataHist.cxx:984
 RooDataHist.cxx:985
 RooDataHist.cxx:986
 RooDataHist.cxx:987
 RooDataHist.cxx:988
 RooDataHist.cxx:989
 RooDataHist.cxx:990
 RooDataHist.cxx:991
 RooDataHist.cxx:992
 RooDataHist.cxx:993
 RooDataHist.cxx:994
 RooDataHist.cxx:995
 RooDataHist.cxx:996
 RooDataHist.cxx:997
 RooDataHist.cxx:998
 RooDataHist.cxx:999
 RooDataHist.cxx:1000
 RooDataHist.cxx:1001
 RooDataHist.cxx:1002
 RooDataHist.cxx:1003
 RooDataHist.cxx:1004
 RooDataHist.cxx:1005
 RooDataHist.cxx:1006
 RooDataHist.cxx:1007
 RooDataHist.cxx:1008
 RooDataHist.cxx:1009
 RooDataHist.cxx:1010
 RooDataHist.cxx:1011
 RooDataHist.cxx:1012
 RooDataHist.cxx:1013
 RooDataHist.cxx:1014
 RooDataHist.cxx:1015
 RooDataHist.cxx:1016
 RooDataHist.cxx:1017
 RooDataHist.cxx:1018
 RooDataHist.cxx:1019
 RooDataHist.cxx:1020
 RooDataHist.cxx:1021
 RooDataHist.cxx:1022
 RooDataHist.cxx:1023
 RooDataHist.cxx:1024
 RooDataHist.cxx:1025
 RooDataHist.cxx:1026
 RooDataHist.cxx:1027
 RooDataHist.cxx:1028
 RooDataHist.cxx:1029
 RooDataHist.cxx:1030
 RooDataHist.cxx:1031
 RooDataHist.cxx:1032
 RooDataHist.cxx:1033
 RooDataHist.cxx:1034
 RooDataHist.cxx:1035
 RooDataHist.cxx:1036
 RooDataHist.cxx:1037
 RooDataHist.cxx:1038
 RooDataHist.cxx:1039
 RooDataHist.cxx:1040
 RooDataHist.cxx:1041
 RooDataHist.cxx:1042
 RooDataHist.cxx:1043
 RooDataHist.cxx:1044
 RooDataHist.cxx:1045
 RooDataHist.cxx:1046
 RooDataHist.cxx:1047
 RooDataHist.cxx:1048
 RooDataHist.cxx:1049
 RooDataHist.cxx:1050
 RooDataHist.cxx:1051
 RooDataHist.cxx:1052
 RooDataHist.cxx:1053
 RooDataHist.cxx:1054
 RooDataHist.cxx:1055
 RooDataHist.cxx:1056
 RooDataHist.cxx:1057
 RooDataHist.cxx:1058
 RooDataHist.cxx:1059
 RooDataHist.cxx:1060
 RooDataHist.cxx:1061
 RooDataHist.cxx:1062
 RooDataHist.cxx:1063
 RooDataHist.cxx:1064
 RooDataHist.cxx:1065
 RooDataHist.cxx:1066
 RooDataHist.cxx:1067
 RooDataHist.cxx:1068
 RooDataHist.cxx:1069
 RooDataHist.cxx:1070
 RooDataHist.cxx:1071
 RooDataHist.cxx:1072
 RooDataHist.cxx:1073
 RooDataHist.cxx:1074
 RooDataHist.cxx:1075
 RooDataHist.cxx:1076
 RooDataHist.cxx:1077
 RooDataHist.cxx:1078
 RooDataHist.cxx:1079
 RooDataHist.cxx:1080
 RooDataHist.cxx:1081
 RooDataHist.cxx:1082
 RooDataHist.cxx:1083
 RooDataHist.cxx:1084
 RooDataHist.cxx:1085
 RooDataHist.cxx:1086
 RooDataHist.cxx:1087
 RooDataHist.cxx:1088
 RooDataHist.cxx:1089
 RooDataHist.cxx:1090
 RooDataHist.cxx:1091
 RooDataHist.cxx:1092
 RooDataHist.cxx:1093
 RooDataHist.cxx:1094
 RooDataHist.cxx:1095
 RooDataHist.cxx:1096
 RooDataHist.cxx:1097
 RooDataHist.cxx:1098
 RooDataHist.cxx:1099
 RooDataHist.cxx:1100
 RooDataHist.cxx:1101
 RooDataHist.cxx:1102
 RooDataHist.cxx:1103
 RooDataHist.cxx:1104
 RooDataHist.cxx:1105
 RooDataHist.cxx:1106
 RooDataHist.cxx:1107
 RooDataHist.cxx:1108
 RooDataHist.cxx:1109
 RooDataHist.cxx:1110
 RooDataHist.cxx:1111
 RooDataHist.cxx:1112
 RooDataHist.cxx:1113
 RooDataHist.cxx:1114
 RooDataHist.cxx:1115
 RooDataHist.cxx:1116
 RooDataHist.cxx:1117
 RooDataHist.cxx:1118
 RooDataHist.cxx:1119
 RooDataHist.cxx:1120
 RooDataHist.cxx:1121
 RooDataHist.cxx:1122
 RooDataHist.cxx:1123
 RooDataHist.cxx:1124
 RooDataHist.cxx:1125
 RooDataHist.cxx:1126
 RooDataHist.cxx:1127
 RooDataHist.cxx:1128
 RooDataHist.cxx:1129
 RooDataHist.cxx:1130
 RooDataHist.cxx:1131
 RooDataHist.cxx:1132
 RooDataHist.cxx:1133
 RooDataHist.cxx:1134
 RooDataHist.cxx:1135
 RooDataHist.cxx:1136
 RooDataHist.cxx:1137
 RooDataHist.cxx:1138
 RooDataHist.cxx:1139
 RooDataHist.cxx:1140
 RooDataHist.cxx:1141
 RooDataHist.cxx:1142
 RooDataHist.cxx:1143
 RooDataHist.cxx:1144
 RooDataHist.cxx:1145
 RooDataHist.cxx:1146
 RooDataHist.cxx:1147
 RooDataHist.cxx:1148
 RooDataHist.cxx:1149
 RooDataHist.cxx:1150
 RooDataHist.cxx:1151
 RooDataHist.cxx:1152
 RooDataHist.cxx:1153
 RooDataHist.cxx:1154
 RooDataHist.cxx:1155
 RooDataHist.cxx:1156
 RooDataHist.cxx:1157
 RooDataHist.cxx:1158
 RooDataHist.cxx:1159
 RooDataHist.cxx:1160
 RooDataHist.cxx:1161
 RooDataHist.cxx:1162
 RooDataHist.cxx:1163
 RooDataHist.cxx:1164
 RooDataHist.cxx:1165
 RooDataHist.cxx:1166
 RooDataHist.cxx:1167
 RooDataHist.cxx:1168
 RooDataHist.cxx:1169
 RooDataHist.cxx:1170
 RooDataHist.cxx:1171
 RooDataHist.cxx:1172
 RooDataHist.cxx:1173
 RooDataHist.cxx:1174
 RooDataHist.cxx:1175
 RooDataHist.cxx:1176
 RooDataHist.cxx:1177
 RooDataHist.cxx:1178
 RooDataHist.cxx:1179
 RooDataHist.cxx:1180
 RooDataHist.cxx:1181
 RooDataHist.cxx:1182
 RooDataHist.cxx:1183
 RooDataHist.cxx:1184
 RooDataHist.cxx:1185
 RooDataHist.cxx:1186
 RooDataHist.cxx:1187
 RooDataHist.cxx:1188
 RooDataHist.cxx:1189
 RooDataHist.cxx:1190
 RooDataHist.cxx:1191
 RooDataHist.cxx:1192
 RooDataHist.cxx:1193
 RooDataHist.cxx:1194
 RooDataHist.cxx:1195
 RooDataHist.cxx:1196
 RooDataHist.cxx:1197
 RooDataHist.cxx:1198
 RooDataHist.cxx:1199
 RooDataHist.cxx:1200
 RooDataHist.cxx:1201
 RooDataHist.cxx:1202
 RooDataHist.cxx:1203
 RooDataHist.cxx:1204
 RooDataHist.cxx:1205
 RooDataHist.cxx:1206
 RooDataHist.cxx:1207
 RooDataHist.cxx:1208
 RooDataHist.cxx:1209
 RooDataHist.cxx:1210
 RooDataHist.cxx:1211
 RooDataHist.cxx:1212
 RooDataHist.cxx:1213
 RooDataHist.cxx:1214
 RooDataHist.cxx:1215
 RooDataHist.cxx:1216
 RooDataHist.cxx:1217
 RooDataHist.cxx:1218
 RooDataHist.cxx:1219
 RooDataHist.cxx:1220
 RooDataHist.cxx:1221
 RooDataHist.cxx:1222
 RooDataHist.cxx:1223
 RooDataHist.cxx:1224
 RooDataHist.cxx:1225
 RooDataHist.cxx:1226
 RooDataHist.cxx:1227
 RooDataHist.cxx:1228
 RooDataHist.cxx:1229
 RooDataHist.cxx:1230
 RooDataHist.cxx:1231
 RooDataHist.cxx:1232
 RooDataHist.cxx:1233
 RooDataHist.cxx:1234
 RooDataHist.cxx:1235
 RooDataHist.cxx:1236
 RooDataHist.cxx:1237
 RooDataHist.cxx:1238
 RooDataHist.cxx:1239
 RooDataHist.cxx:1240
 RooDataHist.cxx:1241
 RooDataHist.cxx:1242
 RooDataHist.cxx:1243
 RooDataHist.cxx:1244
 RooDataHist.cxx:1245
 RooDataHist.cxx:1246
 RooDataHist.cxx:1247
 RooDataHist.cxx:1248
 RooDataHist.cxx:1249
 RooDataHist.cxx:1250
 RooDataHist.cxx:1251
 RooDataHist.cxx:1252
 RooDataHist.cxx:1253
 RooDataHist.cxx:1254
 RooDataHist.cxx:1255
 RooDataHist.cxx:1256
 RooDataHist.cxx:1257
 RooDataHist.cxx:1258
 RooDataHist.cxx:1259
 RooDataHist.cxx:1260
 RooDataHist.cxx:1261
 RooDataHist.cxx:1262
 RooDataHist.cxx:1263
 RooDataHist.cxx:1264
 RooDataHist.cxx:1265
 RooDataHist.cxx:1266
 RooDataHist.cxx:1267
 RooDataHist.cxx:1268
 RooDataHist.cxx:1269
 RooDataHist.cxx:1270
 RooDataHist.cxx:1271
 RooDataHist.cxx:1272
 RooDataHist.cxx:1273
 RooDataHist.cxx:1274
 RooDataHist.cxx:1275
 RooDataHist.cxx:1276
 RooDataHist.cxx:1277
 RooDataHist.cxx:1278
 RooDataHist.cxx:1279
 RooDataHist.cxx:1280
 RooDataHist.cxx:1281
 RooDataHist.cxx:1282
 RooDataHist.cxx:1283
 RooDataHist.cxx:1284
 RooDataHist.cxx:1285
 RooDataHist.cxx:1286
 RooDataHist.cxx:1287
 RooDataHist.cxx:1288
 RooDataHist.cxx:1289
 RooDataHist.cxx:1290
 RooDataHist.cxx:1291
 RooDataHist.cxx:1292
 RooDataHist.cxx:1293
 RooDataHist.cxx:1294
 RooDataHist.cxx:1295
 RooDataHist.cxx:1296
 RooDataHist.cxx:1297
 RooDataHist.cxx:1298
 RooDataHist.cxx:1299
 RooDataHist.cxx:1300
 RooDataHist.cxx:1301
 RooDataHist.cxx:1302
 RooDataHist.cxx:1303
 RooDataHist.cxx:1304
 RooDataHist.cxx:1305
 RooDataHist.cxx:1306
 RooDataHist.cxx:1307
 RooDataHist.cxx:1308
 RooDataHist.cxx:1309
 RooDataHist.cxx:1310
 RooDataHist.cxx:1311
 RooDataHist.cxx:1312
 RooDataHist.cxx:1313
 RooDataHist.cxx:1314
 RooDataHist.cxx:1315
 RooDataHist.cxx:1316
 RooDataHist.cxx:1317
 RooDataHist.cxx:1318
 RooDataHist.cxx:1319
 RooDataHist.cxx:1320
 RooDataHist.cxx:1321
 RooDataHist.cxx:1322
 RooDataHist.cxx:1323
 RooDataHist.cxx:1324
 RooDataHist.cxx:1325
 RooDataHist.cxx:1326
 RooDataHist.cxx:1327
 RooDataHist.cxx:1328
 RooDataHist.cxx:1329
 RooDataHist.cxx:1330
 RooDataHist.cxx:1331
 RooDataHist.cxx:1332
 RooDataHist.cxx:1333
 RooDataHist.cxx:1334
 RooDataHist.cxx:1335
 RooDataHist.cxx:1336
 RooDataHist.cxx:1337
 RooDataHist.cxx:1338
 RooDataHist.cxx:1339
 RooDataHist.cxx:1340
 RooDataHist.cxx:1341
 RooDataHist.cxx:1342
 RooDataHist.cxx:1343
 RooDataHist.cxx:1344
 RooDataHist.cxx:1345
 RooDataHist.cxx:1346
 RooDataHist.cxx:1347
 RooDataHist.cxx:1348
 RooDataHist.cxx:1349
 RooDataHist.cxx:1350
 RooDataHist.cxx:1351
 RooDataHist.cxx:1352
 RooDataHist.cxx:1353
 RooDataHist.cxx:1354
 RooDataHist.cxx:1355
 RooDataHist.cxx:1356
 RooDataHist.cxx:1357
 RooDataHist.cxx:1358
 RooDataHist.cxx:1359
 RooDataHist.cxx:1360
 RooDataHist.cxx:1361
 RooDataHist.cxx:1362
 RooDataHist.cxx:1363
 RooDataHist.cxx:1364
 RooDataHist.cxx:1365
 RooDataHist.cxx:1366
 RooDataHist.cxx:1367
 RooDataHist.cxx:1368
 RooDataHist.cxx:1369
 RooDataHist.cxx:1370
 RooDataHist.cxx:1371
 RooDataHist.cxx:1372
 RooDataHist.cxx:1373
 RooDataHist.cxx:1374
 RooDataHist.cxx:1375
 RooDataHist.cxx:1376
 RooDataHist.cxx:1377
 RooDataHist.cxx:1378
 RooDataHist.cxx:1379
 RooDataHist.cxx:1380
 RooDataHist.cxx:1381
 RooDataHist.cxx:1382
 RooDataHist.cxx:1383
 RooDataHist.cxx:1384
 RooDataHist.cxx:1385
 RooDataHist.cxx:1386
 RooDataHist.cxx:1387
 RooDataHist.cxx:1388
 RooDataHist.cxx:1389
 RooDataHist.cxx:1390
 RooDataHist.cxx:1391
 RooDataHist.cxx:1392
 RooDataHist.cxx:1393
 RooDataHist.cxx:1394
 RooDataHist.cxx:1395
 RooDataHist.cxx:1396
 RooDataHist.cxx:1397
 RooDataHist.cxx:1398
 RooDataHist.cxx:1399
 RooDataHist.cxx:1400
 RooDataHist.cxx:1401
 RooDataHist.cxx:1402
 RooDataHist.cxx:1403
 RooDataHist.cxx:1404
 RooDataHist.cxx:1405
 RooDataHist.cxx:1406
 RooDataHist.cxx:1407
 RooDataHist.cxx:1408
 RooDataHist.cxx:1409
 RooDataHist.cxx:1410
 RooDataHist.cxx:1411
 RooDataHist.cxx:1412
 RooDataHist.cxx:1413
 RooDataHist.cxx:1414
 RooDataHist.cxx:1415
 RooDataHist.cxx:1416
 RooDataHist.cxx:1417
 RooDataHist.cxx:1418
 RooDataHist.cxx:1419
 RooDataHist.cxx:1420
 RooDataHist.cxx:1421
 RooDataHist.cxx:1422
 RooDataHist.cxx:1423
 RooDataHist.cxx:1424
 RooDataHist.cxx:1425
 RooDataHist.cxx:1426
 RooDataHist.cxx:1427
 RooDataHist.cxx:1428
 RooDataHist.cxx:1429
 RooDataHist.cxx:1430
 RooDataHist.cxx:1431
 RooDataHist.cxx:1432
 RooDataHist.cxx:1433
 RooDataHist.cxx:1434
 RooDataHist.cxx:1435
 RooDataHist.cxx:1436
 RooDataHist.cxx:1437
 RooDataHist.cxx:1438
 RooDataHist.cxx:1439
 RooDataHist.cxx:1440
 RooDataHist.cxx:1441
 RooDataHist.cxx:1442
 RooDataHist.cxx:1443
 RooDataHist.cxx:1444
 RooDataHist.cxx:1445
 RooDataHist.cxx:1446
 RooDataHist.cxx:1447
 RooDataHist.cxx:1448
 RooDataHist.cxx:1449
 RooDataHist.cxx:1450
 RooDataHist.cxx:1451
 RooDataHist.cxx:1452
 RooDataHist.cxx:1453
 RooDataHist.cxx:1454
 RooDataHist.cxx:1455
 RooDataHist.cxx:1456
 RooDataHist.cxx:1457
 RooDataHist.cxx:1458
 RooDataHist.cxx:1459
 RooDataHist.cxx:1460
 RooDataHist.cxx:1461
 RooDataHist.cxx:1462
 RooDataHist.cxx:1463
 RooDataHist.cxx:1464
 RooDataHist.cxx:1465
 RooDataHist.cxx:1466
 RooDataHist.cxx:1467
 RooDataHist.cxx:1468
 RooDataHist.cxx:1469
 RooDataHist.cxx:1470
 RooDataHist.cxx:1471
 RooDataHist.cxx:1472
 RooDataHist.cxx:1473
 RooDataHist.cxx:1474
 RooDataHist.cxx:1475
 RooDataHist.cxx:1476
 RooDataHist.cxx:1477
 RooDataHist.cxx:1478
 RooDataHist.cxx:1479
 RooDataHist.cxx:1480
 RooDataHist.cxx:1481
 RooDataHist.cxx:1482
 RooDataHist.cxx:1483
 RooDataHist.cxx:1484
 RooDataHist.cxx:1485
 RooDataHist.cxx:1486
 RooDataHist.cxx:1487
 RooDataHist.cxx:1488
 RooDataHist.cxx:1489
 RooDataHist.cxx:1490
 RooDataHist.cxx:1491
 RooDataHist.cxx:1492
 RooDataHist.cxx:1493
 RooDataHist.cxx:1494
 RooDataHist.cxx:1495
 RooDataHist.cxx:1496
 RooDataHist.cxx:1497
 RooDataHist.cxx:1498
 RooDataHist.cxx:1499
 RooDataHist.cxx:1500
 RooDataHist.cxx:1501
 RooDataHist.cxx:1502
 RooDataHist.cxx:1503
 RooDataHist.cxx:1504
 RooDataHist.cxx:1505
 RooDataHist.cxx:1506
 RooDataHist.cxx:1507
 RooDataHist.cxx:1508
 RooDataHist.cxx:1509
 RooDataHist.cxx:1510
 RooDataHist.cxx:1511
 RooDataHist.cxx:1512
 RooDataHist.cxx:1513
 RooDataHist.cxx:1514
 RooDataHist.cxx:1515
 RooDataHist.cxx:1516
 RooDataHist.cxx:1517
 RooDataHist.cxx:1518
 RooDataHist.cxx:1519
 RooDataHist.cxx:1520
 RooDataHist.cxx:1521
 RooDataHist.cxx:1522
 RooDataHist.cxx:1523
 RooDataHist.cxx:1524
 RooDataHist.cxx:1525
 RooDataHist.cxx:1526
 RooDataHist.cxx:1527
 RooDataHist.cxx:1528
 RooDataHist.cxx:1529
 RooDataHist.cxx:1530
 RooDataHist.cxx:1531
 RooDataHist.cxx:1532
 RooDataHist.cxx:1533
 RooDataHist.cxx:1534
 RooDataHist.cxx:1535
 RooDataHist.cxx:1536
 RooDataHist.cxx:1537
 RooDataHist.cxx:1538
 RooDataHist.cxx:1539
 RooDataHist.cxx:1540
 RooDataHist.cxx:1541
 RooDataHist.cxx:1542
 RooDataHist.cxx:1543
 RooDataHist.cxx:1544
 RooDataHist.cxx:1545
 RooDataHist.cxx:1546
 RooDataHist.cxx:1547
 RooDataHist.cxx:1548
 RooDataHist.cxx:1549
 RooDataHist.cxx:1550
 RooDataHist.cxx:1551
 RooDataHist.cxx:1552
 RooDataHist.cxx:1553
 RooDataHist.cxx:1554
 RooDataHist.cxx:1555
 RooDataHist.cxx:1556
 RooDataHist.cxx:1557
 RooDataHist.cxx:1558
 RooDataHist.cxx:1559
 RooDataHist.cxx:1560
 RooDataHist.cxx:1561
 RooDataHist.cxx:1562
 RooDataHist.cxx:1563
 RooDataHist.cxx:1564
 RooDataHist.cxx:1565
 RooDataHist.cxx:1566
 RooDataHist.cxx:1567
 RooDataHist.cxx:1568
 RooDataHist.cxx:1569
 RooDataHist.cxx:1570
 RooDataHist.cxx:1571
 RooDataHist.cxx:1572
 RooDataHist.cxx:1573
 RooDataHist.cxx:1574
 RooDataHist.cxx:1575
 RooDataHist.cxx:1576
 RooDataHist.cxx:1577
 RooDataHist.cxx:1578
 RooDataHist.cxx:1579
 RooDataHist.cxx:1580
 RooDataHist.cxx:1581
 RooDataHist.cxx:1582
 RooDataHist.cxx:1583
 RooDataHist.cxx:1584
 RooDataHist.cxx:1585
 RooDataHist.cxx:1586
 RooDataHist.cxx:1587
 RooDataHist.cxx:1588
 RooDataHist.cxx:1589
 RooDataHist.cxx:1590
 RooDataHist.cxx:1591
 RooDataHist.cxx:1592
 RooDataHist.cxx:1593
 RooDataHist.cxx:1594
 RooDataHist.cxx:1595
 RooDataHist.cxx:1596
 RooDataHist.cxx:1597
 RooDataHist.cxx:1598
 RooDataHist.cxx:1599
 RooDataHist.cxx:1600
 RooDataHist.cxx:1601
 RooDataHist.cxx:1602
 RooDataHist.cxx:1603
 RooDataHist.cxx:1604
 RooDataHist.cxx:1605
 RooDataHist.cxx:1606
 RooDataHist.cxx:1607
 RooDataHist.cxx:1608
 RooDataHist.cxx:1609
 RooDataHist.cxx:1610
 RooDataHist.cxx:1611
 RooDataHist.cxx:1612
 RooDataHist.cxx:1613
 RooDataHist.cxx:1614
 RooDataHist.cxx:1615
 RooDataHist.cxx:1616
 RooDataHist.cxx:1617
 RooDataHist.cxx:1618
 RooDataHist.cxx:1619
 RooDataHist.cxx:1620
 RooDataHist.cxx:1621
 RooDataHist.cxx:1622
 RooDataHist.cxx:1623
 RooDataHist.cxx:1624
 RooDataHist.cxx:1625
 RooDataHist.cxx:1626
 RooDataHist.cxx:1627
 RooDataHist.cxx:1628
 RooDataHist.cxx:1629
 RooDataHist.cxx:1630
 RooDataHist.cxx:1631
 RooDataHist.cxx:1632
 RooDataHist.cxx:1633
 RooDataHist.cxx:1634
 RooDataHist.cxx:1635
 RooDataHist.cxx:1636
 RooDataHist.cxx:1637
 RooDataHist.cxx:1638
 RooDataHist.cxx:1639
 RooDataHist.cxx:1640
 RooDataHist.cxx:1641
 RooDataHist.cxx:1642
 RooDataHist.cxx:1643
 RooDataHist.cxx:1644
 RooDataHist.cxx:1645
 RooDataHist.cxx:1646
 RooDataHist.cxx:1647
 RooDataHist.cxx:1648
 RooDataHist.cxx:1649
 RooDataHist.cxx:1650
 RooDataHist.cxx:1651
 RooDataHist.cxx:1652
 RooDataHist.cxx:1653
 RooDataHist.cxx:1654
 RooDataHist.cxx:1655
 RooDataHist.cxx:1656
 RooDataHist.cxx:1657
 RooDataHist.cxx:1658
 RooDataHist.cxx:1659
 RooDataHist.cxx:1660
 RooDataHist.cxx:1661
 RooDataHist.cxx:1662
 RooDataHist.cxx:1663
 RooDataHist.cxx:1664
 RooDataHist.cxx:1665
 RooDataHist.cxx:1666
 RooDataHist.cxx:1667
 RooDataHist.cxx:1668
 RooDataHist.cxx:1669
 RooDataHist.cxx:1670
 RooDataHist.cxx:1671
 RooDataHist.cxx:1672
 RooDataHist.cxx:1673
 RooDataHist.cxx:1674
 RooDataHist.cxx:1675
 RooDataHist.cxx:1676
 RooDataHist.cxx:1677
 RooDataHist.cxx:1678
 RooDataHist.cxx:1679
 RooDataHist.cxx:1680
 RooDataHist.cxx:1681
 RooDataHist.cxx:1682
 RooDataHist.cxx:1683
 RooDataHist.cxx:1684
 RooDataHist.cxx:1685
 RooDataHist.cxx:1686
 RooDataHist.cxx:1687
 RooDataHist.cxx:1688
 RooDataHist.cxx:1689
 RooDataHist.cxx:1690
 RooDataHist.cxx:1691
 RooDataHist.cxx:1692
 RooDataHist.cxx:1693
 RooDataHist.cxx:1694
 RooDataHist.cxx:1695
 RooDataHist.cxx:1696
 RooDataHist.cxx:1697
 RooDataHist.cxx:1698
 RooDataHist.cxx:1699
 RooDataHist.cxx:1700
 RooDataHist.cxx:1701
 RooDataHist.cxx:1702
 RooDataHist.cxx:1703
 RooDataHist.cxx:1704
 RooDataHist.cxx:1705
 RooDataHist.cxx:1706
 RooDataHist.cxx:1707
 RooDataHist.cxx:1708
 RooDataHist.cxx:1709
 RooDataHist.cxx:1710
 RooDataHist.cxx:1711
 RooDataHist.cxx:1712
 RooDataHist.cxx:1713
 RooDataHist.cxx:1714
 RooDataHist.cxx:1715
 RooDataHist.cxx:1716
 RooDataHist.cxx:1717
 RooDataHist.cxx:1718
 RooDataHist.cxx:1719
 RooDataHist.cxx:1720
 RooDataHist.cxx:1721
 RooDataHist.cxx:1722
 RooDataHist.cxx:1723
 RooDataHist.cxx:1724
 RooDataHist.cxx:1725
 RooDataHist.cxx:1726
 RooDataHist.cxx:1727
 RooDataHist.cxx:1728
 RooDataHist.cxx:1729
 RooDataHist.cxx:1730
 RooDataHist.cxx:1731
 RooDataHist.cxx:1732
 RooDataHist.cxx:1733
 RooDataHist.cxx:1734
 RooDataHist.cxx:1735
 RooDataHist.cxx:1736
 RooDataHist.cxx:1737
 RooDataHist.cxx:1738
 RooDataHist.cxx:1739
 RooDataHist.cxx:1740
 RooDataHist.cxx:1741
 RooDataHist.cxx:1742
 RooDataHist.cxx:1743
 RooDataHist.cxx:1744
 RooDataHist.cxx:1745
 RooDataHist.cxx:1746
 RooDataHist.cxx:1747
 RooDataHist.cxx:1748
 RooDataHist.cxx:1749
 RooDataHist.cxx:1750
 RooDataHist.cxx:1751
 RooDataHist.cxx:1752
 RooDataHist.cxx:1753
 RooDataHist.cxx:1754
 RooDataHist.cxx:1755
 RooDataHist.cxx:1756
 RooDataHist.cxx:1757
 RooDataHist.cxx:1758
 RooDataHist.cxx:1759
 RooDataHist.cxx:1760
 RooDataHist.cxx:1761
 RooDataHist.cxx:1762
 RooDataHist.cxx:1763
 RooDataHist.cxx:1764
 RooDataHist.cxx:1765
 RooDataHist.cxx:1766
 RooDataHist.cxx:1767
 RooDataHist.cxx:1768
 RooDataHist.cxx:1769
 RooDataHist.cxx:1770
 RooDataHist.cxx:1771
 RooDataHist.cxx:1772
 RooDataHist.cxx:1773
 RooDataHist.cxx:1774
 RooDataHist.cxx:1775
 RooDataHist.cxx:1776
 RooDataHist.cxx:1777
 RooDataHist.cxx:1778
 RooDataHist.cxx:1779
 RooDataHist.cxx:1780
 RooDataHist.cxx:1781
 RooDataHist.cxx:1782
 RooDataHist.cxx:1783
 RooDataHist.cxx:1784
 RooDataHist.cxx:1785
 RooDataHist.cxx:1786
 RooDataHist.cxx:1787
 RooDataHist.cxx:1788
 RooDataHist.cxx:1789
 RooDataHist.cxx:1790
 RooDataHist.cxx:1791
 RooDataHist.cxx:1792
 RooDataHist.cxx:1793
 RooDataHist.cxx:1794
 RooDataHist.cxx:1795
 RooDataHist.cxx:1796
 RooDataHist.cxx:1797
 RooDataHist.cxx:1798
 RooDataHist.cxx:1799
 RooDataHist.cxx:1800
 RooDataHist.cxx:1801
 RooDataHist.cxx:1802
 RooDataHist.cxx:1803
 RooDataHist.cxx:1804
 RooDataHist.cxx:1805
 RooDataHist.cxx:1806
 RooDataHist.cxx:1807
 RooDataHist.cxx:1808
 RooDataHist.cxx:1809
 RooDataHist.cxx:1810
 RooDataHist.cxx:1811
 RooDataHist.cxx:1812
 RooDataHist.cxx:1813
 RooDataHist.cxx:1814
 RooDataHist.cxx:1815
 RooDataHist.cxx:1816
 RooDataHist.cxx:1817
 RooDataHist.cxx:1818
 RooDataHist.cxx:1819
 RooDataHist.cxx:1820
 RooDataHist.cxx:1821
 RooDataHist.cxx:1822
 RooDataHist.cxx:1823
 RooDataHist.cxx:1824
 RooDataHist.cxx:1825
 RooDataHist.cxx:1826
 RooDataHist.cxx:1827
 RooDataHist.cxx:1828
 RooDataHist.cxx:1829
 RooDataHist.cxx:1830
 RooDataHist.cxx:1831
 RooDataHist.cxx:1832
 RooDataHist.cxx:1833
 RooDataHist.cxx:1834
 RooDataHist.cxx:1835
 RooDataHist.cxx:1836
 RooDataHist.cxx:1837
 RooDataHist.cxx:1838
 RooDataHist.cxx:1839
 RooDataHist.cxx:1840
 RooDataHist.cxx:1841
 RooDataHist.cxx:1842
 RooDataHist.cxx:1843
 RooDataHist.cxx:1844
 RooDataHist.cxx:1845
 RooDataHist.cxx:1846
 RooDataHist.cxx:1847
 RooDataHist.cxx:1848
 RooDataHist.cxx:1849
 RooDataHist.cxx:1850
 RooDataHist.cxx:1851
 RooDataHist.cxx:1852
 RooDataHist.cxx:1853
 RooDataHist.cxx:1854
 RooDataHist.cxx:1855
 RooDataHist.cxx:1856
 RooDataHist.cxx:1857
 RooDataHist.cxx:1858
 RooDataHist.cxx:1859
 RooDataHist.cxx:1860
 RooDataHist.cxx:1861
 RooDataHist.cxx:1862
 RooDataHist.cxx:1863
 RooDataHist.cxx:1864
 RooDataHist.cxx:1865
 RooDataHist.cxx:1866
 RooDataHist.cxx:1867
 RooDataHist.cxx:1868
 RooDataHist.cxx:1869
 RooDataHist.cxx:1870
 RooDataHist.cxx:1871
 RooDataHist.cxx:1872
 RooDataHist.cxx:1873
 RooDataHist.cxx:1874
 RooDataHist.cxx:1875
 RooDataHist.cxx:1876
 RooDataHist.cxx:1877
 RooDataHist.cxx:1878
 RooDataHist.cxx:1879
 RooDataHist.cxx:1880
 RooDataHist.cxx:1881
 RooDataHist.cxx:1882
 RooDataHist.cxx:1883
 RooDataHist.cxx:1884
 RooDataHist.cxx:1885
 RooDataHist.cxx:1886
 RooDataHist.cxx:1887
 RooDataHist.cxx:1888
 RooDataHist.cxx:1889
 RooDataHist.cxx:1890
 RooDataHist.cxx:1891
 RooDataHist.cxx:1892
 RooDataHist.cxx:1893
 RooDataHist.cxx:1894
 RooDataHist.cxx:1895
 RooDataHist.cxx:1896
 RooDataHist.cxx:1897
 RooDataHist.cxx:1898
 RooDataHist.cxx:1899
 RooDataHist.cxx:1900
 RooDataHist.cxx:1901
 RooDataHist.cxx:1902
 RooDataHist.cxx:1903
 RooDataHist.cxx:1904
 RooDataHist.cxx:1905
 RooDataHist.cxx:1906
 RooDataHist.cxx:1907
 RooDataHist.cxx:1908
 RooDataHist.cxx:1909
 RooDataHist.cxx:1910
 RooDataHist.cxx:1911
 RooDataHist.cxx:1912
 RooDataHist.cxx:1913
 RooDataHist.cxx:1914
 RooDataHist.cxx:1915
 RooDataHist.cxx:1916
 RooDataHist.cxx:1917
 RooDataHist.cxx:1918
 RooDataHist.cxx:1919
 RooDataHist.cxx:1920
 RooDataHist.cxx:1921
 RooDataHist.cxx:1922
 RooDataHist.cxx:1923
 RooDataHist.cxx:1924
 RooDataHist.cxx:1925
 RooDataHist.cxx:1926
 RooDataHist.cxx:1927
 RooDataHist.cxx:1928
 RooDataHist.cxx:1929
 RooDataHist.cxx:1930
 RooDataHist.cxx:1931
 RooDataHist.cxx:1932
 RooDataHist.cxx:1933
 RooDataHist.cxx:1934
 RooDataHist.cxx:1935
 RooDataHist.cxx:1936
 RooDataHist.cxx:1937
 RooDataHist.cxx:1938
 RooDataHist.cxx:1939
 RooDataHist.cxx:1940
 RooDataHist.cxx:1941
 RooDataHist.cxx:1942
 RooDataHist.cxx:1943
 RooDataHist.cxx:1944
 RooDataHist.cxx:1945
 RooDataHist.cxx:1946
 RooDataHist.cxx:1947
 RooDataHist.cxx:1948
 RooDataHist.cxx:1949
 RooDataHist.cxx:1950
 RooDataHist.cxx:1951
 RooDataHist.cxx:1952
 RooDataHist.cxx:1953
 RooDataHist.cxx:1954
 RooDataHist.cxx:1955
 RooDataHist.cxx:1956
 RooDataHist.cxx:1957
 RooDataHist.cxx:1958
 RooDataHist.cxx:1959
 RooDataHist.cxx:1960
 RooDataHist.cxx:1961
 RooDataHist.cxx:1962
 RooDataHist.cxx:1963
 RooDataHist.cxx:1964
 RooDataHist.cxx:1965
 RooDataHist.cxx:1966
 RooDataHist.cxx:1967
 RooDataHist.cxx:1968
 RooDataHist.cxx:1969
 RooDataHist.cxx:1970
 RooDataHist.cxx:1971
 RooDataHist.cxx:1972
 RooDataHist.cxx:1973
 RooDataHist.cxx:1974
 RooDataHist.cxx:1975
 RooDataHist.cxx:1976
 RooDataHist.cxx:1977
 RooDataHist.cxx:1978
 RooDataHist.cxx:1979
 RooDataHist.cxx:1980
 RooDataHist.cxx:1981
 RooDataHist.cxx:1982
 RooDataHist.cxx:1983
 RooDataHist.cxx:1984
 RooDataHist.cxx:1985
 RooDataHist.cxx:1986
 RooDataHist.cxx:1987
 RooDataHist.cxx:1988
 RooDataHist.cxx:1989
 RooDataHist.cxx:1990
 RooDataHist.cxx:1991
 RooDataHist.cxx:1992
 RooDataHist.cxx:1993
 RooDataHist.cxx:1994
 RooDataHist.cxx:1995
 RooDataHist.cxx:1996
 RooDataHist.cxx:1997
 RooDataHist.cxx:1998
 RooDataHist.cxx:1999
 RooDataHist.cxx:2000
 RooDataHist.cxx:2001
 RooDataHist.cxx:2002
 RooDataHist.cxx:2003
 RooDataHist.cxx:2004
 RooDataHist.cxx:2005
 RooDataHist.cxx:2006
 RooDataHist.cxx:2007
 RooDataHist.cxx:2008
 RooDataHist.cxx:2009
 RooDataHist.cxx:2010
 RooDataHist.cxx:2011
 RooDataHist.cxx:2012
 RooDataHist.cxx:2013
 RooDataHist.cxx:2014
 RooDataHist.cxx:2015
 RooDataHist.cxx:2016
 RooDataHist.cxx:2017
 RooDataHist.cxx:2018
 RooDataHist.cxx:2019
 RooDataHist.cxx:2020
 RooDataHist.cxx:2021
 RooDataHist.cxx:2022
 RooDataHist.cxx:2023
 RooDataHist.cxx:2024
 RooDataHist.cxx:2025
 RooDataHist.cxx:2026
 RooDataHist.cxx:2027
 RooDataHist.cxx:2028
 RooDataHist.cxx:2029
 RooDataHist.cxx:2030
 RooDataHist.cxx:2031
 RooDataHist.cxx:2032
 RooDataHist.cxx:2033
 RooDataHist.cxx:2034
 RooDataHist.cxx:2035
 RooDataHist.cxx:2036
 RooDataHist.cxx:2037
 RooDataHist.cxx:2038
 RooDataHist.cxx:2039
 RooDataHist.cxx:2040
 RooDataHist.cxx:2041
 RooDataHist.cxx:2042
 RooDataHist.cxx:2043
 RooDataHist.cxx:2044
 RooDataHist.cxx:2045
 RooDataHist.cxx:2046
 RooDataHist.cxx:2047
 RooDataHist.cxx:2048
 RooDataHist.cxx:2049
 RooDataHist.cxx:2050
 RooDataHist.cxx:2051
 RooDataHist.cxx:2052
 RooDataHist.cxx:2053
 RooDataHist.cxx:2054
 RooDataHist.cxx:2055
 RooDataHist.cxx:2056
 RooDataHist.cxx:2057
 RooDataHist.cxx:2058
 RooDataHist.cxx:2059
 RooDataHist.cxx:2060
 RooDataHist.cxx:2061
 RooDataHist.cxx:2062
 RooDataHist.cxx:2063
 RooDataHist.cxx:2064
 RooDataHist.cxx:2065
 RooDataHist.cxx:2066
 RooDataHist.cxx:2067
 RooDataHist.cxx:2068
 RooDataHist.cxx:2069
 RooDataHist.cxx:2070
 RooDataHist.cxx:2071
 RooDataHist.cxx:2072
 RooDataHist.cxx:2073
 RooDataHist.cxx:2074
 RooDataHist.cxx:2075
 RooDataHist.cxx:2076
 RooDataHist.cxx:2077
 RooDataHist.cxx:2078
 RooDataHist.cxx:2079
 RooDataHist.cxx:2080
 RooDataHist.cxx:2081
 RooDataHist.cxx:2082
 RooDataHist.cxx:2083
 RooDataHist.cxx:2084
 RooDataHist.cxx:2085
 RooDataHist.cxx:2086
 RooDataHist.cxx:2087
 RooDataHist.cxx:2088
 RooDataHist.cxx:2089
 RooDataHist.cxx:2090
 RooDataHist.cxx:2091
 RooDataHist.cxx:2092
 RooDataHist.cxx:2093
 RooDataHist.cxx:2094
 RooDataHist.cxx:2095
 RooDataHist.cxx:2096
 RooDataHist.cxx:2097
 RooDataHist.cxx:2098
 RooDataHist.cxx:2099
 RooDataHist.cxx:2100
 RooDataHist.cxx:2101
 RooDataHist.cxx:2102
 RooDataHist.cxx:2103
 RooDataHist.cxx:2104
 RooDataHist.cxx:2105
 RooDataHist.cxx:2106
 RooDataHist.cxx:2107
 RooDataHist.cxx:2108
 RooDataHist.cxx:2109
 RooDataHist.cxx:2110
 RooDataHist.cxx:2111
 RooDataHist.cxx:2112
 RooDataHist.cxx:2113
 RooDataHist.cxx:2114
 RooDataHist.cxx:2115
 RooDataHist.cxx:2116
 RooDataHist.cxx:2117
 RooDataHist.cxx:2118
 RooDataHist.cxx:2119
 RooDataHist.cxx:2120
 RooDataHist.cxx:2121
 RooDataHist.cxx:2122
 RooDataHist.cxx:2123
 RooDataHist.cxx:2124
 RooDataHist.cxx:2125
 RooDataHist.cxx:2126
 RooDataHist.cxx:2127
 RooDataHist.cxx:2128
 RooDataHist.cxx:2129
 RooDataHist.cxx:2130
 RooDataHist.cxx:2131
 RooDataHist.cxx:2132
 RooDataHist.cxx:2133
 RooDataHist.cxx:2134
 RooDataHist.cxx:2135
 RooDataHist.cxx:2136
 RooDataHist.cxx:2137
 RooDataHist.cxx:2138
 RooDataHist.cxx:2139
 RooDataHist.cxx:2140
 RooDataHist.cxx:2141
 RooDataHist.cxx:2142
 RooDataHist.cxx:2143
 RooDataHist.cxx:2144
 RooDataHist.cxx:2145
 RooDataHist.cxx:2146
 RooDataHist.cxx:2147
 RooDataHist.cxx:2148
 RooDataHist.cxx:2149
 RooDataHist.cxx:2150
 RooDataHist.cxx:2151
 RooDataHist.cxx:2152
 RooDataHist.cxx:2153
 RooDataHist.cxx:2154
 RooDataHist.cxx:2155
 RooDataHist.cxx:2156
 RooDataHist.cxx:2157
 RooDataHist.cxx:2158
 RooDataHist.cxx:2159
 RooDataHist.cxx:2160
 RooDataHist.cxx:2161
 RooDataHist.cxx:2162
 RooDataHist.cxx:2163
 RooDataHist.cxx:2164
 RooDataHist.cxx:2165
 RooDataHist.cxx:2166
 RooDataHist.cxx:2167
 RooDataHist.cxx:2168
 RooDataHist.cxx:2169
 RooDataHist.cxx:2170
 RooDataHist.cxx:2171
 RooDataHist.cxx:2172
 RooDataHist.cxx:2173
 RooDataHist.cxx:2174
 RooDataHist.cxx:2175
 RooDataHist.cxx:2176
 RooDataHist.cxx:2177
 RooDataHist.cxx:2178
 RooDataHist.cxx:2179
 RooDataHist.cxx:2180
 RooDataHist.cxx:2181
 RooDataHist.cxx:2182
 RooDataHist.cxx:2183
 RooDataHist.cxx:2184
 RooDataHist.cxx:2185
 RooDataHist.cxx:2186
 RooDataHist.cxx:2187
 RooDataHist.cxx:2188
 RooDataHist.cxx:2189
 RooDataHist.cxx:2190
 RooDataHist.cxx:2191
 RooDataHist.cxx:2192
 RooDataHist.cxx:2193
 RooDataHist.cxx:2194
 RooDataHist.cxx:2195
 RooDataHist.cxx:2196
 RooDataHist.cxx:2197
 RooDataHist.cxx:2198