#include "RooFit.h"
#include "TIterator.h"
#include "TIterator.h"
#include "TList.h"
#include "RooRealSumPdf.h"
#include "RooRealProxy.h"
#include "RooPlot.h"
#include "RooRealVar.h"
#include "RooAddGenContext.h"
#include "RooRealConstant.h"
#include "RooRealIntegral.h"
ClassImp(RooRealSumPdf)
;
RooRealSumPdf::RooRealSumPdf(const char *name, const char *title) :
  RooAbsPdf(name,title), 
  _codeReg(10),
  _lastFuncIntSet(0),
  _lastFuncNormSet(0),
  _funcIntList(0),
  _funcNormList(0),
  _haveLastCoef(kFALSE),
  _funcList("funcList","List of functions",this),
  _coefList("coefList","List of coefficients",this)
{
  
  _funcIter   = _funcList.createIterator() ;
  _coefIter  = _coefList.createIterator() ;
}
RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
		     RooAbsReal& func1, RooAbsReal& func2, RooAbsReal& coef1) : 
  RooAbsPdf(name,title),
  _codeReg(10),
  _lastFuncIntSet(0),
  _lastFuncNormSet(0),
  _funcIntList(0),
  _funcNormList(0),
  _haveLastCoef(kFALSE),
  _funcList("funcProxyList","List of functions",this),
  _coefList("coefList","List of coefficients",this)
{
  
  _funcIter  = _funcList.createIterator() ;
  _coefIter = _coefList.createIterator() ;
  _funcList.add(func1) ;  
  _funcList.add(func2) ;
  _coefList.add(coef1) ;
}
RooRealSumPdf::RooRealSumPdf(const char *name, const char *title, const RooArgList& funcList, const RooArgList& coefList) :
  RooAbsPdf(name,title),
  _codeReg(10),
  _lastFuncIntSet(0),
  _lastFuncNormSet(0),
  _funcIntList(0),
  _funcNormList(0),
  _haveLastCoef(kFALSE),
  _funcList("funcProxyList","List of functions",this),
  _coefList("coefList","List of coefficients",this)
{ 
  
  
  
  
  
  if (!(funcList.getSize()==coefList.getSize()+1 || funcList.getSize()==coefList.getSize())) {
    cout << "RooRealSumPdf::RooRealSumPdf(" << GetName() 
	 << ") number of pdfs and coefficients inconsistent, must have Nfunc=Ncoef or Nfunc=Ncoef+1" << endl ;
    assert(0) ;
  }
  _funcIter  = _funcList.createIterator() ;
  _coefIter = _coefList.createIterator() ;
 
  
  TIterator* funcIter = funcList.createIterator() ;
  TIterator* coefIter = coefList.createIterator() ;
  RooAbsReal* func ;
  RooAbsReal* coef ;
  while((coef = (RooAbsReal*)coefIter->Next())) {
    func = (RooAbsReal*) funcIter->Next() ;
    if (!dynamic_cast<RooAbsReal*>(coef)) {
      cout << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal, ignored" << endl ;
      continue ;
    }
    if (!dynamic_cast<RooAbsReal*>(func)) {
      cout << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") func " << func->GetName() << " is not of type RooAbsReal, ignored" << endl ;
      continue ;
    }
    _funcList.add(*func) ;
    _coefList.add(*coef) ;    
  }
  func = (RooAbsReal*) funcIter->Next() ;
  if (func) {
    if (!dynamic_cast<RooAbsReal*>(func)) {
      cout << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") last func " << coef->GetName() << " is not of type RooAbsReal, fatal error" << endl ;
      assert(0) ;
    }
    _funcList.add(*func) ;  
  } else {
    _haveLastCoef = kTRUE ;
  }
  
  delete funcIter ;
  delete coefIter  ;
}
RooRealSumPdf::RooRealSumPdf(const RooRealSumPdf& other, const char* name) :
  RooAbsPdf(other,name),
  _codeReg(other._codeReg),
  _lastFuncIntSet(0),
  _lastFuncNormSet(0),
  _funcIntList(0),
  _funcNormList(0),
  _haveLastCoef(other._haveLastCoef),
  _funcList("funcProxyList",this,other._funcList),
  _coefList("coefList",this,other._coefList)
{
  
  _funcIter  = _funcList.createIterator() ;
  _coefIter = _coefList.createIterator() ;
}
RooRealSumPdf::~RooRealSumPdf()
{
  
  delete _funcIter ;
  delete _coefIter ;
  if (_funcIntList) delete _funcIntList ;
  if (_funcNormList) delete _funcNormList ;		     
}
Double_t RooRealSumPdf::evaluate() const 
{
  
  const RooArgSet* nset = _funcList.nset() ;
  Double_t value(0) ;
  
  _funcIter->Reset() ;
  _coefIter->Reset() ;
  RooAbsReal* coef ;
  RooAbsReal* func ;
      
  
  Double_t lastCoef(1) ;
  while((coef=(RooAbsReal*)_coefIter->Next())) {
    func = (RooAbsReal*)_funcIter->Next() ;
    Double_t coefVal = coef->getVal(nset) ;
    if (coefVal) {
      value += func->getVal(nset)*coefVal ;
      lastCoef -= coef->getVal(nset) ;
    }
  }
  
  if (!_haveLastCoef) {
    
    func = (RooAbsReal*) _funcIter->Next() ;
    value += func->getVal(nset)*lastCoef ;
    
    
    if (lastCoef<0 || lastCoef>1) {
      cout << "RooRealSumPdf::evaluate(" << GetName() 
	   << " WARNING: sum of FUNC coefficients not in range [0-1], value=" 
	   << 1-lastCoef << endl ;
    } 
  }
  return value ;
}
Bool_t RooRealSumPdf::checkObservables(const RooArgSet* nset) const 
{
  
  
  
  
  
  
  Bool_t ret(kFALSE) ;
  _funcIter->Reset() ;
  _coefIter->Reset() ;
  RooAbsReal* coef ;
  RooAbsReal* func ;
  while((coef=(RooAbsReal*)_coefIter->Next())) {
    func = (RooAbsReal*)_funcIter->Next() ;
    if (func->observableOverlaps(nset,*coef)) {
      cout << "RooRealSumPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName() 
	   << " and FUNC " << func->GetName() << " have one or more observables in common" << endl ;
      ret = kTRUE ;
    }
    if (coef->dependsOn(*nset)) {
      cout << "RooRealPdf::checkObservables(" << GetName() << "): ERROR coefficient " << coef->GetName() 
	   << " depends on one or more of the following observables" ; nset->Print("1") ;
      ret = kTRUE ;
    }
  }
  
  return ret ;
}
Int_t RooRealSumPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
					     const RooArgSet* normSet2, const char* ) const 
{
  
  if (allVars.getSize()==0) return 0 ;
  if (_forceNumInt) return 0 ;
  
  RooArgSet* allDeps = getObservables(allVars) ;
  RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
  _funcIter->Reset() ;
  _coefIter->Reset() ;
  analVars.add(*allDeps) ;
  Int_t tmp(0) ;
  Int_t masterCode = _codeReg.store(&tmp,1,allDeps,normSet) + 1 ;
  return masterCode ;
}
Double_t RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* ) const 
{
  
  if (code==0) return getVal(normSet2) ;
  
  RooArgSet *allDeps, *normSet ;
  _codeReg.retrieve(code-1,allDeps,normSet) ;
  syncFuncIntList(allDeps) ;
  if (normSet) syncFuncNormList(normSet) ;
  
  TIterator* funcIntIter = _funcIntList->createIterator() ;
  _coefIter->Reset() ;
  RooAbsReal* coef ;
  RooAbsReal* funcInt ;
  Double_t value(0) ;
  
  Double_t lastCoef(1) ;
  while((coef=(RooAbsReal*)_coefIter->Next())) {
    funcInt = (RooAbsReal*)funcIntIter->Next() ;
    Double_t coefVal = coef->getVal(normSet) ;
    if (coefVal) {
      value += funcInt->getVal()*coefVal ;
      lastCoef -= coef->getVal(normSet) ;
    }
  }
  
  if (!_haveLastCoef) {
    
    funcInt = (RooAbsReal*) funcIntIter->Next() ;
    value += funcInt->getVal()*lastCoef ;
    
    
    
    if (lastCoef<0 || lastCoef>1) {
      cout << "RooRealSumPdf::evaluate(" << GetName() 
	   << " WARNING: sum of FUNC coefficients not in range [0-1], value=" 
	   << 1-lastCoef << endl ;
    } 
  }
  delete funcIntIter ;
  
  Double_t normVal(1) ;
  if (normSet) {
    normVal = 0 ;
    
    RooAbsReal* funcNorm ;
    TIterator* funcNormIter = _funcNormList->createIterator() ;
    _coefIter->Reset() ;
    while((coef=(RooAbsReal*)_coefIter->Next())) {
      funcNorm = (RooAbsReal*)funcNormIter->Next() ;
      Double_t coefVal = coef->getVal(normSet) ;
      if (coefVal) {
	normVal += funcNorm->getVal()*coefVal ;
      }
    }
    
    
    if (!_haveLastCoef) {
      funcNorm = (RooAbsReal*) funcNormIter->Next() ;
      normVal += funcNorm->getVal()*lastCoef ;
      
    }
      
    delete funcNormIter ;      
  }
  return value / normVal;
}
void RooRealSumPdf::syncFuncIntList(const RooArgSet* intSet) const
{
  if (intSet==_lastFuncIntSet) return ;
  _lastFuncIntSet = (RooArgSet*) intSet ;
  if (_funcIntList) delete _funcIntList ;
  
  _funcIntList = new RooArgList ;  
  _funcIter->Reset() ;
  RooAbsReal *func  ;
  while((func=(RooAbsReal*)_funcIter->Next())) {
    RooAbsReal* funcInt = func->createIntegral(*intSet) ;
    _funcIntList->addOwned(*funcInt) ;
  }
}
void RooRealSumPdf::syncFuncNormList(const RooArgSet* normSet) const 
{
  if (normSet==_lastFuncNormSet) return ;
  _lastFuncNormSet = (RooArgSet*) normSet ;
  if (_funcNormList) delete _funcNormList ;
  
  
  _funcNormList = new RooArgList ;  
  RooAbsReal *func  ;
  _funcIter->Reset() ;
  while((func=(RooAbsReal*)_funcIter->Next())) {
    RooAbsReal* funcInt = func->createIntegral(*normSet) ;
    _funcNormList->addOwned(*funcInt) ;
  }
}
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.