#include <sstream>
#include <math.h>
#include <stdexcept>
#include "TMath.h"
#include "TH1.h"
#include "Riostream.h"
#include "Riostream.h"
#include "RooFit.h"
#include "RooStats/HistFactory/ParamHistFunc.h"
#include "RooAbsReal.h"
#include "RooAbsPdf.h"
#include "RooConstVar.h"
#include "RooBinning.h"
#include "RooErrorHandler.h"
#include "RooGaussian.h"
#include "RooHistFunc.h"
#include "RooArgSet.h"
#include "RooNLLVar.h"
#include "RooChi2Var.h"
#include "RooMsgService.h"
#include "RooRealVar.h"
#include "RooArgList.h"
#include "RooWorkspace.h"
#include "RooBinning.h"
ClassImp(ParamHistFunc);
ParamHistFunc::ParamHistFunc() : _numBins(0)
{
  ;
}
ParamHistFunc::ParamHistFunc(const char* name, const char* title, 
			     const RooArgList& vars, const RooArgList& paramSet) :
  RooAbsReal(name, title),
  _dataVars("!dataVars","data Vars",       this),
  _paramSet("!paramSet","bin parameters",  this),
  _numBins(0),
  _dataSet( (std::string(name)+"_dataSet").c_str(), "", vars)
{
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  _numBins = GetNumBins( vars );
  
  addVarSet( vars );
  addParamSet( paramSet );
  
}
ParamHistFunc::ParamHistFunc(const char* name, const char* title, 
			     const RooArgList& vars, const RooArgList& paramSet,
			     const TH1* Hist ) :
  RooAbsReal(name, title),
  
  _dataVars("!dataVars","data Vars",       this),
  _paramSet("!paramSet","bin parameters",  this),
  _numBins(0),
  _dataSet( (std::string(name)+"_dataSet").c_str(), "", vars, Hist)
{
  
  
  
  
  
  
  
  
  
  
  
  
  _numBins = GetNumBins( vars );
  
  addVarSet( vars );
  addParamSet( paramSet );
 
}
Int_t ParamHistFunc::GetNumBins( const RooArgSet& vars ) {
  
  
  
  if( vars.getSize() == 0 ) return 0;
    
  Int_t numBins = 1;
  RooFIter varIter = vars.fwdIterator() ;
  RooAbsArg* comp ;
  while((comp = (RooAbsArg*) varIter.next())) {
    if (!dynamic_cast<RooRealVar*>(comp)) {
      std::cout << "ParamHistFunc::GetNumBins" << vars.GetName() << ") ERROR: component " 
	   << comp->GetName() 
	   << " in vars list is not of type RooRealVar" << std::endl ;
      RooErrorHandler::softAbort() ;
      return -1;
    }
    RooRealVar* var = (RooRealVar*) comp;
    Int_t varNumBins = var->numBins();
    numBins *= varNumBins;
  }
    
  return numBins;
}
ParamHistFunc::ParamHistFunc(const ParamHistFunc& other, const char* name) :
  RooAbsReal(other, name), 
  _dataVars("!dataVars", this, other._dataVars ),
  _paramSet("!paramSet", this, other._paramSet),
  _numBins( other._numBins ),
  _binMap( other._binMap ),
  _dataSet( other._dataSet )
{
  ;
  
  
}
ParamHistFunc::~ParamHistFunc() 
{
  ;
}
Int_t ParamHistFunc::getCurrentBin() const {
  
  
  
  
  
  
  
  Int_t dataSetIndex = _dataSet.getIndex( _dataVars ); 
  return dataSetIndex;
}
RooRealVar& ParamHistFunc::getParameter( Int_t index ) const {
  
  
  
  
  
  Int_t gammaIndex = -1;
  if( _binMap.find( index ) != _binMap.end() ) {
    gammaIndex = _binMap[ index ];
  }
  else {
    std::cout << "Error: ParamHistFunc internal bin index map "
	      << "not properly configured" << std::endl;
    throw -1;
  }
  return (RooRealVar&) _paramSet[gammaIndex];
}
RooRealVar& ParamHistFunc::getParameter() const {
  Int_t index = getCurrentBin();
  return getParameter( index );
}
void ParamHistFunc::setParamConst( Int_t index, Bool_t varConst ) {
  RooRealVar& var = getParameter( index );
  var.setConstant( varConst );
}
void ParamHistFunc::setConstant( bool constant ) {
  for( int i=0; i < numBins(); ++i) {
    setParamConst(i, constant);
  }
}
void ParamHistFunc::setShape( TH1* shape ) {
  
  int num_hist_bins = shape->GetNbinsX()*shape->GetNbinsY()*shape->GetNbinsZ();
  if( num_hist_bins != numBins() ) {
    std::cout << "Error - ParamHistFunc: cannot set Shape of ParamHistFunc: " << GetName()
	      << " using histogram: " << shape->GetName()
	      << ". Bins don't match" << std::endl;
    throw std::runtime_error("setShape");
  }
  Int_t TH1BinNumber = 0;
  for( Int_t i = 0; i < numBins(); ++i) {
    
    TH1BinNumber++;
    
    while( shape->IsBinUnderflow(TH1BinNumber) || shape->IsBinOverflow(TH1BinNumber) ){
      TH1BinNumber++;
    }
    
    RooRealVar& var = dynamic_cast<RooRealVar&>(_paramSet[i]);
    var.setVal( shape->GetBinContent(TH1BinNumber) );
  }
}
RooArgList ParamHistFunc::createParamSet(RooWorkspace& w, const std::string& Prefix, 
					 const RooArgList& vars) {
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  RooArgList paramSet;
  Int_t numVars = vars.getSize();
  Int_t numBins = GetNumBins( vars );
  if( numVars == 0 ) {
    std::cout << "Warning - ParamHistFunc::createParamSet() :"
	 << " No Variables provided.  Not making constraint terms." 
	 << std::endl;
    return paramSet;
  }
  else if( numVars == 1 ) {
 
    
    for( Int_t i = 0; i < numBins; ++i) {
      std::stringstream VarNameStream;
      VarNameStream << Prefix << "_bin_" << i;
      std::string VarName = VarNameStream.str();
      RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 ); 
      
      gamma.setMin( 0.0 );
      gamma.setConstant( false );
      w.import( gamma, RooFit::RecycleConflictNodes() );
      RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName.c_str() );
      paramSet.add( *gamma_wspace );
    }
  }
  else if( numVars == 2 ) {
 
    
    
    std::vector< Int_t > Indices(numVars, 0);
    RooRealVar* varx = (RooRealVar*) vars.at(0);
    RooRealVar* vary = (RooRealVar*) vars.at(1);
    
    
    for( Int_t j = 0; j < vary->numBins(); ++j) {
      for( Int_t i = 0; i < varx->numBins(); ++i) {
	
	
	
	std::stringstream VarNameStream;
	VarNameStream << Prefix << "_bin_" << i << "_" << j;
	std::string VarName = VarNameStream.str();
	RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 ); 
	
	gamma.setMin( 0.0 );
	gamma.setConstant( false );
	  
	w.import( gamma, RooFit::RecycleConflictNodes() );
	RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName.c_str() );
	  
	paramSet.add( *gamma_wspace );
	  
      }
    }
  }
  else if( numVars == 3 ) {
 
    
    
    std::vector< Int_t > Indices(numVars, 0);
    RooRealVar* varx = (RooRealVar*) vars.at(0);
    RooRealVar* vary = (RooRealVar*) vars.at(1);
    RooRealVar* varz = (RooRealVar*) vars.at(2);
    
    
    for( Int_t k = 0; k < varz->numBins(); ++k) {
      for( Int_t j = 0; j < vary->numBins(); ++j) {
	for( Int_t i = 0; i < varx->numBins(); ++i) {
	  
	  
	  
	  
	  std::stringstream VarNameStream;
	  VarNameStream << Prefix << "_bin_" << i << "_" << j << "_" << k;
	  std::string VarName = VarNameStream.str();
	
	  RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 ); 
	  
	  gamma.setMin( 0.0 );
	  gamma.setConstant( false );
	  
	  w.import( gamma, RooFit::RecycleConflictNodes() );
	  RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName.c_str() );
	  
	  paramSet.add( *gamma_wspace );
	
	}
      }
    }
  }
  else {
    std::cout << " Error: ParamHistFunc doesn't support dimensions > 3D " <<  std::endl;
  }
  return paramSet;  
}
RooArgList ParamHistFunc::createParamSet(RooWorkspace& w, const std::string& Prefix, 
					 const RooArgList& vars, 
					 Double_t gamma_min, Double_t gamma_max) {
  
  
  
  
  
  
  
  
  
  
  
  
  
  
  RooArgList params = ParamHistFunc::createParamSet( w, Prefix, vars );
  RooFIter paramIter = params.fwdIterator() ;
  RooAbsArg* comp ;
  while((comp = (RooAbsArg*) paramIter.next())) {
    
    RooRealVar* var = (RooRealVar*) comp;
    var->setMin( gamma_min );
    var->setMax( gamma_max );
  }
  return params;
}
RooArgList ParamHistFunc::createParamSet(const std::string& Prefix, Int_t numBins, 
					 Double_t gamma_min, Double_t gamma_max) {
  
  
  
  
  
  
  
  RooArgList paramSet;
  if( gamma_max <= gamma_min ) {
    std::cout << "Warming: gamma_min <= gamma_max: Using default values (0, 10)" << std::endl;
    gamma_min = 0.0;
    gamma_max = 10.0;
  }
  Double_t gamma_nominal = 1.0;
  if( gamma_nominal < gamma_min ) {
    gamma_nominal = gamma_min;
  }
  if( gamma_nominal > gamma_max ) {
    gamma_nominal = gamma_max;
  }
  
  for( Int_t i = 0; i < numBins; ++i) {
    std::stringstream VarNameStream;
    VarNameStream << Prefix << "_bin_" << i;
    std::string VarName = VarNameStream.str();
    RooRealVar* gamma = new RooRealVar( VarName.c_str(), VarName.c_str(), 
					gamma_nominal, gamma_min, gamma_max );
    gamma->setConstant( false );
    paramSet.add( *gamma );
  }
  return paramSet;
}
Int_t ParamHistFunc::addVarSet( const RooArgList& vars ) {
  
  
  
  
  
  
  
  
  int numVars = 0;
  RooFIter varIter = vars.fwdIterator() ;
  RooAbsArg* comp ;
  while((comp = (RooAbsArg*) varIter.next())) {
    if (!dynamic_cast<RooRealVar*>(comp)) {
      coutE(InputArguments) << "ParamHistFunc::(" << GetName() << ") ERROR: component " 
			    << comp->GetName() << " in variables list is not of type RooRealVar" 
			    << std::endl;
      RooErrorHandler::softAbort() ;
      return 1;
    }
    _dataVars.add( *comp );
    numVars++;
  }
  Int_t numBinsX = 1;
  Int_t numBinsY = 1;
  Int_t numBinsZ = 1;
  if( numVars == 1 ) {
    RooRealVar* varX = (RooRealVar*) _dataVars.at(0);
    numBinsX = varX->numBins();
    numBinsY = 1;
    numBinsZ = 1;
  } else  if( numVars == 2 ) {
    RooRealVar* varX = (RooRealVar*) _dataVars.at(0);
    RooRealVar* varY = (RooRealVar*) _dataVars.at(1);
    numBinsX = varX->numBins();
    numBinsY = varY->numBins();
    numBinsZ = 1;
  } else  if( numVars == 3 ) {
    RooRealVar* varX = (RooRealVar*) _dataVars.at(0);
    RooRealVar* varY = (RooRealVar*) _dataVars.at(1);
    RooRealVar* varZ = (RooRealVar*) _dataVars.at(2);
    numBinsX = varX->numBins();
    numBinsY = varY->numBins();
    numBinsZ = varZ->numBins();
  } else {
    std::cout << "ParamHistFunc() - Only works for 1-3 variables (1d-3d)" << std::endl;
    throw -1;  
  }
  
  
  
  _binMap.clear();
  
  for( Int_t i = 0; i < numBinsX; ++i ) {
    for( Int_t j = 0; j < numBinsY; ++j ) {
      for( Int_t k = 0; k < numBinsZ; ++k ) {
	
	Int_t RooDataSetBin = k + j*numBinsZ + i*numBinsY*numBinsZ; 
	Int_t TH1HistBin    = i + j*numBinsX + k*numBinsX*numBinsY; 
	  
	_binMap[RooDataSetBin] = TH1HistBin;
	
      }
    }
  }
  
  return 0;
}
Int_t ParamHistFunc::addParamSet( const RooArgList& params ) {
  
  
  
  
  
  Int_t numVarBins  = _numBins;
  Int_t numElements = params.getSize();
  if( numVarBins != numElements ) {
    std::cout << "ParamHistFunc::addParamSet - ERROR - " 
	      << "Supplied list of parameters " << params.GetName()
	      << " has " << numElements << " elements but the ParamHistFunc"
	      << GetName() << " has " << numVarBins << " bins."
	      << std::endl;
    return 1;
  }
  
  
  
  
  
  RooFIter paramIter = params.fwdIterator() ;
  RooAbsArg* comp ;
  while((comp = (RooAbsArg*) paramIter.next())) {
    if (!dynamic_cast<RooRealVar*>(comp)) {
      coutE(InputArguments) << "ParamHistFunc::(" << GetName() << ") ERROR: component " 
			    << comp->GetName() << " in paramater list is not of type RooRealVar" 
			    << std::endl;
      RooErrorHandler::softAbort() ;
      return 1;
    }
    _paramSet.add( *comp );
  }
  
  return 0;
}
Double_t ParamHistFunc::evaluate() const 
{
  
  
  RooRealVar* param = (RooRealVar*) &(getParameter());
  Double_t value = param->getVal();
  return value;
  
}
Int_t ParamHistFunc::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
					     const RooArgSet* normSet, 
					     const char* ) const 
{
  
  
  if (allVars.getSize()==0) return 0 ;
  if (_forceNumInt) return 0 ;
  
  analVars.add(allVars) ;
  
  Int_t sterileIdx(-1) ;
  CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx,(const char*)0) ;
  if (cache) {
    return _normIntMgr.lastIndex()+1 ;
  }
  
  
  cache = new CacheElem ;
  
  Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,0) ;
  return code+1 ; 
}
Double_t ParamHistFunc::analyticalIntegralWN(Int_t , const RooArgSet* ,
					     const char* ) const 
{
  
  
  Double_t value(0) ;
  
  
  
  
  RooFIter paramIter = _paramSet.fwdIterator();
  RooRealVar* param = NULL;
  Int_t nominalItr = 0;
  while((param = (RooRealVar*) paramIter.next())) {
    
    Double_t paramVal  = (*param).getVal();
    
    
    _dataSet.get( nominalItr );
    Double_t binVolumeDS  = _dataSet.binVolume(); 
    
    
    value += paramVal*binVolumeDS;
    ++nominalItr;
    
  }
  return value;
}
std::list<Double_t>* ParamHistFunc::plotSamplingHint(RooAbsRealLValue& , Double_t , 
						Double_t ) const
{
  
  
  
  return 0;
  
}
std::list<Double_t>* ParamHistFunc::binBoundaries(RooAbsRealLValue& , Double_t , 
						  Double_t ) const 
{
  
  
  
  return 0;
  
}